ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/UseTheForce/Shapes_FF.cpp
Revision: 1650
Committed: Tue Oct 26 22:24:52 2004 UTC (19 years, 8 months ago) by gezelter
File size: 15005 byte(s)
Log Message:
forcefield refactoring for shapes

File Contents

# Content
1 /*
2 * Shapes_FF.cpp
3 * oopse
4 *
5 * Created by Chris Fennell on 10/20/04.
6 * Copyright 2004 __MyCompanyName__. All rights reserved.
7 *
8 */
9
10 #include <stdlib.h>
11 #include <stdio.h>
12 #include <string.h>
13 #include <map>
14 #include <cmath>
15 #include <iostream>
16
17 #ifdef IS_MPI
18 #include <mpi.h>
19 #endif //is_mpi
20
21 #include "UseTheForce/ForceFields.hpp"
22 #include "primitives/SRI.hpp"
23 #include "utils/simError.h"
24 #include "utils/StringUtils.hpp"
25 #include "io/basic_ifstrstream.hpp"
26 #include "math/RealSphericalHarmonic.hpp"
27 #include "math/SquareMatrix3.hpp"
28 #include "types/ShapeAtomType.hpp"
29 #include "UseTheForce/DarkSide/atype_interface.h"
30 #include "UseTheForce/DarkSide/shapes_interface.h"
31
32 #ifdef IS_MPI
33 #include "UseTheForce/mpiForceField.h"
34 #endif // is_mpi
35
36 using namespace std;
37 using namespace oopse;
38
39 Shapes_FF::~Shapes_FF(){
40
41 #ifdef IS_MPI
42 if( worldRank == 0 ){
43 #endif // is_mpi
44
45 forceFile.close();
46
47 #ifdef IS_MPI
48 }
49 #endif // is_mpi
50 }
51
52
53 void Shapes_FF::calcRcut( void ){
54
55 #ifdef IS_MPI
56 double tempShapesRcut = shapesRcut;
57 MPI_Allreduce( &tempShapesRcut, &shapesRcut, 1, MPI_DOUBLE, MPI_MAX,
58 MPI_COMM_WORLD);
59 #endif //is_mpi
60 entry_plug->setDefaultRcut(shapesRcut);
61 }
62
63
64 void Shapes_FF::initForceField(){
65 initFortran(0);
66 }
67
68
69 void Shapes_FF::readParams( void ){
70
71 char readLine[1024];
72
73 string fileName;
74 string shapeFileName;
75 string tempString;
76
77 char *nameToken;
78 char *delim = " ,;\t\n";
79 int nTokens, i;
80 int nContact = 0;
81 int nRange = 0;
82 int nStrength = 0;
83 int myATID;
84 int isError;
85 string nameString;
86 AtomType* at;
87 DirectionalAtomType* dat;
88 ShapeAtomType* st;
89
90 map<string, AtomType*>::iterator iter;
91
92 // vectors for shape transfer to fortran
93 vector<RealSphericalHarmonic*> tempSHVector;
94 vector<int> contactL;
95 vector<int> contactM;
96 vector<int> contactFunc;
97 vector<double> contactCoeff;
98 vector<int> rangeL;
99 vector<int> rangeM;
100 vector<int> rangeFunc;
101 vector<double> rangeCoeff;
102 vector<int> strengthL;
103 vector<int> strengthM;
104 vector<int> strengthFunc;
105 vector<double> strengthCoeff;
106
107 // generate the force file name
108 fileName = "Shapes.frc";
109
110 // attempt to open the file in the current directory first.
111 forceFile.open( fileName.c_str() );
112
113 if( forceFile == NULL ){
114
115 tempString = ffPath;
116 tempString += "/";
117 tempString += fileName;
118 fileName = tempString;
119
120 forceFile.open( fileName.c_str() );
121
122 if( forceFile == NULL ){
123
124 sprintf( painCave.errMsg,
125 "Error opening the force field parameter file:\n"
126 "\t%s\n"
127 "\tHave you tried setting the FORCE_PARAM_PATH environment "
128 "variable?\n",
129 fileName.c_str() );
130 painCave.severity = OOPSE_ERROR;
131 painCave.isFatal = 1;
132 simError();
133 }
134 }
135
136 // read in the shape types.
137
138 findBegin( forceFile, "ShapeTypes" );
139
140 while( !forceFile.eof() ){
141 forceFile.getline( readLine, sizeof(readLine) );
142
143 // toss comment lines
144 if( readLine[0] != '!' && readLine[0] != '#' ){
145
146 if (isEndLine(readLine)) break;
147
148 nTokens = countTokens(readLine, " ,;\t");
149 if (nTokens != 0) {
150 if (nTokens < 2) {
151 sprintf( painCave.errMsg,
152 "Not enough data on a ShapeTypes line in file: %s\n",
153 fileName.c_str() );
154 painCave.severity = OOPSE_ERROR;
155 painCave.isFatal = 1;
156 simError();
157 }
158
159 nameToken = strtok( readLine, delim );
160 shapeFileName = strtok( NULL, delim );
161
162 // strings are not char arrays!
163 nameString = nameToken;
164
165 // does this AtomType name already exist in the map?
166 iter = atomTypeMap.find(nameString);
167 if (iter == atomTypeMap.end()) {
168 // no, it doesn't, so we may proceed:
169
170 st = new ShapeAtomType();
171
172 st->setName(nameString);
173 myATID = atomTypeMap.size();
174 st->setIdent(myATID);
175 parseShapeFile(shapeFileName, st);
176 st->complete();
177 atomTypeMap.insert(make_pair(nameString, st));
178
179 } else {
180 // atomType map already contained this string (i.e. it was
181 // declared in a previous block, and we just need to add
182 // the shape-specific information for this AtomType:
183
184 st = (ShapeAtomType*)(iter->second);
185 parseShapeFile(shapeFileName, st);
186 }
187 }
188 }
189 }
190
191 #ifdef IS_MPI
192
193 // looks like all the processors have their ShapeType vectors ready...
194 sprintf( checkPointMsg,
195 "Shapes_FF shape objects read successfully." );
196 MPIcheckPoint();
197
198 #endif // is_mpi
199
200 // pack up and send the necessary info to fortran
201 for (iter = atomTypeMap.begin(); iter != atomTypeMap.end(); ++iter){
202
203 at = (AtomType*)(iter->second);
204
205 if (at->isDirectional()) {
206
207 dat = (DirectionalAtomType*)at;
208
209 if (dat->isShape()) {
210
211 st = (ShapeAtomType*)at;
212
213 contactL.clear();
214 contactM.clear();
215 contactFunc.clear();
216 contactCoeff.clear();
217
218 tempSHVector = st->getContactFuncs();
219
220 nContact = tempSHVector.size();
221 for (i=0; i<nContact; i++){
222 contactL.push_back(tempSHVector[i]->getL());
223 contactM.push_back(tempSHVector[i]->getM());
224 contactFunc.push_back(tempSHVector[i]->getFunctionType());
225 contactCoeff.push_back(tempSHVector[i]->getCoefficient());
226 }
227
228 rangeL.clear();
229 rangeM.clear();
230 rangeFunc.clear();
231 rangeCoeff.clear();
232
233 tempSHVector = st->getRangeFuncs();
234
235 nRange = tempSHVector.size();
236 for (i=0; i<nRange; i++){
237 rangeL.push_back(tempSHVector[i]->getL());
238 rangeM.push_back(tempSHVector[i]->getM());
239 rangeFunc.push_back(tempSHVector[i]->getFunctionType());
240 rangeCoeff.push_back(tempSHVector[i]->getCoefficient());
241 }
242
243 strengthL.clear();
244 strengthM.clear();
245 strengthFunc.clear();
246 strengthCoeff.clear();
247
248 tempSHVector = st->getStrengthFuncs();
249
250 nStrength = tempSHVector.size();
251 for (i=0; i<nStrength; i++){
252 strengthL.push_back(tempSHVector[i]->getL());
253 strengthM.push_back(tempSHVector[i]->getM());
254 strengthFunc.push_back(tempSHVector[i]->getFunctionType());
255 strengthCoeff.push_back(tempSHVector[i]->getCoefficient());
256 }
257
258 isError = 0;
259 myATID = at->getIdent();
260 makeShape( &nContact, &contactL[0], &contactM[0], &contactFunc[0],
261 &contactCoeff[0],
262 &nRange, &rangeL[0], &rangeM[0], &rangeFunc[0],
263 &rangeCoeff[0],
264 &nStrength, &strengthL[0], &strengthM[0],
265 &strengthFunc[0], &strengthCoeff[0],
266 &myATID,
267 &isError);
268 if( isError ){
269 sprintf( painCave.errMsg,
270 "Error initializing the \"%s\" shape in fortran\n",
271 (iter->first).c_str() );
272 painCave.isFatal = 1;
273 simError();
274 }
275 }
276 }
277 }
278
279 #ifdef IS_MPI
280 sprintf( checkPointMsg,
281 "Shapes_FF atom structures successfully sent to fortran\n" );
282 MPIcheckPoint();
283 #endif // is_mpi
284
285 }
286
287 void Shapes_FF::cleanMe( void ){
288
289 }
290
291 void Shapes_FF::initializeAtoms( int nAtoms, Atom** the_atoms ){
292
293 int i,j,k;
294 map<string, AtomType*>::iterator iter;
295
296 // initialize the atoms
297 DirectionalAtom* dAtom;
298 AtomType* at;
299 DirectionalAtomType* dat;
300 ShapeAtomType* sat;
301 double sigma;
302 double ji[3];
303 double inertialMat[3][3];
304 Mat3x3d momInt;
305 string myTypeString;
306
307 for( i=0; i<nAtoms; i++ ){
308
309 myTypeString = the_atoms[i]->getType();
310 iter = atomTypeMap.find(myTypeString);
311
312 if (iter == atomTypeMap.end()) {
313 sprintf( painCave.errMsg,
314 "AtomType error, %s not found in force file.\n",
315 the_atoms[i]->getType() );
316 painCave.isFatal = 1;
317 simError();
318 } else {
319
320 at = (AtomType*)(iter->second);
321
322 the_atoms[i]->setMass( at->getMass() );
323 the_atoms[i]->setIdent( at->getIdent() );
324
325 if ( at->isShape() ) {
326
327 sat = (ShapeAtomType*)at;
328 sigma = findLargestContactDistance(sat);
329 if (sigma > bigSigma) bigSigma = sigma;
330
331 }
332
333 the_atoms[i]->setHasCharge(at->isCharge());
334
335 if( at->isDirectional() ){
336
337 dat = (DirectionalAtomType*)at;
338 dAtom = (DirectionalAtom *) the_atoms[i];
339
340 momInt = dat->getI();
341
342 // zero out the moments of inertia matrix
343 for( j=0; j<3; j++ )
344 for( k=0; k<3; k++ )
345 inertialMat[j][k] = momInt(j,k);
346 dAtom->setI( inertialMat );
347
348 ji[0] = 0.0;
349 ji[1] = 0.0;
350 ji[2] = 0.0;
351 dAtom->setJ( ji );
352
353 if (dat->isDipole()) {
354 dAtom->setHasDipole( dat->isDipole() );
355 entry_plug->n_dipoles++;
356 }
357 }
358 }
359 }
360 }
361
362 void Shapes_FF::initializeBonds( int nBonds, Bond** BondArray,
363 bond_pair* the_bonds ){
364
365 if( nBonds ){
366 sprintf( painCave.errMsg,
367 "Shapes_FF does not support bonds.\n" );
368 painCave.isFatal = 1;
369 simError();
370 }
371 }
372
373 void Shapes_FF::initializeBends( int nBends, Bend** bendArray,
374 bend_set* the_bends ){
375
376 if( nBends ){
377 sprintf( painCave.errMsg,
378 "Shapes_FF does not support bends.\n" );
379 painCave.isFatal = 1;
380 simError();
381 }
382 }
383
384 void Shapes_FF::initializeTorsions( int nTorsions, Torsion** torsionArray,
385 torsion_set* the_torsions ){
386
387 if( nTorsions ){
388 sprintf( painCave.errMsg,
389 "Shapes_FF does not support torsions.\n" );
390 painCave.isFatal = 1;
391 simError();
392 }
393 }
394
395 void Shapes_FF::parseShapeFile(string shapeFileName, ShapeAtomType* st){
396 const int MAXLEN = 1024;
397 char inLine[MAXLEN];
398 char *token;
399 char *delim = " ,;\t\n";
400 int nTokens;
401 Mat3x3d momInert;
402 RealSphericalHarmonic* rsh;
403 vector<RealSphericalHarmonic*> functionVector;
404 ifstrstream shapeFile;
405 string tempString;
406
407 shapeFile.open( shapeFileName.c_str() );
408
409 if( shapeFile == NULL ){
410
411 tempString = ffPath;
412 tempString += "/";
413 tempString += shapeFileName;
414 shapeFileName = tempString;
415
416 shapeFile.open( shapeFileName.c_str() );
417
418 if( shapeFile == NULL ){
419
420 sprintf( painCave.errMsg,
421 "Error opening the shape file:\n"
422 "\t%s\n"
423 "\tHave you tried setting the FORCE_PARAM_PATH environment "
424 "variable?\n",
425 shapeFileName.c_str() );
426 painCave.severity = OOPSE_ERROR;
427 painCave.isFatal = 1;
428 simError();
429 }
430 }
431
432 // read in the shape types.
433
434 // first grab the values in the ShapeInfo section
435 findBegin( shapeFile, "ShapeInfo");
436
437 shapeFile.getline(inLine, MAXLEN);
438 while( !shapeFile.eof() ) {
439 // toss comment lines
440 if( inLine[0] != '!' && inLine[0] != '#' ){
441 // end marks section completion
442 if (isEndLine(inLine)) break;
443
444 nTokens = countTokens(inLine, delim);
445 if (nTokens != 0) {
446 if (nTokens < 5) {
447 sprintf( painCave.errMsg,
448 "Not enough data on a ShapeInfo line in file: %s\n",
449 shapeFileName.c_str() );
450 painCave.severity = OOPSE_ERROR;
451 painCave.isFatal = 1;
452 simError();
453
454 token = strtok(inLine, delim);
455 token = strtok(NULL, delim);
456 st->setMass(atof(token));
457 token = strtok(NULL, delim);
458 momInert(0,0) = atof(token);
459 token = strtok(NULL, delim);
460 momInert(1,1) = atof(token);
461 token = strtok(NULL, delim);
462 momInert(2,2) = atof(token);
463 st->setI(momInert);
464 }
465 }
466 }
467 }
468
469 // now grab the contact functions
470 findBegin(shapeFile, "ContactFunctions");
471 functionVector.clear();
472
473 shapeFile.getline(inLine, MAXLEN);
474 while( !shapeFile.eof() ) {
475 // toss comment lines
476 if( inLine[0] != '!' && inLine[0] != '#' ){
477 // end marks section completion
478 if (isEndLine(inLine)) break;
479
480 nTokens = countTokens(inLine, delim);
481 if (nTokens != 0) {
482 if (nTokens < 4) {
483 sprintf( painCave.errMsg,
484 "Not enough data on a ContactFunctions line in file: %s\n",
485 shapeFileName.c_str() );
486 painCave.severity = OOPSE_ERROR;
487 painCave.isFatal = 1;
488 simError();
489
490 // read in a spherical harmonic function
491 token = strtok(inLine, delim);
492 rsh = new RealSphericalHarmonic();
493 rsh->setL(atoi(token));
494 token = strtok(NULL, delim);
495 rsh->setM(atoi(token));
496 token = strtok(NULL, delim);
497 if (!strcasecmp("sin",token))
498 rsh->makeSinFunction();
499 else
500 rsh->makeCosFunction();
501 token = strtok(NULL, delim);
502 rsh->setCoefficient(atof(token));
503
504 functionVector.push_back(rsh);
505 }
506 }
507 }
508 }
509
510 // pass contact functions to ShapeType
511 st->setContactFuncs(functionVector);
512
513 // now grab the range functions
514 findBegin(shapeFile, "RangeFunctions");
515 functionVector.clear();
516
517 shapeFile.getline(inLine, MAXLEN);
518 while( !shapeFile.eof() ) {
519 // toss comment lines
520 if( inLine[0] != '!' && inLine[0] != '#' ){
521 // end marks section completion
522 if (isEndLine(inLine)) break;
523
524 nTokens = countTokens(inLine, delim);
525 if (nTokens != 0) {
526 if (nTokens < 4) {
527 sprintf( painCave.errMsg,
528 "Not enough data on a RangeFunctions line in file: %s\n",
529 shapeFileName.c_str() );
530 painCave.severity = OOPSE_ERROR;
531 painCave.isFatal = 1;
532 simError();
533
534 // read in a spherical harmonic function
535 token = strtok(inLine, delim);
536
537 rsh = new RealSphericalHarmonic();
538 rsh->setL(atoi(token));
539 token = strtok(NULL, delim);
540 rsh->setM(atoi(token));
541 token = strtok(NULL, delim);
542 if (!strcasecmp("sin",token))
543 rsh->makeSinFunction();
544 else
545 rsh->makeCosFunction();
546 token = strtok(NULL, delim);
547 rsh->setCoefficient(atof(token));
548
549 functionVector.push_back(rsh);
550 }
551 }
552 }
553 }
554
555 // pass range functions to ShapeType
556 st->setRangeFuncs(functionVector);
557
558 // finally grab the strength functions
559 findBegin(shapeFile, "StrengthFunctions");
560 functionVector.clear();
561
562 shapeFile.getline(inLine, MAXLEN);
563 while( !shapeFile.eof() ) {
564 // toss comment lines
565 if( inLine[0] != '!' && inLine[0] != '#' ){
566 // end marks section completion
567 if (isEndLine(inLine)) break;
568
569 nTokens = countTokens(inLine, delim);
570 if (nTokens != 0) {
571 if (nTokens < 4) {
572 sprintf( painCave.errMsg,
573 "Not enough data on a StrengthFunctions line in file: %s\n",
574 shapeFileName.c_str() );
575 painCave.severity = OOPSE_ERROR;
576 painCave.isFatal = 1;
577 simError();
578
579 // read in a spherical harmonic function
580 token = strtok(inLine, delim);
581 rsh = new RealSphericalHarmonic();
582 rsh->setL(atoi(token));
583 token = strtok(NULL, delim);
584 rsh->setM(atoi(token));
585 token = strtok(NULL, delim);
586 if (!strcasecmp("sin",token))
587 rsh->makeSinFunction();
588 else
589 rsh->makeCosFunction();
590 token = strtok(NULL, delim);
591 rsh->setCoefficient(atof(token));
592
593 functionVector.push_back(rsh);
594 }
595 }
596 }
597 }
598
599 // pass strength functions to ShapeType
600 st->setStrengthFuncs(functionVector);
601
602 // we're done reading from this file
603 shapeFile.close();
604 }
605
606 double Shapes_FF::findLargestContactDistance(ShapeAtomType* st) {
607 int i, j, nSteps;
608 double theta, thetaStep, thetaMin, costheta;
609 double phi, phiStep;
610 double sigma, bigSigma;
611
612 nSteps = 16;
613
614 thetaStep = M_PI / nSteps;
615 thetaMin = thetaStep / 2.0;
616 phiStep = thetaStep * 2.0;
617 bigSigma = 0.0;
618
619 for (i = 0; i < nSteps; i++) {
620
621 theta = thetaMin + i * thetaStep;
622 costheta = cos(theta);
623
624 for (j = 0; j < nSteps; j++) {
625
626 phi = j*phiStep;
627
628 sigma = st->getContactValueAt(costheta, phi);
629
630 if (sigma > bigSigma) bigSigma = sigma;
631 }
632 }
633
634 return bigSigma;
635 }