1 |
#include <iostream> |
2 |
|
3 |
using namespace std; |
4 |
|
5 |
|
6 |
#include <stdlib.h> |
7 |
|
8 |
#ifdef IS_MPI |
9 |
#include <mpi.h> |
10 |
#endif // is_mpi |
11 |
|
12 |
|
13 |
#ifdef PROFILE |
14 |
#include "mdProfile.hpp" |
15 |
#endif |
16 |
|
17 |
#include "simError.h" |
18 |
#include "ForceFields.hpp" |
19 |
#include "Atom.hpp" |
20 |
#include "fortranWrappers.hpp" |
21 |
|
22 |
|
23 |
void ForceFields::calcRcut( void ){ |
24 |
|
25 |
#ifdef IS_MPI |
26 |
double tempBig = bigSigma; |
27 |
MPI_Allreduce( &tempBig, &bigSigma, 1, MPI_DOUBLE, MPI_MAX, |
28 |
MPI_COMM_WORLD); |
29 |
#endif //is_mpi |
30 |
|
31 |
//calc rCut and rList |
32 |
|
33 |
entry_plug->setDefaultRcut( 2.5 * bigSigma ); |
34 |
|
35 |
} |
36 |
|
37 |
void ForceFields::setRcut( double LJrcut ) { |
38 |
|
39 |
#ifdef IS_MPI |
40 |
double tempBig = bigSigma; |
41 |
MPI_Allreduce( &tempBig, &bigSigma, 1, MPI_DOUBLE, MPI_MAX, |
42 |
MPI_COMM_WORLD); |
43 |
#endif //is_mpi |
44 |
|
45 |
if (LJrcut < 2.5 * bigSigma) { |
46 |
sprintf( painCave.errMsg, |
47 |
"Setting Lennard-Jones cutoff radius to %lf.\n" |
48 |
"\tThis value is smaller than %lf, which is\n" |
49 |
"\t2.5 * bigSigma, where bigSigma is the largest\n" |
50 |
"\tvalue of sigma present in the simulation.\n" |
51 |
"\tThis is potentially a problem since the LJ potential may\n" |
52 |
"\tbe appreciable at this distance. If you don't want the\n" |
53 |
"\tsmaller cutoff, change the LJrcut variable.\n", |
54 |
LJrcut, 2.5*bigSigma); |
55 |
painCave.isFatal = 0; |
56 |
simError(); |
57 |
} else { |
58 |
sprintf( painCave.errMsg, |
59 |
"Setting Lennard-Jones cutoff radius to %lf.\n" |
60 |
"\tThis value is larger than %lf, which is\n" |
61 |
"\t2.5 * bigSigma, where bigSigma is the largest\n" |
62 |
"\tvalue of sigma present in the simulation. This should\n" |
63 |
"\tnot be a problem, but could adversely effect performance.\n", |
64 |
LJrcut, 2.5*bigSigma); |
65 |
painCave.isFatal = 0; |
66 |
simError(); |
67 |
} |
68 |
|
69 |
//calc rCut and rList |
70 |
|
71 |
entry_plug->setDefaultRcut( LJrcut ); |
72 |
} |
73 |
|
74 |
void ForceFields::doForces( int calcPot, int calcStress ){ |
75 |
|
76 |
int i, j, isError; |
77 |
double* frc; |
78 |
double* pos; |
79 |
double* trq; |
80 |
double* A; |
81 |
double* u_l; |
82 |
double* rc; |
83 |
double* massRatio; |
84 |
double factor; |
85 |
SimState* config; |
86 |
|
87 |
Molecule* myMols; |
88 |
Atom** myAtoms; |
89 |
int numAtom; |
90 |
int curIndex; |
91 |
double mtot; |
92 |
int numMol; |
93 |
int numCutoffGroups; |
94 |
CutoffGroup* myCutoffGroup; |
95 |
vector<CutoffGroup*>::iterator iterCutoff; |
96 |
double com[3]; |
97 |
vector<double> rcGroup; |
98 |
|
99 |
short int passedCalcPot = (short int)calcPot; |
100 |
short int passedCalcStress = (short int)calcStress; |
101 |
|
102 |
// forces are zeroed here, before any are accumulated. |
103 |
// NOTE: do not rezero the forces in Fortran. |
104 |
|
105 |
for(i=0; i<entry_plug->n_atoms; i++){ |
106 |
entry_plug->atoms[i]->zeroForces(); |
107 |
} |
108 |
|
109 |
#ifdef PROFILE |
110 |
startProfile(pro7); |
111 |
#endif |
112 |
|
113 |
for(i=0; i<entry_plug->n_mol; i++ ){ |
114 |
// CalcForces in molecules takes care of mapping rigid body coordinates |
115 |
// into atomic coordinates |
116 |
entry_plug->molecules[i].calcForces(); |
117 |
} |
118 |
|
119 |
#ifdef PROFILE |
120 |
endProfile( pro7 ); |
121 |
#endif |
122 |
|
123 |
config = entry_plug->getConfiguration(); |
124 |
|
125 |
frc = config->getFrcArray(); |
126 |
pos = config->getPosArray(); |
127 |
trq = config->getTrqArray(); |
128 |
A = config->getAmatArray(); |
129 |
u_l = config->getUlArray(); |
130 |
|
131 |
if(entry_plug->haveCutoffGroups){ |
132 |
myMols = entry_plug->molecules; |
133 |
numMol = entry_plug->n_mol; |
134 |
for(int i = 0; i < numMol; i++){ |
135 |
|
136 |
numCutoffGroups = myMols[i].getNCutoffGroups(); |
137 |
for(myCutoffGroup =myMols[i].beginCutoffGroup(iterCutoff); myCutoffGroup != NULL; |
138 |
myCutoffGroup =myMols[i].nextCutoffGroup(iterCutoff)){ |
139 |
//get center of mass of the cutoff group |
140 |
myCutoffGroup->getCOM(com); |
141 |
|
142 |
rcGroup.push_back(com[0]); |
143 |
rcGroup.push_back(com[1]); |
144 |
rcGroup.push_back(com[2]); |
145 |
|
146 |
}// end for(myCutoffGroup) |
147 |
|
148 |
}//end for(int i = 0) |
149 |
|
150 |
rc = &rcGroup[0]; |
151 |
} |
152 |
else{ |
153 |
// center of mass of the group is the same as position of the atom if cutoff group does not exist |
154 |
rc = pos; |
155 |
} |
156 |
|
157 |
|
158 |
|
159 |
isError = 0; |
160 |
entry_plug->lrPot = 0.0; |
161 |
|
162 |
for (i=0; i<9; i++) { |
163 |
entry_plug->tau[i] = 0.0; |
164 |
} |
165 |
|
166 |
|
167 |
#ifdef PROFILE |
168 |
startProfile(pro8); |
169 |
#endif |
170 |
|
171 |
fortranForceLoop( pos, |
172 |
rc, |
173 |
A, |
174 |
u_l, |
175 |
frc, |
176 |
trq, |
177 |
entry_plug->tau, |
178 |
&(entry_plug->lrPot), |
179 |
&passedCalcPot, |
180 |
&passedCalcStress, |
181 |
&isError ); |
182 |
|
183 |
|
184 |
#ifdef PROFILE |
185 |
endProfile(pro8); |
186 |
#endif |
187 |
|
188 |
|
189 |
if( isError ){ |
190 |
sprintf( painCave.errMsg, |
191 |
"Error returned from the fortran force calculation.\n" ); |
192 |
painCave.isFatal = 1; |
193 |
simError(); |
194 |
} |
195 |
|
196 |
// scale forces if thermodynamic integration is used |
197 |
if (entry_plug->useSolidThermInt || entry_plug->useLiquidThermInt) { |
198 |
factor = pow(entry_plug->thermIntLambda, entry_plug->thermIntK); |
199 |
for (i=0; i < entry_plug->n_atoms; i++) { |
200 |
for (j=0; j< 3; j++) |
201 |
frc[3*i + j] *= factor; |
202 |
if (entry_plug->atoms[i]->isDirectional()) { |
203 |
for (j=0; j< 3; j++) |
204 |
trq[3*i + j] *= factor; |
205 |
} |
206 |
} |
207 |
entry_plug->vRaw = entry_plug->lrPot; |
208 |
entry_plug->lrPot *= factor; |
209 |
} |
210 |
|
211 |
// collect the atomic forces onto rigid bodies |
212 |
for(i=0; i<entry_plug->n_mol; i++ ){ |
213 |
entry_plug->molecules[i].atoms2rigidBodies(); |
214 |
} |
215 |
|
216 |
// do crystal restraint forces for thermodynamic integration |
217 |
if (entry_plug->useSolidThermInt){ |
218 |
entry_plug->lrPot += entry_plug->restraint->Calc_Restraint_Forces(entry_plug->integrableObjects); |
219 |
entry_plug->vHarm = entry_plug->restraint->getVharm(); |
220 |
} |
221 |
|
222 |
|
223 |
#ifdef IS_MPI |
224 |
sprintf( checkPointMsg, |
225 |
"returned from the force calculation.\n" ); |
226 |
MPIcheckPoint(); |
227 |
#endif // is_mpi |
228 |
|
229 |
|
230 |
} |
231 |
|
232 |
|
233 |
void ForceFields::initFortran(int ljMixPolicy, int useReactionField ){ |
234 |
|
235 |
int isError; |
236 |
|
237 |
isError = 0; |
238 |
initFortranFF( &ljMixPolicy, &useReactionField, &isError ); |
239 |
|
240 |
if(isError){ |
241 |
sprintf( painCave.errMsg, |
242 |
"ForceField error: There was an error initializing the forceField in fortran.\n" ); |
243 |
painCave.isFatal = 1; |
244 |
simError(); |
245 |
} |
246 |
|
247 |
|
248 |
#ifdef IS_MPI |
249 |
sprintf( checkPointMsg, "ForceField successfully initialized the fortran component list.\n" ); |
250 |
MPIcheckPoint(); |
251 |
#endif // is_mpi |
252 |
|
253 |
} |
254 |
|
255 |
|
256 |
void ForceFields::initRestraints(){ |
257 |
int i; |
258 |
// store the initial info. |
259 |
// set the omega values to zero |
260 |
for (i=0; i<entry_plug->integrableObjects.size(); i++) |
261 |
entry_plug->integrableObjects[i]->setZangle( 0.0 ); |
262 |
|
263 |
entry_plug->restraint->Store_Init_Info(entry_plug->integrableObjects); |
264 |
|
265 |
} |
266 |
|
267 |
void ForceFields::dumpzAngle(){ |
268 |
|
269 |
// store the initial info. |
270 |
entry_plug->restraint->Write_zAngle_File(entry_plug->integrableObjects); |
271 |
|
272 |
} |