ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/nonbonded/InteractionManager.cpp
Revision: 1535
Committed: Fri Dec 31 18:31:56 2010 UTC (14 years, 4 months ago) by gezelter
File size: 30181 byte(s)
Log Message:
Well, it compiles and builds, but still has a bus error at runtime.

File Contents

# User Rev Content
1 gezelter 1502 /*
2     * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
4     * 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     *
9     * 1. Redistributions of source code must retain the above copyright
10     * notice, this list of conditions and the following disclaimer.
11     *
12     * 2. Redistributions in binary form must reproduce the above copyright
13     * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31     *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40     */
41    
42     #include "nonbonded/InteractionManager.hpp"
43 gezelter 1535 #include "UseTheForce/doForces_interface.h"
44 gezelter 1502
45     namespace OpenMD {
46 gezelter 1535
47     InteractionManager* InteractionManager::_instance = NULL;
48     SimInfo* InteractionManager::info_ = NULL;
49 gezelter 1502 bool InteractionManager::initialized_ = false;
50 gezelter 1535
51     RealType InteractionManager::rCut_ = 0.0;
52     RealType InteractionManager::rSwitch_ = 0.0;
53     RealType InteractionManager::skinThickness_ = 0.0;
54     RealType InteractionManager::listRadius_ = 0.0;
55     CutoffMethod InteractionManager::cutoffMethod_ = SHIFTED_FORCE;
56     SwitchingFunctionType InteractionManager::sft_ = cubic;
57     RealType InteractionManager::vdwScale_[4] = {0.0, 0.0, 0.0, 0.0};
58     RealType InteractionManager::electrostaticScale_[4] = {0.0, 0.0, 0.0, 0.0};
59    
60 gezelter 1502 map<int, AtomType*> InteractionManager::typeMap_;
61     map<pair<AtomType*, AtomType*>, set<NonBondedInteraction*> > InteractionManager::interactions_;
62 gezelter 1528
63     LJ* InteractionManager::lj_ = new LJ();
64     GB* InteractionManager::gb_ = new GB();
65     Sticky* InteractionManager::sticky_ = new Sticky();
66     Morse* InteractionManager::morse_ = new Morse();
67     EAM* InteractionManager::eam_ = new EAM();
68     SC* InteractionManager::sc_ = new SC();
69     Electrostatic* InteractionManager::electrostatic_ = new Electrostatic();
70 gezelter 1532 MAW* InteractionManager::maw_ = new MAW();
71 gezelter 1530 SwitchingFunction* InteractionManager::switcher_ = new SwitchingFunction();
72 gezelter 1502
73     InteractionManager* InteractionManager::Instance() {
74     if (!_instance) {
75     _instance = new InteractionManager();
76     }
77     return _instance;
78     }
79    
80     void InteractionManager::initialize() {
81    
82 gezelter 1535 ForceField* forceField_ = info_->getForceField();
83    
84 gezelter 1502 lj_->setForceField(forceField_);
85     gb_->setForceField(forceField_);
86     sticky_->setForceField(forceField_);
87     eam_->setForceField(forceField_);
88     sc_->setForceField(forceField_);
89     morse_->setForceField(forceField_);
90     electrostatic_->setForceField(forceField_);
91 gezelter 1532 maw_->setForceField(forceField_);
92 gezelter 1502
93 gezelter 1530 ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
94 gezelter 1535
95 gezelter 1530 // Force fields can set options on how to scale van der Waals and electrostatic
96     // interactions for atoms connected via bonds, bends and torsions
97     // in this case the topological distance between atoms is:
98     // 0 = the atom itself
99     // 1 = bonded together
100     // 2 = connected via a bend
101     // 3 = connected via a torsion
102    
103     vdwScale_[0] = 0.0;
104     vdwScale_[1] = fopts.getvdw12scale();
105     vdwScale_[2] = fopts.getvdw13scale();
106     vdwScale_[3] = fopts.getvdw14scale();
107    
108     electrostaticScale_[0] = 0.0;
109     electrostaticScale_[1] = fopts.getelectrostatic12scale();
110     electrostaticScale_[2] = fopts.getelectrostatic13scale();
111 gezelter 1535 electrostaticScale_[3] = fopts.getelectrostatic14scale();
112 gezelter 1530
113 gezelter 1502 ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes();
114     ForceField::AtomTypeContainer::MapTypeIterator i1, i2;
115     AtomType* atype1;
116     AtomType* atype2;
117     pair<AtomType*, AtomType*> key;
118     pair<set<NonBondedInteraction*>::iterator, bool> ret;
119    
120     for (atype1 = atomTypes->beginType(i1); atype1 != NULL;
121     atype1 = atomTypes->nextType(i1)) {
122    
123     // add it to the map:
124     AtomTypeProperties atp = atype1->getATP();
125    
126     pair<map<int,AtomType*>::iterator,bool> ret;
127     ret = typeMap_.insert( pair<int, AtomType*>(atp.ident, atype1) );
128     if (ret.second == false) {
129     sprintf( painCave.errMsg,
130     "InteractionManager already had a previous entry with ident %d\n",
131     atp.ident);
132     painCave.severity = OPENMD_INFO;
133     painCave.isFatal = 0;
134     simError();
135     }
136     }
137    
138     // Now, iterate over all known types and add to the interaction map:
139    
140     map<int, AtomType*>::iterator it1, it2;
141     for (it1 = typeMap_.begin(); it1 != typeMap_.end(); ++it1) {
142     atype1 = (*it1).second;
143    
144     for( it2 = typeMap_.begin(); it2 != typeMap_.end(); ++it2) {
145     atype2 = (*it2).second;
146    
147     bool vdwExplicit = false;
148     bool metExplicit = false;
149     bool hbExplicit = false;
150    
151     key = make_pair(atype1, atype2);
152    
153     if (atype1->isLennardJones() && atype2->isLennardJones()) {
154     interactions_[key].insert(lj_);
155     }
156     if (atype1->isElectrostatic() && atype2->isElectrostatic() ) {
157     interactions_[key].insert(electrostatic_);
158     }
159     if (atype1->isSticky() && atype2->isSticky() ) {
160     interactions_[key].insert(sticky_);
161     }
162     if (atype1->isStickyPower() && atype2->isStickyPower() ) {
163     interactions_[key].insert(sticky_);
164     }
165     if (atype1->isEAM() && atype2->isEAM() ) {
166     interactions_[key].insert(eam_);
167     }
168     if (atype1->isSC() && atype2->isSC() ) {
169     interactions_[key].insert(sc_);
170     }
171     if (atype1->isGayBerne() && atype2->isGayBerne() ) {
172     interactions_[key].insert(gb_);
173     }
174     if ((atype1->isGayBerne() && atype2->isLennardJones())
175     || (atype1->isLennardJones() && atype2->isGayBerne())) {
176     interactions_[key].insert(gb_);
177     }
178    
179     // look for an explicitly-set non-bonded interaction type using the
180     // two atom types.
181     NonBondedInteractionType* nbiType = forceField_->getNonBondedInteractionType(atype1->getName(), atype2->getName());
182 gezelter 1535
183     if (nbiType != NULL) {
184 gezelter 1502
185 gezelter 1535 if (nbiType->isLennardJones()) {
186     // We found an explicit Lennard-Jones interaction.
187     // override all other vdw entries for this pair of atom types:
188     set<NonBondedInteraction*>::iterator it;
189     for (it = interactions_[key].begin();
190     it != interactions_[key].end(); ++it) {
191     InteractionFamily ifam = (*it)->getFamily();
192     if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
193     }
194     interactions_[key].insert(lj_);
195     vdwExplicit = true;
196 gezelter 1502 }
197 gezelter 1535
198     if (nbiType->isMorse()) {
199     if (vdwExplicit) {
200     sprintf( painCave.errMsg,
201     "InteractionManager::initialize found more than one "
202     "explicit \n"
203     "\tvan der Waals interaction for atom types %s - %s\n",
204     atype1->getName().c_str(), atype2->getName().c_str());
205     painCave.severity = OPENMD_ERROR;
206     painCave.isFatal = 1;
207     simError();
208     }
209     // We found an explicit Morse interaction.
210     // override all other vdw entries for this pair of atom types:
211     set<NonBondedInteraction*>::iterator it;
212     for (it = interactions_[key].begin();
213     it != interactions_[key].end(); ++it) {
214     InteractionFamily ifam = (*it)->getFamily();
215     if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
216     }
217     interactions_[key].insert(morse_);
218     vdwExplicit = true;
219 gezelter 1502 }
220 gezelter 1535
221     if (nbiType->isEAM()) {
222     // We found an explicit EAM interaction.
223     // override all other metallic entries for this pair of atom types:
224     set<NonBondedInteraction*>::iterator it;
225     for (it = interactions_[key].begin();
226     it != interactions_[key].end(); ++it) {
227     InteractionFamily ifam = (*it)->getFamily();
228     if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
229     }
230     interactions_[key].insert(eam_);
231     metExplicit = true;
232 gezelter 1502 }
233 gezelter 1535
234     if (nbiType->isSC()) {
235     if (metExplicit) {
236     sprintf( painCave.errMsg,
237     "InteractionManager::initialize found more than one "
238     "explicit\n"
239     "\tmetallic interaction for atom types %s - %s\n",
240     atype1->getName().c_str(), atype2->getName().c_str());
241     painCave.severity = OPENMD_ERROR;
242     painCave.isFatal = 1;
243     simError();
244     }
245     // We found an explicit Sutton-Chen interaction.
246     // override all other metallic entries for this pair of atom types:
247     set<NonBondedInteraction*>::iterator it;
248     for (it = interactions_[key].begin();
249     it != interactions_[key].end(); ++it) {
250     InteractionFamily ifam = (*it)->getFamily();
251     if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it);
252     }
253     interactions_[key].insert(sc_);
254     metExplicit = true;
255 gezelter 1502 }
256 gezelter 1535
257     if (nbiType->isMAW()) {
258     if (vdwExplicit) {
259     sprintf( painCave.errMsg,
260     "InteractionManager::initialize found more than one "
261     "explicit\n"
262     "\tvan der Waals interaction for atom types %s - %s\n",
263     atype1->getName().c_str(), atype2->getName().c_str());
264     painCave.severity = OPENMD_ERROR;
265     painCave.isFatal = 1;
266     simError();
267     }
268     // We found an explicit MAW interaction.
269     // override all other vdw entries for this pair of atom types:
270     set<NonBondedInteraction*>::iterator it;
271     for (it = interactions_[key].begin();
272     it != interactions_[key].end(); ++it) {
273     InteractionFamily ifam = (*it)->getFamily();
274     if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it);
275     }
276     interactions_[key].insert(maw_);
277     vdwExplicit = true;
278     }
279 gezelter 1502 }
280     }
281     }
282    
283 gezelter 1535
284     // make sure every pair of atom types in this simulation has a
285     // non-bonded interaction:
286    
287     set<AtomType*> simTypes = info_->getSimulatedAtomTypes();
288     set<AtomType*>::iterator it, jt;
289     for (it = simTypes.begin(); it != simTypes.end(); ++it) {
290     atype1 = (*it);
291     for (jt = simTypes.begin(); jt != simTypes.end(); ++jt) {
292     atype2 = (*jt);
293 gezelter 1502 key = make_pair(atype1, atype2);
294    
295     if (interactions_[key].size() == 0) {
296     sprintf( painCave.errMsg,
297     "InteractionManager unable to find an appropriate non-bonded\n"
298     "\tinteraction for atom types %s - %s\n",
299     atype1->getName().c_str(), atype2->getName().c_str());
300     painCave.severity = OPENMD_INFO;
301     painCave.isFatal = 1;
302     simError();
303     }
304     }
305     }
306 gezelter 1535
307     setupCutoffs();
308     setupSwitching();
309     setupNeighborlists();
310     notifyFortranSkinThickness(&skinThickness_);
311    
312     int ljsp = cutoffMethod_ == SHIFTED_POTENTIAL ? 1 : 0;
313     int ljsf = cutoffMethod_ == SHIFTED_FORCE ? 1 : 0;
314     notifyFortranCutoffs(&rCut_, &rSwitch_, &ljsp, &ljsf);
315    
316     int isError;
317     initFortranFF(&isError);
318    
319     initialized_ = true;
320     }
321 gezelter 1502
322 gezelter 1535 /**
323     * setupCutoffs
324     *
325     * Sets the values of cutoffRadius and cutoffMethod
326     *
327     * cutoffRadius : realType
328     * If the cutoffRadius was explicitly set, use that value.
329     * If the cutoffRadius was not explicitly set:
330     * Are there electrostatic atoms? Use 12.0 Angstroms.
331     * No electrostatic atoms? Poll the atom types present in the
332     * simulation for suggested cutoff values (e.g. 2.5 * sigma).
333     * Use the maximum suggested value that was found.
334     *
335     * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, SHIFTED_POTENTIAL)
336     * If cutoffMethod was explicitly set, use that choice.
337     * If cutoffMethod was not explicitly set, use SHIFTED_FORCE
338     */
339     void InteractionManager::setupCutoffs() {
340    
341     Globals* simParams_ = info_->getSimParams();
342    
343     if (simParams_->haveCutoffRadius()) {
344     rCut_ = simParams_->getCutoffRadius();
345     } else {
346     if (info_->usesElectrostaticAtoms()) {
347     sprintf(painCave.errMsg,
348     "InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n"
349     "\tOpenMD will use a default value of 12.0 angstroms"
350     "\tfor the cutoffRadius.\n");
351     painCave.isFatal = 0;
352     painCave.severity = OPENMD_INFO;
353     simError();
354     rCut_ = 12.0;
355     } else {
356     RealType thisCut;
357     set<AtomType*>::iterator i;
358     set<AtomType*> atomTypes;
359     atomTypes = info_->getSimulatedAtomTypes();
360     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
361     thisCut = getSuggestedCutoffRadius((*i));
362     rCut_ = max(thisCut, rCut_);
363     }
364     sprintf(painCave.errMsg,
365     "InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n"
366     "\tOpenMD will use %lf angstroms.\n",
367     rCut_);
368     painCave.isFatal = 0;
369     painCave.severity = OPENMD_INFO;
370     simError();
371     }
372     }
373    
374     map<string, CutoffMethod> stringToCutoffMethod;
375     stringToCutoffMethod["HARD"] = HARD;
376     stringToCutoffMethod["SWITCHED"] = SWITCHED;
377     stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;
378     stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
379    
380     if (simParams_->haveCutoffMethod()) {
381     string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
382     map<string, CutoffMethod>::iterator i;
383     i = stringToCutoffMethod.find(cutMeth);
384     if (i == stringToCutoffMethod.end()) {
385     sprintf(painCave.errMsg,
386     "InteractionManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
387     "\tShould be one of: "
388     "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
389     cutMeth.c_str());
390     painCave.isFatal = 1;
391     painCave.severity = OPENMD_ERROR;
392     simError();
393     } else {
394     cutoffMethod_ = i->second;
395     }
396     } else {
397     sprintf(painCave.errMsg,
398     "InteractionManager::setupCutoffs: No value was set for the cutoffMethod.\n"
399     "\tOpenMD will use SHIFTED_FORCE.\n");
400     painCave.isFatal = 0;
401     painCave.severity = OPENMD_INFO;
402     simError();
403     cutoffMethod_ = SHIFTED_FORCE;
404     }
405     }
406    
407    
408     /**
409     * setupSwitching
410     *
411     * Sets the values of switchingRadius and
412     * If the switchingRadius was explicitly set, use that value (but check it)
413     * If the switchingRadius was not explicitly set: use 0.85 * cutoffRadius_
414     */
415     void InteractionManager::setupSwitching() {
416     Globals* simParams_ = info_->getSimParams();
417    
418     if (simParams_->haveSwitchingRadius()) {
419     rSwitch_ = simParams_->getSwitchingRadius();
420     if (rSwitch_ > rCut_) {
421     sprintf(painCave.errMsg,
422     "InteractionManager::setupSwitching: switchingRadius (%f) is larger than cutoffRadius(%f)\n",
423     rSwitch_, rCut_);
424     painCave.isFatal = 1;
425     painCave.severity = OPENMD_ERROR;
426     simError();
427     }
428     } else {
429     rSwitch_ = 0.85 * rCut_;
430     sprintf(painCave.errMsg,
431     "InteractionManager::setupSwitching: No value was set for the switchingRadius.\n"
432     "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
433     "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
434     painCave.isFatal = 0;
435     painCave.severity = OPENMD_WARNING;
436     simError();
437     }
438    
439     if (simParams_->haveSwitchingFunctionType()) {
440     string funcType = simParams_->getSwitchingFunctionType();
441     toUpper(funcType);
442     if (funcType == "CUBIC") {
443     sft_ = cubic;
444     } else {
445     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
446     sft_ = fifth_order_poly;
447     } else {
448     // throw error
449     sprintf( painCave.errMsg,
450     "InteractionManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
451     "\tswitchingFunctionType must be one of: "
452     "\"cubic\" or \"fifth_order_polynomial\".",
453     funcType.c_str() );
454     painCave.isFatal = 1;
455     painCave.severity = OPENMD_ERROR;
456     simError();
457     }
458     }
459     }
460     }
461    
462     /**
463     * setupNeighborlists
464     *
465     * If the skinThickness was explicitly set, use that value (but check it)
466     * If the skinThickness was not explicitly set: use 1.0 angstroms
467     */
468     void InteractionManager::setupNeighborlists() {
469    
470     Globals* simParams_ = info_->getSimParams();
471    
472     if (simParams_->haveSkinThickness()) {
473     skinThickness_ = simParams_->getSkinThickness();
474     } else {
475     skinThickness_ = 1.0;
476     sprintf(painCave.errMsg,
477     "InteractionManager::setupNeighborlists: No value was set for the skinThickness.\n"
478     "\tOpenMD will use a default value of %f Angstroms\n"
479     "\tfor this simulation\n", skinThickness_);
480     painCave.severity = OPENMD_INFO;
481     painCave.isFatal = 0;
482     simError();
483     }
484    
485     listRadius_ = rCut_ + skinThickness_;
486    
487     }
488    
489    
490 gezelter 1505 void InteractionManager::doPrePair(int *atid1, int *atid2, RealType *rij, RealType *rho_i_at_j, RealType *rho_j_at_i){
491 gezelter 1502
492     if (!initialized_) initialize();
493 gezelter 1504
494     DensityData ddat;
495    
496     ddat.atype1 = typeMap_[*atid1];
497     ddat.atype2 = typeMap_[*atid2];
498     ddat.rij = *rij;
499     ddat.rho_i_at_j = *rho_i_at_j;
500     ddat.rho_j_at_i = *rho_j_at_i;
501    
502     pair<AtomType*, AtomType*> key = make_pair(ddat.atype1, ddat.atype2);
503     set<NonBondedInteraction*>::iterator it;
504    
505     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
506     if ((*it)->getFamily() == METALLIC_FAMILY) {
507     dynamic_cast<MetallicInteraction*>(*it)->calcDensity(ddat);
508     }
509     }
510 gezelter 1502
511     return;
512     }
513 gezelter 1504
514 gezelter 1505 void InteractionManager::doPreForce(int *atid, RealType *rho, RealType *frho, RealType *dfrhodrho){
515 gezelter 1502
516 gezelter 1504 if (!initialized_) initialize();
517    
518     FunctionalData fdat;
519 gezelter 1502
520 gezelter 1504 fdat.atype = typeMap_[*atid];
521     fdat.rho = *rho;
522     fdat.frho = *frho;
523     fdat.dfrhodrho = *dfrhodrho;
524    
525     pair<AtomType*, AtomType*> key = make_pair(fdat.atype, fdat.atype);
526     set<NonBondedInteraction*>::iterator it;
527 gezelter 1502
528 gezelter 1504 for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
529     if ((*it)->getFamily() == METALLIC_FAMILY) {
530     dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(fdat);
531     }
532     }
533    
534 gezelter 1502 return;
535     }
536    
537 gezelter 1530 void InteractionManager::doPair(int *atid1, int *atid2, RealType *d, RealType *r, RealType *r2, RealType *rcut, RealType *sw, int *topoDist, RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1, RealType *eFrame2, RealType *A1, RealType *A2, RealType *t1, RealType *t2, RealType *rho1, RealType *rho2, RealType *dfrho1, RealType *dfrho2, RealType *fshift1, RealType *fshift2){
538 gezelter 1504
539 gezelter 1502 if (!initialized_) initialize();
540 gezelter 1504
541 gezelter 1502 InteractionData idat;
542 gezelter 1504
543 gezelter 1502 idat.atype1 = typeMap_[*atid1];
544     idat.atype2 = typeMap_[*atid2];
545     idat.d = Vector3d(d);
546     idat.rij = *r;
547     idat.r2 = *r2;
548     idat.rcut = *rcut;
549     idat.sw = *sw;
550 gezelter 1530 idat.vdwMult = vdwScale_[*topoDist];
551     idat.electroMult = electrostaticScale_[*topoDist];
552 gezelter 1502 idat.pot = *pot;
553     idat.vpair = *vpair;
554     idat.f1 = Vector3d(f1);
555     idat.eFrame1 = Mat3x3d(eFrame1);
556     idat.eFrame2 = Mat3x3d(eFrame2);
557     idat.A1 = RotMat3x3d(A1);
558     idat.A2 = RotMat3x3d(A2);
559     idat.t1 = Vector3d(t1);
560     idat.t2 = Vector3d(t2);
561     idat.rho1 = *rho1;
562     idat.rho2 = *rho2;
563     idat.dfrho1 = *dfrho1;
564     idat.dfrho2 = *dfrho2;
565     idat.fshift1 = *fshift1;
566     idat.fshift2 = *fshift2;
567    
568     pair<AtomType*, AtomType*> key = make_pair(idat.atype1, idat.atype2);
569     set<NonBondedInteraction*>::iterator it;
570    
571     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
572     (*it)->calcForce(idat);
573    
574     f1[0] = idat.f1.x();
575     f1[1] = idat.f1.y();
576     f1[2] = idat.f1.z();
577    
578     t1[0] = idat.t1.x();
579     t1[1] = idat.t1.y();
580     t1[2] = idat.t1.z();
581    
582     t2[0] = idat.t2.x();
583     t2[1] = idat.t2.y();
584     t2[2] = idat.t2.z();
585    
586     return;
587     }
588    
589 gezelter 1505 void InteractionManager::doSkipCorrection(int *atid1, int *atid2, RealType *d, RealType *r, RealType *skippedCharge1, RealType *skippedCharge2, RealType *sw, RealType *electroMult, RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1, RealType *eFrame2, RealType *t1, RealType *t2){
590 gezelter 1502
591     if (!initialized_) initialize();
592 gezelter 1504
593     SkipCorrectionData skdat;
594    
595     skdat.atype1 = typeMap_[*atid1];
596     skdat.atype2 = typeMap_[*atid2];
597     skdat.d = Vector3d(d);
598     skdat.rij = *r;
599     skdat.skippedCharge1 = *skippedCharge1;
600     skdat.skippedCharge2 = *skippedCharge2;
601     skdat.sw = *sw;
602     skdat.electroMult = *electroMult;
603     skdat.pot = *pot;
604     skdat.vpair = *vpair;
605     skdat.f1 = Vector3d(f1);
606     skdat.eFrame1 = Mat3x3d(eFrame1);
607     skdat.eFrame2 = Mat3x3d(eFrame2);
608     skdat.t1 = Vector3d(t1);
609     skdat.t2 = Vector3d(t2);
610 gezelter 1502
611 gezelter 1504 pair<AtomType*, AtomType*> key = make_pair(skdat.atype1, skdat.atype2);
612     set<NonBondedInteraction*>::iterator it;
613    
614     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
615     if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
616     dynamic_cast<ElectrostaticInteraction*>(*it)->calcSkipCorrection(skdat);
617     }
618     }
619 gezelter 1502
620 gezelter 1504 f1[0] = skdat.f1.x();
621     f1[1] = skdat.f1.y();
622     f1[2] = skdat.f1.z();
623 gezelter 1502
624 gezelter 1504 t1[0] = skdat.t1.x();
625     t1[1] = skdat.t1.y();
626     t1[2] = skdat.t1.z();
627    
628     t2[0] = skdat.t2.x();
629     t2[1] = skdat.t2.y();
630     t2[2] = skdat.t2.z();
631 gezelter 1502
632 gezelter 1504 return;
633 gezelter 1502 }
634    
635 gezelter 1505 void InteractionManager::doSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge, RealType *pot, RealType *t){
636 gezelter 1502
637     if (!initialized_) initialize();
638 gezelter 1504
639     SelfCorrectionData scdat;
640    
641     scdat.atype = typeMap_[*atid];
642     scdat.eFrame = Mat3x3d(eFrame);
643     scdat.skippedCharge = *skippedCharge;
644     scdat.pot = *pot;
645     scdat.t = Vector3d(t);
646 gezelter 1502
647 gezelter 1504 pair<AtomType*, AtomType*> key = make_pair(scdat.atype, scdat.atype);
648     set<NonBondedInteraction*>::iterator it;
649 gezelter 1502
650 gezelter 1504 for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){
651     if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) {
652     dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(scdat);
653     }
654     }
655    
656     t[0] = scdat.t.x();
657     t[1] = scdat.t.y();
658     t[2] = scdat.t.z();
659    
660     return;
661 gezelter 1502 }
662    
663 gezelter 1505
664     RealType InteractionManager::getSuggestedCutoffRadius(int *atid) {
665     if (!initialized_) initialize();
666    
667     AtomType* atype = typeMap_[*atid];
668    
669     pair<AtomType*, AtomType*> key = make_pair(atype, atype);
670     set<NonBondedInteraction*>::iterator it;
671     RealType cutoff = 0.0;
672    
673     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
674     cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));
675     return cutoff;
676     }
677    
678 gezelter 1528 RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) {
679     if (!initialized_) initialize();
680    
681     pair<AtomType*, AtomType*> key = make_pair(atype, atype);
682     set<NonBondedInteraction*>::iterator it;
683     RealType cutoff = 0.0;
684    
685     for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it)
686     cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype));
687     return cutoff;
688     }
689    
690 gezelter 1530
691     void InteractionManager::setSwitch(RealType *rIn, RealType *rOut) {
692     switcher_->setSwitch(*rIn, *rOut);
693     }
694    
695     void InteractionManager::getSwitch(RealType *r2, RealType *sw, RealType *dswdr, RealType *r,
696     int *in_switching_region) {
697     bool isr = switcher_->getSwitch(*r2, *sw, *dswdr, *r);
698     *in_switching_region = (int)isr;
699     }
700    
701 gezelter 1502 } //end namespace OpenMD
702    
703     extern "C" {
704    
705     #define fortranDoPrePair FC_FUNC(do_prepair, DO_PREPAIR)
706     #define fortranDoPreForce FC_FUNC(do_preforce, DO_PREFORCE)
707     #define fortranDoPair FC_FUNC(do_pair, DO_PAIR)
708     #define fortranDoSkipCorrection FC_FUNC(do_skip_correction, DO_SKIP_CORRECTION)
709     #define fortranDoSelfCorrection FC_FUNC(do_self_correction, DO_SELF_CORRECTION)
710     #define fortranGetCutoff FC_FUNC(get_cutoff, GET_CUTOFF)
711 gezelter 1530 #define fortranSetSwitch FC_FUNC(set_switch, SET_SWITCH)
712     #define fortranGetSwitch FC_FUNC(get_switch, GET_SWITCH)
713 gezelter 1502
714     void fortranDoPrePair(int *atid1, int *atid2, RealType *rij,
715     RealType *rho_i_at_j, RealType *rho_j_at_i) {
716    
717 gezelter 1505 return OpenMD::InteractionManager::Instance()->doPrePair(atid1, atid2, rij,
718     rho_i_at_j,
719     rho_j_at_i);
720 gezelter 1502 }
721 gezelter 1528 void fortranDoPreForce(int *atid, RealType *rho, RealType *frho,
722 gezelter 1502 RealType *dfrhodrho) {
723    
724 gezelter 1505 return OpenMD::InteractionManager::Instance()->doPreForce(atid, rho, frho,
725     dfrhodrho);
726 gezelter 1502 }
727    
728     void fortranDoPair(int *atid1, int *atid2, RealType *d, RealType *r,
729 gezelter 1530 RealType *r2, RealType *rcut, RealType *sw, int *topoDist,
730     RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1,
731 gezelter 1504 RealType *eFrame2, RealType *A1, RealType *A2,
732     RealType *t1, RealType *t2, RealType *rho1, RealType *rho2,
733     RealType *dfrho1, RealType *dfrho2, RealType *fshift1,
734     RealType *fshift2){
735 gezelter 1502
736 gezelter 1505 return OpenMD::InteractionManager::Instance()->doPair(atid1, atid2, d, r,
737 gezelter 1530 r2, rcut, sw, topoDist,
738 gezelter 1505 pot, vpair, f1,
739     eFrame1, eFrame2,
740     A1, A2, t1, t2, rho1,
741     rho2, dfrho1, dfrho2,
742     fshift1, fshift2);
743 gezelter 1502 }
744 gezelter 1505
745     void fortranDoSkipCorrection(int *atid1, int *atid2, RealType *d,
746     RealType *r, RealType *skippedCharge1,
747     RealType *skippedCharge2, RealType *sw,
748     RealType *electroMult, RealType *pot,
749 gezelter 1502 RealType *vpair, RealType *f1,
750     RealType *eFrame1, RealType *eFrame2,
751     RealType *t1, RealType *t2){
752    
753 gezelter 1505 return OpenMD::InteractionManager::Instance()->doSkipCorrection(atid1,
754     atid2, d,
755     r,
756     skippedCharge1,
757     skippedCharge2,
758     sw, electroMult, pot,
759     vpair, f1, eFrame1,
760     eFrame2, t1, t2);
761 gezelter 1502 }
762 gezelter 1505
763 gezelter 1502 void fortranDoSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge,
764     RealType *pot, RealType *t) {
765    
766 gezelter 1505 return OpenMD::InteractionManager::Instance()->doSelfCorrection(atid,
767     eFrame,
768     skippedCharge,
769     pot, t);
770 gezelter 1502 }
771 gezelter 1505 RealType fortranGetCutoff(int *atid) {
772     return OpenMD::InteractionManager::Instance()->getSuggestedCutoffRadius(atid);
773 gezelter 1502 }
774 gezelter 1530
775     void fortranGetSwitch(RealType *r2, RealType *sw, RealType *dswdr, RealType *r,
776     int *in_switching_region) {
777    
778     return OpenMD::InteractionManager::Instance()->getSwitch(r2, sw, dswdr, r,
779     in_switching_region);
780     }
781    
782     void fortranSetSwitch(RealType *rIn, RealType *rOut) {
783     return OpenMD::InteractionManager::Instance()->setSwitch(rIn, rOut);
784     }
785    
786 gezelter 1502 }

Properties

Name Value
svn:eol-style native
svn:executable *