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