source: pcp/src/energy.c@ 99bed3

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

Merge branch 'ConcaveHull' of ../espack2 into ConcaveHull

Conflicts:

molecuilder/src/boundary.cpp
molecuilder/src/boundary.hpp
molecuilder/src/builder.cpp
molecuilder/src/linkedcell.cpp
molecuilder/src/linkedcell.hpp
molecuilder/src/vector.cpp
molecuilder/src/vector.hpp
util/src/NanoCreator.c

Basically, this resulted from a lot of conversions two from spaces to one tab, which is my standard indentation. The mess was caused by eclipse auto-indenting. And in espack2:ConcaveHull was the new stuff, so all from ConcaveHull was replaced in case of doubt.
Additionally, vector had ofstream << operator instead ostream << ...

  • Property mode set to 100755
File size: 36.0 KB
Line 
1/** \file energy.c
2 * Energy calculation, especially kinetic and electronic.
3 *
4 * Contains routines for output of calculated values EnergyOutput(),
5 * initialization of the Energy structure InitEnergyArray() or inital calculations
6 * InitPsiEnergyCalculation(),
7 * updating UpdateEnergyArray() and UpdatePsiEnergyCalculation() and of course
8 * calculation of the various parts: Density CalculateDensityEnergy(), ionic
9 * CalculateIonsEnergy(), electronic CalculatePsiEnergy(), kinetic with RiemannTensor
10 * CalculateKineticEnergyUseRT() (calling EnergyGNSP()) and without
11 * CalculateKineticEnergyNoRT() (using EnergyLapNSP()). CalculateGapEnergy() evaluates
12 * would-be addition from unoccupied states to total energy.
13 * Smaller functions such as EnergyAllReduce() gather intermediate results from all
14 * processes and sum up to total value.
15 * \note Some others are stored in excor.c (exchange energies) and pseudo.c (density
16 * energies)
17 *
18 *
19 Project: ParallelCarParrinello
20 \author Jan Hamaekers
21 \date 2000
22
23 File: energy.c
24 $Id: energy.c,v 1.68 2007/02/09 09:13:48 foo Exp $
25*/
26
27#ifdef HAVE_CONFIG_H
28#include <config.h>
29#endif
30
31#include <stdlib.h>
32#include <stdio.h>
33#include <stddef.h>
34#include <math.h>
35#include <string.h>
36
37#include "data.h"
38#include "errors.h"
39#include "helpers.h"
40#include "energy.h"
41#include "riemann.h"
42#include "excor.h"
43#include "myfft.h"
44#include "mymath.h"
45//#include "output.h"
46#include "run.h"
47
48/** Calculates kinetic energy for a given wave function in reciprocal base.
49 * \param *P Problem at hand
50 * \param *Lev LatticeLevel structure
51 * \param *LPsiDatA complex wave function coefficients \f$c_{i,G}\f$ (reciprocal base)
52 * \param *T kinetic eigenvalue, herein result is additionally stored
53 * \param Factor Additional factor for return value (such as occupation number \f$f_i\f$ )
54 * \return \f$E_{kin} = \frac{1}{2} f_i \sum_G |G|^2 |c_{i,G}|^2\f$
55 * \sa used in CalculateKineticEnergyNoRT()
56 */
57static double EnergyGNSP(const struct Problem *P, const struct LatticeLevel *Lev, fftw_complex *LPsiDatA, double *T, const double Factor) {
58 double LocalSP=0.0,PsiSP;
59 int i,s = 0;
60 struct OneGData *G = Lev->GArray;
61 if (Lev->GArray[0].GSq == 0.0) { // only for continuity (see other functions)
62 LocalSP += G[0].GSq*LPsiDatA[0].re*LPsiDatA[0].re;
63 s++;
64 }
65 /// Performs complex product for each reciprocal grid vector times twice this vector's norm.
66 for (i=s; i < Lev->MaxG; i++) {
67 LocalSP += G[i].GSq*2.*(LPsiDatA[i].re*LPsiDatA[i].re+LPsiDatA[i].im*LPsiDatA[i].im); // factor 2 due to gamma point symmetry
68 }
69 /// Gathers partial results by summation from all other coefficients sharing processes.
70 MPI_Allreduce ( &LocalSP, &PsiSP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
71 /// Multiplies by \a Factor and returns result.
72 *T = Factor*PsiSP;
73 return(Factor*PsiSP);
74}
75
76/** Calculates kinetic energy for a given wave function in reciprocal base in case of Riemann tensor use.
77 * \param *P Problem at hand
78 * \param *RT RiemannTensor structure
79 * \param *src replaces reciprocal grid vectors
80 * \param *Lap replaces complex Psi wave function coefficients
81 * \param *T kinetic eigenvalue, herein result is additionally stored
82 * \param Factor Additional factor for return value
83 * \sa used in CalculateKineticEnergyUseRT(), compare with EnergyGNSP()
84 */
85static double EnergyLapNSP(const struct Problem *P, const struct RiemannTensor *RT, fftw_complex *src, fftw_complex *Lap, double *T, const double Factor) {
86 double LocalSP=0.0,PsiSP;
87 int i;
88 int s = 0;
89 struct LatticeLevel *Lev = RT->LevS;
90 if (Lev->GArray[0].GSq == 0.0) {
91 LocalSP += src[0].re*Lap[0].re+src[0].im*Lap[0].im;
92 s++;
93 }
94 for (i=s; i < Lev->MaxG; i++) {
95 LocalSP += 2.*(src[i].re*Lap[i].re+src[i].im*Lap[i].im);
96 }
97 MPI_Allreduce ( &LocalSP, &PsiSP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
98 *T = Factor*PsiSP;
99 return(Factor*PsiSP);
100}
101
102/** Calculcates Kinetic energy without Riemann tensor.
103 * \param *P Problem at hand
104 * \param p local psi wave function number
105 * \sa EnergyGNSP()
106 */
107static void CalculateKineticEnergyNoRT(struct Problem *P, const int p)
108{
109 struct Lattice *Lat = &P->Lat;
110 struct Energy *E = Lat->E;
111 struct RunStruct *R = &P->R;
112 struct LatticeLevel *LevS = R->LevS;
113 struct LevelPsi *LPsi = LevS->LPsi;
114 struct Psis *Psi = &Lat->Psi;
115 const double PsiFactor = Psi->LocalPsiStatus[p].PsiFactor;
116 fftw_complex *psi;
117 psi = LPsi->LocalPsi[p];
118 E->PsiEnergy[KineticEnergy][p] = EnergyGNSP(P, LevS, psi, &Psi->AddData[p].T, 0.5*PsiFactor);
119}
120
121/** Calculcates Kinetic energy with Riemann tensor.
122 * CalculateRiemannLaplaceS() is called and Lattice->RT.RTLaplaceS used instead
123 * of Lattice::Psi.
124 * \param *P Problem at hand
125 * \param p local psi wave function number
126 * \sa EnergyLapNSP()
127 */
128static void CalculateKineticEnergyUseRT(struct Problem *P, const int p)
129{
130 struct Lattice *Lat = &P->Lat;
131 struct Energy *E = Lat->E;
132 struct RunStruct *R = &P->R;
133 struct LatticeLevel *LevS = R->LevS;
134 struct LevelPsi *LPsi = LevS->LPsi;
135 struct Psis *Psi = &Lat->Psi;
136 const double PsiFactor = Psi->LocalPsiStatus[p].PsiFactor;
137 fftw_complex *psi;
138 psi = LPsi->LocalPsi[p];
139 CalculateRiemannLaplaceS(Lat, &Lat->RT, psi, Lat->RT.RTLaplaceS);
140 E->PsiEnergy[KineticEnergy][p] = EnergyLapNSP(P, &Lat->RT, psi, Lat->RT.RTLaplaceS, &Psi->AddData[p].T, 0.5*PsiFactor);
141 //fprintf(stderr,"CalculateKineticEnergyUseRT: PsiEnergy[KineticEnergy] = %lg\n",E->PsiEnergy[KineticEnergy][p]);
142}
143
144/** Sums up calculated energies from all processes, calculating Energy::TotalEnergy.
145 * Via MPI_Allreduce is Energy::AllLocalPsiEnergy summed up and the sum
146 * redistributed back to all processes for the various spin cases. In
147 * case of SpinType#SpinUp or SpinType#SpinDown the "other" energy (e.g. down
148 * for an up process) is via MPI_SendRecv communicated, and lateron
149 * Energy#AllTotalPsiEnergy as the sum of two states generated in each process.<br>
150 * Finally, the Energy::TotalEnergy[0] is the sum of Energy::AllTotalPsiEnergy[],
151 * Energy::AllTotalDensityEnergy[], - Energy::AllTotalIonsEnergy[GaussSelfEnergy]
152 * and Energy::AllTotalIonsEnergy[EwaldEnergy].
153 * Again, as in UpdateEnergyArray(), this is split up completely for PsiTypeTag#UnOccupied
154 * and PsiTypeTag#Occupied.
155 * \param *P Problem at hand
156 */
157void EnergyAllReduce(struct Problem *P)
158{
159 struct Psis *Psi = &P->Lat.Psi;
160 struct RunStruct *R = &P->R;
161 struct Lattice *Lat = &P->Lat;
162 struct Energy *E = P->Lat.E;
163 struct OnePsiElement *OnePsiB;
164 MPI_Status status;
165 int i,j,m;
166 double lambda, overlap, lambda2, overlap2; //, part1, part2, tmp;
167 MPI_Allreduce(E->AllLocalDensityEnergy, E->AllTotalDensityEnergy, MAXALLDENSITYENERGY, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi );
168 switch (Psi->PsiST) {
169 case SpinDouble:
170 MPI_Allreduce(E->AllLocalPsiEnergy, E->AllTotalPsiEnergy, MAXALLPSIENERGY, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_PsiT);
171 break;
172 case SpinUp:
173 MPI_Allreduce(E->AllLocalPsiEnergy, E->AllUpPsiEnergy, MAXALLPSIENERGY, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_PsiT);
174 MPI_Sendrecv(E->AllUpPsiEnergy, MAXALLPSIENERGY, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag1,
175 E->AllDownPsiEnergy, MAXALLPSIENERGY, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag2,
176 P->Par.comm_STInter, &status );
177 for (i=0; i< MAXALLPSIENERGY; i++)
178 E->AllTotalPsiEnergy[i] = E->AllUpPsiEnergy[i] + E->AllDownPsiEnergy[i];
179 break;
180 case SpinDown:
181 MPI_Allreduce(E->AllLocalPsiEnergy, E->AllDownPsiEnergy, MAXALLPSIENERGY, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_PsiT);
182 MPI_Sendrecv(E->AllDownPsiEnergy, MAXALLPSIENERGY, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag2,
183 E->AllUpPsiEnergy, MAXALLPSIENERGY, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag1,
184 P->Par.comm_STInter, &status );
185 for (i=0; i< MAXALLPSIENERGY; i++)
186 E->AllTotalPsiEnergy[i] = E->AllUpPsiEnergy[i] + E->AllDownPsiEnergy[i];
187 break;
188 }
189 switch (R->CurrentMin) {
190 case Occupied:
191 E->TotalEnergy[0] = E->AllTotalPsiEnergy[KineticEnergy];
192 E->TotalEnergy[0] += E->AllTotalPsiEnergy[NonLocalEnergy]
193 ;
194 E->TotalEnergy[0] += E->AllTotalDensityEnergy[CorrelationEnergy]+E->AllTotalDensityEnergy[ExchangeEnergy];
195 E->TotalEnergy[0] += E->AllTotalDensityEnergy[PseudoEnergy];
196 E->TotalEnergy[0] += E->AllTotalDensityEnergy[HartreePotentialEnergy];
197 E->TotalEnergy[0] -= E->AllTotalIonsEnergy[GaussSelfEnergy];
198 E->TotalEnergy[0] += E->AllTotalIonsEnergy[EwaldEnergy];
199// fprintf(stderr,"(%i) Tot %1.5f ---\t Diff %1.5f\t ---\t Kin %1.5f\tNLErg %1.5f ---\t Corr %1.5f\tXC %1.5f\t Pseudo %1.5f\tHP %1.5f\n", P->Par.me,
200// E->TotalEnergy[0],(E->TotalEnergy[1]-E->TotalEnergy[0])/fabs(E->TotalEnergy[0]),
201// E->AllTotalPsiEnergy[KineticEnergy],E->AllTotalPsiEnergy[NonLocalEnergy],
202// E->AllTotalDensityEnergy[CorrelationEnergy],E->AllTotalDensityEnergy[ExchangeEnergy],
203// E->AllTotalDensityEnergy[PseudoEnergy],E->AllTotalDensityEnergy[HartreePotentialEnergy]);
204 break;
205 case UnOccupied:
206 E->TotalEnergy[0] = E->AllTotalPsiEnergy[KineticEnergy]+E->AllTotalPsiEnergy[NonLocalEnergy];
207 E->TotalEnergy[0] += E->AllTotalDensityEnergy[CorrelationEnergy]+E->AllTotalDensityEnergy[ExchangeEnergy];
208 E->TotalEnergy[0] += E->AllTotalDensityEnergy[PseudoEnergy]+E->AllTotalDensityEnergy[HartreePotentialEnergy];
209 break;
210 case Perturbed_P0:
211 case Perturbed_P1:
212 case Perturbed_P2:
213 case Perturbed_RxP0:
214 case Perturbed_RxP1:
215 case Perturbed_RxP2:
216 E->TotalEnergy[0] = E->AllTotalPsiEnergy[Perturbed1_0Energy]+E->AllTotalPsiEnergy[Perturbed0_1Energy];
217 E->parts[0] = E->AllTotalPsiEnergy[Perturbed1_0Energy]+E->AllTotalPsiEnergy[Perturbed0_1Energy];
218 //fprintf(stderr,"(%i) summed single: AllTotalPsiEnergy[total] = %lg + %lg = %lg\n",P->Par.me,E->AllTotalPsiEnergy[Perturbed1_0Energy],E->AllTotalPsiEnergy[Perturbed0_1Energy], E->TotalEnergy[0]);
219 //E->TotalEnergy[0] = 0.;
220 m = -1;
221 E->parts[1] = E->parts[2] = 0.;
222 for (i=0; i < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; i++) { // go through all wave functions: \sum_k
223 OnePsiB = &Psi->AllPsiStatus[i]; // grab OnePsiB
224 if (OnePsiB->PsiType == R->CurrentMin) { // drop all but occupied types
225 m++; // increase m if it is occupied wave function
226 lambda = overlap = 0.;
227 if (OnePsiB->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) { // local?
228 //fprintf(stderr,"(%i) Lambda[%i] = %lg\n",P->Par.me, m,Lat->Psi.AddData[i].Lambda);
229 lambda = Lat->Psi.AddData[OnePsiB->MyLocalNo].Lambda * Lat->Psi.LocalPsiStatus[OnePsiB->MyLocalNo].PsiFactor;
230 for (j=0;j<Psi->NoOfPsis;j++) { // over all Psis: \sum_l
231 //lambda -= Lat->Psi.lambda[j][i - Psi->TypeStartIndex[R->CurrentMin]]*GradSP(P,R->LevS,R->LevS->LPsi->LocalPsi[j + Psi->TypeStartIndex[R->CurrentMin]],R->LevS->LPsi->LocalPsi[i]);
232 overlap -= Lat->Psi.lambda[j][m]*Lat->Psi.Overlap[j][m];
233 }
234 // send results to other processes
235 }
236 //fprintf(stderr,"(%i) before MPI: \tlambda[%i] = %lg\t overlap[%i] = %lg\n",P->Par.me,m,lambda,m,overlap);
237 MPI_Allreduce( &lambda, &lambda2, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_PsiT); // lambda is just local
238 MPI_Allreduce( &overlap, &overlap2, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_PsiT); // we just summed up over local m's!
239 //fprintf(stderr,"(%i) after MPI: \t lambda[%i] = %lg\t overlap[%i] = %lg\n",P->Par.me, OnePsiB->MyGlobalNo,lambda2, m, overlap2);
240 E->TotalEnergy[0] += lambda2;
241 E->TotalEnergy[0] += overlap2;
242 E->parts[1] += lambda2;
243 E->parts[2] += overlap2;
244 }
245 }
246
247 switch(Psi->PsiST) {
248 default:
249 lambda = overlap = 0.;
250 break;
251 case SpinUp:
252 MPI_Sendrecv(&E->parts[1], 1, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag1,
253 &lambda, 1, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag2, P->Par.comm_STInter, &status );
254 E->TotalEnergy[0] += lambda;
255 MPI_Sendrecv(&E->parts[2],1, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag3,
256 &overlap, 1, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag4, P->Par.comm_STInter, &status );
257 E->TotalEnergy[0] += overlap;
258 break;
259 case SpinDown:
260 MPI_Sendrecv(&E->parts[1], 1, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag2,
261 &lambda, 1, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag1, P->Par.comm_STInter, &status );
262 E->TotalEnergy[0] += lambda;
263 MPI_Sendrecv(&E->parts[2],1, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag4,
264 &overlap, 1, MPI_DOUBLE, P->Par.me_comm_ST, EnergyTag3, P->Par.comm_STInter, &status );
265 E->TotalEnergy[0] += overlap;
266 break;
267 }
268 //fprintf(stderr,"(%i) summed single: lambda[total] = %lg\t overlap[total] = %lg\n",P->Par.me,E->parts[1]+lambda, E->parts[2]+overlap);
269 break;
270 default:
271 break;
272 }
273}
274
275/** Initializes Energy Array.
276 * Arrays are malloc'd and set to zero, then InitExchangeCorrelationEnergy()
277 * called.
278 * \param *P Problem at hand
279 */
280void InitEnergyArray(struct Problem *P)
281{
282 struct Lattice *Lat = &P->Lat;
283 struct Energy *E = Lat->E;
284 struct Psis *Psi = &Lat->Psi;
285 int i, type, oldtype = P->R.CurrentMin;
286 for (type=Occupied;type<Extra;type++) {
287 SetCurrentMinState(P, type);
288 for (i=0; i<MAXALLPSIENERGY; i++) {
289 Lat->Energy[type].PsiEnergy[i] = (double *)
290 Malloc(sizeof(double)*Psi->LocalNo,"InitEnergyCalculations: PsiKineticEnergy");
291 SetArrayToDouble0(Lat->Energy[type].PsiEnergy[i], Psi->LocalNo);
292 }
293 SetArrayToDouble0(E->AllLocalDensityEnergy, MAXALLDENSITYENERGY);
294 SetArrayToDouble0(E->AllTotalDensityEnergy, MAXALLDENSITYENERGY);
295 SetArrayToDouble0(E->AllTotalIonsEnergy, MAXALLIONSENERGY);
296 SetArrayToDouble0(E->TotalEnergy, MAXOLD);
297 }
298 InitExchangeCorrelationEnergy(P, &P->ExCo);
299 SetCurrentMinState(P, oldtype);
300}
301
302/** Shifts stored total energy values up to MAXOLD.
303 * Energy::PsiEnergy, Energy::AllLocalDensityEnergy, Energy::AllTotalDensityEnergy
304 * are all set to zero. Energy::TotalEnergy is shifted by one to make
305 * place for next value. This is done either for occupied or for the gap energies,
306 * stated by \a IncType.
307 * \param *P Problem at hand
308 */
309void UpdateEnergyArray(struct Problem *P)
310{
311 struct Lattice *Lat = &P->Lat;
312 struct RunStruct *R = &P->R;
313 struct Energy *E = Lat->E;
314 int i;
315 for (i=0; i<MAXALLPSIENERGY; i++)
316 E->PsiEnergy[i][R->OldActualLocalPsiNo] = 0.0;
317 SetArrayToDouble0(E->AllLocalDensityEnergy, MAXALLDENSITYENERGY);
318 SetArrayToDouble0(E->AllTotalDensityEnergy, MAXALLDENSITYENERGY);
319 for (i=MAXOLD-1; i > 0; i--)
320 Lat->Energy[R->CurrentMin].TotalEnergy[i] = Lat->Energy[R->CurrentMin].TotalEnergy[i-1];
321}
322
323/** Calculates ion energy.
324 * CalculateGaussSelfEnergyNoRT() is called.
325 * \param *P Problem at hand
326 * \warning CalculateIonsUseRT not implemented
327 */
328void CalculateIonsEnergy(struct Problem *P)
329{
330 struct Lattice *Lat = &P->Lat;
331 switch (Lat->RT.ActualUse) {
332 case inactive:
333 case standby:
334 CalculateGaussSelfEnergyNoRT(P, &P->PP, &P->Ion);
335 break;
336 case active:
337 Error(SomeError, "CalculateIonsUseRT: Not implemented");
338 break;
339 }
340}
341
342/** Calculates density energy.
343 * Depending on RT::ActualUse CalculateXCEnergyNoRT(), CalculateSomeEnergyAndHGNoRT()
344 * and CalculateGapEnergy() are called within SpeedMeasure() brackets.
345 * Afterwards, CalculateGradientNoRT() for RunStruct::ActualLocalPsiNo.
346 * \param *P Problem at hand
347 * \param UseInitTime whether SpeedMeasure uses InitLocTime (Init) or LocTime (Update)
348 * \warning CalculateSomeEnergyAndHGUseRT not implemented
349 * \warning call EnergyAllReduce() afterwards!
350 */
351void CalculateDensityEnergy(struct Problem *P, const int UseInitTime)
352{
353 struct Lattice *Lat = &P->Lat;
354 struct RunStruct *R = &P->R;
355 struct LatticeLevel *LevS = R->LevS;
356 switch (Lat->RT.ActualUse) {
357 case inactive:
358 case standby:
359 SpeedMeasure(P, (UseInitTime == 1 ? InitLocTime : LocTime), StartTimeDo);
360 CalculateXCEnergyNoRT(P);
361 CalculateSomeEnergyAndHGNoRT(P);
362 SpeedMeasure(P, (UseInitTime == 1 ? InitLocTime : LocTime), StopTimeDo);
363 CalculateGradientNoRT(P,
364 LevS->LPsi->LocalPsi[R->ActualLocalPsiNo],
365 Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].PsiFactor,
366 &Lat->Psi.AddData[R->ActualLocalPsiNo].Lambda,
367 P->PP.fnl[R->ActualLocalPsiNo],
368 P->Grad.GradientArray[ActualGradient],
369 P->Grad.GradientArray[OldActualGradient],
370 P->Grad.GradientArray[HcGradient]);
371 break;
372 case active:
373 CalculateXCEnergyUseRT(P);
374 /* FIXME - Hier fehlt noch was */
375 Error(SomeError, "CalculateSomeEnergyAndHGUseRT: Not implemented");
376 break;
377 }
378 /* Achtung danach noch EnergyAllReduce aufrufen */
379}
380
381/** Calculation of an additional total energy from unoccupied states.
382 * By including additional (unoccupied) orbitals into the minimisation process,
383 * which also stand orthogonal on all others, however are minimised in the field
384 * of occupied ones (without altering it in a self-consistent manner), a first
385 * approximation to the band gap energy can be obtained. The density part of this
386 * virtual energy is stored in the PsiTypeTag#UnOccupied part of Lattice#Energy.
387 * Note that the band gap is approximated via the eigenvalues of the additional
388 * orbitals. This energy here is purely fictious.
389 * \param *P Problem at hand
390 * \note No RiemannTensor use so far implemented. And mostly copy&paste from
391 * CalculateXCEnergyNoRT() and CalculateSomeEnergyAndHGNoRT().
392 * \bug Understand CoreCorrection part in \f$E_{XC}\f$
393 */
394void CalculateGapEnergy(struct Problem *P) {
395 struct Lattice *Lat = &P->Lat;
396 struct Psis *Psi = &Lat->Psi;
397 struct Energy *E = Lat->E;
398 struct PseudoPot *PP = &P->PP;
399 struct RunStruct *R = &P->R;
400 struct LatticeLevel *Lev0 = R->Lev0;
401 struct LatticeLevel *LevS = R->LevS;
402 struct Density *Dens = Lev0->Dens;
403 struct ExCor *EC = &P->ExCo;
404 struct Ions *I = &P->Ion;
405 double SumEc = 0;
406 double rs, p = 0.0, pUp = 0.0, pDown = 0.0, zeta, rsUp, rsDown, SumEx=0.0;
407 double Factor = R->XCEnergyFactor/Lev0->MaxN;
408 fftw_complex vp,rhoe,rhoe_gap;
409 fftw_complex rp,rhog;
410 double SumEHP,Fac,SumEH,SumEG,SumEPS;
411 int g,s=0,it,Index;
412 int i;
413 fftw_complex *HG = Dens->DensityCArray[HGDensity];
414 SetArrayToDouble0((double *)HG,2*Dens->TotalSize);
415
416 //if (isnan(HG[0].re)) { fprintf(stderr,"WARNGING: CalculateGapEnergy(): HG[%i]_%i = NaN!\n", 0, R->ActualLocalPsiNo); Error(SomeError, "NaN-Fehler!"); }
417 SpeedMeasure(P, GapTime, StartTimeDo);
418 // === Calculate exchange and correlation energy ===
419 for (i = 0; i < Dens->LocalSizeR; i++) { // for each node in radial mesh
420 // put (corecorrected) densities in p, pUp, pDown
421 p = Dens->DensityArray[TotalDensity][i];
422 //if (isnan(p)) { fprintf(stderr,"WARNGING: CalculateGapEnergy(): Dens->DensityArray[TotalUpDensity][%i]_%i = NaN!\n", i, R->ActualLocalPsiNo); Error(SomeError, "NaN-Fehler!"); }
423 switch (Psi->PsiST) {
424 case SpinDouble:
425 pUp = 0.5*p;
426 pDown = 0.5*p;
427 break;
428 case SpinUp:
429 case SpinDown:
430 pUp = Dens->DensityArray[TotalUpDensity][i];
431 pDown = Dens->DensityArray[TotalDownDensity][i];
432 break;
433 default:
434 ;
435 }
436 if ((p < 0) || (pUp < 0) || (pDown < 0)) {
437 /*fprintf(stderr,"index %i pc %g p %g\n",i,Dens->DensityArray[CoreWaveDensity][i],p);*/
438 p = 0.0;
439 pUp = 0.0;
440 pDown = 0.0;
441 }
442 rs = Calcrs(EC,p);
443 rsUp = Calcrs(EC,pUp);
444 rsDown = Calcrs(EC,pDown);
445 zeta = CalcZeta(EC,pUp,pDown);
446 // Calculation with the densities and summation
447 SumEx += CalcSEXr(EC, rsUp, rsDown, pUp, pDown);
448 SumEc += CalcSECr(EC, rs, zeta, p);
449 }
450 E->AllLocalDensityEnergy[CorrelationEnergy] = Factor*SumEc;
451 E->AllLocalDensityEnergy[ExchangeEnergy] = Factor*SumEx;
452
453 // === Calculate electrostatic and local pseudopotential energy ===
454 SumEHP = 0.0;
455 SumEH = 0.0;
456 SumEG = 0.0;
457 SumEPS = 0.0;
458 // calculates energy of local pseudopotential
459 if (Lev0->GArray[0].GSq == 0.0) {
460 Index = Lev0->GArray[0].Index;
461 c_re(vp) = 0.0;
462 c_im(vp) = 0.0;
463 for (it = 0; it < I->Max_Types; it++) {
464 c_re(vp) += (c_re(I->I[it].SFactor[0])*PP->phi_ps_loc[it][0]);
465 c_im(vp) += (c_im(I->I[it].SFactor[0])*PP->phi_ps_loc[it][0]);
466 }
467 //if (isnan(c_re(Dens->DensityCArray[GapDensity][Index]))) { fprintf(stderr,"WARNGING: CalculateGapEnergy(): Dens->DensityArray[GapDensity][%i]_%i.re = NaN!\n", Index, R->ActualLocalPsiNo); Error(SomeError, "NaN-Fehler!"); }
468 SumEPS += (c_re(Dens->DensityCArray[GapDensity][Index])*c_re(vp) +
469 c_im(Dens->DensityCArray[GapDensity][Index])*c_im(vp))*R->HGcFactor;
470 s++;
471 }
472 for (g=s; g < Lev0->MaxG; g++) {
473 Index = Lev0->GArray[g].Index;
474 c_re(vp) = 0.0;
475 c_im(vp) = 0.0;
476 c_re(rp) = 0.0;
477 c_im(rp) = 0.0;
478 Fac = 2.*PI/(Lev0->GArray[g].GSq); // 4*.. -> 2*...: 1/2 here ! (due to lacking dependence of first n^{tot} (r) on gap psis)
479 for (it = 0; it < I->Max_Types; it++) {
480 c_re(vp) += (c_re(I->I[it].SFactor[g])*PP->phi_ps_loc[it][g]);
481 c_im(vp) += (c_im(I->I[it].SFactor[g])*PP->phi_ps_loc[it][g]);
482 c_re(rp) += (c_re(I->I[it].SFactor[g])*PP->FacGauss[it][g]);
483 c_im(rp) += (c_im(I->I[it].SFactor[g])*PP->FacGauss[it][g]);
484 }
485 //if (isnan(c_re(Dens->DensityCArray[TotalDensity][Index]))) { fprintf(stderr,"WARNGING: CalculateGapEnergy(): Dens->DensityArray[TotalDensity][%i]_%i.re = NaN!\n", Index, R->ActualLocalPsiNo); Error(SomeError, "NaN-Fehler!"); }
486 c_re(rhog) = c_re(Dens->DensityCArray[TotalDensity][Index])*R->HGcFactor+c_re(rp);
487 c_im(rhog) = c_im(Dens->DensityCArray[TotalDensity][Index])*R->HGcFactor+c_im(rp);
488 c_re(rhoe) = c_re(Dens->DensityCArray[TotalDensity][Index])*R->HGcFactor;
489 c_im(rhoe) = c_im(Dens->DensityCArray[TotalDensity][Index])*R->HGcFactor;
490 //if (isnan(c_re(Dens->DensityCArray[GapDensity][Index]))) { fprintf(stderr,"WARNGING: CalculateGapEnergy(): Dens->DensityArray[GapDensity][%i]_%i.re = NaN!\n", Index, R->ActualLocalPsiNo); Error(SomeError, "NaN-Fehler!"); }
491 c_re(rhoe_gap) = c_re(Dens->DensityCArray[GapDensity][Index])*R->HGcFactor;
492 c_im(rhoe_gap) = c_im(Dens->DensityCArray[GapDensity][Index])*R->HGcFactor;
493 //c_re(PP->VCoulombc[g]) = Fac*c_re(rhog);
494 //c_im(PP->VCoulombc[g]) = -Fac*c_im(rhog);
495 c_re(HG[Index]) = c_re(vp) + Fac * (c_re(Dens->DensityCArray[TotalDensity][Index])*R->HGcFactor+c_re(rp));
496 c_im(HG[Index]) = c_im(vp) + Fac * (c_im(Dens->DensityCArray[TotalDensity][Index])*R->HGcFactor+c_im(rp));
497 //if (isnan(c_re(HG[Index]))) { fprintf(stderr,"WARNGING: CalculateGapEnergy(): HG[%i]_%i.re = NaN!\n", Index, R->ActualLocalPsiNo); Error(SomeError, "NaN-Fehler!"); }
498 //if (isnan(c_im(HG[Index]))) { fprintf(stderr,"WARNGING: CalculateGapEnergy(): HG[%i]_%i.im = NaN!\n", Index, R->ActualLocalPsiNo); Error(SomeError, "NaN-Fehler!"); }
499 /*if (P->first) */
500 //SumEG += Fac*(c_re(rp)*c_re(rp)+c_im(rp)*c_im(rp)); // Coulomb energy
501 // changed: took out gaussian part! rhog -> rhoe
502 SumEHP += 2*Fac*(c_re(rhoe)*c_re(rhoe_gap)+c_im(rhoe)*c_im(rhoe_gap)); // E_ES first part
503 //SumEH += Fac*(c_re(rhoe)*c_re(rhoe)+c_im(rhoe)*c_im(rhoe));
504 SumEPS += 2.*(c_re(Dens->DensityCArray[GapDensity][Index])*c_re(vp)+
505 c_im(Dens->DensityCArray[GapDensity][Index])*c_im(vp))*R->HGcFactor;
506 }
507 //
508 for (i=0; i<Lev0->MaxDoubleG; i++) {
509 HG[Lev0->DoubleG[2*i+1]].re = HG[Lev0->DoubleG[2*i]].re;
510 HG[Lev0->DoubleG[2*i+1]].im = -HG[Lev0->DoubleG[2*i]].im;
511 }
512 /*if (P->first) */
513 //E->AllLocalDensityEnergy[GaussEnergy] = SumEG*Lat->Volume;
514 E->AllLocalDensityEnergy[PseudoEnergy] = SumEPS*Lat->Volume; // local pseudopotential
515 E->AllLocalDensityEnergy[HartreePotentialEnergy] = SumEHP*Lat->Volume; // Gauss n^Gauss and Hartree n potential
516 //E->AllLocalDensityEnergy[HartreeEnergy] = SumEH*Lat->Volume;
517 SpeedMeasure(P, GapTime, StopTimeDo);
518
519 // === Calculate Gradient ===
520 //if (isnan(P->R.LevS->LPsi->LocalPsi[P->R.ActualLocalPsiNo][0].re)) { fprintf(stderr,"(%i) WARNING in CalculateGapEnergy(): ActualLocalPsi_%i[%i] = NaN!\n", P->Par.me, P->R.ActualLocalPsiNo, 0); Error(SomeError, "NaN-Fehler!"); }
521 CalculateGradientNoRT(P,
522 LevS->LPsi->LocalPsi[R->ActualLocalPsiNo],
523 Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].PsiFactor,
524 &Lat->Psi.AddData[R->ActualLocalPsiNo].Lambda,
525 P->PP.fnl[R->ActualLocalPsiNo],
526 P->Grad.GradientArray[ActualGradient],
527 P->Grad.GradientArray[OldActualGradient],
528 P->Grad.GradientArray[HcGradient]);
529 //if (isnan(P->Grad.GradientArray[ActualGradient][0].re)) { fprintf(stderr,"(%i) WARNING in CalculateGapEnergy(): ActualGradient_%i[%i] = NaN!\n", P->Par.me, P->R.ActualLocalPsiNo, 0); Error(SomeError, "NaN-Fehler!"); }
530}
531
532/** Calculates the energy of a wave function psi.
533 * Calculates kinetc and non local energies, taking care of possible
534 * Riemann tensor usage.
535 * \param *P Problem at hand
536 * \param p number of wave function Psi
537 * \param UseInitTime Whether speed measuring uses InitLocTime (Init) or LocTime (Update)
538 * \sa EnergyGNSP(), EnergyLapNSP()
539 * \warning call EnergyAllReduce() afterwards!
540 */
541void CalculatePsiEnergy(struct Problem *P, const int p, const int UseInitTime)
542{
543 struct Lattice *Lat = &(P->Lat);
544 switch (Lat->RT.ActualUse) {
545 case inactive:
546 case standby:
547 SpeedMeasure(P, (UseInitTime == 1 ? InitLocTime : LocTime), StartTimeDo);
548 CalculateKineticEnergyNoRT(P, p);
549 SpeedMeasure(P, (UseInitTime == 1 ? InitLocTime : LocTime), StopTimeDo);
550 SpeedMeasure(P, (UseInitTime == 1 ? InitNonLocTime : NonLocTime), StartTimeDo);
551 CalculateNonLocalEnergyNoRT(P, p);
552 SpeedMeasure(P, (UseInitTime == 1 ? InitNonLocTime : NonLocTime), StopTimeDo);
553 break;
554 case active:
555 SpeedMeasure(P, (UseInitTime == 1 ? InitLocTime : LocTime), StartTimeDo);
556 CalculateKineticEnergyUseRT(P, p);
557 SpeedMeasure(P, (UseInitTime == 1 ? InitLocTime : LocTime), StopTimeDo);
558 SpeedMeasure(P, (UseInitTime == 1 ? InitNonLocTime : NonLocTime), StartTimeDo);
559 CalculateNonLocalEnergyUseRT(P, p);
560 SpeedMeasure(P, (UseInitTime == 1 ? InitNonLocTime : NonLocTime), StopTimeDo);
561 break;
562 }
563 /* Achtung danach noch EnergyAllReduce aufrufen */
564}
565
566/** Initialization of wave function energy calculation.
567 * For every wave function CalculatePsiEnergy() is called,
568 * Energy::AllLocalPsiEnergy is set to sum of these.
569 * \param *P Problem at hand
570 * \param type minimisation type PsiTypeTag to calculate psi energy for
571 * \warning call EnergyAllReduce() afterwards!
572 */
573void InitPsiEnergyCalculation(struct Problem *P, enum PsiTypeTag type)
574{
575 struct Lattice *Lat = &(P->Lat);
576 int p,i, oldtype = P->R.CurrentMin;
577 SetCurrentMinState(P,type);
578 for (p=0; p < Lat->Psi.LocalNo; p++) {
579 if (P->Lat.Psi.LocalPsiStatus[p].PsiType == P->R.CurrentMin)
580 CalculatePsiEnergy(P, p, 1);
581 }
582 for (i=0; i< MAXALLPSIENERGY; i++) {
583 Lat->E->AllLocalPsiEnergy[i] = 0.0;
584 for (p=0; p < Lat->Psi.LocalNo; p++)
585 if (P->Lat.Psi.LocalPsiStatus[p].PsiType == P->R.CurrentMin)
586 Lat->E->AllLocalPsiEnergy[i] += Lat->E->PsiEnergy[i][p];
587 }
588 SetCurrentMinState(P,oldtype);
589 /* Achtung danach noch EnergyAllReduce aufrufen */
590}
591
592/** Updating wave function energy.
593 * CalculatePsiEnergy() is called for RunStruct::OldActualLocalPsiNo,
594 * Energy::AllLocalPsiEnergy is set to sum of all Psi energies.
595 * \param *P Problem at hand
596 * \warning call EnergyAllReduce() afterwards!
597 */
598void UpdatePsiEnergyCalculation(struct Problem *P)
599{
600 struct Lattice *Lat = &(P->Lat);
601 struct RunStruct *R = &P->R;
602 int p = R->OldActualLocalPsiNo, i;
603 if (P->Lat.Psi.LocalPsiStatus[p].PsiType == P->R.CurrentMin)
604 CalculatePsiEnergy(P, p, 0);
605 for (i=0; i< MAXALLPSIENERGY; i++) {
606 Lat->E->AllLocalPsiEnergy[i] = 0.0;
607 for (p=0; p < P->Lat.Psi.LocalNo; p++)
608 if (P->Lat.Psi.LocalPsiStatus[p].PsiType == P->R.CurrentMin)
609 Lat->E->AllLocalPsiEnergy[i] += Lat->E->PsiEnergy[i][p];
610 }
611 /* Achtung danach noch EnergyAllReduce aufrufen */
612}
613
614/** Output of energies.
615 * Prints either Time, PsiEnergy(Up/Down), DensityEnergy, IonsEnergy, KineticEnergy,
616 * TotalEnergy and Temperature to screen or<br>
617 * Time, TotalEnergy, NonLocalEnergy, CorrelationEnergy, ExchangeEnergy,
618 * PseudoEnergy, HartreePotentialEnergy, GaussSelfEnergy, EwaldEnergy,
619 * KineticEnergy and TotalEnergy+KineticEnergy to the energy file
620 * \param *P Problem at hand
621 * \param doout 1 - print stuff to stderr, 0 - put compact form into energy file
622 */
623void EnergyOutput(struct Problem *P, int doout)
624{
625 struct Lattice *Lat = &P->Lat;
626 struct RunStruct *R = &P->R;
627 struct LatticeLevel *Lev0 = R->Lev0;
628 struct LatticeLevel *LevS = R->LevS;
629 struct Ions *I = &P->Ion;
630 struct Psis *Psi = &Lat->Psi;
631 struct FileData *F = &P->Files;
632 FILE *eout = NULL;
633 FILE *convergence = NULL;
634 int i, it;
635 if (P->Call.out[LeaderOut] && (P->Par.me == 0) && doout) {
636 fprintf(stderr, "Time T: :\t%e\n",P->R.t);
637 fprintf(stderr, "PsiEnergy :");
638 for (i=0; i<MAXALLPSIENERGY; i++) {
639 fprintf(stderr," %e",Lat->Energy[Occupied].AllTotalPsiEnergy[i]);
640 }
641 fprintf(stderr,"\n");
642 if (Psi->PsiST != SpinDouble) {
643 fprintf(stderr, "PsiEnergyUp :");
644 for (i=0; i<MAXALLPSIENERGY; i++) {
645 fprintf(stderr,"\t%e",Lat->Energy[Occupied].AllUpPsiEnergy[i]);
646 }
647 fprintf(stderr,"\n");
648 fprintf(stderr, "PsiEnergyDown:");
649 for (i=0; i<MAXALLPSIENERGY; i++) {
650 fprintf(stderr,"\t%e",Lat->Energy[Occupied].AllDownPsiEnergy[i]);
651 }
652 fprintf(stderr,"\n");
653 }
654 fprintf(stderr, "DensityEnergy:");
655 for (i=0; i<MAXALLDENSITYENERGY; i++) {
656 fprintf(stderr,"\t%e",Lat->Energy[Occupied].AllTotalDensityEnergy[i]);
657 }
658 fprintf(stderr,"\n");
659 fprintf(stderr, "IonsEnergy :");
660 for (i=0; i<MAXALLIONSENERGY; i++) {
661 fprintf(stderr,"\t%e",Lat->Energy[Occupied].AllTotalIonsEnergy[i]);
662 }
663 fprintf(stderr,"\n");
664 fprintf(stderr, "TotalEnergy :\t%e\n",Lat->Energy[Occupied].TotalEnergy[0]);
665 if (R->DoPerturbation)
666 fprintf(stderr, "PerturbedEnergy :\t%e\t%e\t%e\t%e\t%e\t%e\n",Lat->Energy[Perturbed_P0].TotalEnergy[0], Lat->Energy[Perturbed_P1].TotalEnergy[0], Lat->Energy[Perturbed_P2].TotalEnergy[0], Lat->Energy[Perturbed_RxP0].TotalEnergy[0], Lat->Energy[Perturbed_RxP1].TotalEnergy[0], Lat->Energy[Perturbed_RxP2].TotalEnergy[0]);
667 if (R->DoUnOccupied) {
668 fprintf(stderr, "GapEnergy :\t%e\t%e\t%e\t%e\n",Lat->Energy[UnOccupied].AllTotalPsiEnergy[KineticEnergy],Lat->Energy[UnOccupied].AllTotalPsiEnergy[NonLocalEnergy], Lat->Energy[UnOccupied].AllTotalDensityEnergy[PseudoEnergy],Lat->Energy[UnOccupied].AllTotalDensityEnergy[HartreePotentialEnergy]);
669 fprintf(stderr, "TotalGapEnergy :\t%e\n",Lat->Energy[UnOccupied].TotalEnergy[0]);
670 fprintf(stderr, "BandGap :\t%e\n", Lat->Energy[UnOccupied].bandgap);
671 }
672 fprintf(stderr, "IonKinEnergy :\t%e\n",P->Ion.EKin);
673 fprintf(stderr, "Temperature :\t%e\n",P->Ion.ActualTemp);
674 }
675 if (P->Files.MeOutMes == 0 || doout == 0) return;
676 eout = F->EnergyFile;
677 if (eout != NULL) {
678 if (R->DoUnOccupied) {
679 fprintf(eout, "%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\n",
680 P->R.t,
681 Lat->Energy[Occupied].TotalEnergy[0], Lat->Energy[UnOccupied].TotalEnergy[0],
682 Lat->Energy[Occupied].AllTotalPsiEnergy[KineticEnergy], Lat->Energy[Occupied].AllTotalPsiEnergy[NonLocalEnergy],
683 Lat->Energy[UnOccupied].AllTotalPsiEnergy[KineticEnergy]+Lat->Energy[UnOccupied].AllTotalPsiEnergy[NonLocalEnergy],
684 Lat->Energy[Occupied].AllTotalDensityEnergy[CorrelationEnergy], Lat->Energy[Occupied].AllTotalDensityEnergy[ExchangeEnergy],
685 Lat->Energy[Occupied].AllTotalDensityEnergy[PseudoEnergy], Lat->Energy[Occupied].AllTotalDensityEnergy[HartreePotentialEnergy],
686 Lat->Energy[UnOccupied].AllTotalDensityEnergy[PseudoEnergy]+Lat->Energy[UnOccupied].AllTotalDensityEnergy[HartreePotentialEnergy],
687 -Lat->Energy[Occupied].AllTotalIonsEnergy[GaussSelfEnergy],
688 Lat->Energy[Occupied].AllTotalIonsEnergy[EwaldEnergy],
689 P->Ion.EKin,
690 Lat->Energy[Occupied].TotalEnergy[0]+P->Ion.EKin);
691 } else {
692 fprintf(eout, "%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\n",
693 P->R.t,
694 Lat->Energy[Occupied].TotalEnergy[0],
695 Lat->Energy[Occupied].AllTotalPsiEnergy[KineticEnergy], Lat->Energy[Occupied].AllTotalPsiEnergy[NonLocalEnergy],
696 Lat->Energy[Occupied].AllTotalDensityEnergy[CorrelationEnergy], Lat->Energy[Occupied].AllTotalDensityEnergy[ExchangeEnergy],
697 Lat->Energy[Occupied].AllTotalDensityEnergy[PseudoEnergy], Lat->Energy[Occupied].AllTotalDensityEnergy[HartreePotentialEnergy],
698 -Lat->Energy[Occupied].AllTotalIonsEnergy[GaussSelfEnergy],
699 Lat->Energy[Occupied].AllTotalIonsEnergy[EwaldEnergy],
700 P->Ion.EKin,
701 Lat->Energy[Occupied].TotalEnergy[0]+P->Ion.EKin);
702 }
703 fflush(eout);
704 }
705
706 // Output for testing convergence
707 if (!P->Par.me) {
708 if (R->DoPerturbation) {
709 if ((convergence = fopen("pcp.convergence.csv","r")) != NULL) { // only prepend with header if file doesn't exist yet
710 fclose(convergence);
711 convergence = fopen("pcp.convergence.csv","a");
712 } else {
713 convergence = fopen("pcp.convergence.csv","w");
714 if (convergence != NULL) {
715 fprintf(convergence,"Ecut\ta.u.\tMaxG\tN0\tN1\tN2\tChi00\tChi11\tChi22\tP0\tP1\tP2\tRxP0\tRxP1\tRxP2\tTE\t");
716 for (it=0; it < I->Max_Types; it++) { // over all types
717 fprintf(convergence,"Sigma00_%s\tSigma11_%s\tSigma22_%s\t",I->I[it].Symbol,I->I[it].Symbol,I->I[it].Symbol);
718 fprintf(convergence,"Sigma00_%s_rezi\tSigma11_%s_rezi\tSigma22_%s_rezi\t",I->I[it].Symbol,I->I[it].Symbol,I->I[it].Symbol);
719 }
720 fprintf(convergence,"Volume\tN^3\tSawStart\t");
721 fprintf(convergence,"Perturbedrxp-0\tPerturbedrxp-1\tPerturbedrxp-2\t");
722 fprintf(convergence,"Perturbedp-0\tPerturbedp-1\tPerturbedp-2\n");
723 }
724 }
725 if (convergence != NULL) {
726 fprintf(convergence,"%e\t%e\t%i\t%i\t%i\t%i\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t",
727 Lat->ECut,
728 sqrt(Lat->RealBasisSQ[0]), LevS->TotalAllMaxG, Lev0->Plan0.plan->N[0], Lev0->Plan0.plan->N[1], Lev0->Plan0.plan->N[2],
729 I->I[0].chi[0+0*NDIM], I->I[0].chi[1+1*NDIM], I->I[0].chi[2+2*NDIM],
730 Lat->Energy[Perturbed_P0].TotalEnergy[0], Lat->Energy[Perturbed_P1].TotalEnergy[0], Lat->Energy[Perturbed_P2].TotalEnergy[0],
731 Lat->Energy[Perturbed_RxP0].TotalEnergy[0], Lat->Energy[Perturbed_RxP1].TotalEnergy[0],Lat->Energy[Perturbed_RxP2].TotalEnergy[0],
732 Lat->Energy[Occupied].TotalEnergy[0]);
733 for (it=0; it < I->Max_Types; it++) { // over all types
734 //fprintf(convergence,"%e\t%e\t%e\t%e\t%e\t%e\t",
735 fprintf(convergence,"%e\t%e\t%e\t",
736 //I->I[it].sigma[0][0+0*NDIM], I->I[it].sigma[0][1+1*NDIM], I->I[it].sigma[0][2+2*NDIM],
737 I->I[it].sigma_rezi[0][0+0*NDIM], I->I[it].sigma_rezi[0][1+1*NDIM], I->I[it].sigma_rezi[0][2+2*NDIM]);
738 }
739 fprintf(convergence,"%e\t%i\t%e\t",Lat->Volume, Lev0->Plan0.plan->N[0]*Lev0->Plan0.plan->N[1]*Lev0->Plan0.plan->N[2],
740 P->Lat.SawtoothStart);
741 fprintf(convergence,"%e\t%e\t%e\t",Lat->E[Perturbed_RxP0].AllTotalPsiEnergy[Perturbed1_0Energy],Lat->E[Perturbed_RxP1].AllTotalPsiEnergy[Perturbed1_0Energy],Lat->E[Perturbed_RxP2].AllTotalPsiEnergy[Perturbed1_0Energy]);
742 fprintf(convergence,"%e\t%e\t%e\n",Lat->E[Perturbed_P0].AllTotalPsiEnergy[Perturbed1_0Energy],Lat->E[Perturbed_P1].AllTotalPsiEnergy[Perturbed1_0Energy],Lat->E[Perturbed_P2].AllTotalPsiEnergy[Perturbed1_0Energy]);
743 fclose(convergence);
744 }
745 } else {
746 if ((convergence = fopen("pcp.convergence.csv","r")) != NULL) { // only prepend with header if file doesn't exist yet
747 fclose(convergence);
748 convergence = fopen("pcp.convergence.csv","a");
749 } else {
750 convergence = fopen("pcp.convergence.csv","w");
751 if (convergence != NULL)
752 fprintf(convergence,"Ecut\ta.u.\tMaxG\tN0\tN1\tN2\tTE\tVolume\tN^3\tSawStart\n");
753 }
754 if (convergence != NULL) {
755 fprintf(convergence,"%e\t%e\t%i\t%i\t%i\t%i\t%e\t%e\t%i\t%e\n",
756 Lat->ECut,
757 sqrt(Lat->RealBasisSQ[0]), LevS->TotalAllMaxG, Lev0->Plan0.plan->N[0], Lev0->Plan0.plan->N[1], Lev0->Plan0.plan->N[2],
758 Lat->Energy[Occupied].TotalEnergy[0],
759 Lat->Volume, Lev0->Plan0.plan->N[0]*Lev0->Plan0.plan->N[1]*Lev0->Plan0.plan->N[2],
760 P->Lat.SawtoothStart);
761 fclose(convergence);
762 }
763 }
764 }
765}
Note: See TracBrowser for help on using the repository browser.