ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpWriter.cpp
Revision: 419
Committed: Thu Mar 27 15:07:29 2003 UTC (21 years, 3 months ago) by gezelter
File size: 12329 byte(s)
Log Message:
fixes for local->global atom numbering in MPI

File Contents

# User Rev Content
1 mmeineke 377 #include <cstring>
2     #include <iostream>
3     #include <fstream>
4    
5     #ifdef IS_MPI
6     #include <mpi.h>
7 gezelter 416 #include <mpi++.h>
8 mmeineke 377 #include "mpiSimulation.hpp"
9     #define TAKE_THIS_TAG 0
10     #endif //is_mpi
11    
12     #include "ReadWrite.hpp"
13     #include "simError.h"
14    
15     DumpWriter::DumpWriter( SimInfo* the_entry_plug ){
16    
17     entry_plug = the_entry_plug;
18    
19     #ifdef IS_MPI
20     if(worldRank == 0 ){
21     #endif // is_mpi
22    
23    
24    
25     strcpy( outName, entry_plug->sampleName );
26    
27     outFile.open(outName, ios::out | ios::trunc );
28    
29     if( !outFile ){
30    
31     sprintf( painCave.errMsg,
32     "Could not open \"%s\" for dump output.\n",
33     outName);
34     painCave.isFatal = 1;
35     simError();
36     }
37    
38     //outFile.setf( ios::scientific );
39    
40     #ifdef IS_MPI
41     }
42    
43     sprintf( checkPointMsg,
44     "Sucessfully opened output file for dumping.\n");
45     MPIcheckPoint();
46     #endif // is_mpi
47     }
48    
49     DumpWriter::~DumpWriter( ){
50    
51     #ifdef IS_MPI
52     if(worldRank == 0 ){
53     #endif // is_mpi
54    
55     outFile.close();
56    
57     #ifdef IS_MPI
58     }
59     #endif // is_mpi
60     }
61    
62     void DumpWriter::writeDump( double currentTime ){
63    
64     const int BUFFERSIZE = 2000;
65     char tempBuffer[BUFFERSIZE];
66     char writeLine[BUFFERSIZE];
67    
68 gezelter 419 int i, j, which_node, done, game_over, which_atom, local_index;
69 mmeineke 377 double q[4];
70     DirectionalAtom* dAtom;
71     int nAtoms = entry_plug->n_atoms;
72     Atom** atoms = entry_plug->atoms;
73    
74    
75     #ifndef IS_MPI
76    
77     outFile << nAtoms << "\n";
78    
79     outFile << currentTime << "\t"
80     << entry_plug->box_x << "\t"
81     << entry_plug->box_y << "\t"
82     << entry_plug->box_z << "\n";
83    
84     for( i=0; i<nAtoms; i++ ){
85    
86    
87     sprintf( tempBuffer,
88     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
89     atoms[i]->getType(),
90     atoms[i]->getX(),
91     atoms[i]->getY(),
92     atoms[i]->getZ(),
93     atoms[i]->get_vx(),
94     atoms[i]->get_vy(),
95     atoms[i]->get_vz());
96     strcpy( writeLine, tempBuffer );
97    
98     if( atoms[i]->isDirectional() ){
99    
100     dAtom = (DirectionalAtom *)atoms[i];
101     dAtom->getQ( q );
102    
103     sprintf( tempBuffer,
104     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
105     q[0],
106     q[1],
107     q[2],
108     q[3],
109     dAtom->getJx(),
110     dAtom->getJy(),
111     dAtom->getJz());
112     strcat( writeLine, tempBuffer );
113     }
114     else
115     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
116    
117     outFile << writeLine;
118     }
119     outFile.flush();
120    
121     #else // is_mpi
122 gezelter 416
123     MPI::Status istatus;
124     int *AtomToProcMap = mpiSim->getAtomToProcMap();
125 gezelter 415
126 mmeineke 377 // write out header and node 0's coordinates
127 gezelter 415
128 mmeineke 377 if( worldRank == 0 ){
129     outFile << mpiSim->getTotAtoms() << "\n";
130 gezelter 415
131 mmeineke 377 outFile << currentTime << "\t"
132     << entry_plug->box_x << "\t"
133     << entry_plug->box_y << "\t"
134     << entry_plug->box_z << "\n";
135 gezelter 415
136 gezelter 416 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
137 gezelter 417 // Get the Node number which has this atom;
138 mmeineke 377
139 gezelter 415 which_node = AtomToProcMap[i];
140    
141 gezelter 416 if (which_node == mpiSim->getMyNode()) {
142 gezelter 415
143     sprintf( tempBuffer,
144     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
145     atoms[i]->getType(),
146     atoms[i]->getX(),
147     atoms[i]->getY(),
148     atoms[i]->getZ(),
149     atoms[i]->get_vx(),
150     atoms[i]->get_vy(),
151     atoms[i]->get_vz());
152     strcpy( writeLine, tempBuffer );
153 mmeineke 377
154 gezelter 415 if( atoms[i]->isDirectional() ){
155 mmeineke 377
156 gezelter 415 dAtom = (DirectionalAtom *)atoms[i];
157     dAtom->getQ( q );
158 mmeineke 377
159 gezelter 415 sprintf( tempBuffer,
160     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
161     q[0],
162     q[1],
163     q[2],
164     q[3],
165     dAtom->getJx(),
166     dAtom->getJy(),
167     dAtom->getJz());
168     strcat( writeLine, tempBuffer );
169     }
170     else
171     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
172    
173     } else {
174    
175     MPI::COMM_WORLD.Send(&i, 1, MPI_INT, which_node, TAKE_THIS_TAG);
176 gezelter 416 MPI::COMM_WORLD.Recv(writeLine, BUFFERSIZE, MPI_CHAR, which_node,
177 gezelter 415 TAKE_THIS_TAG, istatus);
178 mmeineke 377 }
179 gezelter 415
180 mmeineke 377 outFile << writeLine;
181     }
182    
183 gezelter 415 // kill everyone off:
184 gezelter 416 game_over = -1;
185     for (j = 0; j < mpiSim->getNumberProcessors(); j++) {
186     MPI::COMM_WORLD.Send(&game_over, 1, MPI_INT, j, TAKE_THIS_TAG);
187 mmeineke 377 }
188    
189 gezelter 415 } else {
190    
191     done = 0;
192     while (!done) {
193 gezelter 416 MPI::COMM_WORLD.Recv(&which_atom, 1, MPI_INT, 0,
194 gezelter 419 TAKE_THIS_TAG, istatus);
195    
196 gezelter 416 if (which_atom == -1) {
197 gezelter 415 done=1;
198     continue;
199     } else {
200 gezelter 419 local_index=-1;
201     for (j=0; j < mpiSim->getMyNlocal(); j++) {
202     if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
203     }
204     if (local_index != -1) {
205     //format the line
206     sprintf( tempBuffer,
207     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
208     atoms[local_index]->getType(),
209     atoms[local_index]->getX(),
210     atoms[local_index]->getY(),
211     atoms[local_index]->getZ(),
212     atoms[local_index]->get_vx(),
213     atoms[local_index]->get_vy(),
214     atoms[local_index]->get_vz()); // check here.
215     strcpy( writeLine, tempBuffer );
216 mmeineke 377
217 gezelter 419 if( atoms[local_index]->isDirectional() ){
218    
219     dAtom = (DirectionalAtom *)atoms[local_index];
220     dAtom->getQ( q );
221    
222     sprintf( tempBuffer,
223     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
224     q[0],
225     q[1],
226     q[2],
227     q[3],
228     dAtom->getJx(),
229     dAtom->getJy(),
230     dAtom->getJz());
231     strcat( writeLine, tempBuffer );
232     }
233     else
234     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
235    
236     MPI::COMM_WORLD.Send(writeLine, BUFFERSIZE, MPI_CHAR, 0,
237     TAKE_THIS_TAG);
238     } else {
239     strcpy( writeLine, "ATOM NOT FOUND ON THIS PROCESSOR");
240     MPI::COMM_WORLD.Send(writeLine, BUFFERSIZE, MPI_CHAR, 0,
241     TAKE_THIS_TAG);
242     }
243 mmeineke 377 }
244     }
245 gezelter 415 }
246     outFile.flush();
247     sprintf( checkPointMsg,
248     "Sucessfully took a dump.\n");
249     MPIcheckPoint();
250 mmeineke 377 #endif // is_mpi
251     }
252    
253     void DumpWriter::writeFinal(){
254 gezelter 416
255 mmeineke 377 char finalName[500];
256     ofstream finalOut;
257 gezelter 416
258     const int BUFFERSIZE = 2000;
259     char tempBuffer[BUFFERSIZE];
260     char writeLine[BUFFERSIZE];
261    
262     double q[4];
263     DirectionalAtom* dAtom;
264     int nAtoms = entry_plug->n_atoms;
265     Atom** atoms = entry_plug->atoms;
266 gezelter 419 int i, j, which_node, done, game_over, which_atom, local_index;
267 mmeineke 377
268 gezelter 416
269 mmeineke 377 #ifdef IS_MPI
270     if(worldRank == 0 ){
271     #endif // is_mpi
272    
273     strcpy( finalName, entry_plug->finalName );
274    
275     finalOut.open( finalName, ios::out | ios::trunc );
276     if( !finalOut ){
277     sprintf( painCave.errMsg,
278     "Could not open \"%s\" for final dump output.\n",
279     finalName );
280     painCave.isFatal = 1;
281     simError();
282     }
283    
284     // finalOut.setf( ios::scientific );
285    
286     #ifdef IS_MPI
287     }
288    
289     sprintf(checkPointMsg,"Opened file for final configuration\n");
290     MPIcheckPoint();
291    
292     #endif //is_mpi
293    
294 gezelter 415
295 mmeineke 377 #ifndef IS_MPI
296    
297     finalOut << nAtoms << "\n";
298    
299     finalOut << entry_plug->box_x << "\t"
300     << entry_plug->box_y << "\t"
301     << entry_plug->box_z << "\n";
302 gezelter 416
303 mmeineke 377 for( i=0; i<nAtoms; i++ ){
304    
305     sprintf( tempBuffer,
306     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
307     atoms[i]->getType(),
308     atoms[i]->getX(),
309     atoms[i]->getY(),
310     atoms[i]->getZ(),
311     atoms[i]->get_vx(),
312     atoms[i]->get_vy(),
313     atoms[i]->get_vz());
314     strcpy( writeLine, tempBuffer );
315    
316     if( atoms[i]->isDirectional() ){
317    
318     dAtom = (DirectionalAtom *)atoms[i];
319     dAtom->getQ( q );
320    
321     sprintf( tempBuffer,
322     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
323     q[0],
324     q[1],
325     q[2],
326     q[3],
327     dAtom->getJx(),
328     dAtom->getJy(),
329     dAtom->getJz());
330     strcat( writeLine, tempBuffer );
331     }
332     else
333     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
334    
335     finalOut << writeLine;
336     }
337     finalOut.flush();
338 gezelter 415 finalOut.close();
339 mmeineke 377
340     #else // is_mpi
341 gezelter 415
342 gezelter 416 MPI::Status istatus;
343     int *AtomToProcMap = mpiSim->getAtomToProcMap();
344    
345 mmeineke 377 // write out header and node 0's coordinates
346 gezelter 415
347 mmeineke 377 if( worldRank == 0 ){
348     finalOut << mpiSim->getTotAtoms() << "\n";
349 gezelter 415
350 mmeineke 377 finalOut << entry_plug->box_x << "\t"
351 gezelter 415 << entry_plug->box_y << "\t"
352     << entry_plug->box_z << "\n";
353 mmeineke 377
354 gezelter 416 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
355 gezelter 415 // Get the Node number which has this molecule:
356 mmeineke 377
357 gezelter 415 which_node = AtomToProcMap[i];
358    
359 gezelter 416 if (which_node == mpiSim->getMyNode()) {
360 gezelter 415
361     sprintf( tempBuffer,
362     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
363     atoms[i]->getType(),
364     atoms[i]->getX(),
365     atoms[i]->getY(),
366     atoms[i]->getZ(),
367     atoms[i]->get_vx(),
368     atoms[i]->get_vy(),
369     atoms[i]->get_vz());
370     strcpy( writeLine, tempBuffer );
371 mmeineke 377
372 gezelter 415 if( atoms[i]->isDirectional() ){
373 mmeineke 377
374 gezelter 415 dAtom = (DirectionalAtom *)atoms[i];
375     dAtom->getQ( q );
376 mmeineke 377
377 gezelter 415 sprintf( tempBuffer,
378     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
379     q[0],
380     q[1],
381     q[2],
382     q[3],
383     dAtom->getJx(),
384     dAtom->getJy(),
385     dAtom->getJz());
386     strcat( writeLine, tempBuffer );
387     }
388     else
389     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
390    
391     } else {
392    
393     MPI::COMM_WORLD.Send(&i, 1, MPI_INT, which_node, TAKE_THIS_TAG);
394 gezelter 416 MPI::COMM_WORLD.Recv(writeLine, BUFFERSIZE, MPI_CHAR, which_node,
395 gezelter 415 TAKE_THIS_TAG, istatus);
396 mmeineke 377 }
397 gezelter 415
398 mmeineke 377 finalOut << writeLine;
399     }
400    
401 gezelter 415 // kill everyone off:
402 gezelter 416 game_over = -1;
403     for (j = 0; j < mpiSim->getNumberProcessors(); j++) {
404     MPI::COMM_WORLD.Send(&game_over, 1, MPI_INT, j, TAKE_THIS_TAG);
405 mmeineke 377 }
406    
407 gezelter 415 } else {
408    
409     done = 0;
410     while (!done) {
411 gezelter 416 MPI::COMM_WORLD.Recv(&which_atom, 1, MPI_INT, 0,
412     TAKE_THIS_TAG, istatus);
413 mmeineke 377
414 gezelter 416 if (which_atom == -1) {
415 gezelter 415 done=1;
416     continue;
417     } else {
418 mmeineke 377
419 gezelter 419 local_index=-1;
420     for (j=0; j < mpiSim->getMyNlocal(); j++) {
421     if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
422     }
423     if (local_index != -1) {
424    
425     //format the line
426     sprintf( tempBuffer,
427     "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
428     atoms[local_index]->getType(),
429     atoms[local_index]->getX(),
430     atoms[local_index]->getY(),
431     atoms[local_index]->getZ(),
432     atoms[local_index]->get_vx(),
433     atoms[local_index]->get_vy(),
434     atoms[local_index]->get_vz()); // check here.
435     strcpy( writeLine, tempBuffer );
436 mmeineke 377
437 gezelter 419 if( atoms[local_index]->isDirectional() ){
438    
439     dAtom = (DirectionalAtom *)atoms[local_index];
440     dAtom->getQ( q );
441 mmeineke 377
442 gezelter 419 sprintf( tempBuffer,
443     "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
444     q[0],
445     q[1],
446     q[2],
447     q[3],
448     dAtom->getJx(),
449     dAtom->getJy(),
450     dAtom->getJz());
451     strcat( writeLine, tempBuffer );
452     }
453     else
454     strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
455    
456     MPI::COMM_WORLD.Send(writeLine, BUFFERSIZE, MPI_CHAR, 0,
457     TAKE_THIS_TAG);
458     } else {
459     strcpy( writeLine, "ATOM NOT FOUND ON THIS PROCESSOR");
460     MPI::COMM_WORLD.Send(writeLine, BUFFERSIZE, MPI_CHAR, 0,
461     TAKE_THIS_TAG);
462     }
463 mmeineke 377 }
464     }
465 gezelter 419 }
466 gezelter 415 finalOut.flush();
467     sprintf( checkPointMsg,
468     "Sucessfully took a dump.\n");
469     MPIcheckPoint();
470 mmeineke 377
471 gezelter 415 if( worldRank == 0 ) finalOut.close();
472 mmeineke 377 #endif // is_mpi
473     }