| [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;
 | 
|---|
| [6df7db] | 723 |   struct Ions *I = &P->Ion;
 | 
|---|
| [a0bcf1] | 724 |   //struct FileData *F = &P->Files;
 | 
|---|
 | 725 | //  int i;
 | 
|---|
 | 726 | //  double norm;
 | 
|---|
 | 727 |   //double dEdt0,ddEddt0,HartreeddEddt0,XCddEddt0, d[4], D[4],ConDirHConDir;
 | 
|---|
 | 728 |   struct LatticeLevel *LevS = R->LevS;
 | 
|---|
 | 729 |   int ElementSize = (sizeof(fftw_complex) / sizeof(double));
 | 
|---|
 | 730 |   int iter = 0, status, max_iter=10;
 | 
|---|
 | 731 |   const gsl_min_fminimizer_type *T;
 | 
|---|
 | 732 |   gsl_min_fminimizer *s;
 | 
|---|
 | 733 |   double m, a, b;
 | 
|---|
 | 734 |   double f_m = 0., f_a, f_b;
 | 
|---|
 | 735 |   double dcos, dsin;
 | 
|---|
 | 736 |   int g;
 | 
|---|
 | 737 |   fftw_complex *ConDir = P->Grad.GradientArray[ConDirGradient];
 | 
|---|
 | 738 |   fftw_complex *source = NULL, *oldsource = NULL; 
 | 
|---|
 | 739 |   gsl_function F;
 | 
|---|
 | 740 |   F.function = &fn1;
 | 
|---|
 | 741 |   F.params = (void *) P;
 | 
|---|
 | 742 |   T = gsl_min_fminimizer_brent;
 | 
|---|
 | 743 |   s = gsl_min_fminimizer_alloc (T);
 | 
|---|
 | 744 |   int DoBrent, StartLocalPsiNo;
 | 
|---|
 | 745 |  
 | 
|---|
 | 746 |   ResetBrent(P,Psi);
 | 
|---|
 | 747 |   *Stop = 0;
 | 
|---|
| [6df7db] | 748 |   if ((P->Call.ReadSrcFiles != DoNotParse) && (!I->StructOpt)) {
 | 
|---|
| [a0bcf1] | 749 |     if (!ReadSrcPsiDensity(P,Occupied,1, R->LevSNo)) {  // if file for level exists and desired, read from file
 | 
|---|
| [d6f7f3] | 750 |       P->Call.ReadSrcFiles = DoNotParse; // -r was bogus, remove it, have to start anew
 | 
|---|
| [4931e0] | 751 |       if(P->Call.out[MinOut]) fprintf(stderr,"(%i) Re-initializing, files are missing/corrupted...\n", P->Par.me);
 | 
|---|
| [a0bcf1] | 752 |       InitPsisValue(P, Psi->TypeStartIndex[Occupied], Psi->TypeStartIndex[Occupied+1]);  // initialize perturbed array for this run
 | 
|---|
 | 753 |       ResetGramSchTagType(P, Psi, Occupied, NotOrthogonal); // loaded values are orthonormal
 | 
|---|
 | 754 |     } else {
 | 
|---|
 | 755 |       SpeedMeasure(P, InitSimTime, StartTimeDo);  
 | 
|---|
| [487182] | 756 |       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] | 757 |       ReadSrcPsiDensity(P, Occupied, 0, R->LevSNo);
 | 
|---|
| [487182] | 758 |       //ResetGramSchTagType(P, Psi, Occupied, IsOrthonormal); // loaded values are orthonormal
 | 
|---|
 | 759 |       // 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] | 760 |     }
 | 
|---|
| [487182] | 761 |     SpeedMeasure(P, InitGramSchTime, StartTimeDo);  
 | 
|---|
 | 762 |     GramSch(P, R->LevS, Psi, Orthonormalize);
 | 
|---|
 | 763 |     SpeedMeasure(P, InitGramSchTime, StopTimeDo);  
 | 
|---|
| [a0bcf1] | 764 |     SpeedMeasure(P, InitDensityTime, StartTimeDo);
 | 
|---|
 | 765 |     InitDensityCalculation(P);
 | 
|---|
 | 766 |     SpeedMeasure(P, InitDensityTime, StopTimeDo);
 | 
|---|
 | 767 |     InitPsiEnergyCalculation(P, Occupied);  // go through all orbitals calculating kinetic and non-local
 | 
|---|
 | 768 |     StartLocalPsiNo = R->ActualLocalPsiNo;
 | 
|---|
 | 769 |     do { // otherwise OnePsiElementAddData#Lambda is calculated only for current Psi not for all
 | 
|---|
 | 770 |       CalculateDensityEnergy(P, 0);
 | 
|---|
 | 771 |       UpdateActualPsiNo(P, Occupied);
 | 
|---|
 | 772 |     } while (R->ActualLocalPsiNo != StartLocalPsiNo);
 | 
|---|
 | 773 |     CalculateIonsEnergy(P);
 | 
|---|
 | 774 |     EnergyAllReduce(P);
 | 
|---|
 | 775 |     SpeedMeasure(P, InitSimTime, StopTimeDo);  
 | 
|---|
 | 776 |     R->LevS->Step++;
 | 
|---|
 | 777 |     EnergyOutput(P,0);
 | 
|---|
| [487182] | 778 |   } 
 | 
|---|
| [6df7db] | 779 |   if ((I->StructOpt) || ((P->Call.ReadSrcFiles != DoReadAllSrcDensities) && (P->Call.ReadSrcFiles != DoReadOccupiedSrcDensities))) {  // otherwise minimise oneself
 | 
|---|
| [4931e0] | 780 |     if(P->Call.out[LeaderOut]) fprintf(stderr,"(%i)Beginning minimisation of type %s ...\n", P->Par.me, R->MinimisationName[Occupied]);
 | 
|---|
| [a0bcf1] | 781 |     while (*Stop != 1) {  // loop testing condition over all Psis
 | 
|---|
 | 782 |       // in the following loop, we have two cases:
 | 
|---|
 | 783 |       // 1) still far away and just guessing: Use the normal CalculateNewWave() to improve Psi
 | 
|---|
 | 784 |       // 2) closer (DoBrent=-1): use brent line search instead
 | 
|---|
 | 785 |       // and due to these two cases, we also have two ifs inside each in order to catch stepping from one case
 | 
|---|
 | 786 |       // to the other - due to decreasing DoBrent and/or stepping to the next Psi (which may not yet be DoBrent==1)
 | 
|---|
 | 787 | 
 | 
|---|
 | 788 |       // case 1)
 | 
|---|
 | 789 |       if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent != 1) {
 | 
|---|
 | 790 |         //SetArrayToDouble0((double *)LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo],LevS->MaxG*2);
 | 
|---|
| [487182] | 791 |         if (R->DoBrent == 1) {
 | 
|---|
 | 792 |           memcpy(LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
 | 
|---|
 | 793 |           //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);
 | 
|---|
 | 794 |           f_m = P->Lat.E->TotalEnergy[0]; // grab first value
 | 
|---|
 | 795 |           m = 0.;
 | 
|---|
 | 796 |         }
 | 
|---|
| [a0bcf1] | 797 |         CalculateNewWave(P,NULL);
 | 
|---|
 | 798 |         if ((R->DoBrent == 1) && (fabs(Lat->E->delta[0]) < M_PI/4.)) 
 | 
|---|
 | 799 |           Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent--;
 | 
|---|
 | 800 |         if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent != 1) {
 | 
|---|
 | 801 |           UpdateActualPsiNo(P, Occupied);
 | 
|---|
 | 802 |           UpdateEnergyArray(P);
 | 
|---|
 | 803 |           CalculateEnergy(P); // just to get a sensible delta
 | 
|---|
 | 804 |           if ((R->ActualLocalPsiNo != R->OldActualLocalPsiNo) && (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent == 1)) {
 | 
|---|
 | 805 |              R->OldActualLocalPsiNo = R->ActualLocalPsiNo;
 | 
|---|
 | 806 |             // if we stepped on to a new Psi, which is already down at DoBrent=1 unlike the last one,
 | 
|---|
 | 807 |             // then an up-to-date gradient is missing for the following Brent line search
 | 
|---|
| [487182] | 808 |             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] | 809 |             memcpy(LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
 | 
|---|
 | 810 |             //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);
 | 
|---|
 | 811 |             f_m = P->Lat.E->TotalEnergy[0]; // grab first value
 | 
|---|
 | 812 |             m = 0.;
 | 
|---|
 | 813 |             DoBrent = Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent;
 | 
|---|
 | 814 |             Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent = 2;
 | 
|---|
 | 815 |             CalculateNewWave(P,NULL);
 | 
|---|
 | 816 |             Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent = DoBrent;
 | 
|---|
 | 817 |           }
 | 
|---|
 | 818 |           //fprintf(stderr,"(%i) fnl(%lg) = %lg\n", P->Par.me, m, f_m);
 | 
|---|
 | 819 |         } 
 | 
|---|
 | 820 |       }
 | 
|---|
 | 821 |       
 | 
|---|
 | 822 |       // case 2)
 | 
|---|
 | 823 |       if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent == 1) {
 | 
|---|
 | 824 |         R->PsiStep=R->MaxPsiStep;  // no more fresh gradients from this point for current ActualLocalPsiNo
 | 
|---|
 | 825 |         a = b = 0.5*fabs(Lat->E->delta[0]);
 | 
|---|
 | 826 |         // we have a meaningful first minimum guess from above CalculateNewWave() resp. from end of this if of last step: Lat->E->delta[0]
 | 
|---|
 | 827 |         source = LevS->LPsi->LocalPsi[R->ActualLocalPsiNo];
 | 
|---|
 | 828 |         oldsource = LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo];
 | 
|---|
 | 829 |         //SetArrayToDouble0((double *)source,LevS->MaxG*2);
 | 
|---|
 | 830 |         do { 
 | 
|---|
 | 831 |           a -= fabs(Lat->E->delta[0]) == 0 ? 0.1 : fabs(Lat->E->delta[0]);
 | 
|---|
 | 832 |           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)
 | 
|---|
 | 833 |           dcos = cos(a);
 | 
|---|
 | 834 |           dsin = sin(a);
 | 
|---|
 | 835 |           for (g = 0; g < LevS->MaxG; g++) {  // Here all coefficients are updated for the new found wave function
 | 
|---|
 | 836 |             //if (isnan(ConDir[g].re)) { fprintf(stderr,"WARNGING: CalculateLineSearch(): ConDir_%i(%i) = NaN!\n", R->ActualLocalPsiNo, g); Error(SomeError, "NaN-Fehler!"); }
 | 
|---|
 | 837 |             c_re(source[g]) = c_re(oldsource[g])*dcos + c_re(ConDir[g])*dsin;
 | 
|---|
 | 838 |             c_im(source[g]) = c_im(oldsource[g])*dcos + c_im(ConDir[g])*dsin;
 | 
|---|
 | 839 |           }
 | 
|---|
 | 840 |           CalculateEnergy(P);
 | 
|---|
 | 841 |           f_a = P->Lat.E->TotalEnergy[0]; // grab second value at left border
 | 
|---|
 | 842 |           //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]);
 | 
|---|
 | 843 |         } while (f_a < f_m);
 | 
|---|
 | 844 | 
 | 
|---|
 | 845 |         //SetArrayToDouble0((double *)source,LevS->MaxG*2);
 | 
|---|
 | 846 |         do {
 | 
|---|
 | 847 |           b += fabs(Lat->E->delta[0]) == 0 ? 0.1 : fabs(Lat->E->delta[0]);
 | 
|---|
 | 848 |           if (b > M_PI/2.) b = M_PI/2.;
 | 
|---|
 | 849 |           dcos = cos(b);
 | 
|---|
 | 850 |           dsin = sin(b);
 | 
|---|
 | 851 |           for (g = 0; g < LevS->MaxG; g++) {  // Here all coefficients are updated for the new found wave function
 | 
|---|
 | 852 |             //if (isnan(ConDir[g].re)) { fprintf(stderr,"WARNGING: CalculateLineSearch(): ConDir_%i(%i) = NaN!\n", R->ActualLocalPsiNo, g); Error(SomeError, "NaN-Fehler!"); }
 | 
|---|
 | 853 |             c_re(source[g]) = c_re(oldsource[g])*dcos + c_re(ConDir[g])*dsin;
 | 
|---|
 | 854 |             c_im(source[g]) = c_im(oldsource[g])*dcos + c_im(ConDir[g])*dsin;
 | 
|---|
 | 855 |           }
 | 
|---|
 | 856 |           CalculateEnergy(P);
 | 
|---|
 | 857 |           f_b = P->Lat.E->TotalEnergy[0]; // grab second value at left border
 | 
|---|
 | 858 |           //fprintf(stderr,"(%i) fnl(%lg) = %lg\n", P->Par.me, b, f_b);
 | 
|---|
 | 859 |         } while (f_b < f_m);
 | 
