source: pcp/src/run.c@ 79290f

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

config.h is included in each and every file. After trying to compile on JUMP (with xlc).

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