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