|---|
 | 860 |         
 | 
|---|
 | 861 |         memcpy(source, oldsource, ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
 | 
|---|
 | 862 |         //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);
 | 
|---|
 | 863 |         CalculateEnergy(P);
 | 
|---|
 | 864 | 
 | 
|---|
| [4931e0] | 865 |         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] | 866 |         iter=0; 
 | 
|---|
 | 867 |         gsl_min_fminimizer_set_with_values (s, &F, m, f_m, a, f_a, b, f_b);
 | 
|---|
| [4931e0] | 868 |         if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) using %s method\n",P->Par.me, gsl_min_fminimizer_name (s));
 | 
|---|
 | 869 |         if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) %5s [%9s, %9s] %9s %9s\n",P->Par.me, "iter", "lower", "upper", "min", "err(est)");
 | 
|---|
 | 870 |         if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) %5d [%.7f, %.7f] %.7f %.7f\n",P->Par.me, iter, a, b, m, b - a);
 | 
|---|
| [a0bcf1] | 871 |         do {
 | 
|---|
 | 872 |           iter++;
 | 
|---|
 | 873 |           status = gsl_min_fminimizer_iterate (s);
 | 
|---|
 | 874 |        
 | 
|---|
 | 875 |           m = gsl_min_fminimizer_x_minimum (s);
 | 
|---|
 | 876 |           a = gsl_min_fminimizer_x_lower (s);
 | 
|---|
 | 877 |           b = gsl_min_fminimizer_x_upper (s);
 | 
|---|
 | 878 |        
 | 
|---|
 | 879 |           status = gsl_min_test_interval (a, b, 0.001, 0.0);
 | 
|---|
 | 880 |        
 | 
|---|
 | 881 |           if (status == GSL_SUCCESS)
 | 
|---|
| [4931e0] | 882 |             if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) Converged:\n",P->Par.me);
 | 
|---|
| [a0bcf1] | 883 |        
 | 
|---|
| [4931e0] | 884 |           if(P->Call.out[ValueOut]) fprintf (stderr,"(%i) %5d [%.7f, %.7f] %.7f %.7f\n",P->Par.me,
 | 
|---|
| [a0bcf1] | 885 |             iter, a, b, m, b - a);
 | 
|---|
 | 886 |         } while (status == GSL_CONTINUE && iter < max_iter);
 | 
|---|
 | 887 |         CalculateNewWave(P,&m);
 | 
|---|
 | 888 |         TestGramSch(P,LevS,Psi,Occupied);
 | 
|---|
 | 889 |         UpdateActualPsiNo(P, Occupied); // step on due setting to MaxPsiStep further above
 | 
|---|
 | 890 |         UpdateEnergyArray(P);
 | 
|---|
 | 891 |         CalculateEnergy(P);
 | 
|---|
 | 892 |         //fprintf(stderr,"(%i) Final value for Psi %i: %lg\n", P->Par.me, R->ActualLocalPsiNo, P->Lat.E->TotalEnergy[0]);
 | 
|---|
 | 893 |         R->MinStopStep =  R->ActualMaxMinStopStep; // check stop condition every time
 | 
|---|
 | 894 |         if (*SuperStop != 1)
 | 
|---|
 | 895 |           *SuperStop = CheckCPULIM(P);
 | 
|---|
 | 896 |         *Stop = CalculateMinimumStop(P, *SuperStop);
 | 
|---|
 | 897 |         R->OldActualLocalPsiNo = R->ActualLocalPsiNo;
 | 
|---|
 | 898 |         if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent == 1) { // new wave function means new gradient!
 | 
|---|
 | 899 |           DoBrent = Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent;
 | 
|---|
 | 900 |           Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent = 2;
 | 
|---|
 | 901 |           //SetArrayToDouble0((double *)LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo],LevS->MaxG*2);
 | 
|---|
 | 902 |           memcpy(LevS->LPsi->OldLocalPsi[R->ActualLocalPsiNo], LevS->LPsi->LocalPsi[R->ActualLocalPsiNo], ElementSize*LevS->MaxG*sizeof(double)); // restore old Psi from OldPsi
 | 
|---|
 | 903 |           //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);
 | 
|---|
 | 904 |           f_m = P->Lat.E->TotalEnergy[0]; // grab first value
 | 
|---|
 | 905 |           m = 0.;
 | 
|---|
 | 906 |           CalculateNewWave(P,NULL);
 | 
|---|
 | 907 |           Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent = DoBrent;
 | 
|---|
 | 908 |         }
 | 
|---|
 | 909 |       }
 | 
|---|
 | 910 | 
 | 
|---|
 | 911 |       if (Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].DoBrent != 1) { // otherwise the following checks eliminiate stop=1 from above
 | 
|---|
 | 912 |         if (*SuperStop != 1)
 | 
|---|
 | 913 |           *SuperStop = CheckCPULIM(P);
 | 
|---|
 | 914 |         *Stop = CalculateMinimumStop(P, *SuperStop);
 | 
|---|
 | 915 |       }
 | 
|---|
 | 916 |       /*EnergyOutput(P, Stop);*/
 | 
|---|
 | 917 |       P->Speed.Steps++;
 | 
|---|
 | 918 |       R->LevS->Step++;
 | 
|---|
 | 919 |       /*ControlNativeDensity(P);*/
 | 
|---|
 | 920 |       //fprintf(stderr,"(%i) Stop %i\n",P->Par.me, Stop);
 | 
|---|
 | 921 |     }
 | 
|---|
| [5a538b] | 922 |     if (*SuperStop == 1) OutputVisSrcFiles(P, Occupied); // is now done after localization (ComputeMLWF())
 | 
|---|
| [a0bcf1] | 923 |   }
 | 
|---|
 | 924 |   TestGramSch(P,R->LevS,Psi, Occupied);
 | 
|---|
 | 925 | }
 | 
|---|
 | 926 | 
 | 
|---|
 | 927 | /** Minimisation of the PsiTagType#UnOccupied orbitals in the field of the occupied ones.
 | 
|---|
 | 928 |  * Assuming RunStruct#ActualLocalPsiNo is currenlty still an occupied wave function, we stop onward to the first
 | 
|---|
 | 929 |  * unoccupied and reset RunStruct#OldActualLocalPsiNo. Then it is checked whether CallOptions#ReadSrcFiles is set
 | 
|---|
 | 930 |  * and thus coefficients for the level have to be read from file and afterwards initialized.
 | 
|---|
 | 931 |  * 
 | 
|---|
 | 932 |  * Then follows the main loop, until a stop condition is met:
 | 
|---|
 | 933 |  * -# CalculateNewWave()\n
 | 
|---|
 | 934 |  *    Over a conjugate gradient method the next (minimal) wave function is sought for.
 | 
|---|
 | 935 |  * -# UpdateActualPsiNo()\n
 | 
|---|
 | 936 |  *    Switch local Psi to next one.
 | 
|---|
 | 937 |  * -# UpdateEnergyArray()\n
 | 
|---|
 | 938 |  *    Shift archived energy values to make space for next one.
 | 
|---|
 | 939 |  * -# UpdateDensityCalculation(), SpeedMeasure()'d in DensityTime\n
 | 
|---|
 | 940 |  *    Calculate TotalLocalDensity of LocalPsis and gather results as TotalDensity.
 | 
|---|
 | 941 |  * -# UpdatePsiEnergyCalculation()\n
 | 
|---|
 | 942 |  *    Calculate kinetic and non-local energy contributons from the Psis.
 | 
|---|
 | 943 |  * -# CalculateGapEnergy()\n
 | 
|---|
 | 944 |  *    Calculate Gap energies (Hartreepotential, Pseudo) and the gradient.
 | 
|---|
 | 945 |  * -# EnergyAllReduce()\n
 | 
|---|
 | 946 |  *    Gather PsiEnergy results from all processes and sum up together with all other contributions to TotalEnergy.
 | 
|---|
 | 947 |  * -# CheckCPULIM()\n
 | 
|---|
 | 948 |  *    Check if external signal has been received (e.g. end of time slit on cluster), break operation at next possible moment.
 | 
|---|
 | 949 |  * -# CalculateMinimumStop()\n
 | 
|---|
 | 950 |  *    Evaluates stop condition if desired precision or steps or ... have been reached. Otherwise go to
 | 
|---|
 | 951 |  *    CalculateNewWave().
 | 
|---|
 | 952 |  * 
 | 
|---|
 | 953 |  * Afterwards, the coefficients are written to file by OutputVisSrcFiles() if desired. Orthonormality is tested, we step
 | 
|---|
 | 954 |  * back to the occupied wave functions and the densities are re-initialized. 
 | 
|---|
 | 955 |  * \param *P Problem at hand
 | 
|---|
 | 956 |  * \param *Stop flag to determine if epsilon stop conditions have met
 | 
|---|
 | 957 |  * \param *SuperStop flag to determinte whether external signal's required end of calculations
 | 
|---|
 | 958 |  */
 | 
|---|
 | 959 | static void MinimiseUnoccupied (struct Problem *P, int *Stop, int *SuperStop) {
 | 
|---|
 | 960 |   struct RunStruct *R = &P->R;
 | 
|---|
 | 961 |   struct Lattice *Lat = &P->Lat;
 | 
|---|
 | 962 |   struct Psis *Psi = &Lat->Psi;
 | 
|---|
 | 963 |   int StartLocalPsiNo;
 | 
|---|
 | 964 | 
 | 
|---|
 | 965 |   *Stop = 0;
 | 
|---|
 | 966 |   R->PsiStep = R->MaxPsiStep; // in case it's zero from CalculateForce()
 | 
|---|
 | 967 |   UpdateActualPsiNo(P, UnOccupied); // step on to next unoccupied one
 | 
|---|
 | 968 |   R->OldActualLocalPsiNo = R->ActualLocalPsiNo; // reset, otherwise OldActualLocalPsiNo still points to occupied wave function
 | 
|---|
 | 969 |   UpdateGramSchOldActualPsiNo(P,Psi);
 | 
|---|
| [d6f7f3] | 970 |   if ((P->Call.ReadSrcFiles == DoReadAllSrcDensities) && ReadSrcPsiDensity(P,UnOccupied,1, R->LevSNo)) {
 | 
|---|
| [a0bcf1] | 971 |     SpeedMeasure(P, InitSimTime, StartTimeDo);  
 | 
|---|
| [4931e0] | 972 |     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] | 973 |     ReadSrcPsiDensity(P, UnOccupied, 0, R->LevSNo);
 | 
|---|
| [d6f7f3] | 974 |     if (P->Call.ReadSrcFiles < DoReadAndMinimise) {
 | 
|---|
| [a0bcf1] | 975 |       ResetGramSchTagType(P, Psi, UnOccupied, IsOrthonormal); // loaded values are orthonormal
 | 
|---|
 | 976 |       SpeedMeasure(P, DensityTime, StartTimeDo);
 | 
|---|
 | 977 |       InitDensityCalculation(P);
 | 
|---|
 | 978 |       SpeedMeasure(P, DensityTime, StopTimeDo);
 | 
|---|
 | 979 |       InitPsiEnergyCalculation(P,UnOccupied);  // go through all orbitals calculating kinetic and non-local
 | 
|---|
 | 980 |       //CalculateDensityEnergy(P, 0);
 | 
|---|
 | 981 |       StartLocalPsiNo = R->ActualLocalPsiNo;
 | 
|---|
 | 982 |       do { // otherwise OnePsiElementAddData#Lambda is calculated only for current Psi not for all
 | 
|---|
 | 983 |         CalculateGapEnergy(P);
 | 
|---|
 | 984 |         UpdateActualPsiNo(P, Occupied);
 | 
|---|
 | 985 |       } while (R->ActualLocalPsiNo != StartLocalPsiNo);
 | 
|---|
 | 986 |       EnergyAllReduce(P);
 | 
|---|
 | 987 |     }
 | 
|---|
 | 988 |     SpeedMeasure(P, InitSimTime, StopTimeDo);  
 | 
|---|
 | 989 |   }
 | 
