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

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