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