|---|
| [d6f7f3] | 990 |   if (P->Call.ReadSrcFiles != DoReadAllSrcDensities) {
 | 
|---|
| [a0bcf1] | 991 |     SpeedMeasure(P, InitSimTime, StartTimeDo);  
 | 
|---|
 | 992 |     ResetGramSchTagType(P, Psi, UnOccupied, NotOrthogonal);
 | 
|---|
 | 993 |     SpeedMeasure(P, GramSchTime, StartTimeDo);  
 | 
|---|
 | 994 |     GramSch(P, R->LevS, Psi, Orthonormalize);
 | 
|---|
 | 995 |     SpeedMeasure(P, GramSchTime, StopTimeDo);  
 | 
|---|
 | 996 |     SpeedMeasure(P, InitDensityTime, StartTimeDo);
 | 
|---|
 | 997 |     InitDensityCalculation(P);
 | 
|---|
 | 998 |     SpeedMeasure(P, InitDensityTime, StopTimeDo);
 | 
|---|
 | 999 |     InitPsiEnergyCalculation(P,UnOccupied);  // go through all orbitals calculating kinetic and non-local
 | 
|---|
 | 1000 |     //CalculateDensityEnergy(P, 0);
 | 
|---|
 | 1001 |     CalculateGapEnergy(P);
 | 
|---|
 | 1002 |     EnergyAllReduce(P);
 | 
|---|
 | 1003 |     SpeedMeasure(P, InitSimTime, StopTimeDo);  
 | 
|---|
 | 1004 |     R->LevS->Step++;
 | 
|---|
 | 1005 |     EnergyOutput(P,0);
 | 
|---|
| [4931e0] | 1006 |     if(P->Call.out[LeaderOut]) fprintf(stderr,"(%i)Beginning minimisation of type %s ...\n", P->Par.me, R->MinimisationName[R->CurrentMin]);
 | 
|---|
| [a0bcf1] | 1007 |     while (*Stop != 1) {
 | 
|---|
 | 1008 |       CalculateNewWave(P,NULL);
 | 
|---|
 | 1009 |       UpdateActualPsiNo(P, UnOccupied);
 | 
|---|
 | 1010 |       /* New */
 | 
|---|
 | 1011 |       UpdateEnergyArray(P);
 | 
|---|
 | 1012 |       SpeedMeasure(P, DensityTime, StartTimeDo);
 | 
|---|
 | 1013 |       UpdateDensityCalculation(P);
 | 
|---|
 | 1014 |       SpeedMeasure(P, DensityTime, StopTimeDo);
 | 
|---|
 | 1015 |       UpdatePsiEnergyCalculation(P);
 | 
|---|
 | 1016 |       //CalculateDensityEnergy(P, 0);
 | 
|---|
 | 1017 |       CalculateGapEnergy(P);  //calculates XC, HGDensity, afterwards gradient, where V_xc is added upon HGDensity
 | 
|---|
 | 1018 |       EnergyAllReduce(P);
 | 
|---|
 | 1019 |       if (*SuperStop != 1)
 | 
|---|
 | 1020 |         *SuperStop = CheckCPULIM(P);
 | 
|---|
 | 1021 |       *Stop = CalculateMinimumStop(P, *SuperStop);
 | 
|---|
 | 1022 |       /*EnergyOutput(P, Stop);*/
 | 
|---|
 | 1023 |       P->Speed.Steps++;
 | 
|---|
 | 1024 |       R->LevS->Step++;
 | 
|---|
 | 1025 |       /*ControlNativeDensity(P);*/
 | 
|---|
 | 1026 |     }
 | 
|---|
 | 1027 |     OutputVisSrcFiles(P, UnOccupied);
 | 
|---|
 | 1028 | //    if (!TestReadnWriteSrcDensity(P,UnOccupied))
 | 
|---|
 | 1029 | //      Error(SomeError,"TestReadnWriteSrcDensity failed!");
 | 
|---|
 | 1030 |   }
 | 
|---|
 | 1031 |   TestGramSch(P,R->LevS,Psi, UnOccupied);
 | 
|---|
 | 1032 |   ResetGramSchTagType(P, Psi, UnOccupied, NotUsedToOrtho);
 | 
|---|
 | 1033 |   // re-calculate Occupied values (in preparation for perturbed ones)
 | 
|---|
 | 1034 |   UpdateActualPsiNo(P, Occupied); // step on to next occupied one
 | 
|---|
 | 1035 |   SpeedMeasure(P, DensityTime, StartTimeDo);
 | 
|---|
 | 1036 |   InitDensityCalculation(P);  // re-init densities to occupied standard
 | 
|---|
 | 1037 |   SpeedMeasure(P, DensityTime, StopTimeDo);
 | 
|---|
 | 1038 | //  InitPsiEnergyCalculation(P,Occupied);
 | 
|---|
 | 1039 | //  CalculateDensityEnergy(P, 0);
 | 
|---|
 | 1040 | //  EnergyAllReduce(P);
 | 
|---|
 | 1041 | }
 | 
|---|
 | 1042 | 
 | 
|---|
 | 1043 | 
 | 
|---|
 | 1044 | /** Calculate the forces.
 | 
|---|
 | 1045 |  * From RunStruct::LevSNo downto RunStruct::InitLevSNo the following routines are called in a loop:
 | 
|---|
 | 1046 |  * -# In case of RunStruct#DoSeparated another loop begins for the unoccupied states with some reinitalization
 | 
|---|
 | 1047 |  *    before and afterwards. The loop however is much the same as the one above.
 | 
|---|
 | 1048 |  * -# ChangeToLevUp()\n
 | 
|---|
 | 1049 |  *    Repeat the loop or when the above stop is reached, the level is changed and the loop repeated.
 | 
|---|
 | 1050 |  * 
 | 
|---|
 | 1051 |  * Afterwards comes the actual force and energy calculation by calling:
 | 
|---|
 | 1052 |  * -# ControlNativeDensity()\n
 | 
|---|
 | 1053 |  *    Checks if the density still reproduces particle number.
 | 
|---|
 | 1054 |  * -# CalculateIonLocalForce(), SpeedMeasure()'d in LocFTime\n
 | 
|---|
 | 1055 |  *    Calculale local part of force acting on Ions.
 | 
|---|
 | 1056 |  * -# CalculateIonNonLocalForce(), SpeedMeasure()'d in NonLocFTime\n
 | 
|---|
 | 1057 |  *    Calculale local part of force acting on Ions.
 | 
|---|
 | 1058 |  * -# CalculateEwald()\n
 | 
|---|
 | 1059 |  *    Calculate Ewald force acting on Ions.
 | 
|---|
 | 1060 |  * -# CalculateIonForce()\n
 | 
|---|
 | 1061 |  *    Sum up those three contributions.
 | 
|---|
 | 1062 |  * -# CorrectForces()\n
 | 
|---|
 | 1063 |  *    Shifts center of gravity of all forces for each Ion, so that the cell itself remains at rest.
 | 
|---|
 | 1064 |  * -# GetOuterStop()
 | 
|---|
 | 1065 |  *    Calculates a mean force per Ion.
 | 
|---|
 | 1066 |  * \param *P Problem at hand
 | 
|---|
 | 1067 |  * \return 1 - cpulim received, break operation, 0 - continue as normal
 | 
|---|
 | 1068 |  */
 | 
|---|
 | 1069 | int CalculateForce(struct Problem *P)
 | 
