1 |
gezelter |
246 |
/* |
2 |
|
|
* Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
chrisfen |
153 |
* |
4 |
gezelter |
246 |
* The University of Notre Dame grants you ("Licensee") a |
5 |
|
|
* non-exclusive, royalty free, license to use, modify and |
6 |
|
|
* redistribute this software in source and binary code form, provided |
7 |
|
|
* that the following conditions are met: |
8 |
chrisfen |
153 |
* |
9 |
gezelter |
246 |
* 1. Acknowledgement of the program authors must be made in any |
10 |
|
|
* publication of scientific results based in part on use of the |
11 |
|
|
* program. An acceptable form of acknowledgement is citation of |
12 |
|
|
* the article in which the program was described (Matthew |
13 |
|
|
* A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
|
|
* J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
|
|
* Parallel Simulation Engine for Molecular Dynamics," |
16 |
|
|
* J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
|
|
* |
18 |
|
|
* 2. Redistributions of source code must retain the above copyright |
19 |
|
|
* notice, this list of conditions and the following disclaimer. |
20 |
|
|
* |
21 |
|
|
* 3. Redistributions in binary form must reproduce the above copyright |
22 |
|
|
* notice, this list of conditions and the following disclaimer in the |
23 |
|
|
* documentation and/or other materials provided with the |
24 |
|
|
* distribution. |
25 |
|
|
* |
26 |
|
|
* This software is provided "AS IS," without a warranty of any |
27 |
|
|
* kind. All express or implied conditions, representations and |
28 |
|
|
* warranties, including any implied warranty of merchantability, |
29 |
|
|
* fitness for a particular purpose or non-infringement, are hereby |
30 |
|
|
* excluded. The University of Notre Dame and its licensors shall not |
31 |
|
|
* be liable for any damages suffered by licensee as a result of |
32 |
|
|
* using, modifying or distributing the software or its |
33 |
|
|
* derivatives. In no event will the University of Notre Dame or its |
34 |
|
|
* licensors be liable for any lost revenue, profit or data, or for |
35 |
|
|
* direct, indirect, special, consequential, incidental or punitive |
36 |
|
|
* damages, however caused and regardless of the theory of liability, |
37 |
|
|
* arising out of the use of or inability to use software, even if the |
38 |
|
|
* University of Notre Dame has been advised of the possibility of |
39 |
|
|
* such damages. |
40 |
chrisfen |
153 |
*/ |
41 |
gezelter |
246 |
|
42 |
gezelter |
31 |
#include <stdlib.h> |
43 |
|
|
#include <stdio.h> |
44 |
|
|
#include <string.h> |
45 |
chrisfen |
153 |
#include <map> |
46 |
gezelter |
157 |
#include <cmath> |
47 |
chrisfen |
153 |
#include <iostream> |
48 |
gezelter |
31 |
|
49 |
|
|
#ifdef IS_MPI |
50 |
|
|
#include <mpi.h> |
51 |
|
|
#endif //is_mpi |
52 |
|
|
|
53 |
|
|
#include "UseTheForce/ForceFields.hpp" |
54 |
|
|
#include "primitives/SRI.hpp" |
55 |
|
|
#include "utils/simError.h" |
56 |
gezelter |
157 |
#include "utils/StringUtils.hpp" |
57 |
chrisfen |
153 |
#include "io/basic_ifstrstream.hpp" |
58 |
|
|
#include "math/RealSphericalHarmonic.hpp" |
59 |
|
|
#include "math/SquareMatrix3.hpp" |
60 |
gezelter |
157 |
#include "types/ShapeAtomType.hpp" |
61 |
chrisfen |
153 |
#include "UseTheForce/DarkSide/atype_interface.h" |
62 |
gezelter |
120 |
#include "UseTheForce/DarkSide/shapes_interface.h" |
63 |
gezelter |
31 |
|
64 |
|
|
#ifdef IS_MPI |
65 |
|
|
#include "UseTheForce/mpiForceField.h" |
66 |
|
|
#endif // is_mpi |
67 |
|
|
|
68 |
gezelter |
246 |
|
69 |
gezelter |
157 |
using namespace oopse; |
70 |
gezelter |
31 |
|
71 |
|
|
Shapes_FF::~Shapes_FF(){ |
72 |
|
|
|
73 |
|
|
#ifdef IS_MPI |
74 |
|
|
if( worldRank == 0 ){ |
75 |
|
|
#endif // is_mpi |
76 |
|
|
|
77 |
gezelter |
157 |
forceFile.close(); |
78 |
gezelter |
31 |
|
79 |
|
|
#ifdef IS_MPI |
80 |
|
|
} |
81 |
|
|
#endif // is_mpi |
82 |
|
|
} |
83 |
|
|
|
84 |
|
|
|
85 |
|
|
void Shapes_FF::calcRcut( void ){ |
86 |
|
|
|
87 |
|
|
#ifdef IS_MPI |
88 |
chrisfen |
194 |
double tempShapesRcut = bigContact; |
89 |
gezelter |
31 |
MPI_Allreduce( &tempShapesRcut, &shapesRcut, 1, MPI_DOUBLE, MPI_MAX, |
90 |
|
|
MPI_COMM_WORLD); |
91 |
chrisfen |
194 |
#else |
92 |
|
|
shapesRcut = bigContact; |
93 |
gezelter |
31 |
#endif //is_mpi |
94 |
|
|
entry_plug->setDefaultRcut(shapesRcut); |
95 |
|
|
} |
96 |
|
|
|
97 |
|
|
|
98 |
gezelter |
135 |
void Shapes_FF::initForceField(){ |
99 |
|
|
initFortran(0); |
100 |
gezelter |
31 |
} |
101 |
|
|
|
102 |
|
|
|
103 |
gezelter |
157 |
void Shapes_FF::readParams( void ){ |
104 |
gezelter |
31 |
|
105 |
chrisfen |
153 |
char readLine[1024]; |
106 |
gezelter |
157 |
|
107 |
gezelter |
246 |
std::string fileName; |
108 |
|
|
std::string shapeFileName; |
109 |
|
|
std::string tempString; |
110 |
gezelter |
157 |
|
111 |
chrisfen |
153 |
char *nameToken; |
112 |
|
|
char *delim = " ,;\t\n"; |
113 |
|
|
int nTokens, i; |
114 |
|
|
int nContact = 0; |
115 |
|
|
int nRange = 0; |
116 |
|
|
int nStrength = 0; |
117 |
|
|
int myATID; |
118 |
gezelter |
157 |
int isError; |
119 |
gezelter |
246 |
std::string nameString; |
120 |
gezelter |
157 |
AtomType* at; |
121 |
|
|
DirectionalAtomType* dat; |
122 |
|
|
ShapeAtomType* st; |
123 |
|
|
|
124 |
gezelter |
246 |
std::map<string, AtomType*>::iterator iter; |
125 |
gezelter |
31 |
|
126 |
chrisfen |
153 |
// vectors for shape transfer to fortran |
127 |
gezelter |
246 |
std::vector<RealSphericalHarmonic*> tempSHVector; |
128 |
|
|
std::vector<int> contactL; |
129 |
|
|
std::vector<int> contactM; |
130 |
|
|
std::vector<int> contactFunc; |
131 |
|
|
std::vector<double> contactCoeff; |
132 |
|
|
std::vector<int> rangeL; |
133 |
|
|
std::vector<int> rangeM; |
134 |
|
|
std::vector<int> rangeFunc; |
135 |
|
|
std::vector<double> rangeCoeff; |
136 |
|
|
std::vector<int> strengthL; |
137 |
|
|
std::vector<int> strengthM; |
138 |
|
|
std::vector<int> strengthFunc; |
139 |
|
|
std::vector<double> strengthCoeff; |
140 |
chrisfen |
153 |
|
141 |
|
|
// generate the force file name |
142 |
gezelter |
157 |
fileName = "Shapes.frc"; |
143 |
chrisfen |
153 |
|
144 |
|
|
// attempt to open the file in the current directory first. |
145 |
gezelter |
157 |
forceFile.open( fileName.c_str() ); |
146 |
chrisfen |
153 |
|
147 |
gezelter |
157 |
if( forceFile == NULL ){ |
148 |
gezelter |
31 |
|
149 |
gezelter |
157 |
tempString = ffPath; |
150 |
|
|
tempString += "/"; |
151 |
|
|
tempString += fileName; |
152 |
|
|
fileName = tempString; |
153 |
gezelter |
31 |
|
154 |
gezelter |
157 |
forceFile.open( fileName.c_str() ); |
155 |
gezelter |
31 |
|
156 |
gezelter |
157 |
if( forceFile == NULL ){ |
157 |
chrisfen |
153 |
|
158 |
|
|
sprintf( painCave.errMsg, |
159 |
|
|
"Error opening the force field parameter file:\n" |
160 |
|
|
"\t%s\n" |
161 |
|
|
"\tHave you tried setting the FORCE_PARAM_PATH environment " |
162 |
|
|
"variable?\n", |
163 |
gezelter |
157 |
fileName.c_str() ); |
164 |
chrisfen |
153 |
painCave.severity = OOPSE_ERROR; |
165 |
gezelter |
31 |
painCave.isFatal = 1; |
166 |
|
|
simError(); |
167 |
|
|
} |
168 |
chrisfen |
153 |
} |
169 |
|
|
|
170 |
|
|
// read in the shape types. |
171 |
|
|
|
172 |
gezelter |
157 |
findBegin( forceFile, "ShapeTypes" ); |
173 |
chrisfen |
153 |
|
174 |
gezelter |
157 |
while( !forceFile.eof() ){ |
175 |
|
|
forceFile.getline( readLine, sizeof(readLine) ); |
176 |
gezelter |
31 |
|
177 |
chrisfen |
153 |
// toss comment lines |
178 |
|
|
if( readLine[0] != '!' && readLine[0] != '#' ){ |
179 |
|
|
|
180 |
|
|
if (isEndLine(readLine)) break; |
181 |
|
|
|
182 |
gezelter |
157 |
nTokens = countTokens(readLine, " ,;\t"); |
183 |
chrisfen |
153 |
if (nTokens != 0) { |
184 |
|
|
if (nTokens < 2) { |
185 |
|
|
sprintf( painCave.errMsg, |
186 |
gezelter |
157 |
"Not enough data on a ShapeTypes line in file: %s\n", |
187 |
|
|
fileName.c_str() ); |
188 |
chrisfen |
153 |
painCave.severity = OOPSE_ERROR; |
189 |
|
|
painCave.isFatal = 1; |
190 |
|
|
simError(); |
191 |
|
|
} |
192 |
gezelter |
31 |
|
193 |
chrisfen |
153 |
nameToken = strtok( readLine, delim ); |
194 |
|
|
shapeFileName = strtok( NULL, delim ); |
195 |
gezelter |
31 |
|
196 |
chrisfen |
153 |
// strings are not char arrays! |
197 |
|
|
nameString = nameToken; |
198 |
gezelter |
31 |
|
199 |
chrisfen |
153 |
// does this AtomType name already exist in the map? |
200 |
|
|
iter = atomTypeMap.find(nameString); |
201 |
|
|
if (iter == atomTypeMap.end()) { |
202 |
|
|
// no, it doesn't, so we may proceed: |
203 |
|
|
|
204 |
gezelter |
157 |
st = new ShapeAtomType(); |
205 |
gezelter |
31 |
|
206 |
chrisfen |
153 |
st->setName(nameString); |
207 |
chrisfen |
194 |
myATID = atomTypeMap.size() + 1; |
208 |
chrisfen |
153 |
st->setIdent(myATID); |
209 |
|
|
parseShapeFile(shapeFileName, st); |
210 |
|
|
st->complete(); |
211 |
|
|
atomTypeMap.insert(make_pair(nameString, st)); |
212 |
|
|
|
213 |
|
|
} else { |
214 |
gezelter |
246 |
// atomType map already contained this std::string (i.e. it was |
215 |
chrisfen |
153 |
// declared in a previous block, and we just need to add |
216 |
|
|
// the shape-specific information for this AtomType: |
217 |
gezelter |
31 |
|
218 |
gezelter |
157 |
st = (ShapeAtomType*)(iter->second); |
219 |
chrisfen |
153 |
parseShapeFile(shapeFileName, st); |
220 |
|
|
} |
221 |
|
|
} |
222 |
gezelter |
31 |
} |
223 |
|
|
} |
224 |
|
|
|
225 |
chrisfen |
153 |
#ifdef IS_MPI |
226 |
gezelter |
31 |
|
227 |
chrisfen |
153 |
// looks like all the processors have their ShapeType vectors ready... |
228 |
|
|
sprintf( checkPointMsg, |
229 |
|
|
"Shapes_FF shape objects read successfully." ); |
230 |
|
|
MPIcheckPoint(); |
231 |
gezelter |
31 |
|
232 |
chrisfen |
153 |
#endif // is_mpi |
233 |
gezelter |
31 |
|
234 |
chrisfen |
153 |
// pack up and send the necessary info to fortran |
235 |
|
|
for (iter = atomTypeMap.begin(); iter != atomTypeMap.end(); ++iter){ |
236 |
gezelter |
31 |
|
237 |
gezelter |
157 |
at = (AtomType*)(iter->second); |
238 |
gezelter |
31 |
|
239 |
chrisfen |
153 |
if (at->isDirectional()) { |
240 |
gezelter |
31 |
|
241 |
chrisfen |
153 |
dat = (DirectionalAtomType*)at; |
242 |
gezelter |
31 |
|
243 |
chrisfen |
153 |
if (dat->isShape()) { |
244 |
gezelter |
31 |
|
245 |
chrisfen |
153 |
st = (ShapeAtomType*)at; |
246 |
|
|
|
247 |
|
|
contactL.clear(); |
248 |
|
|
contactM.clear(); |
249 |
|
|
contactFunc.clear(); |
250 |
|
|
contactCoeff.clear(); |
251 |
|
|
|
252 |
|
|
tempSHVector = st->getContactFuncs(); |
253 |
|
|
|
254 |
|
|
nContact = tempSHVector.size(); |
255 |
|
|
for (i=0; i<nContact; i++){ |
256 |
gezelter |
157 |
contactL.push_back(tempSHVector[i]->getL()); |
257 |
|
|
contactM.push_back(tempSHVector[i]->getM()); |
258 |
|
|
contactFunc.push_back(tempSHVector[i]->getFunctionType()); |
259 |
|
|
contactCoeff.push_back(tempSHVector[i]->getCoefficient()); |
260 |
chrisfen |
153 |
} |
261 |
|
|
|
262 |
|
|
rangeL.clear(); |
263 |
|
|
rangeM.clear(); |
264 |
|
|
rangeFunc.clear(); |
265 |
|
|
rangeCoeff.clear(); |
266 |
|
|
|
267 |
|
|
tempSHVector = st->getRangeFuncs(); |
268 |
|
|
|
269 |
|
|
nRange = tempSHVector.size(); |
270 |
|
|
for (i=0; i<nRange; i++){ |
271 |
gezelter |
157 |
rangeL.push_back(tempSHVector[i]->getL()); |
272 |
|
|
rangeM.push_back(tempSHVector[i]->getM()); |
273 |
|
|
rangeFunc.push_back(tempSHVector[i]->getFunctionType()); |
274 |
|
|
rangeCoeff.push_back(tempSHVector[i]->getCoefficient()); |
275 |
chrisfen |
153 |
} |
276 |
|
|
|
277 |
|
|
strengthL.clear(); |
278 |
|
|
strengthM.clear(); |
279 |
|
|
strengthFunc.clear(); |
280 |
|
|
strengthCoeff.clear(); |
281 |
|
|
|
282 |
|
|
tempSHVector = st->getStrengthFuncs(); |
283 |
|
|
|
284 |
|
|
nStrength = tempSHVector.size(); |
285 |
|
|
for (i=0; i<nStrength; i++){ |
286 |
gezelter |
157 |
strengthL.push_back(tempSHVector[i]->getL()); |
287 |
|
|
strengthM.push_back(tempSHVector[i]->getM()); |
288 |
|
|
strengthFunc.push_back(tempSHVector[i]->getFunctionType()); |
289 |
|
|
strengthCoeff.push_back(tempSHVector[i]->getCoefficient()); |
290 |
chrisfen |
153 |
} |
291 |
|
|
|
292 |
|
|
isError = 0; |
293 |
|
|
myATID = at->getIdent(); |
294 |
chrisfen |
194 |
|
295 |
gezelter |
157 |
makeShape( &nContact, &contactL[0], &contactM[0], &contactFunc[0], |
296 |
|
|
&contactCoeff[0], |
297 |
|
|
&nRange, &rangeL[0], &rangeM[0], &rangeFunc[0], |
298 |
|
|
&rangeCoeff[0], |
299 |
|
|
&nStrength, &strengthL[0], &strengthM[0], |
300 |
|
|
&strengthFunc[0], &strengthCoeff[0], |
301 |
chrisfen |
153 |
&myATID, |
302 |
|
|
&isError); |
303 |
chrisfen |
194 |
|
304 |
chrisfen |
153 |
if( isError ){ |
305 |
|
|
sprintf( painCave.errMsg, |
306 |
|
|
"Error initializing the \"%s\" shape in fortran\n", |
307 |
gezelter |
157 |
(iter->first).c_str() ); |
308 |
chrisfen |
153 |
painCave.isFatal = 1; |
309 |
|
|
simError(); |
310 |
|
|
} |
311 |
|
|
} |
312 |
gezelter |
31 |
} |
313 |
|
|
} |
314 |
|
|
|
315 |
chrisfen |
194 |
isError = 0; |
316 |
|
|
completeShapeFF(&isError); |
317 |
|
|
if( isError ){ |
318 |
|
|
sprintf( painCave.errMsg, |
319 |
|
|
"Error completing Shape FF in fortran\n"); |
320 |
|
|
painCave.isFatal = 1; |
321 |
|
|
simError(); |
322 |
|
|
} |
323 |
|
|
|
324 |
gezelter |
31 |
#ifdef IS_MPI |
325 |
|
|
sprintf( checkPointMsg, |
326 |
|
|
"Shapes_FF atom structures successfully sent to fortran\n" ); |
327 |
|
|
MPIcheckPoint(); |
328 |
|
|
#endif // is_mpi |
329 |
|
|
|
330 |
|
|
} |
331 |
|
|
|
332 |
gezelter |
157 |
void Shapes_FF::cleanMe( void ){ |
333 |
|
|
|
334 |
|
|
} |
335 |
|
|
|
336 |
|
|
void Shapes_FF::initializeAtoms( int nAtoms, Atom** the_atoms ){ |
337 |
chrisfen |
153 |
|
338 |
|
|
int i,j,k; |
339 |
gezelter |
246 |
std::map<string, AtomType*>::iterator iter; |
340 |
gezelter |
31 |
|
341 |
|
|
// initialize the atoms |
342 |
chrisfen |
153 |
DirectionalAtom* dAtom; |
343 |
|
|
AtomType* at; |
344 |
|
|
DirectionalAtomType* dat; |
345 |
gezelter |
157 |
ShapeAtomType* sat; |
346 |
chrisfen |
194 |
double longCutoff; |
347 |
chrisfen |
153 |
double ji[3]; |
348 |
|
|
double inertialMat[3][3]; |
349 |
|
|
Mat3x3d momInt; |
350 |
gezelter |
246 |
std::string myTypeString; |
351 |
chrisfen |
153 |
|
352 |
chrisfen |
194 |
bigContact = 0.0; |
353 |
|
|
|
354 |
gezelter |
31 |
for( i=0; i<nAtoms; i++ ){ |
355 |
|
|
|
356 |
gezelter |
246 |
myTypeString = the_atoms[i]->getType().c_str(); |
357 |
chrisfen |
153 |
iter = atomTypeMap.find(myTypeString); |
358 |
|
|
|
359 |
|
|
if (iter == atomTypeMap.end()) { |
360 |
gezelter |
31 |
sprintf( painCave.errMsg, |
361 |
|
|
"AtomType error, %s not found in force file.\n", |
362 |
gezelter |
246 |
the_atoms[i]->getType().c_str() ); |
363 |
gezelter |
31 |
painCave.isFatal = 1; |
364 |
|
|
simError(); |
365 |
chrisfen |
153 |
} else { |
366 |
gezelter |
31 |
|
367 |
chrisfen |
153 |
at = (AtomType*)(iter->second); |
368 |
gezelter |
31 |
|
369 |
chrisfen |
153 |
the_atoms[i]->setMass( at->getMass() ); |
370 |
|
|
the_atoms[i]->setIdent( at->getIdent() ); |
371 |
gezelter |
157 |
|
372 |
|
|
if ( at->isShape() ) { |
373 |
|
|
|
374 |
|
|
sat = (ShapeAtomType*)at; |
375 |
chrisfen |
194 |
longCutoff = findCutoffDistance(sat); |
376 |
|
|
if (longCutoff > bigContact) bigContact = longCutoff; |
377 |
|
|
cout << bigContact << " is the cutoff value\n"; |
378 |
|
|
|
379 |
|
|
entry_plug->useShapes = 1; |
380 |
chrisfen |
153 |
} |
381 |
gezelter |
31 |
|
382 |
chrisfen |
153 |
the_atoms[i]->setHasCharge(at->isCharge()); |
383 |
gezelter |
31 |
|
384 |
chrisfen |
153 |
if( at->isDirectional() ){ |
385 |
gezelter |
31 |
|
386 |
chrisfen |
153 |
dat = (DirectionalAtomType*)at; |
387 |
|
|
dAtom = (DirectionalAtom *) the_atoms[i]; |
388 |
gezelter |
179 |
|
389 |
chrisfen |
153 |
momInt = dat->getI(); |
390 |
gezelter |
31 |
|
391 |
chrisfen |
153 |
// zero out the moments of inertia matrix |
392 |
|
|
for( j=0; j<3; j++ ) |
393 |
|
|
for( k=0; k<3; k++ ) |
394 |
|
|
inertialMat[j][k] = momInt(j,k); |
395 |
|
|
dAtom->setI( inertialMat ); |
396 |
gezelter |
31 |
|
397 |
chrisfen |
153 |
ji[0] = 0.0; |
398 |
|
|
ji[1] = 0.0; |
399 |
|
|
ji[2] = 0.0; |
400 |
|
|
dAtom->setJ( ji ); |
401 |
gezelter |
31 |
|
402 |
chrisfen |
153 |
if (dat->isDipole()) { |
403 |
|
|
dAtom->setHasDipole( dat->isDipole() ); |
404 |
|
|
entry_plug->n_dipoles++; |
405 |
|
|
} |
406 |
|
|
} |
407 |
gezelter |
31 |
} |
408 |
|
|
} |
409 |
|
|
} |
410 |
|
|
|
411 |
chrisfen |
153 |
void Shapes_FF::initializeBonds( int nBonds, Bond** BondArray, |
412 |
|
|
bond_pair* the_bonds ){ |
413 |
gezelter |
31 |
|
414 |
chrisfen |
153 |
if( nBonds ){ |
415 |
gezelter |
31 |
sprintf( painCave.errMsg, |
416 |
chrisfen |
153 |
"Shapes_FF does not support bonds.\n" ); |
417 |
gezelter |
31 |
painCave.isFatal = 1; |
418 |
|
|
simError(); |
419 |
|
|
} |
420 |
chrisfen |
153 |
} |
421 |
gezelter |
31 |
|
422 |
chrisfen |
153 |
void Shapes_FF::initializeBends( int nBends, Bend** bendArray, |
423 |
|
|
bend_set* the_bends ){ |
424 |
|
|
|
425 |
|
|
if( nBends ){ |
426 |
gezelter |
31 |
sprintf( painCave.errMsg, |
427 |
chrisfen |
153 |
"Shapes_FF does not support bends.\n" ); |
428 |
gezelter |
31 |
painCave.isFatal = 1; |
429 |
|
|
simError(); |
430 |
|
|
} |
431 |
chrisfen |
153 |
} |
432 |
gezelter |
31 |
|
433 |
chrisfen |
153 |
void Shapes_FF::initializeTorsions( int nTorsions, Torsion** torsionArray, |
434 |
|
|
torsion_set* the_torsions ){ |
435 |
gezelter |
31 |
|
436 |
chrisfen |
153 |
if( nTorsions ){ |
437 |
gezelter |
31 |
sprintf( painCave.errMsg, |
438 |
chrisfen |
153 |
"Shapes_FF does not support torsions.\n" ); |
439 |
gezelter |
31 |
painCave.isFatal = 1; |
440 |
|
|
simError(); |
441 |
|
|
} |
442 |
chrisfen |
153 |
} |
443 |
gezelter |
31 |
|
444 |
gezelter |
157 |
void Shapes_FF::parseShapeFile(string shapeFileName, ShapeAtomType* st){ |
445 |
chrisfen |
153 |
const int MAXLEN = 1024; |
446 |
|
|
char inLine[MAXLEN]; |
447 |
|
|
char *token; |
448 |
|
|
char *delim = " ,;\t\n"; |
449 |
gezelter |
157 |
int nTokens; |
450 |
chrisfen |
153 |
Mat3x3d momInert; |
451 |
gezelter |
157 |
RealSphericalHarmonic* rsh; |
452 |
gezelter |
246 |
std::vector<RealSphericalHarmonic*> functionVector; |
453 |
gezelter |
157 |
ifstrstream shapeFile; |
454 |
gezelter |
246 |
std::string tempString; |
455 |
gezelter |
31 |
|
456 |
gezelter |
157 |
shapeFile.open( shapeFileName.c_str() ); |
457 |
gezelter |
31 |
|
458 |
gezelter |
157 |
if( shapeFile == NULL ){ |
459 |
|
|
|
460 |
|
|
tempString = ffPath; |
461 |
|
|
tempString += "/"; |
462 |
|
|
tempString += shapeFileName; |
463 |
|
|
shapeFileName = tempString; |
464 |
|
|
|
465 |
|
|
shapeFile.open( shapeFileName.c_str() ); |
466 |
|
|
|
467 |
|
|
if( shapeFile == NULL ){ |
468 |
|
|
|
469 |
|
|
sprintf( painCave.errMsg, |
470 |
|
|
"Error opening the shape file:\n" |
471 |
|
|
"\t%s\n" |
472 |
|
|
"\tHave you tried setting the FORCE_PARAM_PATH environment " |
473 |
|
|
"variable?\n", |
474 |
|
|
shapeFileName.c_str() ); |
475 |
|
|
painCave.severity = OOPSE_ERROR; |
476 |
|
|
painCave.isFatal = 1; |
477 |
|
|
simError(); |
478 |
|
|
} |
479 |
|
|
} |
480 |
|
|
|
481 |
|
|
// read in the shape types. |
482 |
|
|
|
483 |
chrisfen |
153 |
// first grab the values in the ShapeInfo section |
484 |
gezelter |
157 |
findBegin( shapeFile, "ShapeInfo"); |
485 |
gezelter |
31 |
|
486 |
chrisfen |
153 |
shapeFile.getline(inLine, MAXLEN); |
487 |
|
|
while( !shapeFile.eof() ) { |
488 |
|
|
// toss comment lines |
489 |
|
|
if( inLine[0] != '!' && inLine[0] != '#' ){ |
490 |
|
|
// end marks section completion |
491 |
|
|
if (isEndLine(inLine)) break; |
492 |
|
|
|
493 |
gezelter |
157 |
nTokens = countTokens(inLine, delim); |
494 |
chrisfen |
153 |
if (nTokens != 0) { |
495 |
|
|
if (nTokens < 5) { |
496 |
|
|
sprintf( painCave.errMsg, |
497 |
gezelter |
157 |
"Not enough data on a ShapeInfo line in file: %s\n", |
498 |
|
|
shapeFileName.c_str() ); |
499 |
chrisfen |
153 |
painCave.severity = OOPSE_ERROR; |
500 |
|
|
painCave.isFatal = 1; |
501 |
|
|
simError(); |
502 |
chrisfen |
194 |
} else { |
503 |
chrisfen |
153 |
token = strtok(inLine, delim); |
504 |
|
|
token = strtok(NULL, delim); |
505 |
|
|
st->setMass(atof(token)); |
506 |
|
|
token = strtok(NULL, delim); |
507 |
|
|
momInert(0,0) = atof(token); |
508 |
|
|
token = strtok(NULL, delim); |
509 |
|
|
momInert(1,1) = atof(token); |
510 |
|
|
token = strtok(NULL, delim); |
511 |
|
|
momInert(2,2) = atof(token); |
512 |
|
|
st->setI(momInert); |
513 |
|
|
} |
514 |
|
|
} |
515 |
gezelter |
31 |
} |
516 |
gezelter |
179 |
shapeFile.getline(inLine, MAXLEN); |
517 |
gezelter |
31 |
} |
518 |
|
|
|
519 |
chrisfen |
153 |
// now grab the contact functions |
520 |
|
|
findBegin(shapeFile, "ContactFunctions"); |
521 |
|
|
functionVector.clear(); |
522 |
gezelter |
31 |
|
523 |
chrisfen |
153 |
shapeFile.getline(inLine, MAXLEN); |
524 |
|
|
while( !shapeFile.eof() ) { |
525 |
|
|
// toss comment lines |
526 |
|
|
if( inLine[0] != '!' && inLine[0] != '#' ){ |
527 |
|
|
// end marks section completion |
528 |
|
|
if (isEndLine(inLine)) break; |
529 |
gezelter |
157 |
nTokens = countTokens(inLine, delim); |
530 |
chrisfen |
153 |
if (nTokens != 0) { |
531 |
|
|
if (nTokens < 4) { |
532 |
|
|
sprintf( painCave.errMsg, |
533 |
gezelter |
157 |
"Not enough data on a ContactFunctions line in file: %s\n", |
534 |
|
|
shapeFileName.c_str() ); |
535 |
chrisfen |
153 |
painCave.severity = OOPSE_ERROR; |
536 |
|
|
painCave.isFatal = 1; |
537 |
|
|
simError(); |
538 |
chrisfen |
194 |
} else { |
539 |
chrisfen |
153 |
// read in a spherical harmonic function |
540 |
|
|
token = strtok(inLine, delim); |
541 |
gezelter |
157 |
rsh = new RealSphericalHarmonic(); |
542 |
|
|
rsh->setL(atoi(token)); |
543 |
chrisfen |
153 |
token = strtok(NULL, delim); |
544 |
gezelter |
157 |
rsh->setM(atoi(token)); |
545 |
chrisfen |
153 |
token = strtok(NULL, delim); |
546 |
|
|
if (!strcasecmp("sin",token)) |
547 |
gezelter |
157 |
rsh->makeSinFunction(); |
548 |
chrisfen |
153 |
else |
549 |
gezelter |
157 |
rsh->makeCosFunction(); |
550 |
chrisfen |
153 |
token = strtok(NULL, delim); |
551 |
gezelter |
157 |
rsh->setCoefficient(atof(token)); |
552 |
chrisfen |
153 |
|
553 |
gezelter |
157 |
functionVector.push_back(rsh); |
554 |
chrisfen |
153 |
} |
555 |
|
|
} |
556 |
gezelter |
31 |
} |
557 |
gezelter |
179 |
shapeFile.getline(inLine, MAXLEN); |
558 |
gezelter |
31 |
} |
559 |
|
|
|
560 |
chrisfen |
153 |
// pass contact functions to ShapeType |
561 |
chrisfen |
194 |
|
562 |
chrisfen |
153 |
st->setContactFuncs(functionVector); |
563 |
gezelter |
31 |
|
564 |
chrisfen |
153 |
// now grab the range functions |
565 |
|
|
findBegin(shapeFile, "RangeFunctions"); |
566 |
|
|
functionVector.clear(); |
567 |
gezelter |
31 |
|
568 |
chrisfen |
153 |
shapeFile.getline(inLine, MAXLEN); |
569 |
|
|
while( !shapeFile.eof() ) { |
570 |
|
|
// toss comment lines |
571 |
|
|
if( inLine[0] != '!' && inLine[0] != '#' ){ |
572 |
|
|
// end marks section completion |
573 |
|
|
if (isEndLine(inLine)) break; |
574 |
|
|
|
575 |
gezelter |
157 |
nTokens = countTokens(inLine, delim); |
576 |
chrisfen |
153 |
if (nTokens != 0) { |
577 |
|
|
if (nTokens < 4) { |
578 |
|
|
sprintf( painCave.errMsg, |
579 |
gezelter |
157 |
"Not enough data on a RangeFunctions line in file: %s\n", |
580 |
|
|
shapeFileName.c_str() ); |
581 |
chrisfen |
153 |
painCave.severity = OOPSE_ERROR; |
582 |
|
|
painCave.isFatal = 1; |
583 |
|
|
simError(); |
584 |
chrisfen |
194 |
} else { |
585 |
chrisfen |
153 |
|
586 |
|
|
// read in a spherical harmonic function |
587 |
|
|
token = strtok(inLine, delim); |
588 |
gezelter |
157 |
|
589 |
|
|
rsh = new RealSphericalHarmonic(); |
590 |
|
|
rsh->setL(atoi(token)); |
591 |
chrisfen |
153 |
token = strtok(NULL, delim); |
592 |
gezelter |
157 |
rsh->setM(atoi(token)); |
593 |
chrisfen |
153 |
token = strtok(NULL, delim); |
594 |
|
|
if (!strcasecmp("sin",token)) |
595 |
gezelter |
157 |
rsh->makeSinFunction(); |
596 |
chrisfen |
153 |
else |
597 |
gezelter |
157 |
rsh->makeCosFunction(); |
598 |
chrisfen |
153 |
token = strtok(NULL, delim); |
599 |
gezelter |
157 |
rsh->setCoefficient(atof(token)); |
600 |
chrisfen |
153 |
|
601 |
gezelter |
157 |
functionVector.push_back(rsh); |
602 |
chrisfen |
153 |
} |
603 |
|
|
} |
604 |
gezelter |
31 |
} |
605 |
gezelter |
179 |
shapeFile.getline(inLine, MAXLEN); |
606 |
chrisfen |
153 |
} |
607 |
gezelter |
31 |
|
608 |
chrisfen |
153 |
// pass range functions to ShapeType |
609 |
|
|
st->setRangeFuncs(functionVector); |
610 |
gezelter |
31 |
|
611 |
chrisfen |
153 |
// finally grab the strength functions |
612 |
|
|
findBegin(shapeFile, "StrengthFunctions"); |
613 |
|
|
functionVector.clear(); |
614 |
gezelter |
31 |
|
615 |
chrisfen |
153 |
shapeFile.getline(inLine, MAXLEN); |
616 |
|
|
while( !shapeFile.eof() ) { |
617 |
|
|
// toss comment lines |
618 |
|
|
if( inLine[0] != '!' && inLine[0] != '#' ){ |
619 |
|
|
// end marks section completion |
620 |
|
|
if (isEndLine(inLine)) break; |
621 |
|
|
|
622 |
gezelter |
157 |
nTokens = countTokens(inLine, delim); |
623 |
chrisfen |
153 |
if (nTokens != 0) { |
624 |
|
|
if (nTokens < 4) { |
625 |
|
|
sprintf( painCave.errMsg, |
626 |
gezelter |
157 |
"Not enough data on a StrengthFunctions line in file: %s\n", |
627 |
|
|
shapeFileName.c_str() ); |
628 |
chrisfen |
153 |
painCave.severity = OOPSE_ERROR; |
629 |
|
|
painCave.isFatal = 1; |
630 |
|
|
simError(); |
631 |
chrisfen |
194 |
} else { |
632 |
chrisfen |
153 |
|
633 |
|
|
// read in a spherical harmonic function |
634 |
|
|
token = strtok(inLine, delim); |
635 |
gezelter |
157 |
rsh = new RealSphericalHarmonic(); |
636 |
|
|
rsh->setL(atoi(token)); |
637 |
chrisfen |
153 |
token = strtok(NULL, delim); |
638 |
gezelter |
157 |
rsh->setM(atoi(token)); |
639 |
chrisfen |
153 |
token = strtok(NULL, delim); |
640 |
|
|
if (!strcasecmp("sin",token)) |
641 |
gezelter |
157 |
rsh->makeSinFunction(); |
642 |
chrisfen |
153 |
else |
643 |
gezelter |
157 |
rsh->makeCosFunction(); |
644 |
chrisfen |
153 |
token = strtok(NULL, delim); |
645 |
gezelter |
157 |
rsh->setCoefficient(atof(token)); |
646 |
chrisfen |
153 |
|
647 |
gezelter |
157 |
functionVector.push_back(rsh); |
648 |
chrisfen |
153 |
} |
649 |
|
|
} |
650 |
gezelter |
31 |
} |
651 |
gezelter |
179 |
shapeFile.getline(inLine, MAXLEN); |
652 |
gezelter |
31 |
} |
653 |
|
|
|
654 |
chrisfen |
153 |
// pass strength functions to ShapeType |
655 |
|
|
st->setStrengthFuncs(functionVector); |
656 |
|
|
|
657 |
|
|
// we're done reading from this file |
658 |
|
|
shapeFile.close(); |
659 |
gezelter |
31 |
} |
660 |
chrisfen |
153 |
|
661 |
gezelter |
157 |
double Shapes_FF::findLargestContactDistance(ShapeAtomType* st) { |
662 |
|
|
int i, j, nSteps; |
663 |
|
|
double theta, thetaStep, thetaMin, costheta; |
664 |
|
|
double phi, phiStep; |
665 |
chrisfen |
194 |
double sigma, bs; |
666 |
|
|
|
667 |
gezelter |
157 |
nSteps = 16; |
668 |
|
|
|
669 |
|
|
thetaStep = M_PI / nSteps; |
670 |
|
|
thetaMin = thetaStep / 2.0; |
671 |
|
|
phiStep = thetaStep * 2.0; |
672 |
chrisfen |
194 |
bs = 0.0; |
673 |
gezelter |
157 |
|
674 |
|
|
for (i = 0; i < nSteps; i++) { |
675 |
|
|
|
676 |
|
|
theta = thetaMin + i * thetaStep; |
677 |
|
|
costheta = cos(theta); |
678 |
|
|
|
679 |
|
|
for (j = 0; j < nSteps; j++) { |
680 |
|
|
|
681 |
|
|
phi = j*phiStep; |
682 |
|
|
|
683 |
|
|
sigma = st->getContactValueAt(costheta, phi); |
684 |
|
|
|
685 |
chrisfen |
194 |
if (sigma > bs) bs = sigma; |
686 |
gezelter |
157 |
} |
687 |
|
|
} |
688 |
|
|
|
689 |
chrisfen |
194 |
return bs; |
690 |
gezelter |
157 |
} |
691 |
chrisfen |
194 |
|
692 |
|
|
|
693 |
|
|
double Shapes_FF::findCutoffDistance(ShapeAtomType* st) { |
694 |
|
|
int i, j, nSteps; |
695 |
|
|
double theta, thetaStep, thetaMin, costheta; |
696 |
|
|
double phi, phiStep; |
697 |
|
|
double sigma, range; |
698 |
|
|
double bigCut, tempCut; |
699 |
|
|
|
700 |
|
|
nSteps = 16; |
701 |
|
|
|
702 |
|
|
thetaStep = M_PI / nSteps; |
703 |
|
|
thetaMin = thetaStep / 2.0; |
704 |
|
|
phiStep = thetaStep * 2.0; |
705 |
|
|
bigCut = 0.0; |
706 |
|
|
|
707 |
|
|
for (i = 0; i < nSteps; i++) { |
708 |
|
|
|
709 |
|
|
theta = thetaMin + i * thetaStep; |
710 |
|
|
costheta = cos(theta); |
711 |
|
|
|
712 |
|
|
for (j = 0; j < nSteps; j++) { |
713 |
|
|
|
714 |
|
|
phi = j*phiStep; |
715 |
|
|
|
716 |
|
|
sigma = st->getContactValueAt(costheta, phi); |
717 |
|
|
range = st->getRangeValueAt(costheta, phi); |
718 |
|
|
|
719 |
chrisfen |
478 |
// cutoff for a shape is taken to be (1.5*rangeVal + contactVal) |
720 |
chrisfen |
213 |
tempCut = 1.5*range + sigma; |
721 |
chrisfen |
194 |
|
722 |
|
|
if (tempCut > bigCut) bigCut = tempCut; |
723 |
|
|
} |
724 |
|
|
} |
725 |
|
|
|
726 |
|
|
return bigCut; |
727 |
|
|
} |