source: pcp/src/perturbed.c@ 1d77026

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

MinImageConv(): now transform whole vector and does not return single axis value

Due to the switch from hoping for a rectangular simulation box to using Matrix trafo with Real- and ReciBasis, MinImageConv may use these trafos as well and hence may act on the whole vector at once, yielding the minimum image convention in all three directions.
This has impact on numerous functions inside perturbed.c

  • Property mode set to 100644
File size: 217.7 KB
Line 
1/** \file perturbed.c
2 * Perturbation calculation due to external magnetic field.
3 *
4 * Central function is MinimisePerturbed() wherein the actual minimisation of the two different operators with each
5 * three components takes place subsequently. Helpful routines are CalculatePerturbationOperator_P() - which applies a
6 * specified component of p on the current wave function - and CalculatePerturbationOperator_RxP() - which does the
7 * same for the RxP operator.
8 * The actual minimisation loop FindPerturbedMinimum() depends on the same routines also used for the occupied orbitals,
9 * however with a different energy functional and derivatives, evaluated in Calculate1stPerturbedDerivative() and
10 * Calculate2ndPerturbedDerivative(). InitPerturbedEnergyCalculation() calculates the total energy functional
11 * perturbed in second order for all wave functions, UpdatePerturbedEnergyCalculation() just updates the one
12 * for the wave function after it has been minimised during the line search. Both use CalculatePerturbedEnergy() which
13 * evaluates the energy functional (and the gradient if specified).
14 * Finally, FillCurrentDensity() evaluates the current density at a given point in space using the perturbed
15 * wave functions. Afterwards by calling CalculateMagneticSusceptibility() or
16 * CalculateChemicalShieldingByReciprocalCurrentDensity() susceptibility respectively shielding tensor are possible uses
17 * of this current density.
18 *
19 * There are also some test routines: TestCurrent() checks whether the integrated current is zero in each component.
20 * test_fft_symmetry() tests the "pulling out imaginary unit" before fourier transformation on a given wave function.
21 * CheckOrbitalOverlap() outputs the overlap matrix for the wave functions of a given minimisation state, this might
22 * be important for the additional \f$\Delta J{ij}\f$ contribution to the current density, which is non-zero for
23 * non-zero mutual overlap, which is evaluated if FillDeltaCurrentDensity() is called.
24 *
25 * Finally, there are also some smaller routines: truedist() gives the correct relative distance between two points
26 * in the unit cell under periodic boundary conditions with minimum image convention. ApplyTotalHamiltonian() returns
27 * the hamiltonian applied to a given wave function. sawtooth() is a sawtooth implementation which is needed in order
28 * to avoid flipping of position eigenvalues for nodes close to or on the cell boundary. CalculateOverlap()
29 * is used in the energy functional derivatives, keeping an overlap table between perturbed wave functions up to date.
30 * fft_Psi() is very similar to CalculateOneDensityR(), it does the extension of the wave function to the upper level
31 * RunStruct#Lev0 while fouriertransforming it to real space. cross() gives correct indices in evaluating a vector cross
32 * product. AllocCurrentDensity() and DisAllocCurrentDensity() mark the current density arrays as currently being in use or not.
33 *
34 Project: ParallelCarParrinello
35 \author Frederik Heber
36 \date 2006
37
38*/
39
40#include <stdlib.h>
41#include <stdio.h>
42#include <math.h>
43#include <string.h>
44#include <time.h>
45#include <gsl/gsl_matrix.h>
46#include <gsl/gsl_eigen.h>
47#include <gsl/gsl_complex.h>
48#include <gsl/gsl_complex_math.h>
49#include <gsl/gsl_sort_vector.h>
50#include <gsl/gsl_linalg.h>
51#include <gsl/gsl_multimin.h>
52
53#include "data.h"
54#include "density.h"
55#include "energy.h"
56#include "excor.h"
57#include "errors.h"
58#include "grad.h"
59#include "gramsch.h"
60#include "mergesort2.h"
61#include "helpers.h"
62#include "init.h"
63#include "myfft.h"
64#include "mymath.h"
65#include "output.h"
66#include "pcp.h"
67#include "perturbed.h"
68#include "run.h"
69#include "wannier.h"
70
71
72/** Minimisation of the PsiTypeTag#Perturbed_RxP0, PsiTypeTag#Perturbed_P0 and other orbitals.
73 * For each of the above PsiTypeTag we go through the following before the minimisation loop:
74 * -# ResetGramSchTagType() resets current type that is to be minimised to NotOrthogonal.
75 * -# UpdateActualPsiNo() steps on to next perturbed of current PsiTypeTag type.
76 * -# GramSch() orthonormalizes perturbed wave functions.
77 * -# TestGramSch() tests if orthonormality was achieved.
78 * -# InitDensityCalculation() gathers densities from all wave functions (and all processes), within SpeedMeasure() DensityTime.
79 * -# InitPerturbedEnergyCalculation() performs initial calculation of the perturbed energy functional.
80 * -# RunStruct#OldActualLocalPsiNo is set to RunStruct#ActualLocalPsiNo, immediately followed by UpdateGramSchOldActualPsiNo()
81 * to bring info on all processes on par.
82 * -# UpdatePerturbedEnergyCalculation() re-calculates Gradient and GradientTypes#H1cGradient for RunStruct#ActualLocalPsiNo
83 * -# EnergyAllReduce() gathers various energy terms and sums up into Energy#TotalEnergy.
84 *
85 * And during the minimisation loop:
86 * -# FindPerturbedMinimum() performs the gradient conjugation, the line search and wave function update.
87 * -# UpdateActualPsiNo() steps on to the next wave function, orthonormalizing by GramSch() if necessary.
88 * -# UpdateEnergyArray() shifts TotalEnergy values to make space for new one.
89 * -# There is no density update as the energy function does not depend on the changing perturbed density but only on the fixed
90 * unperturbed one.
91 * -# UpdatePerturbedEnergyCalculation() re-calculates the perturbed energy of the changed wave function.
92 * -# EnergyAllReduce() gathers energy terms and sums up.
93 * -# CheckCPULIM() checks if external Stop signal has been given.
94 * -# CalculateMinimumStop() checks whether we have dropped below a certain minimum change during minimisation of total energy.
95 * -# finally step counters LatticeLevel#Step and SpeedStruct#Steps are increased.
96 *
97 * After the minimisation loop:
98 * -# SetGramSchExtraPsi() removes extra Psis from orthogonaliy check.
99 * -# ResetGramSchTagType() sets GramSchToDoType to NotUsedtoOrtho.
100 *
101 * And after all minimisation runs are done:
102 * -# UpdateActualPsiNo() steps back to PsiTypeTag#Occupied type.
103 *
104 * At the end we return to Occupied wave functions.
105 * \param *P at hand
106 * \param *Stop flag to determine if epsilon stop conditions have met
107 * \param *SuperStop flag to determinte whether external signal's required end of calculations
108 */
109void MinimisePerturbed (struct Problem *P, int *Stop, int *SuperStop) {
110 struct RunStruct *R = &P->R;
111 struct Lattice *Lat = &P->Lat;
112 struct Psis *Psi = &Lat->Psi;
113 int type, flag = 0;//,i;
114
115 for (type=Perturbed_P0;type<=Perturbed_RxP2;type++) { // go through each perturbation group separately //
116 *Stop=0; // reset stop flag
117 if(P->Call.out[LeaderOut]) fprintf(stderr,"(%i)Beginning perturbed minimisation of type %s ...\n", P->Par.me, R->MinimisationName[type]);
118 //OutputOrbitalPositions(P, Occupied);
119 R->PsiStep = R->MaxPsiStep; // reset in-Psi-minimisation-counter, so that we really advance to the next wave function
120 UpdateActualPsiNo(P, type); // step on to next perturbed one
121
122 if(P->Call.out[MinOut]) fprintf(stderr, "(%i) Re-initializing perturbed psi array for type %s ", P->Par.me, R->MinimisationName[type]);
123 if (P->Call.ReadSrcFiles && (flag = ReadSrcPsiDensity(P,type,1, R->LevS->LevelNo))) {// in flag store whether stored Psis are readible or not
124 SpeedMeasure(P, InitSimTime, StartTimeDo);
125 if(P->Call.out[MinOut]) fprintf(stderr,"from source file of recent calculation\n");
126 ReadSrcPsiDensity(P,type, 0, R->LevS->LevelNo);
127 ResetGramSchTagType(P, Psi, type, IsOrthogonal); // loaded values are orthonormal
128 SpeedMeasure(P, DensityTime, StartTimeDo);
129 //InitDensityCalculation(P);
130 SpeedMeasure(P, DensityTime, StopTimeDo);
131 R->OldActualLocalPsiNo = R->ActualLocalPsiNo; // needed otherwise called routines in function below crash
132 UpdateGramSchOldActualPsiNo(P,Psi);
133 InitPerturbedEnergyCalculation(P, 1); // go through all orbitals calculate each H^{(0)}-eigenvalue, recalc HGDensity, cause InitDensityCalc zero'd it
134 UpdatePerturbedEnergyCalculation(P); // H1cGradient and Gradient must be current ones
135 EnergyAllReduce(P); // gather energies for minimum search
136 SpeedMeasure(P, InitSimTime, StopTimeDo);
137 }
138 if ((P->Call.ReadSrcFiles != 1) || (!flag)) { // read and don't minimise only if SrcPsi were parsable!
139 SpeedMeasure(P, InitSimTime, StartTimeDo);
140 ResetGramSchTagType(P, Psi, type, NotOrthogonal); // perturbed now shall be orthonormalized
141 if ((P->Call.ReadSrcFiles != 2) || (!flag)) {
142 if (R->LevSNo == Lat->MaxLevel-1) { // is it the starting level? (see InitRunLevel())
143 if(P->Call.out[MinOut]) fprintf(stderr, "randomly.\n");
144 InitPsisValue(P, Psi->TypeStartIndex[type], Psi->TypeStartIndex[type+1]); // initialize perturbed array for this run
145 } else {
146 if(P->Call.out[MinOut]) fprintf(stderr, "from source file of last level.\n");
147 ReadSrcPerturbedPsis(P, type);
148 }
149 }
150 SpeedMeasure(P, InitGramSchTime, StartTimeDo);
151 GramSch(P, R->LevS, Psi, Orthogonalize);
152 SpeedMeasure(P, InitGramSchTime, StopTimeDo);
153 SpeedMeasure(P, InitDensityTime, StartTimeDo);
154 //InitDensityCalculation(P);
155 SpeedMeasure(P, InitDensityTime, StopTimeDo);
156 InitPerturbedEnergyCalculation(P, 1); // go through all orbitals calculate each H^{(0)}-eigenvalue, recalc HGDensity, cause InitDensityCalc zero'd it
157 R->OldActualLocalPsiNo = R->ActualLocalPsiNo; // needed otherwise called routines in function below crash
158 UpdateGramSchOldActualPsiNo(P,Psi);
159 UpdatePerturbedEnergyCalculation(P); // H1cGradient and Gradient must be current ones
160 EnergyAllReduce(P); // gather energies for minimum search
161 SpeedMeasure(P, InitSimTime, StopTimeDo);
162 R->LevS->Step++;
163 EnergyOutput(P,0);
164 while (*Stop != 1) {
165 //debug(P,"FindPerturbedMinimum");
166 FindPerturbedMinimum(P); // find minimum
167 //debug(P,"UpdateActualPsiNo");
168 UpdateActualPsiNo(P, type); // step on to next perturbed Psi
169 //debug(P,"UpdateEnergyArray");
170 UpdateEnergyArray(P); // shift energy values in their array by one
171 //debug(P,"UpdatePerturbedEnergyCalculation");
172 UpdatePerturbedEnergyCalculation(P); // re-calc energies (which is hopefully lower)
173 EnergyAllReduce(P); // gather from all processes and sum up to total energy
174 //ControlNativeDensity(P); // check total density (summed up PertMixed must be zero!)
175 //printf ("(%i,%i,%i)S(%i,%i,%i):\t %5d %10.5f\n",P->Par.my_color_comm_ST,P->Par.me_comm_ST, P->Par.me_comm_ST_PsiT, R->MinStep, R->ActualLocalPsiNo, R->PsiStep, (int)iter, s_multi->f);
176 if (*SuperStop != 1)
177 *SuperStop = CheckCPULIM(P);
178 *Stop = CalculateMinimumStop(P, *SuperStop);
179 P->Speed.Steps++; // step on
180 R->LevS->Step++;
181 }
182 // now release normalization condition and minimize wrt to norm
183 if(P->Call.out[MinOut]) fprintf(stderr,"(%i) Writing %s srcpsi to disk\n", P->Par.me, R->MinimisationName[type]);
184 OutputSrcPsiDensity(P, type);
185// if (!TestReadnWriteSrcDensity(P,type))
186// Error(SomeError,"TestReadnWriteSrcDensity failed!");
187 }
188
189 TestGramSch(P,R->LevS,Psi, type); // functions are orthonormal?
190 // calculate current density summands
191 //if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) Filling current density grid ...\n",P->Par.me);
192 SpeedMeasure(P, CurrDensTime, StartTimeDo);
193 if (*SuperStop != 1) {
194 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
195 R->DoFullCurrent = 1; // set to 1 if it was 2 but Check...() yielded necessity
196 //debug(P,"Filling with Delta j ...");
197 //FillDeltaCurrentDensity(P);
198 }// else
199 //debug(P,"There is no overlap between orbitals.");
200 //debug(P,"Filling with j ...");
201 FillCurrentDensity(P);
202 }
203 SpeedMeasure(P, CurrDensTime, StopTimeDo);
204
205 SetGramSchExtraPsi(P,Psi,NotUsedToOrtho); // remove extra Psis from orthogonality check
206 ResetGramSchTagType(P, Psi, type, NotUsedToOrtho); // remove this group from the check for the next minimisation group as well!
207 }
208 UpdateActualPsiNo(P, Occupied); // step on back to an occupied one
209}
210
211/** Tests overlap matrix between each pair of orbitals for non-diagonal form.
212 * We simply check whether the overlap matrix Psis#lambda has off-diagonal entries greater MYEPSILON or not.
213 * \param *P Problem at hand
214 * \note The routine is meant as atest criteria if \f$\Delta J_[ij]\f$ contribution is necessary, as it is only non-zero if
215 * there is mutual overlap between the two orbitals.
216 */
217int CheckOrbitalOverlap(struct Problem *P)
218{
219 struct Lattice *Lat = &P->Lat;
220 struct Psis *Psi = &Lat->Psi;
221 int i,j;
222 int counter = 0;
223
224 // output matrix
225 if (P->Par.me == 0) fprintf(stderr, "(%i) S_ij =\n", P->Par.me);
226 for (i=0;i<Psi->NoOfPsis;i++) {
227 for (j=0;j<Psi->NoOfPsis;j++) {
228 if (fabs(Psi->lambda[i][j]) > MYEPSILON) counter++;
229 if (P->Par.me == 0) fprintf(stderr, "%e\t", Psi->lambda[i][j]); //Overlap[i][j]
230 }
231 if (P->Par.me == 0) fprintf(stderr, "\n");
232 }
233
234 fprintf(stderr, "(%i) CheckOverlap: %i overlaps found.\t", P->Par.me, counter);
235 if (counter > 0) return (1);
236 else return(0);
237}
238
239/** Initialization of perturbed energy.
240 * For each local wave function of the current minimisation type RunStruct#CurrentMin it is called:
241 * - CalculateNonLocalEnergyNoRT(): for the coefficient-dependent form factors
242 * - CalculatePerturbedEnergy(): for the perturbed energy, yet without gradient calculation
243 * - CalculateOverlap(): for the overlap between the perturbed wave functions of the current RunStruct#CurrentMin state.
244 *
245 * Afterwards for the two types AllPsiEnergyTypes#Perturbed1_0Energy and AllPsiEnergyTypes#Perturbed0_1Energy the
246 * energy contribution from each wave function is added up in Energy#AllLocalPsiEnergy.
247 * \param *P Problem at hand
248 * \param first state whether it is the first (1) or successive call (0), which avoids some initial calculations.
249 * \sa UpdatePerturbedEnergy()
250 * \note Afterwards EnergyAllReduce() must be called.
251 */
252void InitPerturbedEnergyCalculation(struct Problem *P, const int first)
253{
254 struct Lattice *Lat = &(P->Lat);
255 int p,i;
256 const enum PsiTypeTag state = P->R.CurrentMin;
257 for (p=Lat->Psi.TypeStartIndex[state]; p < Lat->Psi.TypeStartIndex[state+1]; p++) {
258 //if (p < 0 || p >= Lat->Psi.LocalNo) Error(SomeError,"InitPerturbedEnergyCalculation: p out of range");
259 //CalculateNonLocalEnergyNoRT(P, p); // recalculating non-local form factors which are coefficient dependent!
260 CalculatePsiEnergy(P,p,1);
261 CalculatePerturbedEnergy(P, p, 0, first);
262 CalculateOverlap(P, p, state);
263 }
264 for (i=0; i<= Perturbed0_1Energy; i++) {
265 Lat->E->AllLocalPsiEnergy[i] = 0.0;
266 for (p=0; p < Lat->Psi.LocalNo; p++)
267 if (P->Lat.Psi.LocalPsiStatus[p].PsiType == state)
268 Lat->E->AllLocalPsiEnergy[i] += Lat->E->PsiEnergy[i][p];
269 }
270}
271
272
273/** Updating of perturbed energy.
274 * For current and former (if not the same) local wave function RunStruct#ActualLocal, RunStruct#OldActualLocalPsiNo it is called:
275 * - CalculateNonLocalEnergyNoRT(): for the form factors
276 * - CalculatePerturbedEnergy(): for the perturbed energy, gradient only for RunStruct#ActualLocal
277 * - CalculatePerturbedOverlap(): for the overlap between the perturbed wave functions
278 *
279 * Afterwards for the two types AllPsiEnergyTypes#Perturbed1_0Energy and AllPsiEnergyTypes#Perturbed0_1Energy the
280 * energy contribution from each wave function is added up in Energy#AllLocalPsiEnergy.
281 * \param *P Problem at hand
282 * \sa CalculatePerturbedEnergy() called from here.
283 * \note Afterwards EnergyAllReduce() must be called.
284 */
285void UpdatePerturbedEnergyCalculation(struct Problem *P)
286{
287 struct Lattice *Lat = &(P->Lat);
288 struct Psis *Psi = &Lat->Psi;
289 struct RunStruct *R = &P->R;
290 const enum PsiTypeTag state = R->CurrentMin;
291 int p = R->ActualLocalPsiNo;
292 const int p_old = R->OldActualLocalPsiNo;
293 int i;
294
295 if (p != p_old) {
296 //if (p_old < 0 || p_old >= Lat->Psi.LocalNo) Error(SomeError,"UpdatePerturbedEnergyCalculation: p_old out of range");
297 //CalculateNonLocalEnergyNoRT(P, p_old);
298 CalculatePsiEnergy(P,p_old,0);
299 CalculatePerturbedEnergy(P, p_old, 0, 0);
300 CalculateOverlap(P, p_old, state);
301 }
302 //if (p < 0 || p >= Lat->Psi.LocalNo) Error(SomeError,"InitPerturbedEnergyCalculation: p out of range");
303 // recalculating non-local form factors which are coefficient dependent!
304 //CalculateNonLocalEnergyNoRT(P,p);
305 CalculatePsiEnergy(P,p,0);
306 CalculatePerturbedEnergy(P, p, 1, 0);
307 CalculateOverlap(P, p, state);
308
309 for (i=0; i<= Perturbed0_1Energy; i++) {
310 Lat->E->AllLocalPsiEnergy[i] = 0.0;
311 for (p=0; p < Psi->LocalNo; p++)
312 if (Psi->LocalPsiStatus[p].PsiType == state)
313 Lat->E->AllLocalPsiEnergy[i] += Lat->E->PsiEnergy[i][p];
314 }
315}
316
317/** Calculates gradient and evaluates second order perturbed energy functional for specific wave function.
318 * The in second order perturbed energy functional reads as follows.
319 * \f[
320 * E^{(2)} = \sum_{kl} \langle \varphi_k^{(1)} | H^{(0)} \delta_{kl} - \lambda_{kl} | \varphi_l^{(1)} \rangle
321 * + \underbrace{\langle \varphi_l^{(0)} | H^{(1)} | \varphi_l^{(1)} \rangle + \langle \varphi_l^{(1)} | H^{(1)} | \varphi_l^{(0)} \rangle}_{2 {\cal R} \langle \varphi_l^{(1)} | H^{(1)} | \varphi_l^{(0)} \rangle}
322 * \f]
323 * And the gradient
324 * \f[
325 * \widetilde{\varphi}_k^{(1)} = - \sum_l ({\cal H}^{(0)} \delta_{kl} - \lambda_{kl} | \varphi_l^{(1)} \rangle + {\cal H}^{(1)} | \varphi_k^{(0)} \rangle
326 * \f]
327 * First, the HGDensity is recalculated if \a first says so - see ApplyTotalHamiltonian().
328 *
329 * Next, we need the perturbation hamiltonian acting on both the respective occupied and current wave function,
330 * see perturbed.c for respective function calls.
331 *
332 * Finally, the scalar product between the wave function and Hc_Gradient yields the eigenvalue of the hamiltonian,
333 * which is summed up over all reciprocal grid vectors and stored in OnePsiElementAddData#Lambda. The Gradient is
334 * the inverse of Hc_Gradient and with the following summation over all perturbed wave functions (MPI exchange of
335 * non-local coefficients) the gradient is computed. Here we need Psis#lambda, which is computed in CalculateHamiltonian().
336 *
337 * Also \f${\cal H}^{(1)} | \varphi_l^{(0)} \rangle\f$ is stored in GradientTypes#H1cGradient.
338 * \param *P Problem at hand, contains RunStruct, Lattice, LatticeLevel RunStruct#LevS
339 * \param l offset of perturbed wave function within Psi#LocalPsiStatus (\f$\varphi_l^{(1)}\f$)
340 * \param DoGradient (1 = yes, 0 = no) whether gradient shall be calculated or not
341 * \param first recaculate HGDensity (1) or not (0)
342 * \note DensityTypes#ActualPsiDensity must be recent for gradient calculation!
343 * \sa CalculateGradientNoRT() - same procedure for evaluation of \f${\cal H}^{(0)}| \varphi_l^{(1)} \rangle\f$
344 * \note without the simplification of \f$2 {\cal R} \langle \varphi_l^{(1)} | H^{(1)} | \varphi_l^{(0)} \rangle\f$ the
345 * calculation would be impossible due to non-local nature of perturbed wave functions. The position operator would
346 * be impossible to apply in a sensible manner.
347 */
348void CalculatePerturbedEnergy(struct Problem *P, const int l, const int DoGradient, const int first)
349{
350 struct Lattice *Lat = &P->Lat;
351 struct Psis *Psi = &Lat->Psi;
352 struct Energy *E = Lat->E;
353 struct PseudoPot *PP = &P->PP;
354 struct RunStruct *R = &P->R;
355 struct LatticeLevel *LevS = R->LevS;
356 const int state = R->CurrentMin;
357 const int l_normal = Psi->TypeStartIndex[Occupied] + (l - Psi->TypeStartIndex[state]); // offset l to \varphi_l^{(0)}
358 const int ActNum = l - Psi->TypeStartIndex[state] + Psi->TypeStartIndex[1] * Psi->LocalPsiStatus[l].my_color_comm_ST_Psi;
359 int g, i, m, j;
360 double lambda, Lambda;
361 double RElambda10, RELambda10;
362 const fftw_complex *source = LevS->LPsi->LocalPsi[l];
363 fftw_complex *grad = P->Grad.GradientArray[ActualGradient];
364 fftw_complex *Hc_grad = P->Grad.GradientArray[HcGradient];
365 fftw_complex *H1c_grad = P->Grad.GradientArray[H1cGradient];
366 fftw_complex *TempPsi_0 = H1c_grad;
367 fftw_complex *varphi_1, *varphi_0;
368 struct OnePsiElement *OnePsiB, *LOnePsiB;
369 fftw_complex *LPsiDatB=NULL;
370 const int ElementSize = (sizeof(fftw_complex) / sizeof(double));
371 int RecvSource;
372 MPI_Status status;
373
374 // ============ Calculate H^(0) psi^(1) =============================
375 //if (Hc_grad != P->Grad.GradientArray[HcGradient]) Error(SomeError,"CalculatePerturbedEnergy: Hc_grad corrupted");
376 SetArrayToDouble0((double *)Hc_grad,2*R->InitLevS->MaxG);
377 ApplyTotalHamiltonian(P,source,Hc_grad, PP->fnl[l], 1, first);
378
379 // ============ ENERGY FUNCTIONAL Evaluation PART 1 ================
380 //if (l_normal < 0 || l_normal >= Psi->LocalNo) Error(SomeError,"CalculatePerturbedEnergy: l_normal out of range");
381 varphi_0 = LevS->LPsi->LocalPsi[l_normal];
382 //if (l < 0 || l >= Psi->LocalNo) Error(SomeError,"CalculatePerturbedEnergy: l out of range");
383 varphi_1 = LevS->LPsi->LocalPsi[l];
384 //if (TempPsi_0 != P->Grad.GradientArray[H1cGradient]) Error(SomeError,"CalculatePerturbedEnergy: TempPsi_0 corrupted");
385 SetArrayToDouble0((double *)TempPsi_0,2*R->InitLevS->MaxG);
386 switch (state) {
387 case Perturbed_P0:
388 CalculatePerturbationOperator_P(P,varphi_0,TempPsi_0,0); // \nabla_0 | \varphi_l^{(0)} \rangle
389 break;
390 case Perturbed_P1:
391 CalculatePerturbationOperator_P(P,varphi_0,TempPsi_0,1); // \nabla_1 | \varphi_l^{(0)} \rangle
392 break;
393 case Perturbed_P2:
394 CalculatePerturbationOperator_P(P,varphi_0,TempPsi_0,2); // \nabla_1 | \varphi_l^{(0)} \rangle
395 break;
396 case Perturbed_RxP0:
397 CalculatePerturbationOperator_RxP(P,varphi_0,TempPsi_0,l_normal,0); // r \times \nabla | \varphi_l^{(0)} \rangle
398 break;
399 case Perturbed_RxP1:
400 CalculatePerturbationOperator_RxP(P,varphi_0,TempPsi_0,l_normal,1); // r \times \nabla | \varphi_l^{(0)} \rangle
401 break;
402 case Perturbed_RxP2:
403 CalculatePerturbationOperator_RxP(P,varphi_0,TempPsi_0,l_normal,2); // r \times \nabla | \varphi_l^{(0)} \rangle
404 break;
405 default:
406 fprintf(stderr,"(%i) CalculatePerturbedEnergy called whilst not within perturbation run: CurrentMin = %i !\n",P->Par.me, R->CurrentMin);
407 break;
408 }
409
410 // ============ GRADIENT and EIGENVALUE Evaluation Part 1==============
411 lambda = 0.0;
412 if ((DoGradient) && (grad != NULL)) {
413 g = 0;
414 if (LevS->GArray[0].GSq == 0.0) {
415 lambda += Hc_grad[0].re*source[0].re;
416 //if (grad != P->Grad.GradientArray[ActualGradient]) Error(SomeError,"CalculatePerturbedEnergy: grad corrupted");
417 grad[0].re = -(Hc_grad[0].re + TempPsi_0[0].re);
418 grad[0].im = -(Hc_grad[0].im + TempPsi_0[0].im);
419 g++;
420 }
421 for (;g<LevS->MaxG;g++) {
422 lambda += 2.*(Hc_grad[g].re*source[g].re + Hc_grad[g].im*source[g].im);
423 //if (grad != P->Grad.GradientArray[ActualGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"CalculatePerturbedEnergy: grad corrupted");
424 grad[g].re = -(Hc_grad[g].re + TempPsi_0[g].re);
425 grad[g].im = -(Hc_grad[g].im + TempPsi_0[g].im);
426 }
427
428 m = -1;
429 for (j=0; j < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; j++) { // go through all wave functions
430 OnePsiB = &Psi->AllPsiStatus[j]; // grab OnePsiB
431 if (OnePsiB->PsiType == state) { // drop all but the ones of current min state
432 m++; // increase m if it is type-specific wave function
433 if (OnePsiB->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // local?
434 LOnePsiB = &Psi->LocalPsiStatus[OnePsiB->MyLocalNo];
435 else
436 LOnePsiB = NULL;
437 if (LOnePsiB == NULL) { // if it's not local ... receive it from respective process into TempPsi
438 RecvSource = OnePsiB->my_color_comm_ST_Psi;
439 MPI_Recv( LevS->LPsi->TempPsi, LevS->MaxG*ElementSize, MPI_DOUBLE, RecvSource, PerturbedTag, P->Par.comm_ST_PsiT, &status );
440 LPsiDatB=LevS->LPsi->TempPsi;
441 } else { // .. otherwise send it to all other processes (Max_me... - 1)
442 for (i=0;i<P->Par.Max_me_comm_ST_PsiT;i++)
443 if (i != OnePsiB->my_color_comm_ST_Psi)
444 MPI_Send( LevS->LPsi->LocalPsi[OnePsiB->MyLocalNo], LevS->MaxG*ElementSize, MPI_DOUBLE, i, PerturbedTag, P->Par.comm_ST_PsiT);
445 LPsiDatB=LevS->LPsi->LocalPsi[OnePsiB->MyLocalNo];
446 } // LPsiDatB is now set to the coefficients of OnePsi either stored or MPI_Received
447
448 g = 0;
449 if (LevS->GArray[0].GSq == 0.0) { // perform the summation
450 //if (grad != P->Grad.GradientArray[ActualGradient]) Error(SomeError,"CalculatePerturbedEnergy: grad corrupted");
451 grad[0].re += Lat->Psi.lambda[ActNum][m]*LPsiDatB[0].re;
452 grad[0].im += Lat->Psi.lambda[ActNum][m]*LPsiDatB[0].im;
453 g++;
454 }
455 for (;g<LevS->MaxG;g++) {
456 //if (grad != P->Grad.GradientArray[ActualGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"CalculatePerturbedEnergy: grad corrupted");
457 grad[g].re += Lat->Psi.lambda[ActNum][m]*LPsiDatB[g].re;
458 grad[g].im += Lat->Psi.lambda[ActNum][m]*LPsiDatB[g].im;
459 }
460 }
461 }
462 } else {
463 lambda = GradSP(P,LevS,Hc_grad,source);
464 }
465 MPI_Allreduce ( &lambda, &Lambda, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
466 //fprintf(stderr,"(%i) Lambda[%i] = %lg\n",P->Par.me, l, Lambda);
467 //if (l < 0 || l >= Psi->LocalNo) Error(SomeError,"CalculatePerturbedEnergy: l out of range");
468 Lat->Psi.AddData[l].Lambda = Lambda;
469
470 // ============ ENERGY FUNCTIONAL Evaluation PART 2 ================
471 // varphi_1 jas negative symmetry, returning TempPsi_0 from CalculatePerturbedOperator also, thus real part of scalar product
472 // "-" due to purely imaginary wave function is on left hand side, thus becomes complex conjugated: i -> -i
473 // (-i goes into pert. op., "-" remains when on right hand side)
474 RElambda10 = GradSP(P,LevS,varphi_1,TempPsi_0) * sqrt(Psi->LocalPsiStatus[l].PsiFactor * Psi->LocalPsiStatus[l_normal].PsiFactor);
475 //RElambda01 = -GradSP(P,LevS,varphi_0,TempPsi_1) * sqrt(Psi->LocalPsiStatus[l].PsiFactor * Psi->LocalPsiStatus[l_normal].PsiFactor);
476
477 MPI_Allreduce ( &RElambda10, &RELambda10, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
478 //MPI_Allreduce ( &RElambda01, &RELambda01, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
479
480 //if (l < 0 || l >= Psi->LocalNo) Error(SomeError,"CalculatePerturbedEnergy: l out of range");
481 E->PsiEnergy[Perturbed1_0Energy][l] = RELambda10;
482 E->PsiEnergy[Perturbed0_1Energy][l] = RELambda10;
483// if (P->Par.me == 0) {
484// fprintf(stderr,"RE.Lambda10[%i-%i] = %lg\t RE.Lambda01[%i-%i] = %lg\n", l, l_normal, RELambda10, l_normal, l, RELambda01);
485// }
486 // GradImSP() is only applicable to a product of wave functions with uneven symmetry!
487 // Otherwise, due to the nature of symmetry, a sum over only half of the coefficients will in most cases not result in zero!
488}
489
490/** Applies \f$H^{(0)}\f$ to a given \a source.
491 * The DensityTypes#HGDensity is computed, the exchange potential added and the
492 * whole multiplied - coefficient by coefficient - with the current wave function, taken from its density coefficients,
493 * on the upper LatticeLevel (RunStruct#Lev0), which (DensityTypes#ActualPsiDensity) is updated beforehand.
494 * After an inverse fft (now G-dependent) the non-local potential is added and
495 * within the reciprocal basis set, the kinetic energy can be evaluated easily.
496 * \param *P Problem at hand
497 * \param *source pointer to source coefficient array, \f$| \varphi(G) \rangle\f$
498 * \param *dest pointer to dest coefficient array,\f$H^{(0)} | \varphi(G) \rangle\f$
499 * \param **fnl pointer to non-local form factor array
500 * \param PsiFactor occupation number of orbital
501 * \param first 1 - Re-calculate DensityTypes#HGDensity, 0 - don't
502 * \sa CalculateConDirHConDir() - same procedure
503 */
504void ApplyTotalHamiltonian(struct Problem *P, const fftw_complex *source, fftw_complex *dest, fftw_complex ***fnl, const double PsiFactor, const int first) {
505 struct Lattice *Lat = &P->Lat;
506 struct RunStruct *R = &P->R;
507 struct LatticeLevel *LevS = R->LevS;
508 struct LatticeLevel *Lev0 = R->Lev0;
509 struct Density *Dens0 = Lev0->Dens;
510 struct fft_plan_3d *plan = Lat->plan;
511 struct PseudoPot *PP = &P->PP;
512 struct Ions *I = &P->Ion;
513 fftw_complex *work = Dens0->DensityCArray[TempDensity];
514 fftw_real *HGcR = Dens0->DensityArray[HGcDensity];
515 fftw_complex *HGcRC = (fftw_complex*)HGcR;
516 fftw_complex *HGC = Dens0->DensityCArray[HGDensity];
517 fftw_real *HGCR = (fftw_real *)HGC;
518 fftw_complex *PsiC = Dens0->DensityCArray[ActualPsiDensity];
519 fftw_real *PsiCR = (fftw_real *)PsiC;
520 //const fftw_complex *dest_bak = dest;
521 int nx,ny,nz,iS,i0;
522 const int Nx = LevS->Plan0.plan->local_nx;
523 const int Ny = LevS->Plan0.plan->N[1];
524 const int Nz = LevS->Plan0.plan->N[2];
525 const int NUpx = LevS->NUp[0];
526 const int NUpy = LevS->NUp[1];
527 const int NUpz = LevS->NUp[2];
528 const double HGcRCFactor = 1./LevS->MaxN;
529 int g, Index, i, it;
530 fftw_complex vp,rp,rhog,TotalPsiDensity;
531 double Fac;
532
533 if (first) {
534 // recalculate HGDensity
535 //if (HGC != Dens0->DensityCArray[HGDensity]) Error(SomeError,"ApplyTotalHamiltonian: HGC corrupted");
536 SetArrayToDouble0((double *)HGC,2*Dens0->TotalSize);
537 g=0;
538 if (Lev0->GArray[0].GSq == 0.0) {
539 Index = Lev0->GArray[0].Index;
540 c_re(vp) = 0.0;
541 c_im(vp) = 0.0;
542 for (it = 0; it < I->Max_Types; it++) {
543 c_re(vp) += (c_re(I->I[it].SFactor[0])*PP->phi_ps_loc[it][0]);
544 c_im(vp) += (c_im(I->I[it].SFactor[0])*PP->phi_ps_loc[it][0]);
545 }
546 //if (HGC != Dens0->DensityCArray[HGDensity] || Index<0 || Index>=Dens0->LocalSizeC) Error(SomeError,"ApplyTotalHamiltonian: HGC corrupted");
547 c_re(HGC[Index]) = c_re(vp);
548 c_re(TotalPsiDensity) = c_re(Dens0->DensityCArray[TotalDensity][Index]);
549 c_im(TotalPsiDensity) = c_im(Dens0->DensityCArray[TotalDensity][Index]);
550
551 g++;
552 }
553 for (; g < Lev0->MaxG; g++) {
554 Index = Lev0->GArray[g].Index;
555 Fac = 4.*PI/(Lev0->GArray[g].GSq);
556 c_re(vp) = 0.0;
557 c_im(vp) = 0.0;
558 c_re(rp) = 0.0;
559 c_im(rp) = 0.0;
560 for (it = 0; it < I->Max_Types; it++) {
561 c_re(vp) += (c_re(I->I[it].SFactor[g])*PP->phi_ps_loc[it][g]);
562 c_im(vp) += (c_im(I->I[it].SFactor[g])*PP->phi_ps_loc[it][g]);
563 c_re(rp) += (c_re(I->I[it].SFactor[g])*PP->FacGauss[it][g]);
564 c_im(rp) += (c_im(I->I[it].SFactor[g])*PP->FacGauss[it][g]);
565 } // rp = n^{Gauss)(G)
566
567 // n^{tot} = n^0 + \lambda n^1 + ...
568 //if (isnan(c_re(Dens0->DensityCArray[TotalDensity][Index]))) { fprintf(stderr,"(%i) WARNING in CalculatePerturbedEnergy(): TotalDensity[%i] = NaN!\n", P->Par.me, Index); Error(SomeError, "NaN-Fehler!"); }
569 c_re(TotalPsiDensity) = c_re(Dens0->DensityCArray[TotalDensity][Index]);
570 c_im(TotalPsiDensity) = c_im(Dens0->DensityCArray[TotalDensity][Index]);
571
572 c_re(rhog) = c_re(TotalPsiDensity)*R->HGcFactor+c_re(rp);
573 c_im(rhog) = c_im(TotalPsiDensity)*R->HGcFactor+c_im(rp);
574 // rhog = n(G) + n^{Gauss}(G), rhoe = n(G)
575 //if (HGC != Dens0->DensityCArray[HGDensity] || Index<0 || Index>=Dens0->LocalSizeC) Error(SomeError,"ApplyTotalHamiltonian: HGC corrupted");
576 c_re(HGC[Index]) = c_re(vp)+Fac*c_re(rhog);
577 c_im(HGC[Index]) = c_im(vp)+Fac*c_im(rhog);
578 }
579 //
580 for (i=0; i<Lev0->MaxDoubleG; i++) {
581 //if (HGC != Dens0->DensityCArray[HGDensity] || Lev0->DoubleG[2*i+1]<0 || Lev0->DoubleG[2*i+1]>Dens0->LocalSizeC || Lev0->DoubleG[2*i]<0 || Lev0->DoubleG[2*i]>Dens0->LocalSizeC) Error(SomeError,"CalculatePerturbedEnergy: grad corrupted");
582 HGC[Lev0->DoubleG[2*i+1]].re = HGC[Lev0->DoubleG[2*i]].re;
583 HGC[Lev0->DoubleG[2*i+1]].im = -HGC[Lev0->DoubleG[2*i]].im;
584 }
585 }
586 // ============ GRADIENT and EIGENVALUE Evaluation Part 1==============
587 // \lambda_l^{(1)} = \langle \varphi_l^{(1)} | H^{(0)} | \varphi_l^{(1)} \rangle and gradient calculation
588 SpeedMeasure(P, LocTime, StartTimeDo);
589 // back-transform HGDensity: (G) -> (R)
590 //if (HGC != Dens0->DensityCArray[HGDensity]) Error(SomeError,"ApplyTotalHamiltonian: HGC corrupted");
591 if (first) fft_3d_complex_to_real(plan, Lev0->LevelNo, FFTNF1, HGC, work);
592 // evaluate exchange potential with this density, add up onto HGCR
593 //if (HGCR != (fftw_real *)Dens0->DensityCArray[HGDensity]) Error(SomeError,"ApplyTotalHamiltonian: HGCR corrupted");
594 if (first) CalculateXCPotentialNoRT(P, HGCR); // add V^{xc} on V^H + V^{ps}
595 // make sure that ActualPsiDensity is recent
596 CalculateOneDensityR(Lat, LevS, Dens0, source, Dens0->DensityArray[ActualDensity], R->FactorDensityR*PsiFactor, 1);
597 for (nx=0;nx<Nx;nx++)
598 for (ny=0;ny<Ny;ny++)
599 for (nz=0;nz<Nz;nz++) {
600 i0 = nz*NUpz+Nz*NUpz*(ny*NUpy+Ny*NUpy*nx*NUpx);
601 iS = nz+Nz*(ny+Ny*nx);
602 //if (HGcR != Dens0->DensityArray[HGcDensity] || iS<0 || iS>=LevS->Dens->LocalSizeR) Error(SomeError,"ApplyTotalHamiltonian: HGC corrupted");
603 HGcR[iS] = HGCR[i0]*PsiCR[i0]; /* Matrix Vector Mult */
604 }
605 // (R) -> (G)
606 //if (HGcRC != (fftw_complex *)Dens0->DensityArray[HGcDensity]) Error(SomeError,"ApplyTotalHamiltonian: HGcRC corrupted");
607 fft_3d_real_to_complex(plan, LevS->LevelNo, FFTNF1, HGcRC, work);
608 SpeedMeasure(P, LocTime, StopTimeDo);
609 /* NonLocalPP */
610 SpeedMeasure(P, NonLocTime, StartTimeDo);
611 //if (dest != dest_bak) Error(SomeError,"ApplyTotalHamiltonian: dest corrupted");
612 CalculateAddNLPot(P, dest, fnl, PsiFactor); // wave function hidden in form factors fnl, also resets Hc_grad beforehand
613 SpeedMeasure(P, NonLocTime, StopTimeDo);
614
615 /* create final vector */
616 for (g=0;g<LevS->MaxG;g++) {
617 Index = LevS->GArray[g].Index; /* FIXME - factoren */
618 //if (dest != dest_bak || g<0 || g>=LevS->MaxG) Error(SomeError,"ApplyTotalHamiltonian: dest corrupted");
619 dest[g].re += PsiFactor*(HGcRC[Index].re*HGcRCFactor + 0.5*LevS->GArray[g].GSq*source[g].re);
620 dest[g].im += PsiFactor*(HGcRC[Index].im*HGcRCFactor + 0.5*LevS->GArray[g].GSq*source[g].im);
621 }
622}
623
624#define stay_above 0.001 //!< value above which the coefficient of the wave function will always remain
625
626/** Finds the minimum of perturbed energy in regards of actual wave function.
627 * The following happens step by step:
628 * -# The Gradient is copied into GradientTypes#GraSchGradient (which is nothing but a pointer to
629 * one array in LPsiDat) and orthonormalized via GramSch() to all occupied wave functions
630 * except to the current perturbed one.
631 * -# Then comes pre-conditioning, analogous to CalculatePreConGrad().
632 * -# The Gradient is projected onto the current perturbed wave function and this is subtracted, i.e.
633 * vector is the conjugate gradient.
634 * -# Finally, Calculate1stPerturbedDerivative() and Calculate2ndPerturbedDerivative() are called and
635 * with these results and the current total energy, CalculateDeltaI() finds the parameter for the one-
636 * dimensional minimisation. The current wave function is set to newly found minimum and approximated
637 * total energy is printed.
638 *
639 * \param *P Problem at hand
640 * \sa CalculateNewWave() and functions therein
641 */
642void FindPerturbedMinimum(struct Problem *P)
643{
644 struct Lattice *Lat = &P->Lat;
645 struct RunStruct *R = &P->R;
646 struct Psis *Psi = &Lat->Psi;
647 struct PseudoPot *PP = &P->PP;
648 struct LatticeLevel *LevS = R->LevS;
649 struct LatticeLevel *Lev0 = R->Lev0;
650 struct Density *Dens = Lev0->Dens;
651 struct Energy *En = Lat->E;
652 struct FileData *F = &P->Files;
653 int g,p,i;
654 int step = R->PsiStep;
655 double *GammaDiv = &Lat->Psi.AddData[R->ActualLocalPsiNo].Gamma;
656 const int ElementSize = (sizeof(fftw_complex) / sizeof(double));
657 fftw_complex *source = LevS->LPsi->LocalPsi[R->ActualLocalPsiNo];
658 fftw_complex *grad = P->Grad.GradientArray[ActualGradient];
659 fftw_complex *GradOrtho = P->Grad.GradientArray[GraSchGradient];
660 fftw_complex *PCgrad = P->Grad.GradientArray[PreConGradient];
661 fftw_complex *PCOrtho = P->Grad.GradientArray[GraSchGradient];
662 fftw_complex *ConDir = P->Grad.GradientArray[ConDirGradient];
663 fftw_complex *ConDir_old = P->Grad.GradientArray[OldConDirGradient];
664 fftw_complex *Ortho = P->Grad.GradientArray[GraSchGradient];
665 const fftw_complex *Hc_grad = P->Grad.GradientArray[HcGradient];
666 const fftw_complex *H1c_grad = P->Grad.GradientArray[H1cGradient];
667 fftw_complex *HConDir = Dens->DensityCArray[ActualDensity];
668 const double PsiFactor = Lat->Psi.LocalPsiStatus[R->ActualLocalPsiNo].PsiFactor;
669 //double Lambda = Lat->Psi.AddData[R->ActualLocalPsiNo].Lambda;
670 double T;
671 double x, K; //, dK;
672 double dS[2], S[2], Gamma, GammaDivOld = *GammaDiv;
673 double LocalSP, PsiSP;
674 double dEdt0, ddEddt0, ConDirHConDir, ConDirConDir;//, sourceHsource;
675 double E0, E, delta;
676 //double E0, E, dE, ddE, delta, dcos, dsin;
677 //double EI, dEI, ddEI, deltaI, dcosI, dsinI;
678 //double HartreeddEddt0, XCddEddt0;
679 double d[4],D[4], Diff;
680 const int Num = Psi->NoOfPsis;
681
682 // ORTHOGONALIZED-GRADIENT
683 for (g=0;g<LevS->MaxG;g++) {
684 //if (GradOrtho != P->Grad.GradientArray[GraSchGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: GradOrtho corrupted");
685 GradOrtho[g].re = grad[g].re; //+Lambda*source[g].re;
686 GradOrtho[g].im = grad[g].im; //+Lambda*source[g].im;
687 }
688 // include the ExtraPsi (which is the GraSchGradient!)
689 SetGramSchExtraPsi(P, Psi, NotOrthogonal);
690 // exclude the minimised Psi
691 SetGramSchActualPsi(P, Psi, NotUsedToOrtho);
692 SpeedMeasure(P, GramSchTime, StartTimeDo);
693 // makes conjugate gradient orthogonal to all other orbits
694 //fprintf(stderr,"CalculateCGGradient: GramSch() for extra orbital\n");
695 GramSch(P, LevS, Psi, Orthogonalize);
696 SpeedMeasure(P, GramSchTime, StopTimeDo);
697 //if (grad != P->Grad.GradientArray[ActualGradient]) Error(SomeError,"FindPerturbedMinimum: grad corrupted");
698 memcpy(grad, GradOrtho, ElementSize*LevS->MaxG*sizeof(double));
699 //memcpy(PCOrtho, GradOrtho, ElementSize*LevS->MaxG*sizeof(double));
700
701 // PRE-CONDITION-GRADIENT
702 //if (fabs(T) < MYEPSILON) T = 1;
703 T = 0.;
704 for (i=0;i<Num;i++)
705 T += Psi->lambda[i][i];
706 for (g=0;g<LevS->MaxG;g++) {
707 x = .5*LevS->GArray[g].GSq;
708 // FIXME: Good way of accessing reciprocal Lev0 Density coefficients on LevS! (not so trivial)
709 //x += sqrt(Dens->DensityCArray[HGDensity][g].re*Dens->DensityCArray[HGDensity][g].re+Dens->DensityCArray[HGDensity][g].im*Dens->DensityCArray[HGDensity][g].im);
710 x -= T/(double)Num;
711 K = x/(x*x+stay_above*stay_above);
712 //if (PCOrtho != P->Grad.GradientArray[GraSchGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: PCOrtho corrupted");
713 c_re(PCOrtho[g]) = K*c_re(grad[g]);
714 c_im(PCOrtho[g]) = K*c_im(grad[g]);
715 }
716 SetGramSchExtraPsi(P, Psi, NotOrthogonal);
717 SpeedMeasure(P, GramSchTime, StartTimeDo);
718 // preconditioned direction is orthogonalized
719 //fprintf(stderr,"CalculatePreConGrad: GramSch() for extra orbital\n");
720 GramSch(P, LevS, Psi, Orthogonalize);
721 SpeedMeasure(P, GramSchTime, StopTimeDo);
722 //if (PCgrad != P->Grad.GradientArray[PreConGradient]) Error(SomeError,"FindPerturbedMinimum: PCgrad corrupted");
723 memcpy(PCgrad, PCOrtho, ElementSize*LevS->MaxG*sizeof(double));
724
725 //debug(P, "Before ConDir");
726 //fprintf(stderr,"|(%i)|^2 = %lg\t |PCgrad|^2 = %lg\t |PCgrad,(%i)| = %lg\n", R->ActualLocalPsiNo, GradSP(P,LevS,source,source),GradSP(P,LevS,PCgrad,PCgrad), R->ActualLocalPsiNo, GradSP(P,LevS,PCgrad,source));
727 // CONJUGATE-GRADIENT
728 LocalSP = GradSP(P, LevS, PCgrad, grad);
729 MPI_Allreduce ( &LocalSP, &PsiSP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
730 *GammaDiv = dS[0] = PsiSP;
731 dS[1] = GammaDivOld;
732 S[0]=dS[0]; S[1]=dS[1];
733 /*MPI_Allreduce ( dS, S, 2, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_PsiT);*/
734 if (step) { // only in later steps is the scalar product used, but always condir stored in oldcondir and Ortho (working gradient)
735 if (fabs(S[1]) < MYEPSILON) fprintf(stderr,"CalculateConDir: S[1] = %lg\n",S[1]);
736 Gamma = S[0]/S[1];
737 if (fabs(S[1]) < MYEPSILON) {
738 if (fabs(S[0]) < MYEPSILON)
739 Gamma = 1.0;
740 else
741 Gamma = 0.0;
742 }
743 for (g=0; g < LevS->MaxG; g++) {
744 //if (ConDir != P->Grad.GradientArray[ConDirGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: ConDir corrupted");
745 c_re(ConDir[g]) = c_re(PCgrad[g]) + Gamma*c_re(ConDir_old[g]);
746 c_im(ConDir[g]) = c_im(PCgrad[g]) + Gamma*c_im(ConDir_old[g]);
747 //if (ConDir_old != P->Grad.GradientArray[OldConDirGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: ConDir_old corrupted");
748 c_re(ConDir_old[g]) = c_re(ConDir[g]);
749 c_im(ConDir_old[g]) = c_im(ConDir[g]);
750 //if (Ortho != P->Grad.GradientArray[GraSchGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: Ortho corrupted");
751 c_re(Ortho[g]) = c_re(ConDir[g]);
752 c_im(Ortho[g]) = c_im(ConDir[g]);
753 }
754 } else {
755 Gamma = 0.0;
756 for (g=0; g < LevS->MaxG; g++) {
757 //if (ConDir != P->Grad.GradientArray[ConDirGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: ConDir corrupted");
758 c_re(ConDir[g]) = c_re(PCgrad[g]);
759 c_im(ConDir[g]) = c_im(PCgrad[g]);
760 //if (ConDir_old != P->Grad.GradientArray[OldConDirGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: ConDir_old corrupted");
761 c_re(ConDir_old[g]) = c_re(ConDir[g]);
762 c_im(ConDir_old[g]) = c_im(ConDir[g]);
763 //if (Ortho != P->Grad.GradientArray[GraSchGradient] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: Ortho corrupted");
764 c_re(Ortho[g]) = c_re(ConDir[g]);
765 c_im(Ortho[g]) = c_im(ConDir[g]);
766 }
767 }
768 // orthonormalize
769 SetGramSchExtraPsi(P, Psi, NotOrthogonal);
770 SpeedMeasure(P, GramSchTime, StartTimeDo);
771 //fprintf(stderr,"CalculateConDir: GramSch() for extra orbital\n");
772 GramSch(P, LevS, Psi, Orthogonalize);
773 SpeedMeasure(P, GramSchTime, StopTimeDo);
774 //if (ConDir != P->Grad.GradientArray[ConDirGradient]) Error(SomeError,"FindPerturbedMinimum: ConDir corrupted");
775 memcpy(ConDir, Ortho, ElementSize*LevS->MaxG*sizeof(double));
776 //debug(P, "Before LineSearch");
777 //fprintf(stderr,"|(%i)|^2 = %lg\t |ConDir|^2 = %lg\t |ConDir,(%i)| = %lg\n", R->ActualLocalPsiNo, GradSP(P,LevS,source,source),GradSP(P,LevS,ConDir,ConDir), R->ActualLocalPsiNo, GradSP(P,LevS,ConDir,source));
778 SetGramSchActualPsi(P, Psi, IsOrthogonal);
779
780 //fprintf(stderr,"(%i) Testing conjugate gradient for Orthogonality ...\n", P->Par.me);
781 //TestForOrth(P,LevS,ConDir);
782
783 // ONE-DIMENSIONAL LINE-SEARCH
784
785 // ========= dE / dt | 0 ============
786 p = Lat->Psi.TypeStartIndex[Occupied] + (R->ActualLocalPsiNo - Lat->Psi.TypeStartIndex[R->CurrentMin]);
787 //if (Hc_grad != P->Grad.GradientArray[HcGradient]) Error(SomeError,"FindPerturbedMinimum: Hc_grad corrupted");
788 //if (H1c_grad != P->Grad.GradientArray[H1cGradient]) Error(SomeError,"FindPerturbedMinimum: H1c_grad corrupted");
789 d[0] = Calculate1stPerturbedDerivative(P, LevS->LPsi->LocalPsi[p], source, ConDir, Hc_grad, H1c_grad);
790 //CalculateConDirHConDir(P, ConDir, PsiFactor, &d[1], &d[2], &d[3]);
791 //if (ConDir != P->Grad.GradientArray[ConDirGradient]) Error(SomeError,"FindPerturbedMinimum: ConDir corrupted");
792 CalculateCDfnl(P, ConDir, PP->CDfnl); // calculate needed non-local form factors
793 //if (HConDir != Dens->DensityCArray[ActualDensity]) Error(SomeError,"FindPerturbedMinimum: HConDir corrupted");
794 SetArrayToDouble0((double *)HConDir,Dens->TotalSize*2);
795 //if (ConDir != P->Grad.GradientArray[ConDirGradient]) Error(SomeError,"FindPerturbedMinimum: ConDir corrupted");
796 ApplyTotalHamiltonian(P,ConDir,HConDir, PP->CDfnl, PsiFactor, 0); // applies H^(0) with total perturbed density!
797 d[1] = GradSP(P,LevS,ConDir,HConDir);
798 d[2] = GradSP(P,LevS,ConDir,ConDir);
799 d[3] = 0.;
800
801 // gather results
802 MPI_Allreduce ( &d, &D, 4, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
803 // ========== ddE / ddt | 0 =========
804 dEdt0 = D[0];
805 for (i=MAXOLD-1; i > 0; i--)
806 En->dEdt0[i] = En->dEdt0[i-1];
807 En->dEdt0[0] = dEdt0;
808 ConDirHConDir = D[1];
809 ConDirConDir = D[2];
810 ddEddt0 = 0.0;
811 //if (ConDir != P->Grad.GradientArray[ConDirGradient]) Error(SomeError,"FindPerturbedMinimum: ConDir corrupted");
812 //if (H1c_grad != P->Grad.GradientArray[H1cGradient]) Error(SomeError,"FindPerturbedMinimum: H1c_grad corrupted");
813 ddEddt0 = Calculate2ndPerturbedDerivative(P, LevS->LPsi->LocalPsi[p], source, ConDir, Lat->Psi.AddData[R->ActualLocalPsiNo].Lambda * Psi->LocalPsiStatus[R->ActualLocalPsiNo].PsiFactor, ConDirHConDir, ConDirConDir);
814
815 for (i=MAXOLD-1; i > 0; i--)
816 En->ddEddt0[i] = En->ddEddt0[i-1];
817 En->ddEddt0[0] = ddEddt0;
818 E0 = En->TotalEnergy[0];
819 // delta
820 //if (isnan(E0)) { fprintf(stderr,"(%i) WARNING in CalculateLineSearch(): E0_%i[%i] = NaN!\n", P->Par.me, i, 0); Error(SomeError, "NaN-Fehler!"); }
821 //if (isnan(dEdt0)) { fprintf(stderr,"(%i) WARNING in CalculateLineSearch(): dEdt0_%i[%i] = NaN!\n", P->Par.me, i, 0); Error(SomeError, "NaN-Fehler!"); }
822 //if (isnan(ddEddt0)) { fprintf(stderr,"(%i) WARNING in CalculateLineSearch(): ddEddt0_%i[%i] = NaN!\n", P->Par.me, i, 0); Error(SomeError, "NaN-Fehler!"); }
823
824 ////deltaI = CalculateDeltaI(E0, dEdt0, ddEddt0,
825 //// &EI, &dEI, &ddEI, &dcosI, &dsinI);
826 ////delta = deltaI; E = EI; dE = dEI; ddE = ddEI; dcos = dcosI; dsin = dsinI;
827 if (ddEddt0 > 0) {
828 delta = - dEdt0/ddEddt0;
829 E = E0 + delta * dEdt0 + delta*delta/2. * ddEddt0;
830 } else {
831 delta = 0.;
832 E = E0;
833 fprintf(stderr,"(%i) Taylor approximation leads not to minimum!\n",P->Par.me);
834 }
835
836 // shift energy delta values
837 for (i=MAXOLD-1; i > 0; i--) {
838 En->delta[i] = En->delta[i-1];
839 En->ATE[i] = En->ATE[i-1];
840 }
841 // store new one
842 En->delta[0] = delta;
843 En->ATE[0] = E;
844 if (En->TotalEnergy[1] != 0.)
845 Diff = fabs(En->TotalEnergy[1] - E0)/(En->TotalEnergy[1] - E0)*fabs((E0 - En->ATE[1])/E0);
846 else
847 Diff = 0.;
848 R->Diffcount += pow(Diff,2);
849
850 // reinstate actual density (only needed for UpdateDensityCalculation) ...
851 //CalculateOneDensityR(Lat, LevS, Dens, source, Dens->DensityArray[ActualDensity], R->FactorDensityR*Psi->LocalPsiStatus[R->ActualLocalPsiNo].PsiFactor, 1);
852 // ... before changing actual local Psi
853 for (g = 0; g < LevS->MaxG; g++) { // Here all coefficients are updated for the new found wave function
854 //if (isnan(ConDir[g].re)) { fprintf(stderr,"WARNGING: CalculateLineSearch(): ConDir_%i(%i) = NaN!\n", R->ActualLocalPsiNo, g); Error(SomeError, "NaN-Fehler!"); }
855 //if (source != LevS->LPsi->LocalPsi[R->ActualLocalPsiNo] || g<0 || g>=LevS->MaxG) Error(SomeError,"FindPerturbedMinimum: source corrupted");
856 ////c_re(source[g]) = c_re(source[g])*dcos + c_re(ConDir[g])*dsin;
857 ////c_im(source[g]) = c_im(source[g])*dcos + c_im(ConDir[g])*dsin;
858 c_re(source[g]) = c_re(source[g]) + c_re(ConDir[g])*delta;
859 c_im(source[g]) = c_im(source[g]) + c_im(ConDir[g])*delta;
860 }
861 if (P->Call.out[StepLeaderOut]) {
862 fprintf(stderr, "(%i,%i,%i)S(%i,%i,%i):\tTE: %e\tATE: %e\t Diff: %e\t --- d: %e\tdEdt0: %e\tddEddt0: %e\n",P->Par.my_color_comm_ST,P->Par.me_comm_ST, P->Par.me_comm_ST_PsiT, R->MinStep, R->ActualLocalPsiNo, R->PsiStep, E0, E, Diff,delta, dEdt0, ddEddt0);
863 //fprintf(stderr, "(%i,%i,%i)S(%i,%i,%i):\tp0: %e p1: %e p2: %e \tATE: %e\t Diff: %e\t --- d: %e\tdEdt0: %e\tddEddt0: %e\n",P->Par.my_color_comm_ST,P->Par.me_comm_ST, P->Par.me_comm_ST_PsiT, R->MinStep, R->ActualLocalPsiNo, R->PsiStep, En->parts[0], En->parts[1], En->parts[2], E, Diff,delta, dEdt0, ddEddt0);
864 }
865 if (P->Par.me == 0) {
866 fprintf(F->MinimisationFile, "%i\t%i\t%i\t%e\t%e\t%e\t%e\t%e\n",R->MinStep, R->ActualLocalPsiNo, R->PsiStep, E0, E, delta, dEdt0, ddEddt0);
867 fflush(F->MinimisationFile);
868 }
869}
870
871/** Applies perturbation operator \f$\nabla_{index}\f$ to \a *source.
872 * As wave functions are stored in the reciprocal basis set, the application is straight-forward,
873 * for every G vector, the by \a index specified component is multiplied with the respective
874 * coefficient. Afterwards, 1/i is applied by flipping real and imaginary components (and an additional minus sign on the new imaginary term).
875 * \param *P Problem at hand
876 * \param *source complex coefficients of wave function \f$\varphi(G)\f$
877 * \param *dest returned complex coefficients of wave function \f$\widehat{p}_{index}|\varphi(G)\f$
878 * \param index_g vectorial index of operator to be applied
879 */
880void CalculatePerturbationOperator_P(struct Problem *P, const fftw_complex *source, fftw_complex *dest, const int index_g)
881{
882 struct RunStruct *R = &P->R;
883 struct LatticeLevel *LevS = R->LevS;
884 //const fftw_complex *dest_bak = dest;
885 int g = 0;
886 if (LevS->GArray[0].GSq == 0.0) {
887 //if (dest != dest_bak) Error(SomeError,"CalculatePerturbationOperator_P: dest corrupted");
888 dest[0].re = LevS->GArray[0].G[index_g]*source[0].im;
889 dest[0].im = -LevS->GArray[0].G[index_g]*source[0].re;
890 g++;
891 }
892 for (;g<LevS->MaxG;g++) {
893 //if (dest != dest_bak || g<0 || g>=LevS->MaxG) Error(SomeError,"CalculatePerturbationOperator_P: g out of range");
894 dest[g].re = LevS->GArray[g].G[index_g]*source[g].im;
895 dest[g].im = -LevS->GArray[g].G[index_g]*source[g].re;
896 }
897 // don't put dest[0].im = 0! Otherwise real parts of perturbed01/10 are not the same anymore!
898}
899
900/** Applies perturbation operator \f$\widehat{r}_{index}\f$ to \a *source.
901 * The \a *source wave function is blown up onto upper level LatticeLevel RunStruct#Lev0, fourier
902 * transformed. Afterwards, for each point on the real mesh the coefficient is multiplied times the real
903 * vector pointing within the cell to the mesh point, yet on LatticeLevel RunStruct#LevS. The new wave
904 * function is inverse fourier transformed and the resulting reciprocal coefficients stored in *dest.
905 * \param *P Problem at hand
906 * \param *source source coefficients
907 * \param *source2 second source coefficients, e.g. in the evaluation of a scalar product
908 * \param *dest destination coefficienta array, is overwrittten!
909 * \param index_r index of real vector.
910 * \param wavenr index of respective PsiTypeTag#Occupied(!) OnePsiElementAddData for the needed Wanner centre of the wave function.
911 */
912void CalculatePerturbationOperator_R(struct Problem *P, const fftw_complex *source, fftw_complex *dest, const fftw_complex *source2, const int index_r, const int wavenr)
913{
914 struct Lattice *Lat = &P->Lat;
915 struct RunStruct *R = &P->R;
916 struct LatticeLevel *Lev0 = R->Lev0;
917 struct LatticeLevel *LevS = R->LevS;
918 struct Density *Dens0 = Lev0->Dens;
919 struct fft_plan_3d *plan = Lat->plan;
920 fftw_complex *TempPsi = Dens0->DensityCArray[Temp2Density];
921 fftw_real *TempPsiR = (fftw_real *) TempPsi;
922 fftw_complex *workC = Dens0->DensityCArray[TempDensity];
923 fftw_complex *PsiC = Dens0->DensityCArray[ActualPsiDensity];
924 fftw_real *PsiCR = (fftw_real *) PsiC;
925 fftw_complex *tempdestRC = (fftw_complex *)Dens0->DensityArray[TempDensity];
926 fftw_complex *posfac, *destsnd, *destrcv;
927 double x[NDIM], X[NDIM], fac[NDIM], Wcentre[NDIM];
928 const int k_normal = Lat->Psi.TypeStartIndex[Occupied] + (wavenr - Lat->Psi.TypeStartIndex[R->CurrentMin]);
929 int n[NDIM], n0, g, Index, pos, iS, i0;
930 int N[NDIM], NUp[NDIM];
931 const int N0 = LevS->Plan0.plan->local_nx;
932 N[0] = LevS->Plan0.plan->N[0];
933 N[1] = LevS->Plan0.plan->N[1];
934 N[2] = LevS->Plan0.plan->N[2];
935 NUp[0] = LevS->NUp[0];
936 NUp[1] = LevS->NUp[1];
937 NUp[2] = LevS->NUp[2];
938 Wcentre[0] = Lat->Psi.AddData[k_normal].WannierCentre[0];
939 Wcentre[1] = Lat->Psi.AddData[k_normal].WannierCentre[1];
940 Wcentre[2] = Lat->Psi.AddData[k_normal].WannierCentre[2];
941 // init pointers and values
942 const int myPE = P->Par.me_comm_ST_Psi;
943 const double FFTFactor = 1./LevS->MaxN;
944 double vector;
945 //double result, Result;
946
947 // blow up source coefficients
948 LockDensityArray(Dens0,TempDensity,real); // tempdestRC
949 LockDensityArray(Dens0,Temp2Density,imag); // TempPsi
950 LockDensityArray(Dens0,ActualPsiDensity,imag); // PsiC
951 //if (tempdestRC != (fftw_complex *)Dens0->DensityArray[TempDensity]) Error(SomeError,"CalculatePerturbationOperator_R: tempdestRC corrupted");
952 SetArrayToDouble0((double *)tempdestRC ,Dens0->TotalSize*2);
953 //if (TempPsi != Dens0->DensityCArray[Temp2Density]) Error(SomeError,"CalculatePerturbationOperator_R: TempPsi corrupted");
954 SetArrayToDouble0((double *)TempPsi ,Dens0->TotalSize*2);
955 //if (PsiC != Dens0->DensityCArray[ActualPsiDensity]) Error(SomeError,"CalculatePerturbationOperator_R: PsiC corrupted");
956 SetArrayToDouble0((double *)PsiC,Dens0->TotalSize*2);
957 for (g=0; g<LevS->MaxG; g++) {
958 Index = LevS->GArray[g].Index;
959 posfac = &LevS->PosFactorUp[LevS->MaxNUp*g];
960 destrcv = &tempdestRC[LevS->MaxNUp*Index];
961 for (pos=0; pos < LevS->MaxNUp; pos++) {
962 //if (destrcv != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->LocalSizeC) Error(SomeError,"CalculatePerturbationOperator_R: destrcv corrupted");
963 destrcv [pos].re = (( source[g].re)*posfac[pos].re-(source[g].im)*posfac[pos].im);
964 destrcv [pos].im = (( source[g].re)*posfac[pos].im+(source[g].im)*posfac[pos].re);
965 }
966 }
967 for (g=0; g<LevS->MaxDoubleG; g++) {
968 destsnd = &tempdestRC [LevS->DoubleG[2*g]*LevS->MaxNUp];
969 destrcv = &tempdestRC [LevS->DoubleG[2*g+1]*LevS->MaxNUp];
970 for (pos=0; pos<LevS->MaxNUp; pos++) {
971 //if (destrcv != &tempdestRC [LevS->DoubleG[2*g+1]*LevS->MaxNUp] || LevS->DoubleG[2*g]*LevS->MaxNUp+pos<0 || LevS->DoubleG[2*g]*LevS->MaxNUp+pos>=Dens0->LocalSizeC|| LevS->DoubleG[2*g+1]*LevS->MaxNUp+pos<0 || LevS->DoubleG[2*g+1]*LevS->MaxNUp+pos>=Dens0->LocalSizeC) Error(SomeError,"CalculatePerturbationOperator_R: destrcv corrupted");
972 destrcv [pos].re = destsnd [pos].re;
973 destrcv [pos].im = -destsnd [pos].im;
974 }
975 }
976 // fourier transform blown up wave function
977 //if (tempdestRC != (fftw_complex *)Dens0->DensityArray[TempDensity]) Error(SomeError,"CalculatePerturbationOperator_R: tempdestRC corrupted");
978 //if (workC != Dens0->DensityCArray[TempDensity]) Error(SomeError,"CalculatePerturbationOperator_R: workC corrupted");
979 fft_3d_complex_to_real(plan,LevS->LevelNo, FFTNFUp, tempdestRC , workC);
980 //if (tempdestRC != (fftw_complex *)Dens0->DensityArray[TempDensity]) Error(SomeError,"CalculatePerturbationOperator_R: tempdestRC corrupted");
981 //if (TempPsiR != (fftw_real *)Dens0->DensityCArray[Temp2Density]) Error(SomeError,"CalculatePerturbationOperator_R: TempPsiR corrupted");
982 DensityRTransformPos(LevS,(fftw_real*)tempdestRC ,TempPsiR );
983 UnLockDensityArray(Dens0,TempDensity,real); // TempdestRC
984
985 //result = 0.;
986 // for every point on the real grid multiply with component of position vector
987 for (n0=0; n0<N0; n0++)
988 for (n[1]=0; n[1]<N[1]; n[1]++)
989 for (n[2]=0; n[2]<N[2]; n[2]++) {
990 n[0] = n0 + N0 * myPE;
991 fac[0] = (double)(n[0])/(double)((N[0]));
992 fac[1] = (double)(n[1])/(double)((N[1]));
993 fac[2] = (double)(n[2])/(double)((N[2]));
994 RMat33Vec3(x,Lat->RealBasis,fac);
995 iS = n[2] + N[2]*(n[1] + N[1]*n0); // mind splitting of x axis due to multiple processes
996 i0 = n[2]*NUp[2]+N[2]*NUp[2]*(n[1]*NUp[1]+N[1]*NUp[1]*n0*NUp[0]);
997 //PsiCR[iS] = ((double)n[0]/(double)N[0]*Lat->RealBasis[0] - fabs(Wcentre[0]))*TempPsiR[i0] - ((double)n[1]/(double)N[1]*Lat->RealBasis[4] - fabs(Wcentre[1]))*TempPsi2R[i0];
998 //fprintf(stderr,"(%i) R[%i] = (%lg,%lg,%lg)\n",P->Par.me, i0, x[0], x[1], x[2]);
999 //else fprintf(stderr,"(%i) WCentre[%i] = %e \n",P->Par.me, index_r, Wcentre[index_r]);
1000 MinImageConv(Lat,x, Wcentre, X);
1001 vector = sawtooth(Lat,X[index_r],index_r);
1002 //vector = 1.;//sin((double)(n[index_r])/(double)((N[index_r]))*2*PI);
1003 PsiCR[iS] = vector * TempPsiR[i0];
1004 //fprintf(stderr,"(%i) vector(%i/%i,%i/%i,%i/%i): %lg\tx[%i] = %e\tWcentre[%i] = %e\tTempPsiR[%i] = %e\tPsiCR[%i] = %e\n",P->Par.me, n[0], N[0], n[1], N[1], n[2], N[2], vector, index_r, x[index_r],index_r, Wcentre[index_r],i0,TempPsiR[i0],iS,PsiCR[iS]);
1005
1006 //truedist(Lat,x[cross(index_r,2)],Wcentre[cross(index_r,2)],cross(index_r,2)) * TempPsiR[i0];
1007 //tmp += truedist(Lat,x[index_r],WCentre[index_r],index_r) * RealPhiR[i0];
1008 //tmp += sawtooth(Lat,truedist(Lat,x[index_r],WCentre[index_r],index_r), index_r)*RealPhiR[i0];
1009 //(Fehler mit falschem Ort ist vor dieser Stelle!): ueber result = RealPhiR[i0] * (x[index_r]) * RealPhiR[i0]; gecheckt
1010 //result += TempPsiR[i0] * PsiCR[iS];
1011 }
1012 UnLockDensityArray(Dens0,Temp2Density,imag); // TempPsi
1013 //MPI_Allreduce( &result, &Result, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
1014 //if (P->Par.me == 0) fprintf(stderr,"(%i) PerturbationOpertator_R: %e\n",P->Par.me, Result/LevS->MaxN);
1015 // inverse fourier transform
1016 fft_3d_real_to_complex(plan,LevS->LevelNo, FFTNF1, PsiC, workC);
1017 //fft_3d_real_to_complex(plan,LevS->LevelNo, FFTNF1, Psi2C, workC);
1018
1019 // copy to destination array
1020 for (g=0; g<LevS->MaxG; g++) {
1021 Index = LevS->GArray[g].Index;
1022 dest[g].re = ( PsiC[Index].re)*FFTFactor;
1023 dest[g].im = ( PsiC[Index].im)*FFTFactor;
1024 }
1025 UnLockDensityArray(Dens0,ActualPsiDensity,imag); //PsiC
1026 //if (LevS->GArray[0].GSq == 0)
1027 // dest[0].im = 0; // imaginary of G=0 is zero
1028}
1029/*
1030{
1031 struct RunStruct *R = &P->R;
1032 struct LatticeLevel *Lev0 = R->Lev0;
1033 struct LatticeLevel *LevS = R->LevS;
1034 struct Lattice *Lat = &P->Lat;
1035 struct fft_plan_3d *plan = Lat->plan;
1036 struct Density *Dens0 = Lev0->Dens;
1037 fftw_complex *tempdestRC = Dens0->DensityCArray[TempDensity];
1038 fftw_real *tempdestR = (fftw_real *) tempdestRC;
1039 fftw_complex *work = Dens0->DensityCArray[Temp2Density];
1040 fftw_complex *PsiC = (fftw_complex *) Dens0->DensityCArray[ActualPsiDensity];;
1041 fftw_real *PsiCR = (fftw_real *) PsiC;
1042 fftw_real *RealPhiR = (fftw_real *) Dens0->DensityArray[Temp2Density];
1043 fftw_complex *posfac, *destsnd, *destrcv;
1044 double x[NDIM], fac[NDIM], WCentre[NDIM];
1045 int n[NDIM], N0, n0, g, Index, pos, iS, i0;
1046
1047 // init pointers and values
1048 int myPE = P->Par.me_comm_ST_Psi;
1049 double FFTFactor = 1./LevS->MaxN;
1050 int N[NDIM], NUp[NDIM];
1051 N[0] = LevS->Plan0.plan->N[0];
1052 N[1] = LevS->Plan0.plan->N[1];
1053 N[2] = LevS->Plan0.plan->N[2];
1054 NUp[0] = LevS->NUp[0];
1055 NUp[1] = LevS->NUp[1];
1056 NUp[2] = LevS->NUp[2];
1057 N0 = LevS->Plan0.plan->local_nx;
1058 wavenr = Lat->Psi.TypeStartIndex[Occupied] + (wavenr - Lat->Psi.TypeStartIndex[R->CurrentMin]);
1059 Wcentre[0] = Lat->Psi.AddData[wavenr].WannierCentre[0];
1060 Wcentre[1] = Lat->Psi.AddData[wavenr].WannierCentre[1];
1061 Wcentre[2] = Lat->Psi.AddData[wavenr].WannierCentre[2];
1062
1063 // blow up source coefficients
1064 SetArrayToDouble0((double *)tempdestRC,Dens0->TotalSize*2);
1065 SetArrayToDouble0((double *)RealPhiR,Dens0->TotalSize*2);
1066 SetArrayToDouble0((double *)PsiC,Dens0->TotalSize*2);
1067 for (g=0; g<LevS->MaxG; g++) {
1068 Index = LevS->GArray[g].Index;
1069 posfac = &LevS->PosFactorUp[LevS->MaxNUp*g];
1070 destrcv = &tempdestRC[LevS->MaxNUp*Index];
1071 for (pos=0; pos<LevS->MaxNUp; pos++) {
1072 destrcv[pos].re = (( source[g].re)*posfac[pos].re-( source[g].im)*posfac[pos].im);
1073 destrcv[pos].im = (( source[g].re)*posfac[pos].im+( source[g].im)*posfac[pos].re);
1074 }
1075 }
1076 for (g=0; g<LevS->MaxDoubleG; g++) {
1077 destsnd = &tempdestRC[LevS->DoubleG[2*g]*LevS->MaxNUp];
1078 destrcv = &tempdestRC[LevS->DoubleG[2*g+1]*LevS->MaxNUp];
1079 for (pos=0; pos<LevS->MaxNUp; pos++) {
1080 destrcv[pos].re = destsnd[pos].re;
1081 destrcv[pos].im = -destsnd[pos].im;
1082 }
1083 }
1084
1085 // fourier transform blown up wave function
1086 fft_3d_complex_to_real(plan,LevS->LevelNo, FFTNFUp, tempdestRC, work);
1087 DensityRTransformPos(LevS,tempdestR,RealPhiR);
1088
1089 //fft_Psi(P,source,RealPhiR,0,0);
1090
1091 // for every point on the real grid multiply with component of position vector
1092 for (n0=0; n0<N0; n0++)
1093 for (n[1]=0; n[1]<N[1]; n[1]++)
1094 for (n[2]=0; n[2]<N[2]; n[2]++) {
1095 n[0] = n0 + N0 * myPE;
1096 fac[0] = (double)(n[0])/(double)((N[0]));
1097 fac[1] = (double)(n[1])/(double)((N[1]));
1098 fac[2] = (double)(n[2])/(double)((N[2]));
1099 RMat33Vec3(x,Lat->RealBasis,fac);
1100 iS = n[2] + N[2]*(n[1] + N[1]*n0); // mind splitting of x axis due to multiple processes
1101 i0 = n[2]*NUp[2]+N[2]*NUp[2]*(n[1]*NUp[1]+N[1]*NUp[1]*n0*NUp[0]);
1102 //PsiCR[iS] = (x[index_r]) * RealPhiR[i0]; //- WCentre[index_r]
1103 PsiCR[iS] = truedist(Lat,x[index_r],WCentre[index_r],index_r) * RealPhiR[i0];
1104 //PsiCR[iS] = truedist(Lat,x[index_r],0.,index_r) * RealPhiR[i0];
1105 //PsiCR[iS] = sawtooth(Lat,truedist(Lat,x[index_r],WCentre[index_r],index_r), index_r)*RealPhiR[i0];
1106 //(Fehler mit falschem Ort ist vor dieser Stelle!): ueber result = RealPhiR[i0] * (x[index_r]) * RealPhiR[i0]; gecheckt
1107 }
1108
1109 // inverse fourier transform
1110 fft_3d_real_to_complex(plan,LevS->LevelNo, FFTNF1, PsiC, work);
1111
1112 // copy to destination array
1113 for (g=0; g<LevS->MaxG; g++) {
1114 Index = LevS->GArray[g].Index;
1115 dest[g].re = ( PsiC[Index].re)*FFTFactor;
1116 dest[g].im = ( PsiC[Index].im)*FFTFactor;
1117 if (LevS->GArray[g].GSq == 0)
1118 dest[g].im = 0; // imaginary of G=0 is zero
1119 }
1120}*/
1121
1122/** Prints the positions of all unperturbed orbitals to screen.
1123 * \param *P Problem at hand
1124 * \param type PsiTypeTag specifying group of orbitals
1125 * \sa CalculatePerturbationOperator_R()
1126 */
1127void OutputOrbitalPositions(struct Problem *P, const enum PsiTypeTag type)
1128{
1129 struct Lattice *Lat = &P->Lat;
1130 struct Psis *Psi = &Lat->Psi;
1131 struct RunStruct *R = &P->R;
1132 struct LatticeLevel *LevS = R->LevS;
1133 fftw_complex *temp = LevS->LPsi->TempPsi;
1134 fftw_complex *source;
1135 int wavenr, index;
1136 double result[NDIM], Result[NDIM];
1137 //double imsult[NDIM], Imsult[NDIM];
1138 double norm[NDIM], Norm[NDIM];
1139 //double imnorm[NDIM], imNorm[NDIM];
1140 double Wcentre[NDIM];
1141
1142 // for every unperturbed wave function
1143 for (wavenr=Psi->TypeStartIndex[type]; wavenr<Psi->TypeStartIndex[type+1]; wavenr++) {
1144 source = LevS->LPsi->LocalPsi[wavenr];
1145 Wcentre[0] = Psi->AddData[wavenr].WannierCentre[0];
1146 Wcentre[1] = Psi->AddData[wavenr].WannierCentre[1];
1147 Wcentre[2] = Psi->AddData[wavenr].WannierCentre[2];
1148 for (index=0; index<NDIM; index++) {
1149 SetArrayToDouble0((double *)temp,2*R->InitLevS->MaxG);
1150 // apply position operator
1151 CalculatePerturbationOperator_R(P,source,temp,source,index, wavenr + Psi->TypeStartIndex[R->CurrentMin]);
1152 // take scalar product
1153 result[index] = GradSP(P,LevS,source,temp);
1154 //imsult[index] = GradImSP(P,LevS,source,temp);
1155 norm[index] = GradSP(P,LevS,source,source);
1156 //imnorm[index] = GradImSP(P,LevS,source,source);
1157 MPI_Allreduce( result, Result, NDIM, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
1158 //MPI_Allreduce( imsult, Imsult, NDIM, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
1159 MPI_Allreduce( norm, Norm, NDIM, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
1160 //MPI_Allreduce( imnorm, imNorm, NDIM, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
1161 }
1162 // print output to stderr
1163 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) Position of Orbital %i: (%e,%e,%e)\n",P->Par.me, wavenr, Result[0]/Norm[0]+Wcentre[0], Result[1]/Norm[1]+Wcentre[1], Result[2]/Norm[2]+Wcentre[2]);
1164 //fprintf(stderr,"(%i) Position of Orbital %i wrt Wannier: (%e,%e,%e)\n",P->Par.me, wavenr, Result[0]/Norm[0], Result[1]/Norm[1], Result[2]/Norm[2]);
1165 //fprintf(stderr,"(%i) with Norm: (%e,%e,%e) + i (%e,%e,%e)\n",P->Par.me, Norm[0], Norm[1], Norm[2], imNorm[0], imNorm[1], imNorm[2]);
1166 //if (P->Par.me == 0) fprintf(stderr,"(%i) Position of Orbital %i: (%e,%e,%e)\n",P->Par.me, wavenr, Result[0]/Norm[0], Result[1]/Norm[1], Result[2]/Norm[2]);
1167 }
1168}
1169
1170#define borderstart 0.9
1171
1172/** Applies perturbation operator \f$(\widehat{r} \times \nabla)_{index}\f$ to \a *source.
1173 * The source is fourier-transformed by transforming it to a density (on the next higher level RunStruct#Lev0)
1174 * and at the same time multiply it with the respective component of the reciprocal G vector - the momentum. This
1175 * is done by callinf fft_Psi(). Thus we get \f$\nabla_k | \varphi (R) \rangle\f$.
1176 *
1177 * Next, we apply the two of three components of the position operator r, which ones stated by cross(), while going
1178 * in a loop through every point of the grid. In order to do this sensibly, truedist() is used to map the coordinates
1179 * onto -L/2...L/2, by subtracting the OneElementPsiAddData#WannierCentre R and wrapping. Also, due to the breaking up
1180 * of the x axis into equally sized chunks for each coefficient sharing process, we need to step only over local
1181 * x-axis grid points, however shift them to the global position when being used as position. In the end, we get
1182 * \f$\epsilon_{index,j,k} (\widehat{r}-R)_j \nabla_k | \varphi (R) \rangle\f$.
1183 *
1184 * One last fft brings the wave function back to reciprocal basis and it is copied to \a *dest.
1185 * \param *P Problem at hand
1186 * \param *source complex coefficients of wave function \f$\varphi(G)\f$
1187 * \param *dest returned complex coefficients of wave function \f$(\widehat{r} \times \widehat{p})_{index}|\varphi(G)\rangle\f$
1188 * \param phi0nr number within LocalPsi of the unperturbed pendant of the given perturbed wavefunction \a *source.
1189 * \param index_rxp index desired of the vector product
1190 * \sa CalculateConDirHConDir() - the procedure of fft and inverse fft is very similar.
1191 */
1192void CalculatePerturbationOperator_RxP(struct Problem *P, const fftw_complex *source, fftw_complex *dest, const int phi0nr, const int index_rxp)
1193
1194{
1195 struct Lattice *Lat = &P->Lat;
1196 struct RunStruct *R = &P->R;
1197 struct LatticeLevel *Lev0 = R->Lev0;
1198 struct LatticeLevel *LevS = R->LevS;
1199 struct Density *Dens0 = Lev0->Dens;
1200 struct fft_plan_3d *plan = Lat->plan;
1201 fftw_complex *TempPsi = Dens0->DensityCArray[Temp2Density];
1202 fftw_real *TempPsiR = (fftw_real *) TempPsi;
1203 fftw_complex *TempPsi2 = (fftw_complex *)Dens0->DensityArray[Temp2Density];
1204 fftw_real *TempPsi2R = (fftw_real *) TempPsi2;
1205 fftw_complex *workC = Dens0->DensityCArray[TempDensity];
1206 fftw_complex *PsiC = Dens0->DensityCArray[ActualPsiDensity];
1207 fftw_real *PsiCR = (fftw_real *) PsiC;
1208 double x[NDIM], X[NDIM], fac[NDIM], Wcentre[NDIM];
1209 int n[NDIM], n0, g, Index, iS, i0; //pos,
1210 int N[NDIM], NUp[NDIM];
1211 const int N0 = LevS->Plan0.plan->local_nx;
1212 N[0] = LevS->Plan0.plan->N[0];
1213 N[1] = LevS->Plan0.plan->N[1];
1214 N[2] = LevS->Plan0.plan->N[2];
1215 NUp[0] = LevS->NUp[0];
1216 NUp[1] = LevS->NUp[1];
1217 NUp[2] = LevS->NUp[2];
1218 Wcentre[0] = Lat->Psi.AddData[phi0nr].WannierCentre[0];
1219 Wcentre[1] = Lat->Psi.AddData[phi0nr].WannierCentre[1];
1220 Wcentre[2] = Lat->Psi.AddData[phi0nr].WannierCentre[2];
1221 // init pointers and values
1222 const int myPE = P->Par.me_comm_ST_Psi;
1223 const double FFTFactor = 1./LevS->MaxN; //
1224// double max[NDIM], max_psi[NDIM];
1225// double max_n[NDIM];
1226 int index[4];
1227// double smooth, wall[NDIM];
1228// for (g=0;g<NDIM;g++) {
1229// max[g] = 0.;
1230// max_psi[g] = 0.;
1231// max_n[g] = -1.;
1232// }
1233
1234 //fprintf(stderr,"(%i) Wannier[%i] (%2.13e, %2.13e, %2.13e)\n", P->Par.me, phi0nr, 10.-Wcentre[0], 10.-Wcentre[1], 10.-Wcentre[2]);
1235 for (g=0;g<4;g++)
1236 index[g] = cross(index_rxp,g);
1237
1238 // blow up source coefficients
1239 LockDensityArray(Dens0,Temp2Density,imag); // TempPsi
1240 LockDensityArray(Dens0,Temp2Density,real); // TempPsi2
1241 LockDensityArray(Dens0,ActualPsiDensity,imag); // PsiC
1242
1243 fft_Psi(P,source,TempPsiR ,index[1],7);
1244 fft_Psi(P,source,TempPsi2R,index[3],7);
1245
1246 //result = 0.;
1247 // for every point on the real grid multiply with component of position vector
1248 for (n0=0; n0<N0; n0++)
1249 for (n[1]=0; n[1]<N[1]; n[1]++)
1250 for (n[2]=0; n[2]<N[2]; n[2]++) {
1251 n[0] = n0 + N0 * myPE;
1252 fac[0] = (double)(n[0])/(double)((N[0]));
1253 fac[1] = (double)(n[1])/(double)((N[1]));
1254 fac[2] = (double)(n[2])/(double)((N[2]));
1255 RMat33Vec3(x,Lat->RealBasis,fac);
1256// fac[0] = (fac[0] > .9) ? fac[0]-0.9 : 0.;
1257// fac[1] = (fac[1] > .9) ? fac[1]-0.9 : 0.;
1258// fac[2] = (fac[2] > .9) ? fac[2]-0.9 : 0.;
1259// RMat33Vec3(wall,Lat->RealBasis,fac);
1260// smooth = exp(wall[0]*wall[0]+wall[1]*wall[1]+wall[2]*wall[2]); // smoothing near the borders of the virtual cell
1261 iS = n[2] + N[2]*(n[1] + N[1]*n0); // mind splitting of x axis due to multiple processes
1262 i0 = n[2]*NUp[2]+N[2]*NUp[2]*(n[1]*NUp[1]+N[1]*NUp[1]*n0*NUp[0]);
1263
1264// if (fabs(truedist(Lat,x[index[1]],Wcentre[index[1]],index[1])) >= borderstart * sqrt(Lat->RealBasisSQ[index[1]])/2.)
1265// if (max[index[1]] < sawtooth(Lat,truedist(Lat,x[index[1]],Wcentre[index[1]],index[1]),index[1]) * TempPsiR [i0]) {
1266// max[index[1]] = sawtooth(Lat,truedist(Lat,x[index[1]],Wcentre[index[1]],index[1]),index[1]) * TempPsiR [i0];
1267// max_psi[index[1]] = TempPsiR [i0];
1268// max_n[index[1]] = truedist(Lat,x[index[1]],Wcentre[index[1]],index[1]);
1269// }
1270//
1271// if (fabs(truedist(Lat,x[index[3]],Wcentre[index[3]],index[3])) >= borderstart * sqrt(Lat->RealBasisSQ[index[3]])/2.)
1272// if (max[index[3]] < sawtooth(Lat,truedist(Lat,x[index[3]],Wcentre[index[3]],index[3]),index[3]) * TempPsiR [i0]) {
1273// max[index[3]] = sawtooth(Lat,truedist(Lat,x[index[3]],Wcentre[index[3]],index[3]),index[3]) * TempPsiR [i0];
1274// max_psi[index[3]] = TempPsiR [i0];
1275// max_n[index[3]] = truedist(Lat,x[index[3]],Wcentre[index[3]],index[3]);
1276// }
1277
1278 MinImageConv(Lat, x, Wcentre, X);
1279 PsiCR[iS] = //vector * TempPsiR[i0];
1280 sawtooth(Lat,X[index[0]],index[0]) * TempPsiR [i0]
1281 -sawtooth(Lat,X[index[2]],index[2]) * TempPsi2R[i0];
1282// ShiftGaugeOrigin(P,X[index[0]],index[0]) * TempPsiR [i0]
1283// -ShiftGaugeOrigin(P,X[index[2]],index[2]) * TempPsi2R[i0];
1284// PsiCR[iS] = (x[index[0]] - Wcentre[index[0]]) * TempPsiR [i0] - (x[index[2]] - Wcentre[index[2]]) * TempPsi2R[i0];
1285 }
1286 //if (P->Par.me == 0) fprintf(stderr,"(%i) PerturbationOpertator_R(xP): %e\n",P->Par.me, Result/LevS->MaxN);
1287 UnLockDensityArray(Dens0,Temp2Density,imag); // TempPsi
1288 UnLockDensityArray(Dens0,Temp2Density,real); // TempPsi2
1289
1290// // print maximum values
1291// fprintf (stderr,"(%i) RxP: Maximum values = (",P->Par.me);
1292// for (g=0;g<NDIM;g++)
1293// fprintf(stderr,"%lg\t", max[g]);
1294// fprintf(stderr,"\b)\t(");
1295// for (g=0;g<NDIM;g++)
1296// fprintf(stderr,"%lg\t", max_psi[g]);
1297// fprintf(stderr,"\b)\t");
1298// fprintf (stderr,"at (");
1299// for (g=0;g<NDIM;g++)
1300// fprintf(stderr,"%lg\t", max_n[g]);
1301// fprintf(stderr,"\b)\n");
1302
1303 // inverse fourier transform
1304 //if (PsiC != Dens0->DensityCArray[ActualPsiDensity]) Error(SomeError,"CalculatePerturbationOperator_RxP: PsiC corrupted");
1305 fft_3d_real_to_complex(plan,LevS->LevelNo, FFTNF1, PsiC, workC);
1306
1307 // copy to destination array
1308 SetArrayToDouble0((double *)dest, 2*R->InitLevS->MaxG);
1309 for (g=0; g<LevS->MaxG; g++) {
1310 Index = LevS->GArray[g].Index;
1311 dest[g].re += ( PsiC[Index].re)*FFTFactor; // factor confirmed, see grad.c:CalculateConDirHConDir()
1312 dest[g].im += ( PsiC[Index].im)*FFTFactor;
1313 //fprintf(stderr,"(%i) PsiC[(%lg,%lg,%lg)] = %lg +i %lg\n", P->Par.me, LevS->GArray[g].G[0], LevS->GArray[g].G[1], LevS->GArray[g].G[2], dest[g].re, dest[g].im);
1314 }
1315 UnLockDensityArray(Dens0,ActualPsiDensity,imag); // PsiC
1316 //if (LevS->GArray[0].GSq == 0.)
1317 //dest[0].im = 0.; // don't do this, see ..._P()
1318}
1319
1320/** Applies perturbation operator \f$-(\nabla \times \widehat{r})_{index}\f$ to \a *source.
1321 * Is analogous to CalculatePerturbationOperator_RxP(), only the order is reversed, first position operator, then
1322 * momentum operator
1323 * \param *P Problem at hand
1324 * \param *source complex coefficients of wave function \f$\varphi(G)\f$
1325 * \param *dest returned complex coefficients of wave function \f$(\widehat{r} \times \widehat{p})_{index}|\varphi(G)\rangle\f$
1326 * \param phi0nr number within LocalPsi of the unperturbed pendant of the given perturbed wavefunction \a *source.
1327 * \param index_pxr index of position operator
1328 * \note Only third component is important due to initial rotiation of cell such that B field is aligned with z axis.
1329 * \sa CalculateConDirHConDir() - the procedure of fft and inverse fft is very similar.
1330 * \bug routine is not tested (but should work), as it offers no advantage over CalculatePerturbationOperator_RxP()
1331 */
1332void CalculatePerturbationOperator_PxR(struct Problem *P, const fftw_complex *source, fftw_complex *dest, const int phi0nr, const int index_pxr)
1333
1334{
1335 struct Lattice *Lat = &P->Lat;
1336 struct RunStruct *R = &P->R;
1337 struct LatticeLevel *Lev0 = R->Lev0;
1338 struct LatticeLevel *LevS = R->LevS;
1339 struct Density *Dens0 = Lev0->Dens;
1340 struct fft_plan_3d *plan = Lat->plan;
1341 fftw_complex *TempPsi = Dens0->DensityCArray[Temp2Density];
1342 fftw_real *TempPsiR = (fftw_real *) TempPsi;
1343 fftw_complex *workC = Dens0->DensityCArray[TempDensity];
1344 fftw_complex *PsiC = Dens0->DensityCArray[ActualPsiDensity];
1345 fftw_real *PsiCR = (fftw_real *) PsiC;
1346 fftw_complex *Psi2C = Dens0->DensityCArray[ActualDensity];
1347 fftw_real *Psi2CR = (fftw_real *) Psi2C;
1348 fftw_complex *tempdestRC = (fftw_complex *)Dens0->DensityArray[Temp2Density];
1349 fftw_complex *posfac, *destsnd, *destrcv;
1350 double x[NDIM], X[NDIM], fac[NDIM], Wcentre[NDIM];
1351 int n[NDIM], n0, g, Index, pos, iS, i0;
1352 int N[NDIM], NUp[NDIM];
1353 const int N0 = LevS->Plan0.plan->local_nx;
1354 N[0] = LevS->Plan0.plan->N[0];
1355 N[1] = LevS->Plan0.plan->N[1];
1356 N[2] = LevS->Plan0.plan->N[2];
1357 NUp[0] = LevS->NUp[0];
1358 NUp[1] = LevS->NUp[1];
1359 NUp[2] = LevS->NUp[2];
1360 Wcentre[0] = Lat->Psi.AddData[phi0nr].WannierCentre[0];
1361 Wcentre[1] = Lat->Psi.AddData[phi0nr].WannierCentre[1];
1362 Wcentre[2] = Lat->Psi.AddData[phi0nr].WannierCentre[2];
1363 // init pointers and values
1364 const int myPE = P->Par.me_comm_ST_Psi;
1365 const double FFTFactor = 1./LevS->MaxN;
1366
1367 // blow up source coefficients
1368 SetArrayToDouble0((double *)tempdestRC ,Dens0->TotalSize*2);
1369 SetArrayToDouble0((double *)TempPsi ,Dens0->TotalSize*2);
1370 SetArrayToDouble0((double *)PsiC,Dens0->TotalSize*2);
1371 SetArrayToDouble0((double *)Psi2C,Dens0->TotalSize*2);
1372 for (g=0; g<LevS->MaxG; g++) {
1373 Index = LevS->GArray[g].Index;
1374 posfac = &LevS->PosFactorUp[LevS->MaxNUp*g];
1375 destrcv = &tempdestRC[LevS->MaxNUp*Index];
1376 for (pos=0; pos < LevS->MaxNUp; pos++) {
1377 destrcv [pos].re = (( source[g].re)*posfac[pos].re-( source[g].im)*posfac[pos].im);
1378 destrcv [pos].im = (( source[g].re)*posfac[pos].im+( source[g].im)*posfac[pos].re);
1379 }
1380 }
1381 for (g=0; g<LevS->MaxDoubleG; g++) {
1382 destsnd = &tempdestRC [LevS->DoubleG[2*g]*LevS->MaxNUp];
1383 destrcv = &tempdestRC [LevS->DoubleG[2*g+1]*LevS->MaxNUp];
1384 for (pos=0; pos<LevS->MaxNUp; pos++) {
1385 destrcv [pos].re = destsnd [pos].re;
1386 destrcv [pos].im = -destsnd [pos].im;
1387 }
1388 }
1389 // fourier transform blown up wave function
1390 fft_3d_complex_to_real(plan,LevS->LevelNo, FFTNFUp, tempdestRC , workC);
1391 DensityRTransformPos(LevS,(fftw_real*)tempdestRC ,TempPsiR );
1392
1393 //fft_Psi(P,source,TempPsiR ,cross(index_pxr,1),7);
1394 //fft_Psi(P,source,TempPsi2R,cross(index_pxr,3),7);
1395
1396 //result = 0.;
1397 // for every point on the real grid multiply with component of position vector
1398 for (n0=0; n0<N0; n0++)
1399 for (n[1]=0; n[1]<N[1]; n[1]++)
1400 for (n[2]=0; n[2]<N[2]; n[2]++) {
1401 n[0] = n0 + N0 * myPE;
1402 fac[0] = (double)(n[0])/(double)((N[0]));
1403 fac[1] = (double)(n[1])/(double)((N[1]));
1404 fac[2] = (double)(n[2])/(double)((N[2]));
1405 RMat33Vec3(x,Lat->RealBasis,fac);
1406 iS = n[2] + N[2]*(n[1] + N[1]*n0); // mind splitting of x axis due to multiple processes
1407 i0 = n[2]*NUp[2]+N[2]*NUp[2]*(n[1]*NUp[1]+N[1]*NUp[1]*n0*NUp[0]);
1408// PsiCR[iS] = sawtooth(Lat,X[cross(index_pxr,1)],cross(index_pxr,1)) * TempPsiR[i0];
1409// Psi2CR[iS] = sawtooth(Lat,X[cross(index_pxr,3)],cross(index_pxr,3)) * TempPsiR[i0];
1410 MinImageConv(Lat,x,Wcentre,X);
1411 PsiCR[iS] = ShiftGaugeOrigin(P,X[cross(index_pxr,1)],cross(index_pxr,1)) * TempPsiR[i0];
1412 Psi2CR[iS] = ShiftGaugeOrigin(P,X[cross(index_pxr,3)],cross(index_pxr,3)) * TempPsiR[i0];
1413 }
1414
1415 // inverse fourier transform
1416 fft_3d_real_to_complex(plan,LevS->LevelNo, FFTNF1, PsiC, workC);
1417 fft_3d_real_to_complex(plan,LevS->LevelNo, FFTNF1, Psi2C, workC);
1418
1419 // copy to destination array
1420 for (g=0; g<LevS->MaxG; g++) {
1421 Index = LevS->GArray[g].Index;
1422 dest[g].re = -LevS->GArray[g].G[cross(index_pxr,0)]*( PsiC[Index].im)*FFTFactor;
1423 dest[g].im = -LevS->GArray[g].G[cross(index_pxr,0)]*(-PsiC[Index].re)*FFTFactor;
1424 dest[g].re -= -LevS->GArray[g].G[cross(index_pxr,2)]*( Psi2C[Index].im)*FFTFactor;
1425 dest[g].im -= -LevS->GArray[g].G[cross(index_pxr,2)]*(-Psi2C[Index].re)*FFTFactor;
1426 }
1427 if (LevS->GArray[0].GSq == 0.)
1428 dest[0].im = 0.; // don't do this, see ..._P()
1429}
1430
1431/** Evaluates first derivative of perturbed energy functional with respect to minimisation parameter \f$\Theta\f$.
1432 * \f[
1433 * \frac{\delta {\cal E}^{(2)}} {\delta \Theta} =
1434 * 2 {\cal R} \langle \widetilde{\varphi}_i^{(1)} | {\cal H}^{(0)} | \varphi_i^{(1)} \rangle
1435 * - \sum_l \lambda_{il} \langle \widetilde{\varphi}_i^{(1)} | \varphi_l^{(1)} \rangle
1436 * - \sum_k \lambda_{ki} \langle \varphi_k^{(1)} | \widetilde{\varphi}_i^{(1)} \rangle
1437 * + 2 {\cal R} \langle \widetilde{\varphi}_i^{(1)} | {\cal H}^{(1)} | \varphi_i^{(0)} \rangle
1438 * \f]
1439 *
1440 * The summation over all Psis has again to be done with an MPI exchange of non-local coefficients, as the conjugate
1441 * directions are not the same in situations where PePGamma > 1 (Psis split up among processes = multiple minimisation)
1442 * \param *P Problem at hand
1443 * \param source0 unperturbed wave function \f$\varphi_l^{(0)}\f$
1444 * \param source perturbed wave function \f$\varphi_l^{(1)} (G)\f$
1445 * \param ConDir normalized conjugate direction \f$\widetilde{\varphi}_l^{(1)} (G)\f$
1446 * \param Hc_grad complex coefficients of \f$H^{(0)} | \varphi_l^{(1)} (G) \rangle\f$, see GradientArray#HcGradient
1447 * \param H1c_grad complex coefficients of \f$H^{(1)} | \varphi_l^{(0)} (G) \rangle\f$, see GradientArray#H1cGradient
1448 * \sa CalculateLineSearch() - used there, \sa CalculateConDirHConDir() - same principles
1449 * \warning The MPI_Allreduce for the scalar product in the end has not been done and must not have been done for given
1450 * parameters yet!
1451 */
1452double Calculate1stPerturbedDerivative(struct Problem *P, const fftw_complex *source0, const fftw_complex *source, const fftw_complex *ConDir, const fftw_complex *Hc_grad, const fftw_complex *H1c_grad)
1453{
1454 struct RunStruct *R = &P->R;
1455 struct Psis *Psi = &P->Lat.Psi;
1456 struct LatticeLevel *LevS = R->LevS;
1457 double result = 0., E0 = 0., Elambda = 0., E1 = 0.;//, E2 = 0.;
1458 int i,m,j;
1459 const int state = R->CurrentMin;
1460 //const int l_normal = R->ActualLocalPsiNo - Psi->TypeStartIndex[state] + Psi->TypeStartIndex[Occupied];
1461 const int ActNum = R->ActualLocalPsiNo - Psi->TypeStartIndex[state] + Psi->TypeStartIndex[1] * Psi->LocalPsiStatus[R->ActualLocalPsiNo].my_color_comm_ST_Psi;
1462 //int l = R->ActualLocalPsiNo;
1463 //int l_normal = Psi->TypeStartIndex[Occupied] + (l - Psi->TypeStartIndex[state]); // offset l to \varphi_l^{(0)}
1464 struct OnePsiElement *OnePsiB, *LOnePsiB;
1465 //fftw_complex *HConGrad = LevS->LPsi->TempPsi;
1466 fftw_complex *LPsiDatB=NULL;
1467 const int ElementSize = (sizeof(fftw_complex) / sizeof(double));
1468 int RecvSource;
1469 MPI_Status status;
1470
1471 //CalculateCDfnl(P,ConDir,PP->CDfnl);
1472 //ApplyTotalHamiltonian(P,ConDir,HConDir, PP->CDfnl, 1, 0);
1473 //E0 = (GradSP(P, LevS, ConDir, Hc_grad) + GradSP(P, LevS, source, HConDir)) * Psi->LocalPsiStatus[R->ActualLocalPsiNo].PsiFactor;
1474 E0 = 2.*GradSP(P, LevS, ConDir, Hc_grad) * Psi->LocalPsiStatus[R->ActualLocalPsiNo].PsiFactor;
1475 result = E0;
1476 //fprintf(stderr,"(%i) 1st: E0 = \t\t%lg\n", P->Par.me, E0);
1477
1478 m = -1;
1479 for (j=0; j < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; j++) { // go through all wave functions
1480 OnePsiB = &Psi->AllPsiStatus[j]; // grab OnePsiB
1481 if (OnePsiB->PsiType == state) { // drop all but the ones of current min state
1482 m++; // increase m if it is type-specific wave function
1483 if (OnePsiB->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // local?
1484 LOnePsiB = &Psi->LocalPsiStatus[OnePsiB->MyLocalNo];
1485 else
1486 LOnePsiB = NULL;
1487 if (LOnePsiB == NULL) { // if it's not local ... receive it from respective process into TempPsi
1488 RecvSource = OnePsiB->my_color_comm_ST_Psi;
1489 MPI_Recv( LevS->LPsi->TempPsi, LevS->MaxG*ElementSize, MPI_DOUBLE, RecvSource, PerturbedTag, P->Par.comm_ST_PsiT, &status );
1490 LPsiDatB=LevS->LPsi->TempPsi;
1491 } else { // .. otherwise send it to all other processes (Max_me... - 1)
1492 for (i=0;i<P->Par.Max_me_comm_ST_PsiT;i++)
1493 if (i != OnePsiB->my_color_comm_ST_Psi)
1494 MPI_Send( LevS->LPsi->LocalPsi[OnePsiB->MyLocalNo], LevS->MaxG*ElementSize, MPI_DOUBLE, i, PerturbedTag, P->Par.comm_ST_PsiT);
1495 LPsiDatB=LevS->LPsi->LocalPsi[OnePsiB->MyLocalNo];
1496 } // LPsiDatB is now set to the coefficients of OnePsi either stored or MPI_Received
1497
1498 Elambda -= 2.*Psi->lambda[ActNum][m]*GradSP(P, LevS, ConDir, LPsiDatB) * OnePsiB->PsiFactor; // lambda is symmetric
1499 }
1500 }
1501 result += Elambda;
1502 //fprintf(stderr,"(%i) 1st: Elambda = \t%lg\n", P->Par.me, Elambda);
1503
1504 E1 = 2.*GradSP(P,LevS,ConDir,H1c_grad) * sqrt(Psi->AllPsiStatus[ActNum].PsiFactor*Psi->LocalPsiStatus[R->ActualLocalPsiNo].PsiFactor);
1505 result += E1;
1506 //fprintf(stderr,"(%i) 1st: E1 = \t\t%lg\n", P->Par.me, E1);
1507
1508 return result;
1509}
1510
1511
1512/** Evaluates second derivative of perturbed energy functional with respect to minimisation parameter \f$\Theta\f$.
1513 * \f[
1514 * \frac{\delta^2 {\cal E}^{(2)}} {\delta \Theta^2} =
1515 * 2 \bigl( \langle \widetilde{\varphi}_l^{(1)} | {\cal H}^{(0)} | \widetilde{\varphi}_l^{(1)} \rangle
1516 * - \langle \varphi_l^{(1)} | {\cal H}^{(0)} | \varphi_l^{(1)} \rangle \bigr )
1517 * + 2 \sum_{i,i \neq l } \lambda_{il} \langle \varphi_i^{(1)} | \varphi_l^{(1)} \rangle
1518 * - 2 {\cal R} \langle \varphi_l^{(1)} | {\cal H}^{(1)} | \varphi_l^{(0)} \rangle
1519 * \f]
1520 *
1521 * The energy eigenvalues of \a ConDir and \a source must be supplied, they can be calculated via CalculateConDirHConDir() and/or
1522 * by the due to CalculatePerturbedEnergy() already present OnePsiElementAddData#Lambda eigenvalue. The summation over the
1523 * unperturbed lambda within the scalar product of perturbed wave functions is evaluated with Psis#lambda and Psis#Overlap.
1524 * Afterwards, the ConDir density is calculated and also the i-th perturbed density to first degree. With these in a sum over
1525 * all real mesh points the exchange-correlation first and second derivatives and also the Hartree potential ones can be calculated
1526 * and summed up.
1527 * \param *P Problem at hand
1528 * \param source0 unperturbed wave function \f$\varphi_l^{(0)}\f$
1529 * \param source wave function \f$\varphi_l^{(1)}\f$
1530 * \param ConDir conjugated direction \f$\widetilde{\varphi}_l^{(1)}\f$
1531 * \param sourceHsource eigenvalue of wave function \f$\langle \varphi_l^{(1)} | H^{(0)} | \varphi_l^{(1)}\rangle\f$
1532 * \param ConDirHConDir perturbed eigenvalue of conjugate direction \f$\langle \widetilde{\varphi}_l^{(1)} | H^{(0)} | \widetilde{\varphi}_l^{(1)}\rangle\f$
1533 * \param ConDirConDir norm of conjugate direction \f$\langle \widetilde{\varphi}_l^{(1)} | \widetilde{\varphi}_l^{(1)}\rangle\f$
1534 * \warning No MPI_AllReduce() takes place, parameters have to be reduced already.
1535 */
1536double Calculate2ndPerturbedDerivative(struct Problem *P, const fftw_complex *source0,const fftw_complex *source, const fftw_complex *ConDir,const double sourceHsource, const double ConDirHConDir, const double ConDirConDir)
1537{
1538 struct RunStruct *R = &P->R;
1539 struct Psis *Psi = &P->Lat.Psi;
1540 //struct Lattice *Lat = &P->Lat;
1541 //struct Energy *E = Lat->E;
1542 double result = 0.;
1543 double Con0 = 0., Elambda = 0.;//, E0 = 0., E1 = 0.;
1544 //int i;
1545 const int state = R->CurrentMin;
1546 //const int l_normal = R->ActualLocalPsiNo - Psi->TypeStartIndex[state] + Psi->TypeStartIndex[Occupied];
1547 const int ActNum = R->ActualLocalPsiNo - Psi->TypeStartIndex[state] + Psi->TypeStartIndex[1] * Psi->LocalPsiStatus[R->ActualLocalPsiNo].my_color_comm_ST_Psi;
1548
1549 Con0 = 2.*ConDirHConDir;
1550 result += Con0;
1551 ////E0 = -2.*sourceHsource;
1552 ////result += E0;
1553 ////E1 = -E->PsiEnergy[Perturbed1_0Energy][R->ActualLocalPsiNo] - E->PsiEnergy[Perturbed0_1Energy][R->ActualLocalPsiNo];
1554 ////result += E1;
1555 //fprintf(stderr,"(%i) 2nd: E1 = \t%lg\n", P->Par.me, E1);
1556
1557 ////for (i=0;i<Lat->Psi.NoOfPsis;i++) {
1558 //// if (i != ActNum) Elambda += Psi->lambda[i][ActNum]*Psi->Overlap[i][ActNum]+ Psi->lambda[ActNum][i]*Psi->Overlap[ActNum][i]; // overlap contains PsiFactor
1559 ////}
1560 ////Elambda = Psi->lambda[ActNum][ActNum]*Psi->Overlap[ActNum][ActNum];
1561 Elambda = 2.*Psi->lambda[ActNum][ActNum]*ConDirConDir;
1562 result -= Elambda;
1563
1564 //fprintf(stderr,"(%i) 2ndPerturbedDerivative: Result = Con0 + E0 + E1 + Elambda + dEdt0_XC + ddEddt0_XC + dEdt0_H + ddEddt0_H = %lg + %lg + %lg + %lg + %lg + %lg + %lg + %lg = %lg\n", P->Par.me, Con0, E0, E1, Elambda, VolumeFactorR*dEdt0_XC, VolumeFactorR*ddEddt0_XC, dEdt0_H, ddEddt0_H, result);
1565
1566 return (result);
1567}
1568
1569/** Returns index of specific component in 3x3 cross product.
1570 * \param i vector product component index, ranging from 0..NDIM
1571 * \param j index specifies which one of the four vectors in x*y - y*x, ranging from 0..3 (0,1 positive sign, 2,3 negative sign)
1572 * \return Component 0..2 of vector to be taken to evaluate a vector product
1573 * \sa crossed() - is the same but vice versa, return value must be specified, \a i is returned.
1574 */
1575inline int cross(int i, int j)
1576{
1577 const int matrix[NDIM*4] = {1,2,2,1,2,0,0,2,0,1,1,0};
1578 if (i>=0 && i<NDIM && j>=0 && j<4)
1579 return (matrix[i*4+j]);
1580 else {
1581 Error(SomeError,"cross: i or j out of range!");
1582 return (0);
1583 }
1584}
1585
1586/** Returns index of resulting vector component in 3x3 cross product.
1587 * In the column specified by the \a j index \a i is looked for and the found row index returned.
1588 * \param i vector component index, ranging from 0..NDIM
1589 * \param j index specifies which one of the four vectors in x*y - y*x, ranging from 0..3 (0,1 positive sign, 2,3 negative sign)
1590 * \return Component 0..2 of resulting vector
1591 * \sa cross() - is the same but vice versa, return value must be specified, \a i is returned.
1592 */
1593inline int crossed(int i, int j)
1594{
1595 const int matrix[NDIM*4] = {1,2,2,1,2,0,0,2,0,1,1,0};
1596 int k;
1597 if (i>=0 && i<NDIM && j>=0 && j<4) {
1598 for (k=0;k<NDIM;k++)
1599 if (matrix[4*k+j] == i) return(k);
1600 Error(SomeError,"crossed: given component not found!");
1601 return(-1);
1602 } else {
1603 Error(SomeError,"crossed: i or j out of range!");
1604 return (-1);
1605 }
1606}
1607
1608#define Nsin 16 //!< should be dependent on MaxG/MaxN per axis!
1609
1610/** Returns sawtooth shaped profile for position operator within cell.
1611 * This is a mapping from -L/2...L/2 (L = length of unit cell derived from Lattice#RealBasisSQ) to -L/2 to L/2 with a smooth transition:
1612 * \f[
1613 * f(x): x \rightarrow \left \{
1614 * \begin{array}{l}
1615 * -\frac{L}{2} \cdot \sin \left ( \frac{x}{0,05\cdot L} \cdot \frac{\pi}{2} \right ), 0<x<0,05\cdot L \\
1616 * (x - 0,05\cdot L) \cdot \frac{10}{9} - \frac{L}{2}, 0,05\cdot L \leq x<0,95\cdot L \\
1617 * \frac{L}{2} \cdot \cos \left ( \frac{x-0,95\cdot L}{0,05\cdot L} \cdot \frac{\pi}{2} \right), 0,95\cdot L<x<L
1618 * \end{array} \right \}
1619 * \f]
1620 * \param *Lat pointer to Lattice structure for Lattice#RealBasisSQ
1621 * \param L parameter x
1622 * \param index component index for Lattice#RealBasisSQ
1623 */
1624inline double sawtooth(struct Lattice *Lat, double L, const int index)
1625{
1626 double axis = sqrt(Lat->RealBasisSQ[index]);
1627 double sawstart = Lat->SawtoothStart;
1628 double sawend = 1. - sawstart;
1629 double sawfactor = (sawstart+sawend)/(sawend-sawstart);
1630 //return(L);
1631
1632 //fprintf(stderr, "sawstart: %e\tsawend: %e\tsawfactor: %e\tL: %e\n", sawstart, sawend, sawfactor, L);
1633 // transform and return (sawtooth profile checked, 04.08.06)
1634 L += axis/2.; // transform to 0 ... L
1635 if (L < (sawstart*axis)) return (-axis/(2*sawfactor)*sin(L/(sawstart*axis)*PI/2.)); // first smooth transition from 0 ... -L/2
1636 if (L > (sawend*axis)) return ( axis/(2*sawfactor)*cos((L-sawend*axis)/(sawstart*axis)*PI/2.)); // second smooth transition from +L/2 ... 0
1637 //fprintf(stderr,"L %e\t sawstart %e\t sawend %e\t sawfactor %e\t axis%e\n", L, sawstart, sawend, sawfactor, axis);
1638 //return ((L - sawstart*axis) - axis/(2*sawfactor)); // area in between scale to -L/2 ... +L/2
1639 return (L - axis/2); // area in between return as it was
1640}
1641
1642/** Shifts the origin of the gauge according to the CSDGT method.
1643 * \f[
1644 * d(r) = r - \sum_{I_s,I_a} (r-R_{I_s,I_a}) exp{(-\alpha_{I_s,I_a}(r-R_{I_s,I_a})^4)}
1645 * \f]
1646 * This trafo is necessary as the current otherweise (CSGT) sensitively depends on the current around
1647 * the core region inadequately/only moderately well approximated by a plane-wave-pseudo-potential-method.
1648 * \param *P Problem at hand, containing Lattice and Ions
1649 * \param r parameter r
1650 * \param index index of the basis vector
1651 * \return \f$d(r)\f$
1652 * \note Continuous Set of Damped Gauge Transformations according to Keith and Bader
1653 */
1654inline double ShiftGaugeOrigin(struct Problem *P, double r, const int index)
1655{
1656 struct Ions *I = &P->Ion;
1657 struct Lattice *Lat = &P->Lat;
1658 double x, tmp;
1659 int is,ia;
1660
1661 // loop over all ions to calculate the sum
1662 x = r;
1663 for (is=0; is < I->Max_Types; is++)
1664 for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
1665 tmp = (r - I->I[is].R[NDIM*ia]);
1666 x -= tmp*exp(- I->I[is].alpha[ia] * tpow(tmp,4));
1667 }
1668
1669 return(sawtooth(Lat,x,index)); // still use sawtooth due to the numerical instability around the border region of the cell
1670}
1671
1672/** Print sawtooth() for each node along one axis.
1673 * \param *P Problem at hand, containing RunStruct, Lattice and LatticeLevel RunStruct#LevS
1674 * \param index index of axis
1675 */
1676void TestSawtooth(struct Problem *P, const int index)
1677{
1678 struct RunStruct *R = &P->R;
1679 struct LatticeLevel *LevS = R->LevS;
1680 struct Lattice *Lat =&P->Lat;
1681 double x[NDIM];
1682 double n[NDIM];
1683 int N[NDIM];
1684 N[0] = LevS->Plan0.plan->N[0];
1685 N[1] = LevS->Plan0.plan->N[1];
1686 N[2] = LevS->Plan0.plan->N[2];
1687
1688 n[0] = n[1] = n[2] = 0.;
1689 for (n[index]=0;n[index]<N[index];n[index]++) {
1690 n[index] = (double)n[index]/(double)N[index] * sqrt(Lat->RealBasisSQ[index]);
1691 //fprintf(stderr,"(%i) x %e\t Axis/2 %e\n",P->Par.me, x, sqrt(Lat->RealBasisSQ[index])/2. );
1692 MinImageConv(Lat, n, Lat->RealBasisCenter, x);
1693 fprintf(stderr,"%e\t%e\n", n[index], sawtooth(Lat,n[index],index));
1694 }
1695}
1696
1697/** Secures minimum image convention between two given points \a R[] and \a r[] within periodic boundary.
1698 * Each distance component within a periodic boundary must always be between -L/2 ... L/2
1699 * \param *Lat pointer to Lattice structure
1700 * \param R[] first vector, NDIM, each must be between 0...L
1701 * \param r[] second vector, NDIM, each must be between 0...L
1702 * \param result[] return vector
1703 */
1704inline void MinImageConv(struct Lattice *Lat, const double R[NDIM], const double r[NDIM], double *result)
1705{
1706 //double axis = Lat->RealBasisQ[index];
1707 double x[NDIM], X[NDIM], Result[NDIM];
1708 int i;
1709
1710 for(i=0;i<NDIM;i++)
1711 result[i] = x[i] = x[i] = 0.;
1712 //fprintf(stderr, "R = (%lg, %lg, %lg), r = (%lg, %lg, %lg)\n", R[0], R[1], R[2], r[0], r[1], r[2]);
1713 RMat33Vec3(X, Lat->ReciBasis, R); // transform both to [0,1]^3
1714 RMat33Vec3(x, Lat->ReciBasis, r);
1715 //fprintf(stderr, "X = (%lg, %lg, %lg), x = (%lg, %lg, %lg)\n", X[0], X[1], X[2], x[0], x[1], x[2]);
1716 for(i=0;i<NDIM;i++) {
1717// if (fabs(X[i]) > 1.)
1718// fprintf(stderr,"X[%i] > 1. : %lg!\n", i, X[i]);
1719// if (fabs(x[i]) > 1.)
1720// fprintf(stderr,"x[%i] > 1. : %lg!\n", i, x[i]);
1721 if (fabs(Result[i] = X[i] - x[i] + 2.*PI) < PI) { }
1722 else if (fabs(Result[i] = X[i] - x[i]) <= PI) { }
1723 else if (fabs(Result[i] = X[i] - x[i] - 2.*PI) < PI) { }
1724 else Error(SomeError, "MinImageConv: None of the three cases applied!");
1725 }
1726 for(i=0;i<NDIM;i++) // ReciBasis is not true inverse, but times 2.*PI
1727 Result[i] /= 2.*PI;
1728 RMat33Vec3(result, Lat->RealBasis, Result);
1729}
1730
1731/** Linear interpolation for coordinate \a R that lies between grid nodes of \a *grid.
1732 * \param *P Problem at hand
1733 * \param *Lat Lattice structure for grid axis
1734 * \param *Lev LatticeLevel structure for grid axis node counts
1735 * \param R[] coordinate vector
1736 * \param *grid grid with fixed nodes
1737 * \return linearly interpolated value of \a *grid for position \a R[NDIM]
1738 */
1739double LinearInterpolationBetweenGrid(struct Problem *P, struct Lattice *Lat, struct LatticeLevel *Lev, double R[NDIM], fftw_real *grid)
1740{
1741 double x[2][NDIM];
1742 const int myPE = P->Par.me_comm_ST_Psi;
1743 int N[NDIM];
1744 const int N0 = Lev->Plan0.plan->local_nx;
1745 N[0] = Lev->Plan0.plan->N[0];
1746 N[1] = Lev->Plan0.plan->N[1];
1747 N[2] = Lev->Plan0.plan->N[2];
1748 int g;
1749 double n[NDIM];
1750 int k[2][NDIM];
1751 double sigma;
1752
1753 RMat33Vec3(n, Lat->ReciBasis, &R[0]); // transform real coordinates to [0,1]^3 vector
1754 for (g=0;g<NDIM;g++) {
1755 // k[i] are right and left nearest neighbour node to true position
1756 k[0][g] = floor(n[g]/(2.*PI)*(double)N[g]); // n[2] is floor grid
1757 k[1][g] = ceil(n[g]/(2.*PI)*(double)N[g]); // n[1] is ceil grid
1758 // x[i] give weights of left and right neighbours (the nearer the true point is to one, the closer its weight to 1)
1759 x[0][g] = (k[1][g] - n[g]/(2.*PI)*(double)N[g]);
1760 x[1][g] = 1. - x[0][g];
1761 //fprintf(stderr,"(%i) n = %lg, n_floor[%i] = %i\tn_ceil[%i] = %i --- x_floor[%i] = %e\tx_ceil[%i] = %e\n",P->Par.me, n[g], g,k[0][g], g,k[1][g], g,x[0][g], g,x[1][g]);
1762 }
1763 sigma = 0.;
1764 for (g=0;g<2;g++) { // interpolate linearly between adjacent grid points per axis
1765 if ((k[g][0] >= N0*myPE) && (k[g][0] < N0*(myPE+1))) {
1766 //fprintf(stderr,"(%i) grid[%i]: sigma = %e\n", P->Par.me, k[g][2]+N[2]*(k[g][1]+N[1]*(k[g][0]-N0*myPE)), sigma);
1767 sigma += (x[g][0]*x[0][1]*x[0][2])*grid[k[0][2]+N[2]*(k[0][1]+N[1]*(k[g][0]-N0*myPE))]*mu0; // if it's local and factor from inverse fft
1768 //fprintf(stderr,"(%i) grid[%i]: sigma += %e * %e \n", P->Par.me, k[g][2]+N[2]*(k[g][1]+N[1]*(k[g][0]-N0*myPE)), (x[g][0]*x[0][1]*x[0][2]), grid[k[0][2]+N[2]*(k[0][1]+N[1]*(k[g][0]-N0*myPE))]*mu0);
1769 sigma += (x[g][0]*x[0][1]*x[1][2])*grid[k[1][2]+N[2]*(k[0][1]+N[1]*(k[g][0]-N0*myPE))]*mu0; // if it's local and factor from inverse fft
1770 //fprintf(stderr,"(%i) grid[%i]: sigma += %e * %e \n", P->Par.me, k[g][2]+N[2]*(k[g][1]+N[1]*(k[g][0]-N0*myPE)), (x[g][0]*x[0][1]*x[1][2]), grid[k[1][2]+N[2]*(k[0][1]+N[1]*(k[g][0]-N0*myPE))]*mu0);
1771 sigma += (x[g][0]*x[1][1]*x[0][2])*grid[k[0][2]+N[2]*(k[1][1]+N[1]*(k[g][0]-N0*myPE))]*mu0; // if it's local and factor from inverse fft
1772 //fprintf(stderr,"(%i) grid[%i]: sigma += %e * %e \n", P->Par.me, k[g][2]+N[2]*(k[g][1]+N[1]*(k[g][0]-N0*myPE)), (x[g][0]*x[1][1]*x[0][2]), grid[k[0][2]+N[2]*(k[1][1]+N[1]*(k[g][0]-N0*myPE))]*mu0);
1773 sigma += (x[g][0]*x[1][1]*x[1][2])*grid[k[1][2]+N[2]*(k[1][1]+N[1]*(k[g][0]-N0*myPE))]*mu0; // if it's local and factor from inverse fft
1774 //fprintf(stderr,"(%i) grid[%i]: sigma += %e * %e \n", P->Par.me, k[g][2]+N[2]*(k[g][1]+N[1]*(k[g][0]-N0*myPE)), (x[g][0]*x[1][1]*x[1][2]), grid[k[1][2]+N[2]*(k[1][1]+N[1]*(k[g][0]-N0*myPE))]*mu0);
1775 }
1776 }
1777 return sigma;
1778}
1779
1780/** Linear Interpolation from all eight corners of the box that singles down to a point on the lower level.
1781 * \param *P Problem at hand
1782 * \param *Lev LatticeLevel structure for node numbers
1783 * \param upperNode Node around which to interpolate
1784 * \param *upperGrid array of grid points
1785 * \return summed up and then averaged octant around \a upperNode
1786 */
1787double LinearPullDownFromUpperLevel(struct Problem *P, struct LatticeLevel *Lev, int upperNode, fftw_real *upperGrid)
1788{
1789 const int N0 = Lev->Plan0.plan->local_nx;
1790 const int N1 = Lev->Plan0.plan->N[1];
1791 const int N2 = Lev->Plan0.plan->N[2];
1792 double lowerGrid = 0.;
1793 int nr=1;
1794 lowerGrid += upperGrid[upperNode];
1795 if (upperNode % N0 != N0-1) {
1796 lowerGrid += upperGrid[upperNode+1];
1797 nr++;
1798 if (upperNode % N1 != N1-1) {
1799 lowerGrid += upperGrid[upperNode + 0 + N2*(1 + N1*1)];
1800 nr++;
1801 if (upperNode % N2 != N2-1) {
1802 lowerGrid += upperGrid[upperNode + 1 + N2*(1 + N1*1)];
1803 nr++;
1804 }
1805 }
1806 if (upperNode % N2 != N2-1) {
1807 lowerGrid += upperGrid[upperNode + 1 + N2*(0 + N1*1)];
1808 nr++;
1809 }
1810 }
1811 if (upperNode % N1 != N1-1) {
1812 lowerGrid += upperGrid[upperNode + 0 + N2*(1 + N1*0)];
1813 nr++;
1814 if (upperNode % N2 != N2-1) {
1815 lowerGrid += upperGrid[upperNode + 1 + N2*(1 + N1*0)];
1816 nr++;
1817 }
1818 }
1819 if (upperNode % N2 != N2-1) {
1820 lowerGrid += upperGrid[upperNode + 1 + N2*(0 + N1*0)];
1821 nr++;
1822 }
1823 return (lowerGrid/(double)nr);
1824}
1825
1826/** Evaluates the 1-stern in order to evaluate the first derivative on the grid.
1827 * \param *P Problem at hand
1828 * \param *Lev Level to interpret the \a *density on
1829 * \param *density array with gridded values
1830 * \param *n 3 vector with indices on the grid
1831 * \param axis axis along which is derived
1832 * \param myPE number of processes who share the density
1833 * \return [+1/2 -1/2] of \a *n
1834 */
1835double FirstDiscreteDerivative(struct Problem *P, struct LatticeLevel *Lev, fftw_real *density, int *n, int axis, int myPE)
1836{
1837 int *N = Lev->Plan0.plan->N; // maximum nodes per axis
1838 const int N0 = Lev->Plan0.plan->local_nx; // special local number due to parallel split up
1839 double ret[NDIM], Ret[NDIM]; // return value local/global
1840 int i;
1841
1842 for (i=0;i<NDIM;i++) {
1843 ret[i] = Ret[i] = 0.;
1844 }
1845 if (((n[0]+1)%N[0] >= N0*myPE) && ((n[0]+1)%N[0] < N0*(myPE+1))) // next cell belongs to this process
1846 ret[0] += 1./2. * (density[n[2]+N[2]*(n[1]+N[1]*(n[0]+1-N0*myPE))]);
1847 if (((n[0]-1)%N[0] >= N0*myPE) && ((n[0]-1)%N[0] < N0*(myPE+1))) // previous cell belongs to this process
1848 ret[0] -= 1./2. * (density[n[2]+N[2]*(n[1]+N[1]*(n[0]-1-N0*myPE))]);
1849 if ((n[0] >= N0*myPE) && (n[0] < N0*(myPE+1))) {
1850 ret[1] += 1./2. * (density[n[2]+N[2]*((n[1]+1)%N[1] + N[1]*(n[0]%N0))]);
1851 ret[1] -= 1./2. * (density[n[2]+N[2]*((n[1]-1)%N[1] + N[1]*(n[0]%N0))]);
1852 }
1853 if ((n[0] >= N0*myPE) && (n[0] < N0*(myPE+1))) {
1854 ret[2] += 1./2. * (density[(n[2]+1)%N[2] + N[2]*(n[1]+N[1]*(n[0]%N0))]);
1855 ret[2] -= 1./2. * (density[(n[2]-1)%N[2] + N[2]*(n[1]+N[1]*(n[0]%N0))]);
1856 }
1857
1858 if (MPI_Allreduce(ret, Ret, 3, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi) != MPI_SUCCESS)
1859 Error(SomeError, "FirstDiscreteDerivative: MPI_Allreduce failure!");
1860
1861 for (i=0;i<NDIM;i++) // transform from node count to [0,1]^3
1862 Ret[i] *= N[i];
1863 RMat33Vec3(ret, P->Lat.ReciBasis, Ret); // this actually divides it by mesh length in real coordinates
1864 //fprintf(stderr, "(%i) sum at (%i,%i,%i) : %lg\n",P->Par.me, n[0],n[1],n[2], ret[axis]);
1865 return ret[axis]; ///(P->Lat.RealBasisQ[axis]/N[axis]);
1866}
1867
1868/** Fouriertransforms given \a source.
1869 * By the use of the symmetry parameter an additional imaginary unit and/or the momentum operator can
1870 * be applied at the same time.
1871 * \param *P Problem at hand
1872 * \param *Psi source array of reciprocal coefficients
1873 * \param *PsiR destination array, becoming filled with real coefficients
1874 * \param index_g component of G vector (only needed for symmetry=4..7)
1875 * \param symmetry 0 - do nothing, 1 - factor by "-1", 2 - factor by "i", 3 - factor by "1/i = -i", from 4 to 7 the same
1876 * but additionally with momentum operator
1877 */
1878void fft_Psi(struct Problem *P, const fftw_complex *Psi, fftw_real *PsiR, const int index_g, const int symmetry)
1879{
1880 struct Lattice *Lat = &P->Lat;
1881 struct RunStruct *R = &P->R;
1882 struct LatticeLevel *Lev0 = R->Lev0;
1883 struct LatticeLevel *LevS = R->LevS;
1884 struct Density *Dens0 = Lev0->Dens;
1885 struct fft_plan_3d *plan = Lat->plan;
1886 fftw_complex *tempdestRC = (fftw_complex *)Dens0->DensityArray[TempDensity];
1887 fftw_complex *work = Dens0->DensityCArray[TempDensity];
1888 fftw_complex *posfac, *destpos, *destRCS, *destRCD;
1889 int i, Index, pos;
1890
1891 LockDensityArray(Dens0,TempDensity,imag); // tempdestRC
1892 SetArrayToDouble0((double *)tempdestRC, Dens0->TotalSize*2);
1893 SetArrayToDouble0((double *)PsiR, Dens0->TotalSize*2);
1894 switch (symmetry) {
1895 case 0:
1896 for (i=0;i<LevS->MaxG;i++) { // incoming is positive, outgoing is positive
1897 Index = LevS->GArray[i].Index;
1898 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
1899 destpos = &tempdestRC[LevS->MaxNUp*Index];
1900 for (pos=0; pos < LevS->MaxNUp; pos++) {
1901 //if (destpos != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destpos corrupted");
1902 destpos[pos].re = (Psi[i].re)*posfac[pos].re-(Psi[i].im)*posfac[pos].im;
1903 destpos[pos].im = (Psi[i].re)*posfac[pos].im+(Psi[i].im)*posfac[pos].re;
1904 }
1905 }
1906 break;
1907 case 1:
1908 for (i=0;i<LevS->MaxG;i++) { // incoming is positive, outgoing is - positive
1909 Index = LevS->GArray[i].Index;
1910 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
1911 destpos = &tempdestRC[LevS->MaxNUp*Index];
1912 for (pos=0; pos < LevS->MaxNUp; pos++) {
1913 //if (destpos != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destpos corrupted");
1914 destpos[pos].re = -((Psi[i].re)*posfac[pos].re-(Psi[i].im)*posfac[pos].im);
1915 destpos[pos].im = -((Psi[i].re)*posfac[pos].im+(Psi[i].im)*posfac[pos].re);
1916 }
1917 }
1918 break;
1919 case 2:
1920 for (i=0;i<LevS->MaxG;i++) { // incoming is positive, outgoing is negative
1921 Index = LevS->GArray[i].Index;
1922 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
1923 destpos = &tempdestRC[LevS->MaxNUp*Index];
1924 for (pos=0; pos < LevS->MaxNUp; pos++) {
1925 //if (destpos != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destpos corrupted");
1926 destpos[pos].re = (-Psi[i].im)*posfac[pos].re-(Psi[i].re)*posfac[pos].im;
1927 destpos[pos].im = (-Psi[i].im)*posfac[pos].im+(Psi[i].re)*posfac[pos].re;
1928 }
1929 }
1930 break;
1931 case 3:
1932 for (i=0;i<LevS->MaxG;i++) { // incoming is negative, outgoing is positive
1933 Index = LevS->GArray[i].Index;
1934 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
1935 destpos = &tempdestRC[LevS->MaxNUp*Index];
1936 for (pos=0; pos < LevS->MaxNUp; pos++) {
1937 //if (destpos != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destpos corrupted");
1938 destpos[pos].re = (Psi[i].im)*posfac[pos].re-(-Psi[i].re)*posfac[pos].im;
1939 destpos[pos].im = (Psi[i].im)*posfac[pos].im+(-Psi[i].re)*posfac[pos].re;
1940 }
1941 }
1942 break;
1943 case 4:
1944 for (i=0;i<LevS->MaxG;i++) { // incoming is positive, outgoing is positive
1945 Index = LevS->GArray[i].Index;
1946 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
1947 destpos = &tempdestRC[LevS->MaxNUp*Index];
1948 for (pos=0; pos < LevS->MaxNUp; pos++) {
1949 //if (destpos != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destpos corrupted");
1950 destpos[pos].re = LevS->GArray[i].G[index_g]*((Psi[i].re)*posfac[pos].re-(Psi[i].im)*posfac[pos].im);
1951 destpos[pos].im = LevS->GArray[i].G[index_g]*((Psi[i].re)*posfac[pos].im+(Psi[i].im)*posfac[pos].re);
1952 }
1953 }
1954 break;
1955 case 5:
1956 for (i=0;i<LevS->MaxG;i++) { // incoming is positive, outgoing is - positive
1957 Index = LevS->GArray[i].Index;
1958 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
1959 destpos = &tempdestRC[LevS->MaxNUp*Index];
1960 for (pos=0; pos < LevS->MaxNUp; pos++) {
1961 //if (destpos != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destpos corrupted");
1962 destpos[pos].re = -LevS->GArray[i].G[index_g]*((Psi[i].re)*posfac[pos].re-(Psi[i].im)*posfac[pos].im);
1963 destpos[pos].im = -LevS->GArray[i].G[index_g]*((Psi[i].re)*posfac[pos].im+(Psi[i].im)*posfac[pos].re);
1964 }
1965 }
1966 break;
1967 case 6:
1968 for (i=0;i<LevS->MaxG;i++) { // incoming is positive, outgoing is negative
1969 Index = LevS->GArray[i].Index;
1970 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
1971 destpos = &tempdestRC[LevS->MaxNUp*Index];
1972 for (pos=0; pos < LevS->MaxNUp; pos++) {
1973 //if (destpos != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destpos corrupted");
1974 destpos[pos].re = LevS->GArray[i].G[index_g]*((-Psi[i].im)*posfac[pos].re-(Psi[i].re)*posfac[pos].im);
1975 destpos[pos].im = LevS->GArray[i].G[index_g]*((-Psi[i].im)*posfac[pos].im+(Psi[i].re)*posfac[pos].re);
1976 }
1977 }
1978 break;
1979 case 7:
1980 for (i=0;i<LevS->MaxG;i++) { // incoming is negative, outgoing is positive
1981 Index = LevS->GArray[i].Index;
1982 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
1983 destpos = &tempdestRC[LevS->MaxNUp*Index];
1984 for (pos=0; pos < LevS->MaxNUp; pos++) {
1985 //if (destpos != &tempdestRC[LevS->MaxNUp*Index] || LevS->MaxNUp*Index+pos<0 || LevS->MaxNUp*Index+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destpos corrupted");
1986 destpos[pos].re = LevS->GArray[i].G[index_g]*((Psi[i].im)*posfac[pos].re-(-Psi[i].re)*posfac[pos].im);
1987 destpos[pos].im = LevS->GArray[i].G[index_g]*((Psi[i].im)*posfac[pos].im+(-Psi[i].re)*posfac[pos].re);
1988 }
1989 }
1990 break;
1991 }
1992 for (i=0; i<LevS->MaxDoubleG; i++) {
1993 destRCS = &tempdestRC[LevS->DoubleG[2*i]*LevS->MaxNUp];
1994 destRCD = &tempdestRC[LevS->DoubleG[2*i+1]*LevS->MaxNUp];
1995 for (pos=0; pos < LevS->MaxNUp; pos++) {
1996 //if (destRCD != &tempdestRC[LevS->DoubleG[2*i+1]*LevS->MaxNUp] || LevS->DoubleG[2*i+1]*LevS->MaxNUp+pos<0 || LevS->DoubleG[2*i+1]*LevS->MaxNUp+pos>=Dens0->TotalSize) Error(SomeError,"fft_Psi: destRCD corrupted");
1997 destRCD[pos].re = destRCS[pos].re;
1998 destRCD[pos].im = -destRCS[pos].im;
1999 }
2000 }
2001 fft_3d_complex_to_real(plan, LevS->LevelNo, FFTNFUp, tempdestRC, work);
2002 DensityRTransformPos(LevS,(fftw_real*)tempdestRC, PsiR);
2003 UnLockDensityArray(Dens0,TempDensity,imag); // tempdestRC
2004}
2005
2006/** Locks all NDIM_NDIM current density arrays
2007 * \param Dens0 Density structure to be locked (in the current parts)
2008 */
2009void AllocCurrentDensity(struct Density *Dens0) {
2010 // real
2011 LockDensityArray(Dens0,CurrentDensity0,real); // CurrentDensity[B_index]
2012 LockDensityArray(Dens0,CurrentDensity1,real); // CurrentDensity[B_index]
2013 LockDensityArray(Dens0,CurrentDensity2,real); // CurrentDensity[B_index]
2014 LockDensityArray(Dens0,CurrentDensity3,real); // CurrentDensity[B_index]
2015 LockDensityArray(Dens0,CurrentDensity4,real); // CurrentDensity[B_index]
2016 LockDensityArray(Dens0,CurrentDensity5,real); // CurrentDensity[B_index]
2017 LockDensityArray(Dens0,CurrentDensity6,real); // CurrentDensity[B_index]
2018 LockDensityArray(Dens0,CurrentDensity7,real); // CurrentDensity[B_index]
2019 LockDensityArray(Dens0,CurrentDensity8,real); // CurrentDensity[B_index]
2020 // imaginary
2021 LockDensityArray(Dens0,CurrentDensity0,imag); // CurrentDensity[B_index]
2022 LockDensityArray(Dens0,CurrentDensity1,imag); // CurrentDensity[B_index]
2023 LockDensityArray(Dens0,CurrentDensity2,imag); // CurrentDensity[B_index]
2024 LockDensityArray(Dens0,CurrentDensity3,imag); // CurrentDensity[B_index]
2025 LockDensityArray(Dens0,CurrentDensity4,imag); // CurrentDensity[B_index]
2026 LockDensityArray(Dens0,CurrentDensity5,imag); // CurrentDensity[B_index]
2027 LockDensityArray(Dens0,CurrentDensity6,imag); // CurrentDensity[B_index]
2028 LockDensityArray(Dens0,CurrentDensity7,imag); // CurrentDensity[B_index]
2029 LockDensityArray(Dens0,CurrentDensity8,imag); // CurrentDensity[B_index]
2030}
2031
2032/** Reset and unlocks all NDIM_NDIM current density arrays
2033 * \param Dens0 Density structure to be unlocked/resetted (in the current parts)
2034 */
2035void DisAllocCurrentDensity(struct Density *Dens0) {
2036 //int i;
2037 // real
2038// for(i=0;i<NDIM*NDIM;i++)
2039// SetArrayToDouble0((double *)Dens0->DensityArray[i], Dens0->TotalSize*2);
2040 UnLockDensityArray(Dens0,CurrentDensity0,real); // CurrentDensity[B_index]
2041 UnLockDensityArray(Dens0,CurrentDensity1,real); // CurrentDensity[B_index]
2042 UnLockDensityArray(Dens0,CurrentDensity2,real); // CurrentDensity[B_index]
2043 UnLockDensityArray(Dens0,CurrentDensity3,real); // CurrentDensity[B_index]
2044 UnLockDensityArray(Dens0,CurrentDensity4,real); // CurrentDensity[B_index]
2045 UnLockDensityArray(Dens0,CurrentDensity5,real); // CurrentDensity[B_index]
2046 UnLockDensityArray(Dens0,CurrentDensity6,real); // CurrentDensity[B_index]
2047 UnLockDensityArray(Dens0,CurrentDensity7,real); // CurrentDensity[B_index]
2048 UnLockDensityArray(Dens0,CurrentDensity8,real); // CurrentDensity[B_index]
2049 // imaginary
2050// for(i=0;i<NDIM*NDIM;i++)
2051// SetArrayToDouble0((double *)Dens0->DensityCArray[i], Dens0->TotalSize*2);
2052 UnLockDensityArray(Dens0,CurrentDensity0,imag); // CurrentDensity[B_index]
2053 UnLockDensityArray(Dens0,CurrentDensity1,imag); // CurrentDensity[B_index]
2054 UnLockDensityArray(Dens0,CurrentDensity2,imag); // CurrentDensity[B_index]
2055 UnLockDensityArray(Dens0,CurrentDensity3,imag); // CurrentDensity[B_index]
2056 UnLockDensityArray(Dens0,CurrentDensity4,imag); // CurrentDensity[B_index]
2057 UnLockDensityArray(Dens0,CurrentDensity5,imag); // CurrentDensity[B_index]
2058 UnLockDensityArray(Dens0,CurrentDensity6,imag); // CurrentDensity[B_index]
2059 UnLockDensityArray(Dens0,CurrentDensity7,imag); // CurrentDensity[B_index]
2060 UnLockDensityArray(Dens0,CurrentDensity8,imag); // CurrentDensity[B_index]
2061}
2062
2063// these defines safe-guard same symmetry for same kind of wave function
2064#define Psi0symmetry 0 // //0 //0 //0 // regard psi0 as real
2065#define Psi1symmetry 0 // //3 //0 //0 // regard psi0 as real
2066#define Psip0symmetry 6 //6 //6 //6 //6 // momentum times "i" due to operation on left hand
2067#define Psip1symmetry 7 //7 //4 //6 //7 // momentum times "-i" as usual (right hand)
2068
2069/** Evaluates the 3x3 current density arrays.
2070 * The formula we want to evaluate is as follows
2071 * \f[
2072 * j_k(r) = \langle \psi_k^{(0)} | \Bigl ( p|r'\rangle\langle r' | + | r' \rangle \langle r' | p \Bigr )
2073 \Bigl [ | \psi_k^{(r\times p )} \rangle - r' \times | \psi_k^{(p)} \rangle \Bigr ] \cdot B.
2074 * \f]
2075 * Most of the DensityTypes-arrays are locked for temporary use. Pointers are set to their
2076 * start address and afterwards the current density arrays locked and reset'ed. Then for every
2077 * unperturbed wave function we do:
2078 * -# FFT unperturbed p-perturbed and rxp-perturbed wave function
2079 * -# FFT wave function with applied momentum operator for all three indices
2080 * -# For each index of the momentum operator:
2081 * -# FFT p-perturbed wave function
2082 * -# For every index of the external field:
2083 * -# FFT rxp-perturbed wave function
2084 * -# Evaluate current density for these momentum index and external field indices
2085 *
2086 * Afterwards the temporary densities are unlocked and the density ones gathered from all Psi-
2087 * sharing processes.
2088 *
2089 * \param *P Problem at hand, containing Lattice and RunStruct
2090 */
2091void FillCurrentDensity(struct Problem *P)
2092{
2093 struct Lattice *Lat = &P->Lat;
2094 struct RunStruct *R = &P->R;
2095 struct Psis *Psi = &Lat->Psi;
2096 struct LatticeLevel *LevS = R->LevS;
2097 struct LatticeLevel *Lev0 = R->Lev0;
2098 struct Density *Dens0 = Lev0->Dens;
2099 fftw_complex *Psi0;
2100 fftw_real *Psi0R, *Psip0R;
2101 fftw_real *CurrentDensity[NDIM*NDIM];
2102 fftw_real *Psi1R;
2103 fftw_real *Psip1R;
2104 fftw_real *tempArray; // intendedly the same
2105 double r_bar[NDIM], x[NDIM], X[NDIM], fac[NDIM];
2106 double Current;//, current;
2107 const double UnitsFactor = 1.; ///LevS->MaxN; // 1/N (from ff-backtransform)
2108 int i, index, B_index;
2109 int k, j, i0;
2110 int n[NDIM], n0;
2111 int *N;
2112 N = Lev0->Plan0.plan->N;
2113 const int N0 = Lev0->Plan0.plan->local_nx;
2114 //int ActNum;
2115 const int myPE = P->Par.me_comm_ST_Psi;
2116 const int type = R->CurrentMin;
2117 MPI_Status status;
2118 int cross_lookup_1[4], cross_lookup_3[4], l_1 = 0, l_3 = 0;
2119 double Factor;//, factor;
2120
2121 //fprintf(stderr,"(%i) FactoR %e\n", P->Par.me, R->FactorDensityR);
2122
2123 // Init values and pointers
2124 if (P->Call.out[PsiOut]) {
2125 fprintf(stderr,"(%i) LockArray: ", P->Par.me);
2126 for(i=0;i<MaxDensityTypes;i++)
2127 fprintf(stderr,"(%i,%i) ",Dens0->LockArray[i],Dens0->LockCArray[i]);
2128 fprintf(stderr,"\n");
2129 }
2130 LockDensityArray(Dens0,Temp2Density,real); // Psi1R
2131 LockDensityArray(Dens0,Temp2Density,imag); // Psip1R and tempArray
2132 LockDensityArray(Dens0,GapDensity,real); // Psi0R
2133 LockDensityArray(Dens0,GapLocalDensity,real); // Psip0R
2134
2135 Psi0R = (fftw_real *)Dens0->DensityArray[GapDensity];
2136 Psip0R = (fftw_real *)Dens0->DensityArray[GapLocalDensity];
2137 Psi1R = (fftw_real *)Dens0->DensityArray[Temp2Density];
2138 tempArray = Psip1R = (fftw_real *)Dens0->DensityCArray[Temp2Density];
2139 SetArrayToDouble0((double *)Psi0R,Dens0->TotalSize*2);
2140 SetArrayToDouble0((double *)Psip0R,Dens0->TotalSize*2);
2141 SetArrayToDouble0((double *)Psi1R,Dens0->TotalSize*2);
2142 SetArrayToDouble0((double *)Psip1R,Dens0->TotalSize*2);
2143
2144 if (P->Call.out[PsiOut]) {
2145 fprintf(stderr,"(%i) LockArray: ", P->Par.me);
2146 for(i=0;i<MaxDensityTypes;i++)
2147 fprintf(stderr,"(%i,%i) ",Dens0->LockArray[i],Dens0->LockCArray[i]);
2148 fprintf(stderr,"\n");
2149 }
2150
2151 // don't put the following stuff into a for loop, they might not be continuous! (preprocessor values: CurrentDensity...)
2152 CurrentDensity[0] = (fftw_real *) Dens0->DensityArray[CurrentDensity0];
2153 CurrentDensity[1] = (fftw_real *) Dens0->DensityArray[CurrentDensity1];
2154 CurrentDensity[2] = (fftw_real *) Dens0->DensityArray[CurrentDensity2];
2155 CurrentDensity[3] = (fftw_real *) Dens0->DensityArray[CurrentDensity3];
2156 CurrentDensity[4] = (fftw_real *) Dens0->DensityArray[CurrentDensity4];
2157 CurrentDensity[5] = (fftw_real *) Dens0->DensityArray[CurrentDensity5];
2158 CurrentDensity[6] = (fftw_real *) Dens0->DensityArray[CurrentDensity6];
2159 CurrentDensity[7] = (fftw_real *) Dens0->DensityArray[CurrentDensity7];
2160 CurrentDensity[8] = (fftw_real *) Dens0->DensityArray[CurrentDensity8];
2161
2162 // initialize the array if it is the first of all six perturbation run
2163 if ((R->DoFullCurrent == 0) && (R->CurrentMin == Perturbed_P0)) { // reset if FillDelta...() hasn't done it before
2164 debug(P,"resetting CurrentDensity...");
2165 for (B_index=0; B_index<NDIM*NDIM; B_index++) // initialize current density array
2166 SetArrayToDouble0((double *)CurrentDensity[B_index],Dens0->TotalSize*2); // DensityArray is fftw_real, no 2*LocalSizeR here!
2167 }
2168
2169 switch(type) { // set j (which is linked to the index from derivation wrt to B^{ext})
2170 case Perturbed_P0:
2171 case Perturbed_P1:
2172 case Perturbed_P2:
2173 j = type - Perturbed_P0;
2174 l_1 = crossed(j,1);
2175 l_3 = crossed(j,3);
2176 for(k=0;k<4;k++) {
2177 cross_lookup_1[k] = cross(l_1,k);
2178 cross_lookup_3[k] = cross(l_3,k);
2179 }
2180 break;
2181 case Perturbed_RxP0:
2182 case Perturbed_RxP1:
2183 case Perturbed_RxP2:
2184 j = type - Perturbed_RxP0;
2185 break;
2186 default:
2187 j = 0;
2188 Error(SomeError,"FillCurrentDensity() called while not in perturbed minimisation!");
2189 break;
2190 }
2191
2192 int wished = -1;
2193 FILE *file = fopen(P->Call.MainParameterFile,"r");
2194 if (!ParseForParameter(0,file,"Orbital",0,1,1,int_type,&wished, 1, optional)) {
2195 if (P->Call.out[ReadOut]) fprintf(stderr,"Desired Orbital missing, using: All!\n");
2196 wished = -1;
2197 } else if (wished != -1) {
2198 if (P->Call.out[ReadOut]) fprintf(stderr,"Desired Orbital is: %i.\n", wished);
2199 } else {
2200 if (P->Call.out[ReadOut]) fprintf(stderr,"Desired Orbital is: All.\n");
2201 }
2202 fclose(file);
2203
2204 // Commence grid filling
2205 for (k=Psi->TypeStartIndex[Occupied];k<Psi->TypeStartIndex[Occupied+1];k++) // every local wave functions adds up its part of the current
2206 if ((k + P->Par.me_comm_ST_PsiT*(Psi->TypeStartIndex[UnOccupied]-Psi->TypeStartIndex[Occupied]) == wished) || (wished == -1)) { // compare with global number
2207 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i)Calculating Current Density Summand of type %s for Psi (%i/%i) ... \n", P->Par.me, R->MinimisationName[type], Psi->LocalPsiStatus[k].MyGlobalNo, k);
2208 //ActNum = k - Psi->TypeStartIndex[Occupied] + Psi->TypeStartIndex[1] * Psi->LocalPsiStatus[k].my_color_comm_ST_Psi; // global number of unperturbed Psi
2209 Psi0 = LevS->LPsi->LocalPsi[k]; // Local unperturbed psi
2210
2211 // now some preemptive ffts for the whole grid
2212 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) Bringing |Psi0> one level up and fftransforming\n", P->Par.me);
2213 fft_Psi(P, Psi0, Psi0R, 0, Psi0symmetry); //0 // 0 //0
2214
2215 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) Bringing |Psi1> one level up and fftransforming\n", P->Par.me);
2216 fft_Psi(P, LevS->LPsi->LocalPsi[Psi->TypeStartIndex[type]+k], Psi1R, 0, Psi1symmetry); //3 //0 //0
2217
2218 for (index=0;index<NDIM;index++) { // for all NDIM components of momentum operator
2219
2220 if ((P->Call.out[StepLeaderOut]) && (!index)) fprintf(stderr,"(%i) Bringing p|Psi0> one level up and fftransforming\n", P->Par.me);
2221 fft_Psi(P, Psi0, Psip0R, index, Psip0symmetry); //6 //6 //6
2222
2223 if ((P->Call.out[StepLeaderOut]) && (!index)) fprintf(stderr,"(%i) Bringing p|Psi1> one level up and fftransforming\n", P->Par.me);
2224 fft_Psi(P, LevS->LPsi->LocalPsi[Psi->TypeStartIndex[type]+k], Psip1R, index, Psip1symmetry); //4 //6 //7
2225
2226 // then for every point on the grid in real space ...
2227
2228 //if (Psi1R != (fftw_real *)Dens0->DensityArray[Temp2Density] || i0<0 || i0>=Dens0->LocalSizeR) Error(SomeError,"fft_Psi: Psi1R corrupted");
2229 //Psi1R[i0] = (Psi1_rxp_R[j])[i0] - (r_bar[cross(j,0)] * (Psi1_p_R[cross(j,1)])[i0] - r_bar[cross(j,2)] * (Psi1_p_R[cross(j,3)])[i0]); //
2230 //if (Psip1R != (fftw_real *)Dens0->DensityCArray[Temp2Density] || i0<0 || i0>=Dens0->LocalSizeR) Error(SomeError,"fft_Psi: Psip1R corrupted");
2231 //Psip1R[i0] = Psi1_rxp_pR[i0] - (r_bar[cross(j,0)] * (Psi1_p_pR[cross(j,1)])[i0] - r_bar[cross(j,2)] * (Psi1_p_pR[cross(j,3)])[i0]); //
2232
2233 switch(type) {
2234 case Perturbed_P0:
2235 case Perturbed_P1:
2236 case Perturbed_P2:
2237/* // evaluate factor to compensate r x normalized phi(r) against normalized phi(rxp)
2238 factor = 0.;
2239 for (n0=0;n0<N0;n0++) // only local points on x axis
2240 for (n[1]=0;n[1]<N[1];n[1]++)
2241 for (n[2]=0;n[2]<N[2];n[2]++) {
2242 i0 = n[2]+N[2]*(n[1]+N[1]*n0);
2243 n[0]=n0 + N0*myPE; // global relative coordinate: due to partitoning of x-axis in PEPGamma>1 case
2244 fac[0] = (double)n[0]/(double)N[0];
2245 fac[1] = (double)n[1]/(double)N[1];
2246 fac[2] = (double)n[2]/(double)N[2];
2247 RMat33Vec3(x, Lat->RealBasis, fac); // relative coordinate times basis matrix gives absolute ones
2248 MinImageConv(Lat, x, Psi->AddData[k].WannierCentre, X)
2249 for (i=0;i<NDIM;i++) // build gauge-translated r_bar evaluation point
2250 r_bar[i] = sawtooth(Lat,X,i);
2251// ShiftGaugeOrigin(P,truedist(Lat, x[i], Psi->AddData[k].WannierCentre[i], i),i);
2252 //truedist(Lat, x[i], Psi->AddData[k].WannierCentre[i], i);
2253 factor += Psi1R[i0] * (r_bar[cross_lookup_1[0]] * Psi1R[i0]);
2254 }
2255 MPI_Allreduce (&factor, &Factor, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
2256 Factor *= R->FactorDensityR; // discrete integration constant
2257 fprintf(stderr,"(%i) normalization factor of Phi^(RxP%i)_{%i} is %lg\n", P->Par.me, type, k, Factor);
2258 Factor = 1./sqrt(fabs(Factor)); //Factor/fabs(Factor) */
2259 Factor = 1.;
2260 for (n0=0;n0<N0;n0++) // only local points on x axis
2261 for (n[1]=0;n[1]<N[1];n[1]++)
2262 for (n[2]=0;n[2]<N[2];n[2]++) {
2263 i0 = n[2]+N[2]*(n[1]+N[1]*n0);
2264 n[0]=n0 + N0*myPE; // global relative coordinate: due to partitoning of x-axis in PEPGamma>1 case
2265 fac[0] = (double)n[0]/(double)N[0];
2266 fac[1] = (double)n[1]/(double)N[1];
2267 fac[2] = (double)n[2]/(double)N[2];
2268 RMat33Vec3(x, Lat->RealBasis, fac); // relative coordinate times basis matrix gives absolute ones
2269 MinImageConv(Lat, x, Psi->AddData[k].WannierCentre, X);
2270 for (i=0;i<NDIM;i++) // build gauge-translated r_bar evaluation point
2271 r_bar[i] = sawtooth(Lat,X[i],i);
2272// ShiftGaugeOrigin(P,X[i],i);
2273 //X[i];
2274 Current = Psip0R[i0] * (r_bar[cross_lookup_1[0]] * Psi1R[i0]);
2275 Current += (Psi0R[i0] * r_bar[cross_lookup_1[0]] * Psip1R[i0]);
2276 Current *= .5 * UnitsFactor * Psi->LocalPsiStatus[k].PsiFactor * R->FactorDensityR; // factor confirmed, see CalculateOneDensityR() and InitDensityCalculation()
2277 ////if (CurrentDensity[index+j*NDIM] != (fftw_real *) Dens0->DensityArray[CurrentDensity0 + index+j*NDIM] || i0<0 || i0>=Dens0->LocalSizeR || (index+j*NDIM)<0 || (index+j*NDIM)>=NDIM*NDIM) Error(SomeError,"FillCurrentDensity: CurrentDensity[] corrupted");
2278 CurrentDensity[index+l_1*NDIM][i0] -= Current; // note: sign of cross product resides in Current itself (here: plus)
2279 Current = - Psip0R[i0] * (r_bar[cross_lookup_3[2]] * Psi1R[i0]);
2280 Current += - (Psi0R[i0] * r_bar[cross_lookup_3[2]] * Psip1R[i0]);
2281 Current *= .5 * UnitsFactor * Psi->LocalPsiStatus[k].PsiFactor * R->FactorDensityR; // factor confirmed, see CalculateOneDensityR() and InitDensityCalculation()
2282 ////if (CurrentDensity[index+j*NDIM] != (fftw_real *) Dens0->DensityArray[CurrentDensity0 + index+j*NDIM] || i0<0 || i0>=Dens0->LocalSizeR || (index+j*NDIM)<0 || (index+j*NDIM)>=NDIM*NDIM) Error(SomeError,"FillCurrentDensity: CurrentDensity[] corrupted");
2283 CurrentDensity[index+l_3*NDIM][i0] -= Current; // note: sign of cross product resides in Current itself (here: minus)
2284 }
2285 break;
2286 case Perturbed_RxP0:
2287 case Perturbed_RxP1:
2288 case Perturbed_RxP2:
2289 for (n0=0;n0<N0;n0++) // only local points on x axis
2290 for (n[1]=0;n[1]<N[1];n[1]++)
2291 for (n[2]=0;n[2]<N[2];n[2]++) {
2292 i0 = n[2]+N[2]*(n[1]+N[1]*n0);
2293 Current = (Psip0R[i0] * Psi1R[i0] + Psi0R[i0] * Psip1R[i0]);
2294 Current *= .5 * UnitsFactor * Psi->LocalPsiStatus[k].PsiFactor * R->FactorDensityR; // factor confirmed, see CalculateOneDensityR() and InitDensityCalculation()
2295 ////if (CurrentDensity[index+j*NDIM] != (fftw_real *) Dens0->DensityArray[CurrentDensity0 + index+j*NDIM] || i0<0 || i0>=Dens0->LocalSizeR || (index+j*NDIM)<0 || (index+j*NDIM)>=NDIM*NDIM) Error(SomeError,"FillCurrentDensity: CurrentDensity[] corrupted");
2296 CurrentDensity[index+j*NDIM][i0] += Current;
2297 }
2298 break;
2299 default:
2300 break;
2301 }
2302 }
2303 //OutputCurrentDensity(P);
2304 }
2305
2306 //debug(P,"Unlocking arrays");
2307 //debug(P,"GapDensity");
2308 UnLockDensityArray(Dens0,GapDensity,real); // Psi0R
2309 //debug(P,"GapLocalDensity");
2310 UnLockDensityArray(Dens0,GapLocalDensity,real); // Psip0R
2311 //debug(P,"Temp2Density");
2312 UnLockDensityArray(Dens0,Temp2Density,real); // Psi1R
2313
2314// if (P->Call.out[StepLeaderOut])
2315// fprintf(stderr,"\n\n");
2316
2317 //debug(P,"MPI operation");
2318 // and in the end gather partial densities from other processes
2319 if (type == Perturbed_RxP2) // exchange all (due to shared wave functions) only after last pertubation run
2320 for (index=0;index<NDIM*NDIM;index++) {
2321 //if (tempArray != (fftw_real *)Dens0->DensityCArray[Temp2Density]) Error(SomeError,"FillCurrentDensity: tempArray corrupted");
2322 //debug(P,"tempArray to zero");
2323 SetArrayToDouble0((double *)tempArray, Dens0->TotalSize*2);
2324 ////if (CurrentDensity[index] != (fftw_real *) Dens0->DensityArray[CurrentDensity0 + index]) Error(SomeError,"FillCurrentDensity: CurrentDensity[] corrupted");
2325 //debug(P,"CurrentDensity exchange");
2326 MPI_Allreduce( CurrentDensity[index], tempArray, Dens0->LocalSizeR, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_PsiT); // gather results from all wave functions ...
2327 switch(Psi->PsiST) { // ... and also from SpinUp/Downs
2328 default:
2329 //debug(P,"CurrentDensity = tempArray");
2330 for (i=0;i<Dens0->LocalSizeR;i++) {
2331 ////if (CurrentDensity[index] != (fftw_real *) Dens0->DensityArray[CurrentDensity0 + index] || i<0 || i>=Dens0->LocalSizeR) Error(SomeError,"FillCurrentDensity: CurrentDensity[] corrupted");
2332 CurrentDensity[index][i] = tempArray[i];
2333 }
2334 break;
2335 case SpinUp:
2336 //debug(P,"CurrentDensity exchange spinup");
2337 MPI_Sendrecv(tempArray, Dens0->LocalSizeR, MPI_DOUBLE, P->Par.me_comm_ST, CurrentTag1,
2338 CurrentDensity[index], Dens0->LocalSizeR, MPI_DOUBLE, P->Par.me_comm_ST, CurrentTag2, P->Par.comm_STInter, &status );
2339 //debug(P,"CurrentDensity += tempArray");
2340 for (i=0;i<Dens0->LocalSizeR;i++) {
2341 ////if (CurrentDensity[index] != (fftw_real *) Dens0->DensityArray[CurrentDensity0 + index] || i<0 || i>=Dens0->LocalSizeR) Error(SomeError,"FillCurrentDensity: CurrentDensity[] corrupted");
2342 CurrentDensity[index][i] += tempArray[i];
2343 }
2344 break;
2345 case SpinDown:
2346 //debug(P,"CurrentDensity exchange spindown");
2347 MPI_Sendrecv(tempArray, Dens0->LocalSizeR, MPI_DOUBLE, P->Par.me_comm_ST, CurrentTag2,
2348 CurrentDensity[index], Dens0->LocalSizeR, MPI_DOUBLE, P->Par.me_comm_ST, CurrentTag1, P->Par.comm_STInter, &status );
2349 //debug(P,"CurrentDensity += tempArray");
2350 for (i=0;i<Dens0->LocalSizeR;i++) {
2351 ////if (CurrentDensity[index] != (fftw_real *) Dens0->DensityArray[CurrentDensity0 + index] || i<0 || i>=Dens0->LocalSizeR) Error(SomeError,"FillCurrentDensity: CurrentDensity[] corrupted");
2352 CurrentDensity[index][i] += tempArray[i];
2353 }
2354 break;
2355 }
2356 }
2357 //debug(P,"Temp2Density");
2358 UnLockDensityArray(Dens0,Temp2Density,imag); // Psip1R and tempArray
2359 //debug(P,"CurrentDensity end");
2360}
2361
2362/** Structure holding Problem at hand and two indices, defining the greens function to be inverted.
2363 */
2364struct params
2365{
2366 struct Problem *P;
2367 int *k;
2368 int *l;
2369 int *iter;
2370 fftw_complex *x_l;
2371};
2372
2373/** Wrapper function to solve G_kl x = b for x.
2374 * \param *x above x
2375 * \param *param additional parameters, here Problem at hand
2376 * \return evaluated to be minimized functional \f$\frac{1}{2}x \cdot Ax - xb\f$ at \a x on return
2377 */
2378static double DeltaCurrent_f(const gsl_vector * x, void * param)
2379{
2380 struct Problem *P = ((struct params *)param)->P;
2381 struct RunStruct *R = &P->R;
2382 struct LatticeLevel *LevS = R->LevS;
2383 struct Psis *Psi = &P->Lat.Psi;
2384 struct PseudoPot *PP = &P->PP;
2385 const double PsiFactor = Psi->AllPsiStatus[*((struct params *)param)->k].PsiFactor;
2386 double result = 0.;
2387 fftw_complex *TempPsi = LevS->LPsi->TempPsi;
2388 fftw_complex *TempPsi2 = LevS->LPsi->TempPsi2;
2389 int u;
2390
2391 //fprintf(stderr,"Evaluating f(%i,%i) for %i-th time\n", *((struct params *)param)->k, *((struct params *)param)->l, *((struct params *)param)->iter);
2392
2393 // extract gsl_vector
2394 for (u=0;u<LevS->MaxG;u++) {
2395 TempPsi[u].re = gsl_vector_get(x, 2*u);
2396 TempPsi[u].im = gsl_vector_get(x, 2*u+1);
2397 }
2398 // generate fnl
2399 CalculateCDfnl(P, TempPsi, PP->CDfnl); // calculate needed non-local form factors
2400 // Apply Hamiltonian to x
2401 ApplyTotalHamiltonian(P,TempPsi,TempPsi2, PP->CDfnl,PsiFactor,0);
2402 // take scalar product to get eigen value
2403 result = .5 * PsiFactor * (((*((struct params *)param)->k == *((struct params *)param)->l ? GradSP(P,LevS,TempPsi,TempPsi2) : 0.) - Psi->lambda[*((struct params *)param)->k][*((struct params *)param)->l])) - GradSP(P,LevS,TempPsi,LevS->LPsi->LocalPsi[*((struct params *)param)->l]);
2404 return result;
2405}
2406
2407/** Wrapper function to solve G_kl x = b for x.
2408 * \param *x above x
2409 * \param *param additional parameters, here Problem at hand
2410 * \param *g gradient vector on return
2411 * \return error code
2412 */
2413static void DeltaCurrent_df(const gsl_vector * x, void * param, gsl_vector * g)
2414{
2415 struct Problem *P = ((struct params *)param)->P;
2416 struct RunStruct *R = &P->R;
2417 struct LatticeLevel *LevS = R->LevS;
2418 struct Psis *Psi = &P->Lat.Psi;
2419 struct PseudoPot *PP = &P->PP;
2420 const double PsiFactor = Psi->AllPsiStatus[*((struct params *)param)->k].PsiFactor;
2421 fftw_complex *TempPsi = LevS->LPsi->TempPsi;
2422 fftw_complex *TempPsi2 = LevS->LPsi->TempPsi2;
2423 fftw_complex *x_l = ((struct params *)param)->x_l;
2424 int u;
2425
2426 //fprintf(stderr,"Evaluating df(%i,%i) for %i-th time\n", *((struct params *)param)->k, *((struct params *)param)->l, *((struct params *)param)->iter);
2427
2428 // extract gsl_vector
2429 for (u=0;u<LevS->MaxG;u++) {
2430 TempPsi[u].re = gsl_vector_get(x, 2*u);
2431 TempPsi[u].im = gsl_vector_get(x, 2*u+1);
2432 }
2433 // generate fnl
2434 CalculateCDfnl(P, TempPsi, PP->CDfnl); // calculate needed non-local form factors
2435 // Apply Hamiltonian to x
2436 ApplyTotalHamiltonian(P,TempPsi,TempPsi2, PP->CDfnl,PsiFactor,0);
2437 // put into returning vector
2438 for (u=0;u<LevS->MaxG;u++) {
2439 gsl_vector_set(g, 2*u, TempPsi2[u].re - x_l[u].re);
2440 gsl_vector_set(g, 2*u+1, TempPsi2[u].im - x_l[u].im);
2441 }
2442}
2443
2444/** Wrapper function to solve G_kl x = b for x.
2445 * \param *x above x
2446 * \param *param additional parameters, here Problem at hand
2447 * \param *f evaluated to be minimized functional \f$\frac{1}{2}x \cdot Ax - xb\f$ at \a x on return
2448 * \param *g gradient vector on return
2449 * \return error code
2450 */
2451static void DeltaCurrent_fdf(const gsl_vector * x, void * param, double * f, gsl_vector * g)
2452{
2453 struct Problem *P = ((struct params *)param)->P;
2454 struct RunStruct *R = &P->R;
2455 struct LatticeLevel *LevS = R->LevS;
2456 struct Psis *Psi = &P->Lat.Psi;
2457 struct PseudoPot *PP = &P->PP;
2458 const double PsiFactor = Psi->AllPsiStatus[*((struct params *)param)->k].PsiFactor;
2459 fftw_complex *TempPsi = LevS->LPsi->TempPsi;
2460 fftw_complex *TempPsi2 = LevS->LPsi->TempPsi2;
2461 fftw_complex *x_l = ((struct params *)param)->x_l;
2462 int u;
2463
2464 //fprintf(stderr,"Evaluating fdf(%i,%i) for %i-th time\n", *((struct params *)param)->k, *((struct params *)param)->l, *((struct params *)param)->iter);
2465
2466 // extract gsl_vector
2467 for (u=0;u<LevS->MaxG;u++) {
2468 TempPsi[u].re = gsl_vector_get(x, 2*u);
2469 TempPsi[u].im = gsl_vector_get(x, 2*u+1);
2470 }
2471 // generate fnl
2472 CalculateCDfnl(P, TempPsi, PP->CDfnl); // calculate needed non-local form factors
2473 // Apply Hamiltonian to x
2474 ApplyTotalHamiltonian(P,TempPsi,TempPsi2, PP->CDfnl,PsiFactor,0);
2475 // put into returning vector
2476 for (u=0;u<LevS->MaxG;u++) {
2477 gsl_vector_set(g, 2*u, TempPsi[u].re - x_l[u].re);
2478 gsl_vector_set(g, 2*u+1, TempPsi[u].im - x_l[u].im);
2479 }
2480
2481 *f = .5 * PsiFactor * (((*((struct params *)param)->k == *((struct params *)param)->l ? GradSP(P,LevS,TempPsi,TempPsi2) : 0.) - Psi->lambda[*((struct params *)param)->k][*((struct params *)param)->l])) - GradSP(P,LevS,TempPsi,LevS->LPsi->LocalPsi[*((struct params *)param)->l]);
2482}
2483
2484/** Evaluates the \f$\Delta j_k(r')\f$ component of the current density.
2485 * \f[
2486 * \Delta j_k(r') = \frac{e}{m} \sum_l \langle \varphi^{(0)}_k | \left ( p |r'\rangle \langle r'| + | r'\rangle\langle r'|p \right ) {\cal G}_{kl} (d_k - d_l) \times p | \varphi^{(1)}_l \rangle \cdot B
2487 * \f]
2488 * \param *P Problem at hand
2489 * \note result has not yet been MPI_Allreduced for ParallelSimulationData#comm_ST_inter or ParallelSimulationData#comm_ST_PsiT groups!
2490 * \warning the routine is checked but does not yet produce sensible results.
2491 */
2492void FillDeltaCurrentDensity(struct Problem *P)
2493{
2494 struct Lattice *Lat = &P->Lat;
2495 struct RunStruct *R = &P->R;
2496 struct Psis *Psi = &Lat->Psi;
2497 struct LatticeLevel *Lev0 = R->Lev0;
2498 struct LatticeLevel *LevS = R->LevS;
2499 struct Density *Dens0 = Lev0->Dens;
2500 int i,j,s;
2501 int k,l,u, in, dex, index,i0;
2502 //const int Num = Psi->NoOfPsis;
2503 int RecvSource;
2504 MPI_Status status;
2505 struct OnePsiElement *OnePsiB, *LOnePsiB, *OnePsiA, *LOnePsiA;
2506 const int ElementSize = (sizeof(fftw_complex) / sizeof(double));
2507 int n[NDIM], n0;
2508 int N[NDIM];
2509 N[0] = Lev0->Plan0.plan->N[0];
2510 N[1] = Lev0->Plan0.plan->N[1];
2511 N[2] = Lev0->Plan0.plan->N[2];
2512 const int N0 = Lev0->Plan0.plan->local_nx;
2513 fftw_complex *LPsiDatB;
2514 fftw_complex *Psi0, *Psi1;
2515 fftw_real *Psi0R, *Psip0R;
2516 fftw_real *Psi1R, *Psip1R;
2517 fftw_complex *x_l = LevS->LPsi->TempPsi;//, **x_l_bak;
2518 fftw_real *CurrentDensity[NDIM*NDIM];
2519 int mem_avail, MEM_avail;
2520 double Current;
2521 double X[NDIM];
2522 const double UnitsFactor = 1.;
2523 int cross_lookup[4];
2524 struct params param;
2525 double factor; // temporary factor in Psi1 pre-evaluation
2526
2527 LockDensityArray(Dens0,GapDensity,real); // Psi0R
2528 LockDensityArray(Dens0,GapLocalDensity,real); // Psip0R
2529 LockDensityArray(Dens0,Temp2Density,imag); // Psi1
2530 LockDensityArray(Dens0,GapUpDensity,real); // Psi1R
2531 LockDensityArray(Dens0,GapDownDensity,real); // Psip1R
2532
2533 CurrentDensity[0] = (fftw_real *) Dens0->DensityArray[CurrentDensity0];
2534 CurrentDensity[1] = (fftw_real *) Dens0->DensityArray[CurrentDensity1];
2535 CurrentDensity[2] = (fftw_real *) Dens0->DensityArray[CurrentDensity2];
2536 CurrentDensity[3] = (fftw_real *) Dens0->DensityArray[CurrentDensity3];
2537 CurrentDensity[4] = (fftw_real *) Dens0->DensityArray[CurrentDensity4];
2538 CurrentDensity[5] = (fftw_real *) Dens0->DensityArray[CurrentDensity5];
2539 CurrentDensity[6] = (fftw_real *) Dens0->DensityArray[CurrentDensity6];
2540 CurrentDensity[7] = (fftw_real *) Dens0->DensityArray[CurrentDensity7];
2541 CurrentDensity[8] = (fftw_real *) Dens0->DensityArray[CurrentDensity8];
2542
2543 Psi0R = (fftw_real *)Dens0->DensityArray[GapDensity];
2544 Psip0R = (fftw_real *)Dens0->DensityArray[GapLocalDensity];
2545 Psi1 = (fftw_complex *) Dens0->DensityCArray[Temp2Density];
2546 Psi1R = (fftw_real *)Dens0->DensityArray[GapUpDensity];
2547 Psip1R = (fftw_real *)Dens0->DensityArray[GapDownDensity];
2548
2549// if (R->CurrentMin == Perturbed_P0)
2550// for (B_index=0; B_index<NDIM*NDIM; B_index++) { // initialize current density array
2551// debug(P,"resetting CurrentDensity...");
2552// SetArrayToDouble0((double *)CurrentDensity[B_index],Dens0->TotalSize*2); // DensityArray is fftw_real, no 2*LocalSizeR here!
2553// }
2554 //if (Psi1 != (fftw_complex *) Dens0->DensityCArray[Temp2Density]) Error(SomeError,"FillDeltaCurrentDensity: Psi1 corrupted");
2555 SetArrayToDouble0((double *)Psi1,2*Dens0->TotalSize);
2556
2557// gsl_vector *x = gsl_vector_alloc(Num);
2558// gsl_matrix *G = gsl_matrix_alloc(Num,Num);
2559// gsl_permutation *p = gsl_permutation_alloc(Num);
2560 //int signum;
2561 // begin of GSL linearer CG solver stuff
2562 int iter, Status;
2563
2564 const gsl_multimin_fdfminimizer_type *T;
2565 gsl_multimin_fdfminimizer *minset;
2566
2567 /* Position of the minimum (1,2). */
2568 //double par[2] = { 1.0, 2.0 };
2569
2570 gsl_vector *x;
2571 gsl_multimin_function_fdf my_func;
2572
2573 param.P = P;
2574 param.k = &k;
2575 param.l = &l;
2576 param.iter = &iter;
2577 param.x_l = x_l;
2578
2579 my_func.f = &DeltaCurrent_f;
2580 my_func.df = &DeltaCurrent_df;
2581 my_func.fdf = &DeltaCurrent_fdf;
2582 my_func.n = 2*LevS->MaxG;
2583 my_func.params = (void *)&param;
2584
2585 T = gsl_multimin_fdfminimizer_conjugate_pr;
2586 minset = gsl_multimin_fdfminimizer_alloc (T, 2*LevS->MaxG);
2587 x = gsl_vector_alloc (2*LevS->MaxG);
2588 // end of GSL CG stuff
2589
2590
2591// // construct G_kl = - (H^{(0)} \delta_{kl} -\langle \varphi^{(0)}_k |H^{(0)}| \varphi^{(0)}_l|rangle)^{-1} = A^{-1}
2592// for (k=0;k<Num;k++)
2593// for (l=0;l<Num;l++)
2594// gsl_matrix_set(G, k, l, k == l ? 0. : Psi->lambda[k][l]);
2595// // and decompose G_kl = L U
2596
2597 mem_avail = MEM_avail = 0;
2598// x_l_bak = x_l = (fftw_complex **) Malloc(sizeof(fftw_complex *)*Num,"FillDeltaCurrentDensity: *x_l");
2599// for (i=0;i<Num;i++) {
2600// x_l[i] = NULL;
2601// x_l[i] = (fftw_complex *) malloc(sizeof(fftw_complex)*LevS->MaxG);
2602// if (x_l[i] == NULL) {
2603// mem_avail = 1; // there was not enough memory for this node
2604// fprintf(stderr,"(%i) FillDeltaCurrentDensity: x_l[%i] ... insufficient memory.\n",P->Par.me,i);
2605// }
2606// }
2607// MPI_Allreduce(&mem_avail,&MEM_avail,1,MPI_INT,MPI_SUM,P->Par.comm_ST); // sum results from all processes
2608
2609 if (MEM_avail != 0) { // means at least node couldn't allocate sufficient memory, skipping...
2610 fprintf(stderr,"(%i) FillDeltaCurrentDensity: x_l[], not enough memory: %i! Skipping FillDeltaCurrentDensity evaluation.", P->Par.me, MEM_avail);
2611 } else {
2612 // sum over k and calculate \Delta j_k(r')
2613 k=-1;
2614 for (i=0; i < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; i++) { // go through all wave functions
2615 //fprintf(stderr,"(%i) GlobalNo: %d\tLocalNo: %d\n", P->Par.me,Psi->AllPsiStatus[i].MyGlobalNo,Psi->AllPsiStatus[i].MyLocalNo);
2616 OnePsiA = &Psi->AllPsiStatus[i]; // grab OnePsiA
2617 if (OnePsiA->PsiType == Occupied) { // drop the extra and perturbed ones
2618 k++;
2619 if (OnePsiA->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // local?
2620 LOnePsiA = &Psi->LocalPsiStatus[OnePsiA->MyLocalNo];
2621 else
2622 LOnePsiA = NULL;
2623 if (LOnePsiA != NULL) {
2624 Psi0=LevS->LPsi->LocalPsi[OnePsiA->MyLocalNo];
2625
2626 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) Bringing |Psi0> one level up and fftransforming\n", P->Par.me);
2627 //if (Psi0R != (fftw_real *)Dens0->DensityArray[GapDensity]) Error(SomeError,"FillDeltaCurrentDensity: Psi0R corrupted");
2628 fft_Psi(P,Psi0,Psi0R, 0, Psi0symmetry); //0 // 0 //0
2629
2630 for (in=0;in<NDIM;in++) { // in is the index from derivation wrt to B^{ext}
2631 l = -1;
2632 for (j=0; j < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; j++) { // go through all wave functions
2633 OnePsiB = &Psi->AllPsiStatus[j]; // grab OnePsiA
2634 if (OnePsiB->PsiType == Occupied)
2635 l++;
2636 if ((OnePsiB != OnePsiA) && (OnePsiB->PsiType == Occupied)) { // drop the same and the extra ones
2637 if (OnePsiB->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // local?
2638 LOnePsiB = &Psi->LocalPsiStatus[OnePsiB->MyLocalNo];
2639 else
2640 LOnePsiB = NULL;
2641 if (LOnePsiB == NULL) { // if it's not local ... receive x from respective process
2642 RecvSource = OnePsiB->my_color_comm_ST_Psi;
2643 MPI_Recv( x_l, LevS->MaxG*ElementSize, MPI_DOUBLE, RecvSource, HamiltonianTag, P->Par.comm_ST_PsiT, &status );
2644 } else { // .. otherwise setup wave function as x ...
2645 // Evaluate cross product: \epsilon_{ijm} (d_k - d_l)_j p_m | \varphi^{(0)} \rangle = b_i ... and
2646 LPsiDatB=LevS->LPsi->LocalPsi[OnePsiB->MyLocalNo];
2647 //LPsiDatx=LevS->LPsi->LocalPsi[OnePsiB->MyLocalNo+Psi->TypeStartIndex[Perturbed_P0]];
2648 //CalculatePerturbationOperator_P(P,LPsiDatB,LPsiDatB_p0,cross(in,1),0);
2649 //CalculatePerturbationOperator_P(P,LPsiDatB,LPsiDatB_p1,cross(in,3),0);
2650 for (dex=0;dex<4;dex++)
2651 cross_lookup[dex] = cross(in,dex);
2652 MinImageConv(Lat,Psi->AddData[LOnePsiA->MyLocalNo].WannierCentre, Psi->AddData[LOnePsiB->MyLocalNo].WannierCentre,X);
2653 for(s=0;s<LevS->MaxG;s++) {
2654 //if (x_l != x_l_bak || s<0 || s>LevS->MaxG) Error(SomeError,"FillDeltaCurrentDensity: x_l[] corrupted");
2655 factor = (X[cross_lookup[0]] * LevS->GArray[s].G[cross_lookup[1]] - X[cross_lookup[2]] * LevS->GArray[s].G[cross_lookup[3]]);
2656 x_l[s].re = factor * (-LPsiDatB[s].im); // switched due to factorization with "-i G"
2657 x_l[s].im = factor * (LPsiDatB[s].re);
2658 }
2659 // ... and send it to all other processes (Max_me... - 1)
2660 for (u=0;u<P->Par.Max_me_comm_ST_PsiT;u++)
2661 if (u != OnePsiB->my_color_comm_ST_Psi)
2662 MPI_Send( x_l, LevS->MaxG*ElementSize, MPI_DOUBLE, u, HamiltonianTag, P->Par.comm_ST_PsiT);
2663 } // x_l row is now filled (either by receiving result or evaluating it on its own)
2664 // Solve Ax = b by minimizing 1/2 xAx -xb (gradient is residual Ax - b) with conjugate gradient polak-ribiere
2665
2666 debug(P,"fill starting point x with values from b");
2667 /* Starting point, x = b */
2668 for (u=0;u<LevS->MaxG;u++) {
2669 gsl_vector_set (x, 2*u, x_l[u].re);
2670 gsl_vector_set (x, 2*u+1, x_l[u].im);
2671 }
2672
2673 gsl_multimin_fdfminimizer_set (minset, &my_func, x, 0.01, 1e-4);
2674
2675 fprintf(stderr,"(%i) Start solving for (%i,%i) and index %i\n",P->Par.me, k,l,in);
2676 // start solving
2677 iter = 0;
2678 do
2679 {
2680 iter++;
2681 Status = gsl_multimin_fdfminimizer_iterate (minset);
2682
2683 if (Status)
2684 break;
2685
2686 Status = gsl_multimin_test_gradient (minset->gradient, 1e-3);
2687
2688 if (Status == GSL_SUCCESS)
2689 fprintf (stderr,"(%i) Minimum found after %i iterations.\n", P->Par.me, iter);
2690
2691 } while (Status == GSL_CONTINUE && iter < 100);
2692
2693 debug(P,"Put solution into Psi1");
2694 // ... and what do we do now? Put solution into Psi1!
2695 for(s=0;s<LevS->MaxG;s++) {
2696 //if (Psi1 != (fftw_complex *) Dens0->DensityCArray[Temp2Density] || s<0 || s>LevS->MaxG) Error(SomeError,"FillDeltaCurrentDensity: Psi1 corrupted");
2697 Psi1[s].re = gsl_vector_get (minset->x, 2*s);
2698 Psi1[s].im = gsl_vector_get (minset->x, 2*s+1);
2699 }
2700
2701 // // Solve A^{-1} b_i = x
2702 // for(s=0;s<LevS->MaxG;s++) {
2703 // // REAL PART
2704 // // retrieve column from gathered matrix
2705 // for(u=0;u<Num;u++)
2706 // gsl_vector_set(x,u,x_l[u][s].re);
2707 //
2708 // // solve: sum_l A_{kl}^(-1) b_l (s) = x_k (s)
2709 // gsl_linalg_LU_svx (G, p, x);
2710 //
2711 // // put solution back into x_l[s]
2712 // for(u=0;u<Num;u++) {
2713 // //if (x_l != x_l_bak || s<0 || s>=LevS->MaxG) Error(SomeError,"FillDeltaCurrentDensity: x_l[] corrupted");
2714 // x_l[u][s].re = gsl_vector_get(x,u);
2715 // }
2716 //
2717 // // IMAGINARY PART
2718 // // retrieve column from gathered matrix
2719 // for(u=0;u<Num;u++)
2720 // gsl_vector_set(x,u,x_l[u][s].im);
2721 //
2722 // // solve: sum_l A_{kl}^(-1) b_l (s) = x_k (s)
2723 // gsl_linalg_LU_svx (G, p, x);
2724 //
2725 // // put solution back into x_l[s]
2726 // for(u=0;u<Num;u++) {
2727 // //if (x_l != x_l_bak || s<0 || s>=LevS->MaxG) Error(SomeError,"FillDeltaCurrentDensity: x_l[] corrupted");
2728 // x_l[u][s].im = gsl_vector_get(x,u);
2729 // }
2730 // } // now we have in x_l a vector similar to "Psi1" which we use to evaluate the current density
2731 //
2732 // // evaluate \Delta J_k ... mind the minus sign from G_kl!
2733 // // fill Psi1
2734 // for(s=0;s<LevS->MaxG;s++) {
2735 // //if (Psi1 != (fftw_complex *) Dens0->DensityCArray[Temp2Density] || s<0 || s>LevS->MaxG) Error(SomeError,"FillDeltaCurrentDensity: Psi1 corrupted");
2736 // Psi1[s].re = x_l[k][s].re;
2737 // Psi1[s].im = x_l[k][s].im;
2738 // }
2739
2740 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) Bringing |Psi1> one level up and fftransforming\n", P->Par.me);
2741 //if (Psi1R != (fftw_real *)Dens0->DensityArray[GapUpDensity]) Error(SomeError,"FillDeltaCurrentDensity: Psi1R corrupted");
2742 fft_Psi(P,Psi1,Psi1R, 0, Psi1symmetry); //2 // 0 //0
2743
2744 for (index=0;index<NDIM;index++) { // for all NDIM components of momentum operator
2745
2746 if ((P->Call.out[StepLeaderOut]) && (!index)) fprintf(stderr,"(%i) Bringing p|Psi0> one level up and fftransforming\n", P->Par.me);
2747 //if (Psip0R != (fftw_real *)Dens0->DensityArray[GapLocalDensity]) Error(SomeError,"FillDeltaCurrentDensity: Psip0R corrupted");
2748 fft_Psi(P,Psi0,Psip0R, index, Psip0symmetry); //6 //6 //6
2749
2750 if ((P->Call.out[StepLeaderOut]) && (!index)) fprintf(stderr,"(%i) Bringing p|Psi1> one level up and fftransforming\n", P->Par.me);
2751 //if (Psip1R != (fftw_real *)Dens0->DensityArray[GapDownDensity]) Error(SomeError,"FillDeltaCurrentDensity: Psip1R corrupted");
2752 fft_Psi(P,Psi1,Psip1R, index, Psip1symmetry); //4 //6 //6
2753
2754 // then for every point on the grid in real space ...
2755 for (n0=0;n0<N0;n0++) // only local points on x axis
2756 for (n[1]=0;n[1]<N[1];n[1]++)
2757 for (n[2]=0;n[2]<N[2];n[2]++) {
2758 i0 = n[2]+N[2]*(n[1]+N[1]*n0);
2759 // and take the product
2760 Current = (Psip0R[i0] * Psi1R[i0] + Psi0R[i0] * Psip1R[i0]);
2761 Current *= 0.5 * UnitsFactor * Psi->AllPsiStatus[OnePsiA->MyGlobalNo].PsiFactor * R->FactorDensityR;
2762 ////if (CurrentDensity[index+in*NDIM] != (fftw_real *) Dens0->DensityArray[CurrentDensity0 + index+in*NDIM]) Error(SomeError,"FillCurrentDensity: CurrentDensity[] corrupted");
2763 //if (i0<0 || i0>=Dens0->LocalSizeR) Error(SomeError,"FillDeltaCurrentDensity: i0 out of range");
2764 //if ((index+in*NDIM)<0 || (index+in*NDIM)>=NDIM*NDIM) Error(SomeError,"FillDeltaCurrentDensity: index out of range");
2765 CurrentDensity[index+in*NDIM][i0] += Current; // minus sign is from G_kl
2766 }
2767 }
2768 }
2769 }
2770 }
2771 }
2772 }
2773 }
2774 }
2775 UnLockDensityArray(Dens0,GapDensity,real); // Psi0R
2776 UnLockDensityArray(Dens0,GapLocalDensity,real); // Psip0R
2777 UnLockDensityArray(Dens0,Temp2Density,imag); // Psi1
2778 UnLockDensityArray(Dens0,GapUpDensity,real); // Psi1R
2779 UnLockDensityArray(Dens0,GapDownDensity,real); // Psip1R
2780// for (i=0;i<Num;i++)
2781// if (x_l[i] != NULL) Free(x_l[i], "FillDeltaCurrentDensity: x_l[i]");
2782// Free(x_l, "FillDeltaCurrentDensity: x_l");
2783 gsl_multimin_fdfminimizer_free (minset);
2784 gsl_vector_free (x);
2785// gsl_matrix_free(G);
2786// gsl_permutation_free(p);
2787// gsl_vector_free(x);
2788}
2789
2790
2791/** Evaluates the overlap integral between \a state wave functions.
2792 * \f[
2793 * S_{kl} = \langle \varphi_k^{(1)} | \varphi_l^{(1)} \rangle
2794 * \f]
2795 * The scalar product is calculated via GradSP(), MPI_Allreduced among comm_ST_Psi and the result
2796 * stored in Psis#Overlap. The rows have to be MPI exchanged, as otherwise processes will add
2797 * to the TotalEnergy overlaps calculated with old wave functions - they have been minimised after
2798 * the product with exchanged coefficients was taken.
2799 * \param *P Problem at hand
2800 * \param l local number of perturbed wave function.
2801 * \param state PsiTypeTag minimisation state of wave functions to be overlapped
2802 */
2803void CalculateOverlap(struct Problem *P, const int l, const enum PsiTypeTag state)
2804{
2805 struct RunStruct *R = &P->R;
2806 struct Lattice *Lat = &(P->Lat);
2807 struct Psis *Psi = &Lat->Psi;
2808 struct LatticeLevel *LevS = R->LevS;
2809 struct OnePsiElement *OnePsiB, *LOnePsiB;
2810 fftw_complex *LPsiDatB=NULL, *LPsiDatA=NULL;
2811 const int ElementSize = (sizeof(fftw_complex) / sizeof(double));
2812 int RecvSource;
2813 MPI_Status status;
2814 int i,j,m,p;
2815 //const int l_normal = l - Psi->TypeStartIndex[state] + Psi->TypeStartIndex[Occupied];
2816 const int ActNum = l - Psi->TypeStartIndex[state] + Psi->TypeStartIndex[1] * Psi->LocalPsiStatus[l].my_color_comm_ST_Psi;
2817 double *sendbuf, *recvbuf;
2818 double tmp,TMP;
2819 const int gsize = P->Par.Max_me_comm_ST_PsiT; //number of processes in PsiT
2820 int p_num; // number of wave functions (for overlap)
2821
2822 // update overlap table after wave function has changed
2823 LPsiDatA = LevS->LPsi->LocalPsi[l];
2824 m = -1; // to access U matrix element (0..Num-1)
2825 for (j=0; j < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; j++) { // go through all wave functions
2826 OnePsiB = &Psi->AllPsiStatus[j]; // grab OnePsiB
2827 if (OnePsiB->PsiType == state) { // drop all but the ones of current min state
2828 m++; // increase m if it is non-extra wave function
2829 if (OnePsiB->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // local?
2830 LOnePsiB = &Psi->LocalPsiStatus[OnePsiB->MyLocalNo];
2831 else
2832 LOnePsiB = NULL;
2833 if (LOnePsiB == NULL) { // if it's not local ... receive it from respective process into TempPsi
2834 RecvSource = OnePsiB->my_color_comm_ST_Psi;
2835 MPI_Recv( LevS->LPsi->TempPsi, LevS->MaxG*ElementSize, MPI_DOUBLE, RecvSource, OverlapTag, P->Par.comm_ST_PsiT, &status );
2836 LPsiDatB=LevS->LPsi->TempPsi;
2837 } else { // .. otherwise send it to all other processes (Max_me... - 1)
2838 for (p=0;p<P->Par.Max_me_comm_ST_PsiT;p++)
2839 if (p != OnePsiB->my_color_comm_ST_Psi)
2840 MPI_Send( LevS->LPsi->LocalPsi[OnePsiB->MyLocalNo], LevS->MaxG*ElementSize, MPI_DOUBLE, p, OverlapTag, P->Par.comm_ST_PsiT);
2841 LPsiDatB=LevS->LPsi->LocalPsi[OnePsiB->MyLocalNo];
2842 } // LPsiDatB is now set to the coefficients of OnePsi either stored or MPI_Received
2843
2844 tmp = GradSP(P, LevS, LPsiDatA, LPsiDatB) * sqrt(Psi->LocalPsiStatus[l].PsiFactor * OnePsiB->PsiFactor);
2845 MPI_Allreduce ( &tmp, &TMP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
2846 //fprintf(stderr,"(%i) Setting Overlap [%i][%i] = %lg\n",P->Par.me, ActNum,m,TMP);
2847 Psi->Overlap[ActNum][m] = TMP; //= Psi->Overlap[m][ActNum]
2848 }
2849 }
2850
2851 // exchange newly calculated rows among PsiT
2852 p_num = (m+1) + 1; // number of Psis: one more due to ActNum
2853 sendbuf = (double *) Malloc(p_num * sizeof(double), "CalculateOverlap: sendbuf");
2854 sendbuf[0] = ActNum; // first entry is the global row number
2855 for (i=1;i<p_num;i++)
2856 sendbuf[i] = Psi->Overlap[ActNum][i-1]; // then follow up each entry of overlap row
2857 recvbuf = (double *) Malloc(gsize * p_num * sizeof(double), "CalculateOverlap: recvbuf");
2858 MPI_Allgather(sendbuf, p_num, MPI_DOUBLE, recvbuf, p_num, MPI_DOUBLE, P->Par.comm_ST_PsiT);
2859 Free(sendbuf, "CalculateOverlap: sendbuf");
2860 for (i=0;i<gsize;i++) {// extract results from other processes out of receiving buffer
2861 m = recvbuf[i*p_num]; // m is ActNum of the process whose results we've just received
2862 //fprintf(stderr,"(%i) Received row %i from process %i\n", P->Par.me, m, i);
2863 for (j=1;j<p_num;j++)
2864 Psi->Overlap[m][j-1] = Psi->Overlap[j-1][m] = recvbuf[i*p_num+j]; // put each entry into correspondent Overlap row
2865 }
2866 Free(recvbuf, "CalculateOverlap: recvbuf");
2867}
2868
2869
2870/** Calculates magnetic susceptibility from known current density.
2871 * The bulk susceptibility tensor can be expressed as a function of the current density.
2872 * \f[
2873 * \chi_{ij} = \frac{\mu_0}{2\Omega} \frac{\delta}{\delta B_i^{ext}} \int_\Omega d^3 r \left (r \times j(r) \right )_j
2874 * \f]
2875 * Thus the integral over real space and subsequent MPI_Allreduce() over results from ParallelSimulationData#comm_ST_Psi is
2876 * straightforward. Tensor is diagonalized afterwards and split into its various sub-tensors of lower rank (e.g., isometric
2877 * value is tensor of rank 0) which are printed to screen and the tensorial elements to file '....chi.csv'
2878 * \param *P Problem at hand
2879 */
2880void CalculateMagneticSusceptibility(struct Problem *P)
2881{
2882 struct RunStruct *R = &P->R;
2883 struct Lattice *Lat = &P->Lat;
2884 struct LatticeLevel *Lev0 = R->Lev0;
2885 struct Density *Dens0 = R->Lev0->Dens;
2886 struct Ions *I = &P->Ion;
2887 fftw_real *CurrentDensity[NDIM*NDIM];
2888 int in, dex, i, i0, n0;
2889 int n[NDIM];
2890 const int N0 = Lev0->Plan0.plan->local_nx;
2891 int N[NDIM];
2892 N[0] = Lev0->Plan0.plan->N[0];
2893 N[1] = Lev0->Plan0.plan->N[1];
2894 N[2] = Lev0->Plan0.plan->N[2];
2895 double chi[NDIM*NDIM],Chi[NDIM*NDIM], x[NDIM], X[NDIM], fac[NDIM];
2896 const double discrete_factor = Lat->Volume/Lev0->MaxN;
2897 const int myPE = P->Par.me_comm_ST_Psi;
2898 double eta, delta_chi, S, A, iso;
2899 int cross_lookup[4];
2900 char filename[256];
2901 FILE *ChiFile;
2902 time_t seconds;
2903 time(&seconds); // get current time
2904
2905 // set pointers onto current density
2906 CurrentDensity[0] = (fftw_real *) Dens0->DensityArray[CurrentDensity0];
2907 CurrentDensity[1] = (fftw_real *) Dens0->DensityArray[CurrentDensity1];
2908 CurrentDensity[2] = (fftw_real *) Dens0->DensityArray[CurrentDensity2];
2909 CurrentDensity[3] = (fftw_real *) Dens0->DensityArray[CurrentDensity3];
2910 CurrentDensity[4] = (fftw_real *) Dens0->DensityArray[CurrentDensity4];
2911 CurrentDensity[5] = (fftw_real *) Dens0->DensityArray[CurrentDensity5];
2912 CurrentDensity[6] = (fftw_real *) Dens0->DensityArray[CurrentDensity6];
2913 CurrentDensity[7] = (fftw_real *) Dens0->DensityArray[CurrentDensity7];
2914 CurrentDensity[8] = (fftw_real *) Dens0->DensityArray[CurrentDensity8];
2915 //for(i=0;i<NDIM;i++) {
2916// field[i] = Dens0->DensityArray[TempDensity+i];
2917 //LockDensityArray(Dens0,TempDensity+i,real);
2918// SetArrayToDouble0((double *)field[i],Dens0->TotalSize*2);
2919 //}
2920 gsl_matrix_complex *H = gsl_matrix_complex_calloc(NDIM,NDIM);
2921
2922
2923 if (P->Call.out[ValueOut]) fprintf(stderr,"(%i) magnetic susceptibility tensor \\Chi_ij = \n",P->Par.me);
2924 for (in=0; in<NDIM; in++) { // index i of integrand vector component
2925 for(dex=0;dex<4;dex++) // initialise cross lookup
2926 cross_lookup[dex] = cross(in,dex);
2927 for (dex=0; dex<NDIM; dex++) { // index j of derivation wrt B field
2928 chi[in+dex*NDIM] = 0.;
2929 // do the integration over real space
2930 for(n0=0;n0<N0;n0++)
2931 for(n[1]=0;n[1]<N[1];n[1]++)
2932 for(n[2]=0;n[2]<N[2];n[2]++) {
2933 n[0]=n0 + N0*myPE; // global relative coordinate: due to partitoning of x-axis in PEPGamma>1 case
2934 fac[0] = (double)(n[0])/(double)N[0];
2935 fac[1] = (double)(n[1])/(double)N[1];
2936 fac[2] = (double)(n[2])/(double)N[2];
2937 RMat33Vec3(x, Lat->RealBasis, fac);
2938 i0 = n[2]+N[2]*(n[1]+N[1]*(n0)); // the index of current density must match LocalSizeR!
2939 MinImageConv(Lat,x, Lat->RealBasisCenter, X);
2940 chi[in+dex*NDIM] += X[cross_lookup[0]] * CurrentDensity[dex*NDIM+cross_lookup[1]][i0]; // x[cross(in,0)], Lat->RealBasisCenter[cross_lookup[0]]
2941 chi[in+dex*NDIM] -= X[cross_lookup[2]] * CurrentDensity[dex*NDIM+cross_lookup[3]][i0]; // x[cross(in,2)], Lat->RealBasisCenter[cross_lookup[2]]
2942// if (in == dex) field[in][i0] =
2943// truedist(Lat,x[cross_lookup[0]], sqrt(Lat->RealBasisSQ[c[0]])/2.,cross_lookup[0]) * CurrentDensity[dex*NDIM+cross_lookup[1]][i0]
2944// - truedist(Lat,x[cross_lookup[2]], sqrt(Lat->RealBasisSQ[c[2]])/2.,cross_lookup[2]) * CurrentDensity[dex*NDIM+cross_lookup[3]][i0];
2945 //fprintf(stderr,"(%i) temporary susceptiblity \\chi[%i][%i] += %e * %e = r[%i] * CurrDens[%i][%i] = %e\n",P->Par.me,in,dex,(double)n[cross_lookup[0]]/(double)N[cross_lookup[0]]*(sqrt(Lat->RealBasisSQ[cross_lookup[0]])),CurrentDensity[dex*NDIM+cross_lookup[1]][i0],cross_lookup[0],dex*NDIM+cross_lookup[1],i0,chi[in*NDIM+dex]);
2946 }
2947 chi[in+dex*NDIM] *= mu0*discrete_factor/(2.*Lat->Volume); // integral factor
2948 chi[in+dex*NDIM] *= (-1625.); // empirical gauge factor ... sigh
2949 MPI_Allreduce ( &chi[in+dex*NDIM], &Chi[in+dex*NDIM], 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // sum "LocalSize to TotalSize"
2950 I->I[0].chi[in+dex*NDIM] = Chi[in+dex*NDIM];
2951 Chi[in+dex*NDIM] *= Lat->Volume*loschmidt_constant; // factor for _molar_ susceptibility
2952 if (P->Call.out[ValueOut]) {
2953 fprintf(stderr,"%e\t", Chi[in+dex*NDIM]);
2954 if (dex == NDIM-1) fprintf(stderr,"\n");
2955 }
2956 }
2957 }
2958 // store symmetrized matrix
2959 for (in=0;in<NDIM;in++)
2960 for (dex=0;dex<NDIM;dex++)
2961 gsl_matrix_complex_set(H,in,dex,gsl_complex_rect((Chi[in+dex*NDIM]+Chi[dex+in*NDIM])/2.,0));
2962 // output tensor to file
2963 if (P->Par.me == 0) {
2964 sprintf(filename, ".chi.L%i.csv", Lev0->LevelNo);
2965 OpenFile(P, &ChiFile, filename, "a", P->Call.out[ReadOut]);
2966 fprintf(ChiFile,"# magnetic susceptibility tensor chi[01,02,03,10,11,12,20,21,22], seed %i, config %s, run on %s", R->Seed, P->Files.default_path, ctime(&seconds));
2967 fprintf(ChiFile,"%lg\t", P->Lat.ECut/(Lat->LevelSizes[0]*Lat->LevelSizes[0]));
2968 for (in=0;in<NDIM*NDIM;in++)
2969 fprintf(ChiFile,"%e\t", Chi[in]);
2970 fprintf(ChiFile,"\n");
2971 fclose(ChiFile);
2972 }
2973 // diagonalize chi
2974 gsl_vector *eval = gsl_vector_alloc(NDIM);
2975 gsl_eigen_herm_workspace *w = gsl_eigen_herm_alloc(NDIM);
2976 gsl_eigen_herm(H, eval, w);
2977 gsl_eigen_herm_free(w);
2978 gsl_sort_vector(eval); // sort eigenvalues
2979 // print eigenvalues
2980 iso = 0;
2981 for (i=0;i<NDIM;i++) {
2982 I->I[0].chi_PAS[i] = gsl_vector_get(eval,i);
2983 iso += Chi[i+i*NDIM]/3.;
2984 }
2985 eta = (gsl_vector_get(eval,1)-gsl_vector_get(eval,0))/(gsl_vector_get(eval,2)-iso);
2986 delta_chi = gsl_vector_get(eval,2) - 0.5*(gsl_vector_get(eval,0)+gsl_vector_get(eval,1));
2987 S = (delta_chi*delta_chi)*(1+1./3.*eta*eta);
2988 A = 0.;
2989 for (i=0;i<NDIM;i++) {
2990 in = cross(i,0);
2991 dex = cross(i,1);
2992 A += pow(-1,i)*pow(0.5*(Chi[in+dex*NDIM]-Chi[dex+in*NDIM]),2);
2993 }
2994 if (P->Call.out[ValueOut]) {
2995 fprintf(stderr,"(%i) converted to Principal Axis System\n==================\nDiagonal entries:", P->Par.me);
2996 for (i=0;i<NDIM;i++)
2997 fprintf(stderr,"\t%lg",gsl_vector_get(eval,i));
2998 fprintf(stderr,"\nsusceptib. : %e\n", iso);
2999 fprintf(stderr,"anisotropy : %e\n", delta_chi);
3000 fprintf(stderr,"asymmetry : %e\n", eta);
3001 fprintf(stderr,"S : %e\n", S);
3002 fprintf(stderr,"A : %e\n", A);
3003 fprintf(stderr,"==================\n");
3004 }
3005 //for(i=0;i<NDIM;i++)
3006 //UnLockDensityArray(Dens0,TempDensity+i,real);
3007 gsl_vector_free(eval);
3008 gsl_matrix_complex_free(H);
3009}
3010
3011/** Fouriertransforms all nine current density components and calculates shielding tensor.
3012 * \f[
3013 * \sigma_{ij} = \left ( \frac{G}{|G|^2} \times J_i(G) \right )_j
3014 * \f]
3015 * The CurrentDensity has to be fouriertransformed to reciprocal subspace in order to be useful, and the final
3016 * product \f$\sigma_{ij}(G)\f$ has to be back-transformed to real space. However, the shielding is the only evaluated
3017 * at the grid points and not where the real ion position is. The shieldings there are interpolated between the eight
3018 * adjacent grid points by a simple linear weighting. Afterwards follows the same analaysis and printout of the rank-2-tensor
3019 * as in the case of CalculateMagneticShielding().
3020 * \param *P Problem at hand
3021 * \note Lots of arrays are used temporarily during the routine for the fft'ed Current density tensor.
3022 * \note MagneticSusceptibility is needed for G=0-component and thus has to be computed beforehand
3023 */
3024void CalculateChemicalShieldingByReciprocalCurrentDensity(struct Problem *P)
3025{
3026 struct RunStruct *R = &P->R;
3027 struct Lattice *Lat = &P->Lat;
3028 struct LatticeLevel *Lev0 = R->Lev0;
3029 struct Ions *I = &P->Ion;
3030 struct Density *Dens0 = Lev0->Dens;
3031 struct OneGData *GArray = Lev0->GArray;
3032 struct fft_plan_3d *plan = Lat->plan;
3033 fftw_real *CurrentDensity[NDIM*NDIM];
3034 fftw_complex *CurrentDensityC[NDIM*NDIM];
3035 fftw_complex *work = (fftw_complex *)Dens0->DensityCArray[TempDensity];
3036 //fftw_complex *sigma_imag = (fftw_complex *)Dens0->DensityCArray[Temp2Density];
3037 //fftw_real *sigma_real = (fftw_real *)sigma_imag;
3038 fftw_complex *sigma_imag[NDIM_NDIM];
3039 fftw_real *sigma_real[NDIM_NDIM];
3040 double sigma,Sigma;
3041 double x[2][NDIM];
3042 int it, ion, in, dex, g, Index, i;
3043 //const double FFTfactor = 1.;///Lev0->MaxN;
3044 int n[2][NDIM];
3045 double eta, delta_sigma, S, A, iso, tmp;
3046 FILE *SigmaFile;
3047 char suffixsigma[255];
3048 const int myPE = P->Par.me_comm_ST_Psi;
3049 int N[NDIM];
3050 int cross_lookup[4]; // cross lookup table
3051 N[0] = Lev0->Plan0.plan->N[0];
3052 N[1] = Lev0->Plan0.plan->N[1];
3053 N[2] = Lev0->Plan0.plan->N[2];
3054 const int N0 = Lev0->Plan0.plan->local_nx;
3055 const double factorDC = R->FactorDensityC;
3056 gsl_matrix_complex *H = gsl_matrix_complex_calloc(NDIM,NDIM);
3057
3058 time_t seconds;
3059 time(&seconds); // get current time
3060
3061 // inverse Fourier transform current densities
3062 CurrentDensityC[0] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity0];
3063 CurrentDensityC[1] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity1];
3064 CurrentDensityC[2] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity2];
3065 CurrentDensityC[3] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity3];
3066 CurrentDensityC[4] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity4];
3067 CurrentDensityC[5] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity5];
3068 CurrentDensityC[6] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity6];
3069 CurrentDensityC[7] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity7];
3070 CurrentDensityC[8] = (fftw_complex *) Dens0->DensityCArray[CurrentDensity8];
3071 // don't put the following stuff into a for loop, they are not continuous! (preprocessor values CurrentDensity.)
3072 CurrentDensity[0] = (fftw_real *) Dens0->DensityArray[CurrentDensity0];
3073 CurrentDensity[1] = (fftw_real *) Dens0->DensityArray[CurrentDensity1];
3074 CurrentDensity[2] = (fftw_real *) Dens0->DensityArray[CurrentDensity2];
3075 CurrentDensity[3] = (fftw_real *) Dens0->DensityArray[CurrentDensity3];
3076 CurrentDensity[4] = (fftw_real *) Dens0->DensityArray[CurrentDensity4];
3077 CurrentDensity[5] = (fftw_real *) Dens0->DensityArray[CurrentDensity5];
3078 CurrentDensity[6] = (fftw_real *) Dens0->DensityArray[CurrentDensity6];
3079 CurrentDensity[7] = (fftw_real *) Dens0->DensityArray[CurrentDensity7];
3080 CurrentDensity[8] = (fftw_real *) Dens0->DensityArray[CurrentDensity8];
3081
3082 // inverse Fourier transform current densities
3083 if (P->Call.out[LeaderOut]) fprintf(stderr,"(%i) Transforming and checking J_{ij} (G=0) = 0 for each i,j ... \n",P->Par.me);
3084 for (in=0;in<NDIM*NDIM;in++) {
3085 CalculateOneDensityC(Lat, R->LevS, Dens0, CurrentDensity[in], CurrentDensityC[in], factorDC);
3086 TestReciprocalCurrent(P, CurrentDensityC[in], GArray, in);
3087 }
3088
3089 // linking pointers to the arrays
3090 for (in=0;in<NDIM*NDIM;in++) {
3091 LockDensityArray(Dens0,in,real); // Psi1R
3092 sigma_imag[in] = (fftw_complex *) Dens0->DensityArray[in];
3093 sigma_real[in] = (fftw_real *) sigma_imag[in];
3094 }
3095
3096 LockDensityArray(Dens0,TempDensity,imag); // work
3097 LockDensityArray(Dens0,Temp2Density,imag); // tempdestRC and field
3098 // go through reciprocal nodes and calculate shielding tensor sigma
3099 for (in=0; in<NDIM; in++) {// index i of vector component in integrand
3100 for(dex=0;dex<4;dex++) // initialise cross lookup
3101 cross_lookup[dex] = cross(in,dex);
3102 for (dex=0; dex<NDIM; dex++) { // index j of B component derivation in current density tensor
3103 //if (tempdestRC != (fftw_complex *)Dens0->DensityCArray[Temp2Density]) Error(SomeError,"CalculateChemicalShieldingByReciprocalCurrentDensity: tempdestRC corrupted");
3104 SetArrayToDouble0((double *)sigma_imag[in+dex*NDIM],Dens0->TotalSize*2);
3105 for (g=0; g < Lev0->MaxG; g++)
3106 if (GArray[g].GSq > MYEPSILON) { // skip due to divisor
3107 Index = GArray[g].Index; // re = im, im = -re due to "i" in formula
3108 //if (tempdestRC != (fftw_complex *)Dens0->DensityCArray[Temp2Density] || Index<0 || Index>=Dens0->LocalSizeC) Error(SomeError,"CalculateChemicalShieldingByReciprocalCurrentDensity: tempdestRC corrupted");
3109 sigma_imag[in+dex*NDIM][Index].re = GArray[g].G[cross_lookup[0]] * (-CurrentDensityC[dex*NDIM+cross_lookup[1]][Index].im)/GArray[g].GSq;//*FFTfactor;
3110 sigma_imag[in+dex*NDIM][Index].re -= GArray[g].G[cross_lookup[2]] * (-CurrentDensityC[dex*NDIM+cross_lookup[3]][Index].im)/GArray[g].GSq;//*FFTfactor;
3111 sigma_imag[in+dex*NDIM][Index].im = GArray[g].G[cross_lookup[0]] * ( CurrentDensityC[dex*NDIM+cross_lookup[1]][Index].re)/GArray[g].GSq;//*FFTfactor;
3112 sigma_imag[in+dex*NDIM][Index].im -= GArray[g].G[cross_lookup[2]] * ( CurrentDensityC[dex*NDIM+cross_lookup[3]][Index].re)/GArray[g].GSq;//*FFTfactor;
3113 } else { // G=0-component stems from magnetic susceptibility
3114 sigma_imag[in+dex*NDIM][GArray[g].Index].re = 2./3.*I->I[0].chi[in+dex*NDIM];//-4.*M_PI*(0.5*I->I[0].chi[0+0*NDIM]+0.5*I->I[0].chi[1+1*NDIM]+2./3.*I->I[0].chi[2+2*NDIM]);
3115 }
3116 for (g=0; g<Lev0->MaxDoubleG; g++) { // apply symmetry
3117 //if (tempdestRC != (fftw_complex *)Dens0->DensityCArray[Temp2Density] || Lev0->DoubleG[2*g+1]<0 || Lev0->DoubleG[2*g+1]>=Dens0->LocalSizeC) Error(SomeError,"CalculateChemicalShieldingByReciprocalCurrentDensity: tempdestRC corrupted");
3118 sigma_imag[in+dex*NDIM][Lev0->DoubleG[2*g+1]].re = sigma_imag[in+dex*NDIM][Lev0->DoubleG[2*g]].re;
3119 sigma_imag[in+dex*NDIM][Lev0->DoubleG[2*g+1]].im = -sigma_imag[in+dex*NDIM][Lev0->DoubleG[2*g]].im;
3120 }
3121 // fourier transformation of sigma
3122 //if (tempdestRC != (fftw_complex *)Dens0->DensityCArray[Temp2Density]) Error(SomeError,"CalculateChemicalShieldingByReciprocalCurrentDensity: tempdestRC corrupted");
3123 fft_3d_complex_to_real(plan, Lev0->LevelNo, FFTNF1, sigma_imag[in+dex*NDIM], work);
3124
3125 for (it=0; it < I->Max_Types; it++) { // integration over all types
3126 for (ion=0; ion < I->I[it].Max_IonsOfType; ion++) { // and each ion of type
3127 // read transformed sigma at core position and MPI_Allreduce
3128 for (g=0;g<NDIM;g++) {
3129 n[0][g] = floor(I->I[it].R[NDIM*ion+g]/sqrt(Lat->RealBasisSQ[g])*(double)N[g]);
3130 n[1][g] = ceil(I->I[it].R[NDIM*ion+g]/sqrt(Lat->RealBasisSQ[g])*(double)N[g]);
3131 x[1][g] = I->I[it].R[NDIM*ion+g]/sqrt(Lat->RealBasisSQ[g])*(double)N[g] - (double)n[0][g];
3132 x[0][g] = 1. - x[1][g];
3133 //fprintf(stderr,"(%i) n_floor[%i] = %i\tn_ceil[%i] = %i --- x_floor[%i] = %e\tx_ceil[%i] = %e\n",P->Par.me, g,n[0][g], g,n[1][g], g,x[0][g], g,x[1][g]);
3134 }
3135 sigma = 0.;
3136 for (i=0;i<2;i++) { // interpolate linearly between adjacent grid points per axis
3137 if ((n[i][0] >= N0*myPE) && (n[i][0] < N0*(myPE+1))) {
3138// fprintf(stderr,"(%i) field[%i]: sigma = %e\n", P->Par.me, n[i][2]+N[2]*(n[i][1]+N[1]*(n[i][0]-N0*myPE)), sigma);
3139 sigma += (x[i][0]*x[0][1]*x[0][2])*sigma_real[in+dex*NDIM][n[0][2]+N[2]*(n[0][1]+N[1]*(n[i][0]-N0*myPE))]*mu0; // if it's local and factor from inverse fft
3140 //fprintf(stderr,"(%i) field[%i]: sigma += %e * %e \n", P->Par.me, n[i][2]+N[2]*(n[i][1]+N[1]*(n[i][0]-N0*myPE)), (x[i][0]*x[0][1]*x[0][2]), field[n[0][2]+N[2]*(n[0][1]+N[1]*(n[i][0]-N0*myPE))]*mu0);
3141 sigma += (x[i][0]*x[0][1]*x[1][2])*sigma_real[in+dex*NDIM][n[1][2]+N[2]*(n[0][1]+N[1]*(n[i][0]-N0*myPE))]*mu0; // if it's local and factor from inverse fft
3142 //fprintf(stderr,"(%i) field[%i]: sigma += %e * %e \n", P->Par.me, n[i][2]+N[2]*(n[i][1]+N[1]*(n[i][0]-N0*myPE)), (x[i][0]*x[0][1]*x[1][2]), field[n[1][2]+N[2]*(n[0][1]+N[1]*(n[i][0]-N0*myPE))]*mu0);
3143 sigma += (x[i][0]*x[1][1]*x[0][2])*sigma_real[in+dex*NDIM][n[0][2]+N[2]*(n[1][1]+N[1]*(n[i][0]-N0*myPE))]*mu0; // if it's local and factor from inverse fft
3144 //fprintf(stderr,"(%i) field[%i]: sigma += %e * %e \n", P->Par.me, n[i][2]+N[2]*(n[i][1]+N[1]*(n[i][0]-N0*myPE)), (x[i][0]*x[1][1]*x[0][2]), field[n[0][2]+N[2]*(n[1][1]+N[1]*(n[i][0]-N0*myPE))]*mu0);
3145 sigma += (x[i][0]*x[1][1]*x[1][2])*sigma_real[in+dex*NDIM][n[1][2]+N[2]*(n[1][1]+N[1]*(n[i][0]-N0*myPE))]*mu0; // if it's local and factor from inverse fft
3146 //fprintf(stderr,"(%i) field[%i]: sigma += %e * %e \n", P->Par.me, n[i][2]+N[2]*(n[i][1]+N[1]*(n[i][0]-N0*myPE)), (x[i][0]*x[1][1]*x[1][2]), field[n[1][2]+N[2]*(n[1][1]+N[1]*(n[i][0]-N0*myPE))]*mu0);
3147 }
3148 }
3149 sigma *= -R->FactorDensityR; // factor from inverse fft? (and its defined as negative proportionaly factor)
3150 MPI_Allreduce ( &sigma, &Sigma, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // sum local to total
3151 I->I[it].sigma_rezi[ion][in+dex*NDIM] = Sigma;
3152 }
3153 }
3154 // fabs() all sigma values, as we need them as a positive density: OutputVis plots them in logarithmic scale and
3155 // thus cannot deal with negative values!
3156 for (i=0; i< Dens0->LocalSizeR; i++)
3157 sigma_real[in+dex*NDIM][i] = fabs(sigma_real[in+dex*NDIM][i]);
3158 }
3159 }
3160 UnLockDensityArray(Dens0,TempDensity,imag); // work
3161 UnLockDensityArray(Dens0,Temp2Density,imag); // tempdestRC and field
3162
3163 // output tensor to file
3164 if (P->Par.me == 0) {
3165 sprintf(&suffixsigma[0], ".sigma_chi_rezi.L%i.csv", Lev0->LevelNo);
3166 OpenFile(P, &SigmaFile, suffixsigma, "a", P->Call.out[ReadOut]);
3167 fprintf(SigmaFile,"# chemical shielding tensor sigma_rezi[01,02,03,10,11,12,20,21,22], seed %i, config %s, run on %s", R->Seed, P->Files.default_path, ctime(&seconds));
3168 fprintf(SigmaFile,"%lg\t", P->Lat.ECut/(Lat->LevelSizes[0]*Lat->LevelSizes[0]));
3169 for (in=0;in<NDIM;in++)
3170 for (dex=0;dex<NDIM;dex++)
3171 fprintf(SigmaFile,"%e\t", GSL_REAL(gsl_matrix_complex_get(H,in,dex)));
3172 fprintf(SigmaFile,"\n");
3173 fclose(SigmaFile);
3174 }
3175
3176 gsl_vector *eval = gsl_vector_alloc(NDIM);
3177 gsl_eigen_herm_workspace *w = gsl_eigen_herm_alloc(NDIM);
3178
3179 for (it=0; it < I->Max_Types; it++) { // integration over all types
3180 for (ion=0; ion < I->I[it].Max_IonsOfType; ion++) { // and each ion of type
3181 if (P->Call.out[ValueOut]) fprintf(stderr,"(%i) Shielding Tensor for Ion %i of element %s \\sigma_ij = \n",P->Par.me, ion, I->I[it].Name);
3182 for (in=0; in<NDIM; in++) { // index i of vector component in integrand
3183 for (dex=0; dex<NDIM; dex++) {// index j of B component derivation in current density tensor
3184 gsl_matrix_complex_set(H,in,dex,gsl_complex_rect((I->I[it].sigma_rezi[ion][in+dex*NDIM]+I->I[it].sigma_rezi[ion][dex+in*NDIM])/2.,0));
3185 if (P->Call.out[ValueOut]) fprintf(stderr,"%e\t", I->I[it].sigma_rezi[ion][in+dex*NDIM]);
3186 }
3187 if (P->Call.out[ValueOut]) fprintf(stderr,"\n");
3188 }
3189 // output tensor to file
3190 if (P->Par.me == 0) {
3191 sprintf(&suffixsigma[0], ".sigma_i%i_%s_rezi.L%i.csv", ion, I->I[it].Symbol, Lev0->LevelNo);
3192 OpenFile(P, &SigmaFile, suffixsigma, "a", P->Call.out[ReadOut]);
3193 fprintf(SigmaFile,"# chemical shielding tensor sigma_rezi[01,02,03,10,11,12,20,21,22], seed %i, config %s, run on %s", R->Seed, P->Files.default_path, ctime(&seconds));
3194 fprintf(SigmaFile,"%lg\t", P->Lat.ECut/(Lat->LevelSizes[0]*Lat->LevelSizes[0]));
3195 for (in=0;in<NDIM;in++)
3196 for (dex=0;dex<NDIM;dex++)
3197 fprintf(SigmaFile,"%e\t", I->I[it].sigma_rezi[ion][in+dex*NDIM]);
3198 fprintf(SigmaFile,"\n");
3199 fclose(SigmaFile);
3200 }
3201 // diagonalize sigma
3202 gsl_eigen_herm(H, eval, w);
3203 gsl_sort_vector(eval); // sort eigenvalues
3204// print eigenvalues
3205// if (P->Call.out[ValueOut]) {
3206// fprintf(stderr,"(%i) diagonal shielding for Ion %i of element %s:", P->Par.me, ion, I->I[it].Name);
3207// for (in=0;in<NDIM;in++)
3208// fprintf(stderr,"\t%lg",gsl_vector_get(eval,in));
3209// fprintf(stderr,"\n\n");
3210// }
3211 iso = 0.;
3212 for (i=0;i<NDIM;i++) {
3213 I->I[it].sigma_rezi_PAS[ion][i] = gsl_vector_get(eval,i);
3214 iso += I->I[it].sigma_rezi[ion][i+i*NDIM]/3.;
3215 }
3216 eta = (gsl_vector_get(eval,1)-gsl_vector_get(eval,0))/(gsl_vector_get(eval,2)-iso);
3217 delta_sigma = gsl_vector_get(eval,2) - 0.5*(gsl_vector_get(eval,0)+gsl_vector_get(eval,1));
3218 S = (delta_sigma*delta_sigma)*(1+1./3.*eta*eta);
3219 A = 0.;
3220 for (i=0;i<NDIM;i++) {
3221 in = cross(i,0);
3222 dex = cross(i,1);
3223 A += pow(-1,i)*pow(0.5*(I->I[it].sigma_rezi[ion][in+dex*NDIM]-I->I[it].sigma_rezi[ion][dex+in*NDIM]),2);
3224 }
3225 if (P->Call.out[ValueOut]) {
3226 fprintf(stderr,"(%i) converted to Principal Axis System\n==================\nDiagonal entries:", P->Par.me);
3227 for (i=0;i<NDIM;i++)
3228 fprintf(stderr,"\t%lg",gsl_vector_get(eval,i));
3229 fprintf(stderr,"\nshielding : %e\n", iso);
3230 fprintf(stderr,"anisotropy : %e\n", delta_sigma);
3231 fprintf(stderr,"asymmetry : %e\n", eta);
3232 fprintf(stderr,"S : %e\n", S);
3233 fprintf(stderr,"A : %e\n", A);
3234 fprintf(stderr,"==================\n");
3235
3236 }
3237 }
3238 }
3239
3240 // Output of magnetic field densities for each direction
3241 for (i=0;i<NDIM*NDIM;i++)
3242 OutputVis(P, sigma_real[i]);
3243 // Diagonalizing the tensor "field" B_ij [r]
3244 fprintf(stderr,"(%i) Diagonalizing B_ij [r] ... \n", P->Par.me);
3245 for (i=0; i< Dens0->LocalSizeR; i++) {
3246 for (in=0; in<NDIM; in++) // index i of vector component in integrand
3247 for (dex=0; dex<NDIM; dex++) { // index j of B component derivation in current density tensor
3248 //fprintf(stderr,"(%i) Setting B_(%i,%i)[%i] ... \n", P->Par.me, in,dex,i);
3249 gsl_matrix_complex_set(H,in,dex,gsl_complex_rect((sigma_real[in+dex*NDIM][i]+sigma_real[dex+in*NDIM][i])/2.,0));
3250 }
3251 gsl_eigen_herm(H, eval, w);
3252 gsl_sort_vector(eval); // sort eigenvalues
3253 for (in=0;in<NDIM;in++)
3254 sigma_real[in][i] = gsl_vector_get(eval,in);
3255 }
3256 // Output of diagonalized magnetic field densities for each direction
3257 for (i=0;i<NDIM;i++)
3258 OutputVis(P, sigma_real[i]);
3259 for (i=0;i<NDIM*NDIM;i++)
3260 UnLockDensityArray(Dens0,i,real); // sigma_imag/real free
3261
3262 gsl_eigen_herm_free(w);
3263 gsl_vector_free(eval);
3264 gsl_matrix_complex_free(H);
3265}
3266
3267
3268/** Calculates the chemical shielding tensor at the positions of the nuclei.
3269 * The chemical shielding tensor at position R is defined as the proportionality factor between the induced and
3270 * the externally applied field.
3271 * \f[
3272 * \sigma_{ij} (R) = \frac{\delta B_j^{ind} (R)}{\delta B_i^{ext}}
3273 * = \frac{\mu_0}{4 \pi} \int d^3 r' \left ( \frac{r'-r}{| r' - r |^3} \times J_i (r') \right )_j
3274 * \f]
3275 * One after another for each nuclear position is the tensor evaluated and the result printed
3276 * to screen. Tensor is diagonalized afterwards.
3277 * \param *P Problem at hand
3278 * \sa CalculateMagneticSusceptibility() - similar calculation, yet without translation to ion centers.
3279 * \warning This routine is out-dated due to being numerically unstable because of the singularity which is not
3280 * considered carefully, recommendend replacement is CalculateChemicalShieldingByReciprocalCurrentDensity().
3281 */
3282void CalculateChemicalShielding(struct Problem *P)
3283{
3284 struct RunStruct *R = &P->R;
3285 struct Lattice *Lat = &P->Lat;
3286 struct LatticeLevel *Lev0 = R->Lev0;
3287 struct Density *Dens0 = R->Lev0->Dens;
3288 struct Ions *I = &P->Ion;
3289 double sigma[NDIM*NDIM],Sigma[NDIM*NDIM];
3290 fftw_real *CurrentDensity[NDIM*NDIM];
3291 int it, ion, in, dex, i0, n[NDIM], n0, i;//, *NUp;
3292 double r[NDIM], fac[NDIM], X[NDIM];
3293 const double discrete_factor = Lat->Volume/Lev0->MaxN;
3294 double eta, delta_sigma, S, A, iso;
3295 const int myPE = P->Par.me_comm_ST_Psi;
3296 int N[NDIM];
3297 N[0] = Lev0->Plan0.plan->N[0];
3298 N[1] = Lev0->Plan0.plan->N[1];
3299 N[2] = Lev0->Plan0.plan->N[2];
3300 const int N0 = Lev0->Plan0.plan->local_nx;
3301 FILE *SigmaFile;
3302 char suffixsigma[255];
3303 time_t seconds;
3304 time(&seconds); // get current time
3305
3306 // set pointers onto current density
3307 CurrentDensity[0] = (fftw_real *) Dens0->DensityArray[CurrentDensity0];
3308 CurrentDensity[1] = (fftw_real *) Dens0->DensityArray[CurrentDensity1];
3309 CurrentDensity[2] = (fftw_real *) Dens0->DensityArray[CurrentDensity2];
3310 CurrentDensity[3] = (fftw_real *) Dens0->DensityArray[CurrentDensity3];
3311 CurrentDensity[4] = (fftw_real *) Dens0->DensityArray[CurrentDensity4];
3312 CurrentDensity[5] = (fftw_real *) Dens0->DensityArray[CurrentDensity5];
3313 CurrentDensity[6] = (fftw_real *) Dens0->DensityArray[CurrentDensity6];
3314 CurrentDensity[7] = (fftw_real *) Dens0->DensityArray[CurrentDensity7];
3315 CurrentDensity[8] = (fftw_real *) Dens0->DensityArray[CurrentDensity8];
3316 gsl_matrix_complex *H = gsl_matrix_complex_calloc(NDIM,NDIM);
3317
3318 for (it=0; it < I->Max_Types; it++) { // integration over all types
3319 for (ion=0; ion < I->I[it].Max_IonsOfType; ion++) { // and each ion of type
3320 if (P->Call.out[ValueOut]) fprintf(stderr,"(%i) Shielding Tensor for Ion %i of element %s \\sigma_ij = \n",P->Par.me, ion, I->I[it].Name);
3321 for (in=0; in<NDIM; in++) {// index i of vector component in integrand
3322 for (dex=0; dex<NDIM; dex++) { // index j of B component derivation in current density tensor
3323 sigma[in+dex*NDIM] = 0.;
3324
3325 for(n0=0;n0<N0;n0++) // do the integration over real space
3326 for(n[1]=0;n[1]<N[1];n[1]++)
3327 for(n[2]=0;n[2]<N[2];n[2]++) {
3328 n[0]=n0 + N0*myPE; // global relative coordinate: due to partitoning of x-axis in PEPGamma>1 case
3329 fac[0] = (double)n[0]/(double)N[0];
3330 fac[1] = (double)n[1]/(double)N[1];
3331 fac[2] = (double)n[2]/(double)N[2];
3332 RMat33Vec3(r, Lat->RealBasis, fac);
3333 MinImageConv(Lat,r, &(I->I[it].R[NDIM*ion]),X);
3334 i0 = n[2]+N[2]*(n[1]+N[1]*(n0)); // the index of current density must match LocalSizeR!
3335 //z = MinImageConv(Lat,r, I->I[it].R[NDIM*ion],in); // "in" always is missing third component in cross product
3336 sigma[in+dex*NDIM] += (X[cross(in,0)] * CurrentDensity[dex*NDIM+cross(in,1)][i0] - X[cross(in,2)] * CurrentDensity[dex*NDIM+cross(in,3)][i0]);
3337 //if (it == 0 && ion == 0) fprintf(stderr,"(%i) moment[%i][%i] += (%e * %e - %e * %e) = %e\n", P->Par.me, in, dex, x,CurrentDensity[dex*NDIM+cross(in,1)][i0],y,CurrentDensity[dex*NDIM+cross(in,3)][i0],moment[in+dex*NDIM]);
3338 }
3339 sigma[in+dex*NDIM] *= -mu0*discrete_factor/(4.*PI); // due to summation instead of integration
3340 MPI_Allreduce ( &sigma[in+dex*NDIM], &Sigma[in+dex*NDIM], 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // sum "LocalSize to TotalSize"
3341 I->I[it].sigma[ion][in+dex*NDIM] = Sigma[in+dex*NDIM];
3342 if (P->Call.out[ValueOut]) fprintf(stderr," %e", Sigma[in+dex*NDIM]);
3343 }
3344 if (P->Call.out[ValueOut]) fprintf(stderr,"\n");
3345 }
3346 // store symmetrized matrix
3347 for (in=0;in<NDIM;in++)
3348 for (dex=0;dex<NDIM;dex++)
3349 gsl_matrix_complex_set(H,in,dex,gsl_complex_rect((Sigma[in+dex*NDIM]+Sigma[dex+in*NDIM])/2.,0));
3350 // output tensor to file
3351 if (P->Par.me == 0) {
3352 sprintf(&suffixsigma[0], ".sigma_i%i_%s.L%i.csv", ion, I->I[it].Symbol, Lev0->LevelNo);
3353 OpenFile(P, &SigmaFile, suffixsigma, "a", P->Call.out[ReadOut]);
3354 fprintf(SigmaFile,"# chemical shielding tensor sigma[01,02,03,10,11,12,20,21,22], seed %i, config %s, run on %s", R->Seed, P->Files.default_path, ctime(&seconds));
3355 fprintf(SigmaFile,"%lg\t", P->Lat.ECut/(Lat->LevelSizes[0]*Lat->LevelSizes[0]));
3356 for (in=0;in<NDIM*NDIM;in++)
3357 fprintf(SigmaFile,"%e\t", Sigma[in]);
3358 fprintf(SigmaFile,"\n");
3359 fclose(SigmaFile);
3360 }
3361 // diagonalize sigma
3362 gsl_vector *eval = gsl_vector_alloc(NDIM);
3363 gsl_eigen_herm_workspace *w = gsl_eigen_herm_alloc(NDIM);
3364 gsl_eigen_herm(H, eval, w);
3365 gsl_eigen_herm_free(w);
3366 gsl_sort_vector(eval); // sort eigenvalues
3367 // print eigenvalues
3368// if (P->Call.out[ValueOut]) {
3369// fprintf(stderr,"(%i) diagonal shielding for Ion %i of element %s:", P->Par.me, ion, I->I[it].Name);
3370// for (in=0;in<NDIM;in++)
3371// fprintf(stderr,"\t%lg",gsl_vector_get(eval,in));
3372// fprintf(stderr,"\n\n");
3373// }
3374 // print eigenvalues
3375 iso = 0;
3376 for (i=0;i<NDIM;i++) {
3377 I->I[it].sigma[ion][i] = gsl_vector_get(eval,i);
3378 iso += Sigma[i+i*NDIM]/3.;
3379 }
3380 eta = (gsl_vector_get(eval,1)-gsl_vector_get(eval,0))/(gsl_vector_get(eval,2)-iso);
3381 delta_sigma = gsl_vector_get(eval,2) - 0.5*(gsl_vector_get(eval,0)+gsl_vector_get(eval,1));
3382 S = (delta_sigma*delta_sigma)*(1+1./3.*eta*eta);
3383 A = 0.;
3384 for (i=0;i<NDIM;i++) {
3385 in = cross(i,0);
3386 dex = cross(i,1);
3387 A += pow(-1,i)*pow(0.5*(Sigma[in+dex*NDIM]-Sigma[dex+in*NDIM]),2);
3388 }
3389 if (P->Call.out[ValueOut]) {
3390 fprintf(stderr,"(%i) converted to Principal Axis System\n==================\nDiagonal entries:", P->Par.me);
3391 for (i=0;i<NDIM;i++)
3392 fprintf(stderr,"\t%lg",gsl_vector_get(eval,i));
3393 fprintf(stderr,"\nshielding : %e\n", iso);
3394 fprintf(stderr,"anisotropy : %e\n", delta_sigma);
3395 fprintf(stderr,"asymmetry : %e\n", eta);
3396 fprintf(stderr,"S : %e\n", S);
3397 fprintf(stderr,"A : %e\n", A);
3398 fprintf(stderr,"==================\n");
3399
3400 }
3401 gsl_vector_free(eval);
3402 }
3403 }
3404
3405 gsl_matrix_complex_free(H);
3406}
3407
3408/** Test if G=0-component of reciprocal current is 0.
3409 * In most cases we do not reach a numerical sensible zero as in MYEPSILON and remain satisfied as long
3410 * as the integrated current density is very small (e.g. compared to single entries in the current density array)
3411 * \param *P Problem at hand
3412 * \param *CurrentC pointer to reciprocal current density
3413 * \param *GArray pointer to array with G vectors
3414 * \param in index of current component
3415 * \sa TestCurrent() these two tests are equivalent and follow by fourier transformation
3416 */
3417void TestReciprocalCurrent(struct Problem *P, const fftw_complex *CurrentC, struct OneGData *GArray, int in)
3418{
3419 double tmp;
3420 tmp = sqrt(CurrentC[0].re*CurrentC[0].re+CurrentC[0].im*CurrentC[0].im);
3421 if ((P->Call.out[LeaderOut]) && (GArray[0].GSq < MYEPSILON)) {
3422 if (in % NDIM == 0) fprintf(stderr,"(%i) ",P->Par.me);
3423 if (tmp > MYEPSILON) {
3424 fprintf(stderr,"J_{%i,%i} = |%e + i%e| < %e ? (%e)\t", in / NDIM, in%NDIM, CurrentC[0].re, CurrentC[0].im, MYEPSILON, tmp - MYEPSILON);
3425 } else {
3426 fprintf(stderr,"J_{%i,%i} ok\t", in / NDIM, in%NDIM);
3427 }
3428 if (in % NDIM == (NDIM-1)) fprintf(stderr,"\n");
3429 }
3430}
3431
3432/** Test if integrated current over cell is 0.
3433 * In most cases we do not reach a numerical sensible zero as in MYEPSILON and remain satisfied as long
3434 * as the integrated current density is very small (e.g. compared to single entries in the current density array)
3435 * \param *P Problem at hand
3436 * \param index index of current component
3437 * \sa CalculateNativeIntDens() for integration of one current tensor component
3438 */
3439 void TestCurrent(struct Problem *P, const int index)
3440{
3441 struct RunStruct *R = &P->R;
3442 struct LatticeLevel *Lev0 = R->Lev0;
3443 struct Density *Dens0 = Lev0->Dens;
3444 fftw_real *CurrentDensity[NDIM*NDIM];
3445 int in;
3446 double result[NDIM*NDIM], res = 0.;
3447
3448 // set pointers onto current density array and get number of grid points in each direction
3449 CurrentDensity[0] = (fftw_real *) Dens0->DensityArray[CurrentDensity0];
3450 CurrentDensity[1] = (fftw_real *) Dens0->DensityArray[CurrentDensity1];
3451 CurrentDensity[2] = (fftw_real *) Dens0->DensityArray[CurrentDensity2];
3452 CurrentDensity[3] = (fftw_real *) Dens0->DensityArray[CurrentDensity3];
3453 CurrentDensity[4] = (fftw_real *) Dens0->DensityArray[CurrentDensity4];
3454 CurrentDensity[5] = (fftw_real *) Dens0->DensityArray[CurrentDensity5];
3455 CurrentDensity[6] = (fftw_real *) Dens0->DensityArray[CurrentDensity6];
3456 CurrentDensity[7] = (fftw_real *) Dens0->DensityArray[CurrentDensity7];
3457 CurrentDensity[8] = (fftw_real *) Dens0->DensityArray[CurrentDensity8];
3458 for(in=0;in<NDIM;in++) {
3459 result[in] = CalculateNativeIntDens(P,Lev0,CurrentDensity[in + NDIM*index],R->FactorDensityR);
3460 res += pow(result[in],2.);
3461 }
3462 res = sqrt(res);
3463 // if greater than 0, complain about it
3464 if ((res > MYEPSILON) && (P->Call.out[LeaderOut]))
3465 fprintf(stderr, "(%i) \\int_\\Omega d^3 r j_%i(r) = (%e,%e,%e), %e > %e!\n",P->Par.me, index, result[0], result[1], result[2], res, MYEPSILON);
3466}
3467
3468/** Testing whether re<->im switches (due to symmetry) confuses fft.
3469 * \param *P Problem at hand
3470 * \param l local wave function number
3471 */
3472void test_fft_symmetry(struct Problem *P, const int l)
3473{
3474 struct Lattice *Lat = &P->Lat;
3475 struct RunStruct *R = &P->R;
3476 struct LatticeLevel *LevS = R->LevS;
3477 struct LatticeLevel *Lev0 = R->Lev0;
3478 struct Density *Dens0 = Lev0->Dens;
3479 struct fft_plan_3d *plan = Lat->plan;
3480 fftw_complex *tempdestRC = (fftw_complex *)Dens0->DensityCArray[Temp2Density];
3481 fftw_complex *work = Dens0->DensityCArray[TempDensity];
3482 fftw_complex *workC = (fftw_complex *)Dens0->DensityArray[TempDensity];
3483 fftw_complex *posfac, *destpos, *destRCS, *destRCD;
3484 fftw_complex *PsiC = Dens0->DensityCArray[ActualPsiDensity];
3485 fftw_real *PsiCR = (fftw_real *) PsiC;
3486 fftw_complex *Psi0 = LevS->LPsi->LocalPsi[l];
3487 fftw_complex *dest = LevS->LPsi->TempPsi;
3488 fftw_real *Psi0R = (fftw_real *)Dens0->DensityArray[Temp2Density];
3489 int i,Index, pos, i0, iS,g; //, NoOfPsis = Psi->TypeStartIndex[UnOccupied] - Psi->TypeStartIndex[Occupied];
3490 int n[NDIM], n0;
3491 const int N0 = LevS->Plan0.plan->local_nx; // we don't want to build global density, but local
3492 int N[NDIM], NUp[NDIM];
3493 N[0] = LevS->Plan0.plan->N[0];
3494 N[1] = LevS->Plan0.plan->N[1];
3495 N[2] = LevS->Plan0.plan->N[2];
3496 NUp[0] = LevS->NUp[0];
3497 NUp[1] = LevS->NUp[1];
3498 NUp[2] = LevS->NUp[2];
3499 //const int k_normal = Lat->Psi.TypeStartIndex[Occupied] + (l - Lat->Psi.TypeStartIndex[R->CurrentMin]);
3500 //const double *Wcentre = Lat->Psi.AddData[k_normal].WannierCentre;
3501 //double x[NDIM], fac[NDIM];
3502 double result1=0., result2=0., result3=0., result4=0.;
3503 double Result1=0., Result2=0., Result3=0., Result4=0.;
3504 const double HGcRCFactor = 1./LevS->MaxN; // factor for inverse fft
3505
3506
3507 // fft to real space
3508 SetArrayToDouble0((double *)tempdestRC, Dens0->TotalSize*2);
3509 SetArrayToDouble0((double *)PsiC, Dens0->TotalSize*2);
3510 for (i=0;i<LevS->MaxG;i++) { // incoming is positive, outgoing is positive
3511 Index = LevS->GArray[i].Index;
3512 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
3513 destpos = &tempdestRC[LevS->MaxNUp*Index];
3514 for (pos=0; pos < LevS->MaxNUp; pos++) {
3515 destpos[pos].re = (Psi0[i].re)*posfac[pos].re-(Psi0[i].im)*posfac[pos].im;
3516 destpos[pos].im = (Psi0[i].re)*posfac[pos].im+(Psi0[i].im)*posfac[pos].re;
3517 //destpos[pos].re = (Psi0[i].im)*posfac[pos].re-(-Psi0[i].re)*posfac[pos].im;
3518 //destpos[pos].im = (Psi0[i].im)*posfac[pos].im+(-Psi0[i].re)*posfac[pos].re;
3519 }
3520 }
3521 for (i=0; i<LevS->MaxDoubleG; i++) {
3522 destRCS = &tempdestRC[LevS->DoubleG[2*i]*LevS->MaxNUp];
3523 destRCD = &tempdestRC[LevS->DoubleG[2*i+1]*LevS->MaxNUp];
3524 for (pos=0; pos < LevS->MaxNUp; pos++) {
3525 destRCD[pos].re = destRCS[pos].re;
3526 destRCD[pos].im = -destRCS[pos].im;
3527 }
3528 }
3529 fft_3d_complex_to_real(plan, LevS->LevelNo, FFTNFUp, tempdestRC, work);
3530 DensityRTransformPos(LevS,(fftw_real*)tempdestRC, Psi0R);
3531
3532 // apply position operator and do first result
3533 for (n0=0;n0<N0;n0++) // only local points on x axis
3534 for (n[1]=0;n[1]<N[1];n[1]++)
3535 for (n[2]=0;n[2]<N[2];n[2]++) {
3536 n[0]=n0 + LevS->Plan0.plan->start_nx; // global relative coordinate: due to partitoning of x-axis in PEPGamma>1 case
3537 i0 = n[2]*NUp[2]+N[2]*NUp[2]*(n[1]*NUp[1]+N[1]*NUp[1]*n0*NUp[0]);
3538 iS = n[2]+N[2]*(n[1]+N[1]*n0);
3539 //x[0] += 1; // shifting expectation value of x coordinate from 0 to 1
3540 PsiCR[iS] = Psi0R[i0]; // truedist(Lat, x[0], Wcentre[0],0) *
3541 result1 += PsiCR[iS] * Psi0R[i0];
3542 }
3543 result1 /= LevS->MaxN; // factor due to discrete integration
3544 MPI_Allreduce ( &result1, &Result1, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // sum "LocalSize to TotalSize"
3545 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) 1st result: %e\n",P->Par.me, Result1);
3546
3547 // fft to reciprocal space and do second result
3548 fft_3d_real_to_complex(plan, LevS->LevelNo, FFTNF1, PsiC, workC);
3549 SetArrayToDouble0((double *)dest, 2*R->InitLevS->MaxG);
3550 for (g=0; g < LevS->MaxG; g++) {
3551 Index = LevS->GArray[g].Index;
3552 dest[g].re = (Psi0[Index].re)*HGcRCFactor;
3553 dest[g].im = (Psi0[Index].im)*HGcRCFactor;
3554 }
3555 result2 = GradSP(P,LevS,Psi0,dest);
3556 MPI_Allreduce ( &result2, &Result2, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // sum "LocalSize to TotalSize"
3557 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) 2nd result: %e\n",P->Par.me, Result2);
3558
3559 // fft again to real space, this time change symmetry
3560 SetArrayToDouble0((double *)tempdestRC, Dens0->TotalSize*2);
3561 SetArrayToDouble0((double *)PsiC, Dens0->TotalSize*2);
3562 for (i=0;i<LevS->MaxG;i++) { // incoming is positive, outgoing is positive
3563 Index = LevS->GArray[i].Index;
3564 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
3565 destpos = &tempdestRC[LevS->MaxNUp*Index];
3566 for (pos=0; pos < LevS->MaxNUp; pos++) {
3567 destpos[pos].re = (Psi0[i].im)*posfac[pos].re-(-Psi0[i].re)*posfac[pos].im;
3568 destpos[pos].im = (Psi0[i].im)*posfac[pos].im+(-Psi0[i].re)*posfac[pos].re;
3569 }
3570 }
3571 for (i=0; i<LevS->MaxDoubleG; i++) {
3572 destRCS = &tempdestRC[LevS->DoubleG[2*i]*LevS->MaxNUp];
3573 destRCD = &tempdestRC[LevS->DoubleG[2*i+1]*LevS->MaxNUp];
3574 for (pos=0; pos < LevS->MaxNUp; pos++) {
3575 destRCD[pos].re = destRCS[pos].re;
3576 destRCD[pos].im = -destRCS[pos].im;
3577 }
3578 }
3579 fft_3d_complex_to_real(plan, LevS->LevelNo, FFTNFUp, tempdestRC, work);
3580 DensityRTransformPos(LevS,(fftw_real*)tempdestRC, Psi0R);
3581
3582 // bring down from Lev0 to LevS
3583 for (n0=0;n0<N0;n0++) // only local points on x axis
3584 for (n[1]=0;n[1]<N[1];n[1]++)
3585 for (n[2]=0;n[2]<N[2];n[2]++) {
3586 i0 = n[2]*NUp[2]+N[2]*NUp[2]*(n[1]*NUp[1]+N[1]*NUp[1]*n0*NUp[0]);
3587 iS = n[2]+N[2]*(n[1]+N[1]*n0);
3588 PsiCR[iS] = Psi0R[i0]; // truedist(Lat, x[0], Wcentre[0],0) *
3589 result3 += PsiCR[iS] * Psi0R[i0];
3590 }
3591 result3 /= LevS->MaxN; // factor due to discrete integration
3592 MPI_Allreduce ( &result3, &Result3, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // sum "LocalSize to TotalSize"
3593 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) 3rd result: %e\n",P->Par.me, Result3);
3594
3595 // fft back to reciprocal space, change symmetry back and do third result
3596 fft_3d_real_to_complex(plan, LevS->LevelNo, FFTNF1, PsiC, workC);
3597 SetArrayToDouble0((double *)dest, 2*R->InitLevS->MaxG);
3598 for (g=0; g < LevS->MaxG; g++) {
3599 Index = LevS->GArray[g].Index;
3600 dest[g].re = (-PsiC[Index].im)*HGcRCFactor;
3601 dest[g].im = ( PsiC[Index].re)*HGcRCFactor;
3602 }
3603 result4 = GradSP(P,LevS,Psi0,dest);
3604 MPI_Allreduce ( &result4, &Result4, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // sum "LocalSize to TotalSize"
3605 if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) 4th result: %e\n",P->Par.me, Result4);
3606}
3607
3608
3609/** Test function to check RxP application.
3610 * Checks applied solution to an analytic for a specific and simple wave function -
3611 * where just one coefficient is unequal to zero.
3612 * \param *P Problem at hand
3613 exp(I b G) - I exp(I b G) b G - exp(I a G) + I exp(I a G) a G
3614 -------------------------------------------------------------
3615 2
3616 G
3617 */
3618void test_rxp(struct Problem *P)
3619{
3620 struct RunStruct *R = &P->R;
3621 struct Lattice *Lat = &P->Lat;
3622 //struct LatticeLevel *Lev0 = R->Lev0;
3623 struct LatticeLevel *LevS = R->LevS;
3624 struct OneGData *GA = LevS->GArray;
3625 //struct Density *Dens0 = Lev0->Dens;
3626 fftw_complex *Psi0 = LevS->LPsi->TempPsi;
3627 fftw_complex *Psi2 = P->Grad.GradientArray[GraSchGradient];
3628 fftw_complex *Psi3 = LevS->LPsi->TempPsi2;
3629 int g, g_bar, i, j, k, k_normal = 0;
3630 double tmp, a,b, G;
3631 //const double *Wcentre = Lat->Psi.AddData[k_normal].WannierCentre;
3632 const double discrete_factor = 1.;//Lat->Volume/LevS->MaxN;
3633 fftw_complex integral;
3634
3635 // reset coefficients
3636 debug (P,"Creating RxP test function.");
3637 SetArrayToDouble0((double *)Psi0,2*R->InitLevS->MaxG);
3638 SetArrayToDouble0((double *)Psi2,2*R->InitLevS->MaxG);
3639
3640 // pick one which becomes non-zero
3641 g = 3;
3642
3643 //for (g=0;g<LevS->MaxG;g++) {
3644 Psi0[g].re = 1.;
3645 Psi0[g].im = 0.;
3646 //}
3647 fprintf(stderr,"(%i) G[%i] = (%e,%e,%e) \n",P->Par.me, g, GA[g].G[0], GA[g].G[1], GA[g].G[2]);
3648 i = 0;
3649
3650 // calculate analytic result
3651 debug (P,"Calculating analytic solution.");
3652 for (g_bar=0;g_bar<LevS->MaxG;g_bar++) {
3653 for (g=0;g<LevS->MaxG;g++) {
3654 if (GA[g].G[i] == GA[g_bar].G[i]) {
3655 j = cross(i,0);
3656 k = cross(i,1);
3657 if (GA[g].G[k] == GA[g_bar].G[k]) {
3658 //b = truedist(Lat, sqrt(Lat->RealBasisSQ[j]), Wcentre[j], j);
3659 b = sqrt(Lat->RealBasisSQ[j]);
3660 //a = truedist(Lat, 0., Wcentre[j], j);
3661 a = 0.;
3662 G = 1; //GA[g].G[k];
3663 if (GA[g].G[j] == GA[g_bar].G[j]) {
3664 Psi2[g_bar].re += G*Psi0[g].re * (.5 * b * b - .5 * a * a) * discrete_factor;
3665 Psi2[g_bar].im += G*Psi0[g].im * (.5 * b * b - .5 * a * a) * discrete_factor;
3666 //if ((G != 0) && ((fabs(Psi0[g].re) > MYEPSILON) || (fabs(Psi0[g].im) > MYEPSILON)))
3667 //fprintf(stderr,"(%i) Psi[%i].re += %e +i %e\n",P->Par.me, g_bar, G*Psi0[g].re * (.5 * b * b - .5 * a * a) * discrete_factor, G*Psi0[g].im * (.5 * b * b - .5 * a * a) * discrete_factor);
3668 } else {
3669 tmp = GA[g].G[j]-GA[g_bar].G[j];
3670 integral.re = (cos(tmp*b)+sin(tmp*b)*b*tmp - cos(tmp*a)-sin(tmp*a)*a*tmp) / (tmp * tmp);
3671 integral.im = (sin(tmp*b)-cos(tmp*b)*b*tmp - sin(tmp*a)+cos(tmp*a)*a*tmp) / (tmp * tmp);
3672 Psi2[g_bar].re += G*(Psi0[g].re*integral.re - Psi0[g].im*integral.im) * discrete_factor;
3673 Psi2[g_bar].im += G*(Psi0[g].re*integral.im + Psi0[g].im*integral.re) * discrete_factor;
3674 //if ((G != 0) && ((fabs(Psi0[g].re) > MYEPSILON) || (fabs(Psi0[g].im) > MYEPSILON)))
3675 //fprintf(stderr,"(%i) Psi[%i].re += %e\tPsi[%i].im += %e \n",P->Par.me, g_bar, G*(Psi0[g].re*integral.re - Psi0[g].im*integral.im) * discrete_factor, g_bar, G*(Psi0[g].re*integral.im + Psi0[g].im*integral.re) * discrete_factor);
3676 }
3677 }
3678 j = cross(i,2);
3679 k = cross(i,3);
3680 if (GA[g].G[k] == GA[g_bar].G[k]) {
3681 //b = truedist(Lat, sqrt(Lat->RealBasisSQ[j]), Wcentre[j], j);
3682 b = sqrt(Lat->RealBasisSQ[j]);
3683 //a = truedist(Lat, 0., Wcentre[j], j);
3684 a = 0.;
3685 G = 1; //GA[g].G[k];
3686 if (GA[g].G[j] == GA[g_bar].G[j]) {
3687 Psi2[g_bar].re += G*Psi0[g].re * (.5 * b * b - .5 * a * a) * discrete_factor;
3688 Psi2[g_bar].im += G*Psi0[g].im * (.5 * b * b - .5 * a * a) * discrete_factor;
3689 //if ((G != 0) && ((fabs(Psi0[g].re) > MYEPSILON) || (fabs(Psi0[g].im) > MYEPSILON)))
3690 //fprintf(stderr,"(%i) Psi[%i].re += %e +i %e\n",P->Par.me, g_bar, G*Psi0[g].re * (.5 * b * b - .5 * a * a) * discrete_factor, G*Psi0[g].im * (.5 * b * b - .5 * a * a) * discrete_factor);
3691 } else {
3692 tmp = GA[g].G[j]-GA[g_bar].G[j];
3693 integral.re = (cos(tmp*b)+sin(tmp*b)*b*tmp - cos(tmp*a)-sin(tmp*a)*a*tmp) / (tmp * tmp);
3694 integral.im = (sin(tmp*b)-cos(tmp*b)*b*tmp - sin(tmp*a)+cos(tmp*a)*a*tmp) / (tmp * tmp);
3695 Psi2[g_bar].re += G*(Psi0[g].re*integral.re - Psi0[g].im*integral.im) * discrete_factor;
3696 Psi2[g_bar].im += G*(Psi0[g].re*integral.im + Psi0[g].im*integral.re) * discrete_factor;
3697 //if ((G != 0) && ((fabs(Psi0[g].re) > MYEPSILON) || (fabs(Psi0[g].im) > MYEPSILON)))
3698 //fprintf(stderr,"(%i) Psi[%i].re += %e\tPsi[%i].im += %e \n",P->Par.me, g_bar, G*(Psi0[g].re*integral.re - Psi0[g].im*integral.im) * discrete_factor, g_bar, G*(Psi0[g].re*integral.im + Psi0[g].im*integral.re) * discrete_factor);
3699 }
3700 }
3701 }
3702 }
3703 }
3704
3705 // apply rxp
3706 debug (P,"Applying RxP to test function.");
3707 CalculatePerturbationOperator_RxP(P,Psi0,Psi3,k_normal,i);
3708
3709 // compare both coefficient arrays
3710 debug(P,"Beginning comparison of analytic and Rxp applied solution.");
3711 for (g=0;g<LevS->MaxG;g++) {
3712 if ((fabs(Psi3[g].re-Psi2[g].re) >= MYEPSILON) || (fabs(Psi3[g].im-Psi2[g].im) >= MYEPSILON))
3713 fprintf(stderr,"(%i) Psi3[%i] = %e +i %e != Psi2[%i] = %e +i %e\n",P->Par.me, g, Psi3[g].re, Psi3[g].im, g, Psi2[g].re, Psi2[g].im);
3714 //else
3715 //fprintf(stderr,"(%i) Psi1[%i] == Psi2[%i] = %e +i %e\n",P->Par.me, g, g, Psi1[g].re, Psi1[g].im);
3716 }
3717 fprintf(stderr,"(%i) <0|1> = <0|r|0> == %e +i %e\n",P->Par.me, GradSP(P,LevS,Psi0,Psi3), GradImSP(P,LevS,Psi0,Psi3));
3718 fprintf(stderr,"(%i) <1|1> = |r|ᅵ == %e +i %e\n",P->Par.me, GradSP(P,LevS,Psi3,Psi3), GradImSP(P,LevS,Psi3,Psi3));
3719 fprintf(stderr,"(%i) <0|0> = %e +i %e\n",P->Par.me, GradSP(P,LevS,Psi0,Psi0), GradImSP(P,LevS,Psi0,Psi0));
3720 fprintf(stderr,"(%i) <0|2> = %e +i %e\n",P->Par.me, GradSP(P,LevS,Psi0,Psi2), GradImSP(P,LevS,Psi0,Psi2));
3721}
3722
3723
3724/** Output of a (X,Y,DX,DY) 2d-vector plot.
3725 * For a printable representation of the induced current two-dimensional vector plots are useful, as three-dimensional
3726 * isospheres are sometimes mis-leading or do not represent the desired flow direction. The routine simply extracts a
3727 * two-dimensional cut orthogonal to one of the lattice axis at a certain node.
3728 * \param *P Problem at hand
3729 * \param B_index direction of B field
3730 * \param n_orth grid node in B_index direction of the plane (the order in which the remaining two coordinate axis
3731 * appear is the same as in a cross product, which is used to determine orthogonality)
3732 */
3733void PlotVectorPlane(struct Problem *P, int B_index, int n_orth)
3734{
3735 struct RunStruct *R = &P->R;
3736 struct LatticeLevel *Lev0 = R->Lev0;
3737 struct Density *Dens0 = Lev0->Dens;
3738 char filename[255];
3739 char *suchpointer;
3740 FILE *PlotFile = NULL;
3741 const int myPE = P->Par.me_comm_ST;
3742 time_t seconds;
3743 fftw_real *CurrentDensity[NDIM*NDIM];
3744 CurrentDensity[0] = (fftw_real *) Dens0->DensityArray[CurrentDensity0];
3745 CurrentDensity[1] = (fftw_real *) Dens0->DensityArray[CurrentDensity1];
3746 CurrentDensity[2] = (fftw_real *) Dens0->DensityArray[CurrentDensity2];
3747 CurrentDensity[3] = (fftw_real *) Dens0->DensityArray[CurrentDensity3];
3748 CurrentDensity[4] = (fftw_real *) Dens0->DensityArray[CurrentDensity4];
3749 CurrentDensity[5] = (fftw_real *) Dens0->DensityArray[CurrentDensity5];
3750 CurrentDensity[6] = (fftw_real *) Dens0->DensityArray[CurrentDensity6];
3751 CurrentDensity[7] = (fftw_real *) Dens0->DensityArray[CurrentDensity7];
3752 CurrentDensity[8] = (fftw_real *) Dens0->DensityArray[CurrentDensity8];
3753 time(&seconds); // get current time
3754
3755 if (!myPE) { // only process 0 writes to file
3756 // open file
3757 sprintf(&filename[0], ".current.L%i.csv", Lev0->LevelNo);
3758 OpenFile(P, &PlotFile, filename, "w", P->Call.out[ReadOut]);
3759 strcpy(filename, ctime(&seconds));
3760 suchpointer = strchr(filename, '\n');
3761 if (suchpointer != NULL)
3762 *suchpointer = '\0';
3763 if (PlotFile != NULL) {
3764 fprintf(PlotFile,"# current vector plot of plane perpendicular to direction e_%i at node %i, seed %i, config %s, run on %s, #cpus %i", B_index, n_orth, R->Seed, P->Files.default_path, filename, P->Par.Max_me_comm_ST_Psi);
3765 fprintf(PlotFile,"\n");
3766 } else { Error(SomeError, "PlotVectorPlane: Opening Plot File"); }
3767 }
3768
3769 // plot density
3770 if (!P->Par.me_comm_ST_PsiT) // only first wave function group as current density of all psis was gathered
3771 PlotRealDensity(P, Lev0, PlotFile, B_index, n_orth, CurrentDensity[B_index*NDIM+cross(B_index,0)], CurrentDensity[B_index*NDIM+cross(B_index,1)]);
3772
3773 if (PlotFile != NULL) {
3774 // close file
3775 fclose(PlotFile);
3776 }
3777}
3778
3779
3780/** Reads psi coefficients of \a type from file and transforms to new level.
3781 * \param *P Problem at hand
3782 * \param type PsiTypeTag of which minimisation group to load from file
3783 * \sa ReadSrcPsiDensity() - reading the coefficients, ChangePsiAndDensToLevUp() - transformation to upper level
3784 */
3785void ReadSrcPerturbedPsis(struct Problem *P, enum PsiTypeTag type)
3786{
3787 struct RunStruct *R = &P->R;
3788 struct Lattice *Lat = &P->Lat;
3789 struct LatticeLevel *Lev0 = &P->Lat.Lev[R->Lev0No+1]; // one level higher than current (ChangeLevUp already occurred)
3790 struct LatticeLevel *LevS = &P->Lat.Lev[R->LevSNo+1];
3791 struct Density *Dens = Lev0->Dens;
3792 struct Psis *Psi = &Lat->Psi;
3793 struct fft_plan_3d *plan = Lat->plan;
3794 fftw_complex *work = (fftw_complex *)Dens->DensityCArray[TempDensity];
3795 fftw_complex *tempdestRC = (fftw_complex *)Dens->DensityArray[TempDensity];
3796 fftw_complex *posfac, *destpos, *destRCS, *destRCD;
3797 fftw_complex *source, *source0;
3798 int Index,i,pos;
3799 double factorC = 1./Lev0->MaxN;
3800 int p,g;
3801
3802 // ================= read coefficients from file to LocalPsi ============
3803 ReadSrcPsiDensity(P, type, 0, R->LevSNo+1);
3804
3805 // ================= transform to upper level ===========================
3806 // for all local Psis do the usual transformation (completing coefficients for all grid vectors, fft, permutation)
3807 LockDensityArray(Dens, TempDensity, real);
3808 LockDensityArray(Dens, TempDensity, imag);
3809 for (p=Psi->LocalNo-1; p >= 0; p--)
3810 if (Psi->LocalPsiStatus[p].PsiType == type) { // only for the desired type
3811 source = LevS->LPsi->LocalPsi[p];
3812 source0 = Lev0->LPsi->LocalPsi[p];
3813 //fprintf(stderr,"(%i) ReadSrcPerturbedPsis: LevSNo %i\t Lev0No %i\tp %i\t source %p\t source0 %p\n", P->Par.me, LevS->LevelNo, Lev0->LevelNo, p, source, source0);
3814 SetArrayToDouble0((double *)tempdestRC, Dens->TotalSize*2);
3815 for (i=0;i<LevS->MaxG;i++) {
3816 Index = LevS->GArray[i].Index;
3817 posfac = &LevS->PosFactorUp[LevS->MaxNUp*i];
3818 destpos = &tempdestRC[LevS->MaxNUp*Index];
3819 //if (isnan(source[i].re)) { fprintf(stderr,"(%i) WARNING in ReadSrcPerturbedPsis(): source_%i[%i] = NaN!\n", P->Par.me, p, i); Error(SomeError, "NaN-Fehler!"); }
3820 for (pos=0; pos < LevS->MaxNUp; pos++) {
3821 destpos[pos].re = source[i].re*posfac[pos].re-source[i].im*posfac[pos].im;
3822 destpos[pos].im = source[i].re*posfac[pos].im+source[i].im*posfac[pos].re;
3823 }
3824 }
3825 for (i=0; i<LevS->MaxDoubleG; i++) {
3826 destRCS = &tempdestRC[LevS->DoubleG[2*i]*LevS->MaxNUp];
3827 destRCD = &tempdestRC[LevS->DoubleG[2*i+1]*LevS->MaxNUp];
3828 for (pos=0; pos < LevS->MaxNUp; pos++) {
3829 destRCD[pos].re = destRCS[pos].re;
3830 destRCD[pos].im = -destRCS[pos].im;
3831 }
3832 }
3833 fft_3d_complex_to_real(plan, LevS->LevelNo, FFTNFUp, tempdestRC, work);
3834 DensityRTransformPos(LevS,(fftw_real*)tempdestRC,(fftw_real *)Dens->DensityCArray[ActualPsiDensity]);
3835 // now we have density in the upper level, fft back to complex and store it as wave function coefficients
3836 fft_3d_real_to_complex(plan, Lev0->LevelNo, FFTNF1, Dens->DensityCArray[ActualPsiDensity], work);
3837 for (g=0; g < Lev0->MaxG; g++) {
3838 Index = Lev0->GArray[g].Index;
3839 source0[g].re = Dens->DensityCArray[ActualPsiDensity][Index].re*factorC;
3840 source0[g].im = Dens->DensityCArray[ActualPsiDensity][Index].im*factorC;
3841 //if (isnan(source0[g].re)) { fprintf(stderr,"(%i) WARNING in ReadSrcPerturbedPsis(): source0_%i[%i] = NaN!\n", P->Par.me, p, g); Error(SomeError, "NaN-Fehler!"); }
3842 }
3843 if (Lev0->GArray[0].GSq == 0.0)
3844 source0[g].im = 0.0;
3845 }
3846 UnLockDensityArray(Dens, TempDensity, real);
3847 UnLockDensityArray(Dens, TempDensity, imag);
3848 // finished.
3849}
3850
3851/** evaluates perturbed energy functional
3852 * \param norm norm of current Psi in functional
3853 * \param *params void-pointer to parameter array
3854 * \return evaluated functional at f(x) with \a norm
3855 */
3856double perturbed_function (double norm, void *params) {
3857 struct Problem *P = (struct Problem *)params;
3858 int i, n = P->R.LevS->MaxG;
3859 double old_norm = GramSchGetNorm2(P,P->R.LevS,P->R.LevS->LPsi->LocalPsi[P->R.ActualLocalPsiNo]);
3860 fftw_complex *currentPsi = P->R.LevS->LPsi->LocalPsi[P->R.ActualLocalPsiNo];
3861 fprintf(stderr,"(%i) perturbed_function: setting norm to %lg ...", P->Par.me, norm);
3862 // set desired norm for current Psi
3863 for (i=0; i< n; i++) {
3864 currentPsi[i].re *= norm/old_norm; // real part
3865 currentPsi[i].im *= norm/old_norm; // imaginary part
3866 }
3867 P->R.PsiStep = 0; // make it not advance to next Psi
3868
3869 //debug(P,"UpdateActualPsiNo");
3870 UpdateActualPsiNo(P, P->R.CurrentMin); // orthogonalize
3871 //debug(P,"UpdateEnergyArray");
3872 UpdateEnergyArray(P); // shift energy values in their array by one
3873 //debug(P,"UpdatePerturbedEnergyCalculation");
3874 UpdatePerturbedEnergyCalculation(P); // re-calc energies (which is hopefully lower)
3875 EnergyAllReduce(P); // gather from all processes and sum up to total energy
3876/*
3877 for (i=0; i< n; i++) {
3878 currentPsi[i].re /= norm/old_norm; // real part
3879 currentPsi[i].im /= norm/old_norm; // imaginary part
3880 }*/
3881
3882 fprintf(stderr,"%lg\n", P->Lat.E->TotalEnergy[0]);
3883 return P->Lat.E->TotalEnergy[0]; // and return evaluated functional
3884}
3885
3886/** evaluates perturbed energy functional.
3887 * \param *x current position in functional
3888 * \param *params void-pointer to parameter array
3889 * \return evaluated functional at f(x)
3890 */
3891double perturbed_f (const gsl_vector *x, void *params) {
3892 struct Problem *P = (struct Problem *)params;
3893 int i, n = P->R.LevS->MaxG*2;
3894 fftw_complex *currentPsi = P->R.LevS->LPsi->LocalPsi[P->R.ActualLocalPsiNo];
3895 //int diff = 0;
3896 //debug(P,"f");
3897 // put x into current Psi
3898 for (i=0; i< n; i+=2) {
3899 //if ((currentPsi[i/2].re != gsl_vector_get (x, i)) || (currentPsi[i/2].im != gsl_vector_get (x, i+1))) diff++;
3900 currentPsi[i/2].re = gsl_vector_get (x, i); // real part
3901 currentPsi[i/2].im = gsl_vector_get (x, i+1); // imaginary part
3902 }
3903 //if (diff) fprintf(stderr,"(%i) %i differences between old and new currentPsi.\n", P->Par.me, diff);
3904 P->R.PsiStep = 0; // make it not advance to next Psi
3905
3906 //debug(P,"UpdateActualPsiNo");
3907 UpdateActualPsiNo(P, P->R.CurrentMin); // orthogonalize
3908 //debug(P,"UpdateEnergyArray");
3909 UpdateEnergyArray(P); // shift energy values in their array by one
3910 //debug(P,"UpdatePerturbedEnergyCalculation");
3911 UpdatePerturbedEnergyCalculation(P); // re-calc energies (which is hopefully lower)
3912 EnergyAllReduce(P); // gather from all processes and sum up to total energy
3913
3914 return P->Lat.E->TotalEnergy[0]; // and return evaluated functional
3915}
3916
3917/** evaluates perturbed energy gradient.
3918 * \param *x current position in functional
3919 * \param *params void-pointer to parameter array
3920 * \param *g array for gradient vector on return
3921 */
3922void perturbed_df (const gsl_vector *x, void *params, gsl_vector *g) {
3923 struct Problem *P = (struct Problem *)params;
3924 int i, n = P->R.LevS->MaxG*2;
3925 fftw_complex *currentPsi = P->R.LevS->LPsi->LocalPsi[P->R.ActualLocalPsiNo];
3926 fftw_complex *gradient = P->Grad.GradientArray[ActualGradient];
3927 //int diff = 0;
3928 //debug(P,"df");
3929 // put x into current Psi
3930 for (i=0; i< n; i+=2) {
3931 //if ((currentPsi[i/2].re != gsl_vector_get (x, i)) || (currentPsi[i/2].im != gsl_vector_get (x, i+1))) diff++;
3932 currentPsi[i/2].re = gsl_vector_get (x, i); // real part
3933 currentPsi[i/2].im = gsl_vector_get (x, i+1); // imaginary part
3934 }
3935 //if (diff) fprintf(stderr,"(%i) %i differences between old and new currentPsi.\n", P->Par.me, diff);
3936 P->R.PsiStep = 0; // make it not advance to next Psi
3937
3938 //debug(P,"UpdateActualPsiNo");
3939 UpdateActualPsiNo(P, P->R.CurrentMin); // orthogonalize
3940 //debug(P,"UpdateEnergyArray");
3941 UpdateEnergyArray(P); // shift energy values in their array by one
3942 //debug(P,"UpdatePerturbedEnergyCalculation");
3943 UpdatePerturbedEnergyCalculation(P); // re-calc energies (which is hopefully lower)
3944 EnergyAllReduce(P); // gather from all processes and sum up to total energy
3945
3946 // checkout gradient
3947 //diff = 0;
3948 for (i=0; i< n; i+=2) {
3949 //if ((-gradient[i/2].re != gsl_vector_get (g, i)) || (-gradient[i/2].im != gsl_vector_get (g, i+1))) diff++;
3950 gsl_vector_set (g, i, -gradient[i/2].re); // real part
3951 gsl_vector_set (g, i+1, -gradient[i/2].im); // imaginary part
3952 }
3953 //if (diff) fprintf(stderr,"(%i) %i differences between old and new gradient.\n", P->Par.me, diff);
3954}
3955
3956/** evaluates perturbed energy functional and gradient.
3957 * \param *x current position in functional
3958 * \param *params void-pointer to parameter array
3959 * \param *f pointer to energy function value on return
3960 * \param *g array for gradient vector on return
3961 */
3962void perturbed_fdf (const gsl_vector *x, void *params, double *f, gsl_vector *g) {
3963 struct Problem *P = (struct Problem *)params;
3964 int i, n = P->R.LevS->MaxG*2;
3965 fftw_complex *currentPsi = P->R.LevS->LPsi->LocalPsi[P->R.ActualLocalPsiNo];
3966 fftw_complex *gradient = P->Grad.GradientArray[ActualGradient];
3967 //int diff = 0;
3968 //debug(P,"fdf");
3969 // put x into current Psi
3970 for (i=0; i< n; i+=2) {
3971 //if ((currentPsi[i/2].re != gsl_vector_get (x, i)) || (currentPsi[i/2].im != gsl_vector_get (x, i+1))) diff++;
3972 currentPsi[i/2].re = gsl_vector_get (x, i); // real part
3973 currentPsi[i/2].im = gsl_vector_get (x, i+1); // imaginary part
3974 }
3975 //if (diff) fprintf(stderr,"(%i) %i differences between old and new currentPsi.\n", P->Par.me, diff);
3976 P->R.PsiStep = 0; // make it not advance to next Psi
3977
3978 //debug(P,"UpdateActualPsiNo");
3979 UpdateActualPsiNo(P, P->R.CurrentMin); // orthogonalize
3980 //debug(P,"UpdateEnergyArray");
3981 UpdateEnergyArray(P); // shift energy values in their array by one
3982 //debug(P,"UpdatePerturbedEnergyCalculation");
3983 UpdatePerturbedEnergyCalculation(P); // re-calc energies (which is hopefully lower)
3984 EnergyAllReduce(P); // gather from all processes and sum up to total energy
3985
3986 // checkout gradient
3987 //diff = 0;
3988 for (i=0; i< n; i+=2) {
3989 //if ((-gradient[i/2].re != gsl_vector_get (g, i)) || (-gradient[i/2].im != gsl_vector_get (g, i+1))) diff++;
3990 gsl_vector_set (g, i, -gradient[i/2].re); // real part
3991 gsl_vector_set (g, i+1, -gradient[i/2].im); // imaginary part
3992 }
3993 //if (diff) fprintf(stderr,"(%i) %i differences between old and new gradient.\n", P->Par.me, diff);
3994
3995 *f = P->Lat.E->TotalEnergy[0]; // and return evaluated functional
3996}
3997
3998/* MinimisePerturbed with all the brent minimisation approach
3999void MinimisePerturbed (struct Problem *P, int *Stop, int *SuperStop) {
4000 struct RunStruct *R = &P->R;
4001 struct Lattice *Lat = &P->Lat;
4002 struct Psis *Psi = &Lat->Psi;
4003 int type;
4004 //int i;
4005
4006 // stuff for GSL minimization
4007 //size_t iter;
4008 //int status, Status
4009 int n = R->LevS->MaxG*2;
4010 const gsl_multimin_fdfminimizer_type *T_multi;
4011 const gsl_min_fminimizer_type *T;
4012 gsl_multimin_fdfminimizer *s_multi;
4013 gsl_min_fminimizer *s;
4014 gsl_vector *x;//, *ss;
4015 gsl_multimin_function_fdf my_func;
4016 gsl_function F;
4017 //fftw_complex *currentPsi;
4018 //double a,b,m, f_m, f_a, f_b;
4019 //double old_norm;
4020
4021 my_func.f = &perturbed_f;
4022 my_func.df = &perturbed_df;
4023 my_func.fdf = &perturbed_fdf;
4024 my_func.n = n;
4025 my_func.params = P;
4026 F.function = &perturbed_function;
4027 F.params = P;
4028
4029 x = gsl_vector_alloc (n);
4030 //ss = gsl_vector_alloc (Psi->NoOfPsis);
4031 T_multi = gsl_multimin_fdfminimizer_vector_bfgs;
4032 s_multi = gsl_multimin_fdfminimizer_alloc (T_multi, n);
4033 T = gsl_min_fminimizer_brent;
4034 s = gsl_min_fminimizer_alloc (T);
4035
4036 for (type=Perturbed_P0;type<=Perturbed_RxP2;type++) { // go through each perturbation group separately //
4037 *Stop=0; // reset stop flag
4038 fprintf(stderr,"(%i)Beginning perturbed minimisation of type %s ...\n", P->Par.me, R->MinimisationName[type]);
4039 //OutputOrbitalPositions(P, Occupied);
4040 R->PsiStep = R->MaxPsiStep; // reset in-Psi-minimisation-counter, so that we really advance to the next wave function
4041 UpdateActualPsiNo(P, type); // step on to next perturbed one
4042 fprintf(stderr, "(%i) Re-initializing perturbed psi array for type %s ", P->Par.me, R->MinimisationName[type]);
4043 if (P->Call.ReadSrcFiles && ReadSrcPsiDensity(P,type,1, R->LevSNo)) {
4044 SpeedMeasure(P, InitSimTime, StartTimeDo);
4045 fprintf(stderr,"from source file of recent calculation\n");
4046 ReadSrcPsiDensity(P,type, 0, R->LevSNo);
4047 ResetGramSchTagType(P, Psi, type, IsOrthogonal); // loaded values are orthonormal
4048 SpeedMeasure(P, DensityTime, StartTimeDo);
4049 //InitDensityCalculation(P);
4050 SpeedMeasure(P, DensityTime, StopTimeDo);
4051 R->OldActualLocalPsiNo = R->ActualLocalPsiNo; // needed otherwise called routines in function below crash
4052 UpdateGramSchOldActualPsiNo(P,Psi);
4053 InitPerturbedEnergyCalculation(P, 1); // go through all orbitals calculate each H^{(0)}-eigenvalue, recalc HGDensity, cause InitDensityCalc zero'd it
4054 UpdatePerturbedEnergyCalculation(P); // H1cGradient and Gradient must be current ones
4055 EnergyAllReduce(P); // gather energies for minimum search
4056 SpeedMeasure(P, InitSimTime, StopTimeDo);
4057 }
4058 if (P->Call.ReadSrcFiles != 1) {
4059 SpeedMeasure(P, InitSimTime, StartTimeDo);
4060 ResetGramSchTagType(P, Psi, type, NotOrthogonal); // perturbed now shall be orthonormalized
4061 if (P->Call.ReadSrcFiles != 2) {
4062 if (R->LevSNo == Lat->MaxLevel-1) { // is it the starting level? (see InitRunLevel())
4063 fprintf(stderr, "randomly.\n");
4064 InitPsisValue(P, Psi->TypeStartIndex[type], Psi->TypeStartIndex[type+1]); // initialize perturbed array for this run
4065 } else {
4066 fprintf(stderr, "from source file of last level.\n");
4067 ReadSrcPerturbedPsis(P, type);
4068 }
4069 }
4070 SpeedMeasure(P, InitGramSchTime, StartTimeDo);
4071 GramSch(P, R->LevS, Psi, Orthogonalize);
4072 SpeedMeasure(P, InitGramSchTime, StopTimeDo);
4073 SpeedMeasure(P, InitDensityTime, StartTimeDo);
4074 //InitDensityCalculation(P);
4075 SpeedMeasure(P, InitDensityTime, StopTimeDo);
4076 InitPerturbedEnergyCalculation(P, 1); // go through all orbitals calculate each H^{(0)}-eigenvalue, recalc HGDensity, cause InitDensityCalc zero'd it
4077 R->OldActualLocalPsiNo = R->ActualLocalPsiNo; // needed otherwise called routines in function below crash
4078 UpdateGramSchOldActualPsiNo(P,Psi);
4079 UpdatePerturbedEnergyCalculation(P); // H1cGradient and Gradient must be current ones
4080 EnergyAllReduce(P); // gather energies for minimum search
4081 SpeedMeasure(P, InitSimTime, StopTimeDo);
4082 R->LevS->Step++;
4083 EnergyOutput(P,0);
4084 while (*Stop != 1) {
4085 // copy current Psi into starting vector
4086 currentPsi = R->LevS->LPsi->LocalPsi[R->ActualLocalPsiNo];
4087 for (i=0; i< n; i+=2) {
4088 gsl_vector_set (x, i, currentPsi[i/2].re); // real part
4089 gsl_vector_set (x, i+1, currentPsi[i/2].im); // imaginary part
4090 }
4091 gsl_multimin_fdfminimizer_set (s_multi, &my_func, x, 0.01, 1e-2);
4092 iter = 0;
4093 status = 0;
4094 do { // look for minimum along current local psi
4095 iter++;
4096 status = gsl_multimin_fdfminimizer_iterate (s_multi);
4097 MPI_Allreduce(&status, &Status, 1, MPI_INT, MPI_MAX, P->Par.comm_ST_Psi);
4098 if (Status)
4099 break;
4100 status = gsl_multimin_test_gradient (s_multi->gradient, 1e-2);
4101 MPI_Allreduce(&status, &Status, 1, MPI_INT, MPI_MAX, P->Par.comm_ST_Psi);
4102 //if (Status == GSL_SUCCESS)
4103 //printf ("Minimum found at:\n");
4104 if (P->Par.me == 0) fprintf (stderr,"(%i,%i,%i)S(%i,%i,%i):\t %5d %10.5f\n",P->Par.my_color_comm_ST,P->Par.me_comm_ST, P->Par.me_comm_ST_PsiT, R->MinStep, R->ActualLocalPsiNo, R->PsiStep, (int)iter, s_multi->f);
4105 //TestGramSch(P,R->LevS,Psi, type); // functions are orthonormal?
4106 } while (Status == GSL_CONTINUE && iter < 3);
4107 // now minimize norm of currentPsi (one-dim)
4108 if (0) {
4109 iter = 0;
4110 status = 0;
4111 m = 1.;
4112 a = MYEPSILON;
4113 b = 100.;
4114 f_a = perturbed_function (a, P);
4115 f_b = perturbed_function (b, P);
4116 f_m = perturbed_function (m, P);
4117 //if ((f_m < f_a) && (f_m < f_b)) {
4118 gsl_min_fminimizer_set (s, &F, m, a, b);
4119 do { // look for minimum along current local psi
4120 iter++;
4121 status = gsl_min_fminimizer_iterate (s);
4122 m = gsl_min_fminimizer_x_minimum (s);
4123 a = gsl_min_fminimizer_x_lower (s);
4124 b = gsl_min_fminimizer_x_upper (s);
4125 status = gsl_min_test_interval (a, b, 0.001, 0.0);
4126 if (status == GSL_SUCCESS)
4127 printf ("Minimum found at:\n");
4128 printf ("%5d [%.7f, %.7f] %.7f %.7f\n",
4129 (int) iter, a, b,
4130 m, b - a);
4131 } while (status == GSL_CONTINUE && iter < 100);
4132 old_norm = GramSchGetNorm2(P,P->R.LevS,P->R.LevS->LPsi->LocalPsi[P->R.ActualLocalPsiNo]);
4133 for (i=0; i< n; i++) {
4134 currentPsi[i].re *= m/old_norm; // real part
4135 currentPsi[i].im *= m/old_norm; // imaginary part
4136 }
4137 } else debug(P,"Norm not minimizable!");
4138 //P->R.PsiStep = P->R.MaxPsiStep; // make it advance to next Psi
4139 FindPerturbedMinimum(P);
4140 //debug(P,"UpdateActualPsiNo");
4141 UpdateActualPsiNo(P, type); // step on to next perturbed Psi
4142 //debug(P,"UpdateEnergyArray");
4143 UpdateEnergyArray(P); // shift energy values in their array by one
4144 //debug(P,"UpdatePerturbedEnergyCalculation");
4145 UpdatePerturbedEnergyCalculation(P); // re-calc energies (which is hopefully lower)
4146 EnergyAllReduce(P); // gather from all processes and sum up to total energy
4147 //ControlNativeDensity(P); // check total density (summed up PertMixed must be zero!)
4148 //printf ("(%i,%i,%i)S(%i,%i,%i):\t %5d %10.5f\n",P->Par.my_color_comm_ST,P->Par.me_comm_ST, P->Par.me_comm_ST_PsiT, R->MinStep, R->ActualLocalPsiNo, R->PsiStep, (int)iter, s_multi->f);
4149 if (*SuperStop != 1)
4150 *SuperStop = CheckCPULIM(P);
4151 *Stop = CalculateMinimumStop(P, *SuperStop);
4152 P->Speed.Steps++; // step on
4153 R->LevS->Step++;
4154 }
4155 // now release normalization condition and minimize wrt to norm
4156 *Stop = 0;
4157 while (*Stop != 1) {
4158 currentPsi = R->LevS->LPsi->LocalPsi[R->ActualLocalPsiNo];
4159 iter = 0;
4160 status = 0;
4161 m = 1.;
4162 a = 0.001;
4163 b = 10.;
4164 f_a = perturbed_function (a, P);
4165 f_b = perturbed_function (b, P);
4166 f_m = perturbed_function (m, P);
4167 if ((f_m < f_a) && (f_m < f_b)) {
4168 gsl_min_fminimizer_set (s, &F, m, a, b);
4169 do { // look for minimum along current local psi
4170 iter++;
4171 status = gsl_min_fminimizer_iterate (s);
4172 m = gsl_min_fminimizer_x_minimum (s);
4173 a = gsl_min_fminimizer_x_lower (s);
4174 b = gsl_min_fminimizer_x_upper (s);
4175 status = gsl_min_test_interval (a, b, 0.001, 0.0);
4176 if (status == GSL_SUCCESS)
4177 printf ("Minimum found at:\n");
4178 printf ("%5d [%.7f, %.7f] %.7f %.7f\n",
4179 (int) iter, a, b,
4180 m, b - a);
4181 } while (status == GSL_CONTINUE && iter < 100);
4182 old_norm = GramSchGetNorm2(P,P->R.LevS,P->R.LevS->LPsi->LocalPsi[P->R.ActualLocalPsiNo]);
4183 for (i=0; i< n; i++) {
4184 currentPsi[i].re *= m/old_norm; // real part
4185 currentPsi[i].im *= m/old_norm; // imaginary part
4186 }
4187 }
4188 P->R.PsiStep = P->R.MaxPsiStep; // make it advance to next Psi
4189 //debug(P,"UpdateActualPsiNo");
4190 UpdateActualPsiNo(P, type); // step on to next perturbed Psi
4191 if (*SuperStop != 1)
4192 *SuperStop = CheckCPULIM(P);
4193 *Stop = CalculateMinimumStop(P, *SuperStop);
4194 P->Speed.Steps++; // step on
4195 R->LevS->Step++;
4196 }
4197 if(P->Call.out[NormalOut]) fprintf(stderr,"(%i) Write %s srcpsi to disk\n", P->Par.me, R->MinimisationName[type]);
4198 OutputSrcPsiDensity(P, type);
4199// if (!TestReadnWriteSrcDensity(P,type))
4200// Error(SomeError,"TestReadnWriteSrcDensity failed!");
4201 }
4202
4203 TestGramSch(P,R->LevS,Psi, type); // functions are orthonormal?
4204 // calculate current density summands
4205 //if (P->Call.out[StepLeaderOut]) fprintf(stderr,"(%i) Filling current density grid ...\n",P->Par.me);
4206 SpeedMeasure(P, CurrDensTime, StartTimeDo);
4207 if (*SuperStop != 1) {
4208 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
4209 R->DoFullCurrent = 1; // set to 1 if it was 2 but Check...() yielded necessity
4210 //debug(P,"Filling with Delta j ...");
4211 //FillDeltaCurrentDensity(P);
4212 }// else
4213 //debug(P,"There is no overlap between orbitals.");
4214 //debug(P,"Filling with j ...");
4215 FillCurrentDensity(P);
4216 }
4217 SpeedMeasure(P, CurrDensTime, StopTimeDo);
4218
4219 SetGramSchExtraPsi(P,Psi,NotUsedToOrtho); // remove extra Psis from orthogonality check
4220 ResetGramSchTagType(P, Psi, type, NotUsedToOrtho); // remove this group from the check for the next minimisation group as well!
4221 }
4222 UpdateActualPsiNo(P, Occupied); // step on back to an occupied one
4223
4224 gsl_multimin_fdfminimizer_free (s_multi);
4225 gsl_min_fminimizer_free (s);
4226 gsl_vector_free (x);
4227 //gsl_vector_free (ss);
4228}
4229*/
Note: See TracBrowser for help on using the repository browser.