|---|
 | 1070 | {
 | 
|---|
 | 1071 |   struct RunStruct *R = &P->R;
 | 
|---|
 | 1072 |   struct Lattice *Lat = &P->Lat;
 | 
|---|
 | 1073 |   struct Psis *Psi = &Lat->Psi;
 | 
|---|
 | 1074 |   struct LatticeLevel *LevS = R->LevS;
 | 
|---|
 | 1075 |   struct FileData *F = &P->Files;
 | 
|---|
 | 1076 |   struct Ions *I = &P->Ion;
 | 
|---|
 | 1077 |   int Stop=0, SuperStop = 0, OuterStop = 0;
 | 
|---|
 | 1078 |   //int i, j;
 | 
|---|
| [5712cb] | 1079 |   SpeedMeasure(P, SimTime, StartTimeDo);
 | 
|---|
 | 1080 |   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
 | 
|---|
 | 1081 |     while ((R->LevSNo > R->InitLevSNo) || (!Stop && R->LevSNo == R->InitLevSNo)) {
 | 
|---|
 | 1082 |       // occupied
 | 
|---|
| [333e84] | 1083 |       R->PsiStep = R->MaxPsiStep; // reset in-Psi-minimisation-counter, so that we really advance to the next wave function
 | 
|---|
| [5712cb] | 1084 |       R->OldActualLocalPsiNo = R->ActualLocalPsiNo; // reset OldActualLocalPsiNo, as it might still point to a perturbed wave function from last level
 | 
|---|
 | 1085 |       UpdateGramSchOldActualPsiNo(P,Psi);
 | 
|---|
 | 1086 |       ControlNativeDensity(P);
 | 
|---|
 | 1087 |       MinimiseOccupied(P, &Stop, &SuperStop);
 | 
|---|
 | 1088 |       if (!I->StructOpt) {
 | 
|---|
| [6df7db] | 1089 |         if (((P->Call.ReadSrcFiles != DoReadAllSrcDensities) && (P->Call.ReadSrcFiles != DoReadOccupiedSrcDensities)) || (!ParseWannierFile(P))) { // only localize and store if they have just been minimised (hence don't come solely from file), otherwise read stored values from file (if possible)
 | 
|---|
| [5712cb] | 1090 |           SpeedMeasure(P, WannierTime, StartTimeDo);
 | 
|---|
 | 1091 |           ComputeMLWF(P);   // localization of orbitals
 | 
|---|
 | 1092 |           SpeedMeasure(P, WannierTime, StopTimeDo);
 | 
|---|
 | 1093 |     //      if (!TestReadnWriteSrcDensity(P,Occupied))
 | 
|---|
 | 1094 |     //        Error(SomeError,"TestReadnWriteSrcDensity failed!");
 | 
|---|
 | 1095 |         }
 | 
|---|
| [6df7db] | 1096 |         // join Wannier orbital to groups with common centres under certain conditions
 | 
|---|
 | 1097 |         //debug (P,"Changing Wannier Centres according to CommonWannier");
 | 
|---|
 | 1098 |         ChangeWannierCentres(P); 
 | 
|---|
 | 1099 |         OutputVisSrcFiles(P, Occupied); // rewrite now localized orbitals
 | 
|---|
 | 1100 |         
 | 
|---|
 | 1101 | 
 | 
|---|
| [a0bcf1] | 1102 | //      // plot psi cuts
 | 
|---|
 | 1103 | //      for (i=0; i < Psi->MaxPsiOfType; i++)  // go through all wave functions (here without the extra ones for each process)
 | 
|---|
 | 1104 | //        if ((Psi->AllPsiStatus[i].PsiType == Occupied) && (Psi->AllPsiStatus[i].my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi))
 | 
|---|
 | 1105 | //          for (j=0;j<NDIM;j++) {
 | 
|---|
 | 1106 | //              //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]);
 | 
|---|
 | 1107 | //              CalculateOneDensityR(Lat, R->LevS, R->Lev0->Dens, R->LevS->LPsi->LocalPsi[Psi->AllPsiStatus[i].MyLocalNo], R->Lev0->Dens->DensityArray[ActualDensity], R->FactorDensityR, 0);
 | 
|---|
 | 1108 | //              PlotSrcPlane(P, j, Lat->Psi.AddData[Psi->AllPsiStatus[i].MyLocalNo].WannierCentre[j], Psi->AllPsiStatus[i].MyGlobalNo, R->Lev0->Dens->DensityArray[ActualDensity]);
 | 
|---|
 | 1109 | //            }
 | 
|---|
 | 1110 |           
 | 
|---|
| [5712cb] | 1111 |         // unoccupied calc 
 | 
|---|
 | 1112 |         if (R->DoUnOccupied) {
 | 
|---|
 | 1113 |           MinimiseUnoccupied(P, &Stop, &SuperStop);
 | 
|---|
| [a0bcf1] | 1114 |         }
 | 
|---|
| [5712cb] | 1115 |         // hamiltonian
 | 
|---|
 | 1116 |         CalculateHamiltonian(P);  // lambda_{kl} needed (and for bandgap after UnOccupied)
 | 
|---|
 | 1117 | 
 | 
|---|
 | 1118 |         //TestSawtooth(P, 0);
 | 
|---|
 | 1119 |         //TestSawtooth(P, 1);
 | 
|---|
 | 1120 |         //TestSawtooth(P, 2);
 | 
|---|
 | 1121 |                 
 | 
|---|
 | 1122 |         // perturbed calc
 | 
|---|
 | 1123 |         if ((R->DoPerturbation)) { // && R->LevSNo <= R->InitLevSNo) {
 | 
|---|
 | 1124 |           AllocCurrentDensity(R->Lev0->Dens);// lock current density arrays
 | 
|---|
 | 1125 |           MinimisePerturbed(P, &Stop, &SuperStop); // herein InitDensityCalculation() is called, thus no need to call it beforehand
 | 
|---|
 | 1126 | 
 | 
|---|
 | 1127 |           SpeedMeasure(P, CurrDensTime, StartTimeDo);
 | 
|---|
 | 1128 |           if (SuperStop != 1) {
 | 
|---|
 | 1129 |             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
 | 
|---|
 | 1130 |               R->DoFullCurrent = 1; // set to 1 if it was 2 but Check...() yielded necessity
 | 
|---|
| [963310a] | 1131 |               //debug(P,"Filling with Delta j ...");
 | 
|---|
| [5712cb] | 1132 |               FillDeltaCurrentDensity(P);
 | 
|---|
| [963310a] | 1133 |             }
 | 
|---|
| [5712cb] | 1134 |           }
 | 
|---|
 | 1135 |           SpeedMeasure(P, CurrDensTime, StopTimeDo);
 | 
|---|
 | 1136 |           TestCurrent(P,0);
 | 
|---|
 | 1137 |           TestCurrent(P,1);
 | 
|---|
 | 1138 |           TestCurrent(P,2);
 | 
|---|
| [963310a] | 1139 |           if (F->DoOutCurr && R->Lev0->LevelNo == 0)  // only output in uppermost level)
 | 
|---|
| [5712cb] | 1140 |             OutputCurrentDensity(P);
 | 
|---|
| [963310a] | 1141 |           if (R->VectorPlane != -1)
 | 
|---|
| [5712cb] | 1142 |             PlotVectorPlane(P,R->VectorPlane,R->VectorCut);
 | 
|---|
 | 1143 |           CalculateMagneticSusceptibility(P);
 | 
|---|
 | 1144 |           debug(P,"Normal calculation of shielding over R-space");
 | 
|---|
| [cc9c36] | 1145 |           CalculateMagneticMoment(P);
 | 
|---|
| [5712cb] | 1146 |           CalculateChemicalShieldingByReciprocalCurrentDensity(P);
 | 
|---|
 | 1147 |           SpeedMeasure(P, CurrDensTime, StopTimeDo);
 | 
|---|
 | 1148 |           DisAllocCurrentDensity(R->Lev0->Dens);    // unlock current density arrays
 | 
|---|
 | 1149 |         }  // end of if perturbation
 | 
|---|
| [a0bcf1] | 1150 |         InitDensityCalculation(P);  // all unperturbed(!) wave functions've "changed" from ComputeMLWF(), thus reinit density
 | 
|---|
| [5712cb] | 1151 |       } else  // end of if StructOpt or MaxOuterStep
 | 
|---|
 | 1152 |                 OutputVisSrcFiles(P, Occupied); // in structopt or MD write for every level
 | 
|---|
| [a0bcf1] | 1153 | 
 | 
|---|
| [5712cb] | 1154 |       if ((!I->StructOpt) && (!R->MaxOuterStep))  // print intermediate levels energy results if we don't do MD or StructOpt
 | 
|---|
| [963310a] | 1155 |             EnergyOutput(P, 1);
 | 
|---|
| [5712cb] | 1156 |       // next level
 | 
|---|
 | 1157 |       ChangeToLevUp(P, &Stop);
 | 
|---|
 | 1158 |       //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!"); }      
 | 
|---|
 | 1159 |       LevS = R->LevS; // re-set pointer that's changed from LevUp
 | 
|---|
 | 1160 |     }
 | 
|---|
 | 1161 |           SpeedMeasure(P, SimTime, StopTimeDo);
 | 
|---|
 | 1162 |           ControlNativeDensity(P);
 | 
|---|
| [2026d4] | 1163 |           TestGramSch(P,LevS,Psi, Occupied);
 | 
|---|
 | 1164 |           // necessary for correct ionic forces ...
 | 
|---|
 | 1165 |                 SpeedMeasure(P, LocFTime, StartTimeDo);
 | 
|---|
 | 1166 |                 CalculateIonLocalForce(P);
 | 
|---|
 | 1167 |                 SpeedMeasure(P, LocFTime, StopTimeDo);
 | 
|---|
 | 1168 |                 SpeedMeasure(P, NonLocFTime, StartTimeDo);
 | 
|---|
 | 1169 |                 CalculateIonNonLocalForce(P);
 | 
|---|
 | 1170 |                 SpeedMeasure(P, NonLocFTime, StopTimeDo);
 | 
|---|
 | 1171 |                 CalculateEwald(P, 1);
 | 
|---|
 | 1172 |                 CalculateIonForce(P);
 | 
|---|
| [a0bcf1] | 1173 |   }
 | 
|---|
| [5712cb] | 1174 |   if (P->Call.ForcesFile != NULL) { // if we parse forces from file, values are written over (in case of DoOutVis)
 | 
|---|
 | 1175 |     fprintf(stderr, "Parsing Forces from file.\n");
 | 
|---|
| [53b5b6] | 1176 |     ParseIonForce(P); 
 | 
|---|
| [5712cb] | 1177 |     //CalculateIonForce(P);
 | 
|---|
| [2026d4] | 1178 |   }
 | 
|---|
| [a0bcf1] | 1179 |   CorrectForces(P);
 | 
|---|
 | 1180 |   // ... on output of densities
 | 
|---|
 | 1181 |   if (F->DoOutOrbitals) { // output of each orbital
 | 
|---|
 | 1182 |     debug(P,"OutputVisAllOrbital");
 | 
|---|
 | 1183 |     OutputVisAllOrbital(P,0,1,Occupied);
 | 
|---|
 | 1184 |   }
 | 
|---|
| [5712cb] | 1185 |   //OutputNorm(P);
 | 
|---|
| [2026d4] | 1186 |   //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] | 1187 |   //OutputVis(P, P->R.Lev0->Dens->DensityArray[TotalDensity]);
 | 
|---|
| [2026d4] | 1188 |   /*TestGramSch(P, R->LevS, &P->Lat.Psi); */
 | 
|---|
| [a0bcf1] | 1189 |   GetOuterStop(P);
 | 
|---|
 | 1190 |   //fprintf(stderr,"(%i) DoubleG: %p, CArray[22]: %p, OldLocalPsi: %p\n", P->Par.me, R->LevS->DoubleG, R->Lev0->Dens->DensityCArray[22], R->LevS->LPsi->OldLocalPsi); 
 | 
|---|
 | 1191 |   if (SuperStop) OuterStop = 1;
 | 
|---|
 | 1192 |   return OuterStop;
 | 
|---|
 | 1193 | }
 | 
|---|
 | 1194 | 
 | 
|---|
| [2026d4] | 1195 | /** Checks whether the given positions \a *v have changed wrt stored in IonData structure.
 | 
|---|
 | 1196 |  * \param *P Problem at hand
 | 
|---|
 | 1197 |  * \param *v gsl_vector storing new positions
 | 
|---|
 | 1198 |  */
 | 
|---|
 | 1199 | int CheckForChangedPositions(struct Problem *P, const gsl_vector *v)
 | 
|---|
 | 1200 | {
 | 
|---|
 | 1201 |   struct Ions *I = &P->Ion;
 | 
|---|
 | 1202 |   int is,ia,k, index=0;
 | 
|---|
 | 1203 |   int diff = 0;
 | 
|---|
 | 1204 |   double *R_ion;
 | 
|---|
 | 1205 |   for (is=0; is < I->Max_Types; is++) // for all elements
 | 
|---|
 | 1206 |     for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {  // for all ions of element
 | 
|---|
 | 1207 |       R_ion = &I->I[is].R[NDIM*ia];
 | 
|---|
 | 1208 |       for (k=0;k<NDIM;k++) { // for all dimensions
 | 
|---|
 | 1209 |         if (fabs(R_ion[k] - gsl_vector_get (v, index++)) > MYEPSILON)
 | 
|---|
 | 1210 |           diff++;
 | 
|---|
 | 1211 |       }      
 | 
|---|
 | 1212 |     }
 | 
|---|
 | 1213 |   return diff;  
 | 
|---|
 | 1214 | }
 | 
|---|
 | 1215 | 
 | 
|---|
 | 1216 | /** Wrapper for CalculateForce() for simplex minimisation of total energy.
 | 
|---|
 | 1217 |  * \param *v vector with degrees of freedom
 | 
|---|
 | 1218 |  * \param *params additional arguments, here Problem at hand
 | 
|---|
 | 1219 |  */
 | 
|---|
 | 1220 | double StructOpt_func(const gsl_vector *v, void *params)
 | 
