ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/EAM_FF.cpp
Revision: 669
Committed: Thu Aug 7 00:47:33 2003 UTC (20 years, 10 months ago) by chuckv
File size: 26095 byte(s)
Log Message:
Bug fixes for eam...

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