1 | /** \file gramsch.c
|
---|
2 | * Gram-Schmidt-Orthonormalisation.
|
---|
3 | * Herein are all the functions necessary to orthogonalize and normalize the wave
|
---|
4 | * functions OnePsiElement, such as initialization FirstInitGramSchData(), norm
|
---|
5 | * GramSchNormalize(), scalar product GramSchSP() and the actual Gram-Schmidt-routine
|
---|
6 | * GramSch(). All depending on the current status of the wave function.
|
---|
7 | *
|
---|
8 | Project: ParallelCarParrinello
|
---|
9 | \author Jan Hamaekers
|
---|
10 | \date 2000
|
---|
11 |
|
---|
12 | File: gramsch.c
|
---|
13 | $Id: gramsch.c,v 1.70.2.1 2007-04-21 12:49:50 foo Exp $
|
---|
14 | */
|
---|
15 |
|
---|
16 | #include <stdlib.h>
|
---|
17 | #include <stdio.h>
|
---|
18 | #include <math.h>
|
---|
19 | #include <string.h>
|
---|
20 |
|
---|
21 | // use double precision fft when we have it
|
---|
22 | #ifdef HAVE_CONFIG_H
|
---|
23 | #include <config.h>
|
---|
24 | #endif
|
---|
25 |
|
---|
26 | #ifdef HAVE_DFFTW_H
|
---|
27 | #include "dfftw.h"
|
---|
28 | #else
|
---|
29 | #include "fftw.h"
|
---|
30 | #endif
|
---|
31 |
|
---|
32 | #include "data.h"
|
---|
33 | #include "helpers.h"
|
---|
34 | #include "errors.h"
|
---|
35 | #include "myfft.h"
|
---|
36 | #include "mymath.h"
|
---|
37 | #include "gramsch.h"
|
---|
38 | #include "run.h"
|
---|
39 | #include "mergesort2.h"
|
---|
40 |
|
---|
41 | /** Deallocates the defined OnePsiElement datatype.
|
---|
42 | */
|
---|
43 | void FreeMPI_OnePsiElement()
|
---|
44 | {
|
---|
45 | MPI_Type_free(&MPI_OnePsiElement);
|
---|
46 | }
|
---|
47 |
|
---|
48 | /** Initialization of Gram-Schmidt-Orthogonalization.
|
---|
49 | * \param *P Problem at hand
|
---|
50 | * \param *Psi wave functions
|
---|
51 | * \sa RemoveEverything()
|
---|
52 | */
|
---|
53 | void FirstInitGramSchData(struct Problem *P, struct Psis *Psi) {
|
---|
54 | int i, type;
|
---|
55 | int GramSchLocalNo = Psi->LocalNo+1;
|
---|
56 | MPI_Datatype type1[10] = { MPI_INT, MPI_INT, MPI_INT, MPI_INT, MPI_INT, MPI_INT, MPI_DOUBLE, MPI_DOUBLE, MPI_INT, MPI_UB}; // type of each OPE array element
|
---|
57 | int blocklen1[10] = { 1, 1, 1, 1, 1, 1, 1, 1, 1, 1}; // block length of each element within the OPE array
|
---|
58 | MPI_Aint base, disp1[10]; // holds adresses in memory
|
---|
59 | struct OnePsiElement OPE[2];
|
---|
60 |
|
---|
61 | /// Create MPI_OnePsiElement, simulacrum of OnePsiElement, enabling exchange of these among the processes
|
---|
62 | // store adresses of its various elements in disp1 array
|
---|
63 | MPI_Address( OPE, disp1);
|
---|
64 | MPI_Address( &OPE[0].my_color_comm_ST_Psi, disp1+1);
|
---|
65 | MPI_Address( &OPE[0].MyLocalNo, disp1+2);
|
---|
66 | MPI_Address( &OPE[0].MyGlobalNo, disp1+3);
|
---|
67 | MPI_Address( &OPE[0].PsiGramSchStatus, disp1+4);
|
---|
68 | MPI_Address( &OPE[0].PsiType, disp1+5);
|
---|
69 | MPI_Address( &OPE[0].PsiFactor, disp1+6);
|
---|
70 | MPI_Address( &OPE[0].PsiReciNorm2, disp1+7);
|
---|
71 | MPI_Address( &OPE[0].DoBrent, disp1+8);
|
---|
72 | MPI_Address( OPE+1, disp1+9);
|
---|
73 | base = disp1[0];
|
---|
74 | for (i=0; i < 10; i++) disp1[i] -= base; // make the adresses of OPE elements relativ to base -> byte displacement of each entry
|
---|
75 | MPI_Type_struct( 10, blocklen1, disp1, type1, &MPI_OnePsiElement); // creates MPI_OnePsiElement as an MPI_struct(ure)
|
---|
76 | MPI_Type_commit( &MPI_OnePsiElement); // commits new data type, now it's usable
|
---|
77 | if (P->Call.out[NormalOut]) fprintf(stderr, "(%i)FirstInitGramSchData\n", P->Par.me);
|
---|
78 |
|
---|
79 | /// Allocates and fills Psis::AllLocalNo (MPI_Allgathered from all other processes).
|
---|
80 | Psi->AllLocalNo = (int *)
|
---|
81 | Malloc(sizeof(int)*P->Par.Max_me_comm_ST_PsiT,"FirstInitGramSchData: Psi->AllLocalNo");
|
---|
82 | MPI_Allgather ( &GramSchLocalNo, 1, MPI_INT,
|
---|
83 | Psi->AllLocalNo, 1, MPI_INT, P->Par.comm_ST_PsiT );
|
---|
84 |
|
---|
85 | /// Calculates from this Psis::MaxPsiOfType.
|
---|
86 | Psi->MaxPsiOfType = 0;
|
---|
87 | for (i=0;i<P->Par.Max_me_comm_ST_PsiT;i++) Psi->MaxPsiOfType += Psi->AllLocalNo[i]-1; // sum up all local (orthogonalizable) Psis in the transposed communicator PsiT
|
---|
88 | fprintf(stderr,"(%i) MaxPsiOfType = %i\n",P->Par.me, Psi->MaxPsiOfType);
|
---|
89 |
|
---|
90 | /// Calculates from this Psis::MaxPsiOfType and at which index this process' Psis start Psis::MyStartNo.
|
---|
91 | Psi->MyStartNo = 0;
|
---|
92 | for (i=0;i<P->Par.me_comm_ST_PsiT;i++) Psi->MyStartNo += Psi->AllLocalNo[i]; // where do my Psis start
|
---|
93 | fprintf(stderr,"(%i) MyStartNo = %i\n",P->Par.me, Psi->MyStartNo);
|
---|
94 |
|
---|
95 | //fprintf(stderr,"(%i) OtherPsiLocalNo %d\n",P->Par.me, RecvCount);
|
---|
96 | /// Allocates arrays Psis::AllPsiStatus, Psis::AllPsiStatusForSort and Psis::LocalPsiStatus (up 'til Extra in PsiTagType)
|
---|
97 | Psi->AllPsiStatus = (struct OnePsiElement *)
|
---|
98 | Malloc(sizeof(struct OnePsiElement)*(Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT),"FirstInitGramSchData: Psi->AllPsiStatus");
|
---|
99 | Psi->AllPsiStatusForSort = (struct OnePsiElement *)
|
---|
100 | Malloc(sizeof(struct OnePsiElement)*(Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT+1),"FirstInitGramSchData: Psi->AllPsiStatusForSort");
|
---|
101 | Psi->LocalPsiStatus = (struct OnePsiElement *)
|
---|
102 | Malloc(sizeof(struct OnePsiElement)*GramSchLocalNo,"FirstInitGramSchData: Psi->LocalPsiStatus");
|
---|
103 | /// Psis::LocalPsiStatus is initialized and distributed among all processes as Psis::AllPsiStatus.
|
---|
104 | for (i=0;i<GramSchLocalNo;i++) {
|
---|
105 | Psi->LocalPsiStatus[i].me_comm_ST_Psi = P->Par.me_comm_ST_Psi;
|
---|
106 | Psi->LocalPsiStatus[i].my_color_comm_ST_Psi = P->Par.my_color_comm_ST_Psi;
|
---|
107 | Psi->LocalPsiStatus[i].MyLocalNo = i;
|
---|
108 | Psi->LocalPsiStatus[i].MyGlobalNo = Psi->MyStartNo + i;
|
---|
109 | Psi->LocalPsiStatus[i].DoBrent = 4;
|
---|
110 | switch (Psi->PsiST) { // set occupation number for the regular local, one extra(!) per process (without current one!) and the additional orbitals (the "latterest" ;) are set to zero of course) (NOTE: extra orbit must always be the very last one (that's why Par->.. - 1)
|
---|
111 | case SpinDouble:
|
---|
112 | for (type=Occupied;type<=Extra;type++) {
|
---|
113 | if (i >= Psi->TypeStartIndex[type] && i < Psi->TypeStartIndex[type+1]) {
|
---|
114 | Psi->LocalPsiStatus[i].PsiType = type;
|
---|
115 | Psi->LocalPsiStatus[i].PsiGramSchStatus = (int)(type != Occupied ? NotUsedToOrtho : NotOrthogonal); // extra or occupied wave function
|
---|
116 | }
|
---|
117 | }
|
---|
118 | if (Psi->LocalPsiStatus[i].PsiType != UnOccupied)
|
---|
119 | Psi->LocalPsiStatus[i].PsiFactor = 2.0;
|
---|
120 | else
|
---|
121 | Psi->LocalPsiStatus[i].PsiFactor = 1.0;
|
---|
122 | break;
|
---|
123 | case SpinUp:
|
---|
124 | for (type=Occupied;type<=Extra;type++) {
|
---|
125 | if (i >= Psi->TypeStartIndex[type] && i < Psi->TypeStartIndex[type+1]) {
|
---|
126 | Psi->LocalPsiStatus[i].PsiType = type;
|
---|
127 | Psi->LocalPsiStatus[i].PsiGramSchStatus = (int)(type != Occupied ? NotUsedToOrtho : NotOrthogonal); // extra or occupied wave function
|
---|
128 | }
|
---|
129 | }
|
---|
130 | Psi->LocalPsiStatus[i].PsiFactor = 1.0;
|
---|
131 | break;
|
---|
132 | case SpinDown:
|
---|
133 | for (type=Occupied;type<=Extra;type++) {
|
---|
134 | if (i >= Psi->TypeStartIndex[type] && i < Psi->TypeStartIndex[type+1]) {
|
---|
135 | Psi->LocalPsiStatus[i].PsiType = type;
|
---|
136 | Psi->LocalPsiStatus[i].PsiGramSchStatus = (int)(type != Occupied ? NotUsedToOrtho : NotOrthogonal); // extra or occupied wave function
|
---|
137 | }
|
---|
138 | }
|
---|
139 | Psi->LocalPsiStatus[i].PsiFactor = 1.0;
|
---|
140 | break;
|
---|
141 | }
|
---|
142 | Psi->LocalPsiStatus[i].PsiReciNorm2 = 0.0;
|
---|
143 | }
|
---|
144 |
|
---|
145 | // Update AllPsiStatus from changed LocalPsiStatus
|
---|
146 | UpdateGramSchAllPsiStatus(P,Psi);
|
---|
147 |
|
---|
148 | /// Psis::TempSendA, Psis::AllActualLocalPsiNo and Psis::AllOldActualLocalPsiNo are allocated, the latter two zeroed.
|
---|
149 | Psi->TempSendA = (int *)
|
---|
150 | Malloc(sizeof(int)*P->Par.Max_me_comm_ST_PsiT,"FirstInitGramSchData: Psi->TempSendA");
|
---|
151 | Psi->AllActualLocalPsiNo = (int *)
|
---|
152 | Malloc(sizeof(int)*P->Par.Max_me_comm_ST_PsiT,"FirstInitGramSchData: Psi->AllActualLocalPsiNo");
|
---|
153 | Psi->AllOldActualLocalPsiNo = (int *)
|
---|
154 | Malloc(sizeof(int)*P->Par.Max_me_comm_ST_PsiT,"FirstInitGramSchData: Psi->AllOldActualLocalPsiNo");
|
---|
155 | for (i=0; i < P->Par.Max_me_comm_ST_PsiT; i++) {
|
---|
156 | Psi->AllActualLocalPsiNo[i] = 0;
|
---|
157 | Psi->AllOldActualLocalPsiNo[i] = 0;
|
---|
158 | }
|
---|
159 | }
|
---|
160 |
|
---|
161 | /** Normalize the coefficients of a given wave function.
|
---|
162 | * Calculates the norm (see GramSchGetNorm2()) and divides each (for all reciprocal grid
|
---|
163 | * vectors) complex coefficient by the norm.
|
---|
164 | * \param *P Problem at hand
|
---|
165 | * \param *Lev LatticeLevel structure
|
---|
166 | * \param *LPsiDat Array of complex wave function coefficients
|
---|
167 | * \param PsiSP If norm already calculated, can be passed on here, otherweise (== 0.0) is calculated
|
---|
168 | * \return Squared norm of wave function
|
---|
169 | */
|
---|
170 | static double GramSchNormalize(const struct Problem *P, struct LatticeLevel *Lev, fftw_complex *LPsiDat, double PsiSP) {
|
---|
171 | double LocalSP=0.0;
|
---|
172 | int i,s = 0;
|
---|
173 | /* Falls PsiSP == 0.0 dann noch SP berechnen */
|
---|
174 | if (PsiSP == 0.0) { // see GramSchGetNorm2()
|
---|
175 | if (Lev->GArray[0].GSq == 0.0) {
|
---|
176 | LocalSP += LPsiDat[0].re*LPsiDat[0].re;
|
---|
177 | s++;
|
---|
178 | }
|
---|
179 | for (i=s; i < Lev->MaxG; i++) {
|
---|
180 | LocalSP += 2*(LPsiDat[i].re*LPsiDat[i].re+LPsiDat[i].im*LPsiDat[i].im);
|
---|
181 | }
|
---|
182 | MPI_Allreduce ( &LocalSP, &PsiSP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
|
---|
183 | }
|
---|
184 | if (PsiSP < MYEPSILON) fprintf(stderr,"GramSchNormalize: PsiSP = %lg\n",PsiSP);
|
---|
185 | PsiSP = sqrt(PsiSP); // take square root
|
---|
186 | for (i=0; i < Lev->MaxG; i++) { // and divide each coefficient by the norm
|
---|
187 | LPsiDat[i].re /= PsiSP;
|
---|
188 | LPsiDat[i].im /= PsiSP;
|
---|
189 | }
|
---|
190 | return(PsiSP);
|
---|
191 | }
|
---|
192 |
|
---|
193 | /** Calculate squared norm of given wave function coefficients.
|
---|
194 | * Go through each node of the reciprocal vector grid, calculate the complex product for this
|
---|
195 | * coefficient and sum up, gathering the results from all processes before return - remember
|
---|
196 | * that the coefficients are - for the parallel calculation of the fft - split up among the
|
---|
197 | * processes.
|
---|
198 | * \param *P Problem at hand
|
---|
199 | * \param *Lev LatticeLevel structure
|
---|
200 | * \param *LPsiDat array over G of complex i-th wave function coefficients \f$c_{i,G}\f$
|
---|
201 | * \return \f$\sum_G c_{i,G} /cdot {c_{i,G}}^{\ast}\f$
|
---|
202 | */
|
---|
203 | double GramSchGetNorm2(const struct Problem *P, struct LatticeLevel *Lev, fftw_complex *LPsiDat) {
|
---|
204 | double LocalSP=0.0, PsiSP=0.0;
|
---|
205 | int i,s = 0;
|
---|
206 | /* Falls PsiSP == 0.0 dann noch SP berechnen */
|
---|
207 | if (LPsiDat != NULL) {
|
---|
208 | if (Lev->GArray[0].GSq == 0.0) {
|
---|
209 | LocalSP += LPsiDat[0].re*LPsiDat[0].re;
|
---|
210 | s++;
|
---|
211 | }
|
---|
212 | for (i=s; i < Lev->MaxG; i++) {
|
---|
213 | LocalSP += 2*(LPsiDat[i].re*LPsiDat[i].re+LPsiDat[i].im*LPsiDat[i].im);
|
---|
214 | }
|
---|
215 | // send local result to all processes and received summed from all into PsiSP
|
---|
216 | MPI_Allreduce ( &LocalSP, &PsiSP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
|
---|
217 | }
|
---|
218 | return(PsiSP);
|
---|
219 | }
|
---|
220 |
|
---|
221 | /** Scalar Product of two arrays of wave function coefficients.
|
---|
222 | * Goes through each reciprocal grid vectors and calculates the complex product
|
---|
223 | * between the two coefficients, summing up, MPI_Allreducing and returning.
|
---|
224 | * (See also GramSchGetNorm2())
|
---|
225 | * \param *P Problem at hand
|
---|
226 | * \param *Lev LatticeLevel structure
|
---|
227 | * \param *LPsiDatA first array of wave function coefficients
|
---|
228 | * \param *LPsiDatB second array of wave function coefficients
|
---|
229 | * \return \f$\sum_G c_{a,G} \cdot c_{b,G}^{\ast}\f$
|
---|
230 | */
|
---|
231 | static double GramSchSP(const struct Problem *P, struct LatticeLevel *Lev, fftw_complex *LPsiDatA, fftw_complex *LPsiDatB) {
|
---|
232 | double LocalSP=0.0,PsiSP;
|
---|
233 | int i,s = 0;
|
---|
234 | if (Lev->GArray[0].GSq == 0.0) {
|
---|
235 | LocalSP += LPsiDatA[0].re*LPsiDatB[0].re;
|
---|
236 | s++;
|
---|
237 | }
|
---|
238 | for (i=s; i < Lev->MaxG; i++) { // go through all nodes and calculate complex scalar product
|
---|
239 | LocalSP += 2*(LPsiDatA[i].re*LPsiDatB[i].re+LPsiDatA[i].im*LPsiDatB[i].im);
|
---|
240 | }
|
---|
241 | // send local result to all processes and received summed from all into PsiSP
|
---|
242 | MPI_Allreduce ( &LocalSP, &PsiSP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi);
|
---|
243 | return(PsiSP);
|
---|
244 | }
|
---|
245 |
|
---|
246 | /** Sort criteria for natueralmergesort(): Returns re-ordered OnePsiElement::PsiGramSchStatus.
|
---|
247 | * The current status in the Gram-Schmidt-Orthonormalization is returned as sort criteria.
|
---|
248 | * \param *a OnePsiElement
|
---|
249 | * \param i i-th wave function
|
---|
250 | * \param *Args
|
---|
251 | * \return integer value for each PsiGramSchStatusType, from IsOrthonormal (0) up to NotOrthogonal(2) and NotUsedToOrtho(3)
|
---|
252 | * \note The enum PsiGramSchStatusType is not simply copied due to a different ordering in the enumeration other than used here.
|
---|
253 | */
|
---|
254 | static double GetKeyOnePsi(void *a, int i, void *Args) {
|
---|
255 | double res=-1;
|
---|
256 | switch ((enum PsiGramSchStatusType)((struct OnePsiElement *)a)[i].PsiGramSchStatus) {
|
---|
257 | case NotOrthogonal:
|
---|
258 | res = 2.;
|
---|
259 | break;
|
---|
260 | case IsOrthogonal:
|
---|
261 | res = 1.;
|
---|
262 | break;
|
---|
263 | case IsOrthonormal:
|
---|
264 | res = 0.;
|
---|
265 | break;
|
---|
266 | case NotUsedToOrtho:
|
---|
267 | res = 100.; // extra before unoccupied and perturbed ones
|
---|
268 | break;
|
---|
269 | }
|
---|
270 | switch (((struct OnePsiElement *)a)[i].PsiType) {
|
---|
271 | case Occupied:
|
---|
272 | res += 0.;
|
---|
273 | break;
|
---|
274 | case UnOccupied:
|
---|
275 | res += 10.;
|
---|
276 | break;
|
---|
277 | case Perturbed_P0:
|
---|
278 | case Perturbed_P1:
|
---|
279 | case Perturbed_P2:
|
---|
280 | case Perturbed_RxP0:
|
---|
281 | case Perturbed_RxP1:
|
---|
282 | case Perturbed_RxP2:
|
---|
283 | res += 20.;
|
---|
284 | break;
|
---|
285 | case Extra:
|
---|
286 | res += 30.;
|
---|
287 | break;
|
---|
288 | }
|
---|
289 | return(res);
|
---|
290 | }
|
---|
291 |
|
---|
292 | /** Sort criteria for natueralmergesort(): Returns the global number of the Psi among all.
|
---|
293 | * \param *a OnePsiElement
|
---|
294 | * \param i i-th wave function
|
---|
295 | * \param *Args unused, for contingency with GetKeyOnePsi()
|
---|
296 | * \return \a i-th OnePsiElement::MyGlobalNo
|
---|
297 | */
|
---|
298 | static double GetKeyOnePsi2(void *a, int i, void *Args) {
|
---|
299 | return(((struct OnePsiElement *)a)[i].MyGlobalNo);
|
---|
300 | }
|
---|
301 |
|
---|
302 | /** Copies wave function OnePsiElement.
|
---|
303 | * Copy each entry in OnePsiElement structure from \a b[j] to \a a[i].
|
---|
304 | * \param *a destination OnePsiElement array
|
---|
305 | * \param i i-th element to be overwritten
|
---|
306 | * \param *b source OnePsiElement array
|
---|
307 | * \param j j-th element's entries to be copied
|
---|
308 | */
|
---|
309 | static void CopyElementOnePsi(void *a, int i, void *b, int j)
|
---|
310 | {
|
---|
311 | ((struct OnePsiElement *)a)[i].me_comm_ST_Psi = ((struct OnePsiElement *)b)[j].me_comm_ST_Psi;
|
---|
312 | ((struct OnePsiElement *)a)[i].my_color_comm_ST_Psi = ((struct OnePsiElement *)b)[j].my_color_comm_ST_Psi;
|
---|
313 | ((struct OnePsiElement *)a)[i].MyLocalNo = ((struct OnePsiElement *)b)[j].MyLocalNo;
|
---|
314 | ((struct OnePsiElement *)a)[i].MyGlobalNo = ((struct OnePsiElement *)b)[j].MyGlobalNo;
|
---|
315 | ((struct OnePsiElement *)a)[i].PsiGramSchStatus = ((struct OnePsiElement *)b)[j].PsiGramSchStatus;
|
---|
316 | ((struct OnePsiElement *)a)[i].PsiType = ((struct OnePsiElement *)b)[j].PsiType;
|
---|
317 | ((struct OnePsiElement *)a)[i].PsiFactor = ((struct OnePsiElement *)b)[j].PsiFactor;
|
---|
318 | ((struct OnePsiElement *)a)[i].PsiReciNorm2 = ((struct OnePsiElement *)b)[j].PsiReciNorm2;
|
---|
319 | }
|
---|
320 |
|
---|
321 |
|
---|
322 | /** Performs Gram-Schmidt-Orthonormalization on all Psis.
|
---|
323 | * Herein the known Gram-Schmidt-Orthogonalization (with subsequent normalization) is implemented in a
|
---|
324 | * parallel way. The problem arises due to the fact that the complex wave function coefficients are not
|
---|
325 | * all accessible from one process, but are shared among them. Thus there are four different cases to
|
---|
326 | * deal with - where O is one orthogonal Psi and P the Psi currently to be orthogonalized:
|
---|
327 | * -# O and P are local\n
|
---|
328 | * The projection is simply calculated via scalar product and subtracted from P.
|
---|
329 | * -# O is local, P not\n
|
---|
330 | * P is received from the respective process and the projetion calculated, noting down this
|
---|
331 | * value for later sending it back to this respective process owning the P coefficients,
|
---|
332 | * who will substract them
|
---|
333 | * -# O is not local, however P is\n
|
---|
334 | * Send the coefficient to every process in need of them and in the end gather projections to
|
---|
335 | * be subtracted from our local P.
|
---|
336 | * -# O and P are not local\n
|
---|
337 | * Nothing to do.
|
---|
338 | *
|
---|
339 | * Afterwards, a division by the norm of the Psi may additionally be called in for. The current status of
|
---|
340 | * a Psi is always noted in OnePsiElement::PsiGramSchStatus.
|
---|
341 | * \param *P Problem at hand
|
---|
342 | * \param *Lev LatticeLevel structure
|
---|
343 | * \param *Psi wave functions structure Psis
|
---|
344 | * \param ToDo states what to do in this function: Orthogonalize or Orthonormalize
|
---|
345 | */
|
---|
346 | void GramSch(struct Problem *P, struct LatticeLevel *Lev, struct Psis *Psi, enum PsiGramSchToDoType ToDo)
|
---|
347 | {
|
---|
348 | int i, j, k, TempRecv, TempSend, RecvSource;
|
---|
349 | //int ResetNo=0;
|
---|
350 | double GlobalSP;
|
---|
351 | struct RunStruct *R = &P->R;
|
---|
352 | struct OnePsiElement *OnePsi = NULL, *LOnePsi = NULL, *ROnePsi = NULL, *RLOnePsi = NULL;
|
---|
353 | int ElementSize = (sizeof(fftw_complex) / sizeof(double));
|
---|
354 | fftw_complex *Temp, *Temp2;
|
---|
355 | int *TempSendA = Psi->TempSendA;
|
---|
356 | MPI_Status status;
|
---|
357 | // Mergesort the wavefunction by their current status from 0 to all plus all extra ones (one for each process)
|
---|
358 | naturalmergesort(Psi->AllPsiStatus,Psi->AllPsiStatusForSort,0,Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT-1,&GetKeyOnePsi,NULL,&CopyElementOnePsi);
|
---|
359 | //fprintf(stderr,"(%i) GramSch: ",P->Par.me);
|
---|
360 | for (i=0; i < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; i++) { // then go through each of the ToDo-order sorted Psis (Each Psi plus an extra one from each process)
|
---|
361 | OnePsi = &Psi->AllPsiStatus[i]; // Mark OnePsi wave function
|
---|
362 | if (OnePsi->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // stored in this process? => L means local
|
---|
363 | LOnePsi = &Psi->LocalPsiStatus[OnePsi->MyLocalNo];
|
---|
364 | else
|
---|
365 | LOnePsi = NULL;
|
---|
366 |
|
---|
367 | switch ((enum PsiGramSchStatusType)OnePsi->PsiGramSchStatus) { // depending on their ToDo-status do ...
|
---|
368 | case NotOrthogonal: // ORTHOGONALIZE!
|
---|
369 | //fprintf(stderr,"(%i) ", Psi->AllPsiStatus[i].MyGlobalNo);
|
---|
370 | //fprintf(stderr,"Orthogonalizing %i, status was: (L)%i\t(A)%i!\n", OnePsi->MyGlobalNo, Psi->LocalPsiStatus[OnePsi->MyLocalNo].PsiGramSchStatus, OnePsi->PsiGramSchStatus);
|
---|
371 | if (LOnePsi != NULL) { // if current Psi is local, copy (reciprocal) complex coefficients from LocalPsi to TempPsi
|
---|
372 | memcpy(Lev->LPsi->TempPsi, Lev->LPsi->LocalPsi[OnePsi->MyLocalNo], ElementSize*Lev->MaxG*sizeof(double));
|
---|
373 | }
|
---|
374 | Temp = Lev->LPsi->TempPsi2; // another complex coefficients array (reciprocal) ...
|
---|
375 | SetArrayToDouble0((double *)Temp, Lev->MaxG*2); // ... which is zeroed
|
---|
376 | TempRecv = 0; // count how often a needed local current Psi has been received (and thus if it has been already)
|
---|
377 | TempSend = 0; // count how often a local current Psi has been sent to other processes
|
---|
378 | for(k=0; k < P->Par.Max_me_comm_ST_PsiT; k++) TempSendA[k] = 0; // zero array counting how often a process has sent its local Psi to others
|
---|
379 |
|
---|
380 | for (j=i-1; j >= 0; j--) { // go through all wave functions from the one before the current downto first one
|
---|
381 | ROnePsi = &Psi->AllPsiStatus[j]; // get the Psi that should be orthogonal to it
|
---|
382 | if (ROnePsi->PsiType <= UnOccupied) { // only orthogonalize against non-perturbed wave functions
|
---|
383 | if (ROnePsi->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // stored in this process? => L means local
|
---|
384 | RLOnePsi = &Psi->LocalPsiStatus[ROnePsi->MyLocalNo];
|
---|
385 | else
|
---|
386 | RLOnePsi = NULL;
|
---|
387 | // if (((OnePsi->PsiType == Extra && (R->CurrentMin <= UnOccupied || ((LOnePsi != NULL && RLOnePsi != NULL) && ROnePsi->MyLocalNo == R->ActualLocalPsiNo)))
|
---|
388 | // || OnePsi->PsiType <= UnOccupied)) {
|
---|
389 | if ((ROnePsi->PsiType == Occupied)
|
---|
390 | || ((ROnePsi->PsiType == UnOccupied) && (OnePsi->PsiType == UnOccupied || OnePsi->PsiType == Extra))
|
---|
391 | || ((LOnePsi != NULL && RLOnePsi != NULL) && ROnePsi->MyLocalNo == R->ActualLocalPsiNo)) {
|
---|
392 | // occupied are orthogonal to occupied
|
---|
393 | // unoccupied are orthogonal to occupied and unoccupied
|
---|
394 | // perturbed are orthogonal to occupied, unoccupied and to their (process-) specific extra
|
---|
395 | // extra is orthogonal dependent on R->CurrentMin (to occupied, occupied&unoccupied, occupied&specific perturbed)
|
---|
396 | if (RLOnePsi != NULL && LOnePsi != NULL) { // if both are stored locally in this process
|
---|
397 | GlobalSP = GramSchSP(P,Lev,Lev->LPsi->LocalPsi[ROnePsi->MyLocalNo],Lev->LPsi->LocalPsi[OnePsi->MyLocalNo]); // scalar product of the two
|
---|
398 | GlobalSP *= RLOnePsi->PsiReciNorm2; // divide by norm
|
---|
399 | Temp = Lev->LPsi->TempPsi;
|
---|
400 | Temp2 = Lev->LPsi->LocalPsi[ROnePsi->MyLocalNo];
|
---|
401 | for(k=0; k < Lev->MaxG; k++) { // orthogonalize it (subtract the projected part, real and imaginary)
|
---|
402 | Temp[k].re -= GlobalSP*Temp2[k].re;
|
---|
403 | Temp[k].im -= GlobalSP*Temp2[k].im;
|
---|
404 | } // the orthogonalized wave function of LocalPsi resides now in Temp!
|
---|
405 | }
|
---|
406 | if (RLOnePsi != NULL && LOnePsi == NULL) { // if the current Psi is not local, the one to which it ought be orthogonal however is
|
---|
407 | /* Recv i and put it to jLocal in TempPsi */
|
---|
408 | if (TempRecv == 0) {
|
---|
409 | MPI_Recv( Lev->LPsi->TempPsi, Lev->MaxG*ElementSize, MPI_DOUBLE, OnePsi->my_color_comm_ST_Psi, GramSchTag1, P->Par.comm_ST_PsiT, &status );
|
---|
410 | TempRecv++;
|
---|
411 | }
|
---|
412 | GlobalSP = GramSchSP(P,Lev,Lev->LPsi->LocalPsi[ROnePsi->MyLocalNo],Lev->LPsi->TempPsi);
|
---|
413 | GlobalSP *= RLOnePsi->PsiReciNorm2;
|
---|
414 | Temp = Lev->LPsi->TempPsi2;
|
---|
415 | Temp2 = Lev->LPsi->LocalPsi[ROnePsi->MyLocalNo];
|
---|
416 | for(k=0; k < Lev->MaxG; k++) {
|
---|
417 | Temp[k].re -= GlobalSP*Temp2[k].re;
|
---|
418 | Temp[k].im -= GlobalSP*Temp2[k].im;
|
---|
419 | } // the negative orthogonal projection resides in Temp (not the local wave function part!)
|
---|
420 | }
|
---|
421 | if (RLOnePsi == NULL && LOnePsi != NULL) { // if the current Psi is local, the one to which it ought be orthogonal yet is not
|
---|
422 | /* Send i to jLocal in TempPsi */
|
---|
423 | if (TempSendA[ROnePsi->my_color_comm_ST_Psi] == 0) { // just send it out to everyone who needs it
|
---|
424 | MPI_Send( Lev->LPsi->LocalPsi[OnePsi->MyLocalNo], Lev->MaxG*ElementSize, MPI_DOUBLE, ROnePsi->my_color_comm_ST_Psi, GramSchTag1, P->Par.comm_ST_PsiT);
|
---|
425 | TempSend++;
|
---|
426 | TempSendA[ROnePsi->my_color_comm_ST_Psi]++; // note that coefficients were sent once more to this process
|
---|
427 | }
|
---|
428 | }
|
---|
429 | }
|
---|
430 | }
|
---|
431 | }
|
---|
432 | if (LOnePsi != NULL) { // holds the current local Psi (TempPsi) and receives results from all other which "ought be orthogonal to this one"
|
---|
433 | /* Hat was in TempPsi und ist local*/
|
---|
434 | for (j=0; j < TempSend; j++) { // each of the recipients before should send something back now
|
---|
435 | MPI_Probe( MPI_ANY_SOURCE, GramSchTag2, P->Par.comm_ST_PsiT, &status );
|
---|
436 | RecvSource = status.MPI_SOURCE;
|
---|
437 | MPI_Recv( Lev->LPsi->TempPsi2, Lev->MaxG*ElementSize, MPI_DOUBLE, RecvSource, GramSchTag2, P->Par.comm_ST_PsiT, &status );
|
---|
438 | Temp2 = Lev->LPsi->TempPsi2;
|
---|
439 | Temp = Lev->LPsi->TempPsi;
|
---|
440 | for(k=0; k < Lev->MaxG; k++) { // sum received projetion onto (temporary) local wave function
|
---|
441 | Temp[k].re += Temp2[k].re;
|
---|
442 | Temp[k].im += Temp2[k].im;
|
---|
443 | }
|
---|
444 | }
|
---|
445 | Temp2 = Lev->LPsi->TempPsi;
|
---|
446 | Temp = Lev->LPsi->LocalPsi[OnePsi->MyLocalNo];
|
---|
447 | memcpy(Temp, Temp2, sizeof(fftw_complex)*Lev->MaxG); // finally copy back onto original one
|
---|
448 | }
|
---|
449 | if (LOnePsi == NULL && TempRecv) { // has calculated a projection to another Psi (TempPsi2) and sends it to the respective (local) owner of the current one
|
---|
450 | /* Hat was in TempPsi2 und ist nicht local*/
|
---|
451 | MPI_Send( Lev->LPsi->TempPsi2, Lev->MaxG*ElementSize, MPI_DOUBLE, OnePsi->my_color_comm_ST_Psi, GramSchTag2, P->Par.comm_ST_PsiT);
|
---|
452 | }
|
---|
453 | /*if (LOnePsi != NULL) { // finally we set the status of our local (multi-projection subtracted) Psi
|
---|
454 | //fprintf(stderr,"Setting L-Status of %i to %i\n",LOnePsi->MyGlobalNo, IsOrthogonal);
|
---|
455 | LOnePsi->PsiGramSchStatus = (int)IsOrthogonal;
|
---|
456 | }
|
---|
457 | fprintf(stderr,"Setting A-Status of %i to %i\n",OnePsi->MyGlobalNo, IsOrthogonal);
|
---|
458 | OnePsi->PsiGramSchStatus = (int)IsOrthogonal; // and also set the status in all processes for this Psi*/
|
---|
459 | // note: There is no break here, normalization will be performed right away!
|
---|
460 | //fprintf(stderr,"-> ");
|
---|
461 | case IsOrthogonal: // NORMALIZE!
|
---|
462 | //fprintf(stderr,"%i ", Psi->AllPsiStatus[i].MyGlobalNo);
|
---|
463 | switch (ToDo) {
|
---|
464 | case Orthonormalize: // ... normalize and store 1 as norm
|
---|
465 | if (LOnePsi != NULL) {
|
---|
466 | //fprintf(stderr,"Setting L-Status of %i to %i\n",LOnePsi->MyLocalNo, IsOrthonormal);
|
---|
467 | LOnePsi->PsiGramSchStatus = (int)IsOrthonormal; // set status and ...
|
---|
468 | GramSchNormalize(P,Lev,Lev->LPsi->LocalPsi[OnePsi->MyLocalNo],0.0); // ... do it
|
---|
469 | /*LOnePsi->PsiReciNorm2 = GramSchGetNorm2(P,Lev,Lev->LPsi->LocalPsi[OnePsi->MyLocalNo]);
|
---|
470 | LOnePsi->PsiReciNorm2 = 1./LOnePsi->PsiReciNorm2;*/
|
---|
471 | LOnePsi->PsiReciNorm2 = 1.;
|
---|
472 | }
|
---|
473 | //fprintf(stderr,"Setting A-Status of %i to %i\n",OnePsi->MyGlobalNo, IsOrthonormal);
|
---|
474 | OnePsi->PsiGramSchStatus = (int)IsOrthonormal;
|
---|
475 | break;
|
---|
476 | case Orthogonalize: // ... calculate norm and store
|
---|
477 | if (LOnePsi != NULL) {
|
---|
478 | //fprintf(stderr,"Setting L-Status of %i to %i\n",LOnePsi->MyLocalNo, IsOrthogonal);
|
---|
479 | LOnePsi->PsiGramSchStatus = (int)IsOrthogonal;
|
---|
480 | LOnePsi->PsiReciNorm2 = GramSchGetNorm2(P,Lev,Lev->LPsi->LocalPsi[OnePsi->MyLocalNo]);
|
---|
481 | if (LOnePsi->PsiReciNorm2 < MYEPSILON) fprintf(stderr,"GramSch: LOnePsi->PsiReciNorm2 Nr. %d = %lg\n",LOnePsi->MyGlobalNo,LOnePsi->PsiReciNorm2);
|
---|
482 | //if (P->Call.out[LeaderOut]) fprintf(stderr,"GramSch: Setting norm of LOnePsi->PsiReciNorm2 Nr. %d = %lg\n",LOnePsi->MyGlobalNo,LOnePsi->PsiReciNorm2);
|
---|
483 | LOnePsi->PsiReciNorm2 = 1./LOnePsi->PsiReciNorm2;
|
---|
484 | }
|
---|
485 | //fprintf(stderr,"Setting A-Status of %i to %i\n",OnePsi->MyGlobalNo, IsOrthogonal);
|
---|
486 | OnePsi->PsiGramSchStatus = (int)IsOrthogonal;
|
---|
487 | break;
|
---|
488 | }
|
---|
489 | break;
|
---|
490 | case IsOrthonormal: // NOTHING TO DO ANY MORE!
|
---|
491 | //fprintf(stderr,"%i ", Psi->AllPsiStatus[i].MyGlobalNo);
|
---|
492 | break;
|
---|
493 | case NotUsedToOrtho:
|
---|
494 | //fprintf(stderr,"[%i] ", Psi->AllPsiStatus[i].MyGlobalNo);
|
---|
495 | break;
|
---|
496 | }
|
---|
497 | }
|
---|
498 | //fprintf(stderr,"\n");
|
---|
499 | /* Reset */
|
---|
500 | naturalmergesort(Psi->AllPsiStatus,Psi->AllPsiStatusForSort,0,Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT-1,&GetKeyOnePsi2,NULL,&CopyElementOnePsi);
|
---|
501 | /* fprintf(stderr,"Setting L-Status of %i to %i\n",Psi->LocalNo, NotUsedToOrtho);
|
---|
502 | Psi->LocalPsiStatus[Psi->LocalNo].PsiGramSchStatus = (int)NotUsedToOrtho;
|
---|
503 | for (i=0; i < P->Par.Max_me_comm_ST_PsiT; i++) {
|
---|
504 | ResetNo += Psi->AllLocalNo[i];
|
---|
505 | OnePsi = &Psi->AllPsiStatus[ResetNo-1];
|
---|
506 | fprintf(stderr,"Setting A-Status of %i to %i\n",OnePsi->MyGlobalNo, NotUsedToOrtho);
|
---|
507 | OnePsi->PsiGramSchStatus = (int)NotUsedToOrtho;
|
---|
508 | }*/
|
---|
509 | }
|
---|
510 |
|
---|
511 | /** Reset status of Gram-Schmidt-Orthogonalization for each and every Psi.
|
---|
512 | * Sets all locally accessible Psis::LocalPsiStatus to PsiGramSchStatusType::NotOrthogonal
|
---|
513 | * and the norm to zero, except the last (extra) and unoccupied ones which are NotUsedToOrtho, then
|
---|
514 | * do the same for all Psis::AllPsiStatus (again exception for extra and unoccupied ones).
|
---|
515 | * \param *P Problem at hand
|
---|
516 | * \param *Psi wave functions structure Psis
|
---|
517 | */
|
---|
518 | void ResetGramSch(const struct Problem *P, struct Psis *Psi)
|
---|
519 | {
|
---|
520 | int i,j, ResetNo=0;
|
---|
521 | struct OnePsiElement *OnePsi = NULL;
|
---|
522 | for (i=0; i < Psi->LocalNo; i++) { // go through all local Psis
|
---|
523 | Psi->LocalPsiStatus[i].PsiGramSchStatus = (Psi->LocalPsiStatus[i].PsiType == Occupied) ? (int)NotOrthogonal : (int)NotUsedToOrtho;
|
---|
524 | //fprintf(stderr,"Setting L-Status of %i to %i\n",i, Psi->LocalPsiStatus[i].PsiGramSchStatus);
|
---|
525 | Psi->LocalPsiStatus[i].PsiReciNorm2 = 0.0;
|
---|
526 | }
|
---|
527 | //fprintf(stderr,"Setting L-Status of %i to %i\n",Psi->LocalNo, NotUsedToOrtho);
|
---|
528 | Psi->LocalPsiStatus[Psi->LocalNo].PsiGramSchStatus = (int)NotUsedToOrtho; // extra wave function
|
---|
529 | Psi->LocalPsiStatus[Psi->LocalNo].PsiReciNorm2 = 0.0;
|
---|
530 | for (i=0; i < P->Par.Max_me_comm_ST_PsiT; i++) {
|
---|
531 | for (j=ResetNo; j < ResetNo+Psi->AllLocalNo[i]-1; j++) {
|
---|
532 | Psi->AllPsiStatus[j].PsiGramSchStatus = (Psi->AllPsiStatus[j].PsiType == Occupied) ? (int)NotOrthogonal : (int)NotUsedToOrtho;
|
---|
533 | //fprintf(stderr,"Setting A-Status of %i to %i\n",j, Psi->AllPsiStatus[j].PsiGramSchStatus);
|
---|
534 | Psi->AllPsiStatus[j].PsiReciNorm2 = 0.0;
|
---|
535 | }
|
---|
536 | ResetNo += Psi->AllLocalNo[i];
|
---|
537 | OnePsi = &Psi->AllPsiStatus[ResetNo-1];
|
---|
538 | //fprintf(stderr,"Setting A-Status of %i to %i\n",ResetNo-1, NotUsedToOrtho);
|
---|
539 | OnePsi->PsiGramSchStatus = (int)NotUsedToOrtho; // extra wave function
|
---|
540 | OnePsi->PsiReciNorm2 = 0.0;
|
---|
541 | }
|
---|
542 | }
|
---|
543 |
|
---|
544 | /** Reset status of Gram-Schmidt-Orthogonalization for each Psi of PsiTagType \a type.
|
---|
545 | * Sets all locally accessible Psis::LocalPsiStatus to PsiGramSchStatusType::NotOrthogonal
|
---|
546 | * and the norm to zero, except the last (extra) and unoccupied ones which are NotUsedToOrtho, then
|
---|
547 | * do the same for all Psis::AllPsiStatus (again exception for extra and unoccupied ones).
|
---|
548 | * \param *P Problem at hand
|
---|
549 | * \param *Psi wave functions structure Psis
|
---|
550 | * \param type PsiTagType of orbitals whose PsiGramSchStatus is to be reset
|
---|
551 | * \param ToDo - set PsiGramSchToDoType for the \a type states to this
|
---|
552 | * \sa ResetGramSch() - same procedure for occupied states
|
---|
553 | */
|
---|
554 | void ResetGramSchTagType(const struct Problem *P, struct Psis *Psi, enum PsiTypeTag type, enum PsiGramSchStatusType ToDo)
|
---|
555 | {
|
---|
556 | int i,j, ResetNo=0;
|
---|
557 | struct OnePsiElement *OnePsi = NULL;
|
---|
558 | for (i=0; i < Psi->LocalNo; i++) { // go through all local Psis
|
---|
559 | if (Psi->LocalPsiStatus[i].PsiType == type) {
|
---|
560 | //fprintf(stderr,"Setting L-Status of %i to %i\n",i, ToDo);
|
---|
561 | Psi->LocalPsiStatus[i].PsiGramSchStatus = ToDo;
|
---|
562 | switch(ToDo) {
|
---|
563 | case NotOrthogonal:
|
---|
564 | Psi->LocalPsiStatus[i].PsiReciNorm2 = 0.0;
|
---|
565 | break;
|
---|
566 | default:
|
---|
567 | break;
|
---|
568 | }
|
---|
569 | }
|
---|
570 | }
|
---|
571 | //fprintf(stderr,"Setting L-Status of %i to %i\n",Psi->LocalNo, NotUsedToOrtho);
|
---|
572 | Psi->LocalPsiStatus[Psi->LocalNo].PsiGramSchStatus = (int)NotUsedToOrtho; // extra wave function
|
---|
573 | Psi->LocalPsiStatus[Psi->LocalNo].PsiReciNorm2 = 0.0;
|
---|
574 | for (i=0; i < P->Par.Max_me_comm_ST_PsiT; i++) {
|
---|
575 | for (j=ResetNo; j < ResetNo+Psi->AllLocalNo[i]-1; j++) {
|
---|
576 | if (Psi->AllPsiStatus[j].PsiType == type) {
|
---|
577 | //fprintf(stderr,"Setting A-Status of %i to %i\n",j, ToDo);
|
---|
578 | Psi->AllPsiStatus[j].PsiGramSchStatus = ToDo;
|
---|
579 | switch(ToDo) {
|
---|
580 | case NotOrthogonal:
|
---|
581 | Psi->AllPsiStatus[j].PsiReciNorm2 = 0.0;
|
---|
582 | break;
|
---|
583 | default:
|
---|
584 | break;
|
---|
585 | }
|
---|
586 | }
|
---|
587 | }
|
---|
588 | ResetNo += Psi->AllLocalNo[i];
|
---|
589 | OnePsi = &Psi->AllPsiStatus[ResetNo-1];
|
---|
590 | //fprintf(stderr,"Setting A-Status of %i to %i\n",ResetNo-1, NotUsedToOrtho);
|
---|
591 | OnePsi->PsiGramSchStatus = (int)NotUsedToOrtho; // extra wave function
|
---|
592 | // OnePsi->PsiReciNorm2 = 0.0;
|
---|
593 | }
|
---|
594 | }
|
---|
595 | /** Set Gram-Schmidt status of the extra Psis::LocalPsiStatus and Psis::AllPsiStatus Psis to \a PsiGramSchStatus.
|
---|
596 | * The number of these "extra" Psis is Psis::AllLocalNo - 1 for each process.
|
---|
597 | * \param *P Problem at hand
|
---|
598 | * \param *Psi wave functions structure Psis
|
---|
599 | * \param PsiGramSchStatus status to be set
|
---|
600 | */
|
---|
601 | void SetGramSchExtraPsi(const struct Problem *P, struct Psis *Psi, enum PsiGramSchStatusType PsiGramSchStatus)
|
---|
602 | {
|
---|
603 | int i, ResetNo=0; // offset to the respective local Psis
|
---|
604 | struct OnePsiElement *OnePsi = NULL;
|
---|
605 | //fprintf(stderr,"Setting L-Status of %i to %i\n",Psi->LocalNo, PsiGramSchStatus);
|
---|
606 | Psi->LocalPsiStatus[Psi->LocalNo].PsiGramSchStatus = (int)PsiGramSchStatus;
|
---|
607 | for (i=0; i < P->Par.Max_me_comm_ST_PsiT; i++) {
|
---|
608 | ResetNo += Psi->AllLocalNo[i];
|
---|
609 | OnePsi = &Psi->AllPsiStatus[ResetNo-1];
|
---|
610 | //fprintf(stderr,"Setting A-Status of %i to %i\n",ResetNo-1, PsiGramSchStatus);
|
---|
611 | OnePsi->PsiGramSchStatus = (int)PsiGramSchStatus;
|
---|
612 | }
|
---|
613 | }
|
---|
614 |
|
---|
615 | /** Set Gram-Schmidt status of the actual Psis::LocalPsiStatus and Psis::AllPsiStatus Psis to \a PsiGramSchStatus.
|
---|
616 | * The number of these "extra" Psis is Psis::AllActualLocalPsiNo for each process.
|
---|
617 | * \param *P Problem at hand
|
---|
618 | * \param *Psi wave functions structure Psis
|
---|
619 | * \param PsiGramSchStatus status to be set
|
---|
620 | */
|
---|
621 | void SetGramSchActualPsi(const struct Problem *P, struct Psis *Psi, enum PsiGramSchStatusType PsiGramSchStatus)
|
---|
622 | {
|
---|
623 | int i, ResetNo=0; // offset to the respective local Psis
|
---|
624 | struct OnePsiElement *OnePsi = NULL;
|
---|
625 | //fprintf(stderr,"Setting L-Status of %i to %i\n",P->R.ActualLocalPsiNo, PsiGramSchStatus);
|
---|
626 | //BUG: Psi->LocalPsiStatus[Psi->LocalNo].PsiGramSchStatus = (int)PsiGramSchStatus;
|
---|
627 | Psi->LocalPsiStatus[P->R.ActualLocalPsiNo].PsiGramSchStatus = (int)PsiGramSchStatus;
|
---|
628 | for (i=0; i < P->Par.Max_me_comm_ST_PsiT; i++) {
|
---|
629 | OnePsi = &Psi->AllPsiStatus[ResetNo+Psi->AllActualLocalPsiNo[i]];
|
---|
630 | //fprintf(stderr,"Setting A-Status of %i to %i\n",ResetNo+Psi->AllActualLocalPsiNo[i], PsiGramSchStatus);
|
---|
631 | OnePsi->PsiGramSchStatus = (int)PsiGramSchStatus;
|
---|
632 | ResetNo += Psi->AllLocalNo[i];
|
---|
633 | }
|
---|
634 | }
|
---|
635 |
|
---|
636 | /** Set Gram-Schmidt status of the former actual Psis::LocalPsiStatus and Psis::AllPsiStatus Psis to \a PsiGramSchStatus.
|
---|
637 | * The former number of these "extra" Psis is Psis::AllOldActualLocalPsiNo for each process.
|
---|
638 | * \param *P Problem at hand
|
---|
639 | * \param *Psi wave functions structure Psis
|
---|
640 | * \param PsiGramSchStatus status to be set
|
---|
641 | */
|
---|
642 | void SetGramSchOldActualPsi(const struct Problem *P, struct Psis *Psi, enum PsiGramSchStatusType PsiGramSchStatus)
|
---|
643 | {
|
---|
644 | int i, ResetNo=0;
|
---|
645 | struct OnePsiElement *OnePsi = NULL;
|
---|
646 | //fprintf(stderr,"Setting L-Status of %i to %i\n",P->R.OldActualLocalPsiNo, PsiGramSchStatus);
|
---|
647 | Psi->LocalPsiStatus[P->R.OldActualLocalPsiNo].PsiGramSchStatus = (int)PsiGramSchStatus;
|
---|
648 | for (i=0; i < P->Par.Max_me_comm_ST_PsiT; i++) {
|
---|
649 | OnePsi = &Psi->AllPsiStatus[ResetNo+Psi->AllOldActualLocalPsiNo[i]];
|
---|
650 | //fprintf(stderr,"Setting A-Status of %i to %i\n",ResetNo+Psi->AllOldActualLocalPsiNo[i], PsiGramSchStatus);
|
---|
651 | OnePsi->PsiGramSchStatus = (int)PsiGramSchStatus;
|
---|
652 | ResetNo += Psi->AllLocalNo[i];
|
---|
653 | }
|
---|
654 | }
|
---|
655 |
|
---|
656 | /** Updates array Psis::AllActualLocalPsiNo from the RunStruct::ActualLocalPsiNo by MPI_Allgather.
|
---|
657 | * \param *P Problem at hand
|
---|
658 | * \param *Psi wave functions structure Psis
|
---|
659 | */
|
---|
660 | void UpdateGramSchActualPsiNo(struct Problem *P, struct Psis *Psi)
|
---|
661 | {
|
---|
662 | struct RunStruct *R = &P->R;
|
---|
663 | MPI_Allgather ( &(R->ActualLocalPsiNo), 1, MPI_INT,
|
---|
664 | Psi->AllActualLocalPsiNo, 1, MPI_INT, P->Par.comm_ST_PsiT );
|
---|
665 | }
|
---|
666 |
|
---|
667 | /** Updates array Psis::AllPsiStatus from the Psis::LocalPsiStatus by MPI_Allgather.
|
---|
668 | * First, calculates number of MPI_OnePsiElement to be received and their displacements in the
|
---|
669 | * global Array Psis::AllPsiStatus, then follows MPI_Allgather and afterwards a printout to screen
|
---|
670 | * if verbose is specified.
|
---|
671 | * \param *P Problem at hand
|
---|
672 | * \param *Psi wave functions structure Psis
|
---|
673 | * \warning Don't use before FirstInitGramSch() due to needed declaration of MPI_OnePsiElement
|
---|
674 | */
|
---|
675 | void UpdateGramSchAllPsiStatus(struct Problem *P, struct Psis *Psi)
|
---|
676 | {
|
---|
677 | int *recvcounts, *displs;
|
---|
678 | int GramSchLocalNo = Psi->LocalNo+1;
|
---|
679 | int MyStartNo = 0;
|
---|
680 | int i;
|
---|
681 |
|
---|
682 | //recvcounts = (int *)Malloc(sizeof(int)*P->Par.Max_me_comm_ST_PsiT,"UpdateGramSchAllPsiStatus: recvcounts");
|
---|
683 | //displs = (int *)Malloc(sizeof(int)*P->Par.Max_me_comm_ST_PsiT,"UpdateGramSchAllPsiStatus: displs");
|
---|
684 | recvcounts = (int *)Malloc(sizeof(int)*P->Par.procs,"UpdateGramSchAllPsiStatus: recvcounts");
|
---|
685 | displs = (int *)Malloc(sizeof(int)*P->Par.procs,"UpdateGramSchAllPsiStatus: displs");
|
---|
686 |
|
---|
687 | for (i=0;i<P->Par.Max_me_comm_ST_PsiT;i++) {
|
---|
688 | recvcounts[i] = Psi->AllLocalNo[i]; // how many Psis should be received
|
---|
689 | displs[i] = MyStartNo; // displacement for these Psis
|
---|
690 | MyStartNo += Psi->AllLocalNo[i]; //
|
---|
691 | }
|
---|
692 | // send all (GramSchLocalNo) own local Psis and gather the AllPsiStatuses of all other processes
|
---|
693 | MPI_Allgatherv ( Psi->LocalPsiStatus, GramSchLocalNo, MPI_OnePsiElement,
|
---|
694 | Psi->AllPsiStatus, recvcounts, displs, MPI_OnePsiElement, P->Par.comm_ST_PsiT );
|
---|
695 |
|
---|
696 | if(P->Call.out[PsiOut])
|
---|
697 | for (i=0;i< MyStartNo;i++)
|
---|
698 | fprintf(stderr,"(%i) MyLocalNo = %i, MyGlobalNo = %i/%i, f = %lg, Type: %i, GramSch: %i, me_comm: %d, my_color_comm: %d \n",P->Par.me, Psi->AllPsiStatus[i].MyLocalNo, i, Psi->AllPsiStatus[i].MyGlobalNo, Psi->AllPsiStatus[i].PsiFactor, Psi->AllPsiStatus[i].PsiType, Psi->AllPsiStatus[i].PsiGramSchStatus, Psi->AllPsiStatus[i].me_comm_ST_Psi, Psi->AllPsiStatus[i].my_color_comm_ST_Psi);
|
---|
699 |
|
---|
700 | Free(recvcounts);
|
---|
701 | Free(displs);
|
---|
702 | }
|
---|
703 |
|
---|
704 | /** Updates array Psis::AllOldActualLocalPsiNo from the RunStruct::OldActualLocalPsiNo by MPI_Allgather.
|
---|
705 | * \param *P Problem at hand
|
---|
706 | * \param *Psi wave functions structure Psis
|
---|
707 | */
|
---|
708 | void UpdateGramSchOldActualPsiNo(struct Problem *P, struct Psis *Psi)
|
---|
709 | {
|
---|
710 | struct RunStruct *R = &P->R;
|
---|
711 | MPI_Allgather ( &(R->OldActualLocalPsiNo), 1, MPI_INT,
|
---|
712 | Psi->AllOldActualLocalPsiNo, 1, MPI_INT, P->Par.comm_ST_PsiT );
|
---|
713 | }
|
---|
714 |
|
---|
715 | #define max_GramSch_iter 3
|
---|
716 |
|
---|
717 | /** Test Gram-Schmidt-Orthogonalization.
|
---|
718 | * Test if all pairs of Psis are orthogonal respectively normalized (scalar product <= 1).
|
---|
719 | * Give output to stderr if not so.
|
---|
720 | * \param *P Problem at hand
|
---|
721 | * \param *Lev LatticeLevel structure
|
---|
722 | * \param *Psi wave functions structure Psis
|
---|
723 | * \param Type2test basically current minimisation type, see RunStruct#CurrentMin
|
---|
724 | */
|
---|
725 | void TestGramSch(struct Problem *P, struct LatticeLevel *Lev, struct Psis *Psi, int Type2test) {
|
---|
726 | double LocalSP=0.0,PsiSP;
|
---|
727 | int i,j,k,s,RecvSource;
|
---|
728 | MPI_Status status;
|
---|
729 | struct OnePsiElement *OnePsiA, *LOnePsiA, *LOnePsiB;
|
---|
730 | int ElementSize = (sizeof(fftw_complex) / sizeof(double));
|
---|
731 | int NotOk; // counts pairs that are not orthogonal
|
---|
732 | int iter = 0;
|
---|
733 | fftw_complex *LPsiDatA, *LPsiDatB;
|
---|
734 |
|
---|
735 | do {
|
---|
736 | NotOk = 0;
|
---|
737 | //fprintf(stderr,"(%i) Testing Orthogonality ... \n", P->Par.me);
|
---|
738 | for (i=0; i < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; i++) { // go through all wave functions (plus the extra one for each process)
|
---|
739 | OnePsiA = &Psi->AllPsiStatus[i]; // grab OnePsiA
|
---|
740 | //fprintf(stderr,"(%i) OnePsiA: Type %d, GlobalNo %d \n", P->Par.me, OnePsiA->PsiType, OnePsiA->MyGlobalNo);
|
---|
741 | if (OnePsiA->PsiGramSchStatus == (int)IsOrthonormal ||
|
---|
742 | OnePsiA->PsiGramSchStatus == (int)IsOrthogonal) { // if it has been orthogonalized
|
---|
743 | //fprintf(stderr,"(%i) ... orthogonal\n", P->Par.me);
|
---|
744 | if (OnePsiA->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // local?
|
---|
745 | LOnePsiA = &Psi->LocalPsiStatus[OnePsiA->MyLocalNo];
|
---|
746 | else
|
---|
747 | LOnePsiA = NULL;
|
---|
748 |
|
---|
749 | if (LOnePsiA == NULL) { // if it's not local ... receive it from respective process into TempPsi
|
---|
750 | RecvSource = OnePsiA->my_color_comm_ST_Psi;
|
---|
751 | MPI_Recv( Lev->LPsi->TempPsi, Lev->MaxG*ElementSize, MPI_DOUBLE, RecvSource, GramSchTag3, P->Par.comm_ST_PsiT, &status );
|
---|
752 | LPsiDatA=Lev->LPsi->TempPsi;
|
---|
753 | } else { // .. otherwise send it to all other processes (Max_me... - 1)
|
---|
754 | for (k=0;k<P->Par.Max_me_comm_ST_PsiT;k++)
|
---|
755 | if (k != OnePsiA->my_color_comm_ST_Psi)
|
---|
756 | MPI_Send( Lev->LPsi->LocalPsi[OnePsiA->MyLocalNo], Lev->MaxG*ElementSize, MPI_DOUBLE, k, GramSchTag3, P->Par.comm_ST_PsiT);
|
---|
757 | LPsiDatA=Lev->LPsi->LocalPsi[OnePsiA->MyLocalNo];
|
---|
758 | } // LPsiDatA is now set to the coefficients of OnePsi either stored or MPI_Received
|
---|
759 |
|
---|
760 | for (j=Psi->TypeStartIndex[Occupied]; j < Psi->TypeStartIndex[Extra]+1; j++) { // for all locally accessible including extra Psis
|
---|
761 | LOnePsiB = &Psi->LocalPsiStatus[j];
|
---|
762 | //if (LOnePsiB->PsiType > UnOccupied || OnePsiA->PsiType > UnOccupied) fprintf(stderr,"(%i) Checking global %i against local %i/%i\n",P->Par.me, OnePsiA->MyGlobalNo, LOnePsiB->MyLocalNo, LOnePsiB->MyGlobalNo);
|
---|
763 | if (LOnePsiB->PsiGramSchStatus == (int)IsOrthonormal ||
|
---|
764 | LOnePsiB->PsiGramSchStatus == (int)IsOrthogonal) { // if it's orthogonalized
|
---|
765 | LPsiDatB=Lev->LPsi->LocalPsi[LOnePsiB->MyLocalNo]; // set LPsiDatB onto it
|
---|
766 | s=0;
|
---|
767 | LocalSP = 0.0;
|
---|
768 | if (Lev->GArray[0].GSq == 0.0) { // calculate scalar product of LPsiDatA and LPsiDatB
|
---|
769 | LocalSP += LPsiDatA[0].re*LPsiDatB[0].re;
|
---|
770 | s++;
|
---|
771 | }
|
---|
772 | for (k=s; k < Lev->MaxG; k++) {
|
---|
773 | LocalSP += 2*(LPsiDatA[k].re*LPsiDatB[k].re+LPsiDatA[k].im*LPsiDatB[k].im);
|
---|
774 | }
|
---|
775 | MPI_Allreduce ( &LocalSP, &PsiSP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // gather by summation results from the group sharing the coefficients
|
---|
776 | //if (P->Call.out[LeaderOut])
|
---|
777 | switch (Type2test) {
|
---|
778 | default:
|
---|
779 | case -1: // test all, checked!
|
---|
780 | if (((LOnePsiB->PsiType <= UnOccupied || (LOnePsiB->MyLocalNo == P->R.ActualLocalPsiNo && OnePsiA->PsiType == Extra)) || (LOnePsiB->MyGlobalNo == OnePsiA->MyGlobalNo))) { // check if it's zero (orthogonal) and give output if wanted
|
---|
781 | if (i == LOnePsiB->MyGlobalNo && LOnePsiB->PsiGramSchStatus == (int)IsOrthonormal) {
|
---|
782 | if (fabs(PsiSP -1.0) >= MYEPSILON) {
|
---|
783 | fprintf(stderr,"(%i)(%i,%i) = %g ?= 1.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,OnePsiA->MyGlobalNo,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP -1.0), MYEPSILON);
|
---|
784 | NotOk++;
|
---|
785 | }
|
---|
786 | } else {
|
---|
787 | if (fabs(PsiSP) >= MYEPSILON && (LOnePsiB != OnePsiA && LOnePsiB->PsiType > UnOccupied)) {
|
---|
788 | fprintf(stderr,"(%i)(%i,%i) = %g ?= 0.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,OnePsiA->MyGlobalNo,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP), MYEPSILON);
|
---|
789 | NotOk++;
|
---|
790 | }
|
---|
791 | }
|
---|
792 | }
|
---|
793 | break;
|
---|
794 | case Occupied: // test unperturbed orthogonality, checked!
|
---|
795 | case UnOccupied:
|
---|
796 | if ((LOnePsiB->PsiType <= UnOccupied) && (OnePsiA->PsiType <= UnOccupied || OnePsiA->PsiType == Extra)) {
|
---|
797 | if (i == LOnePsiB->MyGlobalNo && LOnePsiB->PsiGramSchStatus == (int)IsOrthonormal) {
|
---|
798 | if (fabs(PsiSP -1.0) >= MYEPSILON) {
|
---|
799 | fprintf(stderr,"(%i)(%i,%i) = %g != 1.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,i,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP -1.0), MYEPSILON);
|
---|
800 | NotOk++;
|
---|
801 | } else {
|
---|
802 | //fprintf(stderr,"(%i)(%i,%i) = %g != 1.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,i,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP -1.0), MYEPSILON);
|
---|
803 | }
|
---|
804 | } else {
|
---|
805 | if (fabs(PsiSP) >= MYEPSILON) {
|
---|
806 | fprintf(stderr,"(%i)(%i,%i) = %g != 0.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,i,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP), MYEPSILON);
|
---|
807 | NotOk++;
|
---|
808 | } else {
|
---|
809 | //fprintf(stderr,"(%i)(%i,%i) = %g != 0.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,i,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP), MYEPSILON);
|
---|
810 | }
|
---|
811 | }
|
---|
812 | } else {
|
---|
813 | //fprintf(stderr,"(%i)(%i,%i) Not (Un)Occupied\n",P->Par.my_color_comm_ST,i,LOnePsiB->MyGlobalNo);
|
---|
814 | }
|
---|
815 | break;
|
---|
816 | case Perturbed_P0: // test perturbed orthogonality and normalization of all, checked!
|
---|
817 | case Perturbed_P1:
|
---|
818 | case Perturbed_P2:
|
---|
819 | case Perturbed_RxP0:
|
---|
820 | case Perturbed_RxP1:
|
---|
821 | case Perturbed_RxP2:
|
---|
822 | if ((((LOnePsiB->PsiType <= UnOccupied || LOnePsiB->PsiType == Type2test) && (OnePsiA->PsiType <= UnOccupied || OnePsiA->PsiType == Type2test) && (OnePsiA->PsiType != LOnePsiB->PsiType)) || (LOnePsiB->MyGlobalNo == OnePsiA->MyGlobalNo))) { // check if it's zero (orthogonal) and give output if wanted
|
---|
823 | if (i == LOnePsiB->MyGlobalNo && LOnePsiB->PsiGramSchStatus == (int)IsOrthonormal) {
|
---|
824 | if (fabs(PsiSP -1.0) >= MYEPSILON) {
|
---|
825 | fprintf(stderr,"(%i)(%i,%i) ... is not orthonormal: %g != 1.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,OnePsiA->MyGlobalNo,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP -1.0), MYEPSILON);
|
---|
826 | NotOk++;
|
---|
827 | }
|
---|
828 | } else {
|
---|
829 | if (fabs(PsiSP) >= MYEPSILON && (LOnePsiB->MyGlobalNo != OnePsiA->MyGlobalNo)) {
|
---|
830 | fprintf(stderr,"(%i)(%i,%i) ... is not orthogonal: %g != 0.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,OnePsiA->MyGlobalNo,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP), MYEPSILON);
|
---|
831 | NotOk++;
|
---|
832 | }
|
---|
833 | }
|
---|
834 | }
|
---|
835 | break;
|
---|
836 | case Extra:
|
---|
837 | if (((LOnePsiB->PsiType == Extra) || (LOnePsiB->PsiType == Occupied)) && ((OnePsiA->PsiType == Extra) || (OnePsiA->PsiType == Occupied))) {
|
---|
838 | if (i == LOnePsiB->MyGlobalNo && LOnePsiB->PsiGramSchStatus == (int)IsOrthonormal) {
|
---|
839 | if (fabs(PsiSP -1.0) >= MYEPSILON) {
|
---|
840 | fprintf(stderr,"(%i)(%i,%i) = %g != 1.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,OnePsiA->MyGlobalNo,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP -1.0), MYEPSILON);
|
---|
841 | NotOk++;
|
---|
842 | }
|
---|
843 | } else {
|
---|
844 | if (fabs(PsiSP) >= MYEPSILON) {
|
---|
845 | fprintf(stderr,"(%i)(%i,%i) = %g != 0.0 eps(%g >= %g)\n",P->Par.my_color_comm_ST,OnePsiA->MyGlobalNo,LOnePsiB->MyGlobalNo,PsiSP, fabs(PsiSP), MYEPSILON);
|
---|
846 | NotOk++;
|
---|
847 | }
|
---|
848 | }
|
---|
849 | }
|
---|
850 | break;
|
---|
851 | }
|
---|
852 | }
|
---|
853 | }
|
---|
854 | }
|
---|
855 | }
|
---|
856 | /* if (NotOk != 0) {
|
---|
857 | fprintf(stderr,"(%i) NotOk %i ... re-orthogonalizing type %i for the %ith time\n",P->Par.me, NotOk, Type2test, ++iter);
|
---|
858 | ResetGramSchTagType(P, Psi, Type2test, NotOrthogonal);
|
---|
859 | GramSch(P,Lev,Psi,Orthonormalize);
|
---|
860 | }*/
|
---|
861 | } while ((NotOk != 0) && (iter < max_GramSch_iter));
|
---|
862 | if (P->Call.out[StepLeaderOut]) { // final check if there have been any un-orthogonal pairs
|
---|
863 | if (Type2test != -1) {
|
---|
864 | if (NotOk == 0) {
|
---|
865 | fprintf(stderr,"(%i)TestGramSchm on %s: Ok !\n",P->Par.my_color_comm_ST, P->R.MinimisationName[Type2test]);
|
---|
866 | } else {
|
---|
867 | fprintf(stderr,"(%i)TestGramSchm on %s: There are %i pairs not Ok!\n",P->Par.my_color_comm_ST, P->R.MinimisationName[Type2test], NotOk);
|
---|
868 | //Error(SomeError,"Wave functions not orthonormal as they should be!");
|
---|
869 | }
|
---|
870 | } else {
|
---|
871 | if (NotOk == 0) {
|
---|
872 | fprintf(stderr,"(%i)TestGramSchm on all: Ok !\n",P->Par.my_color_comm_ST);
|
---|
873 | } else {
|
---|
874 | fprintf(stderr,"(%i)TestGramSchm on all: There are %i pairs not Ok!\n",P->Par.my_color_comm_ST,NotOk);
|
---|
875 | //Error(SomeError,"Wave functions not orthonormal as they should be!");
|
---|
876 | }
|
---|
877 | }
|
---|
878 | }
|
---|
879 | }
|
---|
880 |
|
---|
881 |
|
---|
882 | /** Test if a given wave function to all others.
|
---|
883 | * \param *P Problem at hand
|
---|
884 | * \param *Lev LatticeLevel structure
|
---|
885 | * \param *psi pointer to array with wave function coefficients
|
---|
886 | */
|
---|
887 | void TestForOrth(struct Problem *P, struct LatticeLevel *Lev, fftw_complex *psi) {
|
---|
888 | struct Lattice *Lat = &P->Lat;
|
---|
889 | struct Psis *Psi = &Lat->Psi;
|
---|
890 | double LocalSP=0.0,PsiSP;
|
---|
891 | int i,k,RecvSource;
|
---|
892 | MPI_Status status;
|
---|
893 | struct OnePsiElement *OnePsiA, *LOnePsiA;
|
---|
894 | int ElementSize = (sizeof(fftw_complex) / sizeof(double));
|
---|
895 | int NotOk = 0; // counts pairs that are not orthogonal
|
---|
896 | fftw_complex *LPsiDatA;
|
---|
897 | for (i=0; i < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; i++) { // go through all wave functions (plus the extra one for each process)
|
---|
898 | OnePsiA = &Psi->AllPsiStatus[i]; // grab OnePsiA
|
---|
899 | if (OnePsiA->PsiGramSchStatus == (int)IsOrthonormal ||
|
---|
900 | OnePsiA->PsiGramSchStatus == (int)IsOrthogonal) { // if it has been orthogonalized
|
---|
901 | if (OnePsiA->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // local?
|
---|
902 | LOnePsiA = &Psi->LocalPsiStatus[OnePsiA->MyLocalNo];
|
---|
903 | else
|
---|
904 | LOnePsiA = NULL;
|
---|
905 | if (LOnePsiA == NULL) { // if it's not local ... receive it from respective process into TempPsi
|
---|
906 | RecvSource = OnePsiA->my_color_comm_ST_Psi;
|
---|
907 | MPI_Recv( Lev->LPsi->TempPsi, Lev->MaxG*ElementSize, MPI_DOUBLE, RecvSource, GramSchTag3, P->Par.comm_ST_PsiT, &status );
|
---|
908 | LPsiDatA=Lev->LPsi->TempPsi;
|
---|
909 | } else { // .. otherwise send it to all other processes (Max_me... - 1)
|
---|
910 | for (k=0;k<P->Par.Max_me_comm_ST_PsiT;k++)
|
---|
911 | if (k != OnePsiA->my_color_comm_ST_Psi)
|
---|
912 | MPI_Send( Lev->LPsi->LocalPsi[OnePsiA->MyLocalNo], Lev->MaxG*ElementSize, MPI_DOUBLE, k, GramSchTag3, P->Par.comm_ST_PsiT);
|
---|
913 | LPsiDatA=Lev->LPsi->LocalPsi[OnePsiA->MyLocalNo];
|
---|
914 | } // LPsiDatA is now set to the coefficients of OnePsi either stored or MPI_Received
|
---|
915 |
|
---|
916 | LocalSP = 0.0;
|
---|
917 | k=0;
|
---|
918 | if (Lev->GArray[0].GSq == 0.0) { // calculate scalar product of LPsiDatA and LPsiDatB
|
---|
919 | LocalSP += LPsiDatA[0].re*psi[0].re;
|
---|
920 | k++;
|
---|
921 | }
|
---|
922 | for (; k < Lev->MaxG; k++) {
|
---|
923 | LocalSP += 2*(LPsiDatA[k].re*psi[k].re+LPsiDatA[k].im*psi[k].im);
|
---|
924 | }
|
---|
925 | MPI_Allreduce ( &LocalSP, &PsiSP, 1, MPI_DOUBLE, MPI_SUM, P->Par.comm_ST_Psi); // gather by summation results from the group sharing the coefficients
|
---|
926 | if ((fabs(PsiSP -1.0) >= MYEPSILON) && (fabs(PsiSP) >= MYEPSILON)) {
|
---|
927 | fprintf(stderr,"(%i)(%i,Psi) = %g ?= 1.0 or 0.0 eps(%g or %g >= %g)\n",P->Par.my_color_comm_ST,OnePsiA->MyGlobalNo,PsiSP, fabs(PsiSP -1.0), fabs(PsiSP), MYEPSILON);
|
---|
928 | NotOk++;
|
---|
929 | } else
|
---|
930 | fprintf(stderr,"(%i)(%i,Psi) ok.\n",P->Par.my_color_comm_ST,OnePsiA->MyGlobalNo);
|
---|
931 | }
|
---|
932 | }
|
---|
933 | if (P->Call.out[LeaderOut]) { // final check if there have been any un-orthogonal pairs
|
---|
934 | if (NotOk == 0) {
|
---|
935 | fprintf(stderr,"(%i)TestGramSchm: Ok !\n",P->Par.my_color_comm_ST);
|
---|
936 | } else {
|
---|
937 | fprintf(stderr,"(%i)TestGramSchm: There are %i pairs not orthogonal!\n",P->Par.my_color_comm_ST,NotOk);
|
---|
938 | }
|
---|
939 | }
|
---|
940 | }
|
---|