|---|
 | 1221 | {
 | 
|---|
 | 1222 |   struct Problem *P = (struct Problem *)params;
 | 
|---|
 | 1223 |   struct RunStruct *R = &P->R;
 | 
|---|
 | 1224 |   struct Ions *I = &P->Ion;
 | 
|---|
 | 1225 |   struct Energy *E = P->Lat.E;
 | 
|---|
 | 1226 |   int i;
 | 
|---|
 | 1227 |   double *R_ion, *R_old, *R_old_old;//, *FIon;
 | 
|---|
 | 1228 |   //double norm = 0.;
 | 
|---|
 | 1229 |   int is,ia,k,index = 0;
 | 
|---|
 | 1230 |   int OuterStop;
 | 
|---|
 | 1231 |   double diff = 0., tmp;
 | 
|---|
 | 1232 |   debug (P, "StructOpt_func");
 | 
|---|
 | 1233 |   if (CheckForChangedPositions(P,v)) {
 | 
|---|
 | 1234 |     // update ion positions from vector coordinates
 | 
|---|
 | 1235 |     for (is=0; is < I->Max_Types; is++) // for all elements
 | 
|---|
 | 1236 |       for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {  // for all ions of element
 | 
|---|
 | 1237 |         R_ion = &I->I[is].R[NDIM*ia];
 | 
|---|
 | 1238 |         R_old = &I->I[is].R_old[NDIM*ia];
 | 
|---|
 | 1239 |         R_old_old = &I->I[is].R_old_old[NDIM*ia];
 | 
|---|
 | 1240 |         tmp = 0.;
 | 
|---|
 | 1241 |         for (k=0;k<NDIM;k++) { // for all dimensions
 | 
|---|
 | 1242 |           R_old_old[k] = R_old[k];
 | 
|---|
 | 1243 |           R_old[k] = R_ion[k];
 | 
|---|
 | 1244 |           tmp += (R_ion[k]-gsl_vector_get (v, index))*(R_ion[k]-gsl_vector_get (v, index));
 | 
|---|
 | 1245 |           R_ion[k] = gsl_vector_get (v, index++);        
 | 
|---|
 | 1246 |         }
 | 
|---|
 | 1247 |         diff += sqrt(tmp);
 | 
|---|
 | 1248 |       }  
 | 
|---|
 | 1249 |     if (P->Call.out[ValueOut]) fprintf(stderr,"(%i) Summed Difference to former position %lg\n", P->Par.me, diff);
 | 
|---|
 | 1250 |     // recalculate ionic forces (do electronic minimisation)
 | 
|---|
 | 1251 |     R->OuterStep++;
 | 
|---|
 | 1252 |     R->NewRStep++;
 | 
|---|
 | 1253 |     UpdateWaveAfterIonMove(P);
 | 
|---|
 | 1254 |     for (i=MAXOLD-1; i > 0; i--)    // store away old energies
 | 
|---|
 | 1255 |       E->TotalEnergyOuter[i] = E->TotalEnergyOuter[i-1];  
 | 
|---|
 | 1256 |     UpdateToNewWaves(P);
 | 
|---|
 | 1257 |     E->TotalEnergyOuter[0] = E->TotalEnergy[0];
 | 
|---|
 | 1258 |     OuterStop = CalculateForce(P);
 | 
|---|
 | 1259 |     //UpdateIonsU(P);
 | 
|---|
 | 1260 |     //CorrectVelocity(P);
 | 
|---|
 | 1261 |     //CalculateEnergyIonsU(P);
 | 
|---|
 | 1262 |   /*  if ((P->R.ScaleTempStep > 0) && ((R->OuterStep-1) % P->R.ScaleTempStep == 0))
 | 
|---|
 | 1263 |       ScaleTemp(P);*/
 | 
|---|
 | 1264 |     if ((R->OuterStep-1) % P->R.OutSrcStep == 0) 
 | 
|---|
 | 1265 |       OutputVisSrcFiles(P, Occupied);
 | 
|---|
 | 1266 |     if ((R->OuterStep-1) % P->R.OutVisStep == 0) {      
 | 
|---|
 | 1267 |   /*      // recalculate density for the specific wave function ...
 | 
|---|
 | 1268 |         CalculateOneDensityR(Lat, LevS, Dens0, PsiDat, Dens0->DensityArray[ActualDensity], R->FactorDensityR, 0);
 | 
|---|
 | 1269 |         // ... and output (wherein ActualDensity is used instead of TotalDensity)
 | 
|---|
 | 1270 |         OutputVis(P);
 | 
|---|
 | 1271 |         OutputIonForce(P);
 | 
|---|
 | 1272 |         EnergyOutput(P, 1);*/
 | 
|---|
 | 1273 |     }
 | 
|---|
 | 1274 |   }
 | 
|---|
 | 1275 |   if (P->Par.me == 0) fprintf(stderr,"(%i) TE %e\n",P->Par.me, E->TotalEnergy[0]);
 | 
|---|
 | 1276 |   return E->TotalEnergy[0];
 | 
|---|
 | 1277 | }
 | 
|---|
 | 1278 | 
 | 
|---|
| [a0bcf1] | 1279 | /** Wrapper for CalculateForce() for simplex minimisation of ionic forces.
 | 
|---|
 | 1280 |  * \param *v vector with degrees of freedom
 | 
|---|
 | 1281 |  * \param *params additional arguments, here Problem at hand
 | 
|---|
 | 1282 |  */
 | 
|---|
| [f915e1] | 1283 | double StructOpt_f(const gsl_vector *v, void *params)
 | 
|---|
| [a0bcf1] | 1284 | {
 | 
|---|
 | 1285 |   struct Problem *P = (struct Problem *)params;
 | 
|---|
 | 1286 |   struct RunStruct *R = &P->R;
 | 
|---|
 | 1287 |   struct Ions *I = &P->Ion;
 | 
|---|
 | 1288 |   struct Energy *E = P->Lat.E;
 | 
|---|
 | 1289 |   int i;
 | 
|---|
| [27a5bf] | 1290 |   double *R_ion, *R_old, *R_old_old;//, *FIon;
 | 
|---|
 | 1291 |   //double norm = 0.;
 | 
|---|
 | 1292 |   int is,ia,k,index = 0;
 | 
|---|
| [a0bcf1] | 1293 |   int OuterStop;
 | 
|---|
| [27a5bf] | 1294 |   double diff = 0., tmp;
 | 
|---|
 | 1295 |   //debug (P, "StructOpt_f");
 | 
|---|
 | 1296 |   if (CheckForChangedPositions(P,v)) {
 | 
|---|
 | 1297 |     // update ion positions from vector coordinates
 | 
|---|
 | 1298 |     for (is=0; is < I->Max_Types; is++) // for all elements
 | 
|---|
 | 1299 |       for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {  // for all ions of element
 | 
|---|
 | 1300 |         R_ion = &I->I[is].R[NDIM*ia];
 | 
|---|
 | 1301 |         R_old = &I->I[is].R_old[NDIM*ia];
 | 
|---|
 | 1302 |         R_old_old = &I->I[is].R_old_old[NDIM*ia];
 | 
|---|
 | 1303 |         tmp = 0.;
 | 
|---|
 | 1304 |         for (k=0;k<NDIM;k++) { // for all dimensions
 | 
|---|
 | 1305 |           R_old_old[k] = R_old[k];
 | 
|---|
 | 1306 |           R_old[k] = R_ion[k];
 | 
|---|
 | 1307 |           tmp += (R_ion[k]-gsl_vector_get (v, index))*(R_ion[k]-gsl_vector_get (v, index));
 | 
|---|
 | 1308 |           R_ion[k] = gsl_vector_get (v, index++);        
 | 
|---|
 | 1309 |         }
 | 
|---|
 | 1310 |         diff += sqrt(tmp);
 | 
|---|
 | 1311 |       }  
 | 
|---|
 | 1312 |     if (P->Call.out[ValueOut]) fprintf(stderr,"(%i) Summed Difference to former position %lg\n", P->Par.me, diff);
 | 
|---|
 | 1313 |     // recalculate ionic forces (do electronic minimisation)
 | 
|---|
 | 1314 |     //R->OuterStep++;
 | 
|---|
 | 1315 |     R->NewRStep++;
 | 
|---|
 | 1316 |     UpdateWaveAfterIonMove(P);
 | 
|---|
 | 1317 |     for (i=MAXOLD-1; i > 0; i--)    // store away old energies
 | 
|---|
 | 1318 |       E->TotalEnergyOuter[i] = E->TotalEnergyOuter[i-1];  
 | 
|---|
 | 1319 |     UpdateToNewWaves(P);
 | 
|---|
 | 1320 |     E->TotalEnergyOuter[0] = E->TotalEnergy[0];
 | 
|---|
 | 1321 |     OuterStop = CalculateForce(P);
 | 
|---|
 | 1322 |     //UpdateIonsU(P);
 | 
|---|
 | 1323 |     //CorrectVelocity(P);
 | 
|---|
 | 1324 |     //CalculateEnergyIonsU(P);
 | 
|---|
 | 1325 |   /*  if ((P->R.ScaleTempStep > 0) && ((R->OuterStep-1) % P->R.ScaleTempStep == 0))
 | 
|---|
 | 1326 |       ScaleTemp(P);*/
 | 
|---|
 | 1327 |     if ((R->OuterStep-1) % P->R.OutSrcStep == 0) 
 | 
|---|
 | 1328 |       OutputVisSrcFiles(P, Occupied);
 | 
|---|
 | 1329 |     /*if ((R->OuterStep-1) % P->R.OutVisStep == 0) {      
 | 
|---|
 | 1330 |         // recalculate density for the specific wave function ...
 | 
|---|
 | 1331 |         CalculateOneDensityR(Lat, LevS, Dens0, PsiDat, Dens0->DensityArray[ActualDensity], R->FactorDensityR, 0);
 | 
|---|
 | 1332 |         // ... and output (wherein ActualDensity is used instead of TotalDensity)
 | 
|---|
 | 1333 |         OutputVis(P);
 | 
|---|
 | 1334 |         OutputIonForce(P);
 | 
|---|
 | 1335 |         EnergyOutput(P, 1);
 | 
|---|
 | 1336 |     }*/
 | 
|---|
| [a0bcf1] | 1337 |   }
 | 
|---|
| [27a5bf] | 1338 |   GetOuterStop(P);
 | 
|---|
 | 1339 |   //if (P->Call.out[LeaderOut] && (P->Par.me == 0)) fprintf(stderr,"(%i) Absolute Force summed over all Ions %e\n",P->Par.me, norm);
 | 
|---|
 | 1340 |   return R->MeanForce[0];
 | 
|---|
 | 1341 |   //if (P->Call.out[LeaderOut] && (P->Par.me == 0)) fprintf(stderr,"(%i) Struct_optf returning: %lg\n",P->Par.me,E->TotalEnergy[0]); 
 | 
|---|
 | 1342 |   //return E->TotalEnergy[0];
 | 
|---|
| [a0bcf1] | 1343 | }
 | 
|---|
 | 1344 | 
 | 
|---|
| [f915e1] | 1345 | void StructOpt_df(const gsl_vector *v, void *params, gsl_vector *df) 
 | 
|---|
| [a0bcf1] | 1346 | {
 | 
|---|
 | 1347 |   struct Problem *P = (struct Problem *)params;
 | 
|---|
 | 1348 |   struct Ions *I = &P->Ion;
 | 
|---|
 | 1349 |   double *FIon;
 | 
|---|
| [27a5bf] | 1350 |   int is,ia,k, index=0;
 | 
|---|
 | 1351 |   //debug (P, "StructOpt_df");
 | 
|---|
 | 1352 |   // look through coordinate vector if positions have changed sind last StructOpt_f call
 | 
|---|
 | 1353 |   if (CheckForChangedPositions(P,v)) {// if so, recalc to update forces
 | 
|---|
 | 1354 |     debug (P, "Calling StructOpt_f to update");
 | 
|---|
 | 1355 |     StructOpt_f(v, params);
 | 
|---|
 | 1356 |   }
 | 
|---|
| [a0bcf1] | 1357 |   for (is=0; is < I->Max_Types; is++) // for all elements
 | 
|---|
 | 1358 |     for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {  // for all ions of element
 | 
|---|
 | 1359 |       FIon = &I->I[is].FIon[NDIM*ia];
 | 
|---|
 | 1360 |       for (k=0;k<NDIM;k++) { // for all dimensions
 | 
|---|
 | 1361 |         gsl_vector_set (df, index++, FIon[k]);
 | 
|---|
 | 1362 |       }      
 | 
|---|
| [27a5bf] | 1363 |     } 
 | 
|---|
 | 1364 |   if (P->Call.out[LeaderOut] && (P->Par.me == 0)) {
 | 
|---|
 | 1365 |     fprintf(stderr,"(%i) Struct_Optdf returning",P->Par.me); 
 | 
|---|
 | 1366 |     gsl_vector_fprintf(stderr, df, "%lg");
 | 
|---|
 | 1367 |   } 
 | 
|---|
| [a0bcf1] | 1368 | }
 | 
|---|
 | 1369 | 
 | 
|---|
| [f915e1] | 1370 | void StructOpt_fdf (const gsl_vector *x, void *params, double *f, gsl_vector *df)
 | 
|---|
| [a0bcf1] | 1371 | {
 | 
|---|
| [f915e1] | 1372 |   *f = StructOpt_f(x, params);
 | 
|---|
 | 1373 |   StructOpt_df(x, params, df);
 | 
|---|
| [a0bcf1] | 1374 | }
 | 
|---|
 | 1375 | 
 | 
|---|
 | 1376 | 
 | 
|---|
 | 1377 | /** CG implementation for the structure optimization.
 | 
|---|
 | 1378 |  * We follow the example from the GSL manual.
 | 
|---|
 | 1379 |  * \param *P Problem at hand
 | 
|---|
 | 1380 |  */
 | 
|---|
 | 1381 | void UpdateIon_PRCG(struct Problem *P)
 | 
