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/eam_interface.h" |
17 |
|
18 |
//#include "UseTheForce/fortranWrappers.hpp" |
19 |
|
20 |
#ifdef IS_MPI |
21 |
#include "UseTheForce/mpiForceField.h" |
22 |
#endif // is_mpi |
23 |
|
24 |
|
25 |
|
26 |
namespace EAM_NS{ |
27 |
|
28 |
// Declare the structures that will be passed by the parser and MPI |
29 |
|
30 |
typedef struct{ |
31 |
char name[15]; |
32 |
double mass; |
33 |
double lattice_constant; |
34 |
double eam_drho; // The distance between each of the points indexed by rho. |
35 |
double eam_rcut; // The cutoff radius for eam. |
36 |
double eam_dr; // The distance between each of the rho points. |
37 |
int eam_nrho; // Number of points indexed by rho |
38 |
int eam_nr; // The number of points based on r (Both Phi(r) and Rho(r)). |
39 |
int eam_ident; // Atomic number |
40 |
int ident; |
41 |
int last; // 0 -> default |
42 |
// 1 -> in MPI: tells nodes to stop listening |
43 |
} atomStruct; |
44 |
|
45 |
int parseAtom( char *lineBuffer, int lineNum, atomStruct &info, string eamPotFile ); |
46 |
int parseEAM( atomStruct &info, string ffPath, string eamPotFile, |
47 |
double **eam_rvals, double **eam_rhovals, |
48 |
double **eam_Frhovals); |
49 |
#ifdef IS_MPI |
50 |
|
51 |
MPI_Datatype mpiAtomStructType; |
52 |
|
53 |
#endif |
54 |
|
55 |
class LinkedAtomType { |
56 |
public: |
57 |
LinkedAtomType(){ |
58 |
next = NULL; |
59 |
name[0] = '\0'; |
60 |
eam_rvals = NULL; |
61 |
eam_rhovals = NULL; |
62 |
eam_Frhovals = NULL; |
63 |
} |
64 |
|
65 |
~LinkedAtomType(){ |
66 |
if( next != NULL ) delete next; |
67 |
if( eam_rvals != NULL ) delete[] eam_rvals; |
68 |
if( eam_rhovals != NULL ) delete[] eam_rhovals; |
69 |
if( eam_Frhovals != NULL ) delete[] eam_Frhovals; |
70 |
} |
71 |
|
72 |
LinkedAtomType* find(const char* key){ |
73 |
if( !strcmp(name, key) ) return this; |
74 |
if( next != NULL ) return next->find(key); |
75 |
return NULL; |
76 |
} |
77 |
|
78 |
|
79 |
void add( atomStruct &info, double *the_eam_rvals, |
80 |
double *the_eam_rhovals,double *the_eam_Frhovals ){ |
81 |
|
82 |
// check for duplicates |
83 |
|
84 |
if( !strcmp( info.name, name ) ){ |
85 |
sprintf( painCave.errMsg, |
86 |
"Duplicate EAM atom type \"%s\" found in " |
87 |
"the EAM_FF param file./n", |
88 |
name ); |
89 |
painCave.isFatal = 1; |
90 |
simError(); |
91 |
} |
92 |
|
93 |
if( next != NULL ) next->add(info, the_eam_rvals, the_eam_rhovals, the_eam_Frhovals); |
94 |
else{ |
95 |
next = new LinkedAtomType(); |
96 |
strcpy(next->name, info.name); |
97 |
next->mass = info.mass; |
98 |
next->lattice_constant = info.lattice_constant; |
99 |
next->eam_nrho = info.eam_nrho; |
100 |
next->eam_drho = info.eam_drho; |
101 |
next->eam_nr = info.eam_nr; |
102 |
next->eam_dr = info.eam_dr; |
103 |
next->eam_rcut = info.eam_rcut; |
104 |
next->eam_ident = info.eam_ident; |
105 |
next->ident = info.ident; |
106 |
|
107 |
next->eam_rvals = the_eam_rvals; |
108 |
next->eam_rhovals = the_eam_rhovals; |
109 |
next->eam_Frhovals = the_eam_Frhovals; |
110 |
} |
111 |
} |
112 |
|
113 |
|
114 |
#ifdef IS_MPI |
115 |
|
116 |
void duplicate( atomStruct &info ){ |
117 |
strcpy(info.name, name); |
118 |
info.mass = mass; |
119 |
info.lattice_constant = lattice_constant; |
120 |
info.eam_nrho = eam_nrho; |
121 |
info.eam_drho = eam_drho; |
122 |
info.eam_nr = eam_nr; |
123 |
info.eam_dr = eam_dr; |
124 |
info.eam_rcut = eam_rcut; |
125 |
info.eam_ident = eam_ident; |
126 |
info.ident = ident; |
127 |
info.last = 0; |
128 |
} |
129 |
|
130 |
|
131 |
#endif |
132 |
|
133 |
char name[15]; |
134 |
double mass; |
135 |
double lattice_constant; |
136 |
int eam_nrho; // Number of points indexed by rho |
137 |
double eam_drho; // The distance between each of the points indexed by rho. |
138 |
int eam_nr; // The number of points based on r (Both Phi(r) and Rho(r)). |
139 |
double eam_dr; // The distance between each of the rho points. |
140 |
double eam_rcut; // The cutoff radius for eam. |
141 |
|
142 |
double *eam_rvals; // Z of r values |
143 |
double *eam_rhovals; // rho of r values |
144 |
double *eam_Frhovals; // F of rho values |
145 |
int eam_ident; // eam identity (atomic number) |
146 |
int ident; |
147 |
LinkedAtomType* next; |
148 |
}; |
149 |
|
150 |
LinkedAtomType* headAtomType; |
151 |
LinkedAtomType* currentAtomType; |
152 |
|
153 |
} |
154 |
|
155 |
using namespace EAM_NS; |
156 |
|
157 |
//**************************************************************** |
158 |
// begins the actual forcefield stuff. |
159 |
//**************************************************************** |
160 |
|
161 |
EAM_FF::EAM_FF(){ |
162 |
EAM_FF(""); |
163 |
} |
164 |
|
165 |
EAM_FF::EAM_FF(const string &the_variant) : ForceFields(the_variant) { |
166 |
|
167 |
string fileName; |
168 |
string tempString; |
169 |
|
170 |
headAtomType = NULL; |
171 |
currentAtomType = NULL; |
172 |
|
173 |
// Set eamRcut to 0.0 |
174 |
eamRcut = 0.0; |
175 |
|
176 |
#ifdef IS_MPI |
177 |
int i; |
178 |
|
179 |
// ********************************************************************** |
180 |
// Init the atomStruct mpi type |
181 |
|
182 |
atomStruct atomProto; // mpiPrototype |
183 |
int atomBC[3] = {15,5,5}; // block counts |
184 |
MPI_Aint atomDspls[3]; // displacements |
185 |
MPI_Datatype atomMbrTypes[3]; // member mpi types |
186 |
|
187 |
MPI_Address(&atomProto.name, &atomDspls[0]); |
188 |
MPI_Address(&atomProto.mass, &atomDspls[1]); |
189 |
MPI_Address(&atomProto.eam_nrho, &atomDspls[2]); |
190 |
|
191 |
atomMbrTypes[0] = MPI_CHAR; |
192 |
atomMbrTypes[1] = MPI_DOUBLE; |
193 |
atomMbrTypes[2] = MPI_INT; |
194 |
|
195 |
for (i=2; i >= 0; i--) atomDspls[i] -= atomDspls[0]; |
196 |
|
197 |
MPI_Type_struct(3, atomBC, atomDspls, atomMbrTypes, &mpiAtomStructType); |
198 |
MPI_Type_commit(&mpiAtomStructType); |
199 |
|
200 |
// *********************************************************************** |
201 |
|
202 |
if( worldRank == 0 ){ |
203 |
#endif |
204 |
|
205 |
// generate the force file name |
206 |
|
207 |
fileName = "EAM"; |
208 |
|
209 |
if (!the_variant.empty()) { |
210 |
has_variant = 1; |
211 |
fileName += "." + the_variant; |
212 |
|
213 |
sprintf( painCave.errMsg, |
214 |
"Using %s variant of EAM force field.\n", |
215 |
the_variant.c_str() ); |
216 |
painCave.severity = OOPSE_INFO; |
217 |
painCave.isFatal = 0; |
218 |
simError(); |
219 |
} |
220 |
|
221 |
fileName += ".frc"; |
222 |
|
223 |
//fprintf( stderr,"Trying to open %s\n", fileName ); |
224 |
|
225 |
// attempt to open the file in the current directory first. |
226 |
|
227 |
frcFile = fopen( fileName.c_str(), "r" ); |
228 |
|
229 |
if( frcFile == NULL ){ |
230 |
|
231 |
// next see if the force path enviorment variable is set |
232 |
|
233 |
tempString = ffPath + "/" + fileName; |
234 |
fileName = tempString; |
235 |
|
236 |
frcFile = fopen( fileName.c_str(), "r" ); |
237 |
|
238 |
if( frcFile == NULL ){ |
239 |
|
240 |
sprintf( painCave.errMsg, |
241 |
"Error opening the force field parameter file:\n" |
242 |
"\t%s\n" |
243 |
"\tHave you tried setting the FORCE_PARAM_PATH environment " |
244 |
"variable?\n", |
245 |
fileName.c_str() ); |
246 |
painCave.severity = OOPSE_ERROR; |
247 |
painCave.isFatal = 1; |
248 |
simError(); |
249 |
} |
250 |
} |
251 |
|
252 |
|
253 |
#ifdef IS_MPI |
254 |
} |
255 |
|
256 |
sprintf( checkPointMsg, "EAM_FF file opened sucessfully." ); |
257 |
MPIcheckPoint(); |
258 |
|
259 |
#endif // is_mpi |
260 |
} |
261 |
|
262 |
|
263 |
EAM_FF::~EAM_FF(){ |
264 |
|
265 |
if( headAtomType != NULL ) delete headAtomType; |
266 |
|
267 |
#ifdef IS_MPI |
268 |
if( worldRank == 0 ){ |
269 |
#endif // is_mpi |
270 |
|
271 |
fclose( frcFile ); |
272 |
|
273 |
#ifdef IS_MPI |
274 |
} |
275 |
#endif // is_mpi |
276 |
} |
277 |
|
278 |
|
279 |
void EAM_FF::calcRcut( void ){ |
280 |
|
281 |
#ifdef IS_MPI |
282 |
double tempEamRcut = eamRcut; |
283 |
MPI_Allreduce( &tempEamRcut, &eamRcut, 1, MPI_DOUBLE, MPI_MAX, |
284 |
MPI_COMM_WORLD); |
285 |
#endif //is_mpi |
286 |
entry_plug->setDefaultRcut(eamRcut); |
287 |
} |
288 |
|
289 |
|
290 |
void EAM_FF::initForceField( ){ |
291 |
initFortran(0); |
292 |
} |
293 |
|
294 |
void EAM_FF::cleanMe( void ){ |
295 |
|
296 |
#ifdef IS_MPI |
297 |
|
298 |
// keep the linked list in the mpi version |
299 |
|
300 |
#else // is_mpi |
301 |
|
302 |
// delete the linked list in the single processor version |
303 |
|
304 |
if( headAtomType != NULL ) delete headAtomType; |
305 |
|
306 |
#endif // is_mpi |
307 |
} |
308 |
|
309 |
|
310 |
void EAM_FF::readParams( void ){ |
311 |
|
312 |
atomStruct info; |
313 |
info.last = 1; // initialize last to have the last set. |
314 |
// if things go well, last will be set to 0 |
315 |
|
316 |
int identNum; |
317 |
double *eam_rvals; // Z of r values |
318 |
double *eam_rhovals; // rho of r values |
319 |
double *eam_Frhovals; // F of rho values |
320 |
string eamPotFile; |
321 |
|
322 |
|
323 |
|
324 |
bigSigma = 0.0; |
325 |
#ifdef IS_MPI |
326 |
if( worldRank == 0 ){ |
327 |
#endif |
328 |
|
329 |
// read in the atom types. |
330 |
|
331 |
headAtomType = new LinkedAtomType; |
332 |
|
333 |
fastForward( "AtomTypes", "eam atom readParams" ); |
334 |
|
335 |
// we are now at the AtomTypes section. |
336 |
|
337 |
eof_test = fgets( readLine, sizeof(readLine), frcFile ); |
338 |
lineNum++; |
339 |
|
340 |
|
341 |
// read a line, and start parseing out the atom types |
342 |
|
343 |
if( eof_test == NULL ){ |
344 |
sprintf( painCave.errMsg, |
345 |
"Error in reading Atoms from force file at line %d.\n", |
346 |
lineNum ); |
347 |
painCave.isFatal = 1; |
348 |
simError(); |
349 |
} |
350 |
|
351 |
identNum = 1; |
352 |
// stop reading at end of file, or at next section |
353 |
|
354 |
while( readLine[0] != '#' && eof_test != NULL ){ |
355 |
|
356 |
// toss comment lines |
357 |
if( readLine[0] != '!' ){ |
358 |
|
359 |
// the parser returns 0 if the line was blank |
360 |
if( parseAtom( readLine, lineNum, info, eamPotFile ) ){ |
361 |
parseEAM(info, ffPath, eamPotFile, &eam_rvals, |
362 |
&eam_rhovals, &eam_Frhovals); |
363 |
info.ident = identNum; |
364 |
headAtomType->add( info, eam_rvals, |
365 |
eam_rhovals,eam_Frhovals ); |
366 |
identNum++; |
367 |
} |
368 |
} |
369 |
eof_test = fgets( readLine, sizeof(readLine), frcFile ); |
370 |
lineNum++; |
371 |
} |
372 |
|
373 |
|
374 |
|
375 |
#ifdef IS_MPI |
376 |
|
377 |
|
378 |
// send out the linked list to all the other processes |
379 |
|
380 |
sprintf( checkPointMsg, |
381 |
"EAM_FF atom structures read successfully." ); |
382 |
MPIcheckPoint(); |
383 |
|
384 |
currentAtomType = headAtomType->next; //skip the first element who is a place holder. |
385 |
while( currentAtomType != NULL ){ |
386 |
currentAtomType->duplicate( info ); |
387 |
|
388 |
|
389 |
|
390 |
sendFrcStruct( &info, mpiAtomStructType ); |
391 |
|
392 |
// We have to now broadcast the Arrays |
393 |
MPI_Bcast(currentAtomType->eam_rvals, |
394 |
currentAtomType->eam_nr, |
395 |
MPI_DOUBLE,0,MPI_COMM_WORLD); |
396 |
MPI_Bcast(currentAtomType->eam_rhovals, |
397 |
currentAtomType->eam_nr, |
398 |
MPI_DOUBLE,0,MPI_COMM_WORLD); |
399 |
MPI_Bcast(currentAtomType->eam_Frhovals, |
400 |
currentAtomType->eam_nrho, |
401 |
MPI_DOUBLE,0,MPI_COMM_WORLD); |
402 |
|
403 |
sprintf( checkPointMsg, |
404 |
"successfully sent EAM force type: \"%s\"\n", |
405 |
info.name ); |
406 |
MPIcheckPoint(); |
407 |
|
408 |
currentAtomType = currentAtomType->next; |
409 |
} |
410 |
info.last = 1; |
411 |
sendFrcStruct( &info, mpiAtomStructType ); |
412 |
|
413 |
} |
414 |
|
415 |
else{ |
416 |
|
417 |
// listen for node 0 to send out the force params |
418 |
|
419 |
MPIcheckPoint(); |
420 |
|
421 |
headAtomType = new LinkedAtomType; |
422 |
receiveFrcStruct( &info, mpiAtomStructType ); |
423 |
|
424 |
while( !info.last ){ |
425 |
|
426 |
// allocate the arrays |
427 |
|
428 |
eam_rvals = new double[info.eam_nr]; |
429 |
eam_rhovals = new double[info.eam_nr]; |
430 |
eam_Frhovals = new double[info.eam_nrho]; |
431 |
|
432 |
// We have to now broadcast the Arrays |
433 |
MPI_Bcast(eam_rvals, |
434 |
info.eam_nr, |
435 |
MPI_DOUBLE,0,MPI_COMM_WORLD); |
436 |
MPI_Bcast(eam_rhovals, |
437 |
info.eam_nr, |
438 |
MPI_DOUBLE,0,MPI_COMM_WORLD); |
439 |
MPI_Bcast(eam_Frhovals, |
440 |
info.eam_nrho, |
441 |
MPI_DOUBLE,0,MPI_COMM_WORLD); |
442 |
|
443 |
|
444 |
headAtomType->add( info, eam_rvals, eam_rhovals, eam_Frhovals ); |
445 |
|
446 |
MPIcheckPoint(); |
447 |
|
448 |
receiveFrcStruct( &info, mpiAtomStructType ); |
449 |
|
450 |
|
451 |
} |
452 |
} |
453 |
#endif // is_mpi |
454 |
|
455 |
// call new A_types in fortran |
456 |
|
457 |
int isError; |
458 |
|
459 |
// dummy variables |
460 |
|
461 |
currentAtomType = headAtomType->next; |
462 |
while( currentAtomType != NULL ){ |
463 |
|
464 |
if( currentAtomType->name[0] != '\0' ){ |
465 |
|
466 |
AtomType* at = new AtomType(); |
467 |
at->setIdent(currentAtomType->ident); |
468 |
at->setName(currentAtomType->name); |
469 |
at->setEAM(); |
470 |
at->complete(); |
471 |
|
472 |
} |
473 |
|
474 |
currentAtomType = currentAtomType->next; |
475 |
} |
476 |
|
477 |
entry_plug->useLennardJones = 0; |
478 |
entry_plug->useEAM = 1; |
479 |
// Walk down again and send out EAM type |
480 |
currentAtomType = headAtomType->next; |
481 |
while( currentAtomType != NULL ){ |
482 |
|
483 |
if( currentAtomType->name[0] != '\0' ){ |
484 |
isError = 0; |
485 |
|
486 |
newEAMtype( &(currentAtomType->lattice_constant), |
487 |
&(currentAtomType->eam_nrho), |
488 |
&(currentAtomType->eam_drho), |
489 |
&(currentAtomType->eam_nr), |
490 |
&(currentAtomType->eam_dr), |
491 |
&(currentAtomType->eam_rcut), |
492 |
currentAtomType->eam_rvals, |
493 |
currentAtomType->eam_rhovals, |
494 |
currentAtomType->eam_Frhovals, |
495 |
&(currentAtomType->eam_ident), |
496 |
&isError); |
497 |
|
498 |
if( isError ){ |
499 |
sprintf( painCave.errMsg, |
500 |
"Error initializing the \"%s\" atom type in fortran EAM\n", |
501 |
currentAtomType->name ); |
502 |
painCave.isFatal = 1; |
503 |
simError(); |
504 |
} |
505 |
} |
506 |
currentAtomType = currentAtomType->next; |
507 |
} |
508 |
|
509 |
|
510 |
|
511 |
#ifdef IS_MPI |
512 |
sprintf( checkPointMsg, |
513 |
"EAM_FF atom structures successfully sent to fortran\n" ); |
514 |
MPIcheckPoint(); |
515 |
#endif // is_mpi |
516 |
|
517 |
|
518 |
|
519 |
} |
520 |
|
521 |
|
522 |
void EAM_FF::initializeAtoms( int nAtoms, Atom** the_atoms ){ |
523 |
|
524 |
int i; |
525 |
|
526 |
// initialize the atoms |
527 |
|
528 |
for( i=0; i<nAtoms; i++ ){ |
529 |
|
530 |
currentAtomType = headAtomType->find(the_atoms[i]->getType().c_str() ); |
531 |
if( currentAtomType == NULL ){ |
532 |
sprintf( painCave.errMsg, |
533 |
"AtomType error, %s not found in force file.\n", |
534 |
the_atoms[i]->getType().c_str() ); |
535 |
painCave.isFatal = 1; |
536 |
simError(); |
537 |
} |
538 |
|
539 |
the_atoms[i]->setMass( currentAtomType->mass ); |
540 |
the_atoms[i]->setIdent( currentAtomType->ident ); |
541 |
|
542 |
if (eamRcut < currentAtomType->eam_rcut) eamRcut = currentAtomType->eam_rcut; |
543 |
|
544 |
} |
545 |
} |
546 |
|
547 |
void EAM_FF::initializeBonds( int nBonds, Bond** BondArray, |
548 |
bond_pair* the_bonds ){ |
549 |
|
550 |
if( nBonds ){ |
551 |
sprintf( painCave.errMsg, |
552 |
"LJ_FF does not support bonds.\n" ); |
553 |
painCave.isFatal = 1; |
554 |
simError(); |
555 |
} |
556 |
} |
557 |
|
558 |
void EAM_FF::initializeBends( int nBends, Bend** bendArray, |
559 |
bend_set* the_bends ){ |
560 |
|
561 |
if( nBends ){ |
562 |
sprintf( painCave.errMsg, |
563 |
"LJ_FF does not support bends.\n" ); |
564 |
painCave.isFatal = 1; |
565 |
simError(); |
566 |
} |
567 |
} |
568 |
|
569 |
void EAM_FF::initializeTorsions( int nTorsions, Torsion** torsionArray, |
570 |
torsion_set* the_torsions ){ |
571 |
|
572 |
if( nTorsions ){ |
573 |
sprintf( painCave.errMsg, |
574 |
"LJ_FF does not support torsions.\n" ); |
575 |
painCave.isFatal = 1; |
576 |
simError(); |
577 |
} |
578 |
} |
579 |
|
580 |
void EAM_FF::fastForward( char* stopText, char* searchOwner ){ |
581 |
|
582 |
int foundText = 0; |
583 |
char* the_token; |
584 |
|
585 |
rewind( frcFile ); |
586 |
lineNum = 0; |
587 |
|
588 |
eof_test = fgets( readLine, sizeof(readLine), frcFile ); |
589 |
lineNum++; |
590 |
if( eof_test == NULL ){ |
591 |
sprintf( painCave.errMsg, "Error fast forwarding force file for %s: " |
592 |
" file is empty.\n", |
593 |
searchOwner ); |
594 |
painCave.isFatal = 1; |
595 |
simError(); |
596 |
} |
597 |
|
598 |
|
599 |
while( !foundText ){ |
600 |
while( eof_test != NULL && readLine[0] != '#' ){ |
601 |
eof_test = fgets( readLine, sizeof(readLine), frcFile ); |
602 |
lineNum++; |
603 |
} |
604 |
if( eof_test == NULL ){ |
605 |
sprintf( painCave.errMsg, |
606 |
"Error fast forwarding force file for %s at " |
607 |
"line %d: file ended unexpectedly.\n", |
608 |
searchOwner, |
609 |
lineNum ); |
610 |
painCave.isFatal = 1; |
611 |
simError(); |
612 |
} |
613 |
|
614 |
the_token = strtok( readLine, " ,;\t#\n" ); |
615 |
foundText = !strcmp( stopText, the_token ); |
616 |
|
617 |
if( !foundText ){ |
618 |
eof_test = fgets( readLine, sizeof(readLine), frcFile ); |
619 |
lineNum++; |
620 |
|
621 |
if( eof_test == NULL ){ |
622 |
sprintf( painCave.errMsg, |
623 |
"Error fast forwarding force file for %s at " |
624 |
"line %d: file ended unexpectedly.\n", |
625 |
searchOwner, |
626 |
lineNum ); |
627 |
painCave.isFatal = 1; |
628 |
simError(); |
629 |
} |
630 |
} |
631 |
} |
632 |
} |
633 |
|
634 |
|
635 |
|
636 |
int EAM_NS::parseAtom( char *lineBuffer, int lineNum, atomStruct &info, string eamPotFile ){ |
637 |
|
638 |
char* the_token; |
639 |
|
640 |
the_token = strtok( lineBuffer, " \n\t,;" ); |
641 |
if( the_token != NULL ){ |
642 |
|
643 |
strcpy( info.name, the_token ); |
644 |
|
645 |
if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){ |
646 |
sprintf( painCave.errMsg, |
647 |
"Error parseing AtomTypes: line %d\n", lineNum ); |
648 |
painCave.isFatal = 1; |
649 |
simError(); |
650 |
} |
651 |
|
652 |
info.mass = atof( the_token ); |
653 |
|
654 |
if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){ |
655 |
sprintf( painCave.errMsg, |
656 |
"Error parseing AtomTypes: line %d\n", lineNum ); |
657 |
painCave.isFatal = 1; |
658 |
simError(); |
659 |
} |
660 |
|
661 |
eamPotFile = the_token; |
662 |
|
663 |
return 1; |
664 |
} |
665 |
else return 0; |
666 |
} |
667 |
|
668 |
int EAM_NS::parseEAM(atomStruct &info, string ffPath, string eamPotFile, |
669 |
double **eam_rvals, |
670 |
double **eam_rhovals, |
671 |
double **eam_Frhovals){ |
672 |
double* myEam_rvals; |
673 |
double* myEam_rhovals; |
674 |
double* myEam_Frhovals; |
675 |
|
676 |
char* the_token; |
677 |
char* eam_eof_test; |
678 |
FILE *eamFile; |
679 |
const int BUFFERSIZE = 3000; |
680 |
|
681 |
string tempString; |
682 |
int linenumber; |
683 |
int nReadLines; |
684 |
char eam_read_buffer[BUFFERSIZE]; |
685 |
|
686 |
int i,j; |
687 |
|
688 |
linenumber = 0; |
689 |
|
690 |
// Open eam file |
691 |
eamFile = fopen( eamPotFile.c_str(), "r" ); |
692 |
|
693 |
if( eamFile == NULL ){ |
694 |
|
695 |
// next see if the force path enviorment variable is set |
696 |
|
697 |
tempString = ffPath + "/" + eamPotFile; |
698 |
eamPotFile = tempString; |
699 |
|
700 |
eamFile = fopen( eamPotFile.c_str(), "r" ); |
701 |
|
702 |
if( eamFile == NULL ){ |
703 |
|
704 |
sprintf( painCave.errMsg, |
705 |
"Error opening the EAM force parameter file: %s\n" |
706 |
"Have you tried setting the FORCE_PARAM_PATH environment " |
707 |
"variable?\n", |
708 |
eamPotFile.c_str() ); |
709 |
painCave.isFatal = 1; |
710 |
simError(); |
711 |
} |
712 |
} |
713 |
|
714 |
// First line is a comment line, read and toss it.... |
715 |
eam_eof_test = fgets(eam_read_buffer, sizeof(eam_read_buffer),eamFile); |
716 |
linenumber++; |
717 |
if(eam_eof_test == NULL){ |
718 |
sprintf( painCave.errMsg, |
719 |
"error in reading commment in %s\n", eamPotFile.c_str()); |
720 |
painCave.isFatal = 1; |
721 |
simError(); |
722 |
} |
723 |
|
724 |
|
725 |
|
726 |
// The Second line contains atomic number, atomic mass and a lattice constant |
727 |
eam_eof_test = fgets(eam_read_buffer, sizeof(eam_read_buffer),eamFile); |
728 |
linenumber++; |
729 |
if(eam_eof_test == NULL){ |
730 |
sprintf( painCave.errMsg, |
731 |
"error in reading Identifier line in %s\n", eamPotFile.c_str()); |
732 |
painCave.isFatal = 1; |
733 |
simError(); |
734 |
} |
735 |
|
736 |
|
737 |
|
738 |
|
739 |
if ( (the_token = strtok( eam_read_buffer, " \n\t,;")) == NULL){ |
740 |
sprintf( painCave.errMsg, |
741 |
"Error parsing EAM ident line in %s\n", eamPotFile.c_str() ); |
742 |
painCave.isFatal = 1; |
743 |
simError(); |
744 |
} |
745 |
|
746 |
info.eam_ident = atoi( the_token ); |
747 |
|
748 |
if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){ |
749 |
sprintf( painCave.errMsg, |
750 |
"Error parsing EAM mass in %s\n", eamPotFile.c_str() ); |
751 |
painCave.isFatal = 1; |
752 |
simError(); |
753 |
} |
754 |
info.mass = atof( the_token); |
755 |
|
756 |
if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){ |
757 |
sprintf( painCave.errMsg, |
758 |
"Error parsing EAM Lattice Constant %s\n", eamPotFile.c_str() ); |
759 |
painCave.isFatal = 1; |
760 |
simError(); |
761 |
} |
762 |
info.lattice_constant = atof( the_token); |
763 |
|
764 |
// Next line is nrho, drho, nr, dr and rcut |
765 |
eam_eof_test = fgets(eam_read_buffer, sizeof(eam_read_buffer),eamFile); |
766 |
if(eam_eof_test == NULL){ |
767 |
sprintf( painCave.errMsg, |
768 |
"error in reading number of points line in %s\n", |
769 |
eamPotFile.c_str()); |
770 |
painCave.isFatal = 1; |
771 |
simError(); |
772 |
} |
773 |
|
774 |
if ( (the_token = strtok( eam_read_buffer, " \n\t,;")) == NULL){ |
775 |
sprintf( painCave.errMsg, |
776 |
"Error parseing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
777 |
painCave.isFatal = 1; |
778 |
simError(); |
779 |
} |
780 |
|
781 |
info.eam_nrho = atoi( the_token ); |
782 |
|
783 |
if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){ |
784 |
sprintf( painCave.errMsg, |
785 |
"Error parsing EAM drho in %s\n", eamPotFile.c_str() ); |
786 |
painCave.isFatal = 1; |
787 |
simError(); |
788 |
} |
789 |
info.eam_drho = atof( the_token); |
790 |
|
791 |
if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){ |
792 |
sprintf( painCave.errMsg, |
793 |
"Error parsing EAM # r in %s\n", eamPotFile.c_str() ); |
794 |
painCave.isFatal = 1; |
795 |
simError(); |
796 |
} |
797 |
info.eam_nr = atoi( the_token); |
798 |
|
799 |
if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){ |
800 |
sprintf( painCave.errMsg, |
801 |
"Error parsing EAM dr in %s\n", eamPotFile.c_str() ); |
802 |
painCave.isFatal = 1; |
803 |
simError(); |
804 |
} |
805 |
info.eam_dr = atof( the_token); |
806 |
|
807 |
if( ( the_token = strtok( NULL, " \n\t,;" ) ) == NULL ){ |
808 |
sprintf( painCave.errMsg, |
809 |
"Error parsing EAM rcut in %s\n", eamPotFile.c_str() ); |
810 |
painCave.isFatal = 1; |
811 |
simError(); |
812 |
} |
813 |
info.eam_rcut = atof( the_token); |
814 |
|
815 |
|
816 |
|
817 |
|
818 |
|
819 |
// Ok now we have to allocate point arrays and read in number of points |
820 |
// Index the arrays for fortran, starting at 1 |
821 |
myEam_Frhovals = new double[info.eam_nrho]; |
822 |
myEam_rvals = new double[info.eam_nr]; |
823 |
myEam_rhovals = new double[info.eam_nr]; |
824 |
|
825 |
// Parse F of rho vals. |
826 |
|
827 |
// Assume for now that we have a complete number of lines |
828 |
nReadLines = int(info.eam_nrho/5); |
829 |
|
830 |
|
831 |
|
832 |
for (i=0;i<nReadLines;i++){ |
833 |
j = i*5; |
834 |
|
835 |
// Read next line |
836 |
eam_eof_test = fgets(eam_read_buffer, sizeof(eam_read_buffer),eamFile); |
837 |
linenumber++; |
838 |
if(eam_eof_test == NULL){ |
839 |
sprintf( painCave.errMsg, |
840 |
"error in reading EAM file %s at line %d\n", |
841 |
eamPotFile.c_str(), linenumber); |
842 |
painCave.isFatal = 1; |
843 |
simError(); |
844 |
} |
845 |
|
846 |
// Parse 5 values on each line into array |
847 |
// Value 1 |
848 |
if ( (the_token = strtok( eam_read_buffer, " \n\t,;")) == NULL){ |
849 |
sprintf( painCave.errMsg, |
850 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
851 |
painCave.isFatal = 1; |
852 |
simError(); |
853 |
} |
854 |
|
855 |
myEam_Frhovals[j+0] = atof( the_token ); |
856 |
|
857 |
// Value 2 |
858 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
859 |
sprintf( painCave.errMsg, |
860 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
861 |
painCave.isFatal = 1; |
862 |
simError(); |
863 |
} |
864 |
|
865 |
myEam_Frhovals[j+1] = atof( the_token ); |
866 |
|
867 |
// Value 3 |
868 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
869 |
sprintf( painCave.errMsg, |
870 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
871 |
painCave.isFatal = 1; |
872 |
simError(); |
873 |
} |
874 |
|
875 |
myEam_Frhovals[j+2] = atof( the_token ); |
876 |
|
877 |
// Value 4 |
878 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
879 |
sprintf( painCave.errMsg, |
880 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
881 |
painCave.isFatal = 1; |
882 |
simError(); |
883 |
} |
884 |
|
885 |
myEam_Frhovals[j+3] = atof( the_token ); |
886 |
|
887 |
// Value 5 |
888 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
889 |
sprintf( painCave.errMsg, |
890 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
891 |
painCave.isFatal = 1; |
892 |
simError(); |
893 |
} |
894 |
|
895 |
myEam_Frhovals[j+4] = atof( the_token ); |
896 |
|
897 |
} |
898 |
// Parse Z of r vals |
899 |
|
900 |
// Assume for now that we have a complete number of lines |
901 |
nReadLines = int(info.eam_nr/5); |
902 |
|
903 |
for (i=0;i<nReadLines;i++){ |
904 |
j = i*5; |
905 |
|
906 |
// Read next line |
907 |
eam_eof_test = fgets(eam_read_buffer, sizeof(eam_read_buffer),eamFile); |
908 |
linenumber++; |
909 |
if(eam_eof_test == NULL){ |
910 |
sprintf( painCave.errMsg, |
911 |
"error in reading EAM file %s at line %d\n", |
912 |
eamPotFile.c_str(), linenumber); |
913 |
painCave.isFatal = 1; |
914 |
simError(); |
915 |
} |
916 |
|
917 |
// Parse 5 values on each line into array |
918 |
// Value 1 |
919 |
if ( (the_token = strtok( eam_read_buffer, " \n\t,;")) == NULL){ |
920 |
sprintf( painCave.errMsg, |
921 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
922 |
painCave.isFatal = 1; |
923 |
simError(); |
924 |
} |
925 |
|
926 |
myEam_rvals[j+0] = atof( the_token ); |
927 |
|
928 |
// Value 2 |
929 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
930 |
sprintf( painCave.errMsg, |
931 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
932 |
painCave.isFatal = 1; |
933 |
simError(); |
934 |
} |
935 |
|
936 |
myEam_rvals[j+1] = atof( the_token ); |
937 |
|
938 |
// Value 3 |
939 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
940 |
sprintf( painCave.errMsg, |
941 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
942 |
painCave.isFatal = 1; |
943 |
simError(); |
944 |
} |
945 |
|
946 |
myEam_rvals[j+2] = atof( the_token ); |
947 |
|
948 |
// Value 4 |
949 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
950 |
sprintf( painCave.errMsg, |
951 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
952 |
painCave.isFatal = 1; |
953 |
simError(); |
954 |
} |
955 |
|
956 |
myEam_rvals[j+3] = atof( the_token ); |
957 |
|
958 |
// Value 5 |
959 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
960 |
sprintf( painCave.errMsg, |
961 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
962 |
painCave.isFatal = 1; |
963 |
simError(); |
964 |
} |
965 |
|
966 |
myEam_rvals[j+4] = atof( the_token ); |
967 |
|
968 |
} |
969 |
// Parse rho of r vals |
970 |
|
971 |
// Assume for now that we have a complete number of lines |
972 |
|
973 |
for (i=0;i<nReadLines;i++){ |
974 |
j = i*5; |
975 |
|
976 |
// Read next line |
977 |
eam_eof_test = fgets(eam_read_buffer, sizeof(eam_read_buffer),eamFile); |
978 |
linenumber++; |
979 |
if(eam_eof_test == NULL){ |
980 |
sprintf( painCave.errMsg, |
981 |
"error in reading EAM file %s at line %d\n", |
982 |
eamPotFile.c_str(), linenumber); |
983 |
painCave.isFatal = 1; |
984 |
simError(); |
985 |
} |
986 |
|
987 |
// Parse 5 values on each line into array |
988 |
// Value 1 |
989 |
if ( (the_token = strtok( eam_read_buffer, " \n\t,;")) == NULL){ |
990 |
sprintf( painCave.errMsg, |
991 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
992 |
painCave.isFatal = 1; |
993 |
simError(); |
994 |
} |
995 |
|
996 |
myEam_rhovals[j+0] = atof( the_token ); |
997 |
|
998 |
// Value 2 |
999 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
1000 |
sprintf( painCave.errMsg, |
1001 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
1002 |
painCave.isFatal = 1; |
1003 |
simError(); |
1004 |
} |
1005 |
|
1006 |
myEam_rhovals[j+1] = atof( the_token ); |
1007 |
|
1008 |
// Value 3 |
1009 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
1010 |
sprintf( painCave.errMsg, |
1011 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
1012 |
painCave.isFatal = 1; |
1013 |
simError(); |
1014 |
} |
1015 |
|
1016 |
myEam_rhovals[j+2] = atof( the_token ); |
1017 |
|
1018 |
// Value 4 |
1019 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
1020 |
sprintf( painCave.errMsg, |
1021 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
1022 |
painCave.isFatal = 1; |
1023 |
simError(); |
1024 |
} |
1025 |
|
1026 |
myEam_rhovals[j+3] = atof( the_token ); |
1027 |
|
1028 |
// Value 5 |
1029 |
if ( (the_token = strtok( NULL, " \n\t,;")) == NULL){ |
1030 |
sprintf( painCave.errMsg, |
1031 |
"Error parsing EAM nrho: line in %s\n", eamPotFile.c_str() ); |
1032 |
painCave.isFatal = 1; |
1033 |
simError(); |
1034 |
} |
1035 |
|
1036 |
myEam_rhovals[j+4] = atof( the_token ); |
1037 |
|
1038 |
} |
1039 |
*eam_rvals = myEam_rvals; |
1040 |
*eam_rhovals = myEam_rhovals; |
1041 |
*eam_Frhovals = myEam_Frhovals; |
1042 |
|
1043 |
fclose(eamFile); |
1044 |
return 0; |
1045 |
} |