ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/EAM_FF.cpp
Revision: 1097
Committed: Mon Apr 12 20:32:20 2004 UTC (20 years, 2 months ago) by gezelter
File size: 26012 byte(s)
Log Message:
Changes for RigidBody dynamics (Somewhat extensive)

File Contents

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