ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/EAM_FF.cpp
Revision: 1224
Committed: Wed Jun 2 18:27:52 2004 UTC (20 years, 11 months ago) by gezelter
File size: 26077 byte(s)
Log Message:
formatting error messages, dependency fixes

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:\n"
234 "\t%s\n"
235 "\tHave you tried setting the FORCE_PARAM_PATH environment "
236 "variable?\n",
237 fileName );
238 painCave.severity = OOPSE_ERROR;
239 painCave.isFatal = 1;
240 simError();
241 }
242 }
243
244 #ifdef IS_MPI
245 }
246
247 sprintf( checkPointMsg, "EAM_FF file opened sucessfully." );
248 MPIcheckPoint();
249
250 #endif // is_mpi
251 }
252
253
254 EAM_FF::~EAM_FF(){
255
256 if( headAtomType != NULL ) delete headAtomType;
257
258 #ifdef IS_MPI
259 if( worldRank == 0 ){
260 #endif // is_mpi
261
262 fclose( frcFile );
263
264 #ifdef IS_MPI
265 }
266 #endif // is_mpi
267 }
268
269
270 void EAM_FF::calcRcut( void ){
271
272 #ifdef IS_MPI
273 double tempEamRcut = eamRcut;
274 MPI_Allreduce( &tempEamRcut, &eamRcut, 1, MPI_DOUBLE, MPI_MAX,
275 MPI_COMM_WORLD);
276 #endif //is_mpi
277 entry_plug->setDefaultRcut(eamRcut);
278 }
279
280
281 void EAM_FF::initForceField( int ljMixRule ){
282 initFortran( ljMixRule, 0 );
283 }
284
285 void EAM_FF::cleanMe( void ){
286
287 #ifdef IS_MPI
288
289 // keep the linked list in the mpi version
290
291 #else // is_mpi
292
293 // delete the linked list in the single processor version
294
295 if( headAtomType != NULL ) delete headAtomType;
296
297 #endif // is_mpi
298 }
299
300
301 void EAM_FF::readParams( void ){
302
303 atomStruct info;
304 info.last = 1; // initialize last to have the last set.
305 // if things go well, last will be set to 0
306
307 int identNum;
308 double *eam_rvals; // Z of r values
309 double *eam_rhovals; // rho of r values
310 double *eam_Frhovals; // F of rho values
311 char eamPotFile[1000];
312
313
314 bigSigma = 0.0;
315 #ifdef IS_MPI
316 if( worldRank == 0 ){
317 #endif
318
319 // read in the atom types.
320
321 headAtomType = new LinkedAtomType;
322
323 fastForward( "AtomTypes", "eam atom readParams" );
324
325 // we are now at the AtomTypes section.
326
327 eof_test = fgets( readLine, sizeof(readLine), frcFile );
328 lineNum++;
329
330
331 // read a line, and start parseing out the atom types
332
333 if( eof_test == NULL ){
334 sprintf( painCave.errMsg,
335 "Error in reading Atoms from force file at line %d.\n",
336 lineNum );
337 painCave.isFatal = 1;
338 simError();
339 }
340
341 identNum = 1;
342 // stop reading at end of file, or at next section
343
344 while( readLine[0] != '#' && eof_test != NULL ){
345
346 // toss comment lines
347 if( readLine[0] != '!' ){
348
349 // the parser returns 0 if the line was blank
350 if( parseAtom( readLine, lineNum, info, eamPotFile ) ){
351 parseEAM(info,eamPotFile, &eam_rvals,
352 &eam_rhovals, &eam_Frhovals);
353 info.ident = identNum;
354 headAtomType->add( info, eam_rvals,
355 eam_rhovals,eam_Frhovals );
356 identNum++;
357 }
358 }
359 eof_test = fgets( readLine, sizeof(readLine), frcFile );
360 lineNum++;
361 }
362
363
364
365 #ifdef IS_MPI
366
367
368 // send out the linked list to all the other processes
369
370 sprintf( checkPointMsg,
371 "EAM_FF atom structures read successfully." );
372 MPIcheckPoint();
373
374 currentAtomType = headAtomType->next; //skip the first element who is a place holder.
375 while( currentAtomType != NULL ){
376 currentAtomType->duplicate( info );
377
378
379
380 sendFrcStruct( &info, mpiAtomStructType );
381
382 // We have to now broadcast the Arrays
383 MPI_Bcast(currentAtomType->eam_rvals,
384 currentAtomType->eam_nr,
385 MPI_DOUBLE,0,MPI_COMM_WORLD);
386 MPI_Bcast(currentAtomType->eam_rhovals,
387 currentAtomType->eam_nr,
388 MPI_DOUBLE,0,MPI_COMM_WORLD);
389 MPI_Bcast(currentAtomType->eam_Frhovals,
390 currentAtomType->eam_nrho,
391 MPI_DOUBLE,0,MPI_COMM_WORLD);
392
393 sprintf( checkPointMsg,
394 "successfully sent EAM force type: \"%s\"\n",
395 info.name );
396 MPIcheckPoint();
397
398 currentAtomType = currentAtomType->next;
399 }
400 info.last = 1;
401 sendFrcStruct( &info, mpiAtomStructType );
402
403 }
404
405 else{
406
407 // listen for node 0 to send out the force params
408
409 MPIcheckPoint();
410
411 headAtomType = new LinkedAtomType;
412 receiveFrcStruct( &info, mpiAtomStructType );
413
414 while( !info.last ){
415
416 // allocate the arrays
417
418 eam_rvals = new double[info.eam_nr];
419 eam_rhovals = new double[info.eam_nr];
420 eam_Frhovals = new double[info.eam_nrho];
421
422 // We have to now broadcast the Arrays
423 MPI_Bcast(eam_rvals,
424 info.eam_nr,
425 MPI_DOUBLE,0,MPI_COMM_WORLD);
426 MPI_Bcast(eam_rhovals,
427 info.eam_nr,
428 MPI_DOUBLE,0,MPI_COMM_WORLD);
429 MPI_Bcast(eam_Frhovals,
430 info.eam_nrho,
431 MPI_DOUBLE,0,MPI_COMM_WORLD);
432
433
434 headAtomType->add( info, eam_rvals, eam_rhovals, eam_Frhovals );
435
436 MPIcheckPoint();
437
438 receiveFrcStruct( &info, mpiAtomStructType );
439
440
441 }
442 }
443 #endif // is_mpi
444
445 // call new A_types in fortran
446
447 int isError;
448
449 // dummy variables
450 int isLJ = 0;
451 int isDipole = 0;
452 int isSSD = 0;
453 int isGB = 0;
454 int isEAM = 1;
455 int isCharge = 0;
456 double dipole = 0.0;
457 double charge = 0.0;
458 double eamSigma = 0.0;
459 double eamEpslon = 0.0;
460
461 currentAtomType = headAtomType->next;
462 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 &isEAM,
472 &isCharge,
473 &eamEpslon,
474 &eamSigma,
475 &charge,
476 &dipole,
477 &isError );
478 if( isError ){
479 sprintf( painCave.errMsg,
480 "Error initializing the \"%s\" atom type in fortran\n",
481 currentAtomType->name );
482 painCave.isFatal = 1;
483 simError();
484 }
485 }
486 currentAtomType = currentAtomType->next;
487 }
488
489 entry_plug->useLJ = 0;
490 entry_plug->useEAM = 1;
491 // Walk down again and send out EAM type
492 currentAtomType = headAtomType->next;
493 while( currentAtomType != NULL ){
494
495 if( currentAtomType->name[0] != '\0' ){
496 isError = 0;
497
498 newEAMtype( &(currentAtomType->lattice_constant),
499 &(currentAtomType->eam_nrho),
500 &(currentAtomType->eam_drho),
501 &(currentAtomType->eam_nr),
502 &(currentAtomType->eam_dr),
503 &(currentAtomType->eam_rcut),
504 currentAtomType->eam_rvals,
505 currentAtomType->eam_rhovals,
506 currentAtomType->eam_Frhovals,
507 &(currentAtomType->eam_ident),
508 &isError);
509
510 if( isError ){
511 sprintf( painCave.errMsg,
512 "Error initializing the \"%s\" atom type in fortran EAM\n",
513 currentAtomType->name );
514 painCave.isFatal = 1;
515 simError();
516 }
517 }
518 currentAtomType = currentAtomType->next;
519 }
520
521
522
523 #ifdef IS_MPI
524 sprintf( checkPointMsg,
525 "EAM_FF atom structures successfully sent to fortran\n" );
526 MPIcheckPoint();
527 #endif // is_mpi
528
529
530
531 }
532
533
534 void EAM_FF::initializeAtoms( int nAtoms, Atom** the_atoms ){
535
536 int i;
537
538 // initialize the atoms
539
540 for( i=0; i<nAtoms; i++ ){
541
542 currentAtomType = headAtomType->find( the_atoms[i]->getType() );
543 if( currentAtomType == NULL ){
544 sprintf( painCave.errMsg,
545 "AtomType error, %s not found in force file.\n",
546 the_atoms[i]->getType() );
547 painCave.isFatal = 1;
548 simError();
549 }
550
551 the_atoms[i]->setMass( currentAtomType->mass );
552 the_atoms[i]->setIdent( currentAtomType->ident );
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 "variable?\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 }