|---|
 | 1382 | {
 | 
|---|
| [27a5bf] | 1383 |   //struct RunStruct *Run = &P->R;
 | 
|---|
| [a0bcf1] | 1384 |   struct Ions *I = &P->Ion;
 | 
|---|
 | 1385 |   size_t np = NDIM*I->Max_TotalIons;  // d.o.f = number of ions times number of dimensions
 | 
|---|
 | 1386 |   int is, ia, k, index;
 | 
|---|
 | 1387 |   double *R;
 | 
|---|
 | 1388 |  
 | 
|---|
 | 1389 |   const gsl_multimin_fdfminimizer_type *T;
 | 
|---|
 | 1390 |   gsl_multimin_fdfminimizer *s;
 | 
|---|
 | 1391 |   gsl_vector *x;
 | 
|---|
 | 1392 |   gsl_multimin_function_fdf minex_func;
 | 
|---|
 | 1393 | 
 | 
|---|
 | 1394 |   size_t iter = 0;
 | 
|---|
 | 1395 |   int status;
 | 
|---|
 | 1396 | 
 | 
|---|
 | 1397 |   /* Starting point */
 | 
|---|
 | 1398 |   x = gsl_vector_alloc (np);
 | 
|---|
| [27a5bf] | 1399 |   //fprintf(stderr,"(%i) d.o.f. = %i\n", P->Par.me, (int)np);
 | 
|---|
| [a0bcf1] | 1400 | 
 | 
|---|
 | 1401 |   index=0;
 | 
|---|
 | 1402 |   for (is=0; is < I->Max_Types; is++) // for all elements
 | 
|---|
 | 1403 |     for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {  // for all ions of element
 | 
|---|
 | 1404 |       R = &I->I[is].R[NDIM*ia];
 | 
|---|
 | 1405 |       for (k=0;k<NDIM;k++) // for all dimensions
 | 
|---|
 | 1406 |         gsl_vector_set (x, index++, R[k]);      
 | 
|---|
 | 1407 |     }
 | 
|---|
 | 1408 | 
 | 
|---|
 | 1409 |   /* Initialize method and iterate */
 | 
|---|
| [f915e1] | 1410 |   minex_func.f = &StructOpt_f;
 | 
|---|
 | 1411 |   minex_func.df = &StructOpt_df;
 | 
|---|
 | 1412 |   minex_func.fdf = &StructOpt_fdf;
 | 
|---|
| [a0bcf1] | 1413 |   minex_func.n = np;
 | 
|---|
 | 1414 |   minex_func.params = (void *)P;
 | 
|---|
 | 1415 |   
 | 
|---|
 | 1416 |   T = gsl_multimin_fdfminimizer_conjugate_pr;
 | 
|---|
 | 1417 |   s = gsl_multimin_fdfminimizer_alloc (T, np);
 | 
|---|
 | 1418 | 
 | 
|---|
| [27a5bf] | 1419 |   gsl_multimin_fdfminimizer_set (s, &minex_func, x, 0.1, 0.001);
 | 
|---|
| [a0bcf1] | 1420 | 
 | 
|---|
| [27a5bf] | 1421 |   fprintf(stderr,"(%i) Commencing Structure optimization with PRCG: dof %d\n", P->Par.me,(int)np);
 | 
|---|
| [a0bcf1] | 1422 |   do {
 | 
|---|
 | 1423 |     iter++;
 | 
|---|
 | 1424 |     status = gsl_multimin_fdfminimizer_iterate(s);
 | 
|---|
 | 1425 | 
 | 
|---|
 | 1426 |     if (status)
 | 
|---|
 | 1427 |       break;
 | 
|---|
 | 1428 |     
 | 
|---|
| [27a5bf] | 1429 |     status = gsl_multimin_test_gradient (s->gradient, 1e-2);
 | 
|---|
| [a0bcf1] | 1430 |        
 | 
|---|
 | 1431 |     if (status == GSL_SUCCESS)
 | 
|---|
 | 1432 |       if (P->Par.me == 0) fprintf (stderr,"(%i) converged to minimum at\n", P->Par.me);
 | 
|---|
 | 1433 | 
 | 
|---|
| [53b5b6] | 1434 |     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] | 1435 |     if ((P->Call.out[NormalOut]) && (P->Par.me == 0)) fprintf (stderr, "(%i) %5d %10.5f\n", P->Par.me, (int)iter, s->f);
 | 
|---|
| [27a5bf] | 1436 |     //gsl_vector_fprintf(stderr, s->dx, "%lg");
 | 
|---|
 | 1437 |     OutputVis(P, P->R.Lev0->Dens->DensityArray[TotalDensity]);
 | 
|---|
| [774ae8] | 1438 |     OutputIonCoordinates(P, 0);
 | 
|---|
| [27a5bf] | 1439 |     P->R.StructOptStep++;
 | 
|---|
 | 1440 |   } while ((status == GSL_CONTINUE) && (P->R.StructOptStep < P->R.MaxStructOptStep));
 | 
|---|
| [a0bcf1] | 1441 | 
 | 
|---|
 | 1442 |   gsl_vector_free(x);
 | 
|---|
 | 1443 |   gsl_multimin_fdfminimizer_free (s);
 | 
|---|
 | 1444 | }
 | 
|---|
 | 1445 | 
 | 
|---|
| [27a5bf] | 1446 | /** Simplex implementation for the structure optimization.
 | 
|---|
 | 1447 |  * We follow the example from the GSL manual.
 | 
|---|
 | 1448 |  * \param *P Problem at hand
 | 
|---|
 | 1449 |  */
 | 
|---|
 | 1450 | void UpdateIon_Simplex(struct Problem *P)
 | 
|---|
 | 1451 | {
 | 
|---|
 | 1452 |   struct RunStruct *Run = &P->R;
 | 
|---|
 | 1453 |   struct Ions *I = &P->Ion;
 | 
|---|
 | 1454 |   size_t np = NDIM*I->Max_TotalIons;  // d.o.f = number of ions times number of dimensions
 | 
|---|
 | 1455 |   int is, ia, k, index;
 | 
|---|
 | 1456 |   double *R;
 | 
|---|
 | 1457 |  
 | 
|---|
 | 1458 |   const gsl_multimin_fminimizer_type *T;
 | 
|---|
 | 1459 |   gsl_multimin_fminimizer *s;
 | 
|---|
 | 1460 |   gsl_vector *x, *ss;
 | 
|---|
 | 1461 |   gsl_multimin_function minex_func;
 | 
|---|
 | 1462 | 
 | 
|---|
 | 1463 |   size_t iter = 0;
 | 
|---|
 | 1464 |   int status;
 | 
|---|
 | 1465 |   double size;
 | 
|---|
 | 1466 | 
 | 
|---|
 | 1467 |   ss = gsl_vector_alloc (np);
 | 
|---|
 | 1468 |   gsl_vector_set_all(ss, .2);
 | 
|---|
 | 1469 |   /* Starting point */
 | 
|---|
 | 1470 |   x = gsl_vector_alloc (np);
 | 
|---|
 | 1471 |   //fprintf(stderr,"(%i) d.o.f. = %i\n", P->Par.me, (int)np);
 | 
|---|
 | 1472 | 
 | 
|---|
 | 1473 |   index=0;
 | 
|---|
 | 1474 |   for (is=0; is < I->Max_Types; is++) // for all elements
 | 
|---|
 | 1475 |     for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {  // for all ions of element
 | 
|---|
 | 1476 |       R = &I->I[is].R[NDIM*ia];
 | 
|---|
 | 1477 |       for (k=0;k<NDIM;k++) // for all dimensions
 | 
|---|
 | 1478 |         gsl_vector_set (x, index++, R[k]);      
 | 
|---|
 | 1479 |     }
 | 
|---|
 | 1480 | 
 | 
|---|
 | 1481 |   /* Initialize method and iterate */
 | 
|---|
 | 1482 |   minex_func.f = &StructOpt_f;
 | 
|---|
 | 1483 |   minex_func.n = np;
 | 
|---|
 | 1484 |   minex_func.params = (void *)P;
 | 
|---|
 | 1485 |   
 | 
|---|
 | 1486 |   T = gsl_multimin_fminimizer_nmsimplex;
 | 
|---|
 | 1487 |   s = gsl_multimin_fminimizer_alloc (T, np);
 | 
|---|
 | 1488 | 
 | 
|---|
 | 1489 |   gsl_multimin_fminimizer_set (s, &minex_func, x, ss);
 | 
|---|
 | 1490 | 
 | 
|---|
 | 1491 |   fprintf(stderr,"(%i) Commencing Structure optimization with NM simplex: dof %d\n", P->Par.me, (int)np);
 | 
|---|
 | 1492 |   do {
 | 
|---|
 | 1493 |     iter++;
 | 
|---|
 | 1494 |     status = gsl_multimin_fminimizer_iterate(s);
 | 
|---|
 | 1495 | 
 | 
|---|
 | 1496 |     if (status)
 | 
|---|
 | 1497 |       break;
 | 
|---|
 | 1498 |     
 | 
|---|
 | 1499 |     size = gsl_multimin_fminimizer_size (s);
 | 
|---|
 | 1500 |     status = gsl_multimin_test_size (size, 1e-4);
 | 
|---|
 | 1501 |        
 | 
|---|
 | 1502 |     if (status == GSL_SUCCESS)
 | 
|---|
 | 1503 |       if (P->Par.me == 0) fprintf (stderr,"(%i) converged to minimum at\n", P->Par.me);
 | 
|---|
 | 1504 | 
 | 
|---|
 | 1505 |     if (P->Call.out[MinOut]) fprintf(stderr,"(%i) Commencing '%s' step %i ... \n",P->Par.me, gsl_multimin_fminimizer_name(s), P->R.StructOptStep);
 | 
|---|
 | 1506 |     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);
 | 
|---|
 | 1507 |     OutputVis(P, P->R.Lev0->Dens->DensityArray[TotalDensity]);
 | 
|---|
| [774ae8] | 1508 |     OutputIonCoordinates(P, 0);
 | 
|---|
| [27a5bf] | 1509 |     P->R.StructOptStep++;
 | 
|---|
 | 1510 |   } while ((status == GSL_CONTINUE) && (Run->OuterStep < Run->MaxOuterStep));
 | 
|---|
 | 1511 | 
 | 
|---|
 | 1512 |   gsl_vector_free(x);
 | 
|---|
 | 1513 |   gsl_vector_free(ss);
 | 
|---|
 | 1514 |   gsl_multimin_fminimizer_free (s);
 | 
|---|
 | 1515 | }
 | 
|---|
 | 1516 | 
 | 
|---|
| [f915e1] | 1517 | /** Implementation of various thermostats.
 | 
|---|
 | 1518 |  * All these thermostats apply an additional force which has the following forms:
 | 
|---|
 | 1519 |  * -# Woodcock
 | 
|---|
 | 1520 |  *  \f$p_i \rightarrow \sqrt{\frac{T_0}{T}} \cdot p_i\f$
 | 
|---|
 | 1521 |  * -# Gaussian
 | 
|---|
 | 1522 |  *  \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$
 | 
|---|
 | 1523 |  * -# Langevin
 | 
|---|
 | 1524 |  *  \f$p_{i,n} \rightarrow \sqrt{1-\alpha^2} p_{i,0} + \alpha p_r\f$  
 | 
|---|
 | 1525 |  * -# Berendsen
 | 
|---|
 | 1526 |  *  \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$
 | 
|---|
 | 1527 |  * -# Nose-Hoover
 | 
|---|
 | 1528 |  *  \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$
 | 
|---|
 | 1529 |  * These Thermostats either simply rescale the velocities, thus Thermostats() should be called after UpdateIonsU(), and/or
 | 
|---|
 | 1530 |  * have a constraint force acting additionally on the ions. In the latter case, the ion speeds have to be modified
 | 
|---|
 | 1531 |  * belatedly and the constraint force set.
 | 
|---|
 | 1532 |  * \param *P Problem at hand
 | 
|---|
 | 1533 |  * \param i which of the thermostats to take: 0 - none, 1 - Woodcock, 2 - Gaussian, 3 - Langevin, 4 - Berendsen, 5 - Nose-Hoover
 | 
|---|
 | 1534 |  * \sa InitThermostat()
 | 
|---|
 | 1535 |  */
 | 
|---|
 | 1536 | void Thermostats(struct Problem *P, enum thermostats i)
 | 
