ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/EAM_FF.cpp
Revision: 1617
Committed: Wed Oct 20 20:46:20 2004 UTC (19 years, 8 months ago) by chuckv
File size: 26938 byte(s)
Log Message:
Fortran/C++ interface de-obfuscation project (It is a very long story)

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