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