1 | /** \file output.c
|
---|
2 | * Output of forces and energies, Visuals (density and ions for OpenDX).
|
---|
3 | *
|
---|
4 | * Herein all the functions concerning the output of data are gathered:
|
---|
5 | * Initialization of the FileData structure InitOutVisArray() and the files themselves
|
---|
6 | * InitOutputFiles(),
|
---|
7 | * Saving of a calculated state OutputVisSrcFiles() - 1. Psis OutputSrcPsiDensity(), 2. Ions
|
---|
8 | * OutSrcIons() - or retrieving Psis ReadSrcPsiDensity() and Ions ReadSrcIons(),
|
---|
9 | * Preparation (in case of RiemannTensor use) CalculateOutVisDensityPos(), OutVisPosRTransformPosNFRto0(),
|
---|
10 | * Output of visual data (for OpenDx explorer) OutputVis() - 1. OpenDX files CreateDensityOutputGeneral(),
|
---|
11 | * 2. electronic densitiy data OutVisDensity() (uses MPI sending density coefficients OutputOutVisDensity()
|
---|
12 | * and receiving and writing CombineOutVisDensity()), 3. ion data OutVisIons() -
|
---|
13 | * Closing of files CloseOutputFiles().
|
---|
14 | *
|
---|
15 | * There are some more routines: OutputCurrentDensity() outputs the current density for each magnetic field
|
---|
16 | * direction, OutputVisAllOrbital() outputs all orbitals of a certain minimsation type and
|
---|
17 | * TestReadnWriteSrcDensity() checks whether a certain minimisation group can be written and read again
|
---|
18 | * correctly.
|
---|
19 | *
|
---|
20 | * There are some helpers that open files with certain filenames, making extensive usage of all
|
---|
21 | * the suffixes defined in here: OpenFile(), OpenFileNo(), OpenFileNo2(), OpenFileNoNo() and
|
---|
22 | * OpenFileNoPost().
|
---|
23 | *
|
---|
24 | *
|
---|
25 | Project: ParallelCarParrinello
|
---|
26 | \author Jan Hamaekers, Frederik Heber
|
---|
27 | \date 2000, 2006
|
---|
28 |
|
---|
29 | File: output.c
|
---|
30 | $Id: output.c,v 1.51.2.2 2007-04-21 12:55:50 foo Exp $
|
---|
31 | */
|
---|
32 | #include<stdlib.h>
|
---|
33 | #include<stdio.h>
|
---|
34 | #include<string.h>
|
---|
35 | #include<math.h>
|
---|
36 | //#include<sys/time.h>
|
---|
37 | #include <time.h>
|
---|
38 | #include<unistd.h>
|
---|
39 | #include"data.h"
|
---|
40 | #include"density.h"
|
---|
41 | #include"errors.h"
|
---|
42 | #include"gramsch.h"
|
---|
43 | #include"helpers.h"
|
---|
44 | #include "init.h"
|
---|
45 | #include"myfft.h"
|
---|
46 | #include"mymath.h"
|
---|
47 | #include"output.h"
|
---|
48 | #include"pdbformat.h"
|
---|
49 | #include"perturbed.h"
|
---|
50 |
|
---|
51 |
|
---|
52 | /* Konvention: Rueckgabe 0 einer Funktion, bedeutet keinen Fehler (entsprechend exitcode 0) */
|
---|
53 | /* Oeffnet Datei P->mainname+"..." mit what*/
|
---|
54 |
|
---|
55 | /** Opens a file with FileData::mainname+suffix.
|
---|
56 | * Both the path under which the main programme resides and the path to the config file are
|
---|
57 | * tried subsequently.
|
---|
58 | * \param *P Problem at hand
|
---|
59 | * \param **file file pointer array
|
---|
60 | * \param *suffix suffix after mainname
|
---|
61 | * \param *what access parameter for the file: r, w, rw, r+, w+ ...
|
---|
62 | * \param verbose 1 - print status, 0 - don't print anything
|
---|
63 | * \return 0 - could not open file, 1 - file is open
|
---|
64 | */
|
---|
65 | int OpenFile(struct Problem *P, FILE** file, const char* suffix, const char* what, int verbose)
|
---|
66 | {
|
---|
67 | char* name; /* Zu erzeugender Dateiname */
|
---|
68 | name = (char*)
|
---|
69 | Malloc(strlen(P->Files.default_path) + strlen(P->Files.mainname) + strlen(suffix) + 1,"OpenFile");
|
---|
70 | sprintf(name, "%s%s%s", P->Files.default_path, P->Files.mainname, suffix);
|
---|
71 | *file = fopen(name, what);
|
---|
72 | if (*file == NULL) {
|
---|
73 | fprintf(stderr,"(%i) Normal access failed: name %s, what %s\n", P->Par.me, name, what);
|
---|
74 | /* hier default file ausprobieren, falls nur gelesen werden soll! */
|
---|
75 | if(*what == 'r') {
|
---|
76 | name = (char*)
|
---|
77 | Realloc(name,strlen(P->Files.default_path) + strlen(suffix)+1,"OpenFile");
|
---|
78 | sprintf(name, "%s%s", P->Files.default_path,suffix);
|
---|
79 | *file = fopen(name,what);
|
---|
80 | if (*file != NULL) {
|
---|
81 | if (verbose) fprintf(stderr,"(%i) Default file is open: %s\n",P->Par.me,name);
|
---|
82 | Free(name, "OpenFile: name");
|
---|
83 | return(1);
|
---|
84 | } else {
|
---|
85 | fprintf(stderr,"\n(%i)Error: Cannot open neither normal nor default file for reading: %s\n",P->Par.me,name);
|
---|
86 | Free(name, "OpenFile: name");
|
---|
87 | return(0);
|
---|
88 | }
|
---|
89 | } else {
|
---|
90 | fprintf(stderr,"\n(%i)Error: Cannot open normal file for writing: %s\n",P->Par.me,name);
|
---|
91 | Free(name, "OpenFile: name");
|
---|
92 | return(0);
|
---|
93 | }
|
---|
94 | } else {
|
---|
95 | if(verbose) fprintf(stderr,"(%i) File is open: %s\n",P->Par.me,name);
|
---|
96 | Free(name, "OpenFile: name");
|
---|
97 | return(1);
|
---|
98 | }
|
---|
99 | }
|
---|
100 |
|
---|
101 | /** Opens a file with FileData::mainname+suffix+"."+No (2 digits).
|
---|
102 | * Both the path under which the main programme resides and the path to the config file are
|
---|
103 | * tried subsequently.
|
---|
104 | * \param *P Problem at hand
|
---|
105 | * \param **file file pointer array
|
---|
106 | * \param *suffix suffix after mainname
|
---|
107 | * \param No the number with up to two digits
|
---|
108 | * \param *what access parameter for the file: r, w, rw, r+, w+ ...
|
---|
109 | * \param verbose 1 - print status, 0 - don't print anything
|
---|
110 | * \return 0 - could not open file, 1 - file is open
|
---|
111 | */
|
---|
112 | int OpenFileNo2(struct Problem *P, FILE** file, const char* suffix, int No, const char* what, int verbose)
|
---|
113 | {
|
---|
114 | char* name; /* Zu erzeugender Dateiname */
|
---|
115 | name = (char*)
|
---|
116 | Malloc(strlen(P->Files.default_path) + strlen(P->Files.mainname) + strlen(suffix) + 4,"OpenFile");
|
---|
117 | sprintf(name, "%s%s%s.%02i", P->Files.default_path, P->Files.mainname, suffix, No);
|
---|
118 | *file = fopen(name, what);
|
---|
119 | if (*file == NULL) {
|
---|
120 | /* falls nur gelesen wird, auch default_path ausprobieren */
|
---|
121 | if (*what == 'r') {
|
---|
122 | name = (char*)
|
---|
123 | Realloc(name,strlen(P->Files.default_path) + strlen(suffix) + 4,"OpenFileNo2");
|
---|
124 | sprintf(name, "%s%s.%02i", P->Files.default_path, suffix, No);
|
---|
125 | *file = fopen(name, what);
|
---|
126 | if (*file != NULL) {
|
---|
127 | if(verbose) fprintf(stderr,"(%i) Default file is open: %s\n", P->Par.me, name);
|
---|
128 | Free(name, "OpenFileNo2: name");
|
---|
129 | return(1);
|
---|
130 | }
|
---|
131 | }
|
---|
132 | fprintf(stderr,"\n(%i)Error: Cannot open file: %s\n", P->Par.me, name);
|
---|
133 | Free(name, "OpenFileNo2: name");
|
---|
134 | return(0);
|
---|
135 | } else {
|
---|
136 | if(verbose) fprintf(stderr,"(%i) File is open: %s\n",P->Par.me, name);
|
---|
137 | Free(name, "OpenFileNo2: name");
|
---|
138 | return(1);
|
---|
139 | }
|
---|
140 | }
|
---|
141 |
|
---|
142 | /* Oeffnet Datei P->Files.mainname+"...".Nr(4stellig) mit what*/
|
---|
143 | /** Opens a file with FileData::mainname+suffix+"."+No (4 digits).
|
---|
144 | * Both the path under which the main programme resides and the path to the config file are
|
---|
145 | * tried subsequently.
|
---|
146 | * \param *P Problem at hand
|
---|
147 | * \param **file file pointer array
|
---|
148 | * \param *suffix suffix after mainname
|
---|
149 | * \param No the number with up to four digits
|
---|
150 | * \param *what access parameter for the file: r, w, rw, r+, w+ ...
|
---|
151 | * \param verbose 1 - print status, 0 - don't print anything
|
---|
152 | * \return 0 - could not open file, 1 - file is open
|
---|
153 | */
|
---|
154 | int OpenFileNo(struct Problem *P, FILE** file, const char* suffix, int No, const char* what, int verbose)
|
---|
155 | {
|
---|
156 | char* name; /* Zu erzeugender Dateiname */
|
---|
157 | name = (char*)
|
---|
158 | Malloc(strlen(P->Files.default_path) + strlen(P->Files.mainname) + strlen(suffix) + 6,"OpenFileNo");
|
---|
159 | sprintf(name, "%s%s%s.%04i", P->Files.default_path, P->Files.mainname, suffix, No);
|
---|
160 | *file = fopen(name, what);
|
---|
161 | if (*file == NULL) {
|
---|
162 | /* falls nur gelesen wird, auch default_path ausprobieren */
|
---|
163 | if (*what == 'r') {
|
---|
164 | name = (char*)
|
---|
165 | Realloc(name,strlen(P->Files.default_path) + strlen(suffix) + 6,"OpenFileNo");
|
---|
166 | sprintf(name, "%s%s.%04i", P->Files.default_path, suffix, No);
|
---|
167 | *file = fopen(name, what);
|
---|
168 | if (*file != NULL) {
|
---|
169 | if(verbose) fprintf(stderr,"(%i) Default file is open: %s\n", P->Par.me, name);
|
---|
170 | Free(name, "OpenFileNo: name");
|
---|
171 | return(1);
|
---|
172 | }
|
---|
173 | }
|
---|
174 | fprintf(stderr,"\n(%i)Error: Cannot open file: %s\n", P->Par.me, name);
|
---|
175 | Free(name, "OpenFileNo: name");
|
---|
176 | return(0);
|
---|
177 | } else {
|
---|
178 | if(verbose) fprintf(stderr,"(%i) File is open: %s\n", P->Par.me, name);
|
---|
179 | Free(name, "OpenFileNo: name");
|
---|
180 | return(1);
|
---|
181 | }
|
---|
182 | }
|
---|
183 |
|
---|
184 | /* die nachfolgende Routine wird von seq und par benutzt!!! */
|
---|
185 | /* Oeffnet Datei P->Files.mainname+"...".No.postfix mit what*/
|
---|
186 |
|
---|
187 | /** Opens a file with FileData::mainname+suffix+"."+No(4 digits)+postfix.
|
---|
188 | * Only the path under which the main programme resides is tried.
|
---|
189 | * \param *P Problem at hand
|
---|
190 | * \param **file file pointer array
|
---|
191 | * \param *suffix suffix after mainname
|
---|
192 | * \param No the number with up to four digits
|
---|
193 | * \param *postfix post-suffix at the end of filename
|
---|
194 | * \param *what access parameter for the file: r, w, rw, r+, w+ ...
|
---|
195 | * \param verbose 1 - print status, 0 - don't print anything
|
---|
196 | * \return 0 - could not open file, 1 - file is open
|
---|
197 | */
|
---|
198 | int OpenFileNoPost(struct Problem *P, FILE** file, const char* suffix, int No, const char* postfix, const char* what, int verbose)
|
---|
199 | {
|
---|
200 | char* name; /* Zu erzeugender Dateiname */
|
---|
201 | name = (char*)
|
---|
202 | Malloc(strlen(P->Files.default_path) + strlen(P->Files.mainname) + strlen(postfix) + strlen(suffix) + 6,"OpenFileNoPost");
|
---|
203 | sprintf(name, "%s%s%s.%04i%s", P->Files.default_path, P->Files.mainname, suffix, No, postfix);
|
---|
204 | *file = fopen(name, what);
|
---|
205 | if (*file == NULL) {
|
---|
206 | fprintf(stderr,"\n(%i)Error: Cannot open file: %s\n", P->Par.me, name);
|
---|
207 | Free(name, "OpenFileNoPost: name");
|
---|
208 | return(0);
|
---|
209 | } else {
|
---|
210 | if(verbose) fprintf(stderr,"(%i) File is open: %s\n", P->Par.me, name);
|
---|
211 | Free(name, "OpenFileNoPost: name");
|
---|
212 | return(1);
|
---|
213 | }
|
---|
214 | }
|
---|
215 |
|
---|
216 | /** Opens a file with FileData::mainname+suffix+"."+No1(4 digits)+"."+No2(4 digits).
|
---|
217 | * Only the path under which the main programme resides is tried.
|
---|
218 | * \param *P Problem at hand
|
---|
219 | * \param **file file pointer array
|
---|
220 | * \param *suffix suffix after mainname
|
---|
221 | * \param No1 first number with up to four digits
|
---|
222 | * \param No2 second number with up to four digits
|
---|
223 | * \param *what access parameter for the file: r, w, rw, r+, w+ ...
|
---|
224 | * \param verbose 1 - print status, 0 - don't print anything
|
---|
225 | * \return 0 - could not open file, 1 - file is open
|
---|
226 | */
|
---|
227 | int OpenFileNoNo(struct Problem *P, FILE** file, const char* suffix, int No1, int No2, const char* what, int verbose)
|
---|
228 | { /* zuerst suffix; No1: lfd, No2: procId */
|
---|
229 | char *name;
|
---|
230 | name = (char*)
|
---|
231 | Malloc(strlen(P->Files.default_path) + strlen(P->Files.mainname) + strlen(suffix) + 5 + 6 + 1,"OpenFileNoNo");
|
---|
232 | sprintf(name,"%s%s%s.%04i.%04i", P->Files.default_path, P->Files.mainname, suffix, No1, No2);
|
---|
233 | *file = fopen(name, what);
|
---|
234 | if(*file == NULL) {
|
---|
235 | fprintf(stderr,"\n(%i)Error: Cannot open file: %s\n", P->Par.me, name);
|
---|
236 | Free(name, "OpenFileNoNo: name");
|
---|
237 | return(0);
|
---|
238 | }
|
---|
239 | if(verbose) fprintf(stderr,"(%i) File is open: %s\n", P->Par.me, name);
|
---|
240 | Free(name, "OpenFileNoNo: name");
|
---|
241 | return(1);
|
---|
242 | }
|
---|
243 |
|
---|
244 | /** Transformation of wave function of Riemann level to zeroth before output as Vis.
|
---|
245 | * \param *RT RiemannTensor structure, contains RiemannTensor::LatticeLevel
|
---|
246 | * \param *source source wave function array
|
---|
247 | * \param *dest destination wave function array
|
---|
248 | * \param NF dimensional factor (NDIM, generally 3)
|
---|
249 | */
|
---|
250 | static void OutVisPosRTransformPosNFRto0(const struct RiemannTensor *RT, fftw_real *source, fftw_real *dest, const int NF)
|
---|
251 | {
|
---|
252 | struct LatticeLevel *Lev = RT->LevR;
|
---|
253 | int es = Lev->NUp0[2]*NF;
|
---|
254 | unsigned int cpyes = sizeof(fftw_real)*es;
|
---|
255 | int nx=Lev->Plan0.plan->local_nx,ny=Lev->N[1],nz=Lev->N[2],Nx=Lev->NUp0[0],Ny=Lev->NUp0[1];
|
---|
256 | int lx,ly,z,Lx,Ly;
|
---|
257 | for(lx=0; lx < nx; lx++)
|
---|
258 | for(Lx=0; Lx < Nx; Lx++)
|
---|
259 | for(ly=0; ly < ny; ly++)
|
---|
260 | for(Ly=0; Ly < Ny; Ly++)
|
---|
261 | for(z=0; z < nz; z++) {
|
---|
262 | memcpy( &dest[es*(z+nz*(Ly+Ny*(ly+ny*(Lx+Nx*lx))))],
|
---|
263 | &source[es*(Ly+Ny*(Lx+Nx*(z+nz*(ly+ny*lx))))],
|
---|
264 | cpyes);
|
---|
265 | }
|
---|
266 | }
|
---|
267 |
|
---|
268 |
|
---|
269 | /** prints Norm of each wave function to screen.
|
---|
270 | * \param *out output destination (eg. stdout)
|
---|
271 | * \param *P Problem at hand
|
---|
272 | */
|
---|
273 | void OutputNorm (FILE *out, struct Problem *P) {
|
---|
274 | struct Lattice *Lat = &P->Lat;
|
---|
275 | struct Psis *Psi = &Lat->Psi;
|
---|
276 | struct LatticeLevel *Lev = P->R.LevS;
|
---|
277 | struct OnePsiElement *OnePsi, *LOnePsi;
|
---|
278 | int i;
|
---|
279 |
|
---|
280 | for (i=0;i<Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT;i++) {
|
---|
281 | OnePsi =&Psi->AllPsiStatus[i];
|
---|
282 | if (OnePsi->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) { // local one
|
---|
283 | LOnePsi = &Psi->LocalPsiStatus[OnePsi->MyLocalNo];
|
---|
284 | fprintf(out,"(%i) Norm of Psi %i: %e\n", P->Par.me, OnePsi->MyGlobalNo, GramSchGetNorm2(P,Lev,Lev->LPsi->LocalPsi[OnePsi->MyLocalNo]));
|
---|
285 | }
|
---|
286 | }
|
---|
287 | }
|
---|
288 |
|
---|
289 | /** Output of current Psis state of source level RunStruct::LevS's LevelPsi::LocalPsi to file \ref suffixsrcpsidat.
|
---|
290 | * In case of process 0, the doc file is written in a parsable way with minimisation type RunStruct#CurrentMin, level number
|
---|
291 | * LatticeLevel#LevelNo, and number of grid nodes LatticeLevel#N.
|
---|
292 | * The two (three) different Psis::SpinType's are discerned and in case of SpinUpDown separate data files opened.
|
---|
293 | *
|
---|
294 | * Then for every global wave function of the desired type each coefficient of the reciprocal grid is copied into
|
---|
295 | * Density::DensityCArray[TempDensity], fouriertransformed, copied into Density::DensityArray[TempDensity].
|
---|
296 | *
|
---|
297 | * As the file output is only handled by process 0, all other coefficient shares are sent within the Psi group to which
|
---|
298 | * the current global wave function belongs to process 0 of the Psi group and from there on to process 0 of all via
|
---|
299 | * MPI.
|
---|
300 | * \param *P Problem at hand
|
---|
301 | * \param type Current minimisation state
|
---|
302 | * \return 0 - file written, 1 - unable to open files for writing
|
---|
303 | * \note This serves as a backup file, when the process is terminated and one would like to restart it at the
|
---|
304 | * current calculation lateron, see ReadSrcPsiDensity(). Note also that it is not necessary to specify the
|
---|
305 | * same number of processes on later restart, any number may be used under the condition that the number of
|
---|
306 | * grid nodes match and that there 2 sharing wave functions in case of SpinUpDown.
|
---|
307 | * \sa ReadSrcPsiDensity() - same for Reading the coefficients, TestReadnWriteSrcDensity() - checks both routines against
|
---|
308 | * each other
|
---|
309 | */
|
---|
310 | int OutputSrcPsiDensity(struct Problem *P, enum PsiTypeTag type)
|
---|
311 | {
|
---|
312 | int i,j,k, Index, zahl, owner;
|
---|
313 | struct Lattice *Lat = &P->Lat;
|
---|
314 | //struct RunStruct *R = &P->R;
|
---|
315 | struct fft_plan_3d *plan = Lat->plan;
|
---|
316 | struct LatticeLevel *LevS = P->R.LevS;
|
---|
317 | struct Psis *Psi = &Lat->Psi;
|
---|
318 | fftw_complex *work = (fftw_complex *)LevS->Dens->DensityArray[TempDensity];
|
---|
319 | double *destpsiR = (double *)LevS->Dens->DensityArray[TempDensity];
|
---|
320 | fftw_real *srcpsiR = (fftw_real *)LevS->Dens->DensityCArray[TempDensity];
|
---|
321 | fftw_complex *srcpsiC = (fftw_complex *)LevS->Dens->DensityCArray[TempDensity];
|
---|
322 | FILE *SrcPsiData, *SrcPsiDoc;
|
---|
323 | char suffixdat[255], suffixdoc[255];
|
---|
324 | MPI_Status status;
|
---|
325 | struct OnePsiElement *OnePsiA, *LOnePsiA;
|
---|
326 | fftw_complex *LPsiDatA;
|
---|
327 | int Num = 0, colorNo = 0;
|
---|
328 | int Sent = 0, sent = 0;
|
---|
329 |
|
---|
330 | SpeedMeasure(P,ReadnWriteTime,StartTimeDo);
|
---|
331 | sprintf(suffixdat, ".%.254s.L%i", P->R.MinimisationName[type], LevS->LevelNo);
|
---|
332 | strncpy (suffixdoc, suffixdat, 255);
|
---|
333 | // for the various spin cases, output the doc-file if it's process 0
|
---|
334 | if (!(P->Par.me_comm_ST)) {
|
---|
335 | switch (Lat->Psi.PsiST) {
|
---|
336 | case SpinDouble:
|
---|
337 | colorNo = 0;
|
---|
338 | strncat (suffixdat, suffixsrcpsidat, 255-strlen(suffixdat));
|
---|
339 | strncat (suffixdoc, suffixsrcpsidoc, 255-strlen(suffixdoc));
|
---|
340 | Num = Lat->Psi.GlobalNo[PsiMaxNoDouble];
|
---|
341 | break;
|
---|
342 | case SpinUp:
|
---|
343 | colorNo = 0;
|
---|
344 | strncat (suffixdat, suffixsrcpsiupdat, 255-strlen(suffixdat));
|
---|
345 | strncat (suffixdoc, suffixsrcpsiupdoc, 255-strlen(suffixdoc));
|
---|
346 | Num = Lat->Psi.GlobalNo[PsiMaxNoUp];
|
---|
347 | break;
|
---|
348 | case SpinDown:
|
---|
349 | colorNo = 1;
|
---|
350 | strncat (suffixdat, suffixsrcpsidowndat, 255-strlen(suffixdat));
|
---|
351 | strncat (suffixdoc, suffixsrcpsidowndoc, 255-strlen(suffixdoc));
|
---|
352 | Num = Lat->Psi.GlobalNo[PsiMaxNoDown];
|
---|
353 | break;
|
---|
354 | }
|
---|
355 | if (!OpenFileNo(P, &SrcPsiData, suffixdat, colorNo, "wb",P->Call.out[ReadOut])) // open SourcePsiData as write binary
|
---|
356 | fprintf(stderr,"(%i) Error opening file with suffix %s for writing!\n",P->Par.me, suffixdat);
|
---|
357 | if (!(P->Par.my_color_comm_ST_Psi)) { // if we are process 0
|
---|
358 | if (!OpenFile(P, &SrcPsiDoc, suffixdoc, "w",P->Call.out[ReadOut])) // open the (text) doc file
|
---|
359 | fprintf(stderr,"(%i) Error opening file with suffix %s for writing!\n",P->Par.me, suffixdoc);
|
---|
360 | fprintf(SrcPsiDoc, "Mintype\t%i\n", (int)type);
|
---|
361 | fprintf(SrcPsiDoc, "LevelNo\t%i\n", LevS->LevelNo);
|
---|
362 | fprintf(SrcPsiDoc, "GridNodes\t%i\t%i\t%i\n", LevS->N[0], LevS->N[1], LevS->N[2]);
|
---|
363 | fprintf(SrcPsiDoc, "PsiNo\t%i\t%i\n", Num, P->Lat.Psi.GlobalNo[PsiMaxAdd]);
|
---|
364 | fprintf(SrcPsiDoc, "Epsilon\t%lg\t%lg\n", P->R.RelEpsTotalEnergy, P->R.RelEpsKineticEnergy);
|
---|
365 | for (i = 0; i < P->Par.Max_my_color_comm_ST_Psi; i++) {
|
---|
366 | fprintf(SrcPsiDoc, "\t%i", Lat->Psi.RealAllLocalNo[i]); // print number of local Psis
|
---|
367 | }
|
---|
368 | fprintf(SrcPsiDoc, "\n");
|
---|
369 | fclose(SrcPsiDoc);
|
---|
370 | }
|
---|
371 | }
|
---|
372 |
|
---|
373 | // send/receive around and write share of coefficient array of each wave function
|
---|
374 | MPI_Allreduce(&sent, &Sent, 1, MPI_INT, MPI_SUM, P->Par.comm_ST); // catch all at the starter line
|
---|
375 | if (P->Call.out[PsiOut]) fprintf(stderr,"(%i) me (%i/%i) \t Psi (%i/%i)\t PsiT (%i/%i)\n", P->Par.me, P->Par.me_comm_ST, P->Par.Max_me_comm_ST, P->Par.me_comm_ST_Psi, P->Par.Max_me_comm_ST_Psi, P->Par.me_comm_ST_PsiT, P->Par.Max_me_comm_ST_PsiT);
|
---|
376 | k = -1; // k is global PsiNo counter for the desired group
|
---|
377 | for (j=0; j < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; j++) { // go through all wave functions (plus the extra one for each process)
|
---|
378 | OnePsiA = &Psi->AllPsiStatus[j]; // grab OnePsiA
|
---|
379 | if (OnePsiA->PsiType == type) { // only take desired minimisation group
|
---|
380 | k++;
|
---|
381 | owner = 0; // notes down in process 0 of each psi group the owner of the next coefficient share
|
---|
382 | //fprintf(stderr,"(%i) ST_Psi: OnePsiA %i\tP->Par.me %i\n", P->Par.me,OnePsiA->my_color_comm_ST_Psi,P->Par.my_color_comm_ST_Psi);
|
---|
383 | if (OnePsiA->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) { // Belongs to my Psi group?
|
---|
384 | LOnePsiA = &Psi->LocalPsiStatus[OnePsiA->MyLocalNo];
|
---|
385 | LPsiDatA = LevS->LPsi->LocalPsi[OnePsiA->MyLocalNo];
|
---|
386 | SetArrayToDouble0((double *)srcpsiR, LevS->Dens->TotalSize*2); // zero DensityCArray[TempDensity]
|
---|
387 | for (i=0;i<LevS->MaxG;i++) { // for each every unique G grid vector
|
---|
388 | Index = LevS->GArray[i].Index;
|
---|
389 | srcpsiC[Index].re = LPsiDatA[i].re; // copy real value
|
---|
390 | srcpsiC[Index].im = LPsiDatA[i].im; // copy imaginary value
|
---|
391 | }
|
---|
392 | for (i=0; i<LevS->MaxDoubleG; i++) { // also for every doubly appearing G vector (symmetry savings)
|
---|
393 | srcpsiC[LevS->DoubleG[2*i+1]].re = srcpsiC[LevS->DoubleG[2*i]].re;
|
---|
394 | srcpsiC[LevS->DoubleG[2*i+1]].im = -srcpsiC[LevS->DoubleG[2*i]].im;
|
---|
395 | }
|
---|
396 | // do an fft transform from complex to real on these srcPsiR
|
---|
397 | fft_3d_complex_to_real(plan, LevS->LevelNo, FFTNF1, srcpsiC, work);
|
---|
398 |
|
---|
399 | for (i=0; i < LevS->Dens->LocalSizeR; i++)
|
---|
400 | destpsiR[i] = (double)srcpsiR[i];
|
---|
401 | } else
|
---|
402 | LOnePsiA = NULL;
|
---|
403 |
|
---|
404 | if (P->Par.me_comm_ST == 0) { // if we are process 0 of all, only we may access the file
|
---|
405 | for (i=0; i<P->Par.Max_me_comm_ST_Psi;i++) { // for each share of the coefficient in the PsiGroup
|
---|
406 | if (LOnePsiA == NULL) { // if it's not local, receive coefficients from correct PsiGroup (process 0 within that group)
|
---|
407 | if (MPI_Recv( destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, OnePsiA->my_color_comm_ST_Psi, OutputSrcPsiTag, P->Par.comm_ST_PsiT, &status ) != MPI_SUCCESS)
|
---|
408 | Error(SomeError, "OutputSrcPsiDensity: MPI_Recv of loaded coefficients failed!");
|
---|
409 | MPI_Get_count(&status, MPI_DOUBLE, &zahl);
|
---|
410 | if (zahl != LevS->Dens->LocalSizeR) // check number of elements
|
---|
411 | fprintf(stderr,"(%i)OutputSrcPsiDensity: MPI_Recv of loaded coefficients of GlobalNo %i, owner %i failed: Too few coefficients - %i instead of %i!\n", P->Par.me, k, i, zahl, LevS->Dens->LocalSizeR);
|
---|
412 | //else
|
---|
413 | //fprintf(stderr,"(%i)OutputSrcPsiDensity: MPI_Recv of loaded coefficients of GlobalNo %i, owner %i succeeded!\n", P->Par.me, k, i);
|
---|
414 | } else { // if it's local ...
|
---|
415 | if (i != 0) { // but share of array not for us, receive from owner process within Psi group
|
---|
416 | if (MPI_Recv( destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, i, OutputSrcPsiTag, P->Par.comm_ST_Psi, &status ) != MPI_SUCCESS)
|
---|
417 | Error(SomeError, "OutputSrcPsiDensity: MPI_Recv of loaded coefficients failed!");
|
---|
418 | MPI_Get_count(&status, MPI_DOUBLE, &zahl);
|
---|
419 | if (zahl != LevS->Dens->LocalSizeR) // check number of elements
|
---|
420 | fprintf(stderr,"(%i)OutputSrcPsiDensity: MPI_Recv of loaded coefficients of GlobalNo %i, owner %i failed: Too few coefficients - %i instead of %i!\n", P->Par.me, k, i, zahl, LevS->Dens->LocalSizeR);
|
---|
421 | //else
|
---|
422 | //fprintf(stderr,"(%i)OutputSrcPsiDensity: MPI_Recv of loaded coefficients of GlobalNo %i, owner %i succeeded!\n", P->Par.me, k, i);
|
---|
423 | } // otherwise it was our share already
|
---|
424 | }
|
---|
425 | // store the final share on disc
|
---|
426 | if ((zahl = fwrite(destpsiR, sizeof(double), (size_t)(LevS->Dens->LocalSizeR), SrcPsiData)) != (size_t)(LevS->Dens->LocalSizeR)) {
|
---|
427 | fclose(SrcPsiData);
|
---|
428 | //if (P->Par.me == 0)
|
---|
429 | fprintf(stderr, "(%i)OutputSrcPsiDensity: only %i bytes written instead of expected %i\n", P->Par.me, zahl, LevS->Dens->LocalSizeR);
|
---|
430 | Error(SomeError,"OutputSrcPsiDensity: fwrite Error");
|
---|
431 | }
|
---|
432 | }
|
---|
433 | } else { // if we are not process 0 of all, we are but a deliverer
|
---|
434 | if (LOnePsiA != NULL) { // send if it's local
|
---|
435 | if (P->Par.me_comm_ST_Psi == 0) { // if we are process 0 in the group, send final share to our process 0
|
---|
436 | for (owner = 0; owner < P->Par.Max_me_comm_ST_Psi; owner++) { // for all processes of our Psi group
|
---|
437 | if (owner != 0) { // still not our share of coefficients, receive from owner in our Psi group (increasing owner)
|
---|
438 | if (MPI_Recv( destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, owner, OutputSrcPsiTag, P->Par.comm_ST_Psi, &status ) != MPI_SUCCESS)
|
---|
439 | Error(SomeError, "OutputSrcPsiDensity: MPI_Recv of loaded coefficients failed!");
|
---|
440 | MPI_Get_count(&status, MPI_DOUBLE, &zahl);
|
---|
441 | if (zahl != LevS->Dens->LocalSizeR) // check number of elements
|
---|
442 | fprintf(stderr,"(%i)OutputSrcPsiDensity: MPI_Recv of loaded coefficients of GlobalNo %i, owner %i failed: Too few coefficients - %i instead of %i!\n", P->Par.me, k, owner, zahl, LevS->Dens->LocalSizeR);
|
---|
443 | //else
|
---|
444 | //fprintf(stderr,"(%i)OutputSrcPsiDensity: MPI_Recv of loaded coefficients of GlobalNo %i, owner %i succeeded!\n", P->Par.me, k, owner);
|
---|
445 | } else sent++; // only count sent if it was our share
|
---|
446 |
|
---|
447 | if (MPI_Send(destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, 0, OutputSrcPsiTag, P->Par.comm_ST_PsiT) != MPI_SUCCESS)
|
---|
448 | Error(SomeError, "OutputSrcPsiDensity: MPI_Send of loaded coefficients failed!");
|
---|
449 | //else
|
---|
450 | //fprintf(stderr,"(%i)OutputSrcPsiDensity: MPI_Send to process %i in PsiT group of loaded coefficients GlobalNo %i succeeded!\n", P->Par.me, 0, k);
|
---|
451 | }
|
---|
452 | } else {
|
---|
453 | sent++;
|
---|
454 | if (MPI_Send(destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, 0, OutputSrcPsiTag, P->Par.comm_ST_Psi) != MPI_SUCCESS)
|
---|
455 | Error(SomeError, "OutputSrcPsiDensity: MPI_Send of loaded coefficients failed!");
|
---|
456 | //else
|
---|
457 | //fprintf(stderr,"(%i)OutputSrcPsiDensity: MPI_Send to process %i in Psi group of loaded coefficients GlobalNo %i succeeded!\n", P->Par.me, 0, k);
|
---|
458 | }
|
---|
459 | }
|
---|
460 | // otherwise we don't have anything to do with this
|
---|
461 | }
|
---|
462 | }
|
---|
463 | }
|
---|
464 | MPI_Allreduce(&sent, &Sent, 1, MPI_INT, MPI_SUM, P->Par.comm_ST); // catch all again at finish
|
---|
465 | if (P->Call.out[PsiOut]) fprintf(stderr,"(%i) Out of %i shares %i had to be sent in total, %i from this process alone.\n", P->Par.me, P->Par.Max_me_comm_ST_Psi*Psi->NoOfPsis, Sent, sent);
|
---|
466 | if (!(P->Par.me_comm_ST))
|
---|
467 | fclose(SrcPsiData);
|
---|
468 | SpeedMeasure(P,ReadnWriteTime,StopTimeDo);
|
---|
469 |
|
---|
470 | return 0;
|
---|
471 | }
|
---|
472 |
|
---|
473 | /** Tests whether writing and successant reading of coefficient array is working correctly.
|
---|
474 | * The local wave function array is written to a disc (\sa OutputSrcPsiDensity()), the by \a wavenr
|
---|
475 | * specified coefficient array copied to OldPsiDat and afterwards read again via ReadSrcPsiDensity().
|
---|
476 | * Comparison per process of each local coefficient shows incorrect read or writes.
|
---|
477 | * \param *P Problem at hand
|
---|
478 | * \param type minimisation type array to test for read and write
|
---|
479 | * \return 1 - successful, 0 - test failed
|
---|
480 | */
|
---|
481 | int TestReadnWriteSrcDensity(struct Problem *P, enum PsiTypeTag type)
|
---|
482 | {
|
---|
483 | struct RunStruct *R = &P->R;
|
---|
484 | struct LatticeLevel *LevS = R->LevS;
|
---|
485 | struct Lattice *Lat = &P->Lat;
|
---|
486 | struct Psis *Psi = &Lat->Psi;
|
---|
487 | int i,k;
|
---|
488 | fftw_complex *destpsiC, *srcpsiC;
|
---|
489 |
|
---|
490 | //fprintf(stderr,"(%i)TestReadnWriteSrcDensity\n",P->Par.me);
|
---|
491 | // write whole array of type to disc
|
---|
492 | OutputSrcPsiDensity(P,type);
|
---|
493 | debug(P,"array written");
|
---|
494 |
|
---|
495 | // copy specified array to OldPsiDat
|
---|
496 | for (k=0; k < Lat->Psi.LocalNo; k++) // for every local wave function of type, copy coefficients
|
---|
497 | if (Psi->LocalPsiStatus[k].PsiType == type) { // ... yet only for given type
|
---|
498 | srcpsiC = LevS->LPsi->LocalPsi[k];
|
---|
499 | destpsiC = LevS->LPsi->OldLocalPsi[k - Psi->TypeStartIndex[type]];
|
---|
500 | for (i=0;i<LevS->MaxG;i++) { // for each every unique G grid vector
|
---|
501 | destpsiC[i].re = srcpsiC[i].re; // copy real value
|
---|
502 | destpsiC[i].im = srcpsiC[i].im; // copy imaginary value
|
---|
503 | }
|
---|
504 | }
|
---|
505 | debug(P,"array copied");
|
---|
506 |
|
---|
507 | // read whole array again
|
---|
508 | if (!ReadSrcPsiDensity(P,type,0,R->LevSNo))
|
---|
509 | return 0;
|
---|
510 | debug(P,"array read");
|
---|
511 |
|
---|
512 | // compare with copied array
|
---|
513 | for (k=0; k < Lat->Psi.LocalNo; k++) // for every local wave function of type, compare coefficients
|
---|
514 | if (Psi->LocalPsiStatus[k].PsiType == type) { // ... yet only for given type
|
---|
515 | srcpsiC = LevS->LPsi->LocalPsi[k];
|
---|
516 | destpsiC = LevS->LPsi->OldLocalPsi[k - Psi->TypeStartIndex[type]];
|
---|
517 | for (i=0;i<LevS->MaxG;i++) // for each every unique G grid vector
|
---|
518 | if ((fabs(destpsiC[i].re - srcpsiC[i].re) >= MYEPSILON) ||(fabs(destpsiC[i].im - srcpsiC[i].im) >= MYEPSILON)) {
|
---|
519 | fprintf(stderr,"(%i)TestReadnWriteSrcDensity: First difference at index %i - %lg+i%lg against loaded %lg+i%lg\n",P->Par.me, i, srcpsiC[i].re, srcpsiC[i].im,destpsiC[i].re,destpsiC[i].im);
|
---|
520 | return 0;
|
---|
521 | }
|
---|
522 | }
|
---|
523 | debug(P,"array compared");
|
---|
524 | fprintf(stderr,"(%i)TestReadnWriteSrcDensity: OK!\n",P->Par.me);
|
---|
525 | return 1;
|
---|
526 | }
|
---|
527 |
|
---|
528 |
|
---|
529 | /** Read Psis state from an earlier run.
|
---|
530 | * The doc file is opened, mesh sizes LatticeLevel::N[], global number of Psis read and checked against the known
|
---|
531 | * values in the inital level RunStruct::LevS.
|
---|
532 | * Note, only process 0 handles the files, all read coefficients are transfered to their respective owners via MPI
|
---|
533 | * afterwards. Here, process 0 of a certain Psi group is used as a transferer for the coefficients of the other processes
|
---|
534 | * in this Psi group. He receives them all from process 0 and sends them onward accordingly. The complete set of
|
---|
535 | * coefficients on the real grid for one wave function in the Psi group are transformed into complex wave function
|
---|
536 | * coefficients by the usual fft procedure (see ChangeToLevUp()).
|
---|
537 | *
|
---|
538 | * \param *P Problem at hand
|
---|
539 | * \param type minimisation type to read
|
---|
540 | * \param test whether to just test for presence of files (1) or really read them (0)
|
---|
541 | * \param LevSNo level number to be read
|
---|
542 | * \note This is the counterpart to OutputSrcPsiDensity().
|
---|
543 | * \return 1 - success, 0 - failure
|
---|
544 | * \note It is not necessary to specify the same number of processes on later restart, any number may be used under
|
---|
545 | * the condition that the number of grid nodes match and that there at least 2 processes sharing wave functions
|
---|
546 | * in case of SpinUpDown.
|
---|
547 | * \sa OutputSrcPsiDensity() - same for writing the coefficients, TestReadnWriteSrcDensity() - checks both routines against
|
---|
548 | * each other
|
---|
549 | */
|
---|
550 | int ReadSrcPsiDensity(struct Problem *P, enum PsiTypeTag type, int test, int LevSNo)
|
---|
551 | {
|
---|
552 | int i, j, k, Index, owner;
|
---|
553 | struct RunStruct *R = &P->R;
|
---|
554 | struct Lattice *Lat = &P->Lat;
|
---|
555 | struct fft_plan_3d *plan = Lat->plan;
|
---|
556 | struct LatticeLevel *LevS = R->LevS; // keep open for LevelNo read from file
|
---|
557 | struct Psis *Psi = &Lat->Psi;
|
---|
558 | //struct Energy *E = Lat->E;
|
---|
559 | fftw_complex *work;
|
---|
560 | double *destpsiR;
|
---|
561 | fftw_real *srcpsiR;
|
---|
562 | fftw_complex *srcpsiC;
|
---|
563 | FILE *SrcPsiData, *SrcPsiDoc;
|
---|
564 | int N[NDIM], GlobalNo[2];
|
---|
565 | int LevelNo, readnr=0;
|
---|
566 | int zahl, signal = test ? 1 : 2; // 0 - ok, 1 - test failed, 2 - throw Error
|
---|
567 | char suffixdat[255], suffixdoc[255];
|
---|
568 | int read_type, Num = 0, colorNo = 0;
|
---|
569 | char spin[20];
|
---|
570 | double Eps[2];
|
---|
571 | MPI_Status status;
|
---|
572 | struct OnePsiElement *OnePsiA, *LOnePsiA;
|
---|
573 | int Recv=0, recv=0;
|
---|
574 |
|
---|
575 | SpeedMeasure(P,ReadnWriteTime,StartTimeDo);
|
---|
576 | sprintf(suffixdat, ".%.254s.L%i", P->R.MinimisationName[type], LevSNo);
|
---|
577 | strncpy (suffixdoc, suffixdat, 255);
|
---|
578 | // Depending on Psis::SpinType the source psi doc file is opened and header written
|
---|
579 | switch (Lat->Psi.PsiST) {
|
---|
580 | case SpinDouble:
|
---|
581 | colorNo = 0;
|
---|
582 | strncat (suffixdat, suffixsrcpsidat, 255-strlen(suffixdat));
|
---|
583 | strncat (suffixdoc, suffixsrcpsidoc, 255-strlen(suffixdoc));
|
---|
584 | strncpy (spin, "GlobalNoSpinDouble", 20);
|
---|
585 | Num = Lat->Psi.GlobalNo[PsiMaxNoDouble];
|
---|
586 | break;
|
---|
587 | case SpinUp:
|
---|
588 | colorNo = 0;
|
---|
589 | strncat (suffixdat, suffixsrcpsiupdat, 255-strlen(suffixdat));
|
---|
590 | strncat (suffixdoc, suffixsrcpsiupdoc, 255-strlen(suffixdoc));
|
---|
591 | strncpy (spin, "GlobalNoSpinUp", 20);
|
---|
592 | Num = Lat->Psi.GlobalNo[PsiMaxNoUp];
|
---|
593 | break;
|
---|
594 | case SpinDown:
|
---|
595 | colorNo = 1;
|
---|
596 | strncat (suffixdat, suffixsrcpsidowndat, 255-strlen(suffixdat));
|
---|
597 | strncat (suffixdoc, suffixsrcpsidowndoc, 255-strlen(suffixdoc));
|
---|
598 | strncpy (spin, "GlobalNoSpinDown", 20);
|
---|
599 | Num = Lat->Psi.GlobalNo[PsiMaxNoDown];
|
---|
600 | break;
|
---|
601 | }
|
---|
602 | // open doc file ...
|
---|
603 | if (!(P->Par.me_comm_ST)) {
|
---|
604 | if (!OpenFile(P, &SrcPsiDoc, suffixdoc, "r", test ? 0 : P->Call.out[ReadOut])) { // open doc file
|
---|
605 | debug(P,"ReadSrcPsiDensity: doc file pointer NULL\n");
|
---|
606 | if (test) {
|
---|
607 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
608 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
609 | return 0;
|
---|
610 | }
|
---|
611 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
612 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
613 | Error(SomeError,"ReadSrcPsiDensity: cannot open doc file!");
|
---|
614 | }
|
---|
615 | // ... and parse critical ...
|
---|
616 | readnr += ParseForParameter(0,SrcPsiDoc,"Mintype",0,1,1,int_type,(int *)&read_type, 1, test ? optional : critical);
|
---|
617 | readnr += ParseForParameter(0,SrcPsiDoc,"LevelNo",0,1,1,int_type,&LevelNo,1, test ? optional : critical);
|
---|
618 | readnr += 3*ParseForParameter(0,SrcPsiDoc,"GridNodes",0,3,1,row_int,&N[0], 1, test ? optional : critical);
|
---|
619 | readnr += 2*ParseForParameter(0,SrcPsiDoc,"PsiNo",0,2,1,row_int,&GlobalNo[0], 1, test ? optional : critical);
|
---|
620 | // and optional items ...
|
---|
621 | if (ParseForParameter(0,SrcPsiDoc,"Epsilon",0,2,1,row_double,&Eps[0],1,optional))
|
---|
622 | if ((P->Call.ReadSrcFiles == 1) && ((Eps[1] < R->RelEpsKineticEnergy) || (Eps[0] < R->RelEpsTotalEnergy))) {
|
---|
623 | //fprintf(stderr,"(%i) Eps %lg %lg\tRelEps %lg %lg\n", P->Par.me, Eps[0], Eps[1], R->RelEpsTotalEnergy, R->RelEpsKineticEnergy);
|
---|
624 | fprintf(stderr,"(%i) NOTE: Doing minimization after source file parsing due to smaller specified epsilon stop conditions.\n",P->Par.me);
|
---|
625 | P->Call.ReadSrcFiles = 2; // do minimisation even after loading
|
---|
626 | }
|
---|
627 | if (readnr != 7) { // check number of items read
|
---|
628 | debug(P, "ReadSrcPsiDensity: too few doc items in file\n");
|
---|
629 | fclose(SrcPsiDoc);
|
---|
630 | if (test) {
|
---|
631 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
632 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
633 | return 0;
|
---|
634 | }
|
---|
635 | fprintf(stderr,"ReadSrcPsiDensity: Only %i items read!\n",readnr);
|
---|
636 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
637 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
638 | Error(SomeError, "ReadSrcPsiDensity: read error");
|
---|
639 | }
|
---|
640 | // check if levels match
|
---|
641 | if (LevSNo != LevelNo) {
|
---|
642 | debug(P,"ReadSrcPsiDensity: mismatching levels\n");
|
---|
643 | fclose(SrcPsiDoc);
|
---|
644 | if (test) {
|
---|
645 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
646 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
647 | return 0;
|
---|
648 | }
|
---|
649 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
650 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
651 | Error(SomeError,"ReadSrcPsiDensity: Mismatching levels!");
|
---|
652 | } else {
|
---|
653 | LevS = &P->Lat.Lev[LevelNo];
|
---|
654 | }
|
---|
655 |
|
---|
656 | // check if systems in memory and file match
|
---|
657 | if ((read_type != R->CurrentMin) || (N[0] != LevS->N[0]) || (N[1] != LevS->N[1]) || (N[2] != LevS->N[2]) || (GlobalNo[0] != Num) || (GlobalNo[1] != Lat->Psi.GlobalNo[PsiMaxAdd])) { // check read system
|
---|
658 | debug(P,"ReadSrcPsiDensity: srcpsi file does not fit to system\n");
|
---|
659 | fclose(SrcPsiDoc);
|
---|
660 | if (test) {
|
---|
661 | fprintf(stderr,"(%i) Min\t N(x,y,z)\tPsiNo+AddNo\n file: %s\t %i %i %i\t %i + %i\nsystem: %s\t %d %d %d\t %d + %d\n",P->Par.me, R->MinimisationName[read_type], N[0], N[1], N[2], GlobalNo[0], GlobalNo[1], R->MinimisationName[R->CurrentMin], LevS->N[0] , LevS->N[1], LevS->N[2], Lat->Psi.GlobalNo[PsiMaxNoDouble], Lat->Psi.GlobalNo[PsiMaxAdd]);
|
---|
662 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
663 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
664 | return 0;
|
---|
665 | }
|
---|
666 | fprintf(stderr,"ReadSrcPsiDensity: Type %i != CurrentMin %i || N[0] %i != %i || N[1] %i != %i || N[2] %i != %i || %s %i + %i != %i + %i\n", read_type, R->CurrentMin, N[0], LevS->N[0], N[1], LevS->N[1], N[2], LevS->N[2], spin, GlobalNo[0], GlobalNo[1], Num, P->Lat.Psi.GlobalNo[PsiMaxAdd]);
|
---|
667 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
668 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
669 | Error(SomeError,"ReadSrcPsiDensity: srcpsi file does not fit to system");
|
---|
670 | }
|
---|
671 | signal = 0; // everything went alright, signal ok
|
---|
672 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
673 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
674 | } else { // others wait for signal from root process
|
---|
675 | if (MPI_Bcast(&signal,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
676 | Error(SomeError,"ReadSrcPsiDensity: Bcast of signal failed\n");
|
---|
677 | if (signal == 1)
|
---|
678 | return 0;
|
---|
679 | else if (signal == 2)
|
---|
680 | Error(SomeError, "ReadSrcPsiDensity: Something went utterly wrong, see root process");
|
---|
681 | else if (P->Call.out[PsiOut])
|
---|
682 | fprintf(stderr,"(%i) ReadSrcPsiDensity: Everything went alright so far\n", P->Par.me);
|
---|
683 | }
|
---|
684 | if (MPI_Bcast(&LevelNo,1,MPI_INT,0,P->Par.comm_ST) != MPI_SUCCESS)
|
---|
685 | Error(SomeError,"ReadSrcPsiDensity: Bcast of LevelNo failed\n");
|
---|
686 | LevS = &P->Lat.Lev[LevelNo];
|
---|
687 | //if (!test) fprintf(stderr,"(%i) LevelSNo %i\n", P->Par.me, LevS->LevelNo);
|
---|
688 |
|
---|
689 | if (!test) {
|
---|
690 | // set some pointers for work to follow
|
---|
691 | work = (fftw_complex *)LevS->Dens->DensityArray[TempDensity];
|
---|
692 | destpsiR = (double *)LevS->Dens->DensityArray[TempDensity];
|
---|
693 | srcpsiR = (fftw_real *)LevS->Dens->DensityCArray[TempDensity];
|
---|
694 | srcpsiC = (fftw_complex *)LevS->Dens->DensityCArray[TempDensity];
|
---|
695 |
|
---|
696 | // read share of coefficient array for each wave function and send/receive around
|
---|
697 | owner = 0;
|
---|
698 | MPI_Allreduce (&recv, &Recv, 1, MPI_INT, MPI_SUM, P->Par.comm_ST);
|
---|
699 | //fprintf(stderr,"(%i) me (%i/%i) \t Psi (%i/%i)\t PsiT (%i/%i)\n", P->Par.me, P->Par.me_comm_ST, P->Par.Max_me_comm_ST, P->Par.me_comm_ST_Psi, P->Par.Max_me_comm_ST_Psi, P->Par.me_comm_ST_PsiT, P->Par.Max_me_comm_ST_PsiT);
|
---|
700 | k = -1; // k is global PsiNo counter for the desired group
|
---|
701 | for (j=0; j < Psi->MaxPsiOfType+P->Par.Max_me_comm_ST_PsiT; j++) { // go through all wave functions (plus the extra one for each process)
|
---|
702 | OnePsiA = &Psi->AllPsiStatus[j]; // grab OnePsiA
|
---|
703 | if (OnePsiA->PsiType == type) { // only take desired minimisation group
|
---|
704 | k++;
|
---|
705 | //fprintf(stderr,"(%i) ST_Psi: OnePsiA %i\tP->Par.my_color_comm_ST_Psi %i\n", P->Par.me,OnePsiA->my_color_comm_ST_Psi,P->Par.my_color_comm_ST_Psi);
|
---|
706 | if (OnePsiA->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // Belongs to my Psi group?
|
---|
707 | LOnePsiA = &Psi->LocalPsiStatus[OnePsiA->MyLocalNo];
|
---|
708 | else
|
---|
709 | LOnePsiA = NULL;
|
---|
710 |
|
---|
711 | if (P->Par.me_comm_ST == 0) { // if we are process 0 of all, we may access file
|
---|
712 | if (!OpenFileNo(P, &SrcPsiData, suffixdat, colorNo, "r", test ? 0 : P->Call.out[ReadOut])) {
|
---|
713 | Error(SomeError,"ReadSrcPsiDensity: cannot open data file");
|
---|
714 | }
|
---|
715 | for (i=P->Par.Max_me_comm_ST_Psi-1; i>=0;i--) { // load coefficients
|
---|
716 | fseek( SrcPsiData, (long)((k*N[0]*N[1]*N[2]+i*((long)LevS->Dens->LocalSizeR))*sizeof(double)), SEEK_SET); // seek to beginning of this process' coefficients
|
---|
717 | // readin
|
---|
718 | if ((zahl = fread(destpsiR, sizeof(double), (size_t)(LevS->Dens->LocalSizeR), SrcPsiData)) != (size_t)LevS->Dens->LocalSizeR) {
|
---|
719 | fclose(SrcPsiData);
|
---|
720 | fprintf(stderr, "(%i)ReadSrcPsiDensity: only %i bytes read instead of expected %i\n", P->Par.me, zahl, LevS->Dens->LocalSizeR);
|
---|
721 | Error(SomeError,"ReadSrcPsiDensity: fread Error");
|
---|
722 | }
|
---|
723 | if (LOnePsiA == NULL) { // if it's not local, send away coefficients to correct PsiGroup (process 0 within that group)
|
---|
724 | if (MPI_Send(destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, OnePsiA->my_color_comm_ST_Psi, ReadSrcPsiTag, P->Par.comm_ST_PsiT) != MPI_SUCCESS)
|
---|
725 | Error(SomeError, "ReadSrcPsiDensity: MPI_Send of loaded coefficients failed!");
|
---|
726 | //else
|
---|
727 | //fprintf(stderr,"(%i)ReadSrcPsiDensity: MPI_Send to process %i of loaded coefficients GlobalNo %i, owner %i succeeded!\n", P->Par.me, OnePsiA->my_color_comm_ST_Psi, k, i);
|
---|
728 | } else { // if it's local ...
|
---|
729 | if (i != 0) { // but share of array not for us, send to owner process within Psi group
|
---|
730 | if (MPI_Send(destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, i, ReadSrcPsiTag, P->Par.comm_ST_Psi) != MPI_SUCCESS)
|
---|
731 | Error(SomeError, "ReadSrcPsiDensity: MPI_Send within Psi group of loaded coefficients failed!");
|
---|
732 | //else
|
---|
733 | //fprintf(stderr,"(%i)ReadSrcPsiDensity: MPI_Send to process %i within Psi group of loaded coefficients GlobalNo %i succeeded!\n", P->Par.me, i, k);
|
---|
734 | } // otherwise it was our share already
|
---|
735 | }
|
---|
736 | }
|
---|
737 | } else {
|
---|
738 | if (LOnePsiA != NULL) { // receive
|
---|
739 | if (P->Par.me_comm_ST_Psi == 0) { // if we are process 0 in the group, receive share from process 0 of all
|
---|
740 | for (owner = P->Par.Max_me_comm_ST_Psi -1; owner >=0; owner--) { // for all processes of our Psi group
|
---|
741 | if (MPI_Recv(destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, 0, ReadSrcPsiTag, P->Par.comm_ST_PsiT, &status) != MPI_SUCCESS)
|
---|
742 | Error(SomeError, "ReadSrcPsiDensity: MPI_Recv of loaded coefficients failed!");
|
---|
743 | MPI_Get_count(&status, MPI_DOUBLE, &zahl);
|
---|
744 | if (zahl != LevS->Dens->LocalSizeR) // check number of elements
|
---|
745 | fprintf(stderr,"(%i)ReadSrcPsiDensity: MPI_Recv from process 0 of loaded coefficients of GlobalNo %i, owner %i failed: Too few coefficients - %i instead of %i!\n", P->Par.me, k, owner, zahl, LevS->Dens->LocalSizeR);
|
---|
746 | //else
|
---|
747 | //fprintf(stderr,"(%i)ReadSrcPsiDensity: MPI_Recv from process 0 of loaded coefficients of GlobalNo %i, owner %i succeeded!\n", P->Par.me, k, owner);
|
---|
748 |
|
---|
749 | if (owner != 0) { // not our share of coefficients, send to owner in our Psi group
|
---|
750 | if (MPI_Send(destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, owner, ReadSrcPsiTag, P->Par.comm_ST_Psi) != MPI_SUCCESS)
|
---|
751 | Error(SomeError, "ReadSrcPsiDensity: MPI_Send within Psi group of loaded coefficients failed!");
|
---|
752 | //else
|
---|
753 | //fprintf(stderr,"(%i)ReadSrcPsiDensity: MPI_Send to process %i within Psi group of loaded coefficients GlobalNo %i succeeded!\n", P->Par.me, owner, k);
|
---|
754 | } else recv++;
|
---|
755 | }
|
---|
756 | // otherwise it's our share!
|
---|
757 | } else { // our share within Psi Group not belonging to process 0 of all
|
---|
758 | recv++;
|
---|
759 | if (MPI_Recv(destpsiR, LevS->Dens->LocalSizeR, MPI_DOUBLE, 0, ReadSrcPsiTag, P->Par.comm_ST_Psi, &status) != MPI_SUCCESS)
|
---|
760 | Error(SomeError, "ReadSrcPsiDensity: MPI_Recv of loaded coefficients failed!");
|
---|
761 | MPI_Get_count(&status, MPI_DOUBLE, &zahl);
|
---|
762 | if (zahl != LevS->Dens->LocalSizeR) // check number of elements
|
---|
763 | fprintf(stderr,"(%i)ReadSrcPsiDensity: MPI_Recv of loaded coefficients of GlobalNo %i, owner %i failed: Too few coefficients - %i instead of %i!\n", P->Par.me, k, P->Par.me_comm_ST_Psi, zahl, LevS->Dens->LocalSizeR);
|
---|
764 | //else
|
---|
765 | //fprintf(stderr,"(%i)ReadSrcPsiDensity: MPI_Recv of loaded coefficients of GlobalNo %i, owner %i succeeded!\n", P->Par.me, k, P->Par.me_comm_ST_Psi);
|
---|
766 | }
|
---|
767 | }
|
---|
768 | // otherwise we don't have anything to do with this
|
---|
769 | }
|
---|
770 |
|
---|
771 | if (LOnePsiA != NULL) {
|
---|
772 | SetArrayToDouble0((double *)srcpsiR, LevS->Dens->TotalSize*2);
|
---|
773 | for (i=0; i < LevS->Dens->LocalSizeR; i++) // copy dest to src
|
---|
774 | srcpsiR[i] = (fftw_real)destpsiR[i];
|
---|
775 |
|
---|
776 | fft_3d_real_to_complex(plan, LevS->LevelNo, FFTNF1, srcpsiC, work); // fft transform
|
---|
777 | //if (P->Call.out[PsiOut])
|
---|
778 | //fprintf(stderr,"(%i) LevSNo %i\t LocalPsi %p\n", P->Par.me, LevS->LevelNo, LevS->LPsi->LocalPsi[LOnePsiA->MyLocalNo]);
|
---|
779 | for (i=0;i<LevS->MaxG;i++) { // and copy into wave functions coefficients
|
---|
780 | Index = LevS->GArray[i].Index;
|
---|
781 | LevS->LPsi->LocalPsi[LOnePsiA->MyLocalNo][i].re = srcpsiC[Index].re/LevS->MaxN;
|
---|
782 | LevS->LPsi->LocalPsi[LOnePsiA->MyLocalNo][i].im = srcpsiC[Index].im/LevS->MaxN;
|
---|
783 | }
|
---|
784 | }
|
---|
785 | if ((P->Par.me_comm_ST == 0) && (SrcPsiData != NULL)) fclose(SrcPsiData);
|
---|
786 | }
|
---|
787 | }
|
---|
788 | MPI_Allreduce (&recv, &Recv, 1, MPI_INT, MPI_SUM, P->Par.comm_ST);
|
---|
789 | if (P->Call.out[PsiOut]) fprintf(stderr,"(%i) Out of %i shares %i had to be received in total, %i from this process alone.\n", P->Par.me, P->Par.Max_me_comm_ST_Psi*Psi->NoOfPsis, Recv, recv);
|
---|
790 | SpeedMeasure(P,ReadnWriteTime,StopTimeDo);
|
---|
791 | }
|
---|
792 | return 1; // everything went well till the end
|
---|
793 | }
|
---|
794 |
|
---|
795 | /** Creates the density \ref suffixdensdoc and \ref suffixdensdx files for OpenDx.
|
---|
796 | * Opens \ref suffixdensdoc, fills (pos&data file name, byte order, max mesh points, matrix alignment, steps)
|
---|
797 | * and closes it.
|
---|
798 | * Opens \ref suffixdensdx, then for every FileData::OutVisStep a describing structure for DX is written and
|
---|
799 | * the file closed again.
|
---|
800 | * \param *P Problem at hand
|
---|
801 | * \param me my process number in communicator Psi (0 - do nothing, else - do)
|
---|
802 | */
|
---|
803 | static void CreateDensityOutputGeneral(struct Problem *P, const int me)
|
---|
804 | {
|
---|
805 | FILE *DensityDoc, *DensityDx;
|
---|
806 | char *posname, *datname, *suffix;
|
---|
807 | struct LatticeLevel *Lev = &P->Lat.Lev[STANDARTLEVEL];
|
---|
808 | unsigned int i, MaxPoints, N[NDIM];
|
---|
809 | double *RB = P->Lat.RealBasis;
|
---|
810 | if (me) return;
|
---|
811 | N[0] = Lev->N[0]*Lev->NUp[0];
|
---|
812 | N[1] = Lev->N[1]*Lev->NUp[1];
|
---|
813 | N[2] = Lev->N[2]*Lev->NUp[2];
|
---|
814 | MaxPoints = (N[0]+1)*(N[1]+1)*(N[2]+1);
|
---|
815 | posname = (char*)
|
---|
816 | Malloc(strlen(P->Files.mainname) + strlen(suffixdenspos) + 1,"OpenFile");
|
---|
817 | sprintf(posname, "%s%s", P->Files.mainname, suffixdenspos);
|
---|
818 | datname = (char*)
|
---|
819 | Malloc(strlen(P->Files.mainname) + strlen(suffixdensdat) + 1,"OpenFile");
|
---|
820 | sprintf(datname, "%s%s", P->Files.mainname, suffixdensdat);
|
---|
821 | // write doc file
|
---|
822 | suffix = (char *)
|
---|
823 | Malloc(strlen(suffixdensdoc) + 3 + 1,"CreateDensityOutputGeneral: suffix");
|
---|
824 | sprintf(suffix, ".L%i%s", Lev->LevelNo, suffixdensdoc);
|
---|
825 | OpenFile(P, &DensityDoc, suffix, "w",P->Call.out[ReadOut]);
|
---|
826 | fprintf(DensityDoc,"DensityPositions file = %s.####\n", posname);
|
---|
827 | fprintf(DensityDoc,"DensityData file = %s.####\n", datname);
|
---|
828 | fprintf(DensityDoc,"format = ieee float (Bytes %lu) %s\n",(unsigned long) sizeof(float),msb);
|
---|
829 | fprintf(DensityDoc,"points = %i\n", MaxPoints);
|
---|
830 | fprintf(DensityDoc,"majority = row\n");
|
---|
831 | fprintf(DensityDoc,"TimeSeries = %i\n",P->Files.OutVisStep+1);
|
---|
832 | fclose(DensityDoc);
|
---|
833 | Free(suffix, "CreateDensityOutputGeneral: suffix");
|
---|
834 | // write DX file
|
---|
835 | suffix = (char *)
|
---|
836 | Malloc(strlen(suffixdensdx) + 3 + 1,"CreateDensityOutputGeneral: suffix");
|
---|
837 | sprintf(suffix, ".L%i%s", Lev->LevelNo, suffixdensdx);
|
---|
838 | OpenFile(P, &DensityDx, suffix, "w",P->Call.out[ReadOut]);
|
---|
839 | for (i=0; i < (unsigned int)P->Files.OutVisStep+1; i++) { // for every OutVis step
|
---|
840 | if (i==0) {
|
---|
841 | fprintf(DensityDx,"object \"gridcon\" class gridconnections counts %i %i %i\n\n",(N[0]+1),(N[1]+1),(N[2]+1));
|
---|
842 | if (P->Files.OutputPosType[i] != active)
|
---|
843 | fprintf(DensityDx, "object \"posdens\" class gridpositions counts %i %i %i\norigin 0.0 0.0 0.0\ndelta %f %f %f\ndelta %f %f %f\ndelta %f %f %f\n\n",
|
---|
844 | (N[0]+1),(N[1]+1),(N[2]+1),
|
---|
845 | (float)(RB[0]/N[0]),(float)(RB[1]/N[0]),(float)RB[2]/N[0],
|
---|
846 | (float)(RB[3]/N[1]),(float)(RB[4]/N[1]),(float)RB[5]/N[1],
|
---|
847 | (float)(RB[6]/N[2]),(float)(RB[7]/N[2]),(float)RB[8]/N[2]);
|
---|
848 | }
|
---|
849 | if (P->Files.OutputPosType[i] == active) {
|
---|
850 | fprintf(DensityDx,
|
---|
851 | "object \"pos.%04u\" class array type float rank 1 shape 3 items %i %s binary\n",i,MaxPoints,msb);
|
---|
852 | fprintf(DensityDx,"data file %s.%04u,0\n",posname,i);
|
---|
853 | }
|
---|
854 | fprintf(DensityDx,"attribute \"dep\" string \"positions\"\n");
|
---|
855 | fprintf(DensityDx,"# %lu - %lu Bytes\n\n",MaxPoints*i*(unsigned long)sizeof(float)*NDIM,MaxPoints*(i+1)*(unsigned long)sizeof(float)*NDIM-1);
|
---|
856 |
|
---|
857 | fprintf(DensityDx,"object \"dat.%04u\" class array type float rank 0 items %i %s binary\n",i,MaxPoints,msb);
|
---|
858 | fprintf(DensityDx,"data file %s.%04u,0\n",datname,i);
|
---|
859 | fprintf(DensityDx,"attribute \"dep\" string \"positions\"\n");
|
---|
860 | fprintf(DensityDx,"# %lu - %lu Bytes\n\n",MaxPoints*i*(unsigned long)sizeof(float),MaxPoints*(i+1)*(unsigned long)sizeof(float)-1);
|
---|
861 |
|
---|
862 | fprintf(DensityDx,"object \"obj.%04u\" class field\n",i);
|
---|
863 | if (P->Files.OutputPosType[i] == active)
|
---|
864 | fprintf(DensityDx,"component \"positions\" \"pos.%04i\"\n",i);
|
---|
865 | if (P->Files.OutputPosType[i] != active)
|
---|
866 | fprintf(DensityDx,"component \"positions\" \"posdens\"\n");
|
---|
867 | fprintf(DensityDx,"component \"connections\" \"gridcon\"\n");
|
---|
868 | fprintf(DensityDx,"component \"data\" \"dat.%04i\"\n",i);
|
---|
869 | }
|
---|
870 | fprintf(DensityDx,"\nobject \"series\" class series\n");
|
---|
871 | for (i=0; i < (unsigned int)P->Files.OutVisStep+1; i++)
|
---|
872 | fprintf(DensityDx,"member %i \"obj.%04u\" position %f\n",i,i,(float)i);
|
---|
873 | fprintf(DensityDx,"end\n");
|
---|
874 | fclose(DensityDx);
|
---|
875 | Free(suffix, "CreateDensityOutputGeneral: suffix");
|
---|
876 |
|
---|
877 | Free(posname, "CreateDensityOutputGeneral: posname");
|
---|
878 | Free(datname, "CreateDensityOutputGeneral: datname");
|
---|
879 | }
|
---|
880 |
|
---|
881 | /** Calculates the OutVis density of the RiemannTensor level.
|
---|
882 | * The usual pattern arises when a density is fftransformed:
|
---|
883 | * -# over all grid vectors up to MaxG
|
---|
884 | * -# over all doubly grid vectors up to MaxDoubleG
|
---|
885 | * -# call to fft_3d_complex_to_real()
|
---|
886 | *
|
---|
887 | * In this case here followed by call to OutVisPosRTransformPosNFRto0() and finally FileData::work
|
---|
888 | * is copied to FileData::PosR.
|
---|
889 | * \param *Lat Lattice structure, containing Lattice::plan and LatticeLevel
|
---|
890 | * \param *F FileData structure, containing FileData::PosC, FileData::PosR, FileData::work, FileData::Totalsize, FileData::LocalSizeR
|
---|
891 | */
|
---|
892 | static void CalculateOutVisDensityPos(struct Lattice *Lat, struct FileData *F/*, const double FactorC_R, const double FactorR_C*/)
|
---|
893 | {
|
---|
894 | struct fft_plan_3d *plan = Lat->plan;
|
---|
895 | struct RiemannTensor *RT = &Lat->RT;
|
---|
896 | struct LatticeLevel *LevR = RT->LevR;
|
---|
897 | fftw_complex *destC = F->PosC;
|
---|
898 | fftw_real *destR = F->PosR;
|
---|
899 | fftw_complex *work = F->work;
|
---|
900 | fftw_real *workR = (fftw_real*)work;
|
---|
901 | fftw_complex *PosFactor = F->PosFactor;
|
---|
902 | fftw_complex *posfac, *destpos, *destRCS, *destRCD;
|
---|
903 | fftw_complex *coeff;
|
---|
904 | fftw_complex source;
|
---|
905 | int i, Index, pos, n;
|
---|
906 | int NF = NDIM, MaxNUp = F->MaxNUp, TotalSize = F->TotalSize, LocalSizeR = F->LocalSizeR;
|
---|
907 | SetArrayToDouble0((double *)destC, TotalSize*2);
|
---|
908 | for (i=0;i < LevR->MaxG;i++) {
|
---|
909 | Index = LevR->GArray[i].Index;
|
---|
910 | posfac = &PosFactor[MaxNUp*i];
|
---|
911 | destpos = &destC[MaxNUp*Index*NF];
|
---|
912 | coeff = &RT->Coeff[i];
|
---|
913 | for (pos=0; pos < MaxNUp; pos++) {
|
---|
914 | for (n=0; n < NF; n++) {
|
---|
915 | source.re = coeff[n].re;
|
---|
916 | source.im = coeff[n].im;
|
---|
917 | destpos[n+NF*pos].re = source.re*posfac[pos].re-source.im*posfac[pos].im;
|
---|
918 | destpos[n+NF*pos].im = source.re*posfac[pos].im+source.im*posfac[pos].re;
|
---|
919 | }
|
---|
920 | }
|
---|
921 | }
|
---|
922 | for (i=0; i < LevR->MaxDoubleG; i++) {
|
---|
923 | destRCS = &destC[LevR->DoubleG[2*i]*MaxNUp*NF];
|
---|
924 | destRCD = &destC[LevR->DoubleG[2*i+1]*MaxNUp*NF];
|
---|
925 | for (pos=0; pos < MaxNUp; pos++) {
|
---|
926 | for (n=0; n < NF; n++) {
|
---|
927 | destRCD[n+NF*pos].re = destRCS[n+NF*pos].re;
|
---|
928 | destRCD[n+NF*pos].im = -destRCS[n+NF*pos].im;
|
---|
929 | }
|
---|
930 | }
|
---|
931 | }
|
---|
932 | fft_3d_complex_to_real(plan, LevR->LevelNo, FFTNFRVecUp0, destC, work);
|
---|
933 | OutVisPosRTransformPosNFRto0(RT, destR, workR, NF);
|
---|
934 | memcpy(destR,workR,sizeof(fftw_real)*LocalSizeR);
|
---|
935 | }
|
---|
936 |
|
---|
937 | /** Prepare Density::DensityArray for output.
|
---|
938 | * Into FileData::work subsequently each node (all z, all y, all x) is written as \f$\log(1+x)\f$,
|
---|
939 | * where x is Density::DensityArray[TotalDensity]. In the end result is send to process 0 (yet not
|
---|
940 | * received here, see CombineOutVisArray()). In case of RiemannTensor use, some more complex calculations
|
---|
941 | * are made: FileData::PosR is used, the coefficient offset'ed to the current node and the log taken there.
|
---|
942 | * \param *P Problem at hand
|
---|
943 | * \param myPE this ranks process in the Psi communcator ParallelSimulationData::me_comm_ST_Psi
|
---|
944 | * \param *srcdens Pointer to DensityArray which is to be displayed
|
---|
945 | */
|
---|
946 | static void OutputOutVisDensity(struct Problem *P, const int myPE, fftw_real *srcdens)
|
---|
947 | {
|
---|
948 | int N[NDIM], n[NDIM], pos[NDIM];
|
---|
949 | int destpos = 0;
|
---|
950 | double fac[NDIM], posd[NDIM];
|
---|
951 | float posf[NDIM+1];
|
---|
952 | struct Lattice *Lat = &P->Lat;
|
---|
953 | struct LatticeLevel *Lev0 = &Lat->Lev[0];
|
---|
954 | fftw_real *srcpos = P->Files.PosR;
|
---|
955 | //fftw_real *srcdens = Lev0->Dens->DensityArray[ActualDensity]; //[TotalDensity]; trick to display single density
|
---|
956 | float *dest = (float *)P->Files.work;
|
---|
957 | int Nx = Lev0->Plan0.plan->N[0];
|
---|
958 | int i;
|
---|
959 | double min, max;
|
---|
960 |
|
---|
961 | N[0] = Lev0->Plan0.plan->local_nx;
|
---|
962 | N[1] = Lev0->Plan0.plan->N[1];
|
---|
963 | N[2] = Lev0->Plan0.plan->N[2];
|
---|
964 |
|
---|
965 | max = min = srcdens[0];
|
---|
966 | for (i=1;i<P->R.Lev0->Dens->LocalSizeR;i++) {
|
---|
967 | if (srcdens[i] < min) min = srcdens[i];
|
---|
968 | if (srcdens[i] > max) max = srcdens[i];
|
---|
969 | }
|
---|
970 | if (P->Call.out[PsiOut]) fprintf(stderr,"(%i)OutputOutVisDensity: min %e\tmax %e\n",P->Par.me, min, max);
|
---|
971 |
|
---|
972 | // go through all nodes
|
---|
973 | for (n[0]=0; n[0] < N[0]; n[0]++) {
|
---|
974 | pos[0] = (n[0] == N[0] ? 0 : n[0]);
|
---|
975 | for (n[1]=0; n[1] <= N[1]; n[1]++) {
|
---|
976 | pos[1] = (n[1] == N[1] ? 0 : n[1]);
|
---|
977 | for (n[2]=0; n[2] <= N[2]; n[2]++) {
|
---|
978 | pos[2] = (n[2] == N[2] ? 0 : n[2]);
|
---|
979 | // depending on RiemannTensor use, fill FileData::work
|
---|
980 | switch (Lat->RT.ActualUse) {
|
---|
981 | case inactive:
|
---|
982 | case standby:
|
---|
983 | if ((srcdens[pos[2]+N[2]*(pos[1]+N[1]*pos[0])]) > 0.)
|
---|
984 | dest[destpos] = log(1.0+(srcdens[pos[2]+N[2]*(pos[1]+N[1]*pos[0])]));
|
---|
985 | else
|
---|
986 | dest[destpos] = 0.;
|
---|
987 | destpos++;
|
---|
988 | break;
|
---|
989 | case active:
|
---|
990 | posf[0] = srcpos[0+NDIM*(pos[2]+N[2]*(pos[1]+N[1]*pos[0]))];
|
---|
991 | posf[1] = srcpos[1+NDIM*(pos[2]+N[2]*(pos[1]+N[1]*pos[0]))];
|
---|
992 | posf[2] = srcpos[2+NDIM*(pos[2]+N[2]*(pos[1]+N[1]*pos[0]))];
|
---|
993 | fac[0] = ((n[0]+N[0]*myPE)/(double)Nx);
|
---|
994 | fac[1] = (n[1]/(double)N[1]);
|
---|
995 | fac[2] = (n[2]/(double)N[2]);
|
---|
996 | RMat33Vec3(posd, Lat->RealBasis, fac);
|
---|
997 | posf[0] += posd[0];
|
---|
998 | posf[1] += posd[1];
|
---|
999 | posf[2] += posd[2];
|
---|
1000 | if ((srcdens[pos[2]+N[2]*(pos[1]+N[1]*pos[0])]) > 0.)
|
---|
1001 | posf[3] = log(1.0+(srcdens[pos[2]+N[2]*(pos[1]+N[1]*pos[0])]));
|
---|
1002 | else
|
---|
1003 | posf[3] = 0.;
|
---|
1004 | dest[destpos+0] = posf[0];
|
---|
1005 | dest[destpos+1] = posf[1];
|
---|
1006 | dest[destpos+2] = posf[2];
|
---|
1007 | dest[destpos+3] = posf[3];
|
---|
1008 | destpos += 4;
|
---|
1009 | break;
|
---|
1010 | }
|
---|
1011 | }
|
---|
1012 | }
|
---|
1013 | }
|
---|
1014 | if (myPE) MPI_Send(dest, destpos, MPI_FLOAT, 0, OutputDensTag, P->Par.comm_ST_Psi);
|
---|
1015 | }
|
---|
1016 |
|
---|
1017 | /** Combines prepared electronic Psis density and output to file.
|
---|
1018 | * If we are process 0, open file suffixdensdat (only when RiemannTensor is used) and suffixdenspos, receive
|
---|
1019 | * FileData::work logarithmic coefficients sent by the other processes in OutputOutVisDensity(), go through all
|
---|
1020 | * nodes and save the coefficient to file - again depending on RiemannTensor use - followed by FileData::PosTemp
|
---|
1021 | * (for y and z nodes), close file(s).
|
---|
1022 | * \param *P Problem at hand
|
---|
1023 | * \param me this ranks process in the Psi communcator ParallelSimulationData::me_comm_ST_Psi
|
---|
1024 | * \param Maxme number of processes in this Psi communcator ParallelSimulationData::Max_me_comm_ST_Psi
|
---|
1025 | */
|
---|
1026 | static void CombineOutVisDensity(struct Problem *P, const int me, const int Maxme)
|
---|
1027 | {
|
---|
1028 | int i,n[NDIM], N[NDIM];
|
---|
1029 | float posf[NDIM+1];
|
---|
1030 | float *source = (float *)P->Files.work;
|
---|
1031 | double posd[NDIM], fac[NDIM];
|
---|
1032 | float *Temp = (float *)P->Files.PosTemp;
|
---|
1033 | struct Lattice *Lat = &P->Lat;
|
---|
1034 | struct LatticeLevel *Lev0 = &Lat->Lev[0];
|
---|
1035 | float step = P->Files.OutVisStep;
|
---|
1036 | int No=0, destpos;
|
---|
1037 | FILE *DensityData, *DensityPos;
|
---|
1038 | int Nx = Lev0->Plan0.plan->N[0]+1;
|
---|
1039 | MPI_Status status;
|
---|
1040 | if (me) return; // if we are process 0!
|
---|
1041 | N[0] = Lev0->Plan0.plan->local_nx;
|
---|
1042 | N[1] = Lev0->Plan0.plan->N[1]+1;
|
---|
1043 | N[2] = Lev0->Plan0.plan->N[2]+1;
|
---|
1044 | // Open respective file depending on RiemannTensor use
|
---|
1045 | switch (Lat->RT.ActualUse) {
|
---|
1046 | case active:
|
---|
1047 | OpenFileNo(P, &DensityPos, suffixdenspos, (int)step, "wb",P->Call.out[ReadOut]);
|
---|
1048 | case inactive:
|
---|
1049 | case standby:
|
---|
1050 | OpenFileNo(P, &DensityData, suffixdensdat, (int)step, "wb",P->Call.out[ReadOut]);
|
---|
1051 | break;
|
---|
1052 | }
|
---|
1053 | // for all processes in the communicator
|
---|
1054 | for (i=0; i< Maxme; i++) {
|
---|
1055 | if (i) { // if process != 0, receive from this process
|
---|
1056 | /* MPI_Probe( i, OutputDensTag, P->Par.comm_ST_Psi, &status );*/
|
---|
1057 | switch (Lat->RT.ActualUse) {
|
---|
1058 | case inactive:
|
---|
1059 | case standby:
|
---|
1060 | MPI_Recv( source, N[0]*N[1]*N[2], MPI_FLOAT, i, OutputDensTag, P->Par.comm_ST_Psi, &status );
|
---|
1061 | break;
|
---|
1062 | case active:
|
---|
1063 | MPI_Recv( source, N[0]*N[1]*N[2]*4, MPI_FLOAT, i, OutputDensTag, P->Par.comm_ST_Psi, &status );
|
---|
1064 | break;
|
---|
1065 | }
|
---|
1066 | }
|
---|
1067 | destpos = 0;
|
---|
1068 | // go through all nodes and save the coefficient to file DensityData, depending on RiemannTensor
|
---|
1069 | for (n[0]=0; n[0] < N[0]; n[0]++) {
|
---|
1070 | for (n[1]=0; n[1] < N[1]; n[1]++) {
|
---|
1071 | for (n[2]=0; n[2] < N[2]; n[2]++) {
|
---|
1072 | switch (Lat->RT.ActualUse) {
|
---|
1073 | case inactive:
|
---|
1074 | case standby:
|
---|
1075 | posf[3] = source[destpos];
|
---|
1076 | destpos++;
|
---|
1077 | (void)fwrite(&posf[3], sizeof(float), (size_t)(1), DensityData);
|
---|
1078 | No++;
|
---|
1079 | if (i==0 && n[0] == 0)
|
---|
1080 | Temp[(n[2]+N[2]*(n[1]+N[1]*n[0]))] = posf[3];
|
---|
1081 | break;
|
---|
1082 | case active:
|
---|
1083 | posf[0] = source[destpos+0];
|
---|
1084 | posf[1] = source[destpos+1];
|
---|
1085 | posf[2] = source[destpos+2];
|
---|
1086 | posf[3] = source[destpos+3];
|
---|
1087 | destpos += 4;
|
---|
1088 | (void)fwrite(posf, sizeof(float), (size_t)(NDIM), DensityPos);
|
---|
1089 | (void)fwrite(&posf[3], sizeof(float), (size_t)(1), DensityData);
|
---|
1090 | No++;
|
---|
1091 | if (i==0 && n[0] == 0) {
|
---|
1092 | fac[0] = ((n[0]+N[0]*i)/(double)(Nx-1));
|
---|
1093 | fac[1] = (n[1]/(double)(N[1]-1));
|
---|
1094 | fac[2] = (n[2]/(double)(N[2]-1));
|
---|
1095 | RMat33Vec3(posd, Lat->RealBasis, fac);
|
---|
1096 | posf[0] -= posd[0];
|
---|
1097 | posf[1] -= posd[1];
|
---|
1098 | posf[2] -= posd[2];
|
---|
1099 | fac[0] = ((Nx-1)/(double)(Nx-1));
|
---|
1100 | fac[1] = (n[1]/(double)(N[1]-1));
|
---|
1101 | fac[2] = (n[2]/(double)(N[2]-1));
|
---|
1102 | RMat33Vec3(posd, Lat->RealBasis, fac);
|
---|
1103 | posf[0] += posd[0];
|
---|
1104 | posf[1] += posd[1];
|
---|
1105 | posf[2] += posd[2];
|
---|
1106 | Temp[0+(NDIM+1)*(n[2]+N[2]*(n[1]+N[1]*n[0]))] = posf[0];
|
---|
1107 | Temp[1+(NDIM+1)*(n[2]+N[2]*(n[1]+N[1]*n[0]))] = posf[1];
|
---|
1108 | Temp[2+(NDIM+1)*(n[2]+N[2]*(n[1]+N[1]*n[0]))] = posf[2];
|
---|
1109 | Temp[3+(NDIM+1)*(n[2]+N[2]*(n[1]+N[1]*n[0]))] = posf[3];
|
---|
1110 | }
|
---|
1111 | break;
|
---|
1112 | }
|
---|
1113 | }
|
---|
1114 | }
|
---|
1115 | }
|
---|
1116 | }
|
---|
1117 | n[0] = N[0];
|
---|
1118 | for (n[1]=0; n[1] < N[1]; n[1]++) {
|
---|
1119 | for (n[2]=0; n[2] < N[2]; n[2]++) {
|
---|
1120 | switch (Lat->RT.ActualUse) {
|
---|
1121 | case inactive:
|
---|
1122 | case standby:
|
---|
1123 | (void)fwrite(&Temp[n[2]+N[2]*(n[1])], sizeof(float), (size_t)(1), DensityData);
|
---|
1124 | No++;
|
---|
1125 | break;
|
---|
1126 | case active:
|
---|
1127 | (void)fwrite(&Temp[(NDIM+1)*(n[2]+N[2]*(n[1]))], sizeof(float), (size_t)(NDIM), DensityPos);
|
---|
1128 | (void)fwrite(&Temp[(NDIM+1)*(n[2]+N[2]*(n[1]))+3], sizeof(float), (size_t)(1), DensityData);
|
---|
1129 | No++;
|
---|
1130 | break;
|
---|
1131 | }
|
---|
1132 | }
|
---|
1133 | }
|
---|
1134 | if (No != Nx*N[1]*N[2]) Error(SomeError,"CombineOutVisDensity: No != points");
|
---|
1135 | switch (Lat->RT.ActualUse) {
|
---|
1136 | case active:
|
---|
1137 | fclose(DensityPos);
|
---|
1138 | case inactive:
|
---|
1139 | case standby:
|
---|
1140 | fclose(DensityData);
|
---|
1141 | break;
|
---|
1142 | }
|
---|
1143 | }
|
---|
1144 |
|
---|
1145 | /** Main output electronic Psis density for OpenDX.
|
---|
1146 | * If FileData::MeOutVis is set, calls OutputOutVisDensity() followed by CombineOutVisDensity().
|
---|
1147 | * Beforehand CalculateOutVisDensityPos() is called if RiemannTensor is used.
|
---|
1148 | * \param *P Problem at hand
|
---|
1149 | * \param *src_dens Pointer to DensityArray which is to be displayed
|
---|
1150 | */
|
---|
1151 | static void OutVisDensity(struct Problem *P, fftw_real *src_dens)
|
---|
1152 | {
|
---|
1153 | if (!P->Files.MeOutVis) return;
|
---|
1154 | if (P->Lat.RT.ActualUse == active) CalculateOutVisDensityPos(&P->Lat, &P->Files/*, P->Lat.FactorDensityR, P->Lat.FactorDensityC*/);
|
---|
1155 | OutputOutVisDensity(P, P->Par.me_comm_ST_Psi, src_dens);
|
---|
1156 | /* Achtung hier: P->Files.work (RT->TempC, Dens->DensityCArray[TempDensity]) fuer myPE == 0 nicht veraendern !!! */
|
---|
1157 | CombineOutVisDensity(P, P->Par.me_comm_ST_Psi, P->Par.Max_me_comm_ST_Psi);
|
---|
1158 | }
|
---|
1159 |
|
---|
1160 | /** Opening and Initializing of output measurement files.
|
---|
1161 | * If this is process 0, opens and writes top line of FileData::ForcesFile, FileData::EnergyFile.
|
---|
1162 | * and sets FileData::MeOutVis and FileData::MeOutMes (if output desired) to 1, otherwise 0.
|
---|
1163 | * \param *P Problem at hand
|
---|
1164 | */
|
---|
1165 | void InitOutputFiles(struct Problem *P)
|
---|
1166 | {
|
---|
1167 | struct FileData *F = &P->Files;
|
---|
1168 | F->ForcesFile = NULL;
|
---|
1169 | F->EnergyFile = NULL;
|
---|
1170 | F->HamiltonianFile = NULL;
|
---|
1171 | F->MinimisationFile = NULL;
|
---|
1172 | F->SpreadFile = NULL;
|
---|
1173 | F->ReciSpreadFile = NULL;
|
---|
1174 | F->TemperatureFile = NULL;
|
---|
1175 | // process 0 ?
|
---|
1176 | F->MeOutVis = ((P->Par.my_color_comm_ST == 0 && P->Par.my_color_comm_ST_Psi == 0 && F->DoOutVis) ? 1 : 0);
|
---|
1177 | F->MeOutCurr = ((P->Par.my_color_comm_ST == 0 && P->Par.my_color_comm_ST_Psi == 0 && F->DoOutCurr) ? 1 : 0);
|
---|
1178 | F->MeOutMes = ((P->Par.me == 0 && F->DoOutMes) ? 1 : 0);
|
---|
1179 | OpenFile(P, &F->HamiltonianFile, suffixhamiltonianall, "w",P->Call.out[ReadOut]);
|
---|
1180 |
|
---|
1181 | if (!F->MeOutMes) return;
|
---|
1182 | OpenFile(P, &F->ForcesFile, suffixforcesall, "w",P->Call.out[ReadOut]);
|
---|
1183 | if (F->ForcesFile == NULL) fprintf(stderr,"Error opening ForcesFile\n");
|
---|
1184 | // write header of forces file
|
---|
1185 | fprintf(F->ForcesFile, "%s\t%s\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t%s\t%s\t%s\t\t%s\t\t%s\n",
|
---|
1186 | "Type", "No",
|
---|
1187 | "Pos0", "Pos1", "Pos2",
|
---|
1188 | "Total0", "Total1", "Total2",
|
---|
1189 | "Local0", "Local1", "Local2",
|
---|
1190 | "NLocal0", "NLocal1", "NLocal2",
|
---|
1191 | "Magnetic0", "Magnetic1", "Magnetic2",
|
---|
1192 | "Ewald0", "Ewald1", "Ewald2");
|
---|
1193 | OpenFile(P, &F->EnergyFile, suffixenergyall, "w",P->Call.out[ReadOut]);
|
---|
1194 | if (F->EnergyFile == NULL) fprintf(stderr,"Error opening EnergyFile\n");
|
---|
1195 | // write header of energy file
|
---|
1196 | if (P->R.DoUnOccupied) {
|
---|
1197 | fprintf(F->EnergyFile, "%s\t\t%s\t\t%s\t%s\t\t%s\t%s\t\t%s\t%s\t%s\t\t%s\t\t%s\t%s\t\t%s\t\t%s\t\t%s\n",
|
---|
1198 | "Time",
|
---|
1199 | "Total",
|
---|
1200 | "Total+Gap",
|
---|
1201 | "Kinetic", "NonLocal",
|
---|
1202 | "GapPsi",
|
---|
1203 | "Correlation", "Exchange",
|
---|
1204 | "Pseudo", "Hartree",
|
---|
1205 | "GapDensity",
|
---|
1206 | "-Gauss",
|
---|
1207 | "Ewald",
|
---|
1208 | "IonKin",
|
---|
1209 | "ETotal");
|
---|
1210 | } else {
|
---|
1211 | fprintf(F->EnergyFile, "%s\t\t%s\t\t%s\t\t%s\t%s\t%s\t%s\t\t%s\t\t%s\t\t%s\t\t%s\t\t%s\n",
|
---|
1212 | "Time",
|
---|
1213 | "Total",
|
---|
1214 | "Kinetic", "NonLocal",
|
---|
1215 | "Correlation", "Exchange",
|
---|
1216 | "Pseudo", "Hartree",
|
---|
1217 | "-Gauss",
|
---|
1218 | "Ewald",
|
---|
1219 | "IonKin",
|
---|
1220 | "ETotal");
|
---|
1221 | }
|
---|
1222 | OpenFile(P, &F->MinimisationFile, suffixminall, "w",P->Call.out[ReadOut]);
|
---|
1223 | if (F->MinimisationFile == NULL) fprintf(stderr,"Error opening MinimsationFile\n");
|
---|
1224 | // write header of minimsation file
|
---|
1225 | fprintf(F->MinimisationFile, "%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\n","Step", "Psi", "PsiStep", "TE", "ATE", "delta", "dEdt0", "ddEdtt0");
|
---|
1226 |
|
---|
1227 | OpenFile(P, &F->TemperatureFile, suffixtempall, "w",P->Call.out[ReadOut]);
|
---|
1228 | if (F->TemperatureFile == NULL) fprintf(stderr,"Error opening TemperatureFile\n");
|
---|
1229 | // write header of minimsation file
|
---|
1230 | }
|
---|
1231 |
|
---|
1232 | /** Closes all output measurement files.
|
---|
1233 | * Closes FileData::ForcesFile and FileData::EnergyFile
|
---|
1234 | * \param *P Problem at hand
|
---|
1235 | */
|
---|
1236 | void CloseOutputFiles(struct Problem *P)
|
---|
1237 | {
|
---|
1238 | struct FileData *F = &P->Files;
|
---|
1239 | if (!F->MeOutMes) return;
|
---|
1240 | if (F->ForcesFile != NULL) fclose(F->ForcesFile); // only they've been opened (thus not pointing to NULL)
|
---|
1241 | if (F->EnergyFile != NULL) fclose(F->EnergyFile);
|
---|
1242 | if (F->HamiltonianFile != NULL) fclose(F->HamiltonianFile);
|
---|
1243 | if (F->MinimisationFile != NULL) fclose(F->MinimisationFile);
|
---|
1244 | if (F->SpreadFile != NULL) fclose(F->SpreadFile);
|
---|
1245 | if (F->ReciSpreadFile != NULL) fclose(F->ReciSpreadFile);
|
---|
1246 | if (F->TemperatureFile != NULL) fclose(F->TemperatureFile);
|
---|
1247 | }
|
---|
1248 |
|
---|
1249 | /** Initialization of Problem::FileData structure for visual output.
|
---|
1250 | * If this is process 0 (and OutVis desired), allocate memory for FileData::PosTemp, FileData::work,
|
---|
1251 | * set the entries of FileData all to their corresponding values from RiemannTensor,
|
---|
1252 | * FileData::OutVisStep to zero.
|
---|
1253 | * \param *P Problem at hand
|
---|
1254 | */
|
---|
1255 | void InitOutVisArray(struct Problem *P)
|
---|
1256 | {
|
---|
1257 | struct FileData *F = &P->Files;
|
---|
1258 | struct Lattice *Lat = &P->Lat;
|
---|
1259 | struct RiemannTensor *RT = &Lat->RT;
|
---|
1260 | struct LatticeLevel *Lev0 = &Lat->Lev[0];
|
---|
1261 | F->OutputPosType = NULL;
|
---|
1262 | if (!F->MeOutVis) return;
|
---|
1263 | F->OutVisStep = 0;
|
---|
1264 | F->PosTemp = (fftw_complex *)
|
---|
1265 | Malloc(sizeof(float)*(Lev0->Plan0.plan->N[1]+1)*(Lev0->Plan0.plan->N[2]+1)*
|
---|
1266 | ((Lat->RT.Use != UseRT ? 0 : NDIM)+1), "InitOutVisArray: TempC");
|
---|
1267 | F->work = (fftw_complex *)Lev0->Dens->DensityCArray[TempDensity];
|
---|
1268 | if (Lat->RT.Use != UseRT) return;
|
---|
1269 | F->TotalSize = RT->TotalSize[RTAIRT]/NDIM;
|
---|
1270 | F->LocalSizeR = RT->LocalSizeR[RTAiGcg];
|
---|
1271 | F->LocalSizeC = RT->LocalSizeC[RTAiGcg];
|
---|
1272 | F->MaxNUp = RT->MaxNUp[RTPFRto0];
|
---|
1273 | F->PosC = RT->DensityC[RTAiGcg];
|
---|
1274 | F->PosR = (fftw_real *)F->PosC;
|
---|
1275 | F->work = RT->TempC;
|
---|
1276 | F->PosFactor = RT->PosFactor[RTPFRto0];
|
---|
1277 | }
|
---|
1278 |
|
---|
1279 | static const char suffixionfor[] = ".ions.force"; //!< Suffix for ion forces file
|
---|
1280 | static const char suffixionZ[] = ".ions.datZ"; //!< Suffix for ion datZ file
|
---|
1281 | static const char suffixionpos[] = ".ions.pos"; //!< Suffix for ion positions file
|
---|
1282 | static const char suffixiondx[] = ".ions.dx"; //!< Suffix for ions dx file
|
---|
1283 | static const char suffixiondoc[] = ".ions.doc"; //!< Suffix for ions doc file
|
---|
1284 | static const char suffixsrciondoc[] = ".srcion.doc"; //!< Suffix for state ions doc file
|
---|
1285 | static const char suffixsrciondat[] = ".srcion.data"; //!< Suffix for state ions data file
|
---|
1286 |
|
---|
1287 | /** Output current ions state.
|
---|
1288 | * If this is process0, open file suffixsrciondoc for writing, output Ions::Max_Types and
|
---|
1289 | * Ions::Max_IonsOfType of each type - each in a new line - closes it.
|
---|
1290 | * Then opens suffixsrciondat for binary writing, outputs Lattice:RealBasis vectors and
|
---|
1291 | * position IonType::R and speed IonType::U, closes it.
|
---|
1292 | * \param *P Problem at hand
|
---|
1293 | * \note This is the ionic counterpart to the elecontric OutputSrcPsiDensity(), storing a so far made
|
---|
1294 | * calculation to file.
|
---|
1295 | */
|
---|
1296 | static void OutSrcIons(struct Problem *P)
|
---|
1297 | {
|
---|
1298 | struct Ions *I = &P->Ion;
|
---|
1299 | double *U, *pos;
|
---|
1300 | double data[2*NDIM];
|
---|
1301 | int is,ia,i;
|
---|
1302 | FILE *SrcIonDoc, *SrcIonData;
|
---|
1303 |
|
---|
1304 | if (!(P->Par.me == 0)) return;
|
---|
1305 |
|
---|
1306 | // output of ion types and numbers per type
|
---|
1307 | OpenFile(P, &SrcIonDoc, suffixsrciondoc, "w",P->Call.out[ReadOut]);
|
---|
1308 | fprintf(SrcIonDoc,"%i\n", I->Max_Types);
|
---|
1309 | for (is=0; is < I->Max_Types; is++)
|
---|
1310 | fprintf(SrcIonDoc,"%i\n", I->I[is].Max_IonsOfType);
|
---|
1311 | fclose(SrcIonDoc);
|
---|
1312 |
|
---|
1313 | OpenFile(P, &SrcIonData, suffixsrciondat, "wb",P->Call.out[ReadOut]);
|
---|
1314 | (void)fwrite(P->Lat.RealBasis, sizeof(double), (size_t)(NDIM_NDIM), SrcIonData);
|
---|
1315 | for (is=0; is < I->Max_Types; is++) {
|
---|
1316 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
---|
1317 | U = &I->I[is].U[NDIM*ia];
|
---|
1318 | pos = &I->I[is].R[NDIM*ia];
|
---|
1319 | for (i=0;i<NDIM;i++) {
|
---|
1320 | data[i] = pos[i];
|
---|
1321 | data[i+NDIM] = U[i];
|
---|
1322 | }
|
---|
1323 | (void)fwrite(&data, sizeof(double), (size_t)(2*NDIM), SrcIonData);
|
---|
1324 | }
|
---|
1325 | }
|
---|
1326 | fclose(SrcIonData);
|
---|
1327 | }
|
---|
1328 |
|
---|
1329 | /** Read ions state from a file.
|
---|
1330 | * Reads the suffixsrciondoc file and checks it against the current state in Ions regarding
|
---|
1331 | * IonType::MaxTypes and IonType::Max_IonsOfType, closes it.
|
---|
1332 | * Afterwards, opens suffixsrciondat for binary reading, retrieves the basis checking it against
|
---|
1333 | * Problem::Lattice::RealBasis. If ok, reads position IonType::R and speed IonType::U, closes it.
|
---|
1334 | * \param *P Problem at hand
|
---|
1335 | * \return 1 - successful, 0 - failure
|
---|
1336 | * \note This is the ionic counterpart to the electronic ReadSrcPsiDensity(), see also OutSrcIons().
|
---|
1337 | */
|
---|
1338 | int ReadSrcIons(struct Problem *P)
|
---|
1339 | {
|
---|
1340 | struct Ions *I = &P->Ion;
|
---|
1341 | double *U, *pos;
|
---|
1342 | double data[2*NDIM];
|
---|
1343 | int is,ia,i;
|
---|
1344 | int Max_Types;
|
---|
1345 | int *Max_IonsOfType = NULL;
|
---|
1346 | double RealBasis[NDIM_NDIM];
|
---|
1347 | FILE *SrcIonDoc, *SrcIonData;
|
---|
1348 | // read the doc file and check
|
---|
1349 | if (OpenFile(P, &SrcIonDoc, suffixsrciondoc, "r",P->Call.out[ReadOut])) {
|
---|
1350 | if (fscanf(SrcIonDoc,"%i", &Max_Types) != 1)
|
---|
1351 | //Error(SomeError, "ReadSrcIons: read error");
|
---|
1352 | return 0;
|
---|
1353 | if (Max_Types != I->Max_Types)
|
---|
1354 | //Error(SomeError, "ReadSrcIons: srcion file does not fit to system, MaxTypes");
|
---|
1355 | return 0;
|
---|
1356 | Max_IonsOfType = (int *) Malloc(Max_Types*sizeof(int), "ReadSrcIons: Max_IonsOfType");
|
---|
1357 | for (is=0; is < Max_Types; is++) {
|
---|
1358 | if (fscanf(SrcIonDoc,"%i", &Max_IonsOfType[is]) != 1)
|
---|
1359 | //Error(SomeError, "ReadSrcIons: read error");
|
---|
1360 | return 0;
|
---|
1361 | if (Max_IonsOfType[is] != I->I[is].Max_IonsOfType)
|
---|
1362 | //Error(SomeError, "ReadSrcIons: srcion file does not fit to system, Max_IonsOfType");
|
---|
1363 | return 0;
|
---|
1364 | }
|
---|
1365 | fclose(SrcIonDoc);
|
---|
1366 | // read basis, then positions and speeds of ions
|
---|
1367 | if (OpenFile(P, &SrcIonData, suffixsrciondat, "rb",P->Call.out[ReadOut])) {
|
---|
1368 | if (fread(RealBasis, sizeof(double), (size_t)(NDIM_NDIM), SrcIonData) != NDIM_NDIM)
|
---|
1369 | //Error(SomeError, "ReadSrcIons: read error");
|
---|
1370 | return 0;
|
---|
1371 | for (i=0; i < NDIM_NDIM; i++)
|
---|
1372 | if (RealBasis[i] != P->Lat.RealBasis[i])
|
---|
1373 | //Error(SomeError, "ReadSrcIons: srcion file does not fit to system, RealBasis");
|
---|
1374 | return 0;
|
---|
1375 | for (is=0; is < I->Max_Types; is++) {
|
---|
1376 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
---|
1377 | if (fread(&data, sizeof(double), (size_t)(2*NDIM), SrcIonData) != 2*NDIM)
|
---|
1378 | //Error(SomeError, "ReadSrcIons: read error");
|
---|
1379 | return 0;
|
---|
1380 | U = &I->I[is].U[NDIM*ia];
|
---|
1381 | pos = &I->I[is].R[NDIM*ia];
|
---|
1382 | for (i=0;i<NDIM;i++) {
|
---|
1383 | pos[i] = data[i];
|
---|
1384 | U[i] = data[i+NDIM];
|
---|
1385 | }
|
---|
1386 | }
|
---|
1387 | }
|
---|
1388 | fclose(SrcIonData);
|
---|
1389 | }
|
---|
1390 | if (Max_IonsOfType != NULL) Free(Max_IonsOfType, "Max_IonsOfType");
|
---|
1391 | }
|
---|
1392 | return 1;
|
---|
1393 | }
|
---|
1394 |
|
---|
1395 | /** Output of ion doc, dx, forces and positions file for OpenDX.
|
---|
1396 | * If this is process 0,
|
---|
1397 | * open, fill and close IonDoc file suffixiondoc,
|
---|
1398 | * open, fill for each FileData::OutVisStep and close IonDX file suffixiondx
|
---|
1399 | * for every
|
---|
1400 | * open suffixionfor, suffixionpos (and suffixionZ in case of only FileData::OutVisStep), fill
|
---|
1401 | * them with the ion forces IonType::FIon and positions IonType::R of each type and each ion per type,
|
---|
1402 | * close them all.
|
---|
1403 | * \param *P Problem at hand
|
---|
1404 | */
|
---|
1405 | static void OutVisIons(struct Problem *P)
|
---|
1406 | {
|
---|
1407 | struct Ions *I = &P->Ion;
|
---|
1408 | struct FileData *F = &P->Files;
|
---|
1409 | int i,is,ia;
|
---|
1410 | double *fion, *pos;
|
---|
1411 | float data[6]; // holds temporarily twice NDIM values as write buffer
|
---|
1412 | int Z;
|
---|
1413 | char *datnamef, *datnameZ, *posname;
|
---|
1414 | FILE *IonsDataF, *IonsDataZ, *IonsPos, *IonsDoc, *IonsDx;
|
---|
1415 | if (!P->Files.MeOutVis && P->Par.me == 0) return;
|
---|
1416 | // generate file names
|
---|
1417 | datnamef = (char*)
|
---|
1418 | malloc(strlen(P->Files.mainname)+strlen(suffixionfor) + 1);
|
---|
1419 | sprintf(datnamef, "%s%s", P->Files.mainname, suffixionfor);
|
---|
1420 | datnameZ = (char*)
|
---|
1421 | malloc(strlen(P->Files.mainname)+strlen(suffixionZ) + 1);
|
---|
1422 | sprintf(datnameZ, "%s%s", P->Files.mainname, suffixionZ);
|
---|
1423 | posname = (char*)
|
---|
1424 | malloc(strlen(P->Files.mainname)+strlen(suffixionpos) + 1);
|
---|
1425 | sprintf(posname, "%s%s", P->Files.mainname, suffixionpos);
|
---|
1426 | // open, fill and close doc file
|
---|
1427 | if (OpenFile(P, &IonsDoc, suffixiondoc, "w",P->Call.out[ReadOut])) {
|
---|
1428 | fprintf(IonsDoc,"IonsPos file = %s.####\n", posname);
|
---|
1429 | fprintf(IonsDoc,"IonsForce file = %s.####\n", datnamef);
|
---|
1430 | fprintf(IonsDoc,"format = ieee float (Bytes %lu) %s = Force(3)\n",(unsigned long) sizeof(float),msb);
|
---|
1431 | fprintf(IonsDoc,"IonsZ file = %s.####\n", datnameZ);
|
---|
1432 | fprintf(IonsDoc,"format = int (Bytes %lu) %s = Z(1)\n",(unsigned long) sizeof(int),msb);
|
---|
1433 | fprintf(IonsDoc,"points = %i\n", I->Max_TotalIons);
|
---|
1434 | fprintf(IonsDoc,"majority = row\n");
|
---|
1435 | fprintf(IonsDoc,"TimeSeries = %i\n",F->OutVisStep+1);
|
---|
1436 | fclose(IonsDoc);
|
---|
1437 | }
|
---|
1438 | // open dx file and fill it with each output step, close it
|
---|
1439 | if (OpenFile(P, &IonsDx, suffixiondx, "w",P->Call.out[ReadOut])) {
|
---|
1440 | for (i=0; i < F->OutVisStep+1; i++) {
|
---|
1441 | if (i==0) {
|
---|
1442 | /* fprintf(IonsDx,"object \"ioncon\" class array type int rank 1 shape 2 items 0 data follows\nattribute \"element type\" string \"lines\"\nattribute \"ref\" string \"positions\"\n\n"); */
|
---|
1443 | fprintf(IonsDx,"object \"iondatZ\" class array type int rank 0 items %i %s binary\n",I->Max_TotalIons,msb);
|
---|
1444 | fprintf(IonsDx,"data file %s,0\n",datnameZ);
|
---|
1445 | fprintf(IonsDx,"attribute \"dep\" string \"positions\"\n\n");
|
---|
1446 | }
|
---|
1447 |
|
---|
1448 | fprintf(IonsDx,"object \"ionpos.%04i\" class array type float rank 1 shape 3 items %i %s binary\n",i,I->Max_TotalIons,msb);
|
---|
1449 | fprintf(IonsDx,"data file %s.%04i,0\n\n",posname,i);
|
---|
1450 |
|
---|
1451 | fprintf(IonsDx,"object \"iondatF.%04i\" class array type float rank 1 shape 3 items %i %s binary\n",i,I->Max_TotalIons,msb);
|
---|
1452 | fprintf(IonsDx,"data file %s.%04i,0\n",datnamef,i);
|
---|
1453 | fprintf(IonsDx,"attribute \"dep\" string \"positions\"\n\n");
|
---|
1454 |
|
---|
1455 | fprintf(IonsDx,"object \"ionobjF.%04i\" class field\n",i);
|
---|
1456 | fprintf(IonsDx,"component \"positions\" \"ionpos.%04i\"\n",i);
|
---|
1457 | fprintf(IonsDx,"component \"data\" \"iondatF.%04i\"\n",i);
|
---|
1458 | /*fprintf(IonsDx,"component \"connections\" \"ioncon\"\n\n");*/
|
---|
1459 |
|
---|
1460 | fprintf(IonsDx,"object \"ionobjZ.%04i\" class field\n",i);
|
---|
1461 | fprintf(IonsDx,"component \"positions\" \"ionpos.%04i\"\n",i);
|
---|
1462 | fprintf(IonsDx,"component \"data\" \"iondatZ\"\n");
|
---|
1463 | /* fprintf(IonsDx,"component \"connections\" \"ioncon\"\n\n");*/
|
---|
1464 | }
|
---|
1465 | fprintf(IonsDx,"\nobject \"ionseriesF\" class series\n");
|
---|
1466 | for (i=0; i < F->OutVisStep+1; i++)
|
---|
1467 | fprintf(IonsDx,"member %i \"ionobjF.%04i\" position %f\n",i,i,(float)i);
|
---|
1468 | fprintf(IonsDx,"\nobject \"ionseriesZ\" class series\n");
|
---|
1469 | for (i=0; i < F->OutVisStep+1; i++)
|
---|
1470 | fprintf(IonsDx,"member %i \"ionobjZ.%04i\" position %f\n",i,i,(float)i);
|
---|
1471 | fprintf(IonsDx,"end\n");
|
---|
1472 | fclose(IonsDx);
|
---|
1473 | }
|
---|
1474 | Free(datnamef, "OutVisIons: datnamef");
|
---|
1475 | Free(datnameZ, "OutVisIons: datnameZ");
|
---|
1476 | Free(posname, "OutVisIons: posname");
|
---|
1477 | // open IonForces, IonZ and IonPosition file, write forces respectively positions for each ion of each type, close them
|
---|
1478 | if (OpenFileNo(P, &IonsDataF, suffixionfor, F->OutVisStep, "wb",P->Call.out[ReadOut])) {
|
---|
1479 | if (F->OutVisStep == 0)
|
---|
1480 | OpenFile(P, &IonsDataZ, suffixionZ, "wb",P->Call.out[ReadOut]);
|
---|
1481 | if (OpenFileNo(P, &IonsPos, suffixionpos, F->OutVisStep, "wb",P->Call.out[ReadOut])) {
|
---|
1482 | for (is=0; is < I->Max_Types; is++) {
|
---|
1483 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
---|
1484 | fion = &I->I[is].FIon[NDIM*ia];
|
---|
1485 | pos = &I->I[is].R[NDIM*ia];
|
---|
1486 | for (i=0;i<3;i++) {
|
---|
1487 | data[i+3] = fion[i];
|
---|
1488 | data[i] = pos[i];
|
---|
1489 | }
|
---|
1490 | Z = I->I[is].Z;
|
---|
1491 | if (fwrite(&data[0],sizeof(float), (size_t)(3),IonsPos) != 3)
|
---|
1492 | Error(FileOpenParams, "Error writing ion positions!");
|
---|
1493 | if (F->OutVisStep == 0) (void)fwrite(&Z,sizeof(int), (size_t)(1),IonsDataZ);
|
---|
1494 | if (fwrite(&data[3],sizeof(float), (size_t)(3),IonsDataF) != 3)
|
---|
1495 | Error(FileOpenParams, "Error writing ionic forces!");
|
---|
1496 | }
|
---|
1497 | }
|
---|
1498 | fclose(IonsPos);
|
---|
1499 | }
|
---|
1500 | if (F->OutVisStep == 0)
|
---|
1501 | fclose(IonsDataZ);
|
---|
1502 | fclose(IonsDataF);
|
---|
1503 | }
|
---|
1504 | }
|
---|
1505 |
|
---|
1506 | /** Output electronic density and ion state files with so far made calculations.
|
---|
1507 | * If CallOptions::WriteSrcFiles is set, OutputSrcPsiDensity() and OutSrcIons() are called.
|
---|
1508 | * \param *P Problem at hand
|
---|
1509 | * \param type which PsiTypeTag should be put to file
|
---|
1510 | */
|
---|
1511 | void OutputVisSrcFiles(struct Problem *P, enum PsiTypeTag type)
|
---|
1512 | {
|
---|
1513 | if (P->Call.WriteSrcFiles) {
|
---|
1514 | if(P->Call.out[NormalOut]) fprintf(stderr,"(%i) Writing %s srcpsi to disk\n", P->Par.me, P->R.MinimisationName[type]);
|
---|
1515 | OutputSrcPsiDensity(P, type);
|
---|
1516 | if(P->Call.out[NormalOut]) fprintf(stderr,"(%i) Writing srcion to disk\n", P->Par.me);
|
---|
1517 | OutSrcIons(P);
|
---|
1518 | }
|
---|
1519 | // if (!P->Files.MeOutVis) return;
|
---|
1520 | // P->Files.OutputPosType =
|
---|
1521 | // Realloc(P->Files.OutputPosType,sizeof(enum ModeType)*(P->Files.OutVisStep+1),"OutputVis");
|
---|
1522 | // P->Files.OutputPosType[P->Files.OutVisStep] = P->Lat.RT.ActualUse;
|
---|
1523 | // CreateDensityOutputGeneral(P, P->Par.me_comm_ST_Psi);
|
---|
1524 | // OutVisDensity(P);
|
---|
1525 | // OutVisIons(P);
|
---|
1526 | // if(P->Call.out[NormalOut]) fprintf(stderr,"(%i) Written OutVisStep %i to disk\n", P->Par.me, P->Files.OutVisStep);
|
---|
1527 | // /* P->Files.OutVisStep++; Genau ebend nicht hochzaehlen - wird immer ueberschrieben */
|
---|
1528 | }
|
---|
1529 |
|
---|
1530 | /** Main output total electronic density and ion data for OpenDX.
|
---|
1531 | * Calls subsequently preparing CreateDensityOutputGeneral(), then output of electronic
|
---|
1532 | * densities OutVisDensity() and ion data OutVisIons(), increasing finally FileData::OutVisStep.
|
---|
1533 | * \param *P Problem at hand
|
---|
1534 | * \param *srcdens Pointer to DensityArray which is to be displayed
|
---|
1535 | * \note Output is made only RunStruct::OutVisStep steps and if FileData::MeOutVis is set.
|
---|
1536 | */
|
---|
1537 | void OutputVis(struct Problem *P, fftw_real *srcdens)
|
---|
1538 | {
|
---|
1539 | if (!P->Files.MeOutVis) return;
|
---|
1540 | P->Files.OutputPosType = (enum ModeType *) Realloc(P->Files.OutputPosType,sizeof(enum ModeType)*(P->Files.OutVisStep+1),"OutputVis");
|
---|
1541 | P->Files.OutputPosType[P->Files.OutVisStep] = P->Lat.RT.ActualUse;
|
---|
1542 |
|
---|
1543 | CreateDensityOutputGeneral(P, P->Par.me_comm_ST_Psi);
|
---|
1544 | OutVisDensity(P, srcdens);
|
---|
1545 | OutVisIons(P);
|
---|
1546 | if(P->Call.out[MinOut]) fprintf(stderr,"(%i) Written OutVisStep %i to disk\n", P->Par.me, P->Files.OutVisStep);
|
---|
1547 | P->Files.OutVisStep++;
|
---|
1548 | }
|
---|
1549 |
|
---|
1550 | /** Output of each current density component for OpenDX.
|
---|
1551 | * \param *P Problem at hand
|
---|
1552 | * \note Output is made only RunStruct::OutVisStep steps and if FileData::MeOutVis is set.
|
---|
1553 | */
|
---|
1554 | void OutputCurrentDensity(struct Problem *P)
|
---|
1555 | {
|
---|
1556 | int index, i, r;
|
---|
1557 | fftw_real *density = P->R.Lev0->Dens->DensityArray[ActualDensity];
|
---|
1558 | fftw_real *CurrentDensity[NDIM*NDIM];
|
---|
1559 | if (!P->Files.MeOutCurr) return;
|
---|
1560 |
|
---|
1561 | P->Files.OutputPosType = (enum ModeType *) Realloc(P->Files.OutputPosType,sizeof(enum ModeType)*(P->Files.OutVisStep+(1)*NDIM),"OutputVis");
|
---|
1562 | for (i=0;i<(1)*NDIM;i++)
|
---|
1563 | P->Files.OutputPosType[P->Files.OutVisStep+i] = P->Lat.RT.ActualUse;
|
---|
1564 | if(P->Call.out[PsiOut]) fprintf(stderr,"(%i) OutVisStep %i, OutputPosType %p\n",P->Par.me, P->Files.OutVisStep, P->Files.OutputPosType);
|
---|
1565 |
|
---|
1566 | // due to preprocessor values we can't put the following stuff into a loop
|
---|
1567 | CurrentDensity[0] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity0];
|
---|
1568 | CurrentDensity[1] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity1];
|
---|
1569 | CurrentDensity[2] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity2];
|
---|
1570 | CurrentDensity[3] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity3];
|
---|
1571 | CurrentDensity[4] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity4];
|
---|
1572 | CurrentDensity[5] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity5];
|
---|
1573 | CurrentDensity[6] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity6];
|
---|
1574 | CurrentDensity[7] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity7];
|
---|
1575 | CurrentDensity[8] = (fftw_real *) P->R.Lev0->Dens->DensityArray[CurrentDensity8];
|
---|
1576 |
|
---|
1577 | // output current density, not vector component
|
---|
1578 | for (index=0;index<NDIM;index++) {
|
---|
1579 | // evaluate euclidian norm for given B component
|
---|
1580 | SetArrayToDouble0((double *)density,P->R.Lev0->Dens->TotalSize*2); // reset
|
---|
1581 | for (r=0;r<P->R.Lev0->Dens->LocalSizeR;r++) {
|
---|
1582 | for (i=0;i<NDIM;i++)
|
---|
1583 | density[r] += CurrentDensity[i + index*NDIM][r]*CurrentDensity[i + index*NDIM][r];
|
---|
1584 | density[r] = sqrt(density[r]);
|
---|
1585 | }
|
---|
1586 | // output
|
---|
1587 | CreateDensityOutputGeneral(P, P->Par.me_comm_ST_Psi);
|
---|
1588 | OutVisDensity(P, density);
|
---|
1589 | OutVisIons(P);
|
---|
1590 | if(P->Call.out[NormalOut]) fprintf(stderr,"(%i) Written OutVisStep %i to disk\n", P->Par.me, P->Files.OutVisStep);
|
---|
1591 | P->Files.OutVisStep++;
|
---|
1592 | }
|
---|
1593 | }
|
---|
1594 |
|
---|
1595 | /** Output each orbital in a certain step order for OpenDX.
|
---|
1596 | * \param *P Problem at hand
|
---|
1597 | * \param offset from which step do we start
|
---|
1598 | * \param increment by which increment do we advance step-wise
|
---|
1599 | * \param type Only PsiTypeTag orbitals are displayed
|
---|
1600 | */
|
---|
1601 | void OutputVisAllOrbital(struct Problem *P, int offset, int increment, enum PsiTypeTag type) {
|
---|
1602 | struct Lattice *Lat = &P->Lat;
|
---|
1603 | struct Psis *Psi = &Lat->Psi;
|
---|
1604 | struct RunStruct *R = &P->R;
|
---|
1605 | struct LatticeLevel *LevS = R->LevS;
|
---|
1606 | struct LatticeLevel *Lev0 = R->Lev0;
|
---|
1607 | struct Density *Dens0 = Lev0->Dens;
|
---|
1608 | struct OnePsiElement *OnePsiA, *LOnePsiA;
|
---|
1609 | MPI_Status status;
|
---|
1610 | int ElementSize = (sizeof(fftw_complex) / sizeof(double));
|
---|
1611 | fftw_complex *LPsiDatA;
|
---|
1612 | int i, p, RecvSource;
|
---|
1613 | int Num = Psi->NoOfPsis;
|
---|
1614 |
|
---|
1615 | if (!P->Files.MeOutVis) return;
|
---|
1616 |
|
---|
1617 | P->Files.OutputPosType = (enum ModeType *) Realloc(P->Files.OutputPosType,sizeof(enum ModeType)*(P->Files.OutVisStep+(Num)),"OutputVis");
|
---|
1618 |
|
---|
1619 | P->Files.OutVisStep += offset;
|
---|
1620 | P->Files.OutputPosType[P->Files.OutVisStep] = P->Lat.RT.ActualUse;
|
---|
1621 | 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)
|
---|
1622 | OnePsiA = &Psi->AllPsiStatus[i]; // grab the desired OnePsiA
|
---|
1623 | if (OnePsiA->PsiType == type) { // drop if extra one
|
---|
1624 | if (OnePsiA->my_color_comm_ST_Psi == P->Par.my_color_comm_ST_Psi) // local?
|
---|
1625 | LOnePsiA = &Psi->LocalPsiStatus[OnePsiA->MyLocalNo];
|
---|
1626 | else
|
---|
1627 | LOnePsiA = NULL;
|
---|
1628 | if (LOnePsiA == NULL) { // if it's not local ... receive it from respective process into TempPsi
|
---|
1629 | RecvSource = OnePsiA->my_color_comm_ST_Psi;
|
---|
1630 | MPI_Recv( LevS->LPsi->TempPsi, LevS->MaxG*ElementSize, MPI_DOUBLE, RecvSource, WannierTag1, P->Par.comm_ST_PsiT, &status );
|
---|
1631 | LPsiDatA=LevS->LPsi->TempPsi;
|
---|
1632 | } else { // .. otherwise send it to all other processes (Max_me... - 1)
|
---|
1633 | for (p=0;p<P->Par.Max_me_comm_ST_PsiT;p++)
|
---|
1634 | if (p != OnePsiA->my_color_comm_ST_Psi)
|
---|
1635 | MPI_Send( LevS->LPsi->LocalPsi[OnePsiA->MyLocalNo], LevS->MaxG*ElementSize, MPI_DOUBLE, p, WannierTag1, P->Par.comm_ST_PsiT);
|
---|
1636 | LPsiDatA=LevS->LPsi->LocalPsi[OnePsiA->MyLocalNo];
|
---|
1637 | } // LPsiDatA is now set to the coefficients of OnePsi either stored or MPI_Received
|
---|
1638 |
|
---|
1639 | P->Files.OutputPosType[P->Files.OutVisStep] = P->Lat.RT.ActualUse;
|
---|
1640 | // recalculate density for the specific wave function ...
|
---|
1641 | CalculateOneDensityR(Lat, LevS, Dens0, LPsiDatA, Dens0->DensityArray[ActualDensity], R->FactorDensityR, 0);
|
---|
1642 | // ... and output (wherein ActualDensity is used instead of TotalDensity)
|
---|
1643 | //OutputVis(P);
|
---|
1644 | CreateDensityOutputGeneral(P, P->Par.me_comm_ST_Psi);
|
---|
1645 | OutVisDensity(P, Dens0->DensityArray[ActualDensity]);
|
---|
1646 | OutVisIons(P);
|
---|
1647 | if(P->Call.out[NormalOut]) fprintf(stderr,"(%i) Written OutVisStep %i to disk\n", P->Par.me, P->Files.OutVisStep);
|
---|
1648 | P->Files.OutVisStep+=increment;
|
---|
1649 | P->Files.OutputPosType[P->Files.OutVisStep] = P->Lat.RT.ActualUse;
|
---|
1650 | }
|
---|
1651 | }
|
---|
1652 | }
|
---|
1653 |
|
---|
1654 | /** Read source files containing stored calculations.
|
---|
1655 | * Calls ReadSrcPsiDensity() and ReadSrcIons().
|
---|
1656 | * \param *P Problem at hand
|
---|
1657 | */
|
---|
1658 | void ReadSrcFiles(struct Problem *P)
|
---|
1659 | {
|
---|
1660 | if (P->Call.out[NormalOut]) fprintf(stderr, "(%i)ReadSrcPsiDensity\n", P->Par.me);
|
---|
1661 | ReadSrcPsiDensity(P, Occupied, 0, P->R.LevSNo);
|
---|
1662 | ReadSrcPsiDensity(P, UnOccupied, 0, P->R.LevSNo);
|
---|
1663 | if (P->Call.out[NormalOut]) fprintf(stderr, "(%i)ReadSrcIons\n", P->Par.me);
|
---|
1664 | ReadSrcIons(P);
|
---|
1665 | }
|
---|
1666 |
|
---|
1667 | /** Plots a cut plane of the real density of one wave function.
|
---|
1668 | * \param *P Problem at hand
|
---|
1669 | * \param index index of axis (vector orthogonal to plane)
|
---|
1670 | * \param node node specifying where to cut at the given axis
|
---|
1671 | * \param wavenr global number of wave function
|
---|
1672 | * \param *density density array to plot
|
---|
1673 | * \sa PlotVectorPlane() - very similar
|
---|
1674 | */
|
---|
1675 | void PlotSrcPlane(struct Problem *P, int index, double node, int wavenr, fftw_real *density)
|
---|
1676 | {
|
---|
1677 | struct RunStruct *R = &P->R;
|
---|
1678 | struct Lattice *Lat = &P->Lat;
|
---|
1679 | struct LatticeLevel *Lev0 = R->Lev0;
|
---|
1680 | const int myPE = P->Par.me_comm_ST_Psi;
|
---|
1681 | char filename[255], spin[12];
|
---|
1682 | char *suchpointer;
|
---|
1683 | FILE *PlotFile = NULL;
|
---|
1684 | time_t seconds;
|
---|
1685 |
|
---|
1686 |
|
---|
1687 | // open file
|
---|
1688 | time(&seconds); // get current time
|
---|
1689 |
|
---|
1690 | switch (Lat->Psi.PsiST) {
|
---|
1691 | case SpinDouble:
|
---|
1692 | sprintf(&filename[0], ".psi%i_cut%i.csv", wavenr, index);
|
---|
1693 | strncat(spin,"SpinDouble",12);
|
---|
1694 | break;
|
---|
1695 | case SpinUp:
|
---|
1696 | sprintf(&filename[0], ".psi%i_cut%i_up.csv", wavenr, index);
|
---|
1697 | strncat(spin,"SpinUp",12);
|
---|
1698 | break;
|
---|
1699 | case SpinDown:
|
---|
1700 | sprintf(&filename[0], ".psi%i_cut%i_down.csv", wavenr, index);
|
---|
1701 | strncat(spin,"SpinDown",12);
|
---|
1702 | break;
|
---|
1703 | }
|
---|
1704 |
|
---|
1705 | if (!myPE) { // only process 0 writes to file
|
---|
1706 | OpenFile(P, &PlotFile, filename, "w", P->Call.out[ReadOut]);
|
---|
1707 | strcpy(filename, ctime(&seconds));
|
---|
1708 | suchpointer = strchr(filename, '\n');
|
---|
1709 | if (suchpointer != NULL)
|
---|
1710 | *suchpointer = '\0';
|
---|
1711 | if (PlotFile != NULL) {
|
---|
1712 | fprintf(PlotFile,"# Psi %i, type %s (real density) plot of plane perpendicular to direction e_%i at node %lg, seed %i, config %s, run on %s, #cpus %i", wavenr, spin, index, node, R->Seed, P->Files.default_path, filename, P->Par.Max_me_comm_ST_Psi);
|
---|
1713 | fprintf(PlotFile,"\n");
|
---|
1714 | } else { Error(SomeError, "PlotSrcPlane: Opening Plot File"); }
|
---|
1715 | }
|
---|
1716 |
|
---|
1717 | // plot density
|
---|
1718 | PlotRealDensity(P, Lev0, PlotFile, index, node, density, density);
|
---|
1719 |
|
---|
1720 | if (PlotFile != NULL) {
|
---|
1721 | // close file
|
---|
1722 | fclose(PlotFile);
|
---|
1723 | }
|
---|
1724 | }
|
---|
1725 |
|
---|
1726 | /** plots a cut plane of a given 3d real density.
|
---|
1727 | * \param *P Problem at hand, contains pointer to Lattice structure
|
---|
1728 | * \param *Lev LatticeLevel of the real density
|
---|
1729 | * \param PlotFile file pointer (already open and valid)
|
---|
1730 | * \param index index of lattice axis
|
---|
1731 | * \param n_orth position on lattice axis where to cut
|
---|
1732 | * \param *density1 first real density array
|
---|
1733 | * \param *density2 second real density array (point to \a *density1 if not needed)
|
---|
1734 | */
|
---|
1735 | void PlotRealDensity(struct Problem *P, struct LatticeLevel *Lev, FILE *PlotFile, int index, double n_orth, fftw_real *density1, fftw_real *density2)
|
---|
1736 | {
|
---|
1737 | struct Lattice *Lat = &P->Lat;
|
---|
1738 | int n[NDIM], n0;
|
---|
1739 | int N[NDIM];
|
---|
1740 | N[0] = Lev->Plan0.plan->N[0];
|
---|
1741 | N[1] = Lev->Plan0.plan->N[1];
|
---|
1742 | N[2] = Lev->Plan0.plan->N[2];
|
---|
1743 | const int N0 = Lev->Plan0.plan->local_nx;
|
---|
1744 | const int myPE = P->Par.me_comm_ST_Psi;
|
---|
1745 | double fac[NDIM], x[NDIM];
|
---|
1746 | int i0, i = 0;
|
---|
1747 | int PE, zahl;
|
---|
1748 | double *buffer;
|
---|
1749 | MPI_Status status;
|
---|
1750 | int sizes[P->Par.Max_me_comm_ST_Psi], c0, c1;
|
---|
1751 | double nodes[NDIM], node[NDIM];
|
---|
1752 |
|
---|
1753 | for(i=0;i<NDIM;i++) {
|
---|
1754 | nodes[i] = (i == index) ? n_orth : 0.;
|
---|
1755 | node[i] = 0.;
|
---|
1756 | }
|
---|
1757 | RMat33Vec3(node, Lat->ReciBasis, nodes); // transform cartesian coordinates into cell coordinates [0,1]^3
|
---|
1758 | for(i=0;i<NDIM;i++) // now N^3 within node range of discrete grid
|
---|
1759 | node[i] = (int)(node[i]*N[i]/(2.*PI));
|
---|
1760 | fprintf(stderr,"(%i) n_orth %lg, index %i converted to plane offset vector (%lg, %lg, %lg).\n", P->Par.me, n_orth, index, node[0], node[1], node[2]);
|
---|
1761 |
|
---|
1762 | switch (index) {
|
---|
1763 | case 0:
|
---|
1764 | zahl = 4*N[1]*N[2];
|
---|
1765 | break;
|
---|
1766 | case 1:
|
---|
1767 | zahl = 4*N0*N[2];
|
---|
1768 | break;
|
---|
1769 | case 2:
|
---|
1770 | zahl = 4*N0*N[1];
|
---|
1771 | break;
|
---|
1772 | }
|
---|
1773 | fprintf(stderr,"(%i) buffer size %i\n", P->Par.me, zahl);
|
---|
1774 | buffer = Malloc(sizeof(double)*zahl,"PlotRealDensity: buffer");
|
---|
1775 |
|
---|
1776 | c0 = cross(index,0);
|
---|
1777 | c1 = cross(index,1);
|
---|
1778 | // then for every point on the grid in real space ...
|
---|
1779 | i=0;
|
---|
1780 | for (n0=0;n0<N0;n0++) // only local points on x axis
|
---|
1781 | for (n[1]=0;n[1]<N[1];n[1]++)
|
---|
1782 | for (n[2]=0;n[2]<N[2];n[2]++) {
|
---|
1783 | n[0]=n0 + N0*myPE; // global relative coordinate: due to partitoning of x-axis in PEPGamma>1 case
|
---|
1784 | if (n[index] == (int)node[index]) { // only on the correct plane orthogonal to desired axis and at desired node ...
|
---|
1785 | fac[0] = (double)n[0]/(double)N[0];
|
---|
1786 | fac[1] = (double)n[1]/(double)N[1];
|
---|
1787 | fac[2] = (double)n[2]/(double)N[2];
|
---|
1788 | RMat33Vec3(x, Lat->RealBasis, fac); // relative coordinate times basis matrix gives absolute ones
|
---|
1789 | i0 = n[2]+N[2]*(n[1]+N[1]*n0); // index to local density array
|
---|
1790 |
|
---|
1791 | buffer[i++] = x[c0]; // fill buffer
|
---|
1792 | buffer[i++] = x[c1];
|
---|
1793 | buffer[i++] = density1[i0];
|
---|
1794 | buffer[i++] = density2[i0];
|
---|
1795 | if (i > zahl) Error(SomeError, "PlotRealDensity: buffer too small!");
|
---|
1796 | }
|
---|
1797 | }
|
---|
1798 | // exchange sizes of each buffer
|
---|
1799 | MPI_Allgather(&i, 1, MPI_INT, sizes, 1, MPI_INT, P->Par.comm_ST_Psi);
|
---|
1800 | if (myPE == 0) {
|
---|
1801 | for (PE=0; PE < P->Par.Max_me_comm_ST_Psi; PE++) {
|
---|
1802 | if (PE != 0) {
|
---|
1803 | // receive them
|
---|
1804 | if (MPI_Recv(buffer, sizes[PE], MPI_DOUBLE, PE, PlotRealDensityTag, P->Par.comm_ST_Psi, &status) != MPI_SUCCESS)
|
---|
1805 | Error(SomeError, "PlotRealDensity: MPI_Recv failure!");
|
---|
1806 | MPI_Get_count(&status, MPI_DOUBLE, &zahl);
|
---|
1807 | if (zahl != sizes[PE])
|
---|
1808 | Error(SomeError, "PlotRealDensity: received unexpected amount of elements!");
|
---|
1809 | }
|
---|
1810 | //write them: local one (still in buffer) and received ones
|
---|
1811 | for (i0 = 0; i0 < sizes[PE];) {
|
---|
1812 | fprintf(PlotFile,"%e", buffer[(i0)++]);
|
---|
1813 | if ((i0 % 4) == 0) {
|
---|
1814 | fprintf(PlotFile,"\n");
|
---|
1815 | } else {
|
---|
1816 | fprintf(PlotFile,"\t");
|
---|
1817 | }
|
---|
1818 | }
|
---|
1819 | }
|
---|
1820 | } else { // send them
|
---|
1821 | if (MPI_Send(buffer, i, MPI_DOUBLE, 0, PlotRealDensityTag, P->Par.comm_ST_Psi) != MPI_SUCCESS)
|
---|
1822 | Error(SomeError, "PlotRealDensity: MPI_Send failure!");
|
---|
1823 | }
|
---|
1824 | Free(buffer, "PlotRealDensity: buffer");
|
---|
1825 | }
|
---|