ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/ForceFields.cpp
Revision: 1167
Committed: Wed May 12 16:38:45 2004 UTC (20 years, 11 months ago) by tim
File size: 5225 byte(s)
Log Message:
get rid of rc and massratio from simState, creat cutoff group forevery atom
which does not belong to
cutoff group defined at mdl file

File Contents

# Content
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 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
157
158 isError = 0;
159 entry_plug->lrPot = 0.0;
160
161 for (i=0; i<9; i++) {
162 entry_plug->tau[i] = 0.0;
163 }
164
165
166 #ifdef PROFILE
167 startProfile(pro8);
168 #endif
169
170 fortranForceLoop( pos,
171 rc,
172 A,
173 u_l,
174 frc,
175 trq,
176 entry_plug->tau,
177 &(entry_plug->lrPot),
178 &passedCalcPot,
179 &passedCalcStress,
180 &isError );
181
182
183 #ifdef PROFILE
184 endProfile(pro8);
185 #endif
186
187
188 if( isError ){
189 sprintf( painCave.errMsg,
190 "Error returned from the fortran force calculation.\n" );
191 painCave.isFatal = 1;
192 simError();
193 }
194
195 for(i=0; i<entry_plug->n_mol; i++ ){
196 entry_plug->molecules[i].atoms2rigidBodies();
197 }
198
199
200 #ifdef IS_MPI
201 sprintf( checkPointMsg,
202 "returned from the force calculation.\n" );
203 MPIcheckPoint();
204 #endif // is_mpi
205
206
207 }
208
209
210 void ForceFields::initFortran(int ljMixPolicy, int useReactionField ){
211
212 int isError;
213
214 isError = 0;
215 initFortranFF( &ljMixPolicy, &useReactionField, &isError );
216
217 if(isError){
218 sprintf( painCave.errMsg,
219 "ForceField error: There was an error initializing the forceField in fortran.\n" );
220 painCave.isFatal = 1;
221 simError();
222 }
223
224
225 #ifdef IS_MPI
226 sprintf( checkPointMsg, "ForceField successfully initialized the fortran component list.\n" );
227 MPIcheckPoint();
228 #endif // is_mpi
229
230 }