ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-4/src/UseTheForce/LJFF.cpp
Revision: 1702
Committed: Wed Nov 3 18:00:36 2004 UTC (19 years, 8 months ago) by tim
File size: 12004 byte(s)
Log Message:
ForceFiled get compiled. Still a long way to go ......

File Contents

# Content
1 #include <stdlib.h>
2 #include <stdio.h>
3 #include <string.h>
4
5 #include <iostream>
6 using namespace std;
7
8 #ifdef IS_MPI
9 #include <mpi.h>
10 #endif //is_mpi
11
12 #include "UseTheForce/ForceFields.hpp"
13 #include "primitives/SRI.hpp"
14 #include "utils/simError.h"
15 #include "types/AtomType.hpp"
16 #include "UseTheForce/DarkSide/lj_interface.h"
17
18
19 #ifdef IS_MPI
20 #include "UseTheForce/mpiForceField.h"
21 #endif // is_mpi
22
23
24
25 namespace LJ_NS{
26
27 // Declare the structures that will be passed by the parser and MPI
28
29 typedef struct{
30 char name[15];
31 double mass;
32 double epslon;
33 double sigma;
34 int ident;
35 int last; // 0 -> default
36 // 1 -> in MPI: tells nodes to stop listening
37 } atomStruct;
38
39 int parseAtom( char *lineBuffer, int lineNum, atomStruct &info );
40
41 #ifdef IS_MPI
42
43 MPI_Datatype mpiAtomStructType;
44
45 #endif
46
47 class LinkedAtomType {
48 public:
49 LinkedAtomType(){
50 next = NULL;
51 name[0] = '\0';
52 }
53 ~LinkedAtomType(){ if( next != NULL ) delete next; }
54
55 LinkedAtomType* find(const char* key){
56 if( !strcmp(name, key) ) return this;
57 if( next != NULL ) return next->find(key);
58 return NULL;
59 }
60
61
62 void add( atomStruct &info ){
63
64 // check for duplicates
65
66 if( !strcmp( info.name, name ) ){
67 sprintf( painCave.errMsg,
68 "Duplicate LJ atom type \"%s\" found in "
69 "the LJFF param file./n",
70 name );
71 painCave.isFatal = 1;
72 simError();
73 }
74
75 if( next != NULL ) next->add(info);
76 else{
77 next = new LinkedAtomType();
78 strcpy(next->name, info.name);
79 next->mass = info.mass;
80 next->epslon = info.epslon;
81 next->sigma = info.sigma;
82 next->ident = info.ident;
83 }
84 }
85
86
87 #ifdef IS_MPI
88
89 void duplicate( atomStruct &info ){
90 strcpy(info.name, name);
91 info.mass = mass;
92 info.epslon = epslon;
93 info.sigma = sigma;
94 info.ident = ident;
95 info.last = 0;
96 }
97
98
99 #endif
100
101 char name[15];
102 double mass;
103 double epslon;
104 double sigma;
105 int ident;
106 LinkedAtomType* next;
107 };
108
109 LinkedAtomType* headAtomType;
110 LinkedAtomType* currentAtomType;
111
112 }
113
114 using namespace LJ_NS;
115
116 //****************************************************************
117 // begins the actual forcefield stuff.
118 //****************************************************************
119
120
121 LJFF::LJFF() : ForceFields() {
122
123 string fileName;
124 string tempString;
125
126 headAtomType = NULL;
127 currentAtomType = NULL;
128
129 // do the funtion wrapping
130 // wrapMeFF( this );
131
132 #ifdef IS_MPI
133 int i;
134
135 // **********************************************************************
136 // Init the atomStruct mpi type
137
138 atomStruct atomProto; // mpiPrototype
139 int atomBC[3] = {15,3,2}; // block counts
140 MPI_Aint atomDspls[3]; // displacements
141 MPI_Datatype atomMbrTypes[3]; // member mpi types
142
143 MPI_Address(&atomProto.name, &atomDspls[0]);
144 MPI_Address(&atomProto.mass, &atomDspls[1]);
145 MPI_Address(&atomProto.ident, &atomDspls[2]);
146
147 atomMbrTypes[0] = MPI_CHAR;
148 atomMbrTypes[1] = MPI_DOUBLE;
149 atomMbrTypes[2] = MPI_INT;
150
151 for (i=2; i >= 0; i--) atomDspls[i] -= atomDspls[0];
152
153 MPI_Type_struct(3, atomBC, atomDspls, atomMbrTypes, &mpiAtomStructType);
154 MPI_Type_commit(&mpiAtomStructType);
155
156 // ***********************************************************************
157
158 if( worldRank == 0 ){
159 #endif
160
161 // generate the force file name
162
163 fileName = "LJFF.frc";
164 // fprintf( stderr,"Trying to open %s\n", fileName );
165
166 // attempt to open the file in the current directory first.
167
168 frcFile = fopen( fileName.c_str(), "r" );
169
170 if( frcFile == NULL ){
171
172 // next see if the force path enviorment variable is set
173
174 tempString = ffPath + "/" + fileName;
175 fileName = tempString;
176
177 frcFile = fopen( fileName.c_str(), "r" );
178
179 if( frcFile == NULL ){
180
181 sprintf( painCave.errMsg,
182 "Error opening the force field parameter file:\n"
183 "\t%s\n"
184 "\tHave you tried setting the FORCE_PARAM_PATH environment "
185 "variable?\n",
186 fileName.c_str() );
187 painCave.severity = OOPSE_ERROR;
188 painCave.isFatal = 1;
189 simError();
190 }
191 }
192
193 #ifdef IS_MPI
194 }
195
196 sprintf( checkPointMsg, "LJFF file opened sucessfully." );
197 MPIcheckPoint();
198
199 #endif // is_mpi
200 }
201
202
203 LJFF::~LJFF(){
204
205 if( headAtomType != NULL ) delete headAtomType;
206
207 #ifdef IS_MPI
208 if( worldRank == 0 ){
209 #endif // is_mpi
210
211 fclose( frcFile );
212
213 #ifdef IS_MPI
214 }
215 #endif // is_mpi
216 }
217
218 void LJFF::initForceField(){
219
220 initFortran(0);
221 }
222
223 void LJFF::cleanMe( void ){
224
225 #ifdef IS_MPI
226
227 // keep the linked list in the mpi version
228
229 #else // is_mpi
230
231 // delete the linked list in the single processor version
232
233 if( headAtomType != NULL ) delete headAtomType;
234
235 #endif // is_mpi
236 }
237
238 void LJFF::readParams( void ){
239
240 atomStruct info;
241 info.last = 1; // initialize last to have the last set.
242 // if things go well, last will be set to 0
243
244 int identNum;
245 int isError;
246
247 bigSigma = 0.0;
248 #ifdef IS_MPI
249 if( worldRank == 0 ){
250 #endif
251
252 // read in the atom types.
253
254 headAtomType = new LinkedAtomType;
255
256 fastForward( "AtomTypes", "initializeAtoms" );
257
258 // we are now at the AtomTypes section.
259
260 eof_test = fgets( readLine, sizeof(readLine), frcFile );
261 lineNum++;
262
263
264 // read a line, and start parseing out the atom types
265
266 if( eof_test == NULL ){
267 sprintf( painCave.errMsg,
268 "Error in reading Atoms from force file at line %d.\n",
269 lineNum );
270 painCave.isFatal = 1;
271 simError();
272 }
273
274 identNum = 1;
275 // stop reading at end of file, or at next section
276 while( readLine[0] != '#' && eof_test != NULL ){
277
278 // toss comment lines
279 if( readLine[0] != '!' ){
280
281 // the parser returns 0 if the line was blank
282 if( parseAtom( readLine, lineNum, info ) ){
283 info.ident = identNum;
284 headAtomType->add( info );;
285 identNum++;
286 }
287 }
288 eof_test = fgets( readLine, sizeof(readLine), frcFile );
289 lineNum++;
290 }
291
292 #ifdef IS_MPI
293
294 // send out the linked list to all the other processes
295
296 sprintf( checkPointMsg,
297 "LJFF atom structures read successfully." );
298 MPIcheckPoint();
299
300 currentAtomType = headAtomType->next; //skip the first element who is a place holder.
301 while( currentAtomType != NULL ){
302 currentAtomType->duplicate( info );
303
304
305
306 sendFrcStruct( &info, mpiAtomStructType );
307
308 sprintf( checkPointMsg,
309 "successfully sent lJ force type: \"%s\"\n",
310 info.name );
311 MPIcheckPoint();
312
313 currentAtomType = currentAtomType->next;
314 }
315 info.last = 1;
316 sendFrcStruct( &info, mpiAtomStructType );
317
318 }
319
320 else{
321
322 // listen for node 0 to send out the force params
323
324 MPIcheckPoint();
325
326 headAtomType = new LinkedAtomType;
327 receiveFrcStruct( &info, mpiAtomStructType );
328
329 while( !info.last ){
330
331
332
333 headAtomType->add( info );
334
335 MPIcheckPoint();
336
337 receiveFrcStruct( &info, mpiAtomStructType );
338 }
339 }
340 #endif // is_mpi
341
342 currentAtomType = headAtomType;
343 while( currentAtomType != NULL ){
344
345 if( currentAtomType->name[0] != '\0' ){
346
347 AtomType* at = new AtomType();
348 at->setIdent(currentAtomType->ident);
349 printf ("currentName = %s\n", currentAtomType->name);
350 at->setName(currentAtomType->name);
351 printf("Did setName\n");
352 at->setLennardJones();
353 printf("Did setLennardJones\n");
354 at->complete();
355 printf("Did complete\n");
356
357 }
358 currentAtomType = currentAtomType->next;
359 }
360
361 currentAtomType = headAtomType;
362 while( currentAtomType != NULL ){
363
364 if( currentAtomType->name[0] != '\0' ){
365 isError = 0;
366 newLJtype( &(currentAtomType->ident), &(currentAtomType->sigma),
367 &(currentAtomType->epslon), &isError);
368 if( isError ){
369 sprintf( painCave.errMsg,
370 "Error initializing the \"%s\" LJ type in fortran\n",
371 currentAtomType->name );
372 painCave.isFatal = 1;
373 simError();
374 }
375 }
376
377 currentAtomType = currentAtomType->next;
378 }
379
380 entry_plug->useLennardJones = 1;
381
382 #ifdef IS_MPI
383 sprintf( checkPointMsg,
384 "LJFF atom structures successfully sent to fortran\n" );
385 MPIcheckPoint();
386 #endif // is_mpi
387
388 }
389
390 void LJFF::initializeAtoms( int nAtoms, Atom** the_atoms ){
391
392 int i;
393
394 // initialize the atoms
395
396
397 for( i=0; i<nAtoms; i++ ){
398
399 currentAtomType = headAtomType->find(the_atoms[i]->getType().c_str() );
400 if( currentAtomType == NULL ){
401 sprintf( painCave.errMsg,
402 "AtomType error, %s not found in force file.\n",
403 the_atoms[i]->getType().c_str() );
404 painCave.isFatal = 1;
405 simError();
406 }
407
408 the_atoms[i]->setMass( currentAtomType->mass );
409 the_atoms[i]->setIdent( currentAtomType->ident );
410
411 if( bigSigma < currentAtomType->sigma ) bigSigma = currentAtomType->sigma;
412 }
413 }
414
415 void LJFF::initializeBonds( int nBonds, Bond** BondArray,
416 bond_pair* the_bonds ){
417
418 if( nBonds ){
419 sprintf( painCave.errMsg,
420 "LJFF does not support bonds.\n" );
421 painCave.isFatal = 1;
422 simError();
423 }
424 }
425
426 void LJFF::initializeBends( int nBends, Bend** bendArray,
427 bend_set* the_bends ){
428
429 if( nBends ){
430 sprintf( painCave.errMsg,
431 "LJFF does not support bends.\n" );
432 painCave.isFatal = 1;
433 simError();
434 }
435 }
436
437 void LJFF::initializeTorsions( int nTorsions, Torsion** torsionArray,
438 torsion_set* the_torsions ){
439
440 if( nTorsions ){
441 sprintf( painCave.errMsg,
442 "LJFF does not support torsions.\n" );
443 painCave.isFatal = 1;
444 simError();
445 }
446 }
447
448 void LJFF::fastForward( char* stopText, char* searchOwner ){
449
450 int foundText = 0;
451 char* the_token;
452
453 rewind( frcFile );
454 lineNum = 0;
455
456 eof_test = fgets( readLine, sizeof(readLine), frcFile );
457 lineNum++;
458 if( eof_test == NULL ){
459 sprintf( painCave.errMsg, "Error fast forwarding force file for %s: "
460 " file is empty.\n",
461 searchOwner );
462 painCave.isFatal = 1;
463 simError();
464 }
465
466
467 while( !foundText ){
468 while( eof_test != NULL && readLine[0] != '#' ){
469 eof_test = fgets( readLine, sizeof(readLine), frcFile );
470 lineNum++;
471 }
472 if( eof_test == NULL ){
473 sprintf( painCave.errMsg,
474 "Error fast forwarding force file for %s at "
475 "line %d: file ended unexpectedly.\n",
476 searchOwner,
477 lineNum );
478 painCave.isFatal = 1;
479 simError();
480 }
481
482 the_token = strtok( readLine, " ,;\t#\n" );
483 foundText = !strcmp( stopText, the_token );
484
485 if( !foundText ){
486 eof_test = fgets( readLine, sizeof(readLine), frcFile );
487 lineNum++;
488
489 if( eof_test == NULL ){
490 sprintf( painCave.errMsg,
491 "Error fast forwarding force file for %s at "
492 "line %d: file ended unexpectedly.\n",
493 searchOwner,
494 lineNum );
495 painCave.isFatal = 1;
496 simError();
497 }
498 }
499 }
500 }
501
502
503
504 int LJ_NS::parseAtom( char *lineBuffer, int lineNum, atomStruct &info ){
505
506 char* the_token;
507
508 the_token = strtok( lineBuffer, " \n\t,;" );
509 if( the_token != NULL ){
510
511 strcpy( info.name, the_token );
512
513 if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){
514 sprintf( painCave.errMsg,
515 "Error parseing AtomTypes: line %d\n", lineNum );
516 painCave.isFatal = 1;
517 simError();
518 }
519
520 info.mass = atof( the_token );
521
522 if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){
523 sprintf( painCave.errMsg,
524 "Error parseing AtomTypes: line %d\n", lineNum );
525 painCave.isFatal = 1;
526 simError();
527 }
528
529 info.epslon = atof( the_token );
530
531 if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){
532 sprintf( painCave.errMsg,
533 "Error parseing AtomTypes: line %d\n", lineNum );
534 painCave.isFatal = 1;
535 simError();
536 }
537
538 info.sigma = atof( the_token );
539
540 return 1;
541 }
542 else return 0;
543 }
544
545