ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/UseTheForce/EAM_FF.cpp
Revision: 1654
Committed: Wed Oct 27 02:16:34 2004 UTC (19 years, 8 months ago) by gezelter
File size: 25921 byte(s)
Log Message:
more char* -> string conversion

File Contents

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