ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/EAM_FF.cpp
Revision: 627
Committed: Wed Jul 16 21:49:59 2003 UTC (20 years, 11 months ago) by chuckv
File size: 25125 byte(s)
Log Message:
More up to date version of EAM_FF

File Contents

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