1 |
/* |
2 |
* Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
* |
4 |
* The University of Notre Dame grants you ("Licensee") a |
5 |
* non-exclusive, royalty free, license to use, modify and |
6 |
* redistribute this software in source and binary code form, provided |
7 |
* that the following conditions are met: |
8 |
* |
9 |
* 1. Acknowledgement of the program authors must be made in any |
10 |
* publication of scientific results based in part on use of the |
11 |
* program. An acceptable form of acknowledgement is citation of |
12 |
* the article in which the program was described (Matthew |
13 |
* A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
* J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
* Parallel Simulation Engine for Molecular Dynamics," |
16 |
* J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
* |
18 |
* 2. Redistributions of source code must retain the above copyright |
19 |
* notice, this list of conditions and the following disclaimer. |
20 |
* |
21 |
* 3. Redistributions in binary form must reproduce the above copyright |
22 |
* notice, this list of conditions and the following disclaimer in the |
23 |
* documentation and/or other materials provided with the |
24 |
* distribution. |
25 |
* |
26 |
* This software is provided "AS IS," without a warranty of any |
27 |
* kind. All express or implied conditions, representations and |
28 |
* warranties, including any implied warranty of merchantability, |
29 |
* fitness for a particular purpose or non-infringement, are hereby |
30 |
* excluded. The University of Notre Dame and its licensors shall not |
31 |
* be liable for any damages suffered by licensee as a result of |
32 |
* using, modifying or distributing the software or its |
33 |
* derivatives. In no event will the University of Notre Dame or its |
34 |
* licensors be liable for any lost revenue, profit or data, or for |
35 |
* direct, indirect, special, consequential, incidental or punitive |
36 |
* damages, however caused and regardless of the theory of liability, |
37 |
* arising out of the use of or inability to use software, even if the |
38 |
* University of Notre Dame has been advised of the possibility of |
39 |
* such damages. |
40 |
*/ |
41 |
|
42 |
#define _LARGEFILE_SOURCE64 |
43 |
#define _FILE_OFFSET_BITS 64 |
44 |
|
45 |
#include <sys/types.h> |
46 |
#include <sys/stat.h> |
47 |
|
48 |
#include <iostream> |
49 |
#include <math.h> |
50 |
|
51 |
#include <stdio.h> |
52 |
#include <stdlib.h> |
53 |
#include <string.h> |
54 |
|
55 |
#include "primitives/Molecule.hpp" |
56 |
#include "utils/MemoryUtils.hpp" |
57 |
#include "utils/StringTokenizer.hpp" |
58 |
#include "io/RestReader.hpp" |
59 |
#include "utils/simError.h" |
60 |
|
61 |
#ifdef IS_MPI |
62 |
#include <mpi.h> |
63 |
#define TAKE_THIS_TAG_CHAR 0 |
64 |
#define TAKE_THIS_TAG_INT 1 |
65 |
#define TAKE_THIS_TAG_DOUBLE 2 |
66 |
#endif // is_mpi |
67 |
|
68 |
namespace oopse { |
69 |
|
70 |
RestReader::RestReader( SimInfo* info ) : info_(info){ |
71 |
|
72 |
idealName = "idealCrystal.in"; |
73 |
|
74 |
isScanned = false; |
75 |
|
76 |
#ifdef IS_MPI |
77 |
if (worldRank == 0) { |
78 |
#endif |
79 |
|
80 |
inIdealFile = fopen(idealName, "r"); |
81 |
if(inIdealFile == NULL){ |
82 |
sprintf(painCave.errMsg, |
83 |
"RestReader: Cannot open file: %s\n", idealName); |
84 |
painCave.isFatal = 1; |
85 |
simError(); |
86 |
} |
87 |
|
88 |
inIdealFileName = idealName; |
89 |
#ifdef IS_MPI |
90 |
} |
91 |
strcpy( checkPointMsg, |
92 |
"File \"idealCrystal.in\" opened successfully for reading." ); |
93 |
MPIcheckPoint(); |
94 |
#endif |
95 |
return; |
96 |
} |
97 |
|
98 |
RestReader :: ~RestReader( ){ |
99 |
#ifdef IS_MPI |
100 |
if (worldRank == 0) { |
101 |
#endif |
102 |
int error; |
103 |
error = fclose( inIdealFile ); |
104 |
|
105 |
if( error ){ |
106 |
sprintf( painCave.errMsg, |
107 |
"Error closing %s\n", inIdealFileName.c_str()); |
108 |
simError(); |
109 |
} |
110 |
|
111 |
MemoryUtils::deletePointers(framePos); |
112 |
|
113 |
#ifdef IS_MPI |
114 |
} |
115 |
strcpy( checkPointMsg, "Restraint file closed successfully." ); |
116 |
MPIcheckPoint(); |
117 |
#endif |
118 |
|
119 |
return; |
120 |
} |
121 |
|
122 |
|
123 |
void RestReader :: readIdealCrystal(){ |
124 |
|
125 |
int i; |
126 |
unsigned int j; |
127 |
|
128 |
#ifdef IS_MPI |
129 |
int done, which_node, which_atom; // loop counter |
130 |
#endif //is_mpi |
131 |
|
132 |
const int BUFFERSIZE = 2000; // size of the read buffer |
133 |
int nTotObjs; // the number of atoms |
134 |
char read_buffer[BUFFERSIZE]; //the line buffer for reading |
135 |
|
136 |
char *eof_test; // ptr to see when we reach the end of the file |
137 |
char *parseErr; |
138 |
|
139 |
std::vector<StuntDouble*> integrableObjects; |
140 |
|
141 |
Molecule* mol; |
142 |
StuntDouble* integrableObject; |
143 |
SimInfo::MoleculeIterator mi; |
144 |
Molecule::IntegrableObjectIterator ii; |
145 |
|
146 |
#ifndef IS_MPI |
147 |
|
148 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inIdealFile); |
149 |
if( eof_test == NULL ){ |
150 |
sprintf( painCave.errMsg, |
151 |
"RestraintReader error: error reading 1st line of \"%s\"\n", |
152 |
inIdealFileName.c_str() ); |
153 |
painCave.isFatal = 1; |
154 |
simError(); |
155 |
} |
156 |
|
157 |
nTotObjs = atoi( read_buffer ); |
158 |
|
159 |
if( nTotObjs != info_->getNGlobalIntegrableObjects() ){ |
160 |
sprintf( painCave.errMsg, |
161 |
"RestraintReader error. %s nIntegrable, %d, " |
162 |
"does not match the meta-data file's nIntegrable, %d.\n", |
163 |
inIdealFileName.c_str(), nTotObjs, |
164 |
info_->getNGlobalIntegrableObjects()); |
165 |
painCave.isFatal = 1; |
166 |
simError(); |
167 |
} |
168 |
|
169 |
// skip over the comment line |
170 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inIdealFile); |
171 |
if(eof_test == NULL){ |
172 |
sprintf( painCave.errMsg, |
173 |
"error in reading commment in %s\n", inIdealFileName.c_str()); |
174 |
painCave.isFatal = 1; |
175 |
simError(); |
176 |
} |
177 |
|
178 |
// parse the ideal crystal lines |
179 |
/* |
180 |
* Note: we assume that there is a one-to-one correspondence between |
181 |
* integrable objects and lines in the idealCrystal.in file. Thermodynamic |
182 |
* integration is only supported for simple rigid bodies. |
183 |
*/ |
184 |
for (mol = info_->beginMolecule(mi); mol != NULL; |
185 |
mol = info_->nextMolecule(mi)) { |
186 |
|
187 |
for (integrableObject = mol->beginIntegrableObject(ii); |
188 |
integrableObject != NULL; |
189 |
integrableObject = mol->nextIntegrableObject(ii)) { |
190 |
|
191 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inIdealFile); |
192 |
if(eof_test == NULL){ |
193 |
sprintf(painCave.errMsg, |
194 |
"RestReader Error: error in reading file %s\n" |
195 |
"natoms = %d; index = %d\n" |
196 |
"error reading the line from the file.\n", |
197 |
inIdealFileName.c_str(), nTotObjs, |
198 |
integrableObject->getGlobalIndex() ); |
199 |
painCave.isFatal = 1; |
200 |
simError(); |
201 |
} |
202 |
|
203 |
parseErr = parseIdealLine( read_buffer, integrableObject); |
204 |
if( parseErr != NULL ){ |
205 |
strcpy( painCave.errMsg, parseErr ); |
206 |
painCave.isFatal = 1; |
207 |
simError(); |
208 |
} |
209 |
} |
210 |
} |
211 |
|
212 |
// MPI Section of code.......... |
213 |
#else //IS_MPI |
214 |
|
215 |
// first thing first, suspend fatalities. |
216 |
painCave.isEventLoop = 1; |
217 |
|
218 |
int masterNode = 0; |
219 |
int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone |
220 |
|
221 |
MPI_Status istatus; |
222 |
int nCurObj; |
223 |
int nitems; |
224 |
int haveError; |
225 |
|
226 |
nTotObjs = info_->getNGlobalIntegrableObjects(); |
227 |
haveError = 0; |
228 |
|
229 |
if (worldRank == masterNode) { |
230 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inIdealFile); |
231 |
if( eof_test == NULL ){ |
232 |
sprintf( painCave.errMsg, |
233 |
"Error reading 1st line of %s \n ",inIdealFileName.c_str()); |
234 |
painCave.isFatal = 1; |
235 |
simError(); |
236 |
} |
237 |
|
238 |
nitems = atoi( read_buffer ); |
239 |
|
240 |
// Check to see that the number of integrable objects in the |
241 |
// intial configuration file is the same as derived from the |
242 |
// meta-data file. |
243 |
if( nTotObjs != nitems){ |
244 |
sprintf( painCave.errMsg, |
245 |
"RestraintReader Error. %s nIntegrable, %d, " |
246 |
"does not match the meta-data file's nIntegrable, %d.\n", |
247 |
inIdealFileName.c_str(), nTotObjs, |
248 |
info_->getNGlobalIntegrableObjects()); |
249 |
painCave.isFatal = 1; |
250 |
simError(); |
251 |
} |
252 |
|
253 |
// skip over the comment line |
254 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inIdealFile); |
255 |
if(eof_test == NULL){ |
256 |
sprintf( painCave.errMsg, |
257 |
"error in reading commment in %s\n", inIdealFileName.c_str()); |
258 |
painCave.isFatal = 1; |
259 |
simError(); |
260 |
} |
261 |
|
262 |
for (i=0 ; i < info_->getNGlobalMolecules(); i++) { |
263 |
int which_node = info_->getMolToProc(i); |
264 |
|
265 |
if(which_node == masterNode){ |
266 |
//molecules belong to master node |
267 |
|
268 |
mol = info_->getMoleculeByGlobalIndex(i); |
269 |
|
270 |
if(mol == NULL) { |
271 |
sprintf(painCave.errMsg, |
272 |
"RestReader Error: Molecule not found on node %d!\n", |
273 |
worldRank); |
274 |
painCave.isFatal = 1; |
275 |
simError(); |
276 |
} |
277 |
|
278 |
for (integrableObject = mol->beginIntegrableObject(ii); |
279 |
integrableObject != NULL; |
280 |
integrableObject = mol->nextIntegrableObject(ii)){ |
281 |
|
282 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inIdealFile); |
283 |
|
284 |
if(eof_test == NULL){ |
285 |
sprintf(painCave.errMsg, |
286 |
"RestReader Error: error in reading file %s\n" |
287 |
"natoms = %d; index = %d\n" |
288 |
"error reading the line from the file.\n", |
289 |
inIdealFileName.c_str(), nTotObjs, i ); |
290 |
painCave.isFatal = 1; |
291 |
simError(); |
292 |
} |
293 |
|
294 |
parseIdealLine(read_buffer, integrableObjects[j]); |
295 |
} |
296 |
} else { |
297 |
//molecule belongs to slave nodes |
298 |
|
299 |
MPI_Recv(&nCurObj, 1, MPI_INT, which_node, |
300 |
TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus); |
301 |
|
302 |
for(j=0; j < nCurObj; j++){ |
303 |
|
304 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inIdealFile); |
305 |
if(eof_test == NULL){ |
306 |
sprintf(painCave.errMsg, |
307 |
"RestReader Error: error in reading file %s\n" |
308 |
"natoms = %d; index = %d\n" |
309 |
"error reading the line from the file.\n", |
310 |
inIdealFileName.c_str(), nTotObjs, i ); |
311 |
painCave.isFatal = 1; |
312 |
simError(); |
313 |
} |
314 |
|
315 |
MPI_Send(read_buffer, BUFFERSIZE, MPI_CHAR, which_node, |
316 |
TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD); |
317 |
} |
318 |
} |
319 |
} |
320 |
} else { |
321 |
//actions taken at slave nodes |
322 |
for (i=0 ; i < info_->getNGlobalMolecules(); i++) { |
323 |
int which_node = info_->getMolToProc(i); |
324 |
|
325 |
if(which_node == worldRank){ |
326 |
//molecule with global index i belongs to this processor |
327 |
|
328 |
mol = info_->getMoleculeByGlobalIndex(i); |
329 |
|
330 |
if(mol == NULL) { |
331 |
sprintf(painCave.errMsg, |
332 |
"RestReader Error: molecule not found on node %d\n", |
333 |
worldRank); |
334 |
painCave.isFatal = 1; |
335 |
simError(); |
336 |
} |
337 |
|
338 |
nCurObj = mol->getNIntegrableObjects(); |
339 |
|
340 |
MPI_Send(&nCurObj, 1, MPI_INT, masterNode, |
341 |
TAKE_THIS_TAG_INT, MPI_COMM_WORLD); |
342 |
|
343 |
for (integrableObject = mol->beginIntegrableObject(ii); |
344 |
integrableObject != NULL; |
345 |
integrableObject = mol->nextIntegrableObject(ii)){ |
346 |
|
347 |
MPI_Recv(read_buffer, BUFFERSIZE, MPI_CHAR, masterNode, |
348 |
TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD, &istatus); |
349 |
|
350 |
parseErr = parseIdealLine(read_buffer, integrableObject); |
351 |
|
352 |
if( parseErr != NULL ){ |
353 |
strcpy( painCave.errMsg, parseErr ); |
354 |
simError(); |
355 |
} |
356 |
} |
357 |
} |
358 |
} |
359 |
} |
360 |
#endif |
361 |
} |
362 |
|
363 |
char* RestReader::parseIdealLine(char* readLine, StuntDouble* sd){ |
364 |
|
365 |
char *foo; // the pointer to the current string token |
366 |
|
367 |
double pos[3]; // position place holders |
368 |
double q[4]; // the quaternions |
369 |
double RfromQ[3][3]; // the rotation matrix |
370 |
double normalize; // to normalize the reference unit vector |
371 |
double uX, uY, uZ; // reference unit vector place holders |
372 |
double uselessToken; |
373 |
StringTokenizer tokenizer(readLine); |
374 |
int nTokens; |
375 |
|
376 |
nTokens = tokenizer.countTokens(); |
377 |
|
378 |
if (nTokens < 14) { |
379 |
sprintf(painCave.errMsg, |
380 |
"RestReader Error: Not enough Tokens.\n"); |
381 |
painCave.isFatal = 1; |
382 |
simError(); |
383 |
} |
384 |
|
385 |
std::string name = tokenizer.nextToken(); |
386 |
|
387 |
if (name != sd->getType()) { |
388 |
|
389 |
sprintf(painCave.errMsg, |
390 |
"RestReader Error: Atom type [%s] in %s does not " |
391 |
"match Atom Type [%s] in .md file.\n", |
392 |
name.c_str(), inIdealFileName.c_str(), |
393 |
sd->getType().c_str()); |
394 |
painCave.isFatal = 1; |
395 |
simError(); |
396 |
} |
397 |
|
398 |
pos[0] = tokenizer.nextTokenAsDouble(); |
399 |
pos[1] = tokenizer.nextTokenAsDouble(); |
400 |
pos[2] = tokenizer.nextTokenAsDouble(); |
401 |
|
402 |
// store the positions in the stuntdouble as generic data doubles |
403 |
DoubleGenericData* refPosX = new DoubleGenericData(); |
404 |
refPosX->setID("refPosX"); |
405 |
refPosX->setData(pos[0]); |
406 |
sd->addProperty(refPosX); |
407 |
|
408 |
DoubleGenericData* refPosY = new DoubleGenericData(); |
409 |
refPosY->setID("refPosY"); |
410 |
refPosY->setData(pos[1]); |
411 |
sd->addProperty(refPosY); |
412 |
|
413 |
DoubleGenericData* refPosZ = new DoubleGenericData(); |
414 |
refPosZ->setID("refPosZ"); |
415 |
refPosZ->setData(pos[2]); |
416 |
sd->addProperty(refPosZ); |
417 |
|
418 |
// we don't need the velocities |
419 |
uselessToken = tokenizer.nextTokenAsDouble(); |
420 |
uselessToken = tokenizer.nextTokenAsDouble(); |
421 |
uselessToken = tokenizer.nextTokenAsDouble(); |
422 |
|
423 |
if (sd->isDirectional()) { |
424 |
|
425 |
q[0] = tokenizer.nextTokenAsDouble(); |
426 |
q[1] = tokenizer.nextTokenAsDouble(); |
427 |
q[2] = tokenizer.nextTokenAsDouble(); |
428 |
q[3] = tokenizer.nextTokenAsDouble(); |
429 |
|
430 |
// now build the rotation matrix and find the unit vectors |
431 |
RfromQ[0][0] = q[0]*q[0] + q[1]*q[1] - q[2]*q[2] - q[3]*q[3]; |
432 |
RfromQ[0][1] = 2*(q[1]*q[2] + q[0]*q[3]); |
433 |
RfromQ[0][2] = 2*(q[1]*q[3] - q[0]*q[2]); |
434 |
RfromQ[1][0] = 2*(q[1]*q[2] - q[0]*q[3]); |
435 |
RfromQ[1][1] = q[0]*q[0] - q[1]*q[1] + q[2]*q[2] - q[3]*q[3]; |
436 |
RfromQ[1][2] = 2*(q[2]*q[3] + q[0]*q[1]); |
437 |
RfromQ[2][0] = 2*(q[1]*q[3] + q[0]*q[2]); |
438 |
RfromQ[2][1] = 2*(q[2]*q[3] - q[0]*q[1]); |
439 |
RfromQ[2][2] = q[0]*q[0] - q[1]*q[1] - q[2]*q[2] + q[3]*q[3]; |
440 |
|
441 |
normalize = sqrt(RfromQ[2][0]*RfromQ[2][0] + RfromQ[2][1]*RfromQ[2][1] |
442 |
+ RfromQ[2][2]*RfromQ[2][2]); |
443 |
uX = RfromQ[2][0]/normalize; |
444 |
uY = RfromQ[2][1]/normalize; |
445 |
uZ = RfromQ[2][2]/normalize; |
446 |
|
447 |
// store reference unit vectors as generic data in the stuntdouble |
448 |
DoubleGenericData* refVectorX = new DoubleGenericData(); |
449 |
refVectorX->setID("refVectorX"); |
450 |
refVectorX->setData(uX); |
451 |
sd->addProperty(refVectorX); |
452 |
|
453 |
DoubleGenericData* refVectorY = new DoubleGenericData(); |
454 |
refVectorY->setID("refVectorY"); |
455 |
refVectorY->setData(uY); |
456 |
sd->addProperty(refVectorY); |
457 |
|
458 |
DoubleGenericData* refVectorZ = new DoubleGenericData(); |
459 |
refVectorZ->setID("refVectorZ"); |
460 |
refVectorZ->setData(uZ); |
461 |
sd->addProperty(refVectorZ); |
462 |
} |
463 |
|
464 |
// we don't need the angular velocities, so let's exit the line |
465 |
return NULL; |
466 |
} |
467 |
|
468 |
void RestReader::readZangle(){ |
469 |
|
470 |
int i; |
471 |
unsigned int j; |
472 |
int isPresent; |
473 |
|
474 |
Molecule* mol; |
475 |
StuntDouble* integrableObject; |
476 |
SimInfo::MoleculeIterator mi; |
477 |
Molecule::IntegrableObjectIterator ii; |
478 |
|
479 |
#ifdef IS_MPI |
480 |
int done, which_node, which_atom; // loop counter |
481 |
int nProc; |
482 |
MPI_Status istatus; |
483 |
#endif //is_mpi |
484 |
|
485 |
const int BUFFERSIZE = 2000; // size of the read buffer |
486 |
int nTotObjs; // the number of atoms |
487 |
char read_buffer[BUFFERSIZE]; //the line buffer for reading |
488 |
|
489 |
char *eof_test; // ptr to see when we reach the end of the file |
490 |
char *parseErr; |
491 |
|
492 |
std::vector<StuntDouble*> vecParticles; |
493 |
std::vector<double> tempZangs; |
494 |
|
495 |
inAngFileName = info_->getRestFileName(); |
496 |
|
497 |
inAngFileName += "0"; |
498 |
|
499 |
// open the omega value file for reading |
500 |
#ifdef IS_MPI |
501 |
if (worldRank == 0) { |
502 |
#endif |
503 |
isPresent = 1; |
504 |
inAngFile = fopen(inAngFileName.c_str(), "r"); |
505 |
if(!inAngFile){ |
506 |
sprintf(painCave.errMsg, |
507 |
"Restraints Warning: %s file is not present\n" |
508 |
"\tAll omega values will be initialized to zero. If the\n" |
509 |
"\tsimulation is starting from the idealCrystal.in reference\n" |
510 |
"\tconfiguration, this is the desired action. If this is not\n" |
511 |
"\tthe case, the energy calculations will be incorrect.\n", |
512 |
inAngFileName.c_str()); |
513 |
painCave.severity = OOPSE_WARNING; |
514 |
painCave.isFatal = 0; |
515 |
simError(); |
516 |
isPresent = 0; |
517 |
} |
518 |
|
519 |
#ifdef IS_MPI |
520 |
// let the other nodes know the status of the file search |
521 |
MPI_Bcast(&isPresent, 1, MPI_INT, 0, MPI_COMM_WORLD); |
522 |
#endif // is_mpi |
523 |
|
524 |
if (!isPresent) { |
525 |
zeroZangle(); |
526 |
return; |
527 |
} |
528 |
|
529 |
#ifdef IS_MPI |
530 |
} |
531 |
|
532 |
// listen to node 0 to see if we should exit this function |
533 |
if (worldRank != 0) { |
534 |
MPI_Bcast(&isPresent, 1, MPI_INT, 0, MPI_COMM_WORLD); |
535 |
if (!isPresent) { |
536 |
zeroZangle(); |
537 |
return; |
538 |
} |
539 |
} |
540 |
|
541 |
strcpy( checkPointMsg, "zAngle file opened successfully for reading." ); |
542 |
MPIcheckPoint(); |
543 |
#endif |
544 |
|
545 |
#ifndef IS_MPI |
546 |
|
547 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inAngFile); |
548 |
if( eof_test == NULL ){ |
549 |
sprintf( painCave.errMsg, |
550 |
"RestraintReader error: error reading 1st line of \"%s\"\n", |
551 |
inAngFileName.c_str() ); |
552 |
painCave.isFatal = 1; |
553 |
simError(); |
554 |
} |
555 |
|
556 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inAngFile); |
557 |
while ( eof_test != NULL ) { |
558 |
// check for and ignore blank lines |
559 |
if ( read_buffer != NULL ) |
560 |
tempZangs.push_back( atof(read_buffer) ); |
561 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inAngFile); |
562 |
} |
563 |
|
564 |
nTotObjs = info_->getNGlobalIntegrableObjects(); |
565 |
|
566 |
if( nTotObjs != tempZangs.size() ){ |
567 |
sprintf( painCave.errMsg, |
568 |
"RestraintReader zAngle reading error. %s nIntegrable, %d, " |
569 |
"does not match the meta-data file's nIntegrable, %d.\n", |
570 |
inAngFileName.c_str(), tempZangs.size(), nTotObjs ); |
571 |
painCave.isFatal = 1; |
572 |
simError(); |
573 |
} |
574 |
|
575 |
// load the zAngles into the integrable objects |
576 |
i = 0; |
577 |
for (mol = info_->beginMolecule(mi); mol != NULL; |
578 |
mol = info_->nextMolecule(mi)) { |
579 |
|
580 |
for (integrableObject = mol->beginIntegrableObject(ii); |
581 |
integrableObject != NULL; |
582 |
integrableObject = mol->nextIntegrableObject(ii)) { |
583 |
|
584 |
integrableObject->setZangle(tempZangs[i]); |
585 |
i++; |
586 |
} |
587 |
} |
588 |
|
589 |
// MPI Section of code.......... |
590 |
#else //IS_MPI |
591 |
|
592 |
// first thing first, suspend fatalities. |
593 |
painCave.isEventLoop = 1; |
594 |
|
595 |
int masterNode = 0; |
596 |
int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone |
597 |
int haveError; |
598 |
int index; |
599 |
|
600 |
int nCurObj; |
601 |
double angleTranfer; |
602 |
|
603 |
nTotObjs = info_->getNGlobalIntegrableObjects(); |
604 |
haveError = 0; |
605 |
|
606 |
if (worldRank == masterNode) { |
607 |
|
608 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inAngFile); |
609 |
if( eof_test == NULL ){ |
610 |
sprintf( painCave.errMsg, |
611 |
"Error reading 1st line of %s \n ",inAngFileName.c_str()); |
612 |
haveError = 1; |
613 |
simError(); |
614 |
} |
615 |
|
616 |
// let node 0 load the temporary angle vector |
617 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inAngFile); |
618 |
while ( eof_test != NULL ) { |
619 |
// check for and ignore blank lines |
620 |
if ( read_buffer != NULL ) |
621 |
tempZangs.push_back( atof(read_buffer) ); |
622 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inAngFile); |
623 |
} |
624 |
|
625 |
// Check to see that the number of integrable objects in the |
626 |
// intial configuration file is the same as derived from the |
627 |
// meta-data file. |
628 |
if( nTotObjs != tempZangs.size() ){ |
629 |
sprintf( painCave.errMsg, |
630 |
"RestraintReader zAngle reading Error. %s nIntegrable, %d, " |
631 |
"does not match the meta-data file's nIntegrable, %d.\n", |
632 |
inAngFileName.c_str(), tempZangs.size(), nTotObjs); |
633 |
haveError= 1; |
634 |
simError(); |
635 |
} |
636 |
|
637 |
// At this point, node 0 has a tempZangs vector completed, and |
638 |
// everyone else has nada |
639 |
index = 0; |
640 |
|
641 |
for (i=0 ; i < info_->getNGlobalMolecules(); i++) { |
642 |
// Get the Node number which has this atom |
643 |
which_node = info_->getMolToProc(i); |
644 |
|
645 |
if (worldRank == masterNode) { |
646 |
mol = info_->getMoleculeByGlobalIndex(i); |
647 |
|
648 |
if(mol == NULL) { |
649 |
strcpy(painCave.errMsg, "Molecule not found on node 0!"); |
650 |
haveError = 1; |
651 |
simError(); |
652 |
} |
653 |
|
654 |
for (integrableObject = mol->beginIntegrableObject(ii); |
655 |
integrableObject != NULL; |
656 |
integrableObject = mol->nextIntegrableObject(ii)){ |
657 |
|
658 |
integrableObject->setZangle(tempZangs[index]); |
659 |
index++; |
660 |
} |
661 |
|
662 |
} else { |
663 |
// I am MASTER OF THE UNIVERSE, but I don't own this molecule |
664 |
|
665 |
MPI_Recv(&nCurObj, 1, MPI_INT, which_node, |
666 |
TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus); |
667 |
|
668 |
for(j=0; j < nCurObj; j++){ |
669 |
angleTransfer = tempZangs[index]; |
670 |
MPI_Send(&angleTransfer, 1, MPI_DOUBLE, which_node, |
671 |
TAKE_THIS_TAG_DOUBLE, MPI_COMM_WORLD); |
672 |
index++; |
673 |
} |
674 |
|
675 |
} |
676 |
} |
677 |
} else { |
678 |
// I am SLAVE TO THE MASTER |
679 |
for (i=0 ; i < info_->getNGlobalMolecules(); i++) { |
680 |
int which_node = info_->getMolToProc(i); |
681 |
|
682 |
if (which_node == worldRank) { |
683 |
|
684 |
// BUT I OWN THIS MOLECULE!!! |
685 |
|
686 |
mol = info_->getMoleculeByGlobalIndex(i); |
687 |
|
688 |
if(mol == NULL) { |
689 |
sprintf(painCave.errMsg, |
690 |
"RestReader Error: molecule not found on node %d\n", |
691 |
worldRank); |
692 |
painCave.isFatal = 1; |
693 |
simError(); |
694 |
} |
695 |
|
696 |
nCurObj = mol->getNIntegrableObjects(); |
697 |
|
698 |
MPI_Send(&nCurObj, 1, MPI_INT, 0, |
699 |
TAKE_THIS_TAG_INT, MPI_COMM_WORLD); |
700 |
|
701 |
for (integrableObject = mol->beginIntegrableObject(ii); |
702 |
integrableObject != NULL; |
703 |
integrableObject = mol->nextIntegrableObject(ii)){ |
704 |
|
705 |
MPI_Recv(&angleTransfer, 1, MPI_DOUBLE, 0, |
706 |
TAKE_THIS_TAG_DOUBLE, MPI_COMM_WORLD, &istatus); |
707 |
|
708 |
integrableObject->setZangle(angleTransfer); |
709 |
} |
710 |
} |
711 |
} |
712 |
} |
713 |
#endif |
714 |
} |
715 |
|
716 |
void RestReader :: zeroZangle(){ |
717 |
|
718 |
int i; |
719 |
unsigned int j; |
720 |
int nTotObjs; // the number of atoms |
721 |
|
722 |
Molecule* mol; |
723 |
StuntDouble* integrableObject; |
724 |
SimInfo::MoleculeIterator mi; |
725 |
Molecule::IntegrableObjectIterator ii; |
726 |
|
727 |
std::vector<StuntDouble*> vecParticles; |
728 |
|
729 |
#ifndef IS_MPI |
730 |
// set all zAngles to 0.0 |
731 |
for (mol = info_->beginMolecule(mi); mol != NULL; |
732 |
mol = info_->nextMolecule(mi)) |
733 |
|
734 |
for (integrableObject = mol->beginIntegrableObject(ii); |
735 |
integrableObject != NULL; |
736 |
integrableObject = mol->nextIntegrableObject(ii)) |
737 |
integrableObject->setZangle( 0.0 ); |
738 |
|
739 |
|
740 |
// MPI Section of code.......... |
741 |
#else //IS_MPI |
742 |
|
743 |
// first thing first, suspend fatalities. |
744 |
painCave.isEventLoop = 1; |
745 |
|
746 |
int masterNode = 0; |
747 |
int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone |
748 |
int haveError; |
749 |
int which_node; |
750 |
|
751 |
MPI_Status istatus; |
752 |
|
753 |
int nCurObj; |
754 |
double angleTranfer; |
755 |
|
756 |
nTotObjs = info_->getNGlobalIntegrableObjects(); |
757 |
haveError = 0; |
758 |
if (worldRank == masterNode) { |
759 |
|
760 |
for (i=0 ; i < info_->getNGlobalMolecules(); i++) { |
761 |
// Get the Node number which has this atom |
762 |
which_node = info_->getMolToProc(i); |
763 |
|
764 |
// let's let node 0 pass out constant values to all the processors |
765 |
if (worldRank == masterNode) { |
766 |
mol = info_->getMoleculeByGlobalIndex(i); |
767 |
|
768 |
if(mol == NULL) { |
769 |
strcpy(painCave.errMsg, "Molecule not found on node 0!"); |
770 |
haveError = 1; |
771 |
simError(); |
772 |
} |
773 |
|
774 |
for (integrableObject = mol->beginIntegrableObject(ii); |
775 |
integrableObject != NULL; |
776 |
integrableObject = mol->nextIntegrableObject(ii)){ |
777 |
|
778 |
integrableObject->setZangle( 0.0 ); |
779 |
|
780 |
} |
781 |
|
782 |
} else { |
783 |
// I am MASTER OF THE UNIVERSE, but I don't own this molecule |
784 |
|
785 |
MPI_Recv(&nCurObj, 1, MPI_INT, which_node, |
786 |
TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus); |
787 |
|
788 |
for(j=0; j < nCurObj; j++){ |
789 |
angleTransfer = 0.0; |
790 |
MPI_Send(&angleTransfer, 1, MPI_DOUBLE, which_node, |
791 |
TAKE_THIS_TAG_DOUBLE, MPI_COMM_WORLD); |
792 |
|
793 |
} |
794 |
} |
795 |
} |
796 |
} else { |
797 |
// I am SLAVE TO THE MASTER |
798 |
for (i=0 ; i < info_->getNGlobalMolecules(); i++) { |
799 |
int which_node = info_->getMolToProc(i); |
800 |
|
801 |
if (which_node == worldRank) { |
802 |
|
803 |
// BUT I OWN THIS MOLECULE!!! |
804 |
mol = info_->getMoleculeByGlobalIndex(i); |
805 |
|
806 |
if(mol == NULL) { |
807 |
sprintf(painCave.errMsg, |
808 |
"RestReader Error: molecule not found on node %d\n", |
809 |
worldRank); |
810 |
painCave.isFatal = 1; |
811 |
simError(); |
812 |
} |
813 |
|
814 |
nCurObj = mol->getNIntegrableObjects(); |
815 |
|
816 |
MPI_Send(&nCurObj, 1, MPI_INT, 0, |
817 |
TAKE_THIS_TAG_INT, MPI_COMM_WORLD); |
818 |
|
819 |
for (integrableObject = mol->beginIntegrableObject(ii); |
820 |
integrableObject != NULL; |
821 |
integrableObject = mol->nextIntegrableObject(ii)){ |
822 |
|
823 |
MPI_Recv(&angleTransfer, 1, MPI_DOUBLE, 0, |
824 |
TAKE_THIS_TAG_DOUBLE, MPI_COMM_WORLD, &istatus); |
825 |
vecParticles[j]->setZangle(angleTransfer); |
826 |
} |
827 |
} |
828 |
} |
829 |
} |
830 |
#endif |
831 |
} |
832 |
|
833 |
} // end namespace oopse |