|---|
 | 1537 | {
 | 
|---|
 | 1538 |   struct FileData *Files = &P->Files;
 | 
|---|
 | 1539 |   struct Ions *I = &P->Ion;
 | 
|---|
 | 1540 |   int is, ia, d;
 | 
|---|
 | 1541 |   double *U;
 | 
|---|
 | 1542 |   double a, ekin = 0.;
 | 
|---|
 | 1543 |   double E = 0., F = 0.;
 | 
|---|
 | 1544 |   double delta_alpha = 0.;
 | 
|---|
 | 1545 |   const int delta_t = P->R.delta_t;
 | 
|---|
 | 1546 |   double ScaleTempFactor;
 | 
|---|
 | 1547 |   double sigma;
 | 
|---|
 | 1548 |   gsl_rng * r;
 | 
|---|
 | 1549 |   const gsl_rng_type * T;
 | 
|---|
 | 1550 |   
 | 
|---|
 | 1551 |   // calculate current temperature
 | 
|---|
 | 1552 |   CalculateEnergyIonsU(P); // Temperature now in I->ActualTemp
 | 
|---|
 | 1553 |   ScaleTempFactor = P->R.TargetTemp/I->ActualTemp;
 | 
|---|
 | 1554 |   //if ((P->Par.me == 0) && (I->ActualTemp < MYEPSILON)) fprintf(stderr,"Thermostat: (1) I->ActualTemp = %lg",I->ActualTemp);
 | 
|---|
 | 1555 |   if (Files->MeOutMes) fprintf(Files->TemperatureFile, "%d\t%lg",P->R.OuterStep, I->ActualTemp);
 | 
|---|
 | 1556 |     
 | 
|---|
 | 1557 |   // differentating between the various thermostats
 | 
|---|
 | 1558 |   switch(i) {
 | 
|---|
 | 1559 |      case None:
 | 
|---|
 | 1560 |       debug(P, "Applying no thermostat...");
 | 
|---|
 | 1561 |       break;
 | 
|---|
 | 1562 |      case Woodcock:
 | 
|---|
 | 1563 |       if ((P->R.ScaleTempStep > 0) && ((P->R.OuterStep-1) % P->R.ScaleTempStep == 0)) {
 | 
|---|
 | 1564 |         debug(P, "Applying Woodcock thermostat...");
 | 
|---|
 | 1565 |         for (is=0; is < I->Max_Types; is++) {
 | 
|---|
 | 1566 |           a = 0.5*I->I[is].IonMass;
 | 
|---|
 | 1567 |           for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
 | 
|---|
 | 1568 |             U = &I->I[is].U[NDIM*ia];
 | 
|---|
 | 1569 |             if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
 | 
|---|
 | 1570 |               for (d=0; d<NDIM; d++) {
 | 
|---|
 | 1571 |                 U[d] *= sqrt(ScaleTempFactor);
 | 
|---|
 | 1572 |                 ekin += 0.5*I->I[is].IonMass * U[d]*U[d];
 | 
|---|
 | 1573 |               }
 | 
|---|
 | 1574 |           }
 | 
|---|
 | 1575 |         }
 | 
|---|
 | 1576 |       }
 | 
|---|
 | 1577 |       break;
 | 
|---|
 | 1578 |      case Gaussian:
 | 
|---|
 | 1579 |       debug(P, "Applying Gaussian thermostat...");
 | 
|---|
 | 1580 |       for (is=0; is < I->Max_Types; is++) { // sum up constraint constant
 | 
|---|
 | 1581 |         for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
 | 
|---|
 | 1582 |           U = &I->I[is].U[NDIM*ia];
 | 
|---|
 | 1583 |           if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
 | 
|---|
 | 1584 |             for (d=0; d<NDIM; d++) {
 | 
|---|
 | 1585 |               F += U[d] * I->I[is].FIon[d+NDIM*ia]; 
 | 
|---|
 | 1586 |               E += U[d]*U[d]*I->I[is].IonMass;
 | 
|---|
 | 1587 |             }
 | 
|---|
 | 1588 |         }
 | 
|---|
 | 1589 |       }
 | 
|---|
 | 1590 |       if (P->Call.out[ValueOut]) fprintf(stderr, "(%i) Gaussian Least Constraint constant is %lg\n", P->Par.me, F/E);
 | 
|---|
 | 1591 |       for (is=0; is < I->Max_Types; is++) { // apply constraint constant on constraint force and on velocities
 | 
|---|
 | 1592 |         for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
 | 
|---|
 | 1593 |           U = &I->I[is].U[NDIM*ia];
 | 
|---|
 | 1594 |           if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
 | 
|---|
 | 1595 |             for (d=0; d<NDIM; d++) {
 | 
|---|
 | 1596 |               I->I[is].FConstraint[d+NDIM*ia] = (F/E) * (U[d]*I->I[is].IonMass);
 | 
|---|
 | 1597 |               U[d] += delta_t/I->I[is].IonMass * (I->I[is].FConstraint[d+NDIM*ia]);
 | 
|---|
 | 1598 |               ekin += 0.5*I->I[is].IonMass * U[d]*U[d];
 | 
|---|
 | 1599 |             }
 | 
|---|
 | 1600 |         }
 | 
|---|
 | 1601 |       }
 | 
|---|
 | 1602 |       break;
 | 
|---|
 | 1603 |      case Langevin:
 | 
|---|
 | 1604 |       debug(P, "Applying Langevin thermostat...");
 | 
|---|
 | 1605 |       // init random number generator
 | 
|---|
 | 1606 |       gsl_rng_env_setup();
 | 
|---|
 | 1607 |       T = gsl_rng_default;
 | 
|---|
 | 1608 |       r = gsl_rng_alloc (T);
 | 
|---|
 | 1609 |       // Go through each ion
 | 
|---|
 | 1610 |       for (is=0; is < I->Max_Types; is++) {
 | 
|---|
 | 1611 |         sigma  = sqrt(P->R.TargetTemp/I->I[is].IonMass); // sigma = (k_b T)/m (Hartree/atomicmass = atomiclength/atomictime)
 | 
|---|
 | 1612 |         for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
 | 
|---|
 | 1613 |           U = &I->I[is].U[NDIM*ia];
 | 
|---|
 | 1614 |           // throw a dice to determine whether it gets hit by a heat bath particle
 | 
|---|
 | 1615 |           if (((((rand()/(double)RAND_MAX))*P->R.TempFrequency) < 1.)) {  // (I->I[is].IMT[ia] == MoveIon) &&  even FixedIon moves, only not by other's forces
 | 
|---|
 | 1616 |             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]));
 | 
|---|
 | 1617 |             // pick three random numbers from a Boltzmann distribution around the desired temperature T for each momenta axis
 | 
|---|
 | 1618 |             for (d=0; d<NDIM; d++) {
 | 
|---|
 | 1619 |               U[d] = gsl_ran_gaussian (r, sigma);
 | 
|---|
 | 1620 |             }
 | 
|---|
 | 1621 |             if (P->Par.me == 0) fprintf(stderr,"%lg\n", sqrt(U[0]*U[0]+U[1]*U[1]+U[2]*U[2]));
 | 
|---|
 | 1622 |           }
 | 
|---|
 | 1623 |           for (d=0; d<NDIM; d++)
 | 
|---|
 | 1624 |             ekin += 0.5*I->I[is].IonMass * U[d]*U[d];
 | 
|---|
 | 1625 |         }
 | 
|---|
 | 1626 |       }
 | 
|---|
 | 1627 |       break;
 | 
|---|
 | 1628 |      case Berendsen:
 | 
|---|
 | 1629 |       debug(P, "Applying Berendsen-VanGunsteren thermostat...");
 | 
|---|
 | 1630 |       for (is=0; is < I->Max_Types; is++) {
 | 
|---|
 | 1631 |         for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
 | 
|---|
 | 1632 |           U = &I->I[is].U[NDIM*ia];
 | 
|---|
 | 1633 |           if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
 | 
|---|
 | 1634 |             for (d=0; d<NDIM; d++) {
 | 
|---|
 | 1635 |               U[d] *= sqrt(1+(P->R.delta_t/P->R.TempFrequency)*(ScaleTempFactor-1));
 | 
|---|
 | 1636 |               ekin += 0.5*I->I[is].IonMass * U[d]*U[d];
 | 
|---|
 | 1637 |             }
 | 
|---|
 | 1638 |         }
 | 
|---|
 | 1639 |       }
 | 
|---|
 | 1640 |       break;
 | 
|---|
 | 1641 |      case NoseHoover:
 | 
|---|
 | 1642 |       debug(P, "Applying Nose-Hoover thermostat...");
 | 
|---|
 | 1643 |       // dynamically evolve alpha (the additional degree of freedom)
 | 
|---|
 | 1644 |       delta_alpha = 0.;
 | 
|---|
 | 1645 |       for (is=0; is < I->Max_Types; is++) { // sum up constraint constant
 | 
|---|
 | 1646 |         for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
 | 
|---|
 | 1647 |           U = &I->I[is].U[NDIM*ia];
 | 
|---|
 | 1648 |           if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
 | 
|---|
 | 1649 |             for (d=0; d<NDIM; d++) {
 | 
|---|
 | 1650 |               delta_alpha += U[d]*U[d]*I->I[is].IonMass;
 | 
|---|
 | 1651 |             }
 | 
|---|
 | 1652 |         }
 | 
|---|
 | 1653 |       }
 | 
|---|
 | 1654 |       delta_alpha = (delta_alpha - (3.*I->Max_TotalIons+1.) * P->R.TargetTemp)/(P->R.HooverMass*Units2Electronmass);
 | 
|---|
 | 1655 |       P->R.alpha += delta_alpha*delta_t;
 | 
|---|
 | 1656 |       if (P->Par.me == 0) fprintf(stderr,"(%i) alpha = %lg * %i = %lg\n", P->Par.me, delta_alpha, delta_t, P->R.alpha);
 | 
|---|
 | 1657 |       // apply updated alpha as additional force
 | 
|---|
 | 1658 |       for (is=0; is < I->Max_Types; is++) { // apply constraint constant on constraint force and on velocities
 | 
|---|
 | 1659 |         for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
 | 
|---|
 | 1660 |           U = &I->I[is].U[NDIM*ia];
 | 
|---|
 | 1661 |           if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
 | 
|---|
 | 1662 |             for (d=0; d<NDIM; d++) {
 | 
|---|
 | 1663 |               I->I[is].FConstraint[d+NDIM*ia] = - P->R.alpha * (U[d] * I->I[is].IonMass);
 | 
|---|
 | 1664 |               U[d] += delta_t/I->I[is].IonMass * (I->I[is].FConstraint[d+NDIM*ia]);
 | 
|---|
 | 1665 |               ekin += (0.5*I->I[is].IonMass) * U[d]*U[d];
 | 
|---|
 | 1666 |             }
 | 
|---|
 | 1667 |         }
 | 
|---|
 | 1668 |       }
 | 
|---|
 | 1669 |       break;
 | 
|---|
 | 1670 |   }    
 | 
|---|
 | 1671 |   I->EKin = ekin;
 | 
|---|
 | 1672 |   I->ActualTemp = (2./(3.*I->Max_TotalIons)*I->EKin);    
 | 
|---|
 | 1673 |   //if ((P->Par.me == 0) && (I->ActualTemp < MYEPSILON)) fprintf(stderr,"Thermostat: (2) I->ActualTemp = %lg",I->ActualTemp);
 | 
|---|
 | 1674 |   if (Files->MeOutMes) { fprintf(Files->TemperatureFile, "\t%lg\n", I->ActualTemp); fflush(Files->TemperatureFile); }
 | 
|---|
 | 1675 | }
 | 
|---|
 | 1676 | 
 | 
