ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/EAM_FF.cpp
Revision: 941
Committed: Tue Jan 13 23:01:43 2004 UTC (20 years, 5 months ago) by gezelter
File size: 26096 byte(s)
Log Message:
Changes for adding direct charge-charge interactions (with switching function)

File Contents

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