ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/io/RestReader.cpp
Revision: 2115
Committed: Fri Mar 11 00:43:19 2005 UTC (19 years, 3 months ago) by chrisfen
File size: 25935 byte(s)
Log Message:
fixed a bug in MPI restraints

File Contents

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