source: pcp/src/run.c@ f915e1

Last change on this file since f915e1 was f915e1, checked in by Frederik Heber <heber@…>, 17 years ago

function Thermostats() added

  • Property mode set to 100644
File size: 73.3 KB
Line 
1/** \file run.c
2 * Initialization of levels and calculation super-functions.
3 *
4 * Most important functions herein are CalculateForce() and CalculateMD(), which calls various
5 * functions in order to perfom the Molecular Dynamics simulation. MinimiseOccupied() and MinimiseUnOccupied()
6 * call various functions to perform the actual minimisation for the occupied and unoccupied wave functions.
7 * CalculateMinimumStop() evaluates the stop condition for desired precision or step count (or external signals).
8 *
9 * Minor functions are ChangeToLevUp() (which takes the calculation to the next finer level),
10 * UpdateActualPsiNo() (next Psi is minimized and an additional orthonormalization takes place) and UpdateToNewWaves()
11 * (which reinitializes density calculations after the wave functions have changed due to the ionic motion).
12 * OccupyByFermi() re-occupies orbitals according to Fermi distribution if calculated with additional orbitals.
13 * InitRun() and InitRunLevel() prepare the RunStruct with starting values. UpdateIon_PRCG() implements a CG
14 * algorithm to minimize the ionic forces and thus optimize the structure.
15 *
16 *
17 Project: ParallelCarParrinello
18 \author Jan Hamaekers
19 \date 2000
20
21 File: run.c
22 $Id: run.c,v 1.101.2.2 2007-04-21 13:01:13 foo Exp $
23*/
24
25#include <signal.h>
26#include <stdlib.h>
27#include <stdio.h>
28#include <string.h>
29#include <math.h>
30#include <gsl/gsl_multimin.h>
31#include <gsl/gsl_vector.h>
32#include <gsl/gsl_errno.h>
33#include <gsl/gsl_math.h>
34#include <gsl/gsl_min.h>
35#include <gsl/gsl_randist.h>
36#include "mpi.h"
37#include "data.h"
38#include "errors.h"
39#include "helpers.h"
40#include "init.h"
41#include "opt.h"
42#include "myfft.h"
43#include "gramsch.h"
44#include "output.h"
45#include "energy.h"
46#include "density.h"
47#include "ions.h"
48#include "run.h"
49#include "riemann.h"
50#include "mymath.h"
51#include "pcp.h"
52#include "perturbed.h"
53#include "wannier.h"
54
55
56/** Initialization of the (initial) zero and simulation levels in RunStruct structure.
57 * RunStruct::InitLevS is set onto the STANDARTLEVEL in Lattice::Lev[], RunStruct::InitLev0 on
58 * level 0, RunStruct::LevS onto Lattice::MaxLevel-1 (maximum level) and RunStruct::Lev0 onto
59 * Lattice::MaxLevel-2 (one below).
60 * In case of RiemannTensor use an additional Riemann level is intertwined.
61 * \param *P Problem at hand
62 */
63void InitRunLevel(struct Problem *P)
64{
65 struct Lattice *Lat = &P->Lat;
66 struct RunStruct *R = &P->R;
67 struct RiemannTensor *RT = &Lat->RT;
68 int d,i;
69
70 switch (Lat->RT.Use) {
71 case UseNotRT:
72 R->InitLevSNo = STANDARTLEVEL;
73 R->InitLev0No = 0;
74 R->InitLevS = &P->Lat.Lev[R->InitLevSNo];
75 R->InitLev0 = &P->Lat.Lev[R->InitLev0No];
76 R->LevSNo = Lat->MaxLevel-1;
77 R->Lev0No = Lat->MaxLevel-2;
78 R->LevS = &P->Lat.Lev[R->LevSNo];
79 R->Lev0 = &P->Lat.Lev[R->Lev0No];
80 break;
81 case UseRT:
82 R->InitLevSNo = STANDARTLEVEL;
83 R->InitLev0No = 0;
84 R->InitLevS = &P->Lat.Lev[R->InitLevSNo];
85 R->InitLev0 = &P->Lat.Lev[R->InitLev0No];
86
87 /* R->LevSNo = Lat->MaxLevel-1;
88 R->Lev0No = Lat->MaxLevel-2;*/
89 R->LevSNo = Lat->MaxLevel-2;
90 R->Lev0No = Lat->MaxLevel-3;
91
92 R->LevRNo = P->Lat.RT.RiemannLevel;
93 R->LevRSNo = STANDARTLEVEL;
94 R->LevR0No = 0;
95 R->LevS = &P->Lat.Lev[R->LevSNo];
96 R->Lev0 = &P->Lat.Lev[R->Lev0No];
97 R->LevR = &P->Lat.Lev[R->LevRNo];
98 R->LevRS = &P->Lat.Lev[R->LevRSNo];
99 R->LevR0 = &P->Lat.Lev[R->LevR0No];
100 for (d=0; d<NDIM; d++) {
101 RT->NUpLevRS[d] = 1;
102 for (i=R->LevRNo-1; i >= R->LevRSNo; i--)
103 RT->NUpLevRS[d] *= Lat->LevelSizes[i];
104 RT->NUpLevR0[d] = 1;
105 for (i=R->LevRNo-1; i >= R->LevR0No; i--)
106 RT->NUpLevR0[d] *= Lat->LevelSizes[i];
107 }
108 break;
109 }
110}
111
112
113/** Initialization of RunStruct structure.
114 * Most of the actual entries in the RunStruct are set to their starter no-nonsense
115 * values (init if LatticeLevel is not STANDARTLEVEL otherwise normal max): FactorDensity,
116 * all Steps, XCEnergyFactor and HGcFactor, current and archived energie values are zeroed.
117 * \param *P problem at hand
118 */
119void InitRun(struct Problem *P)
120{
121 struct Lattice *Lat = &P->Lat;
122 struct RunStruct *R = &P->R;
123 struct Psis *Psi = &Lat->Psi;
124 int i,j;
125
126#ifndef SHORTSPEED
127 R->MaxMinStepFactor = Psi->AllMaxLocalNo;
128#else
129 R->MaxMinStepFactor = SHORTSPEED;
130#endif
131 if (R->LevSNo == STANDARTLEVEL) {
132 R->ActualMaxMinStep = R->MaxMinStep*R->MaxPsiStep*R->MaxMinStepFactor;
133 R->ActualRelEpsTotalEnergy = R->RelEpsTotalEnergy;
134 R->ActualRelEpsKineticEnergy = R->RelEpsKineticEnergy;
135 R->ActualMaxMinStopStep = R->MaxMinStopStep;
136 R->ActualMaxMinGapStopStep = R->MaxMinGapStopStep;
137 } else {
138 R->ActualMaxMinStep = R->MaxInitMinStep*R->MaxPsiStep*R->MaxMinStepFactor;
139 R->ActualRelEpsTotalEnergy = R->InitRelEpsTotalEnergy;
140 R->ActualRelEpsKineticEnergy = R->InitRelEpsKineticEnergy;
141 R->ActualMaxMinStopStep = R->InitMaxMinStopStep;
142 R->ActualMaxMinGapStopStep = R->InitMaxMinGapStopStep;
143 }
144
145 R->FactorDensityR = 1./Lat->Volume;
146 R->FactorDensityC = Lat->Volume;
147
148 R->OldActualLocalPsiNo = R->ActualLocalPsiNo = 0;
149 R->UseOldPsi = 1;
150 R->MinStep = 0;
151 R->PsiStep = 0;
152 R->AlphaStep = 0;
153 R->DoCalcCGGauss = 0;
154 R->CurrentMin = Occupied;
155
156 R->MinStopStep = 0;
157
158 R->ScanPotential = 0; // in order to deactivate, simply set this to 0
159 R->ScanAtStep = 6; // must not be set to same as ScanPotential (then gradient is never calculated)
160 R->ScanDelta = 0.01; // step size on advance
161 R->ScanFlag = 0; // flag telling that we are scanning
162
163 //R->DoBrent = 0; // InitRun() occurs after ReadParameters(), thus this deactivates DoBrent line search
164
165 /* R->MaxOuterStep = 1;
166 R->MeanForceEps = 0.0;*/
167
168 R->NewRStep = 1;
169 /* Factor */
170 R->XCEnergyFactor = 1.0/R->FactorDensityR;
171 R->HGcFactor = 1.0/Lat->Volume;
172
173 /* Sollte auch geaendert werden */
174 /*Grad->GradientArray[GraSchGradient] = LevS->LPsi->LocalPsi[Psi->LocalNo];*/
175
176 for (j=Occupied;j<Extra;j++)
177 for (i=0; i < RUNMAXOLD; i++) {
178 R->TE[j][i] = 0;
179 R->KE[j][i] = 0;
180 }
181
182 R->MinimisationName = (char **) Malloc((perturbations+3)*(sizeof(char *)), "InitRun: *MinimisationName");
183 for (j=Occupied;j<=Extra;j++)
184 R->MinimisationName[j] = (char *) MallocString(6*(sizeof(char)), "InitRun: MinimisationName[]");
185 strncpy(R->MinimisationName[0],"Occ",6);
186 strncpy(R->MinimisationName[1],"UnOcc",6);
187 strncpy(R->MinimisationName[2],"P0",6);
188 strncpy(R->MinimisationName[3],"P1",6);
189 strncpy(R->MinimisationName[4],"P2",6);
190 strncpy(R->MinimisationName[5],"RxP0",6);
191 strncpy(R->MinimisationName[6],"RxP1",6);
192 strncpy(R->MinimisationName[7],"RxP2",6);
193 strncpy(R->MinimisationName[8],"Extra",6);
194}
195
196/** Re-occupy orbitals according to Fermi (bottom-up energy-wise).
197 * All OnePsiElementAddData#PsiFactor's are set to zero. \a electrons is set to Psi#Use-dependent
198 * Psis#GlobalNo.
199 * Then we go through OnePsiElementAddData#Lambda, find biggest, put one or two electrons into
200 * its PsiFactor, withdraw from \a electrons. Go on as long as there are \a electrons left.
201 * \param *P Problem at hand
202 */
203void OccupyByFermi(struct Problem *P) {
204 struct Lattice *Lat = &P->Lat;
205 struct Psis *Psi = &Lat->Psi;
206 int i, index, electrons = 0;
207 double lambda, electronsperorbit;
208
209 for (i=0; i< Psi->LocalNo; i++) {// set all PsiFactors to zero
210 Psi->LocalPsiStatus[i].PsiFactor = 0.0;
211 Psi->LocalPsiStatus[i].PsiType = UnOccupied;
212 //Psi->LocalPsiStatus[i].PsiGramSchStatus = (R->DoSeparated) ? NotUsedToOrtho : NotOrthogonal;
213 }
214
215 electronsperorbit = (Psi->Use == UseSpinUpDown) ? 1 : 2;
216 switch (Psi->PsiST) { // how many electrons may we re-distribute
217 case SpinDouble:
218 electrons = Psi->GlobalNo[PsiMaxNoDouble];
219 break;
220 case SpinUp:
221 electrons = Psi->GlobalNo[PsiMaxNoUp];
222 break;
223 case SpinDown:
224 electrons = Psi->GlobalNo[PsiMaxNoDown];
225 break;
226 }
227 while (electrons > 0) {
228 lambda = 0.0;
229 index = 0;
230 for (i=0; i< Psi->LocalNo; i++) // seek biggest unoccupied one
231 if ((lambda < Psi->AddData[i].Lambda) && (Psi->LocalPsiStatus[i].PsiFactor == 0.0)) {
232 index = i;
233 lambda = Psi->AddData[i].Lambda;
234 }
235 Psi->LocalPsiStatus[index].PsiFactor = electronsperorbit; // occupy state
236 Psi->LocalPsiStatus[index].PsiType = Occupied;
237 electrons--; // one electron less
238 }
239 for (i=0; i< Psi->LocalNo; i++) // set all PsiFactors to zero
240 if (Psi->LocalPsiStatus[i].PsiType == UnOccupied) Psi->LocalPsiStatus[i].PsiFactor = 1.0;
241
242 SpeedMeasure(P, DensityTime, StartTimeDo);
243 UpdateDensityCalculation(P);
244 SpeedMeasure(P, DensityTime, StopTimeDo);
245 InitPsiEnergyCalculation(P,Occupied); // goes through all orbitals calculating kinetic and non-local
246 CalculateDensityEnergy(P, 0);
247 EnergyAllReduce(P);
248// SetCurrentMinState(P,UnOccupied);
249// InitPsiEnergyCalculation(P,UnOccupied); /* STANDARTLEVEL */
250// CalculateGapEnergy(P); /* STANDARTLEVEL */
251// EnergyAllReduce(P);
252// SetCurrentMinState(P,Occupied);
253}
254
255/** Use next local Psi: Update RunStruct::ActualLocalPsiNo.
256 * Increases OnePsiElementAddData::Step, RunStruct::MinStep and RunStruct::PsiStep.
257 * RunStruct::OldActualLocalPsiNo is set to current one and this distributed
258 * (UpdateGramSchOldActualPsiNo()) among process.
259 * Afterwards RunStruct::ActualLocalPsiNo is increased (modulo Psis::LocalNo of
260 * this process) and again distributed (UpdateGramSchActualPsiNo()).
261 * Due to change in the GramSchmidt-Status, GramSch() is called for Orthonormalization.
262 * \param *P Problem at hand#
263 * \param IncType skip types PsiTypeTag#UnOccupied or PsiTypeTag#Occupied we only want next(thus we can handily advance only through either type)
264 */
265void UpdateActualPsiNo(struct Problem *P, enum PsiTypeTag IncType)
266{
267 struct RunStruct *R = &P->R;
268 if (R->CurrentMin != IncType) {
269 SetCurrentMinState(P,IncType);
270 R->PsiStep = R->MaxPsiStep; // force step to next Psi
271 }
272 P->Lat.Psi.AddData[R->ActualLocalPsiNo].Step++;
273 R->MinStep++;
274 R->PsiStep++;
275 if (R->OldActualLocalPsiNo != R->ActualLocalPsiNo) {
276 R->OldActualLocalPsiNo = R->ActualLocalPsiNo;
277 UpdateGramSchOldActualPsiNo(P, &P->Lat.Psi);
278 }
279 if (R->PsiStep >= R->MaxPsiStep) {
280 R->PsiStep=0;
281 do {
282 R->ActualLocalPsiNo++;
283 R->ActualLocalPsiNo %= P->Lat.Psi.LocalNo;
284 } while (P->Lat.Psi.AllPsiStatus[R->ActualLocalPsiNo].PsiType != IncType);
285 UpdateGramSchActualPsiNo(P, &P->Lat.Psi);
286 //fprintf(stderr,"(%i) ActualLocalNo: %d\n", P->Par.me, R->ActualLocalPsiNo);
287 }
288 if ((R->UseAddGramSch == 1 && (R->OldActualLocalPsiNo != R->ActualLocalPsiNo || P->Lat.Psi.NoOfPsis == 1)) || R->UseAddGramSch == 2) {
289 if (P->Lat.Psi.LocalPsiStatus[R->OldActualLocalPsiNo].PsiGramSchStatus != NotUsedToOrtho) // don't reset by accident last psi of former minimisation run
290 SetGramSchOldActualPsi(P, &P->Lat.Psi, NotOrthogonal);
291 SpeedMeasure(P, GramSchTime, StartTimeDo);
292 if (R->CurrentMin <= UnOccupied)
293 GramSch(P, R->LevS, &P->Lat.Psi, Orthonormalize);
294 else
295 GramSch(P, R->LevS, &P->Lat.Psi, Orthogonalize); //Orthogonalize
296 SpeedMeasure(P, GramSchTime, StopTimeDo);
297 }
298}
299
300/** Resets all OnePsiElement#DoBrent.\
301 * \param *P Problem at hand
302 * \param *Psi pointer to wave functions
303 */
304void ResetBrent(struct Problem *P, struct Psis *Psi) {
305 int i;
306 for (i=0; i< Psi->LocalNo; i++) {// set all PsiFactors to zero
307 //fprintf(stderr,"(%i) DoBrent[%i] = %i\n", P->Par.me, i, Psi->LocalPsiStatus[i].DoBrent);
308 if (Psi->LocalPsiStatus[i].PsiType == Occupied) Psi->LocalPsiStatus[i].DoBrent = 4;
309 }
310}
311
312/** Sets current minimisation state.
313 * Stores given \a state in RunStruct#CurrentMin and sets pointer Lattice#E accordingly.
314 * \param *P Problem at hand
315 * \param state given PsiTypeTag state
316 */
317void SetCurrentMinState(struct Problem *P, enum PsiTypeTag state) {
318 P->R.CurrentMin = state;
319 P->R.TotalEnergy = &(P->R.TE[state][0]);
320 P->R.KineticEnergy = &(P->R.KE[state][0]);
321 P->R.ActualRelTotalEnergy = &(P->R.ActualRelTE[state][0]);
322 P->R.ActualRelKineticEnergy = &(P->R.ActualRelKE[state][0]);
323 P->Lat.E = &(P->Lat.Energy[state]);
324}
325/*{
326 struct RunStruct *R = &P->R;
327 struct Lattice *Lat = &P->Lat;
328 struct Psis *Psi = &Lat->Psi;
329 P->Lat.Psi.AddData[R->ActualLocalPsiNo].Step++;
330 R->MinStep++;
331 R->PsiStep++;
332 if (R->OldActualLocalPsiNo != R->ActualLocalPsiNo) { // remember old actual local number
333 R->OldActualLocalPsiNo = R->ActualLocalPsiNo;
334 UpdateGramSchOldActualPsiNo(P, &P->Lat.Psi);
335 }
336 if (R->PsiStep >= R->MaxPsiStep) { // done enough minimisation steps for this orbital?
337 R->PsiStep=0;
338 do { // step on as long as we are still on a SkipType orbital
339 R->ActualLocalPsiNo++;
340 R->ActualLocalPsiNo %= P->Lat.Psi.LocalNo;
341 } while ((P->Lat.Psi.LocalPsiStatus[R->ActualLocalPsiNo].PsiType == SkipType));
342 UpdateGramSchActualPsiNo(P, &P->Lat.Psi);
343 if (R->UseAddGramSch >= 1) {
344 SetGramSchOldActualPsi(P,Psi,NotOrthogonal);
345 // setze von OldActual bis bla auf nicht orthogonal
346 GramSch(P, R->LevS, &P->Lat.Psi, Orthonormalize);
347 }
348 } else if (R->UseAddGramSch == 2) {
349 SetGramSchOldActualPsi(P, &P->Lat.Psi, NotOrthogonal);
350 //if (SkipType == UnOccupied)
351 //ResetGramSch(P,Psi);
352 //fprintf(stderr,"UpdateActualPsiNo: GramSch() for %i\n",R->OldActualLocalPsiNo);
353 GramSch(P, R->LevS, &P->Lat.Psi, Orthonormalize);
354 }
355}*/
356
357/** Upgrades the calculation to the next finer level.
358 * If we are below the initial level,
359 * ChangePsiAndDensToLevUp() prepares density and Psi coefficients.
360 * Then the level change is made as RunStruct::LevSNo and RunStruct::Lev0No are decreased.
361 * The RunStruct::OldActualLocalPsi is set to current one and both are distributed
362 * (UpdateGramSchActualPsiNo(), UpdateGramSchOldActualPsiNo()).
363 * The PseudoPot'entials adopt the level up by calling ChangePseudoToLevUp().
364 * Now we are prepared to reset Energy::PsiEnergy and local and total density energy and
365 * recalculate them: InitPsiEnergyCalculation(), CalculateDensityEnergy() and CalculateIonsEnergy().
366 * Results are gathered EnergyAllReduce() and the output made EnergyOutput().
367 * Finally, the stop condition are reset for the new level (depending if it's the STANDARTLEVEL or
368 * not).
369 * \param *P Problem at hand
370 * \param *Stop is set to zero if we are below or equal to init level (see CalculateForce())
371 * \sa UpdateToNewWaves() very similar in the procedure, only the update of the Psis and density
372 * (ChangePsiAndDensToLevUp()) is already made there.
373 * \bug Missing TotalEnergy shifting for other PsiTypeTag's!
374 */
375static void ChangeToLevUp(struct Problem *P, int *Stop)
376{
377 struct RunStruct *R = &P->R;
378 struct Lattice *Lat = &P->Lat;
379 struct Psis *Psi = &Lat->Psi;
380 struct Energy *E = Lat->E;
381 struct RiemannTensor *RT = &Lat->RT;
382 int i;
383 if (R->LevSNo <= R->InitLevSNo) {
384 if (P->Call.out[LeaderOut] && (P->Par.me == 0))
385 fprintf(stderr, "(%i) ChangeLevUp: LevSNo(%i) <= InitLevSNo(%i)\n", P->Par.me, R->LevSNo, R->InitLevSNo);
386 *Stop = 1;
387 return;
388 }
389 if (P->Call.out[LeaderOut] && (P->Par.me == 0))
390 fprintf(stderr, "(0) ChangeLevUp: LevSNo(%i) InitLevSNo(%i)\n", R->LevSNo, R->InitLevSNo);
391 *Stop = 0;
392 P->Speed.LevUpSteps++;
393 SpeedMeasure(P, SimTime, StopTimeDo);
394 SpeedMeasure(P, InitSimTime, StartTimeDo);
395 SpeedMeasure(P, InitDensityTime, StartTimeDo);
396 ChangePsiAndDensToLevUp(P);
397 SpeedMeasure(P, InitDensityTime, StopTimeDo);
398 R->LevSNo--;
399 R->Lev0No--;
400 if (RT->ActualUse == standby && R->LevSNo == STANDARTLEVEL) {
401 P->Lat.RT.ActualUse = active;
402 CalculateRiemannTensorData(P);
403 Error(SomeError, "Calculate RT: Not further implemented");
404 }
405 R->LevS = &P->Lat.Lev[R->LevSNo];
406 R->Lev0 = &P->Lat.Lev[R->Lev0No];
407 R->OldActualLocalPsiNo = R->ActualLocalPsiNo;
408 UpdateGramSchActualPsiNo(P, &P->Lat.Psi);
409 UpdateGramSchOldActualPsiNo(P, &P->Lat.Psi);
410 ResetBrent(P, &P->Lat.Psi);
411 R->PsiStep=0;
412 R->MinStep=0;
413 P->Grad.GradientArray[GraSchGradient] = R->LevS->LPsi->LocalPsi[Psi->LocalNo];
414 ChangePseudoToLevUp(P);
415 for (i=0; i<MAXALLPSIENERGY; i++)
416 SetArrayToDouble0(E->PsiEnergy[i], Psi->LocalNo);
417 SetArrayToDouble0(E->AllLocalDensityEnergy, MAXALLDENSITYENERGY);
418 SetArrayToDouble0(E->AllTotalDensityEnergy, MAXALLDENSITYENERGY);
419 for (i=MAXOLD-1; i > 0; i--) {
420 E->TotalEnergy[i] = E->TotalEnergy[i-1];
421 Lat->Energy[UnOccupied].TotalEnergy[i] = Lat->Energy[UnOccupied].TotalEnergy[i-1];
422 }
423 InitPsiEnergyCalculation(P,Occupied);
424 CalculateDensityEnergy(P,1);
425 CalculateIonsEnergy(P);
426 EnergyAllReduce(P);
427/* SetCurrentMinState(P,UnOccupied);
428 InitPsiEnergyCalculation(P,UnOccupied);
429 CalculateGapEnergy(P);
430 EnergyAllReduce(P);
431 SetCurrentMinState(P,Occupied);*/
432 EnergyOutput(P,0);
433 if (R->LevSNo == STANDARTLEVEL) {
434 R->ActualMaxMinStep = R->MaxMinStep*R->MaxPsiStep*R->MaxMinStepFactor;
435 R->ActualRelEpsTotalEnergy = R->RelEpsTotalEnergy;
436 R->ActualRelEpsKineticEnergy = R->RelEpsKineticEnergy;
437 R->ActualMaxMinStopStep = R->MaxMinStopStep;
438 R->ActualMaxMinGapStopStep = R->MaxMinGapStopStep;
439 } else {
440 R->ActualMaxMinStep = R->MaxInitMinStep*R->MaxPsiStep*R->MaxMinStepFactor;
441 R->ActualRelEpsTotalEnergy = R->InitRelEpsTotalEnergy;
442 R->ActualRelEpsKineticEnergy = R->InitRelEpsKineticEnergy;
443 R->ActualMaxMinStopStep = R->InitMaxMinStopStep;
444 R->ActualMaxMinGapStopStep = R->InitMaxMinGapStopStep;
445 }
446 R->MinStopStep = 0;
447 SpeedMeasure(P, InitSimTime, StopTimeDo);
448 SpeedMeasure(P, SimTime, StartTimeDo);
449 if (P->Call.out[LeaderOut] && (P->Par.me == 0))
450 fprintf(stderr, "(0) ChangeLevUp: LevSNo(%i) InitLevSNo(%i) Done\n", R->LevSNo, R->InitLevSNo);
451}
452
453/** Updates after the wave functions have changed (e.g.\ Ion moved).
454 * Old and current RunStruct::ActualLocalPsiNo are zeroed and distributed among the processes.
455 * InitDensityCalculation() is called, afterwards the pseudo potentials update to the new
456 * wave functions UpdatePseudoToNewWaves().
457 * Energy::AllLocalDensityEnergy, Energy::AllTotalDensityEnergy, Energy::AllTotalIonsEnergy and
458 * Energy::PsiEnergy[i] are set to zero.
459 * We are set to recalculate all of the following energies: Psis InitPsiEnergyCalculation(), density
460 * CalculateDensityEnergy(), ionic CalculateIonsEnergy() and ewald CalculateEwald().
461 * Results are gathered from all processes EnergyAllReduce() and EnergyOutput() is called.
462 * Finally, the various conditons in the RunStruct for stopping the calculation are set: number of
463 * minimisation steps, relative total or kinetic energy change or how often stop condition was
464 * evaluated.
465 * \param *P Problem at hand
466 */
467static void UpdateToNewWaves(struct Problem *P)
468{
469 struct RunStruct *R = &P->R;
470 struct Lattice *Lat = &P->Lat;
471 struct Psis *Psi = &Lat->Psi;
472 struct Energy *E = Lat->E;
473 int i,type;
474 R->OldActualLocalPsiNo = R->ActualLocalPsiNo = 0;
475 //if (isnan((double)R->LevS->LPsi->LocalPsi[R->OldActualLocalPsiNo][0].re)) { fprintf(stderr,"(%i) WARNING in UpdateGramSchActualPsiNo(): LPsi->LocalPsi[%i]_[%i] = NaN!\n", P->Par.me, R->OldActualLocalPsiNo, 0); Error(SomeError, "NaN-Fehler!"); }
476 UpdateGramSchActualPsiNo(P, &P->Lat.Psi);
477 UpdateGramSchOldActualPsiNo(P, &P->Lat.Psi);
478 R->PsiStep=0;
479 R->MinStep=0;
480 SpeedMeasure(P, InitDensityTime, StartTimeDo);
481 //if (isnan((double)R->LevS->LPsi->LocalPsi[R->OldActualLocalPsiNo][0].re)) { fprintf(stderr,"(%i) WARNING in Update../InitDensityCalculation(): LPsi->LocalPsi[%i]_[%i] = NaN!\n", P->Par.me, R->OldActualLocalPsiNo, 0); Error(SomeError, "NaN-Fehler!"); }
482 InitDensityCalculation(P);
483 SpeedMeasure(P, InitDensityTime, StopTimeDo);
484 UpdatePseudoToNewWaves(P);
485 for (i=0; i<MAXALLPSIENERGY; i++)
486 SetArrayToDouble0(E->PsiEnergy[i], Psi->LocalNo);
487 SetArrayToDouble0(E->AllLocalDensityEnergy, MAXALLDENSITYENERGY);
488 SetArrayToDouble0(E->AllTotalDensityEnergy, MAXALLDENSITYENERGY);
489 SetArrayToDouble0(E->AllTotalIonsEnergy, MAXALLIONSENERGY);
490 InitPsiEnergyCalculation(P,Occupied);
491 CalculateDensityEnergy(P,1);
492 CalculateIonsEnergy(P);
493 CalculateEwald(P, 0);
494 EnergyAllReduce(P);
495 if (R->DoUnOccupied) {
496 SetCurrentMinState(P,UnOccupied);
497 InitPsiEnergyCalculation(P,UnOccupied); /* STANDARTLEVEL */
498 CalculateGapEnergy(P); /* STANDARTLEVEL */
499 EnergyAllReduce(P);
500 }
501 if (R->DoPerturbation)
502 for(type=Perturbed_P0;type <=Perturbed_RxP2;type++) {
503 SetCurrentMinState(P,type);
504 InitPerturbedEnergyCalculation(P,1); /* STANDARTLEVEL */
505 EnergyAllReduce(P);
506 }
507 SetCurrentMinState(P,Occupied);
508 E->TotalEnergyOuter[0] = E->TotalEnergy[0];
509 EnergyOutput(P,0);
510 R->ActualMaxMinStep = R->MaxMinStep*R->MaxPsiStep*R->MaxMinStepFactor;
511 R->ActualRelEpsTotalEnergy = R->RelEpsTotalEnergy;
512 R->ActualRelEpsKineticEnergy = R->RelEpsKineticEnergy;
513 R->ActualMaxMinStopStep = R->MaxMinStopStep;
514 R->ActualMaxMinGapStopStep = R->MaxMinGapStopStep;
515 R->MinStopStep = 0;
516}
517
518/** Evaluates the stop condition and returns 0 or 1 for occupied states.
519 * Stop is set when:
520 * - SuperStop at best possible point (e.g.\ LevelChange): RunStruct::PsiStep == 0 && SuperStop == 1
521 * - RunStruct::PsiStep && RunStruct::MinStopStep modulo RunStruct::ActualMaxMinStopStep == 0
522 * - To many minimisation steps: RunStruct::MinStep > RunStruct::ActualMaxMinStopStep
523 * - below relative rate of change:
524 * - Remember old values: Shift all RunStruct::TotalEnergy and RunStruct::KineticEnergy by
525 * one and transfer current one from Energy::TotalEnergy and Energy::AllTotalPsiEnergy[KineticEnergy].
526 * - if more than one minimisation step was made, calculate the relative changes of total
527 * energy and kinetic energy and store them in RunStruct::ActualRelTotalEnergy and
528 * RunStruct::ActualRelKineticEnergy and check them against the sought for minimum
529 * values RunStruct::ActualRelEpsTotalEnergy and RunStruct::ActualRelEpsKineticEnergy.
530 * - if RunStruct::PsiStep is zero (default), increase RunStruct::MinStopStep
531 * \param *P Problem at hand
532 * \param SuperStop 1 - external signal: ceasing calculation, 0 - no signal
533 * \return Stop: 1 - stop, 0 - continue
534 */
535int CalculateMinimumStop(struct Problem *P, int SuperStop)
536{
537 int Stop = 0, i;
538 struct RunStruct *R = &P->R;
539 struct Energy *E = P->Lat.E;
540 if (R->PsiStep == 0 && SuperStop) Stop = 1;
541 if (R->PsiStep == 0 && ((R->MinStopStep % R->ActualMaxMinStopStep == 0 && R->CurrentMin != UnOccupied) || (R->MinStopStep % R->ActualMaxMinGapStopStep == 0 && R->CurrentMin == UnOccupied))) {
542// if (R->MinStep >= R->ActualMaxMinStep) {
543// Stop = 1;
544// fprintf(stderr,"(%i) MinStep %i >= %i MaxMinStep.\n", P->Par.me, R->MinStep, R->ActualMaxMinStep);
545// }
546 for (i=RUNMAXOLD-1; i > 0; i--) {
547 R->TotalEnergy[i] = R->TotalEnergy[i-1];
548 R->KineticEnergy[i] = R->KineticEnergy[i-1];
549 }
550 R->TotalEnergy[0] = E->TotalEnergy[0];
551 R->KineticEnergy[0] = E->AllTotalPsiEnergy[KineticEnergy];
552 if (R->MinStopStep) {
553 //if (R->TotalEnergy[1] < MYEPSILON) fprintf(stderr,"CalculateMinimumStop: R->TotalEnergy[1] = %lg\n",R->TotalEnergy[1]);
554 R->ActualRelTotalEnergy[0] = fabs((R->TotalEnergy[0]-R->TotalEnergy[1])/R->TotalEnergy[1]);
555 //if (R->KineticEnergy[1] < MYEPSILON) fprintf(stderr,"CalculateMinimumStop: R->KineticEnergy[1] = %lg\n",R->KineticEnergy[1]);
556 //if (R->CurrentMin < Perturbed_P0)
557 R->ActualRelKineticEnergy[0] = fabs((R->KineticEnergy[0]-R->KineticEnergy[1])/R->KineticEnergy[1]);
558 //else R->ActualRelKineticEnergy[0] = 0.;
559 if (P->Call.out[LeaderOut] && (P->Par.me == 0))
560 switch (R->CurrentMin) {
561 default:
562 fprintf(stderr, "ARelTE: %e\tARelKE: %e\n", R->ActualRelTotalEnergy[0], R->ActualRelKineticEnergy[0]);
563 break;
564 case UnOccupied:
565 fprintf(stderr, "ARelTGE: %e\tARelKGE: %e\n", R->ActualRelTotalEnergy[0], R->ActualRelKineticEnergy[0]);
566 break;
567 }
568 if ((R->ActualRelTotalEnergy[0] < R->ActualRelEpsTotalEnergy) &&
569 (R->ActualRelKineticEnergy[0] < R->ActualRelEpsKineticEnergy))
570 Stop = 1;
571 }
572 }
573 if (R->PsiStep == 0)
574 R->MinStopStep++;
575 if (P->Call.WriteSrcFiles == 2)
576 OutputVisSrcFiles(P, R->CurrentMin);
577 return(Stop);
578}
579
580/** Evaluates the stop condition and returns 0 or 1 for gap energies.
581 * Stop is set when:
582 * - SuperStop at best possible point (e.g.\ LevelChange): RunStruct::PsiStep == 0 && SuperStop == 1
583 * - RunStruct::PsiStep && RunStruct::MinStopStep modulo RunStruct::ActualMaxMinStopStep == 0
584 * - To many minimisation steps: RunStruct::MinStep > RunStruct::ActualMaxMinStopStep
585 * - below relative rate of change:
586 * - Remember old values: Shift all RunStruct::TotalEnergy and RunStruct::KineticEnergy by
587 * one and transfer current one from Energy::TotalEnergy and Energy::AllTotalPsiEnergy[KineticEnergy].
588 * - if more than one minimisation step was made, calculate the relative changes of total
589 * energy and kinetic energy and store them in RunStruct::ActualRelTotalEnergy and
590 * RunStruct::ActualRelKineticEnergy and check them against the sought for minimum
591 * values RunStruct::ActualRelEpsTotalEnergy and RunStruct::ActualRelEpsKineticEnergy.
592 * - if RunStruct::PsiStep is zero (default), increase RunStruct::MinStopStep
593 * \param *P Problem at hand
594 * \param SuperStop 1 - external signal: ceasing calculation, 0 - no signal
595 * \return Stop: 1 - stop, 0 - continue
596 * \sa CalculateMinimumStop() - same procedure for occupied states
597 *//*
598static double CalculateGapStop(struct Problem *P, int SuperStop)
599{
600 int Stop = 0, i;
601 struct RunStruct *R = &P->R;
602 struct Lattice *Lat = &P->Lat;
603 struct Energy *E = P->Lat.E;
604 if (R->PsiStep == 0 && SuperStop) Stop = 1;
605 if (R->PsiStep == 0 && (R->MinStopStep % R->ActualMaxMinGapStopStep) == 0) {
606 if (R->MinStep >= R->ActualMaxMinStep) Stop = 1;
607 for (i=RUNMAXOLD-1; i > 0; i--) {
608 R->TotalGapEnergy[i] = R->TotalGapEnergy[i-1];
609 R->KineticGapEnergy[i] = R->KineticGapEnergy[i-1];
610 }
611 R->TotalGapEnergy[0] = Lat->Energy[UnOccupied].TotalEnergy[0];
612 R->KineticGapEnergy[0] = E->AllTotalPsiEnergy[GapPsiEnergy];
613 if (R->MinStopStep) {
614 if (R->TotalGapEnergy[1] < MYEPSILON) fprintf(stderr,"CalculateMinimumStop: R->TotalGapEnergy[1] = %lg\n",R->TotalGapEnergy[1]);
615 R->ActualRelTotalGapEnergy[0] = fabs((R->TotalGapEnergy[0]-R->TotalGapEnergy[1])/R->TotalGapEnergy[1]);
616 if (R->KineticGapEnergy[1] < MYEPSILON) fprintf(stderr,"CalculateMinimumStop: R->KineticGapEnergy[1] = %lg\n",R->KineticGapEnergy[1]);
617 R->ActualRelKineticGapEnergy[0] = fabs((R->KineticGapEnergy[0]-R->KineticGapEnergy[1])/R->KineticGapEnergy[1]);
618 if (P->Call.out[LeaderOut] && (P->Par.me == 0))
619 fprintf(stderr, "(%i) -------------------------> ARelTGE: %e\tARelKGE: %e\n", P->Par.me, R->ActualRelTotalGapEnergy[0], R->ActualRelKineticGapEnergy[0]);
620 if ((R->ActualRelTotalGapEnergy[0] < R->ActualRelEpsTotalGapEnergy) &&
621 (R->ActualRelKineticGapEnergy[0] < R->ActualRelEpsKineticGapEnergy))
622 Stop = 1;
623 }
624 }
625 if (R->PsiStep == 0)
626 R->MinStopStep++;
627
628 return(Stop);
629}*/
630
631#define StepTolerance 1e-4
632
633static void CalculateEnergy(struct Problem *P) {
634 SpeedMeasure(P, DensityTime, StartTimeDo);
635 UpdateDensityCalculation(P);
636 SpeedMeasure(P, DensityTime, StopTimeDo);
637 UpdatePsiEnergyCalculation(P);
638 CalculateDensityEnergy(P, 0);
639 //CalculateIonsEnergy(P);
640 EnergyAllReduce(P);
641}
642
643/** Energy functional depending on one parameter \a x (for a certain Psi in a certain conjugate direction).
644 * \param x parameter for the which the function must be minimized
645 * \param *params additional params
646 * \return total energy if Psi is changed according to the given parameter
647 */
648static double fn1 (double x, void * params) {
649 struct Problem *P = (struct Problem *)(params);
650 struct RunStruct *R = &P->R;
651 struct Lattice *Lat = &P->Lat;
652 struct LatticeLevel *LevS = R->LevS;
653 int ElementSize = (sizeof(fftw_complex) / sizeof(double));
654 int i=R->ActualLocalPsiNo;
655 double ret;
656
657 //fprintf(stderr,"(%i) Evaluating fnl at %lg ...\n",P->Par.me, x);
658 //TestForOrth(P,R->LevS,P->Grad.GradientArray[GraSchGradient]);
659 CalculateNewWave(P, &x); // also stores Psi to oldPsi
660 //TestGramSch(P,R->LevS,&P->Lat.Psi,Occupied);
661 //fprintf(stderr,"(%i) Testing for orthogonality of %i against ...\n",P->Par.me, R->ActualLocalPsiNo);
662 //TestForOrth(P, LevS, LevS->LPsi->LocalPsi[R->ActualLocalPsiNo]);
663 //UpdateActualPsiNo(P, Occupied);
664 //UpdateEnergyArray(P);
665 CalculateEnergy(P);
666 ret = Lat->E->TotalEnergy[0];
667 memcpy(LevS->LPsi->LocalPsi[i], LevS->LPsi->OldLocalPsi[i], ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
668 //fprintf(stderr,"(%i) Psi %i at %p retrieved from OldPsi at %p: Old[0] %lg+i%lg\n", P->Par.me, R->ActualLocalPsiNo, LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].re, LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].im);
669 CalculateEnergy(P);
670 //fprintf(stderr,"(%i) fnl(%lg) = %lg\n", P->Par.me, x, ret);
671 return ret;
672}
673
674#ifdef HAVE_INLINE
675inline void flip(double *a, double *b) {
676#else
677void flip(double *a, double *b) {
678#endif
679 double tmp = *a;
680 *a = *b;
681 *b = tmp;
682}
683
684
685/** Minimise PsiType#Occupied orbitals.
686 * It is checked whether CallOptions#ReadSrcFiles is set and thus coefficients for the level have to be
687 * read from file and afterwards initialized.
688 *
689 * Then follows the main loop, until a stop condition is met:
690 * -# CalculateNewWave()\n
691 * Over a conjugate gradient method the next (minimal) wave function is sought for.
692 * -# UpdateActualPsiNo()\n
693 * Switch local Psi to next one.
694 * -# UpdateEnergyArray()\n
695 * Shift archived energy values to make space for next one.
696 * -# UpdateDensityCalculation(), SpeedMeasure()'d in DensityTime\n
697 * Calculate TotalLocalDensity of LocalPsis and gather results as TotalDensity.
698 * -# UpdatePsiEnergyCalculation()\n
699 * Calculate kinetic and non-local energy contributons from the Psis.
700 * -# CalculateDensityEnergy()\n
701 * Calculate remaining energy contributions from the Density and adds \f$V_xc\f$ onto DensityTypes#HGDensity.
702 * -# CalculateIonsEnergy()\n
703 * Calculate the Gauss self energy of the Ions.
704 * -# EnergyAllReduce()\n
705 * Gather PsiEnergy results from all processes and sum up together with all other contributions to TotalEnergy.
706 * -# CheckCPULIM()\n
707 * Check if external signal has been received (e.g. end of time slit on cluster), break operation at next possible moment.
708 * -# CalculateMinimumStop()\n
709 * Evaluates stop condition if desired precision or steps or ... have been reached. Otherwise go to
710 * CalculateNewWave().
711 *
712 * Before return orthonormality is tested.
713 * \param *P Problem at hand
714 * \param *Stop flag to determine if epsilon stop conditions have met
715 * \param *SuperStop flag to determinte whether external signal's required end of calculations
716 */
717static void MinimiseOccupied(struct Problem *P, int *Stop, int *SuperStop)
718{
719 struct RunStruct *R = &P->R;
720 struct Lattice *Lat = &P->Lat;
721 struct Psis *Psi = &Lat->Psi;
722 //struct FileData *F = &P->Files;
723// int i;
724// double norm;
725 //double dEdt0,ddEddt0,HartreeddEddt0,XCddEddt0, d[4], D[4],ConDirHConDir;
726 struct LatticeLevel *LevS = R->LevS;
727 int ElementSize = (sizeof(fftw_complex) / sizeof(double));
728 int iter = 0, status, max_iter=10;
729 const gsl_min_fminimizer_type *T;
730 gsl_min_fminimizer *s;
731 double m, a, b;
732 double f_m = 0., f_a, f_b;
733 double dcos, dsin;
734 int g;
735 fftw_complex *ConDir = P->Grad.GradientArray[ConDirGradient];
736 fftw_complex *source = NULL, *oldsource = NULL;
737 gsl_function F;
738 F.function = &fn1;
739 F.params = (void *) P;
740 T = gsl_min_fminimizer_brent;
741 s = gsl_min_fminimizer_alloc (T);
742 int DoBrent, StartLocalPsiNo;
743
744 ResetBrent(P,Psi);
745 *Stop = 0;
746 if (P->Call.ReadSrcFiles) {
747 if (!ReadSrcPsiDensity(P,Occupied,1, R->LevSNo)) { // if file for level exists and desired, read from file
748 P->Call.ReadSrcFiles = 0; // -r was bogus, remove it, have to start anew
749 if(P->Call.out[MinOut]) fprintf(stderr,"(%i) Re-initializing, files are missing/corrupted...\n", P->Par.me);
750 InitPsisValue(P, Psi->TypeStartIndex[Occupied], Psi->TypeStartIndex[Occupied+1]); // initialize perturbed array for this run
751 ResetGramSchTagType(P, Psi, Occupied, NotOrthogonal); // loaded values are orthonormal
752 SpeedMeasure(P, InitGramSchTime, StartTimeDo);
753 GramSch(P, R->LevS, Psi, Orthonormalize);
754 SpeedMeasure(P, InitGramSchTime, StopTimeDo);
755 } else {
756 SpeedMeasure(P, InitSimTime, StartTimeDo);
757 if(P->Call.out[MinOut]) fprintf(stderr,"(%i) Reading from file...\n", P->Par.me);
758 ReadSrcPsiDensity(P, Occupied, 0, R->LevSNo);
759 ResetGramSchTagType(P, Psi, Occupied, IsOrthonormal); // loaded values are orthonormal
760 }
761 SpeedMeasure(P, InitDensityTime, StartTimeDo);
762 InitDensityCalculation(P);
763 SpeedMeasure(P, InitDensityTime, StopTimeDo);
764 InitPsiEnergyCalculation(P, Occupied); // go through all orbitals calculating kinetic and non-local
765 StartLocalPsiNo = R->ActualLocalPsiNo;
766 do { // otherwise OnePsiElementAddData#Lambda is calculated only for current Psi not for all
767 CalculateDensityEnergy(P, 0);
768 UpdateActualPsiNo(P, Occupied);
769 } while (R->ActualLocalPsiNo != StartLocalPsiNo);
770 CalculateIonsEnergy(P);
771 EnergyAllReduce(P);
772 SpeedMeasure(P, InitSimTime, StopTimeDo);
773 R->LevS->Step++;
774 EnergyOutput(P,0);
775 }
776 if (P->Call.ReadSrcFiles != 1) { // otherwise minimise oneself
777 if(P->Call.out[LeaderOut]) fprintf(stderr,"(%i)Beginning minimisation of type %s ...\n", P->Par.me, R->MinimisationName[Occupied]);
778 while (*Stop != 1) { // loop testing condition over all Psis
779 // in the following loop, we have two cases:
780 // 1) still far away and just guessing: Use the normal CalculateNewWave() to improve Psi
781 // 2) closer (DoBrent=-1): use brent line search instead
782 // and due to these two cases, we also have two ifs inside each in order to catch stepping from one case
783 // to the other - due to decreasing DoBrent and/or stepping to the next Psi (which may not yet be DoBrent==1)
784
785 // case 1)
786 if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent != 1) {
787 //SetArrayToDouble0((double *)LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo],LevS->MaxG*2);
788 memcpy(LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
789 //fprintf(stderr,"(%i) Psi %i at %p stored in OldPsi at %p: Old[0] %lg+i%lg\n", P->Par.me, R->ActualLocalPsiNo, LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].re, LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].im);
790 f_m = P->Lat.E->TotalEnergy[0]; // grab first value
791 m = 0.;
792 CalculateNewWave(P,NULL);
793 if ((R->DoBrent == 1) && (fabs(Lat->E->delta[0]) < M_PI/4.))
794 Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent--;
795 if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent != 1) {
796 UpdateActualPsiNo(P, Occupied);
797 UpdateEnergyArray(P);
798 CalculateEnergy(P); // just to get a sensible delta
799 if ((R->ActualLocalPsiNo != R->OldActualLocalPsiNo) && (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent == 1)) {
800 R->OldActualLocalPsiNo = R->ActualLocalPsiNo;
801 // if we stepped on to a new Psi, which is already down at DoBrent=1 unlike the last one,
802 // then an up-to-date gradient is missing for the following Brent line search
803 if(P->Call.out[LeaderOut]) fprintf(stderr,"(%i) We stepped on to a new Psi, which is already in the Brent regime ...re-calc delta\n", P->Par.me);
804 memcpy(LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
805 //fprintf(stderr,"(%i) Psi %i at %p stored in OldPsi at %p: Old[0] %lg+i%lg\n", P->Par.me, R->ActualLocalPsiNo, LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].re, LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].im);
806 f_m = P->Lat.E->TotalEnergy[0]; // grab first value
807 m = 0.;
808 DoBrent = Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent;
809 Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent = 2;
810 CalculateNewWave(P,NULL);
811 Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent = DoBrent;
812 }
813 //fprintf(stderr,"(%i) fnl(%lg) = %lg\n", P->Par.me, m, f_m);
814 }
815 }
816
817 // case 2)
818 if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent == 1) {
819 R->PsiStep=R->MaxPsiStep; // no more fresh gradients from this point for current ActualLocalPsiNo
820 a = b = 0.5*fabs(Lat->E->delta[0]);
821 // we have a meaningful first minimum guess from above CalculateNewWave() resp. from end of this if of last step: Lat->E->delta[0]
822 source = LevS->LPsi->LocalPsi[R->ActualLocalPsiNo];
823 oldsource = LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo];
824 //SetArrayToDouble0((double *)source,LevS->MaxG*2);
825 do {
826 a -= fabs(Lat->E->delta[0]) == 0 ? 0.1 : fabs(Lat->E->delta[0]);
827 if (a < -M_PI/2.) a = -M_PI/2.;// for this to work we need the pre-estimation which leads us into a nice regime (without gradient being the better _initial_ guess for a Psi)
828 dcos = cos(a);
829 dsin = sin(a);
830 for (g = 0; g < LevS->MaxG; g++) { // Here all coefficients are updated for the new found wave function
831 //if (isnan(ConDir[g].re)) { fprintf(stderr,"WARNGING: CalculateLineSearch(): ConDir_%i(%i) = NaN!\n", R->ActualLocalPsiNo, g); Error(SomeError, "NaN-Fehler!"); }
832 c_re(source[g]) = c_re(oldsource[g])*dcos + c_re(ConDir[g])*dsin;
833 c_im(source[g]) = c_im(oldsource[g])*dcos + c_im(ConDir[g])*dsin;
834 }
835 CalculateEnergy(P);
836 f_a = P->Lat.E->TotalEnergy[0]; // grab second value at left border
837 //fprintf(stderr,"(%i) fnl(%lg) = %lg, Check ConDir[0] = %lg+i%lg, source[0] = %lg+i%lg, oldsource[0] = %lg+i%lg, TotDens[0] = %lg\n", P->Par.me, a, f_a, ConDir[0].re, ConDir[0].im, source[0].re, source[0].im, oldsource[0].re, oldsource[0].im, R->Lev0->Dens->DensityArray[TotalDensity][0]);
838 } while (f_a < f_m);
839
840 //SetArrayToDouble0((double *)source,LevS->MaxG*2);
841 do {
842 b += fabs(Lat->E->delta[0]) == 0 ? 0.1 : fabs(Lat->E->delta[0]);
843 if (b > M_PI/2.) b = M_PI/2.;
844 dcos = cos(b);
845 dsin = sin(b);
846 for (g = 0; g < LevS->MaxG; g++) { // Here all coefficients are updated for the new found wave function
847 //if (isnan(ConDir[g].re)) { fprintf(stderr,"WARNGING: CalculateLineSearch(): ConDir_%i(%i) = NaN!\n", R->ActualLocalPsiNo, g); Error(SomeError, "NaN-Fehler!"); }
848 c_re(source[g]) = c_re(oldsource[g])*dcos + c_re(ConDir[g])*dsin;
849 c_im(source[g]) = c_im(oldsource[g])*dcos + c_im(ConDir[g])*dsin;
850 }
851 CalculateEnergy(P);
852 f_b = P->Lat.E->TotalEnergy[0]; // grab second value at left border
853 //fprintf(stderr,"(%i) fnl(%lg) = %lg\n", P->Par.me, b, f_b);
854 } while (f_b < f_m);
855
856 memcpy(source, oldsource, ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
857 //fprintf(stderr,"(%i) Psi %i at %p retrieved from OldPsi at %p: Old[0] %lg+i%lg\n", P->Par.me, R->ActualLocalPsiNo, LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].re, LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].im);
858 CalculateEnergy(P);
859
860 if(P->Call.out[ValueOut]) fprintf(stderr,"(%i) Preparing brent with f(a) (%lg,%lg)\t f(b) (%lg,%lg)\t f(m) (%lg,%lg) ...\n", P->Par.me,a,f_a,b,f_b,m,f_m);
861 iter=0;
862 gsl_min_fminimizer_set_with_values (s, &F, m, f_m, a, f_a, b, f_b);
863 if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) using %s method\n",P->Par.me, gsl_min_fminimizer_name (s));
864 if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) %5s [%9s, %9s] %9s %9s\n",P->Par.me, "iter", "lower", "upper", "min", "err(est)");
865 if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) %5d [%.7f, %.7f] %.7f %.7f\n",P->Par.me, iter, a, b, m, b - a);
866 do {
867 iter++;
868 status = gsl_min_fminimizer_iterate (s);
869
870 m = gsl_min_fminimizer_x_minimum (s);
871 a = gsl_min_fminimizer_x_lower (s);
872 b = gsl_min_fminimizer_x_upper (s);
873
874 status = gsl_min_test_interval (a, b, 0.001, 0.0);
875
876 if (status == GSL_SUCCESS)
877 if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) Converged:\n",P->Par.me);
878
879 if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) %5d [%.7f, %.7f] %.7f %.7f\n",P->Par.me,
880 iter, a, b, m, b - a);
881 } while (status == GSL_CONTINUE && iter < max_iter);
882 CalculateNewWave(P,&m);
883 TestGramSch(P,LevS,Psi,Occupied);
884 UpdateActualPsiNo(P, Occupied); // step on due setting to MaxPsiStep further above
885 UpdateEnergyArray(P);
886 CalculateEnergy(P);
887 //fprintf(stderr,"(%i) Final value for Psi %i: %lg\n", P->Par.me, R->ActualLocalPsiNo, P->Lat.E->TotalEnergy[0]);
888 R->MinStopStep = R->ActualMaxMinStopStep; // check stop condition every time
889 if (*SuperStop != 1)
890 *SuperStop = CheckCPULIM(P);
891 *Stop = CalculateMinimumStop(P, *SuperStop);
892 R->OldActualLocalPsiNo = R->ActualLocalPsiNo;
893 if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent == 1) { // new wave function means new gradient!
894 DoBrent = Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent;
895 Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent = 2;
896 //SetArrayToDouble0((double *)LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo],LevS->MaxG*2);
897 memcpy(LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
898 //fprintf(stderr,"(%i) Psi %i at %p stored in OldPsi at %p: Old[0] %lg+i%lg\n", P->Par.me, R->ActualLocalPsiNo, LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].re, LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo][0].im);
899 f_m = P->Lat.E->TotalEnergy[0]; // grab first value
900 m = 0.;
901 CalculateNewWave(P,NULL);
902 Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent = DoBrent;
903 }
904 }
905
906 if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent != 1) { // otherwise the following checks eliminiate stop=1 from above
907 if (*SuperStop != 1)
908 *SuperStop = CheckCPULIM(P);
909 *Stop = CalculateMinimumStop(P, *SuperStop);
910 }
911 /*EnergyOutput(P, Stop);*/
912 P->Speed.Steps++;
913 R->LevS->Step++;
914 /*ControlNativeDensity(P);*/
915 //fprintf(stderr,"(%i) Stop %i\n",P->Par.me, Stop);
916 }
917 //OutputVisSrcFiles(P, Occupied); // is now done after localization (ComputeMLWF())
918 }
919 TestGramSch(P,R->LevS,Psi, Occupied);
920}
921
922/** Minimisation of the PsiTagType#UnOccupied orbitals in the field of the occupied ones.
923 * Assuming RunStruct#ActualLocalPsiNo is currenlty still an occupied wave function, we stop onward to the first
924 * unoccupied and reset RunStruct#OldActualLocalPsiNo. Then it is checked whether CallOptions#ReadSrcFiles is set
925 * and thus coefficients for the level have to be read from file and afterwards initialized.
926 *
927 * Then follows the main loop, until a stop condition is met:
928 * -# CalculateNewWave()\n
929 * Over a conjugate gradient method the next (minimal) wave function is sought for.
930 * -# UpdateActualPsiNo()\n
931 * Switch local Psi to next one.
932 * -# UpdateEnergyArray()\n
933 * Shift archived energy values to make space for next one.
934 * -# UpdateDensityCalculation(), SpeedMeasure()'d in DensityTime\n
935 * Calculate TotalLocalDensity of LocalPsis and gather results as TotalDensity.
936 * -# UpdatePsiEnergyCalculation()\n
937 * Calculate kinetic and non-local energy contributons from the Psis.
938 * -# CalculateGapEnergy()\n
939 * Calculate Gap energies (Hartreepotential, Pseudo) and the gradient.
940 * -# EnergyAllReduce()\n
941 * Gather PsiEnergy results from all processes and sum up together with all other contributions to TotalEnergy.
942 * -# CheckCPULIM()\n
943 * Check if external signal has been received (e.g. end of time slit on cluster), break operation at next possible moment.
944 * -# CalculateMinimumStop()\n
945 * Evaluates stop condition if desired precision or steps or ... have been reached. Otherwise go to
946 * CalculateNewWave().
947 *
948 * Afterwards, the coefficients are written to file by OutputVisSrcFiles() if desired. Orthonormality is tested, we step
949 * back to the occupied wave functions and the densities are re-initialized.
950 * \param *P Problem at hand
951 * \param *Stop flag to determine if epsilon stop conditions have met
952 * \param *SuperStop flag to determinte whether external signal's required end of calculations
953 */
954static void MinimiseUnoccupied (struct Problem *P, int *Stop, int *SuperStop) {
955 struct RunStruct *R = &P->R;
956 struct Lattice *Lat = &P->Lat;
957 struct Psis *Psi = &Lat->Psi;
958 int StartLocalPsiNo;
959
960 *Stop = 0;
961 R->PsiStep = R->MaxPsiStep; // in case it's zero from CalculateForce()
962 UpdateActualPsiNo(P, UnOccupied); // step on to next unoccupied one
963 R->OldActualLocalPsiNo = R->ActualLocalPsiNo; // reset, otherwise OldActualLocalPsiNo still points to occupied wave function
964 UpdateGramSchOldActualPsiNo(P,Psi);
965 if (P->Call.ReadSrcFiles && ReadSrcPsiDensity(P,UnOccupied,1, R->LevSNo)) {
966 SpeedMeasure(P, InitSimTime, StartTimeDo);
967 if(P->Call.out[MinOut]) fprintf(stderr, "(%i) Re-initializing %s psi array from source file of recent calculation\n", P->Par.me, R->MinimisationName[R->CurrentMin]);
968 ReadSrcPsiDensity(P, UnOccupied, 0, R->LevSNo);
969 if (P->Call.ReadSrcFiles != 2) {
970 ResetGramSchTagType(P, Psi, UnOccupied, IsOrthonormal); // loaded values are orthonormal
971 SpeedMeasure(P, DensityTime, StartTimeDo);
972 InitDensityCalculation(P);
973 SpeedMeasure(P, DensityTime, StopTimeDo);
974 InitPsiEnergyCalculation(P,UnOccupied); // go through all orbitals calculating kinetic and non-local
975 //CalculateDensityEnergy(P, 0);
976 StartLocalPsiNo = R->ActualLocalPsiNo;
977 do { // otherwise OnePsiElementAddData#Lambda is calculated only for current Psi not for all
978 CalculateGapEnergy(P);
979 UpdateActualPsiNo(P, Occupied);
980 } while (R->ActualLocalPsiNo != StartLocalPsiNo);
981 EnergyAllReduce(P);
982 }
983 SpeedMeasure(P, InitSimTime, StopTimeDo);
984 }
985 if (P->Call.ReadSrcFiles != 1) {
986 SpeedMeasure(P, InitSimTime, StartTimeDo);
987 ResetGramSchTagType(P, Psi, UnOccupied, NotOrthogonal);
988 SpeedMeasure(P, GramSchTime, StartTimeDo);
989 GramSch(P, R->LevS, Psi, Orthonormalize);
990 SpeedMeasure(P, GramSchTime, StopTimeDo);
991 SpeedMeasure(P, InitDensityTime, StartTimeDo);
992 InitDensityCalculation(P);
993 SpeedMeasure(P, InitDensityTime, StopTimeDo);
994 InitPsiEnergyCalculation(P,UnOccupied); // go through all orbitals calculating kinetic and non-local
995 //CalculateDensityEnergy(P, 0);
996 CalculateGapEnergy(P);
997 EnergyAllReduce(P);
998 SpeedMeasure(P, InitSimTime, StopTimeDo);
999 R->LevS->Step++;
1000 EnergyOutput(P,0);
1001 if(P->Call.out[LeaderOut]) fprintf(stderr,"(%i)Beginning minimisation of type %s ...\n", P->Par.me, R->MinimisationName[R->CurrentMin]);
1002 while (*Stop != 1) {
1003 CalculateNewWave(P,NULL);
1004 UpdateActualPsiNo(P, UnOccupied);
1005 /* New */
1006 UpdateEnergyArray(P);
1007 SpeedMeasure(P, DensityTime, StartTimeDo);
1008 UpdateDensityCalculation(P);
1009 SpeedMeasure(P, DensityTime, StopTimeDo);
1010 UpdatePsiEnergyCalculation(P);
1011 //CalculateDensityEnergy(P, 0);
1012 CalculateGapEnergy(P); //calculates XC, HGDensity, afterwards gradient, where V_xc is added upon HGDensity
1013 EnergyAllReduce(P);
1014 if (*SuperStop != 1)
1015 *SuperStop = CheckCPULIM(P);
1016 *Stop = CalculateMinimumStop(P, *SuperStop);
1017 /*EnergyOutput(P, Stop);*/
1018 P->Speed.Steps++;
1019 R->LevS->Step++;
1020 /*ControlNativeDensity(P);*/
1021 }
1022 OutputVisSrcFiles(P, UnOccupied);
1023// if (!TestReadnWriteSrcDensity(P,UnOccupied))
1024// Error(SomeError,"TestReadnWriteSrcDensity failed!");
1025 }
1026 TestGramSch(P,R->LevS,Psi, UnOccupied);
1027 ResetGramSchTagType(P, Psi, UnOccupied, NotUsedToOrtho);
1028 // re-calculate Occupied values (in preparation for perturbed ones)
1029 UpdateActualPsiNo(P, Occupied); // step on to next occupied one
1030 SpeedMeasure(P, DensityTime, StartTimeDo);
1031 InitDensityCalculation(P); // re-init densities to occupied standard
1032 SpeedMeasure(P, DensityTime, StopTimeDo);
1033// InitPsiEnergyCalculation(P,Occupied);
1034// CalculateDensityEnergy(P, 0);
1035// EnergyAllReduce(P);
1036}
1037
1038
1039/** Calculate the forces.
1040 * From RunStruct::LevSNo downto RunStruct::InitLevSNo the following routines are called in a loop:
1041 * -# In case of RunStruct#DoSeparated another loop begins for the unoccupied states with some reinitalization
1042 * before and afterwards. The loop however is much the same as the one above.
1043 * -# ChangeToLevUp()\n
1044 * Repeat the loop or when the above stop is reached, the level is changed and the loop repeated.
1045 *
1046 * Afterwards comes the actual force and energy calculation by calling:
1047 * -# ControlNativeDensity()\n
1048 * Checks if the density still reproduces particle number.
1049 * -# CalculateIonLocalForce(), SpeedMeasure()'d in LocFTime\n
1050 * Calculale local part of force acting on Ions.
1051 * -# CalculateIonNonLocalForce(), SpeedMeasure()'d in NonLocFTime\n
1052 * Calculale local part of force acting on Ions.
1053 * -# CalculateEwald()\n
1054 * Calculate Ewald force acting on Ions.
1055 * -# CalculateIonForce()\n
1056 * Sum up those three contributions.
1057 * -# CorrectForces()\n
1058 * Shifts center of gravity of all forces for each Ion, so that the cell itself remains at rest.
1059 * -# GetOuterStop()
1060 * Calculates a mean force per Ion.
1061 * \param *P Problem at hand
1062 * \return 1 - cpulim received, break operation, 0 - continue as normal
1063 */
1064int CalculateForce(struct Problem *P)
1065{
1066 struct RunStruct *R = &P->R;
1067 struct Lattice *Lat = &P->Lat;
1068 struct Psis *Psi = &Lat->Psi;
1069 struct LatticeLevel *LevS = R->LevS;
1070 struct FileData *F = &P->Files;
1071 struct Ions *I = &P->Ion;
1072 int Stop=0, SuperStop = 0, OuterStop = 0;
1073 //int i, j;
1074 while ((R->LevSNo > R->InitLevSNo) || (!Stop && R->LevSNo == R->InitLevSNo)) {
1075 // occupied
1076 R->PsiStep = R->MaxPsiStep; // reset in-Psi-minimisation-counter, so that we really advance to the next wave function
1077 R->OldActualLocalPsiNo = R->ActualLocalPsiNo; // reset OldActualLocalPsiNo, as it might still point to a perturbed wave function from last level
1078 UpdateGramSchOldActualPsiNo(P,Psi);
1079 MinimiseOccupied(P, &Stop, &SuperStop);
1080 if (!I->StructOpt) {
1081 if ((P->Call.ReadSrcFiles != 1) || (!ParseWannierFile(P))) { // only localize and store if they have just been minimised (hence don't come solely from file), otherwise read stored values from file (if possible)
1082 SpeedMeasure(P, WannierTime, StartTimeDo);
1083 ComputeMLWF(P); // localization of orbitals
1084 SpeedMeasure(P, WannierTime, StopTimeDo);
1085 OutputVisSrcFiles(P, Occupied); // rewrite now localized orbitals
1086 // if (!TestReadnWriteSrcDensity(P,Occupied))
1087 // Error(SomeError,"TestReadnWriteSrcDensity failed!");
1088 }
1089
1090// // plot psi cuts
1091// for (i=0; i < Psi->MaxPsiOfType; i++) // go through all wave functions (here without the extra ones for each process)
1092// if ((Psi->AllPsiStatus[i].PsiType == Occupied) && (Psi->AllPsiStatus[i].my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi))
1093// for (j=0;j<NDIM;j++) {
1094// //fprintf(stderr,"(%i) Plotting Psi %i/%i cut axis %i at coordinate %lg \n",P->Par.me, i, Psi->AllPsiStatus[i].MyGlobalNo, j, Lat->Psi.AddData[Psi->AllPsiStatus[i].MyLocalNo].WannierCentre[j]);
1095// CalculateOneDensityR(Lat, R->LevS, R->Lev0->Dens, R->LevS->LPsi->LocalPsi[Psi->AllPsiStatus[i].MyLocalNo], R->Lev0->Dens->DensityArray[ActualDensity], R->FactorDensityR, 0);
1096// PlotSrcPlane(P, j, Lat->Psi.AddData[Psi->AllPsiStatus[i].MyLocalNo].WannierCentre[j], Psi->AllPsiStatus[i].MyGlobalNo, R->Lev0->Dens->DensityArray[ActualDensity]);
1097// }
1098
1099 // unoccupied calc
1100 if (R->DoUnOccupied) {
1101 MinimiseUnoccupied(P, &Stop, &SuperStop);
1102 }
1103 // hamiltonian
1104 CalculateHamiltonian(P); // lambda_{kl} needed (and for bandgap after UnOccupied)
1105
1106 // perturbed calc
1107 if ((R->DoPerturbation)) { // && R->LevSNo <= R->InitLevSNo) {
1108 AllocCurrentDensity(R->Lev0->Dens);// lock current density arrays
1109 MinimisePerturbed(P, &Stop, &SuperStop); // herein InitDensityCalculation() is called, thus no need to call it beforehand
1110
1111 SpeedMeasure(P, CurrDensTime, StartTimeDo);
1112 if (SuperStop != 1) {
1113 if ((R->DoFullCurrent == 1) || ((R->DoFullCurrent == 2) && (CheckOrbitalOverlap(P) == 1))) { //test to check whether orbitals have mutual overlap and thus \\DeltaJ_{xc} must not be dropped
1114 R->DoFullCurrent = 1; // set to 1 if it was 2 but Check...() yielded necessity
1115 debug(P,"Filling with Delta j ...");
1116 FillDeltaCurrentDensity(P);
1117 }// else
1118 //debug(P,"There is no overlap between orbitals.");
1119 //debug(P,"Filling with j ...");
1120 //FillCurrentDensity(P);
1121 }
1122 SpeedMeasure(P, CurrDensTime, StopTimeDo);
1123 TestCurrent(P,0);
1124 TestCurrent(P,1);
1125 TestCurrent(P,2);
1126 if (F->DoOutCurr) {
1127 debug(P,"OutputCurrentDensity");
1128 OutputCurrentDensity(P);
1129 }
1130 if (R->VectorPlane != -1) {
1131 debug(P,"PlotVectorPlane");
1132 PlotVectorPlane(P,R->VectorPlane,R->VectorCut);
1133 }
1134 fprintf(stderr,"(%i) ECut [L%i] = %e Ht\n", P->Par.me, R->Lev0->LevelNo, pow(2*M_PI*M_PI/Lat->Volume*R->Lev0->TotalAllMaxG, 2./3.));
1135 debug(P,"Calculation of magnetic susceptibility");
1136 CalculateMagneticSusceptibility(P);
1137 debug(P,"Normal calculation of shielding over R-space");
1138 CalculateChemicalShielding(P);
1139 debug(P,"Reciprocal calculation of shielding over G-space");
1140 CalculateChemicalShieldingByReciprocalCurrentDensity(P);
1141 SpeedMeasure(P, CurrDensTime, StopTimeDo);
1142 DisAllocCurrentDensity(R->Lev0->Dens); // unlock current density arrays
1143 } else {
1144 InitDensityCalculation(P); // all unperturbed(!) wave functions've "changed" from ComputeMLWF(), thus reinit density
1145 }
1146 //fprintf(stderr,"(%i) DoubleG: %p, CArray[22]: %p, OldLocalPsi: %p\n", P->Par.me, R->LevS->DoubleG, R->Lev0->Dens->DensityCArray[22], R->LevS->LPsi->OldLocalPsi);
1147 }
1148
1149// if (!I->StructOpt && R->DoPerturbation) {
1150// InitDensityCalculation(P); // most of the density array were used during FillCurrentDensity(), thus reinit density
1151// }
1152
1153 // next level
1154 ChangeToLevUp(P, &Stop);
1155 //if (isnan(LevS->LPsi->LocalPsi[R->ActualLocalPsiNo][0].re)) { fprintf(stderr,"(%i) WARNING in ChangeToLevUp(): LPsi->LocalPsi[%i]_[%i] = NaN!\n", P->Par.me, R->ActualLocalPsiNo, 0); Error(SomeError, "NaN-Fehler!"); }
1156 LevS = R->LevS; // re-set pointer that's changed from LevUp
1157 }
1158 //fprintf(stderr,"(%i) DoubleG: %p, CArray[22]: %p, OldLocalPsi: %p\n", P->Par.me, R->LevS->DoubleG, R->Lev0->Dens->DensityCArray[22], R->LevS->LPsi->OldLocalPsi);
1159 // necessary for correct ionic forces ...
1160 SpeedMeasure(P, LocFTime, StartTimeDo);
1161 CalculateIonLocalForce(P);
1162 SpeedMeasure(P, LocFTime, StopTimeDo);
1163 SpeedMeasure(P, NonLocFTime, StartTimeDo);
1164 CalculateIonNonLocalForce(P);
1165 SpeedMeasure(P, NonLocFTime, StopTimeDo);
1166 CalculateEwald(P, 0);
1167 CalculateIonForce(P);
1168 CorrectForces(P);
1169 // ... on output of densities
1170 if (F->DoOutOrbitals) { // output of each orbital
1171 debug(P,"OutputVisAllOrbital");
1172 OutputVisAllOrbital(P,0,1,Occupied);
1173 }
1174
1175 OutputNorm(stderr, P);
1176 //fprintf(stderr,"(%i) DoubleG: %p, CArray[22]: %p, OldLocalPsi: %p\n", P->Par.me, R->LevS->DoubleG, R->Lev0->Dens->DensityCArray[22], R->LevS->LPsi->OldLocalPsi);
1177 OutputVis(P, P->R.Lev0->Dens->DensityArray[TotalDensity]);
1178 TestGramSch(P,LevS,Psi, -1);
1179 SpeedMeasure(P, SimTime, StopTimeDo);
1180 /*TestGramSch(P, R->LevS, &P->Lat.Psi); */
1181 ControlNativeDensity(P);
1182 SpeedMeasure(P, LocFTime, StartTimeDo);
1183 CalculateIonLocalForce(P);
1184 SpeedMeasure(P, LocFTime, StopTimeDo);
1185 SpeedMeasure(P, NonLocFTime, StartTimeDo);
1186 CalculateIonNonLocalForce(P);
1187 SpeedMeasure(P, NonLocFTime, StopTimeDo);
1188 CalculateEwald(P, 0);
1189 CalculateIonForce(P);
1190 CorrectForces(P);
1191 GetOuterStop(P);
1192 //fprintf(stderr,"(%i) DoubleG: %p, CArray[22]: %p, OldLocalPsi: %p\n", P->Par.me, R->LevS->DoubleG, R->Lev0->Dens->DensityCArray[22], R->LevS->LPsi->OldLocalPsi);
1193 if (SuperStop) OuterStop = 1;
1194 return OuterStop;
1195}
1196
1197/** Wrapper for CalculateForce() for simplex minimisation of ionic forces.
1198 * \param *v vector with degrees of freedom
1199 * \param *params additional arguments, here Problem at hand
1200 */
1201double StructOpt_f(const gsl_vector *v, void *params)
1202{
1203 struct Problem *P = (struct Problem *)params;
1204 struct RunStruct *R = &P->R;
1205 struct Ions *I = &P->Ion;
1206 struct Energy *E = P->Lat.E;
1207 int i;
1208 double *R_ion, *R_old, *R_old_old, *FIon;
1209 double norm = 0.;
1210 int is,ia,k,index;
1211 int OuterStop;
1212 // update ion positions from vector coordinates
1213 index=0;
1214 for (is=0; is < I->Max_Types; is++) // for all elements
1215 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) { // for all ions of element
1216 R_ion = &I->I[is].R[NDIM*ia];
1217 R_old = &I->I[is].R_old[NDIM*ia];
1218 R_old_old = &I->I[is].R_old_old[NDIM*ia];
1219 for (k=0;k<NDIM;k++) { // for all dimensions
1220 R_old_old[k] = R_old[k];
1221 R_old[k] = R_ion[k];
1222 R_ion[k] = gsl_vector_get (v, index++);
1223 }
1224 }
1225 // recalculate ionic forces (do electronic minimisation)
1226 R->OuterStep++;
1227 if (P->Call.out[NormalOut]) fprintf(stderr,"(%i) Commencing Fletcher-Reeves step %i ... \n",P->Par.me, R->OuterStep);
1228 R->NewRStep++;
1229 OutputIonCoordinates(P);
1230 UpdateWaveAfterIonMove(P);
1231 for (i=MAXOLD-1; i > 0; i--) // store away old energies
1232 E->TotalEnergyOuter[i] = E->TotalEnergyOuter[i-1];
1233 UpdateToNewWaves(P);
1234 E->TotalEnergyOuter[0] = E->TotalEnergy[0];
1235 OuterStop = CalculateForce(P);
1236 UpdateIonsU(P);
1237 CorrectVelocity(P);
1238 CalculateEnergyIonsU(P);
1239/* if ((P->R.ScaleTempStep > 0) && ((R->OuterStep-1) % P->R.ScaleTempStep == 0))
1240 ScaleTemp(P);*/
1241 if ((R->OuterStep-1) % P->R.OutSrcStep == 0)
1242 OutputVisSrcFiles(P, Occupied);
1243 if ((R->OuterStep-1) % P->R.OutVisStep == 0) {
1244/* // recalculate density for the specific wave function ...
1245 CalculateOneDensityR(Lat, LevS, Dens0, PsiDat, Dens0->DensityArray[ActualDensity], R->FactorDensityR, 0);
1246 // ... and output (wherein ActualDensity is used instead of TotalDensity)
1247 OutputVis(P);
1248 OutputIonForce(P);
1249 EnergyOutput(P, 1);*/
1250 }
1251 // sum up mean force
1252 for (is=0; is < I->Max_Types; is++)
1253 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1254 FIon = &I->I[is].FIon[NDIM*ia];
1255 norm += sqrt(RSP3(FIon,FIon));
1256 }
1257 if (P->Par.me == 0) fprintf(stderr,"(%i) Mean Force over all Ions %e\n",P->Par.me, norm);
1258 return norm;
1259}
1260
1261void StructOpt_df(const gsl_vector *v, void *params, gsl_vector *df)
1262{
1263 struct Problem *P = (struct Problem *)params;
1264 struct Ions *I = &P->Ion;
1265 double *FIon;
1266 int is,ia,k, index=0;
1267 for (is=0; is < I->Max_Types; is++) // for all elements
1268 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) { // for all ions of element
1269 FIon = &I->I[is].FIon[NDIM*ia];
1270 for (k=0;k<NDIM;k++) { // for all dimensions
1271 gsl_vector_set (df, index++, FIon[k]);
1272 }
1273 }
1274}
1275
1276void StructOpt_fdf (const gsl_vector *x, void *params, double *f, gsl_vector *df)
1277{
1278 *f = StructOpt_f(x, params);
1279 StructOpt_df(x, params, df);
1280}
1281
1282
1283/** CG implementation for the structure optimization.
1284 * We follow the example from the GSL manual.
1285 * \param *P Problem at hand
1286 */
1287void UpdateIon_PRCG(struct Problem *P)
1288{
1289 struct RunStruct *Run = &P->R;
1290 struct Ions *I = &P->Ion;
1291 size_t np = NDIM*I->Max_TotalIons; // d.o.f = number of ions times number of dimensions
1292 int is, ia, k, index;
1293 double *R;
1294
1295 const gsl_multimin_fdfminimizer_type *T;
1296 gsl_multimin_fdfminimizer *s;
1297 gsl_vector *x;
1298 gsl_multimin_function_fdf minex_func;
1299
1300 size_t iter = 0;
1301 int status;
1302
1303 /* Starting point */
1304 x = gsl_vector_alloc (np);
1305
1306 index=0;
1307 for (is=0; is < I->Max_Types; is++) // for all elements
1308 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) { // for all ions of element
1309 R = &I->I[is].R[NDIM*ia];
1310 for (k=0;k<NDIM;k++) // for all dimensions
1311 gsl_vector_set (x, index++, R[k]);
1312 }
1313
1314 /* Initialize method and iterate */
1315 minex_func.f = &StructOpt_f;
1316 minex_func.df = &StructOpt_df;
1317 minex_func.fdf = &StructOpt_fdf;
1318 minex_func.n = np;
1319 minex_func.params = (void *)P;
1320
1321 T = gsl_multimin_fdfminimizer_conjugate_pr;
1322 s = gsl_multimin_fdfminimizer_alloc (T, np);
1323
1324 gsl_multimin_fdfminimizer_set (s, &minex_func, x, 0.1, 1e-4);
1325
1326 do {
1327 iter++;
1328 status = gsl_multimin_fdfminimizer_iterate(s);
1329
1330 if (status)
1331 break;
1332
1333 status = gsl_multimin_test_gradient (s->gradient, 1e-4);
1334
1335 if (status == GSL_SUCCESS)
1336 if (P->Par.me == 0) fprintf (stderr,"(%i) converged to minimum at\n", P->Par.me);
1337
1338 //if (P->Call.out[NormalOut]) fprintf(stderr,"(%i) Commencing '%s' step %i ... \n",P->Par.me, gsl_multimin_fdfminimizer_name(s), P->R.StructOptStep);
1339 if ((P->Call.out[NormalOut]) && (P->Par.me == 0)) fprintf (stderr, "(%i) %5d %10.5f\n", P->Par.me, (int)iter, s->f);
1340 } while (status == GSL_CONTINUE && iter < Run->MaxOuterStep);
1341
1342 gsl_vector_free(x);
1343 gsl_multimin_fdfminimizer_free (s);
1344}
1345
1346/** Implementation of various thermostats.
1347 * All these thermostats apply an additional force which has the following forms:
1348 * -# Woodcock
1349 * \f$p_i \rightarrow \sqrt{\frac{T_0}{T}} \cdot p_i\f$
1350 * -# Gaussian
1351 * \f$ \frac{ \sum_i \frac{p_i}{m_i} \frac{\partial V}{\partial q_i}} {\sum_i \frac{p^2_i}{m_i}} \cdot p_i\f$
1352 * -# Langevin
1353 * \f$p_{i,n} \rightarrow \sqrt{1-\alpha^2} p_{i,0} + \alpha p_r\f$
1354 * -# Berendsen
1355 * \f$p_i \rightarrow \left [ 1+ \frac{\delta t}{\tau_T} \left ( \frac{T_0}{T} \right ) \right ]^{\frac{1}{2}} \cdot p_i\f$
1356 * -# Nose-Hoover
1357 * \f$\zeta p_i \f$ with \f$\frac{\partial \zeta}{\partial t} = \frac{1}{M_s} \left ( \sum^N_{i=1} \frac{p_i^2}{m_i} - g k_B T \right )\f$
1358 * These Thermostats either simply rescale the velocities, thus Thermostats() should be called after UpdateIonsU(), and/or
1359 * have a constraint force acting additionally on the ions. In the latter case, the ion speeds have to be modified
1360 * belatedly and the constraint force set.
1361 * \param *P Problem at hand
1362 * \param i which of the thermostats to take: 0 - none, 1 - Woodcock, 2 - Gaussian, 3 - Langevin, 4 - Berendsen, 5 - Nose-Hoover
1363 * \sa InitThermostat()
1364 */
1365void Thermostats(struct Problem *P, enum thermostats i)
1366{
1367 struct FileData *Files = &P->Files;
1368 struct Ions *I = &P->Ion;
1369 int is, ia, d;
1370 double *U;
1371 double a, ekin = 0.;
1372 double E = 0., F = 0.;
1373 double delta_alpha = 0.;
1374 const int delta_t = P->R.delta_t;
1375 double ScaleTempFactor;
1376 double sigma;
1377 gsl_rng * r;
1378 const gsl_rng_type * T;
1379
1380 // calculate current temperature
1381 CalculateEnergyIonsU(P); // Temperature now in I->ActualTemp
1382 ScaleTempFactor = P->R.TargetTemp/I->ActualTemp;
1383 //if ((P->Par.me == 0) && (I->ActualTemp < MYEPSILON)) fprintf(stderr,"Thermostat: (1) I->ActualTemp = %lg",I->ActualTemp);
1384 if (Files->MeOutMes) fprintf(Files->TemperatureFile, "%d\t%lg",P->R.OuterStep, I->ActualTemp);
1385
1386 // differentating between the various thermostats
1387 switch(i) {
1388 case None:
1389 debug(P, "Applying no thermostat...");
1390 break;
1391 case Woodcock:
1392 if ((P->R.ScaleTempStep > 0) && ((P->R.OuterStep-1) % P->R.ScaleTempStep == 0)) {
1393 debug(P, "Applying Woodcock thermostat...");
1394 for (is=0; is < I->Max_Types; is++) {
1395 a = 0.5*I->I[is].IonMass;
1396 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1397 U = &I->I[is].U[NDIM*ia];
1398 if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
1399 for (d=0; d<NDIM; d++) {
1400 U[d] *= sqrt(ScaleTempFactor);
1401 ekin += 0.5*I->I[is].IonMass * U[d]*U[d];
1402 }
1403 }
1404 }
1405 }
1406 break;
1407 case Gaussian:
1408 debug(P, "Applying Gaussian thermostat...");
1409 for (is=0; is < I->Max_Types; is++) { // sum up constraint constant
1410 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1411 U = &I->I[is].U[NDIM*ia];
1412 if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
1413 for (d=0; d<NDIM; d++) {
1414 F += U[d] * I->I[is].FIon[d+NDIM*ia];
1415 E += U[d]*U[d]*I->I[is].IonMass;
1416 }
1417 }
1418 }
1419 if (P->Call.out[ValueOut]) fprintf(stderr, "(%i) Gaussian Least Constraint constant is %lg\n", P->Par.me, F/E);
1420 for (is=0; is < I->Max_Types; is++) { // apply constraint constant on constraint force and on velocities
1421 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1422 U = &I->I[is].U[NDIM*ia];
1423 if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
1424 for (d=0; d<NDIM; d++) {
1425 I->I[is].FConstraint[d+NDIM*ia] = (F/E) * (U[d]*I->I[is].IonMass);
1426 U[d] += delta_t/I->I[is].IonMass * (I->I[is].FConstraint[d+NDIM*ia]);
1427 ekin += 0.5*I->I[is].IonMass * U[d]*U[d];
1428 }
1429 }
1430 }
1431 break;
1432 case Langevin:
1433 debug(P, "Applying Langevin thermostat...");
1434 // init random number generator
1435 gsl_rng_env_setup();
1436 T = gsl_rng_default;
1437 r = gsl_rng_alloc (T);
1438 // Go through each ion
1439 for (is=0; is < I->Max_Types; is++) {
1440 sigma = sqrt(P->R.TargetTemp/I->I[is].IonMass); // sigma = (k_b T)/m (Hartree/atomicmass = atomiclength/atomictime)
1441 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1442 U = &I->I[is].U[NDIM*ia];
1443 // throw a dice to determine whether it gets hit by a heat bath particle
1444 if (((((rand()/(double)RAND_MAX))*P->R.TempFrequency) < 1.)) { // (I->I[is].IMT[ia] == MoveIon) && even FixedIon moves, only not by other's forces
1445 if (P->Par.me == 0) fprintf(stderr,"(%i) Particle %i,%i was hit (sigma %lg): %lg -> ", P->Par.me, is, ia, sigma, sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]));
1446 // pick three random numbers from a Boltzmann distribution around the desired temperature T for each momenta axis
1447 for (d=0; d<NDIM; d++) {
1448 U[d] = gsl_ran_gaussian (r, sigma);
1449 }
1450 if (P->Par.me == 0) fprintf(stderr,"%lg\n", sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]));
1451 }
1452 for (d=0; d<NDIM; d++)
1453 ekin += 0.5*I->I[is].IonMass * U[d]*U[d];
1454 }
1455 }
1456 break;
1457 case Berendsen:
1458 debug(P, "Applying Berendsen-VanGunsteren thermostat...");
1459 for (is=0; is < I->Max_Types; is++) {
1460 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1461 U = &I->I[is].U[NDIM*ia];
1462 if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
1463 for (d=0; d<NDIM; d++) {
1464 U[d] *= sqrt(1+(P->R.delta_t/P->R.TempFrequency)*(ScaleTempFactor-1));
1465 ekin += 0.5*I->I[is].IonMass * U[d]*U[d];
1466 }
1467 }
1468 }
1469 break;
1470 case NoseHoover:
1471 debug(P, "Applying Nose-Hoover thermostat...");
1472 // dynamically evolve alpha (the additional degree of freedom)
1473 delta_alpha = 0.;
1474 for (is=0; is < I->Max_Types; is++) { // sum up constraint constant
1475 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1476 U = &I->I[is].U[NDIM*ia];
1477 if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
1478 for (d=0; d<NDIM; d++) {
1479 delta_alpha += U[d]*U[d]*I->I[is].IonMass;
1480 }
1481 }
1482 }
1483 delta_alpha = (delta_alpha - (3.*I->Max_TotalIons+1.) * P->R.TargetTemp)/(P->R.HooverMass*Units2Electronmass);
1484 P->R.alpha += delta_alpha*delta_t;
1485 if (P->Par.me == 0) fprintf(stderr,"(%i) alpha = %lg * %i = %lg\n", P->Par.me, delta_alpha, delta_t, P->R.alpha);
1486 // apply updated alpha as additional force
1487 for (is=0; is < I->Max_Types; is++) { // apply constraint constant on constraint force and on velocities
1488 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1489 U = &I->I[is].U[NDIM*ia];
1490 if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
1491 for (d=0; d<NDIM; d++) {
1492 I->I[is].FConstraint[d+NDIM*ia] = - P->R.alpha * (U[d] * I->I[is].IonMass);
1493 U[d] += delta_t/I->I[is].IonMass * (I->I[is].FConstraint[d+NDIM*ia]);
1494 ekin += (0.5*I->I[is].IonMass) * U[d]*U[d];
1495 }
1496 }
1497 }
1498 break;
1499 }
1500 I->EKin = ekin;
1501 I->ActualTemp = (2./(3.*I->Max_TotalIons)*I->EKin);
1502 //if ((P->Par.me == 0) && (I->ActualTemp < MYEPSILON)) fprintf(stderr,"Thermostat: (2) I->ActualTemp = %lg",I->ActualTemp);
1503 if (Files->MeOutMes) { fprintf(Files->TemperatureFile, "\t%lg\n", I->ActualTemp); fflush(Files->TemperatureFile); }
1504}
1505
1506/** Does the Molecular Dynamics Calculations.
1507 * All of the following is SpeedMeasure()'d in SimTime.
1508 * Initialization by calling:
1509 * -# CorrectVelocity()\n
1510 * Shifts center of gravity of Ions momenta, so that the cell itself remains at rest.
1511 * -# CalculateEnergyIonsU(), SpeedMeasure()'d in TimeTypes#InitSimTime\n
1512 * Calculates kinetic energy of "movable" Ions.
1513 * -# CalculateForce()\n
1514 * Does the minimisation, calculates densities, then energies and finally the forces.
1515 * -# OutputVisSrcFiles()\n
1516 * If desired, so-far made calculations are stored to file for later restarting.
1517 * -# OutputIonForce()\n
1518 * Write ion forces to file.
1519 * -# EnergyOutput()\n
1520 * Write calculated energies to screen or file.
1521 *
1522 * The simulation phase begins:
1523 * -# UpdateIonsR()\n
1524 * Move Ions according to the calculated force.
1525 * -# UpdateWaveAfterIonMove()\n
1526 * Update wave functions by averaging LocalPsi coefficients after the Ions have been shifted.
1527 * -# UpdateToNewWaves()\n
1528 * Update after wave functions have changed.
1529 * -# CalculateForce()\n
1530 * Does the minimisation, calculates densities, then energies and finally the forces.
1531 * -# UpdateIonsU()\n
1532 * Change ion's velocities according to the calculated acting force.
1533 * -# CorrectVelocity()\n
1534 * Shifts center of gravity of Ions momenta, so that the cell itself remains at rest.
1535 * -# CalculateEnergyIonsU()\n
1536 * Calculates kinetic energy of "movable" Ions.
1537 * -# ScaleTemp()\n
1538 * The temperature is scaled, so the systems energy remains constant (they must not gain momenta out of nothing)
1539 * -# OutputVisSrcFiles()\n
1540 * If desired, so-far made calculations are stored to file for later restarting.
1541 * -# OutputVis()\n
1542 * Visulization data for OpenDX is written at certain steps if desired.
1543 * -# OutputIonForce()\n
1544 * Write ion forces to file.
1545 * -# EnergyOutput()\n
1546 * Write calculated energies to screen or file.
1547 *
1548 * After the ground state is found:
1549 * -# CalculateUnOccupied()\n
1550 * Energies of unoccupied orbitals - that have been left out completely so far -
1551 * are calculated.
1552 * -# TestGramSch()\n
1553 * Test if orbitals are still orthogonal.
1554 * -# CalculateHamiltonian()\n
1555 * Construct Hamiltonian and calculate Eigenvalues.
1556 * -# ComputeMLWF()\n
1557 * Localize orbital wave functions.
1558 *
1559 * \param *P Problem at hand
1560 */
1561void CalculateMD(struct Problem *P)
1562{
1563 struct RunStruct *R = &P->R;
1564 struct Ions *I = &P->Ion;
1565 int OuterStop = 0;
1566 SpeedMeasure(P, SimTime, StartTimeDo);
1567 SpeedMeasure(P, InitSimTime, StartTimeDo);
1568 R->OuterStep = 0;
1569 CorrectVelocity(P);
1570 CalculateEnergyIonsU(P);
1571 OuterStop = CalculateForce(P);
1572 R->OuterStep++;
1573 P->Speed.InitSteps++;
1574 SpeedMeasure(P, InitSimTime, StopTimeDo);
1575 OutputIonForce(P);
1576 EnergyOutput(P, 1);
1577 if (R->MaxOuterStep > 0) {
1578 debug(P,"Commencing Fletcher-Reeves minimisation on ionic structure ...");
1579 UpdateIon_PRCG(P);
1580 }
1581 if (I->StructOpt && !OuterStop) {
1582 I->StructOpt = 0;
1583 OuterStop = CalculateForce(P);
1584 }
1585 /* while (!OuterStop && R->OuterStep <= R->MaxOuterStep) {
1586 R->OuterStep++;
1587 if (P->Call.out[NormalOut]) fprintf(stderr,"(%i) Commencing MD steps %i ... \n",P->Par.me, R->OuterStep);
1588 P->R.t += P->R.delta_t; // increase current time by delta_t
1589 R->NewRStep++;
1590 if (P->Ion.StructOpt == 1) {
1591 UpdateIons(P);
1592 OutputIonCoordinates(P);
1593 } else {
1594 UpdateIonsR(P);
1595 }
1596 UpdateWaveAfterIonMove(P);
1597 for (i=MAXOLD-1; i > 0; i--) // store away old energies
1598 E->TotalEnergyOuter[i] = E->TotalEnergyOuter[i-1];
1599 UpdateToNewWaves(P);
1600 E->TotalEnergyOuter[0] = E->TotalEnergy[0];
1601 OuterStop = CalculateForce(P);
1602 UpdateIonsU(P);
1603 CorrectVelocity(P);
1604 CalculateEnergyIonsU(P);
1605 if ((P->R.ScaleTempStep > 0) && ((R->OuterStep-1) % P->R.ScaleTempStep == 0))
1606 ScaleTemp(P);
1607 if ((R->OuterStep-1) % P->R.OutSrcStep == 0)
1608 OutputVisSrcFiles(P, Occupied);
1609 if ((R->OuterStep-1) % P->R.OutVisStep == 0) {
1610 // recalculate density for the specific wave function ...
1611 //CalculateOneDensityR(Lat, LevS, Dens0, PsiDat, Dens0->DensityArray[ActualDensity], R->FactorDensityR, 0);
1612 // ... and output (wherein ActualDensity is used instead of TotalDensity)
1613 //OutputVis(P);
1614 //OutputIonForce(P);
1615 //EnergyOutput(P, 1);
1616 }
1617 }*/
1618 SpeedMeasure(P, SimTime, StopTimeDo);
1619 // hack to output each orbital before and after spread-minimisation
1620/* if (P->Files.MeOutVis) P->Files.OutputPosType = (enum ModeType *) Realloc(P->Files.OutputPosType,sizeof(enum ModeType)*(P->Files.OutVisStep+P->Lat.Psi.MaxPsiOfType*2),"OutputVis");
1621 OutputVisAllOrbital(P, 0, 2, Occupied);
1622 CalculateHamiltonian(P);
1623 if (P->Files.MeOutVis) P->Files.OutVisStep -= (P->Lat.Psi.MaxPsiOfType)*2;
1624 OutputVisAllOrbital(P, 1, 2, Occupied);*/
1625 CloseOutputFiles(P);
1626}
Note: See TracBrowser for help on using the repository browser.