|---|
| [a0bcf1] | 1677 | /** Does the Molecular Dynamics Calculations.
 | 
|---|
 | 1678 |  * All of the following is SpeedMeasure()'d in SimTime.
 | 
|---|
 | 1679 |  * Initialization by calling:
 | 
|---|
 | 1680 |  * -# CorrectVelocity()\n
 | 
|---|
 | 1681 |  *    Shifts center of gravity of Ions momenta, so that the cell itself remains at rest.
 | 
|---|
 | 1682 |  * -# CalculateEnergyIonsU(), SpeedMeasure()'d in TimeTypes#InitSimTime\n
 | 
|---|
 | 1683 |  *    Calculates kinetic energy of "movable" Ions.
 | 
|---|
 | 1684 |  * -# CalculateForce()\n
 | 
|---|
 | 1685 |  *    Does the minimisation, calculates densities, then energies and finally the forces.
 | 
|---|
 | 1686 |  * -# OutputVisSrcFiles()\n
 | 
|---|
 | 1687 |  *    If desired, so-far made calculations are stored to file for later restarting.
 | 
|---|
 | 1688 |  * -# OutputIonForce()\n
 | 
|---|
 | 1689 |  *    Write ion forces to file.
 | 
|---|
 | 1690 |  * -# EnergyOutput()\n
 | 
|---|
 | 1691 |  *    Write calculated energies to screen or file.
 | 
|---|
 | 1692 |  * 
 | 
|---|
 | 1693 |  * The simulation phase begins: 
 | 
|---|
 | 1694 |  * -# UpdateIonsR()\n
 | 
|---|
 | 1695 |  *    Move Ions according to the calculated force.
 | 
|---|
 | 1696 |  * -# UpdateWaveAfterIonMove()\n
 | 
|---|
 | 1697 |  *    Update wave functions by averaging LocalPsi coefficients after the Ions have been shifted.
 | 
|---|
 | 1698 |  * -# UpdateToNewWaves()\n
 | 
|---|
 | 1699 |  *    Update after wave functions have changed.
 | 
|---|
 | 1700 |  * -# CalculateForce()\n
 | 
|---|
 | 1701 |  *    Does the minimisation, calculates densities, then energies and finally the forces.
 | 
|---|
 | 1702 |  * -# UpdateIonsU()\n
 | 
|---|
 | 1703 |  *    Change ion's velocities according to the calculated acting force.
 | 
|---|
 | 1704 |  * -# CorrectVelocity()\n
 | 
|---|
 | 1705 |  *    Shifts center of gravity of Ions momenta, so that the cell itself remains at rest.
 | 
|---|
 | 1706 |  * -# CalculateEnergyIonsU()\n
 | 
|---|
 | 1707 |  *    Calculates kinetic energy of "movable" Ions.
 | 
|---|
 | 1708 |  * -# ScaleTemp()\n
 | 
|---|
 | 1709 |  *    The temperature is scaled, so the systems energy remains constant (they must not gain momenta out of nothing)
 | 
|---|
 | 1710 |  * -# OutputVisSrcFiles()\n
 | 
|---|
 | 1711 |  *    If desired, so-far made calculations are stored to file for later restarting.
 | 
|---|
 | 1712 |  * -# OutputVis()\n
 | 
|---|
 | 1713 |  *    Visulization data for OpenDX is written at certain steps if desired.
 | 
|---|
 | 1714 |  * -# OutputIonForce()\n
 | 
|---|
 | 1715 |  *    Write ion forces to file.
 | 
|---|
 | 1716 |  * -# EnergyOutput()\n
 | 
|---|
 | 1717 |  *    Write calculated energies to screen or file.
 | 
|---|
 | 1718 |  * 
 | 
|---|
 | 1719 |  * After the ground state is found:
 | 
|---|
 | 1720 |  * -# CalculateUnOccupied()\n
 | 
|---|
 | 1721 |  *    Energies of unoccupied orbitals - that have been left out completely so far -
 | 
|---|
 | 1722 |  *    are calculated.
 | 
|---|
 | 1723 |  * -# TestGramSch()\n
 | 
|---|
 | 1724 |  *    Test if orbitals are still orthogonal.
 | 
|---|
 | 1725 |  * -# CalculateHamiltonian()\n
 | 
|---|
 | 1726 |  *    Construct Hamiltonian and calculate Eigenvalues.
 | 
|---|
 | 1727 |  * -# ComputeMLWF()\n
 | 
|---|
 | 1728 |  *    Localize orbital wave functions.
 | 
|---|
 | 1729 |  * 
 | 
|---|
 | 1730 |  * \param *P Problem at hand
 | 
|---|
 | 1731 |  */
 | 
|---|
 | 1732 | void CalculateMD(struct Problem *P)
 | 
|---|
 | 1733 | {
 | 
|---|
 | 1734 |   struct RunStruct *R = &P->R;
 | 
|---|
 | 1735 |   struct Ions *I = &P->Ion;
 | 
|---|
| [27a5bf] | 1736 |   struct Energy *E = P->Lat.E;
 | 
|---|
| [a0bcf1] | 1737 |   int OuterStop = 0;
 | 
|---|
| [27a5bf] | 1738 |   int i;
 | 
|---|
 | 1739 |   
 | 
|---|
| [a0bcf1] | 1740 |   SpeedMeasure(P, SimTime, StartTimeDo);
 | 
|---|
| [27a5bf] | 1741 |   // initial calculations (bring density on BO surface and output start energies, coordinates, densities, ...)
 | 
|---|
| [a0bcf1] | 1742 |   SpeedMeasure(P, InitSimTime, StartTimeDo);
 | 
|---|
 | 1743 |   R->OuterStep = 0;
 | 
|---|
 | 1744 |   CorrectVelocity(P);
 | 
|---|
 | 1745 |   CalculateEnergyIonsU(P);
 | 
|---|
 | 1746 |   OuterStop = CalculateForce(P);
 | 
|---|
| [27a5bf] | 1747 |   //R->OuterStep++;
 | 
|---|
| [a0bcf1] | 1748 |   P->Speed.InitSteps++;
 | 
|---|
 | 1749 |   SpeedMeasure(P, InitSimTime, StopTimeDo);
 | 
|---|
| [27a5bf] | 1750 | 
 | 
|---|
| [774ae8] | 1751 |   OutputIonCoordinates(P, 1);
 | 
|---|
| [27a5bf] | 1752 |   OutputVis(P, P->R.Lev0->Dens->DensityArray[TotalDensity]);
 | 
|---|
| [a0bcf1] | 1753 |   OutputIonForce(P);
 | 
|---|
 | 1754 |   EnergyOutput(P, 1);
 | 
|---|
| [27a5bf] | 1755 |   
 | 
|---|
 | 1756 |   // if desired perform beforehand a structure relaxation/optimization
 | 
|---|
 | 1757 |   if (I->StructOpt) {
 | 
|---|
 | 1758 |     debug(P,"Commencing minimisation on ionic structure ...");
 | 
|---|
 | 1759 |     R->StructOptStep = 0;
 | 
|---|
 | 1760 |     //UpdateIon_PRCG(P);
 | 
|---|
 | 1761 |     //UpdateIon_Simplex(P);
 | 
|---|
 | 1762 |     while ((R->MeanForce[0] > 1e-4) && (R->StructOptStep < R->MaxStructOptStep)) {
 | 
|---|
 | 1763 |       R->StructOptStep++; 
 | 
|---|
| [774ae8] | 1764 |       OutputIonCoordinates(P, 1);
 | 
|---|
| [27a5bf] | 1765 |       UpdateIons(P); 
 | 
|---|
 | 1766 |       UpdateWaveAfterIonMove(P);
 | 
|---|
 | 1767 |       for (i=MAXOLD-1; i > 0; i--)    // store away old energies
 | 
|---|
 | 1768 |         E->TotalEnergyOuter[i] = E->TotalEnergyOuter[i-1];  
 | 
|---|
 | 1769 |       UpdateToNewWaves(P);
 | 
|---|
 | 1770 |       E->TotalEnergyOuter[0] = E->TotalEnergy[0]; 
 | 
|---|
 | 1771 |       OuterStop = CalculateForce(P);
 | 
|---|
 | 1772 |       CalculateEnergyIonsU(P);
 | 
|---|
 | 1773 |       if ((R->StructOptStep-1) % P->R.OutSrcStep == 0) 
 | 
|---|
 | 1774 |         OutputVisSrcFiles(P, Occupied);
 | 
|---|
 | 1775 |       if ((R->StructOptStep-1) % P->R.OutVisStep == 0) {      
 | 
|---|
 | 1776 |         OutputVis(P, P->R.Lev0->Dens->DensityArray[TotalDensity]);
 | 
|---|
 | 1777 |         OutputIonForce(P);
 | 
|---|
 | 1778 |         EnergyOutput(P, 1);
 | 
|---|
 | 1779 |       }
 | 
|---|
 | 1780 |       if (P->Par.me == 0) fprintf(stderr,"(%i) Mean force is %lg\n", P->Par.me, R->MeanForce[0]); 
 | 
|---|
 | 1781 |     }
 | 
|---|
| [774ae8] | 1782 |     OutputIonCoordinates(P, 1);
 | 
|---|
| [a0bcf1] | 1783 |   }
 | 
|---|
| [ae0078] | 1784 |   if ((I->StructOpt) && (!OuterStop) && (R->DoPerturbation)) { // do one last perturbation if desired calculation
 | 
|---|
| [a0bcf1] | 1785 |     I->StructOpt = 0;
 | 
|---|
 | 1786 |     OuterStop = CalculateForce(P);
 | 
|---|
 | 1787 |   }
 | 
|---|
| [27a5bf] | 1788 |   
 | 
|---|
 | 1789 |   // and now begin with the molecular dynamics simulation
 | 
|---|
 | 1790 |   debug(P,"Commencing MD simulation ...");
 | 
|---|
 | 1791 |   while (!OuterStop && R->OuterStep < R->MaxOuterStep) {
 | 
|---|
| [a0bcf1] | 1792 |     R->OuterStep++;
 | 
|---|
| [27a5bf] | 1793 |     if (P->Par.me == 0) {
 | 
|---|
 | 1794 |       if (R->OuterStep > 1) fprintf(stderr,"\b\b\b\b\b\b\b\b\b\b\b\b");
 | 
|---|
 | 1795 |       fprintf(stderr,"Time: %f fs\r", R->t*Atomictime2Femtoseconds);
 | 
|---|
 | 1796 |       fflush(stderr);
 | 
|---|
| [a0bcf1] | 1797 |     }
 | 
|---|
 | 1798 |     OuterStop = CalculateForce(P);
 | 
|---|
| [27a5bf] | 1799 |     P->R.t += P->R.delta_t;   // increase current time by delta_t
 | 
|---|
 | 1800 |     R->NewRStep++;
 | 
|---|
 | 1801 | 
 | 
|---|
| [a0bcf1] | 1802 |     UpdateIonsU(P);
 | 
|---|
 | 1803 |     CorrectVelocity(P);
 | 
|---|
| [27a5bf] | 1804 |     Thermostats(P, I->Thermostat);
 | 
|---|
| [a0bcf1] | 1805 |     CalculateEnergyIonsU(P);
 | 
|---|
| [27a5bf] | 1806 | 
 | 
|---|
 | 1807 |     UpdateIonsR(P);
 | 
|---|
| [774ae8] | 1808 |     OutputIonCoordinates(P, 1);
 | 
|---|
| [27a5bf] | 1809 | 
 | 
|---|
 | 1810 |     UpdateWaveAfterIonMove(P);
 | 
|---|
 | 1811 |     for (i=MAXOLD-1; i > 0; i--)    // store away old energies
 | 
|---|
 | 1812 |       E->TotalEnergyOuter[i] = E->TotalEnergyOuter[i-1];  
 | 
|---|
 | 1813 |     UpdateToNewWaves(P);
 | 
|---|
 | 1814 |     E->TotalEnergyOuter[0] = E->TotalEnergy[0]; 
 | 
|---|
 | 1815 |     //if ((P->R.ScaleTempStep > 0) && ((R->OuterStep-1) % P->R.ScaleTempStep == 0))
 | 
|---|
 | 1816 |     //  ScaleTemp(P);
 | 
|---|
| [a0bcf1] | 1817 |     if ((R->OuterStep-1) % P->R.OutSrcStep == 0) 
 | 
|---|
 | 1818 |       OutputVisSrcFiles(P, Occupied);
 | 
|---|
 | 1819 |     if ((R->OuterStep-1) % P->R.OutVisStep == 0) {      
 | 
|---|
| [27a5bf] | 1820 |       OutputVis(P, P->R.Lev0->Dens->DensityArray[TotalDensity]);
 | 
|---|
 | 1821 |       OutputIonForce(P);
 | 
|---|
 | 1822 |       EnergyOutput(P, 1);
 | 
|---|
| [a0bcf1] | 1823 |     }
 | 
|---|
| [27a5bf] | 1824 |     ResetForces(P);
 | 
|---|
 | 1825 |   }
 | 
|---|
| [a0bcf1] | 1826 |   SpeedMeasure(P, SimTime, StopTimeDo);
 | 
|---|
 | 1827 |   CloseOutputFiles(P);
 | 
|---|
 | 1828 | }
 | 
|---|