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