00001
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 #include "time.h"
00018 #include <math.h>
00019 #include "charm++.h"
00020
00021 #include "SimParameters.h"
00022 #include "HomePatch.h"
00023 #include "AtomMap.h"
00024 #include "Node.h"
00025 #include "PatchMap.inl"
00026 #include "main.h"
00027 #include "ProxyMgr.decl.h"
00028 #include "ProxyMgr.h"
00029 #include "Migration.h"
00030 #include "Molecule.h"
00031 #include "PatchMgr.h"
00032 #include "Sequencer.h"
00033 #include "LdbCoordinator.h"
00034 #include "Settle.h"
00035 #include "ReductionMgr.h"
00036 #include "Sync.h"
00037 #include "Random.h"
00038 #include "Priorities.h"
00039 #include "ComputeNonbondedUtil.h"
00040 #include "ComputeGBIS.inl"
00041 #include "Priorities.h"
00042 #include "SortAtoms.h"
00043
00044
00045 #define TINY 1.0e-20;
00046 #define MAXHGS 10
00047 #define MIN_DEBUG_LEVEL 2
00048
00049 #include "Debug.h"
00050
00051 #include <vector>
00052 #include <algorithm>
00053 using namespace std;
00054
00055 typedef int HGArrayInt[MAXHGS];
00056 typedef BigReal HGArrayBigReal[MAXHGS];
00057 typedef zVector HGArrayVector[MAXHGS];
00058 typedef BigReal HGMatrixBigReal[MAXHGS][MAXHGS];
00059 typedef zVector HGMatrixVector[MAXHGS][MAXHGS];
00060
00061 int average(CompAtom *qtilde,const HGArrayVector &q,BigReal *lambda,const int n,const int m, const HGArrayBigReal &imass, const HGArrayBigReal &length2, const HGArrayInt &ial, const HGArrayInt &ibl, const HGArrayVector &refab, const BigReal tolf, const int ntrial);
00062
00063 void mollify(CompAtom *qtilde,const HGArrayVector &q0,const BigReal *lambda, HGArrayVector &force,const int n, const int m, const HGArrayBigReal &imass,const HGArrayInt &ial,const HGArrayInt &ibl,const HGArrayVector &refab);
00064
00065
00066
00067 #if NAMD_SeparateWaters != 0
00068
00069
00070
00071
00072
00073
00074
00075
00076
00077
00078 #define IS_HYDROGEN_GROUP_WATER(hgs, mass) \
00079 ((hgs >= 3) && ((mass >= 14.0) && (mass <= 18.0)))
00080
00081 #endif
00082
00083
00084 HomePatch::HomePatch(PatchID pd, int atomCnt) : Patch(pd)
00085
00086 #if NAMD_SeparateWaters != 0
00087 ,tempAtom()
00088 #endif
00089 {
00090 settle_initialized = 0;
00091
00092
00093 numGBISP1Arrived = 0;
00094 numGBISP2Arrived = 0;
00095 numGBISP3Arrived = 0;
00096 phase1BoxClosedCalled = false;
00097 phase2BoxClosedCalled = false;
00098 phase3BoxClosedCalled = false;
00099
00100 min.x = PatchMap::Object()->min_a(patchID);
00101 min.y = PatchMap::Object()->min_b(patchID);
00102 min.z = PatchMap::Object()->min_c(patchID);
00103 max.x = PatchMap::Object()->max_a(patchID);
00104 max.y = PatchMap::Object()->max_b(patchID);
00105 max.z = PatchMap::Object()->max_c(patchID);
00106 center = 0.5*(min+max);
00107
00108 int aAway = PatchMap::Object()->numaway_a();
00109 if ( PatchMap::Object()->periodic_a() ||
00110 PatchMap::Object()->gridsize_a() > aAway + 1 ) {
00111 aAwayDist = (max.x - min.x) * aAway;
00112 } else {
00113 aAwayDist = Node::Object()->simParameters->patchDimension;
00114 }
00115 int bAway = PatchMap::Object()->numaway_b();
00116 if ( PatchMap::Object()->periodic_b() ||
00117 PatchMap::Object()->gridsize_b() > bAway + 1 ) {
00118 bAwayDist = (max.y - min.y) * bAway;
00119 } else {
00120 bAwayDist = Node::Object()->simParameters->patchDimension;
00121 }
00122 int cAway = PatchMap::Object()->numaway_c();
00123 if ( PatchMap::Object()->periodic_c() ||
00124 PatchMap::Object()->gridsize_c() > cAway + 1 ) {
00125 cAwayDist = (max.z - min.z) * cAway;
00126 } else {
00127 cAwayDist = Node::Object()->simParameters->patchDimension;
00128 }
00129
00130 migrationSuspended = false;
00131 allMigrationIn = false;
00132 marginViolations = 0;
00133 patchMapRead = 0;
00134
00135 inMigration = false;
00136 numMlBuf = 0;
00137 flags.sequence = -1;
00138 flags.maxForceUsed = -1;
00139
00140 numAtoms = atomCnt;
00141 replacementForces = 0;
00142
00143 SimParameters *simParams = Node::Object()->simParameters;
00144 doPairlistCheck_newTolerance =
00145 0.5 * ( simParams->pairlistDist - simParams->cutoff );
00146
00147
00148 numFixedAtoms = 0;
00149
00150
00151
00152
00153
00154
00155 #ifdef NODEAWARE_PROXY_SPANNINGTREE
00156 ptnTree.resize(0);
00157
00158
00159 #else
00160 child = new int[proxySpanDim];
00161 nChild = 0;
00162 #endif
00163
00164 #if CMK_PERSISTENT_COMM && USE_PERSISTENT_TREE
00165 nphs = 0;
00166 localphs = NULL;
00167 isProxyChanged = 0;
00168 #endif
00169
00170
00171
00172 #if NAMD_SeparateWaters != 0
00173
00174
00175 tempAtom.resize(numAtoms);
00176 numWaterAtoms = -1;
00177
00178 #endif
00179
00180 if (simParams->watmodel == WAT_TIP4) init_tip4();
00181 else if (simParams->watmodel == WAT_SWM4) init_swm4();
00182
00183 isNewProxyAdded = 0;
00184 }
00185
00186 HomePatch::HomePatch(PatchID pd, FullAtomList al) : Patch(pd), atom(al)
00187
00188 #if NAMD_SeparateWaters != 0
00189 ,tempAtom()
00190 #endif
00191 {
00192 settle_initialized = 0;
00193
00194 numGBISP1Arrived = 0;
00195 numGBISP2Arrived = 0;
00196 numGBISP3Arrived = 0;
00197 phase1BoxClosedCalled = false;
00198 phase2BoxClosedCalled = false;
00199 phase3BoxClosedCalled = false;
00200
00201 min.x = PatchMap::Object()->min_a(patchID);
00202 min.y = PatchMap::Object()->min_b(patchID);
00203 min.z = PatchMap::Object()->min_c(patchID);
00204 max.x = PatchMap::Object()->max_a(patchID);
00205 max.y = PatchMap::Object()->max_b(patchID);
00206 max.z = PatchMap::Object()->max_c(patchID);
00207 center = 0.5*(min+max);
00208
00209 int aAway = PatchMap::Object()->numaway_a();
00210 if ( PatchMap::Object()->periodic_a() ||
00211 PatchMap::Object()->gridsize_a() > aAway + 1 ) {
00212 aAwayDist = (max.x - min.x) * aAway;
00213 } else {
00214 aAwayDist = Node::Object()->simParameters->patchDimension;
00215 }
00216 int bAway = PatchMap::Object()->numaway_b();
00217 if ( PatchMap::Object()->periodic_b() ||
00218 PatchMap::Object()->gridsize_b() > bAway + 1 ) {
00219 bAwayDist = (max.y - min.y) * bAway;
00220 } else {
00221 bAwayDist = Node::Object()->simParameters->patchDimension;
00222 }
00223 int cAway = PatchMap::Object()->numaway_c();
00224 if ( PatchMap::Object()->periodic_c() ||
00225 PatchMap::Object()->gridsize_c() > cAway + 1 ) {
00226 cAwayDist = (max.z - min.z) * cAway;
00227 } else {
00228 cAwayDist = Node::Object()->simParameters->patchDimension;
00229 }
00230
00231 migrationSuspended = false;
00232 allMigrationIn = false;
00233 marginViolations = 0;
00234 patchMapRead = 0;
00235
00236 inMigration = false;
00237 numMlBuf = 0;
00238 flags.sequence = -1;
00239 flags.maxForceUsed = -1;
00240
00241 numAtoms = atom.size();
00242 replacementForces = 0;
00243
00244 SimParameters *simParams = Node::Object()->simParameters;
00245 doPairlistCheck_newTolerance =
00246 0.5 * ( simParams->pairlistDist - simParams->cutoff );
00247
00248
00249 numFixedAtoms = 0;
00250 if ( simParams->fixedAtomsOn ) {
00251 for ( int i = 0; i < numAtoms; ++i ) {
00252 numFixedAtoms += ( atom[i].atomFixed ? 1 : 0 );
00253 }
00254 }
00255
00256 #ifdef NODEAWARE_PROXY_SPANNINGTREE
00257 ptnTree.resize(0);
00258
00259
00260 #else
00261 child = new int[proxySpanDim];
00262 nChild = 0;
00263 #endif
00264
00265 #if CMK_PERSISTENT_COMM && USE_PERSISTENT_TREE
00266 nphs = 0;
00267 localphs = NULL;
00268 isProxyChanged = 0;
00269 #endif
00270
00271
00272
00273 #if NAMD_SeparateWaters != 0
00274
00275
00276 tempAtom.resize(numAtoms);
00277 numWaterAtoms = -1;
00278
00279
00280 separateAtoms();
00281
00282 #endif
00283
00284
00285 if (simParams->watmodel == WAT_TIP4) init_tip4();
00286 else if (simParams->watmodel == WAT_SWM4) init_swm4();
00287
00288 isNewProxyAdded = 0;
00289 }
00290
00291 void HomePatch::write_tip4_props() {
00292 printf("Writing r_om and r_ohc: %f | %f\n", r_om, r_ohc);
00293 }
00294
00295 void HomePatch::init_tip4() {
00296
00297 Molecule *mol = Node::Object()->molecule;
00298 r_om = mol->r_om;
00299 r_ohc = mol->r_ohc;
00300 }
00301
00302
00303 void ::HomePatch::init_swm4() {
00304
00305 Molecule *mol = Node::Object()->molecule;
00306 r_om = mol->r_om;
00307 r_ohc = mol->r_ohc;
00308 }
00309
00310
00311 void HomePatch::reinitAtoms(FullAtomList al) {
00312 atomMapper->unregisterIDsFullAtom(atom.begin(),atom.end());
00313
00314 atom = al;
00315 numAtoms = atom.size();
00316
00317
00318 #if NAMD_SeparateWaters != 0
00319
00320
00321 numWaterAtoms = -1;
00322
00323
00324 separateAtoms();
00325
00326 #endif
00327 }
00328
00329
00330 void HomePatch::useSequencer(Sequencer *sequencerPtr)
00331 { sequencer=sequencerPtr; }
00332
00333
00334 void HomePatch::runSequencer(void)
00335 { sequencer->run(); }
00336
00337 void HomePatch::readPatchMap() {
00338
00339 PatchMap *p = PatchMap::Object();
00340 PatchID nnPatchID[PatchMap::MaxOneAway];
00341
00342 patchMigrationCounter = numNeighbors
00343 = PatchMap::Object()->oneAwayNeighbors(patchID, nnPatchID);
00344 DebugM( 1, "NumNeighbors for pid " <<patchID<<" is "<< numNeighbors << "\n");
00345 int n;
00346 for (n=0; n<numNeighbors; n++) {
00347 realInfo[n].destNodeID = p->node(realInfo[n].destPatchID = nnPatchID[n]);
00348 DebugM( 1, " nnPatchID=" <<nnPatchID[n]<<" nnNodeID="<< realInfo[n].destNodeID<< "\n");
00349 realInfo[n].mList.resize(0);
00350 }
00351
00352
00353 for (int i=0; i<3; i++)
00354 for (int j=0; j<3; j++)
00355 for (int k=0; k<3; k++)
00356 {
00357 int pid = p->pid(p->index_a(patchID)+i-1,
00358 p->index_b(patchID)+j-1, p->index_c(patchID)+k-1);
00359 if (pid < 0) {
00360 DebugM(5, "ERROR, for patchID " << patchID <<" I got neigh pid = " << pid << "\n");
00361 }
00362 if (pid == patchID && ! (
00363 ( (i-1) && p->periodic_a() ) ||
00364 ( (j-1) && p->periodic_b() ) ||
00365 ( (k-1) && p->periodic_c() ) )) {
00366 mInfo[i][j][k] = NULL;
00367 }
00368 else {
00369
00370
00371 for (n = 0; n<numNeighbors; n++) {
00372 if (pid == realInfo[n].destPatchID) {
00373 mInfo[i][j][k] = &realInfo[n];
00374 break;
00375 }
00376 }
00377 if (n == numNeighbors) {
00378 DebugM(4,"BAD News, I could not find PID " << pid << "\n");
00379 }
00380 }
00381 }
00382
00383 DebugM(1,"Patch("<<patchID<<") # of neighbors = " << numNeighbors << "\n");
00384 }
00385
00386 HomePatch::~HomePatch()
00387 {
00388 atomMapper->unregisterIDsFullAtom(atom.begin(),atom.end());
00389 #ifdef NODEAWARE_PROXY_SPANNINGTREE
00390 ptnTree.resize(0);
00391 #ifdef USE_NODEPATCHMGR
00392 delete [] nodeChildren;
00393 #endif
00394 #endif
00395 delete [] child;
00396 }
00397
00398
00399 void HomePatch::boxClosed(int box) {
00400
00401 if (box == 5) {
00402 phase1BoxClosedCalled = true;
00403 if (!psiSumBox.isOpen() && numGBISP1Arrived == proxy.size()) {
00404 if (flags.doGBIS && flags.doNonbonded) {
00405 sequencer->awaken();
00406 }
00407 } else {
00408
00409 }
00410 } else if (box == 6) {
00411
00412 } else if (box == 7) {
00413
00414 } else if (box == 8) {
00415 phase2BoxClosedCalled = true;
00416
00417
00418 if (!dEdaSumBox.isOpen() && numGBISP2Arrived == proxy.size()) {
00419 if (flags.doGBIS && flags.doNonbonded) {
00420 sequencer->awaken();
00421 }
00422 } else {
00423
00424 }
00425 } else if (box == 9) {
00426
00427 } else if (box == 10) {
00428
00429 } else {
00430
00431 }
00432
00433
00434 if ( ! --boxesOpen ) {
00435 if ( replacementForces ) {
00436 for ( int i = 0; i < numAtoms; ++i ) {
00437 if ( replacementForces[i].replace ) {
00438 for ( int j = 0; j < Results::maxNumForces; ++j ) { f[j][i] = 0; }
00439 f[Results::normal][i] = replacementForces[i].force;
00440 }
00441 }
00442 replacementForces = 0;
00443 }
00444 DebugM(1,patchID << ": " << CthSelf() << " awakening sequencer "
00445 << sequencer->thread << "(" << patchID << ") @" << CmiTimer() << "\n");
00446
00447
00448 phase3BoxClosedCalled = true;
00449 if (flags.doGBIS) {
00450 if (flags.doNonbonded) {
00451 sequencer->awaken();
00452 } else {
00453 if (numGBISP1Arrived == proxy.size() &&
00454 numGBISP2Arrived == proxy.size() &&
00455 numGBISP3Arrived == proxy.size()) {
00456 sequencer->awaken();
00457 }
00458 }
00459 } else {
00460 sequencer->awaken();
00461 }
00462 } else {
00463 DebugM(1,patchID << ": " << boxesOpen << " boxes left to close.\n");
00464 }
00465 }
00466
00467 void HomePatch::registerProxy(RegisterProxyMsg *msg) {
00468 DebugM(4, "registerProxy("<<patchID<<") - adding node " <<msg->node<<"\n");
00469 proxy.add(msg->node);
00470 forceBox.clientAdd();
00471
00472 isNewProxyAdded = 1;
00473 #if CMK_PERSISTENT_COMM && USE_PERSISTENT_TREE
00474 isProxyChanged = 1;
00475 #endif
00476
00477 Random((patchID + 37) * 137).reorder(proxy.begin(),proxy.size());
00478 delete msg;
00479 }
00480
00481 void HomePatch::unregisterProxy(UnregisterProxyMsg *msg) {
00482 #if CMK_PERSISTENT_COMM && USE_PERSISTENT_TREE
00483 isProxyChanged = 1;
00484 #endif
00485 int n = msg->node;
00486 NodeID *pe = proxy.begin();
00487 for ( ; *pe != n; ++pe );
00488 forceBox.clientRemove();
00489 proxy.del(pe - proxy.begin());
00490 delete msg;
00491 }
00492
00493 #if USE_TOPOMAP && USE_SPANNING_TREE
00494
00495 int HomePatch::findSubroots(int dim, int* subroots, int psize, int* pidscopy){
00496 int nChild = 0;
00497 int cones[6][proxySpanDim*proxySpanDim+proxySpanDim];
00498 int conesizes[6] = {0,0,0,0,0,0};
00499 int conecounters[6] = {0,0,0,0,0,0};
00500 int childcounter = 0;
00501 nChild = (psize>proxySpanDim)?proxySpanDim:psize;
00502 TopoManager tmgr;
00503 for(int i=0;i<psize;i++){
00504 int cone = tmgr.getConeNumberForRank(pidscopy[i]);
00505 cones[cone][conesizes[cone]++] = pidscopy[i];
00506 }
00507
00508 while(childcounter<nChild){
00509 for(int i=0;i<6;i++){
00510 if(conecounters[i]<conesizes[i]){
00511 subroots[childcounter++] = cones[i][conecounters[i]++];
00512 }
00513 }
00514 }
00515 for(int i=nChild;i<proxySpanDim;i++)
00516 subroots[i] = -1;
00517 return nChild;
00518 }
00519 #endif // USE_TOPOMAP
00520
00521 static int compDistance(const void *a, const void *b)
00522 {
00523 int d1 = abs(*(int *)a - CkMyPe());
00524 int d2 = abs(*(int *)b - CkMyPe());
00525 if (d1 < d2)
00526 return -1;
00527 else if (d1 == d2)
00528 return 0;
00529 else
00530 return 1;
00531 }
00532
00533 void HomePatch::sendProxies()
00534 {
00535 #if USE_NODEPATCHMGR
00536 CProxy_NodeProxyMgr pm(CkpvAccess(BOCclass_group).nodeProxyMgr);
00537 NodeProxyMgr *npm = pm[CkMyNode()].ckLocalBranch();
00538 npm->sendProxyList(patchID, proxy.begin(), proxy.size());
00539 #else
00540 ProxyMgr::Object()->sendProxies(patchID, proxy.begin(), proxy.size());
00541 #endif
00542 }
00543
00544 #ifdef NODEAWARE_PROXY_SPANNINGTREE
00545 void HomePatch::buildNodeAwareSpanningTree(void){
00546 #if defined(PROCTRACE_DEBUG) && defined(NAST_DEBUG)
00547 DebugFileTrace *dft = DebugFileTrace::Object();
00548 dft->openTrace();
00549 dft->writeTrace("HomePatch[%d] has %d proxy on proc[%d] node[%d]\n", patchID, proxy.size(), CkMyPe(), CkMyNode());
00550 dft->writeTrace("Proxies are: ");
00551 for(int i=0; i<proxy.size(); i++) dft->writeTrace("%d(%d), ", proxy[i], CkNodeOf(proxy[i]));
00552 dft->writeTrace("\n");
00553 dft->closeTrace();
00554 #endif
00555
00556
00557 if(! proxy.size()) {
00558
00559
00560
00561 return;
00562 }
00563 ProxyMgr::buildSinglePatchNodeAwareSpanningTree(patchID, proxy, ptnTree);
00564
00565
00566
00567 setupChildrenFromProxySpanningTree();
00568
00569 sendNodeAwareSpanningTree();
00570 }
00571
00572 void HomePatch::setupChildrenFromProxySpanningTree(){
00573 #if CMK_PERSISTENT_COMM && USE_PERSISTENT_TREE
00574 isProxyChanged = 1;
00575 #endif
00576 if(ptnTree.size()==0) {
00577 nChild = 0;
00578 delete [] child;
00579 child = NULL;
00580 #ifdef USE_NODEPATCHMGR
00581 numNodeChild = 0;
00582 delete [] nodeChildren;
00583 nodeChildren = NULL;
00584 #endif
00585 return;
00586 }
00587 proxyTreeNode *rootnode = &ptnTree.item(0);
00588 CmiAssert(rootnode->peIDs[0] == CkMyPe());
00589
00590
00591
00592 int internalChild = rootnode->numPes-1;
00593 int externalChild = ptnTree.size()-1;
00594 externalChild = (proxySpanDim>externalChild)?externalChild:proxySpanDim;
00595 int internalSlots = proxySpanDim-externalChild;
00596 if(internalChild>0){
00597 if(internalSlots==0) {
00598
00599 internalChild = 1;
00600 }else{
00601 internalChild = (internalSlots>internalChild)?internalChild:internalSlots;
00602 }
00603 }
00604
00605 nChild = externalChild+internalChild;
00606 CmiAssert(nChild>0);
00607
00608
00609 delete [] child;
00610 child = new int[nChild];
00611
00612 for(int i=0; i<externalChild; i++) {
00613 child[i] = ptnTree.item(i+1).peIDs[0];
00614 }
00615 for(int i=externalChild, j=1; i<nChild; i++, j++) {
00616 child[i] = rootnode->peIDs[j];
00617 }
00618
00619 #ifdef USE_NODEPATCHMGR
00620
00621
00622 CProxy_NodeProxyMgr pm(CkpvAccess(BOCclass_group).nodeProxyMgr);
00623 NodeProxyMgr *npm = pm[CkMyNode()].ckLocalBranch();
00624 if(rootnode->numPes==1){
00625 npm->registerPatch(patchID, 0, NULL);
00626 }
00627 else{
00628 npm->registerPatch(patchID, rootnode->numPes-1, &rootnode->peIDs[1]);
00629 }
00630
00631
00632 numNodeChild = externalChild;
00633 if(internalChild) numNodeChild++;
00634 delete [] nodeChildren;
00635 nodeChildren = new int[numNodeChild];
00636 for(int i=0; i<externalChild; i++) {
00637 nodeChildren[i] = ptnTree.item(i+1).nodeID;
00638 }
00639
00640
00641 if(internalChild)
00642 nodeChildren[numNodeChild-1] = rootnode->nodeID;
00643 #endif
00644
00645 #if defined(PROCTRACE_DEBUG) && defined(NAST_DEBUG)
00646 DebugFileTrace *dft = DebugFileTrace::Object();
00647 dft->openTrace();
00648 dft->writeTrace("HomePatch[%d] has %d children: ", patchID, nChild);
00649 for(int i=0; i<nChild; i++)
00650 dft->writeTrace("%d ", child[i]);
00651 dft->writeTrace("\n");
00652 dft->closeTrace();
00653 #endif
00654 }
00655 #endif
00656
00657 #ifdef NODEAWARE_PROXY_SPANNINGTREE
00658
00659 void HomePatch::recvNodeAwareSpanningTree(ProxyNodeAwareSpanningTreeMsg *msg){
00660
00661 int treesize = msg->numNodesWithProxies;
00662 ptnTree.resize(treesize);
00663 int *pAllPes = msg->allPes;
00664 for(int i=0; i<treesize; i++) {
00665 proxyTreeNode *oneNode = &ptnTree.item(i);
00666 delete [] oneNode->peIDs;
00667 oneNode->numPes = msg->numPesOfNode[i];
00668 oneNode->nodeID = CkNodeOf(*pAllPes);
00669 oneNode->peIDs = new int[oneNode->numPes];
00670 for(int j=0; j<oneNode->numPes; j++) {
00671 oneNode->peIDs[j] = *pAllPes;
00672 pAllPes++;
00673 }
00674 }
00675
00676 setupChildrenFromProxySpanningTree();
00677
00678 sendNodeAwareSpanningTree();
00679 }
00680
00681 void HomePatch::sendNodeAwareSpanningTree(){
00682 if(ptnTree.size()==0) return;
00683 ProxyNodeAwareSpanningTreeMsg *msg =
00684 ProxyNodeAwareSpanningTreeMsg::getANewMsg(patchID, CkMyPe(), ptnTree.begin(), ptnTree.size());
00685
00686 #if defined(PROCTRACE_DEBUG) && defined(NAST_DEBUG)
00687 msg->printOut("HP::sendST");
00688 #endif
00689
00690 CmiAssert(CkMyPe() == msg->allPes[0]);
00691 ProxyMgr::Object()->sendNodeAwareSpanningTree(msg);
00692
00693 }
00694 #else
00695 void HomePatch::recvNodeAwareSpanningTree(ProxyNodeAwareSpanningTreeMsg *msg){}
00696 void HomePatch::sendNodeAwareSpanningTree(){}
00697 #endif
00698
00699 #ifndef NODEAWARE_PROXY_SPANNINGTREE
00700
00701 void HomePatch::recvSpanningTree(int *t, int n)
00702 {
00703 int i;
00704 nChild=0;
00705 tree.resize(n);
00706 for (i=0; i<n; i++) {
00707 tree[i] = t[i];
00708 }
00709
00710 for (i=1; i<=proxySpanDim; i++) {
00711 if (tree.size() <= i) break;
00712 child[i-1] = tree[i];
00713 nChild++;
00714 }
00715
00716 #if CMK_PERSISTENT_COMM && USE_PERSISTENT_TREE
00717 isProxyChanged = 1;
00718 #endif
00719
00720
00721 sendSpanningTree();
00722 }
00723
00724 void HomePatch::sendSpanningTree()
00725 {
00726 if (tree.size() == 0) return;
00727 ProxySpanningTreeMsg *msg = new ProxySpanningTreeMsg;
00728 msg->patch = patchID;
00729 msg->node = CkMyPe();
00730 msg->tree = tree;
00731 ProxyMgr::Object()->sendSpanningTree(msg);
00732 }
00733 #else
00734 void HomePatch::recvSpanningTree(int *t, int n){}
00735 void HomePatch::sendSpanningTree(){}
00736 #endif
00737
00738 #ifndef NODEAWARE_PROXY_SPANNINGTREE
00739 void HomePatch::buildSpanningTree(void)
00740 {
00741 nChild = 0;
00742 int psize = proxy.size();
00743 if (psize == 0) return;
00744 NodeIDList oldtree = tree;
00745 int oldsize = oldtree.size();
00746 tree.resize(psize + 1);
00747 tree.setall(-1);
00748 tree[0] = CkMyPe();
00749 int s=1, e=psize+1;
00750 NodeIDList::iterator pli;
00751 int patchNodesLast =
00752 ( PatchMap::Object()->numNodesWithPatches() < ( 0.7 * CkNumPes() ) );
00753 int nNonPatch = 0;
00754 #if 1
00755
00756 for ( pli = proxy.begin(); pli != proxy.end(); ++pli )
00757 {
00758 int oldindex = oldtree.find(*pli);
00759 if (oldindex != -1 && oldindex < psize) {
00760 tree[oldindex] = *pli;
00761 }
00762 }
00763 s=1; e=psize;
00764 for ( pli = proxy.begin(); pli != proxy.end(); ++pli )
00765 {
00766 if (tree.find(*pli) != -1) continue;
00767 if ( patchNodesLast && PatchMap::Object()->numPatchesOnNode(*pli) ) {
00768 while (tree[e] != -1) { e--; if (e==-1) e = psize; }
00769 tree[e] = *pli;
00770 } else {
00771 while (tree[s] != -1) { s++; if (s==psize+1) s = 1; }
00772 tree[s] = *pli;
00773 nNonPatch++;
00774 }
00775 }
00776 #if 1
00777 if (oldsize==0 && nNonPatch) {
00778
00779 qsort(&tree[1], nNonPatch, sizeof(int), compDistance);
00780 }
00781 #endif
00782
00783
00784
00785 #if USE_TOPOMAP && USE_SPANNING_TREE
00786
00787
00788 int *treecopy = new int [psize];
00789 int subroots[proxySpanDim];
00790 int subsizes[proxySpanDim];
00791 int subtrees[proxySpanDim][proxySpanDim];
00792 int idxes[proxySpanDim];
00793 int i = 0;
00794
00795 for(i=0;i<proxySpanDim;i++){
00796 subsizes[i] = 0;
00797 idxes[i] = i;
00798 }
00799
00800 for(i=0;i<psize;i++){
00801 treecopy[i] = tree[i+1];
00802 }
00803
00804 TopoManager tmgr;
00805 tmgr.sortRanksByHops(treecopy,nNonPatch);
00806 tmgr.sortRanksByHops(treecopy+nNonPatch,
00807 psize-nNonPatch);
00808
00809
00810 nChild = findSubroots(proxySpanDim,subroots,psize,treecopy);
00811 delete [] treecopy;
00812
00813 for(int i=1;i<psize+1;i++){
00814 int isSubroot=0;
00815 for(int j=0;j<nChild;j++)
00816 if(tree[i]==subroots[j]){
00817 isSubroot=1;
00818 break;
00819 }
00820 if(isSubroot) continue;
00821
00822 int bAdded = 0;
00823 tmgr.sortIndexByHops(tree[i], subroots,
00824 idxes, proxySpanDim);
00825 for(int j=0;j<proxySpanDim;j++){
00826 if(subsizes[idxes[j]]<proxySpanDim){
00827 subtrees[idxes[j]][(subsizes[idxes[j]])++] = tree[i];
00828 bAdded = 1;
00829 break;
00830 }
00831 }
00832 if( psize > proxySpanDim && ! bAdded ) {
00833 NAMD_bug("HomePatch BGL Spanning Tree error: Couldn't find subtree for leaf\n");
00834 }
00835 }
00836
00837 #else
00838
00839 for (int i=1; i<=proxySpanDim; i++) {
00840 if (tree.size() <= i) break;
00841 child[i-1] = tree[i];
00842 nChild++;
00843 }
00844 #endif
00845 #endif
00846
00847 #if 0
00848
00849 CkPrintf("[%d] Spanning tree for %d with %d children %d nNonPatch %d\n", CkMyPe(), patchID, psize, nNonPatch);
00850 for (int i=0; i<psize+1; i++) {
00851 CkPrintf("%d ", tree[i]);
00852 }
00853 CkPrintf("\n");
00854 #endif
00855
00856 sendSpanningTree();
00857 }
00858 #endif
00859
00860
00861 void HomePatch::receiveResults(ProxyResultVarsizeMsg *msg) {
00862
00863 numGBISP3Arrived++;
00864 DebugM(4, "patchID("<<patchID<<") receiveRes() nodeID("<<msg->node<<")\n");
00865 int n = msg->node;
00866 Results *r = forceBox.clientOpen();
00867
00868 char *iszeroPtr = msg->isZero;
00869 Force *msgFPtr = msg->forceArr;
00870
00871 for ( int k = 0; k < Results::maxNumForces; ++k )
00872 {
00873 Force *rfPtr = r->f[k];
00874 for(int i=0; i<msg->flLen[k]; i++, rfPtr++, iszeroPtr++) {
00875 if((*iszeroPtr)!=1) {
00876 *rfPtr += *msgFPtr;
00877 msgFPtr++;
00878 }
00879 }
00880 }
00881 forceBox.clientClose();
00882 delete msg;
00883 }
00884
00885 void HomePatch::receiveResults(ProxyResultMsg *msg) {
00886 numGBISP3Arrived++;
00887 DebugM(4, "patchID("<<patchID<<") receiveRes() nodeID("<<msg->node<<")\n");
00888 int n = msg->node;
00889 Results *r = forceBox.clientOpen();
00890 for ( int k = 0; k < Results::maxNumForces; ++k )
00891 {
00892 Force *f = r->f[k];
00893 register ForceList::iterator f_i, f_e;
00894 f_i = msg->forceList[k].begin();
00895 f_e = msg->forceList[k].end();
00896 for ( ; f_i != f_e; ++f_i, ++f ) *f += *f_i;
00897 }
00898 forceBox.clientClose();
00899 delete msg;
00900 }
00901
00902 void HomePatch::receiveResults(ProxyCombinedResultRawMsg* msg)
00903 {
00904 numGBISP3Arrived++;
00905 DebugM(4, "patchID("<<patchID<<") receiveRes() #nodes("<<msg->nodeSize<<")\n");
00906 Results *r = forceBox.clientOpen(msg->nodeSize);
00907 register char* isNonZero = msg->isForceNonZero;
00908 register Force* f_i = msg->forceArr;
00909 for ( int k = 0; k < Results::maxNumForces; ++k )
00910 {
00911 Force *f = r->f[k];
00912 int nf = msg->flLen[k];
00913 #ifdef ARCH_POWERPC
00914 #pragma disjoint (*f_i, *f)
00915 #endif
00916 for (int count = 0; count < nf; count++) {
00917 if(*isNonZero){
00918 f[count].x += f_i->x;
00919 f[count].y += f_i->y;
00920 f[count].z += f_i->z;
00921 f_i++;
00922 }
00923 isNonZero++;
00924 }
00925 }
00926 forceBox.clientClose(msg->nodeSize);
00927
00928 delete msg;
00929 }
00930
00931
00932 void HomePatch::positionsReady(int doMigration)
00933 {
00934 SimParameters *simParams = Node::Object()->simParameters;
00935
00936 flags.sequence++;
00937
00938 if (!patchMapRead) { readPatchMap(); }
00939
00940 if (numNeighbors && ! simParams->staticAtomAssignment) {
00941 if (doMigration) {
00942 doAtomMigration();
00943 } else {
00944 doMarginCheck();
00945 }
00946 }
00947
00948 #ifdef NAMD_CUDA
00949 if ( doMigration ) {
00950 int n = numAtoms;
00951 FullAtom *a_i = atom.begin();
00952 int *ao = new int[n];
00953 int nfree;
00954 if ( simParams->fixedAtomsOn && ! simParams->fixedAtomsForces ) {
00955 int k = 0;
00956 int k2 = n;
00957 for ( int j=0; j<n; ++j ) {
00958
00959 if ( a_i[j].atomFixed ) ao[--k2] = j;
00960 else ao[k++] = j;
00961 }
00962 nfree = k;
00963 } else {
00964 nfree = n;
00965 for ( int j=0; j<n; ++j ) {
00966 ao[j] = j;
00967 }
00968 }
00969
00970 sortAtomsForCUDA(ao,a_i,nfree,n);
00971
00972 for ( int i=0; i<n; ++i ) {
00973 a_i[i].sortOrder = ao[i];
00974 }
00975 delete [] ao;
00976 }
00977
00978 {
00979 const double charge_scaling = sqrt(COULOMB * ComputeNonbondedUtil::scaling * ComputeNonbondedUtil::dielectric_1);
00980 const Vector ucenter = lattice.unscale(center);
00981 const BigReal ucenter_x = ucenter.x;
00982 const BigReal ucenter_y = ucenter.y;
00983 const BigReal ucenter_z = ucenter.z;
00984 const int n = numAtoms;
00985 cudaAtomList.resize(n);
00986 CudaAtom *ac = cudaAtomPtr = cudaAtomList.begin();
00987 const FullAtom *a = atom.begin();
00988 for ( int k=0; k<n; ++k ) {
00989 int j = a[k].sortOrder;
00990 ac[k].x = a[j].position.x - ucenter_x;
00991 ac[k].y = a[j].position.y - ucenter_y;
00992 ac[k].z = a[j].position.z - ucenter_z;
00993 ac[k].q = charge_scaling * a[j].charge;
00994 }
00995 }
00996 #endif
00997
00998 doMigration = (doMigration && numNeighbors) || ! patchMapRead;
00999
01000
01001 doGroupSizeCheck();
01002
01003
01004 p.resize(numAtoms);
01005 CompAtom *p_i = p.begin();
01006 pExt.resize(numAtoms);
01007 CompAtomExt *pExt_i = pExt.begin();
01008 FullAtom *a_i = atom.begin();
01009 int i; int n = numAtoms;
01010 for ( i=0; i<n; ++i ) {
01011 p_i[i] = a_i[i];
01012 pExt_i[i] = a_i[i];
01013 }
01014
01015
01016 doPairlistCheck();
01017
01018 if (flags.doMolly) mollyAverage();
01019
01020 if (flags.doLoweAndersen) loweAndersenVelocities();
01021
01022
01023 if (flags.doGBIS) {
01024
01025 numGBISP1Arrived = 0;
01026 phase1BoxClosedCalled = false;
01027 numGBISP2Arrived = 0;
01028 phase2BoxClosedCalled = false;
01029 numGBISP3Arrived = 0;
01030 phase3BoxClosedCalled = false;
01031 if (doMigration || isNewProxyAdded)
01032 setGBISIntrinsicRadii();
01033 }
01034
01035 if (flags.doLCPO) {
01036 if (doMigration || isNewProxyAdded) {
01037 setLcpoType();
01038 }
01039 }
01040
01041
01042 NodeIDList::iterator pli;
01043 int *pids = NULL;
01044 int pidsPreAllocated = 1;
01045 int npid;
01046 if (proxySendSpanning == 0) {
01047 npid = proxy.size();
01048 pids = new int[npid];
01049 pidsPreAllocated = 0;
01050 int *pidi = pids;
01051 int *pide = pids + proxy.size();
01052 int patchNodesLast =
01053 ( PatchMap::Object()->numNodesWithPatches() < ( 0.7 * CkNumPes() ) );
01054 for ( pli = proxy.begin(); pli != proxy.end(); ++pli )
01055 {
01056 if ( patchNodesLast && PatchMap::Object()->numPatchesOnNode(*pli) ) {
01057 *(--pide) = *pli;
01058 } else {
01059 *(pidi++) = *pli;
01060 }
01061 }
01062 }
01063 else {
01064 #ifdef NODEAWARE_PROXY_SPANNINGTREE
01065 #ifdef USE_NODEPATCHMGR
01066 npid = numNodeChild;
01067 pids = nodeChildren;
01068 #else
01069 npid = nChild;
01070 pids = child;
01071 #endif
01072 #else
01073 npid = nChild;
01074 pidsPreAllocated = 0;
01075 pids = new int[proxySpanDim];
01076 for (int i=0; i<nChild; i++) pids[i] = child[i];
01077 #endif
01078 }
01079 if (npid) {
01080 int seq = flags.sequence;
01081 int priority = PROXY_DATA_PRIORITY + PATCH_PRIORITY(patchID);
01082
01083 int pdMsgPLLen = p.size();
01084 int pdMsgAvgPLLen = 0;
01085 if(flags.doMolly) {
01086 pdMsgAvgPLLen = p_avg.size();
01087 }
01088
01089 int pdMsgVLLen = 0;
01090 if (flags.doLoweAndersen) {
01091 pdMsgVLLen = v.size();
01092 }
01093
01094
01095 int intRadLen = 0;
01096 if (flags.doGBIS && (doMigration || isNewProxyAdded)) {
01097 intRadLen = numAtoms * 2;
01098 }
01099
01100
01101 int lcpoTypeLen = 0;
01102 if (flags.doLCPO && (doMigration || isNewProxyAdded)) {
01103 lcpoTypeLen = numAtoms;
01104 }
01105
01106 int pdMsgPLExtLen = 0;
01107 if(doMigration || isNewProxyAdded) {
01108 pdMsgPLExtLen = pExt.size();
01109 }
01110
01111 int cudaAtomLen = 0;
01112 #ifdef NAMD_CUDA
01113 cudaAtomLen = numAtoms;
01114 #endif
01115
01116 ProxyDataMsg *nmsg = new (pdMsgPLLen, pdMsgAvgPLLen, pdMsgVLLen, intRadLen,
01117 lcpoTypeLen, pdMsgPLExtLen, cudaAtomLen, PRIORITY_SIZE) ProxyDataMsg;
01118
01119 SET_PRIORITY(nmsg,seq,priority);
01120 nmsg->patch = patchID;
01121 nmsg->flags = flags;
01122 nmsg->plLen = pdMsgPLLen;
01123
01124 memcpy(nmsg->positionList, p.begin(), sizeof(CompAtom)*pdMsgPLLen);
01125 nmsg->avgPlLen = pdMsgAvgPLLen;
01126 if(flags.doMolly) {
01127 memcpy(nmsg->avgPositionList, p_avg.begin(), sizeof(CompAtom)*pdMsgAvgPLLen);
01128 }
01129
01130 nmsg->vlLen = pdMsgVLLen;
01131 if (flags.doLoweAndersen) {
01132 memcpy(nmsg->velocityList, v.begin(), sizeof(CompAtom)*pdMsgVLLen);
01133 }
01134
01135
01136 if (flags.doGBIS && (doMigration || isNewProxyAdded)) {
01137 for (int i = 0; i < numAtoms * 2; i++) {
01138 nmsg->intRadList[i] = intRad[i];
01139 }
01140 }
01141
01142 if (flags.doLCPO && (doMigration || isNewProxyAdded)) {
01143 for (int i = 0; i < numAtoms; i++) {
01144 nmsg->lcpoTypeList[i] = lcpoType[i];
01145 }
01146 }
01147
01148 nmsg->plExtLen = pdMsgPLExtLen;
01149 if(doMigration || isNewProxyAdded){
01150 memcpy(nmsg->positionExtList, pExt.begin(), sizeof(CompAtomExt)*pdMsgPLExtLen);
01151 }
01152
01153 #ifdef NAMD_CUDA
01154 memcpy(nmsg->cudaAtomList, cudaAtomPtr, sizeof(CudaAtom)*cudaAtomLen);
01155 #endif
01156
01157 #if NAMD_SeparateWaters != 0
01158
01159 nmsg->numWaterAtoms = numWaterAtoms;
01160 #endif
01161
01162 #if defined(NODEAWARE_PROXY_SPANNINGTREE) && defined(USE_NODEPATCHMGR) && (CMK_SMP) && defined(NAMDSRC_IMMQD_HACK)
01163 nmsg->isFromImmMsgCall = 0;
01164 #endif
01165
01166 #if defined(PROCTRACE_DEBUG) && defined(NAST_DEBUG)
01167 DebugFileTrace *dft = DebugFileTrace::Object();
01168 dft->openTrace();
01169 dft->writeTrace("HP::posReady: for HomePatch[%d], sending proxy msg to: ", patchID);
01170 for(int i=0; i<npid; i++) {
01171 dft->writeTrace("%d ", pids[i]);
01172 }
01173 dft->writeTrace("\n");
01174 dft->closeTrace();
01175 #endif
01176
01177 #if CMK_PERSISTENT_COMM && USE_PERSISTENT_TREE
01178 if (isProxyChanged || localphs == NULL)
01179 {
01180
01181
01182 if (nphs) {
01183 for (int i=0; i<nphs; i++) {
01184 CmiDestoryPersistent(localphs[i]);
01185 }
01186 delete [] localphs;
01187 }
01188 localphs = new PersistentHandle[npid];
01189 int persist_size = sizeof(envelope) + sizeof(ProxyDataMsg) + sizeof(CompAtom)*(pdMsgPLLen+pdMsgAvgPLLen+pdMsgVLLen) + intRadLen*sizeof(Real) + lcpoTypeLen*sizeof(int) + sizeof(CompAtomExt)*pdMsgPLExtLen + sizeof(CudaAtom)*cudaAtomLen + PRIORITY_SIZE/8 + 2048;
01190 for (int i=0; i<npid; i++) {
01191 #if defined(NODEAWARE_PROXY_SPANNINGTREE) && defined(USE_NODEPATCHMGR)
01192 if (proxySendSpanning)
01193 localphs[i] = CmiCreateNodePersistent(pids[i], persist_size);
01194 else
01195 #endif
01196 localphs[i] = CmiCreatePersistent(pids[i], persist_size);
01197 }
01198 nphs = npid;
01199 }
01200 CmiAssert(nphs == npid && localphs != NULL);
01201 CmiUsePersistentHandle(localphs, nphs);
01202 #endif
01203 if(doMigration || isNewProxyAdded) {
01204 ProxyMgr::Object()->sendProxyAll(nmsg,npid,pids);
01205 }else{
01206 ProxyMgr::Object()->sendProxyData(nmsg,npid,pids);
01207 }
01208 #if CMK_PERSISTENT_COMM && USE_PERSISTENT_TREE
01209 CmiUsePersistentHandle(NULL, 0);
01210 #endif
01211 isNewProxyAdded = 0;
01212 }
01213 isProxyChanged = 0;
01214 if(!pidsPreAllocated) delete [] pids;
01215 DebugM(4, "patchID("<<patchID<<") doing positions Ready\n");
01216
01217 #ifdef REMOVE_PROXYDATAMSG_EXTRACOPY
01218 positionPtrBegin = p.begin();
01219 positionPtrEnd = p.end();
01220 #endif
01221
01222 if(flags.doMolly) {
01223 avgPositionPtrBegin = p_avg.begin();
01224 avgPositionPtrEnd = p_avg.end();
01225 }
01226
01227 if (flags.doLoweAndersen) {
01228 velocityPtrBegin = v.begin();
01229 velocityPtrEnd = v.end();
01230 }
01231
01232 Patch::positionsReady(doMigration);
01233
01234 patchMapRead = 1;
01235
01236
01237 Sync::Object()->PatchReady();
01238 }
01239
01240 void HomePatch::replaceForces(ExtForce *f)
01241 {
01242 replacementForces = f;
01243 }
01244
01245 void HomePatch::saveForce(const int ftag)
01246 {
01247 f_saved[ftag].resize(numAtoms);
01248 for ( int i = 0; i < numAtoms; ++i )
01249 {
01250 f_saved[ftag][i] = f[ftag][i];
01251 }
01252 }
01253
01254
01255 #undef DEBUG_REDISTRIB_FORCE
01256 #undef DEBUG_REDISTRIB_FORCE_VERBOSE
01257
01258
01259
01260
01261
01262
01263
01264
01265
01266
01267
01268
01269
01270
01271
01272
01273
01274
01275
01276
01277
01278
01279
01280
01281
01282
01283 #define FIX_FOR_WATER
01284 void HomePatch::redistrib_lp_force(
01285 Vector& fi, Vector& fj, Vector& fk, Vector& fl,
01286 const Vector& ri, const Vector& rj, const Vector& rk, const Vector& rl,
01287 Tensor *virial, int midpt) {
01288 #ifdef DEBUG_REDISTRIB_FORCE
01289 Vector foldnet, toldnet;
01290 foldnet = fi + fj + fk + fl;
01291 toldnet = cross(ri,fi) + cross(rj,fj) + cross(rk,fk) + cross(rl,fl);
01292 #endif
01293 Vector fja(0), fka(0), fla(0);
01294
01295 Vector r = ri - rj;
01296 BigReal invr2 = r.rlength();
01297 invr2 *= invr2;
01298 BigReal fdot = (fi*r) * invr2;
01299 Vector fr = r * fdot;
01300
01301 fja += fr;
01302
01303 Vector s, t;
01304 if (midpt) {
01305 s = rj - 0.5*(rk + rl);
01306 t = 0.5*(rk - rl);
01307 }
01308 else {
01309 s = rj - rk;
01310 t = rk - rl;
01311 }
01312 BigReal invs2 = s.rlength();
01313 invs2 *= invs2;
01314
01315 Vector p = cross(r,s);
01316 #if !defined(FIX_FOR_WATER)
01317 BigReal invp = p.rlength();
01318 #else
01319 BigReal p2 = p.length2();
01320 #endif
01321
01322 Vector q = cross(s,t);
01323 BigReal invq = q.rlength();
01324
01325 #if !defined(FIX_FOR_WATER)
01326 BigReal fpdot = (fi*p) * invp;
01327 Vector fp = p * fpdot;
01328 Vector ft = fi - fr - fp;
01329 #else
01330 BigReal fpdot;
01331 Vector fp, ft;
01332 if (p2 < 1e-6) {
01333 fpdot = 0;
01334 fp = 0;
01335 ft = fi - fr;
01336 }
01337 else {
01338 fpdot = (fi*p) / sqrt(p2);
01339 fp = p * fpdot;
01340 ft = fi - fr - fp;
01341 }
01342 #endif
01343
01344 fja += ft;
01345 Vector v = cross(r,ft);
01346 ft = cross(s,v) * invs2;
01347 fja -= ft;
01348
01349 if (midpt) {
01350 fka += 0.5 * ft;
01351 fla += 0.5 * ft;
01352 }
01353 else {
01354 fka += ft;
01355 }
01356
01357 BigReal srdot = (s*r) * invs2;
01358 Vector rr = r - s*srdot;
01359 BigReal rrdot = rr.length();
01360 BigReal stdot = (s*t) * invs2;
01361 Vector tt = t - s*stdot;
01362 BigReal invtt = tt.rlength();
01363 BigReal fact = rrdot*fpdot*invtt*invq;
01364 Vector fq = q * fact;
01365
01366 fla += fq;
01367 fja += fp*(1+srdot) + fq*stdot;
01368
01369 ft = fq*(1+stdot) + fp*srdot;
01370
01371 if (midpt) {
01372 fka += -0.5*ft;
01373 fla += -0.5*ft;
01374 }
01375 else {
01376 fka -= ft;
01377 }
01378
01379 if (virial) {
01380 Tensor va = outer(fja,rj);
01381 va += outer(fka,rk);
01382 va += outer(fla,rl);
01383 va -= outer(fi,ri);
01384 *virial += va;
01385 }
01386
01387 fi = 0;
01388 fj += fja;
01389 fk += fka;
01390 fl += fla;
01391
01392 #ifdef DEBUG_REDISTRIB_FORCE
01393 #define TOL_REDISTRIB 1e-4
01394 Vector fnewnet, tnewnet;
01395 fnewnet = fi + fj + fk + fl;
01396 tnewnet = cross(ri,fi) + cross(rj,fj) + cross(rk,fk) + cross(rl,fl);
01397 Vector fdiff = fnewnet - foldnet;
01398 Vector tdiff = tnewnet - toldnet;
01399 if (fdiff.length2() > TOL_REDISTRIB*TOL_REDISTRIB) {
01400 printf("Error: force redistribution for water exceeded tolerance: "
01401 "fdiff=(%f, %f, %f)\n", fdiff.x, fdiff.y, fdiff.z);
01402 }
01403 if (tdiff.length2() > TOL_REDISTRIB*TOL_REDISTRIB) {
01404 printf("Error: torque redistribution for water exceeded tolerance: "
01405 "tdiff=(%f, %f, %f)\n", tdiff.x, tdiff.y, tdiff.z);
01406 }
01407 #endif
01408 }
01409
01410
01411
01412
01413
01414
01415
01416
01417
01418
01419
01420 void HomePatch::redistrib_lp_water_force(
01421 Vector& f_ox, Vector& f_h1, Vector& f_h2, Vector& f_lp,
01422 const Vector& p_ox, const Vector& p_h1, const Vector& p_h2,
01423 const Vector& p_lp, Tensor *virial) {
01424
01425 #ifdef DEBUG_REDISTRIB_FORCE
01426
01427
01428
01429 Vector totforce, tottorque;
01430 totforce = f_ox + f_h1 + f_h2 + f_lp;
01431 tottorque = cross(f_ox, p_ox) + cross(f_h1, p_h1) + cross(f_h2, p_h2);
01432
01433
01434 tottorque += cross(f_lp, p_lp);
01435
01436
01437 #endif
01438
01439
01440 Vector fad_ox(0), fad_h(0);
01441
01442
01443 Vector r_ox_lp = p_lp - p_ox;
01444 BigReal invlen2_r_ox_lp = r_ox_lp.rlength();
01445 invlen2_r_ox_lp *= invlen2_r_ox_lp;
01446 BigReal rad_factor = (f_lp * r_ox_lp) * invlen2_r_ox_lp;
01447 Vector f_rad = r_ox_lp * rad_factor;
01448
01449 fad_ox += f_rad;
01450
01451
01452 Vector r_hcom_ox = p_ox - ( (p_h1 + p_h2) * 0.5 );
01453 Vector r_h2_h1_2 = (p_h1 - p_h2) * 0.5;
01454
01455
01456
01457
01458
01459
01460
01461
01462
01463
01464
01465
01466 Vector f_ang = f_lp - f_rad;
01467
01468
01469 BigReal len_r_ox_lp = r_ox_lp.length();
01470 BigReal invlen_r_hcom_ox = r_hcom_ox.rlength();
01471 BigReal oxcomp = (r_hcom_ox.length() - len_r_ox_lp) * invlen_r_hcom_ox;
01472 BigReal hydcomp = 0.5 * len_r_ox_lp * invlen_r_hcom_ox;
01473
01474 fad_ox += (f_ang * oxcomp);
01475 fad_h += (f_ang * hydcomp);
01476
01477
01478 if (virial) {
01479 Tensor vir = outer(fad_ox, p_ox);
01480 vir += outer(fad_h, p_h1);
01481 vir += outer(fad_h, p_h2);
01482 vir -= outer(f_lp, p_lp);
01483 *virial += vir;
01484 }
01485
01486
01487 f_lp = 0;
01488 f_ox += fad_ox;
01489 f_h1 += fad_h;
01490 f_h2 += fad_h;
01491
01492 #ifdef DEBUG_REDISTRIB_FORCE
01493
01494 Vector newforce, newtorque;
01495 newforce = f_ox + f_h1 + f_h2;
01496 newtorque = cross(f_ox, p_ox) + cross(f_h1, p_h1) + cross(f_h2, p_h2);
01497 Vector fdiff = newforce - totforce;
01498 Vector tdiff = newtorque - tottorque;
01499 BigReal error = fdiff.length();
01500 if (error > 0.0001) {
01501 printf("Error: Force redistribution for water "
01502 "exceeded force tolerance: error=%f\n", error);
01503 }
01504 #ifdef DEBUG_REDISTRIB_FORCE_VERBOSE
01505 printf("Error in net force: %f\n", error);
01506 #endif
01507
01508 error = tdiff.length();
01509 if (error > 0.0001) {
01510 printf("Error: Force redistribution for water "
01511 "exceeded torque tolerance: error=%f\n", error);
01512 }
01513 #ifdef DEBUG_REDISTRIB_FORCE_VERBOSE
01514 printf("Error in net torque: %f\n", error);
01515 #endif
01516 #endif
01517 }
01518
01519 void HomePatch::reposition_lonepair(
01520 Vector& ri, const Vector& rj, const Vector& rk, const Vector& rl,
01521 Real distance, Real angle, Real dihedral)
01522 {
01523 if ( (rj-rk).length2() > 100. || (rj-rl).length2() > 100. ) {
01524 iout << iWARN << "Large distance between lonepair reference atoms: ("
01525 << rj << ") (" << rk << ") (" << rl << ")\n" << endi;
01526 }
01527 BigReal r, t, p, cst, snt, csp, snp, invlen;
01528 Vector v, w, a, b, c;
01529
01530 if (distance >= 0) {
01531 v = rk;
01532 r = distance;
01533 }
01534 else {
01535 v = 0.5*(rk + rl);
01536 r = -distance;
01537 }
01538
01539 t = angle;
01540 p = dihedral;
01541 cst = cos(t);
01542 snt = sin(t);
01543 csp = cos(p);
01544 snp = sin(p);
01545 a = v - rj;
01546 b = rl - v;
01547 invlen = a.rlength();
01548 a *= invlen;
01549 c = cross(b, a);
01550 invlen = c.rlength();
01551 c *= invlen;
01552 b = cross(a, c);
01553 w.x = r*cst;
01554 w.y = r*snt*csp;
01555 w.z = r*snt*snp;
01556 ri.x = rj.x + w.x*a.x + w.y*b.x + w.z*c.x;
01557 ri.y = rj.y + w.x*a.y + w.y*b.y + w.z*c.y;
01558 ri.z = rj.z + w.x*a.z + w.y*b.z + w.z*c.z;
01559 }
01560
01561 void HomePatch::reposition_all_lonepairs(void) {
01562
01563 for (int i=0; i < numAtoms; i++) {
01564 if (atom[i].mass < 0.01) {
01565
01566 AtomID aid = atom[i].id;
01567 Lphost *lph = Node::Object()->molecule->get_lphost(aid);
01568 if (lph == NULL) {
01569 char errmsg[512];
01570 sprintf(errmsg, "reposition lone pairs: "
01571 "no Lphost exists for LP %d\n", aid);
01572 NAMD_die(errmsg);
01573 }
01574 LocalID j = AtomMap::Object()->localID(lph->atom2);
01575 LocalID k = AtomMap::Object()->localID(lph->atom3);
01576 LocalID l = AtomMap::Object()->localID(lph->atom4);
01577 if (j.pid != patchID || k.pid != patchID || l.pid != patchID) {
01578 char errmsg[512];
01579 sprintf(errmsg, "reposition lone pairs: "
01580 "LP %d has some Lphost atom off patch\n", aid);
01581 NAMD_die(errmsg);
01582 }
01583
01584 reposition_lonepair(atom[i].position, atom[j.index].position,
01585 atom[k.index].position, atom[l.index].position,
01586 lph->distance, lph->angle, lph->dihedral);
01587 }
01588 }
01589 }
01590
01591 void HomePatch::swm4_omrepos(Vector *ref, Vector *pos, Vector *vel,
01592 BigReal invdt) {
01593
01594
01595
01596 pos[2] = pos[0] + (0.5 * (pos[3] + pos[4]) - pos[0]) * (r_om / r_ohc);
01597
01598
01599 if (invdt != 0) {
01600 vel[2] = (pos[2] - ref[2]) * invdt;
01601 }
01602
01603 }
01604
01605 void HomePatch::tip4_omrepos(Vector* ref, Vector* pos, Vector* vel, BigReal invdt) {
01606
01607
01608
01609
01610
01611
01612
01613
01614
01615
01616
01617
01618 pos[3] = pos[0] + (0.5 * (pos[1] + pos[2]) - pos[0]) * (r_om / r_ohc);
01619
01620
01621
01622 if (invdt != 0) {
01623 vel[3] = (pos[3] - ref[3]) * invdt;
01624 }
01625
01626
01627 return;
01628 }
01629
01630 void HomePatch::redistrib_lonepair_forces(const int ftag, Tensor *virial) {
01631
01632 ForceList *f_mod = f;
01633 for (int i = 0; i < numAtoms; i++) {
01634 if (atom[i].mass < 0.01) {
01635
01636 AtomID aid = atom[i].id;
01637 Lphost *lph = Node::Object()->molecule->get_lphost(aid);
01638 if (lph == NULL) {
01639 char errmsg[512];
01640 sprintf(errmsg, "redistrib lone pair forces: "
01641 "no Lphost exists for LP %d\n", aid);
01642 NAMD_die(errmsg);
01643 }
01644 LocalID j = AtomMap::Object()->localID(lph->atom2);
01645 LocalID k = AtomMap::Object()->localID(lph->atom3);
01646 LocalID l = AtomMap::Object()->localID(lph->atom4);
01647 if (j.pid != patchID || k.pid != patchID || l.pid != patchID) {
01648 char errmsg[512];
01649 sprintf(errmsg, "redistrib lone pair forces: "
01650 "LP %d has some Lphost atom off patch\n", aid);
01651 NAMD_die(errmsg);
01652 }
01653
01654 int midpt = (lph->distance < 0);
01655 redistrib_lp_force(f_mod[ftag][i], f_mod[ftag][j.index],
01656 f_mod[ftag][k.index], f_mod[ftag][l.index],
01657 atom[i].position, atom[j.index].position,
01658 atom[k.index].position, atom[l.index].position, virial, midpt);
01659 }
01660 }
01661 }
01662
01663 void HomePatch::redistrib_swm4_forces(const int ftag, Tensor *virial) {
01664
01665
01666 ForceList *f_mod = f;
01667 for (int i = 0; i < numAtoms; i++) {
01668 if (atom[i].mass < 0.01) {
01669
01670 redistrib_lp_water_force(f_mod[ftag][i-2], f_mod[ftag][i+1],
01671 f_mod[ftag][i+2], f_mod[ftag][i],
01672 atom[i-2].position, atom[i+1].position,
01673 atom[i+2].position, atom[i].position, virial);
01674 }
01675 }
01676 }
01677
01678 void HomePatch::redistrib_tip4p_forces(const int ftag, Tensor* virial) {
01679
01680
01681
01682
01683 ForceList *f_mod =f;
01684 for (int i=0; i<numAtoms; i++) {
01685 if (atom[i].mass < 0.01) {
01686
01687 redistrib_lp_water_force(f_mod[ftag][i-3], f_mod[ftag][i-2],
01688 f_mod[ftag][i-1], f_mod[ftag][i],
01689 atom[i-3].position, atom[i-2].position,
01690 atom[i-1].position, atom[i].position, virial);
01691 }
01692 }
01693 }
01694
01695
01696 void HomePatch::addForceToMomentum(const BigReal timestep, const int ftag,
01697 const int useSaved)
01698 {
01699 SimParameters *simParams = Node::Object()->simParameters;
01700 const BigReal dt = timestep / TIMEFACTOR;
01701 ForceList *f_use = (useSaved ? f_saved : f);
01702
01703 if ( simParams->fixedAtomsOn ) {
01704 for ( int i = 0; i < numAtoms; ++i ) {
01705 if ( atom[i].atomFixed ) {
01706 atom[i].velocity = 0;
01707 } else {
01708 BigReal recip_val = ( atom[i].mass > 0. ? dt * namd_reciprocal( atom[i].mass ) : 0.);
01709 atom[i].velocity += f_use[ftag][i] * recip_val;
01710 }
01711 }
01712 } else {
01713 FullAtom *atom_arr = atom.begin();
01714 const Force *force_arr = f_use[ftag].const_begin();
01715 #ifdef ARCH_POWERPC
01716 #pragma disjoint (*force_arr, *atom_arr)
01717 #endif
01718 for ( int i = 0; i < numAtoms; ++i ) {
01719 if (atom[i].mass == 0.) continue;
01720 BigReal recip_val = ( atom[i].mass > 0. ? dt * namd_reciprocal( atom[i].mass ) : 0.);
01721
01722 atom_arr[i].velocity.x += force_arr[i].x * recip_val;
01723 atom_arr[i].velocity.y += force_arr[i].y * recip_val;
01724 atom_arr[i].velocity.z += force_arr[i].z * recip_val;
01725 }
01726 }
01727 }
01728
01729 void HomePatch::addVelocityToPosition(const BigReal timestep)
01730 {
01731 SimParameters *simParams = Node::Object()->simParameters;
01732 const BigReal dt = timestep / TIMEFACTOR;
01733 if ( simParams->fixedAtomsOn ) {
01734 for ( int i = 0; i < numAtoms; ++i ) {
01735 if ( ! atom[i].atomFixed ) atom[i].position += atom[i].velocity * dt;
01736 }
01737 } else {
01738 FullAtom *atom_arr = atom.begin();
01739 for ( int i = 0; i < numAtoms; ++i ) {
01740 atom_arr[i].position.x += atom_arr[i].velocity.x * dt;
01741 atom_arr[i].position.y += atom_arr[i].velocity.y * dt;
01742 atom_arr[i].position.z += atom_arr[i].velocity.z * dt;
01743 }
01744 }
01745 }
01746
01747
01748 int HomePatch::rattle1(const BigReal timestep, Tensor *virial,
01749 SubmitReduction *ppreduction)
01750 {
01751 Molecule *mol = Node::Object()->molecule;
01752 SimParameters *simParams = Node::Object()->simParameters;
01753 const int fixedAtomsOn = simParams->fixedAtomsOn;
01754 const int useSettle = simParams->useSettle;
01755 const BigReal dt = timestep / TIMEFACTOR;
01756 const BigReal invdt = (dt == 0.) ? 0. : 1.0 / dt;
01757 BigReal tol2 = 2.0 * simParams->rigidTol;
01758 int maxiter = simParams->rigidIter;
01759 int dieOnError = simParams->rigidDie;
01760 int i, iter;
01761 BigReal dsq[10], tmp;
01762 int ial[10], ibl[10];
01763 Vector ref[10];
01764 Vector refab[10];
01765 Vector pos[10];
01766 Vector vel[10];
01767 Vector netdp[10];
01768 BigReal rmass[10];
01769 int fixed[10];
01770 Tensor wc;
01771 BigReal idz, zmin;
01772 int nslabs;
01773
01774
01775
01776
01777
01778
01779 if ( ! settle_initialized ) {
01780 for ( int ig = 0; ig < numAtoms; ig += atom[ig].hydrogenGroupSize ) {
01781
01782 if (atom[ig].rigidBondLength > 0) {
01783 int oatm;
01784 if (simParams->watmodel == WAT_SWM4) {
01785 oatm = ig+3;
01786
01787
01788 }
01789 else {
01790 oatm = ig+1;
01791
01792 if (atom[ig].mass < 0.5 || atom[ig+1].mass < 0.5) {
01793 oatm += 1;
01794 }
01795 }
01796
01797
01798 settle1init(atom[ig].mass, atom[oatm].mass,
01799 atom[ig].rigidBondLength,
01800 atom[oatm].rigidBondLength,
01801 settle_mOrmT, settle_mHrmT, settle_ra,
01802 settle_rb, settle_rc, settle_rra);
01803 settle_initialized = 1;
01804 break;
01805 }
01806 }
01807 }
01808
01809 if (ppreduction) {
01810 nslabs = simParams->pressureProfileSlabs;
01811 idz = nslabs/lattice.c().z;
01812 zmin = lattice.origin().z - 0.5*lattice.c().z;
01813 }
01814
01815
01816 int wathgsize = 3;
01817 int watmodel = simParams->watmodel;
01818 if (watmodel == WAT_TIP4) wathgsize = 4;
01819 else if (watmodel == WAT_SWM4) wathgsize = 5;
01820
01821 for ( int ig = 0; ig < numAtoms; ig += atom[ig].hydrogenGroupSize ) {
01822 int hgs = atom[ig].hydrogenGroupSize;
01823 if ( hgs == 1 ) continue;
01824
01825 for ( i = 0; i < hgs; ++i ) {
01826 ref[i] = atom[ig+i].position;
01827 pos[i] = atom[ig+i].position;
01828 vel[i] = atom[ig+i].velocity;
01829 rmass[i] = (atom[ig+i].mass > 0. ? 1. / atom[ig+i].mass : 0.);
01830
01831 fixed[i] = ( fixedAtomsOn && atom[ig+i].atomFixed );
01832
01833
01834 if ( fixed[i] ) rmass[i] = 0.; else pos[i] += vel[i] * dt;
01835 }
01836 int icnt = 0;
01837 if ( ( tmp = atom[ig].rigidBondLength ) > 0 ) {
01838 if (hgs != wathgsize) {
01839 NAMD_bug("Hydrogen group error caught in rattle1().");
01840 }
01841
01842
01843 if (useSettle && !atom[ig].groupFixed) {
01844 if (simParams->watmodel == WAT_SWM4) {
01845
01846
01847
01848 Vector lp_ref = ref[2];
01849 Vector lp_pos = pos[2];
01850 Vector lp_vel = vel[2];
01851 ref[2] = ref[0];
01852 pos[2] = pos[0];
01853 vel[2] = vel[0];
01854 settle1(ref+2, pos+2, vel+2, invdt,
01855 settle_mOrmT, settle_mHrmT, settle_ra,
01856 settle_rb, settle_rc, settle_rra);
01857 ref[0] = ref[2];
01858 pos[0] = pos[2];
01859 vel[0] = vel[2];
01860 ref[2] = lp_ref;
01861 pos[2] = lp_pos;
01862 vel[2] = lp_vel;
01863
01864 swm4_omrepos(ref, pos, vel, invdt);
01865 }
01866 else {
01867 settle1(ref, pos, vel, invdt,
01868 settle_mOrmT, settle_mHrmT, settle_ra,
01869 settle_rb, settle_rc, settle_rra);
01870 if (simParams->watmodel == WAT_TIP4) {
01871 tip4_omrepos(ref, pos, vel, invdt);
01872 }
01873 }
01874
01875
01876
01877 int ppoffset, partition;
01878 if ( invdt == 0 ) for ( i = 0; i < wathgsize; ++i ) {
01879 atom[ig+i].position = pos[i];
01880 } else if ( virial == 0 ) for ( i = 0; i < wathgsize; ++i ) {
01881 atom[ig+i].velocity = vel[i];
01882 } else for ( i = 0; i < wathgsize; ++i ) {
01883 Force df = (vel[i] - atom[ig+i].velocity) * ( atom[ig+i].mass * invdt );
01884 Tensor vir = outer(df, ref[i]);
01885 wc += vir;
01886 f[Results::normal][ig+i] += df;
01887 atom[ig+i].velocity = vel[i];
01888 if (ppreduction) {
01889
01890
01891 if (!i) {
01892 BigReal z = pos[i].z;
01893 partition = atom[ig].partition;
01894 int slab = (int)floor((z-zmin)*idz);
01895 if (slab < 0) slab += nslabs;
01896 else if (slab >= nslabs) slab -= nslabs;
01897 ppoffset = 3*(slab + nslabs*partition);
01898 }
01899 ppreduction->item(ppoffset ) += vir.xx;
01900 ppreduction->item(ppoffset+1) += vir.yy;
01901 ppreduction->item(ppoffset+2) += vir.zz;
01902 }
01903 }
01904 continue;
01905 }
01906 if ( !(fixed[1] && fixed[2]) ) {
01907 dsq[icnt] = tmp * tmp; ial[icnt] = 1; ibl[icnt] = 2; ++icnt;
01908 }
01909 }
01910 for ( i = 1; i < hgs; ++i ) {
01911 if ( ( tmp = atom[ig+i].rigidBondLength ) > 0 ) {
01912 if ( !(fixed[0] && fixed[i]) ) {
01913 dsq[icnt] = tmp * tmp; ial[icnt] = 0; ibl[icnt] = i; ++icnt;
01914 }
01915 }
01916 }
01917 if ( icnt == 0 ) continue;
01918 for ( i = 0; i < icnt; ++i ) {
01919 refab[i] = ref[ial[i]] - ref[ibl[i]];
01920 }
01921 for ( i = 0; i < hgs; ++i ) {
01922 netdp[i] = 0.;
01923 }
01924 int done;
01925 int consFailure;
01926 for ( iter = 0; iter < maxiter; ++iter ) {
01927
01928 done = 1;
01929 consFailure = 0;
01930 for ( i = 0; i < icnt; ++i ) {
01931 int a = ial[i]; int b = ibl[i];
01932 Vector pab = pos[a] - pos[b];
01933 BigReal pabsq = pab.x*pab.x + pab.y*pab.y + pab.z*pab.z;
01934 BigReal rabsq = dsq[i];
01935 BigReal diffsq = rabsq - pabsq;
01936 if ( fabs(diffsq) > (rabsq * tol2) ) {
01937 Vector &rab = refab[i];
01938 BigReal rpab = rab.x*pab.x + rab.y*pab.y + rab.z*pab.z;
01939 if ( rpab < ( rabsq * 1.0e-6 ) ) {
01940 done = 0;
01941 consFailure = 1;
01942 continue;
01943 }
01944 BigReal rma = rmass[a];
01945 BigReal rmb = rmass[b];
01946 BigReal gab = diffsq / ( 2.0 * ( rma + rmb ) * rpab );
01947 Vector dp = rab * gab;
01948 pos[a] += rma * dp;
01949 pos[b] -= rmb * dp;
01950 if ( invdt != 0. ) {
01951 dp *= invdt;
01952 netdp[a] += dp;
01953 netdp[b] -= dp;
01954 }
01955 done = 0;
01956 }
01957 }
01958 if ( done ) break;
01959 }
01960
01961 if ( consFailure ) {
01962 if ( dieOnError ) {
01963 iout << iERROR << "Constraint failure in RATTLE algorithm for atom "
01964 << (atom[ig].id + 1) << "!\n" << endi;
01965 return -1;
01966 } else {
01967 iout << iWARN << "Constraint failure in RATTLE algorithm for atom "
01968 << (atom[ig].id + 1) << "!\n" << endi;
01969 }
01970 } else if ( ! done ) {
01971 if ( dieOnError ) {
01972 iout << iERROR << "Exceeded RATTLE iteration limit for atom "
01973 << (atom[ig].id + 1) << "!\n" << endi;
01974 return -1;
01975 } else {
01976 iout << iWARN << "Exceeded RATTLE iteration limit for atom "
01977 << (atom[ig].id + 1) << "!\n" << endi;
01978 }
01979 }
01980
01981
01982 int ppoffset, partition;
01983 if ( invdt == 0 ) for ( i = 0; i < hgs; ++i ) {
01984 atom[ig+i].position = pos[i];
01985 } else if ( virial == 0 ) for ( i = 0; i < hgs; ++i ) {
01986 atom[ig+i].velocity = vel[i] + rmass[i] * netdp[i];
01987 } else for ( i = 0; i < hgs; ++i ) {
01988 Force df = netdp[i] * invdt;
01989 Tensor vir = outer(df, ref[i]);
01990 wc += vir;
01991 f[Results::normal][ig+i] += df;
01992 atom[ig+i].velocity = vel[i] + rmass[i] * netdp[i];
01993 if (ppreduction) {
01994 if (!i) {
01995 BigReal z = pos[i].z;
01996 int partition = atom[ig].partition;
01997 int slab = (int)floor((z-zmin)*idz);
01998 if (slab < 0) slab += nslabs;
01999 else if (slab >= nslabs) slab -= nslabs;
02000 ppoffset = 3*(slab + nslabs*partition);
02001 }
02002 ppreduction->item(ppoffset ) += vir.xx;
02003 ppreduction->item(ppoffset+1) += vir.yy;
02004 ppreduction->item(ppoffset+2) += vir.zz;
02005 }
02006 }
02007 }
02008 if ( dt && virial ) *virial += wc;
02009
02010 return 0;
02011 }
02012
02013
02014 void HomePatch::rattle2(const BigReal timestep, Tensor *virial)
02015 {
02016 Molecule *mol = Node::Object()->molecule;
02017 SimParameters *simParams = Node::Object()->simParameters;
02018 const int fixedAtomsOn = simParams->fixedAtomsOn;
02019 const int useSettle = simParams->useSettle;
02020 const BigReal dt = timestep / TIMEFACTOR;
02021 Tensor wc;
02022 BigReal tol = simParams->rigidTol;
02023 int maxiter = simParams->rigidIter;
02024 int dieOnError = simParams->rigidDie;
02025 int i, iter;
02026 BigReal dsqi[10], tmp;
02027 int ial[10], ibl[10];
02028 Vector ref[10];
02029 Vector refab[10];
02030 Vector vel[10];
02031 BigReal rmass[10];
02032 BigReal redmass[10];
02033 int fixed[10];
02034
02035
02036 int wathgsize = 3;
02037 if (simParams->watmodel == WAT_TIP4) wathgsize = 4;
02038 else if (simParams->watmodel == WAT_SWM4) wathgsize = 5;
02039
02040
02041 for ( int ig = 0; ig < numAtoms; ig += atom[ig].hydrogenGroupSize ) {
02042
02043 int hgs = atom[ig].hydrogenGroupSize;
02044 if ( hgs == 1 ) continue;
02045
02046 for ( i = 0; i < hgs; ++i ) {
02047 ref[i] = atom[ig+i].position;
02048 vel[i] = atom[ig+i].velocity;
02049 rmass[i] = atom[ig+i].mass > 0. ? 1. / atom[ig+i].mass : 0.;
02050 fixed[i] = ( fixedAtomsOn && atom[ig+i].atomFixed );
02051 if ( fixed[i] ) rmass[i] = 0.;
02052 }
02053 int icnt = 0;
02054 if ( ( tmp = atom[ig].rigidBondLength ) > 0 ) {
02055 if ( wathgsize != 4 ) {
02056 NAMD_bug("Hydrogen group error caught in rattle2().");
02057 }
02058
02059 if (useSettle && !fixed[0] && !fixed[1] && !fixed[2]) {
02060 settle2(atom[ig].mass, atom[ig+1].mass, ref, vel, dt, virial);
02061 for (i=0; i<3; i++) {
02062 atom[ig+i].velocity = vel[i];
02063 }
02064 continue;
02065 }
02066 if ( !(fixed[1] && fixed[2]) ) {
02067 redmass[icnt] = 1. / (rmass[1] + rmass[2]);
02068 dsqi[icnt] = 1. / (tmp * tmp); ial[icnt] = 1; ibl[icnt] = 2; ++icnt;
02069 }
02070 }
02071
02072 for ( i = 1; i < hgs; ++i ) {
02073 if ( ( tmp = atom[ig+i].rigidBondLength ) > 0 ) {
02074 if ( !(fixed[0] && fixed[i]) ) {
02075 redmass[icnt] = 1. / (rmass[0] + rmass[i]);
02076 dsqi[icnt] = 1. / (tmp * tmp); ial[icnt] = 0;
02077 ibl[icnt] = i; ++icnt;
02078 }
02079 }
02080 }
02081 if ( icnt == 0 ) continue;
02082
02083 for ( i = 0; i < icnt; ++i ) {
02084 refab[i] = ref[ial[i]] - ref[ibl[i]];
02085 }
02086
02087 int done;
02088 for ( iter = 0; iter < maxiter; ++iter ) {
02089 done = 1;
02090 for ( i = 0; i < icnt; ++i ) {
02091 int a = ial[i]; int b = ibl[i];
02092 Vector vab = vel[a] - vel[b];
02093 Vector &rab = refab[i];
02094 BigReal rabsqi = dsqi[i];
02095 BigReal rvab = rab.x*vab.x + rab.y*vab.y + rab.z*vab.z;
02096 if ( (fabs(rvab) * dt * rabsqi) > tol ) {
02097 Vector dp = rab * (-rvab * redmass[i] * rabsqi);
02098 wc += outer(dp,rab);
02099 vel[a] += rmass[a] * dp;
02100 vel[b] -= rmass[b] * dp;
02101 done = 0;
02102 }
02103 }
02104 if ( done ) break;
02105
02106 }
02107 if ( ! done ) {
02108 if ( dieOnError ) {
02109 NAMD_die("Exceeded maximum number of iterations in rattle2().");
02110 } else {
02111 iout << iWARN <<
02112 "Exceeded maximum number of iterations in rattle2().\n" << endi;
02113 }
02114 }
02115
02116 for ( i = 0; i < hgs; ++i ) {
02117 atom[ig+i].velocity = vel[i];
02118 }
02119 }
02120
02121
02122 *virial += wc / ( 0.5 * dt );
02123
02124 }
02125
02126
02127
02128 void HomePatch::loweAndersenVelocities()
02129 {
02130 DebugM(2, "loweAndersenVelocities\n");
02131 Molecule *mol = Node::Object()->molecule;
02132 SimParameters *simParams = Node::Object()->simParameters;
02133 v.resize(numAtoms);
02134 for (int i = 0; i < numAtoms; ++i) {
02135
02136
02137 v[i].position = atom[i].velocity;
02138 v[i].charge = atom[i].mass;
02139 }
02140 DebugM(2, "loweAndersenVelocities\n");
02141 }
02142
02143 void HomePatch::loweAndersenFinish()
02144 {
02145 DebugM(2, "loweAndersenFinish\n");
02146 v.resize(0);
02147 }
02148
02149
02150
02151 void HomePatch::setLcpoType() {
02152 Molecule *mol = Node::Object()->molecule;
02153 const int *lcpoParamType = mol->getLcpoParamType();
02154
02155 lcpoType.resize(numAtoms);
02156 for (int i = 0; i < numAtoms; i++) {
02157 lcpoType[i] = lcpoParamType[pExt[i].id];
02158 }
02159 }
02160
02161
02162 void HomePatch::setGBISIntrinsicRadii() {
02163 intRad.resize(numAtoms*2);
02164 intRad.setall(0);
02165 Molecule *mol = Node::Object()->molecule;
02166 SimParameters *simParams = Node::Object()->simParameters;
02167 Real offset = simParams->coulomb_radius_offset;
02168 for (int i = 0; i < numAtoms; i++) {
02169 Real rad = MassToRadius(atom[i].mass);
02170 Real screen = MassToScreen(atom[i].mass);
02171 intRad[2*i+0] = rad - offset;
02172 intRad[2*i+1] = screen*(rad - offset);
02173 }
02174 }
02175
02176
02177 void HomePatch::gbisComputeAfterP1() {
02178
02179 SimParameters *simParams = Node::Object()->simParameters;
02180 BigReal alphaMax = simParams->alpha_max;
02181 BigReal delta = simParams->gbis_delta;
02182 BigReal beta = simParams->gbis_beta;
02183 BigReal gamma = simParams->gbis_gamma;
02184 BigReal coulomb_radius_offset = simParams->coulomb_radius_offset;
02185
02186 BigReal rhoi;
02187 BigReal rhoi0;
02188
02189 for (int i = 0; i < numAtoms; i++) {
02190 rhoi0 = intRad[2*i];
02191 rhoi = rhoi0+coulomb_radius_offset;
02192 psiFin[i] += psiSum[i];
02193 psiFin[i] *= rhoi0;
02194 bornRad[i]=1/(1/rhoi0-1/rhoi*tanh(psiFin[i]*(delta+psiFin[i]*(-beta+gamma*psiFin[i]))));
02195 bornRad[i] = (bornRad[i] > alphaMax) ? alphaMax : bornRad[i];
02196 #ifdef PRINT_COMP
02197 CkPrintf("BORNRAD(%04i)[%04i] = % .4e\n",flags.sequence,pExt[i].id,bornRad[i]);
02198 #endif
02199 }
02200
02201 gbisP2Ready();
02202 }
02203
02204
02205 void HomePatch::gbisComputeAfterP2() {
02206
02207 SimParameters *simParams = Node::Object()->simParameters;
02208 BigReal delta = simParams->gbis_delta;
02209 BigReal beta = simParams->gbis_beta;
02210 BigReal gamma = simParams->gbis_gamma;
02211 BigReal epsilon_s = simParams->solvent_dielectric;
02212 BigReal epsilon_p = simParams->dielectric;
02213 BigReal epsilon_s_i = 1/simParams->solvent_dielectric;
02214 BigReal epsilon_p_i = 1/simParams->dielectric;
02215 BigReal coulomb_radius_offset = simParams->coulomb_radius_offset;
02216 BigReal kappa = simParams->kappa;
02217 BigReal fij, expkappa, Dij, dEdai, dedasum;
02218 BigReal rhoi, rhoi0, psii, nbetapsi;
02219 BigReal gammapsi2, tanhi, daidr;
02220 for (int i = 0; i < numAtoms; i++) {
02221
02222 dHdrPrefix[i] += dEdaSum[i];
02223 fij = bornRad[i];
02224 expkappa = exp(-kappa*fij);
02225 Dij = epsilon_p_i - expkappa*epsilon_s_i;
02226
02227 dEdai = -0.5*COULOMB*atom[i].charge*atom[i].charge
02228 *(kappa*epsilon_s_i*expkappa-Dij/fij)/bornRad[i];
02229 dHdrPrefix[i] += dEdai;
02230 dedasum = dHdrPrefix[i];
02231
02232 rhoi0 = intRad[2*i];
02233 rhoi = rhoi0+coulomb_radius_offset;
02234 psii = psiFin[i];
02235 nbetapsi = -beta*psii;
02236 gammapsi2 = gamma*psii*psii;
02237 tanhi = tanh(psii*(delta+nbetapsi+gammapsi2));
02238 daidr = bornRad[i]*bornRad[i]*rhoi0/rhoi*(1-tanhi*tanhi)
02239 * (delta+nbetapsi+nbetapsi+gammapsi2+gammapsi2+gammapsi2);
02240 dHdrPrefix[i] *= daidr;
02241 #ifdef PRINT_COMP
02242 CkPrintf("DHDR(%04i)[%04i] = % .4e\n",flags.sequence,pExt[i].id,dHdrPrefix[i]);
02243 #endif
02244 }
02245 gbisP3Ready();
02246 }
02247
02248
02249 void HomePatch::gbisP2Ready() {
02250 if (proxy.size() > 0) {
02251 CProxy_ProxyMgr cp(CkpvAccess(BOCclass_group).proxyMgr);
02252 for (int i = 0; i < proxy.size(); i++) {
02253 int node = proxy[i];
02254 ProxyGBISP2DataMsg *msg=new(numAtoms,PRIORITY_SIZE) ProxyGBISP2DataMsg;
02255 msg->patch = patchID;
02256 msg->origPe = CkMyPe();
02257 memcpy(msg->bornRad,bornRad.begin(),numAtoms*sizeof(Real));
02258 msg->destPe = node;
02259 int seq = flags.sequence;
02260 int priority = GB2_PROXY_DATA_PRIORITY + PATCH_PRIORITY(patchID);
02261 SET_PRIORITY(msg,seq,priority);
02262 cp[node].recvData(msg);
02263 }
02264 }
02265 Patch::gbisP2Ready();
02266 }
02267
02268
02269 void HomePatch::gbisP3Ready() {
02270 if (proxy.size() > 0) {
02271 CProxy_ProxyMgr cp(CkpvAccess(BOCclass_group).proxyMgr);
02272
02273 int msgAtoms = (flags.doFullElectrostatics) ? numAtoms : 0;
02274 for (int i = 0; i < proxy.size(); i++) {
02275 int node = proxy[i];
02276 ProxyGBISP3DataMsg *msg = new(msgAtoms,PRIORITY_SIZE) ProxyGBISP3DataMsg;
02277 msg->patch = patchID;
02278 msg->dHdrPrefixLen = msgAtoms;
02279 msg->origPe = CkMyPe();
02280 memcpy(msg->dHdrPrefix, dHdrPrefix.begin(), msgAtoms*sizeof(Real));
02281 msg->destPe = node;
02282 int seq = flags.sequence;
02283 int priority = GB3_PROXY_DATA_PRIORITY + PATCH_PRIORITY(patchID);
02284 SET_PRIORITY(msg,seq,priority);
02285 cp[node].recvData(msg);
02286 }
02287 }
02288 Patch::gbisP3Ready();
02289 }
02290
02291
02292 void HomePatch::receiveResult(ProxyGBISP1ResultMsg *msg) {
02293 ++numGBISP1Arrived;
02294 for ( int i = 0; i < msg->psiSumLen; ++i ) {
02295 psiFin[i] += msg->psiSum[i];
02296 }
02297 delete msg;
02298
02299 if (flags.doNonbonded) {
02300
02301 if (phase1BoxClosedCalled == true &&
02302 numGBISP1Arrived==proxy.size() ) {
02303 sequencer->awaken();
02304 }
02305 } else {
02306
02307 if (boxesOpen == 0 &&
02308 numGBISP1Arrived == proxy.size() &&
02309 numGBISP2Arrived == proxy.size() &&
02310 numGBISP3Arrived == proxy.size()) {
02311 sequencer->awaken();
02312 }
02313 }
02314 }
02315
02316
02317 void HomePatch::receiveResult(ProxyGBISP2ResultMsg *msg) {
02318 ++numGBISP2Arrived;
02319
02320 for ( int i = 0; i < msg->dEdaSumLen; ++i ) {
02321 dHdrPrefix[i] += msg->dEdaSum[i];
02322 }
02323 delete msg;
02324
02325 if (flags.doNonbonded) {
02326
02327 if (phase2BoxClosedCalled == true &&
02328 numGBISP2Arrived==proxy.size() ) {
02329 sequencer->awaken();
02330 }
02331 } else {
02332
02333 if (boxesOpen == 0 &&
02334 numGBISP1Arrived == proxy.size() &&
02335 numGBISP2Arrived == proxy.size() &&
02336 numGBISP3Arrived == proxy.size()) {
02337 sequencer->awaken();
02338 }
02339 }
02340 }
02341
02342
02343 void HomePatch::mollyAverage()
02344 {
02345 Molecule *mol = Node::Object()->molecule;
02346 SimParameters *simParams = Node::Object()->simParameters;
02347 BigReal tol = simParams->mollyTol;
02348 int maxiter = simParams->mollyIter;
02349 int i, iter;
02350 HGArrayBigReal dsq;
02351 BigReal tmp;
02352 HGArrayInt ial, ibl;
02353 HGArrayVector ref;
02354 HGArrayVector refab;
02355 HGArrayBigReal rmass;
02356 BigReal *lambda;
02357 CompAtom *avg;
02358 int numLambdas = 0;
02359 HGArrayInt fixed;
02360
02361
02362 p_avg.resize(numAtoms);
02363 for ( i=0; i<numAtoms; ++i ) p_avg[i] = p[i];
02364
02365 for ( int ig = 0; ig < numAtoms; ig += atom[ig].hydrogenGroupSize ) {
02366 int hgs = atom[ig].hydrogenGroupSize;
02367 if ( hgs == 1 ) continue;
02368 for ( i = 0; i < hgs; ++i ) {
02369 ref[i] = atom[ig+i].position;
02370 rmass[i] = 1. / atom[ig+i].mass;
02371 fixed[i] = ( simParams->fixedAtomsOn && atom[ig+i].atomFixed );
02372 if ( fixed[i] ) rmass[i] = 0.;
02373 }
02374 avg = &(p_avg[ig]);
02375 int icnt = 0;
02376
02377 if ( ( tmp = atom[ig].rigidBondLength ) > 0 ) {
02378 if ( hgs != 3 ) {
02379 NAMD_die("Hydrogen group error caught in mollyAverage(). It's a bug!\n");
02380 }
02381 if ( !(fixed[1] && fixed[2]) ) {
02382 dsq[icnt] = tmp * tmp; ial[icnt] = 1; ibl[icnt] = 2; ++icnt;
02383 }
02384 }
02385 for ( i = 1; i < hgs; ++i ) {
02386 if ( ( tmp = atom[ig+i].rigidBondLength ) > 0 ) {
02387 if ( !(fixed[0] && fixed[i]) ) {
02388 dsq[icnt] = tmp * tmp; ial[icnt] = 0; ibl[icnt] = i; ++icnt;
02389 }
02390 }
02391 }
02392 if ( icnt == 0 ) continue;
02393 numLambdas += icnt;
02394 molly_lambda.resize(numLambdas);
02395 lambda = &(molly_lambda[numLambdas - icnt]);
02396 for ( i = 0; i < icnt; ++i ) {
02397 refab[i] = ref[ial[i]] - ref[ibl[i]];
02398 }
02399
02400 iter=average(avg,ref,lambda,hgs,icnt,rmass,dsq,ial,ibl,refab,tol,maxiter);
02401 if ( iter == maxiter ) {
02402 iout << iWARN << "Exceeded maximum number of iterations in mollyAverage().\n"<<endi;
02403 }
02404 }
02405
02406
02407
02408
02409
02410
02411
02412
02413 }
02414
02415
02416
02417 void HomePatch::mollyMollify(Tensor *virial)
02418 {
02419 Molecule *mol = Node::Object()->molecule;
02420 SimParameters *simParams = Node::Object()->simParameters;
02421 Tensor wc;
02422 int i;
02423 HGArrayInt ial, ibl;
02424 HGArrayVector ref;
02425 CompAtom *avg;
02426 HGArrayVector refab;
02427 HGArrayVector force;
02428 HGArrayBigReal rmass;
02429 BigReal *lambda;
02430 int numLambdas = 0;
02431 HGArrayInt fixed;
02432
02433 for ( int ig = 0; ig < numAtoms; ig += atom[ig].hydrogenGroupSize ) {
02434 int hgs = atom[ig].hydrogenGroupSize;
02435 if (hgs == 1 ) continue;
02436 for ( i = 0; i < hgs; ++i ) {
02437 ref[i] = atom[ig+i].position;
02438 force[i] = f[Results::slow][ig+i];
02439 rmass[i] = 1. / atom[ig+i].mass;
02440 fixed[i] = ( simParams->fixedAtomsOn && atom[ig+i].atomFixed );
02441 if ( fixed[i] ) rmass[i] = 0.;
02442 }
02443 int icnt = 0;
02444
02445 if ( atom[ig].rigidBondLength > 0 ) {
02446 if ( hgs != 3 ) {
02447 NAMD_die("Hydrogen group error caught in mollyMollify(). It's a bug!\n");
02448 }
02449 if ( !(fixed[1] && fixed[2]) ) {
02450 ial[icnt] = 1; ibl[icnt] = 2; ++icnt;
02451 }
02452 }
02453 for ( i = 1; i < hgs; ++i ) {
02454 if ( atom[ig+i].rigidBondLength > 0 ) {
02455 if ( !(fixed[0] && fixed[i]) ) {
02456 ial[icnt] = 0; ibl[icnt] = i; ++icnt;
02457 }
02458 }
02459 }
02460
02461 if ( icnt == 0 ) continue;
02462 lambda = &(molly_lambda[numLambdas]);
02463 numLambdas += icnt;
02464 for ( i = 0; i < icnt; ++i ) {
02465 refab[i] = ref[ial[i]] - ref[ibl[i]];
02466 }
02467 avg = &(p_avg[ig]);
02468 mollify(avg,ref,lambda,force,hgs,icnt,rmass,ial,ibl,refab);
02469
02470 for ( i = 0; i < hgs; ++i ) {
02471 wc += outer(force[i]-f[Results::slow][ig+i],ref[i]);
02472 f[Results::slow][ig+i] = force[i];
02473 }
02474 }
02475
02476 *virial += wc;
02477 p_avg.resize(0);
02478 }
02479
02480 void HomePatch::checkpoint(void) {
02481 FullAtomList tmp_a(&atom); checkpoint_atom = tmp_a;
02482 checkpoint_lattice = lattice;
02483
02484
02485 #if NAMD_SeparateWaters != 0
02486 checkpoint_numWaterAtoms = numWaterAtoms;
02487 #endif
02488 }
02489
02490 void HomePatch::revert(void) {
02491 atomMapper->unregisterIDsFullAtom(atom.begin(),atom.end());
02492
02493 FullAtomList tmp_a(&checkpoint_atom); atom = tmp_a;
02494 numAtoms = atom.size();
02495 lattice = checkpoint_lattice;
02496
02497
02498 #if NAMD_SeparateWaters != 0
02499 numWaterAtoms = checkpoint_numWaterAtoms;
02500 #endif
02501 }
02502
02503 void HomePatch::submitLoadStats(int timestep)
02504 {
02505 LdbCoordinator::Object()->patchLoad(patchID,numAtoms,timestep);
02506 }
02507
02508
02509 void HomePatch::doPairlistCheck()
02510 {
02511 SimParameters *simParams = Node::Object()->simParameters;
02512
02513 if ( numAtoms == 0 || ! flags.usePairlists ) {
02514 flags.pairlistTolerance = 0.;
02515 flags.maxAtomMovement = 99999.;
02516 return;
02517 }
02518
02519 int i; int n = numAtoms;
02520 CompAtom *p_i = p.begin();
02521
02522 if ( flags.savePairlists ) {
02523 flags.pairlistTolerance = doPairlistCheck_newTolerance;
02524 flags.maxAtomMovement = 0.;
02525 doPairlistCheck_newTolerance *= (1. - simParams->pairlistShrink);
02526 doPairlistCheck_lattice = lattice;
02527 doPairlistCheck_positions.resize(numAtoms);
02528 CompAtom *psave_i = doPairlistCheck_positions.begin();
02529 for ( i=0; i<n; ++i ) { psave_i[i] = p_i[i]; }
02530 return;
02531 }
02532
02533 Lattice &lattice_old = doPairlistCheck_lattice;
02534 Position center_cur = lattice.unscale(center);
02535 Position center_old = lattice_old.unscale(center);
02536 Vector center_delta = center_cur - center_old;
02537
02538
02539 BigReal max_cd = 0.;
02540 for ( i=0; i<2; ++i ) {
02541 for ( int j=0; j<2; ++j ) {
02542 for ( int k=0; k<2; ++k ) {
02543 ScaledPosition corner( i ? min.x : max.x ,
02544 j ? min.y : max.y ,
02545 k ? min.z : max.z );
02546 Vector corner_delta =
02547 lattice.unscale(corner) - lattice_old.unscale(corner);
02548 corner_delta -= center_delta;
02549 BigReal cd = corner_delta.length2();
02550 if ( cd > max_cd ) max_cd = cd;
02551 }
02552 }
02553 }
02554 max_cd = sqrt(max_cd);
02555
02556
02557 BigReal max_pd = 0.;
02558 CompAtom *p_old_i = doPairlistCheck_positions.begin();
02559 for ( i=0; i<n; ++i ) {
02560 Vector p_delta = p_i[i].position - p_old_i[i].position;
02561 p_delta -= center_delta;
02562 BigReal pd = p_delta.length2();
02563 if ( pd > max_pd ) max_pd = pd;
02564 }
02565 max_pd = sqrt(max_pd);
02566
02567 BigReal max_tol = max_pd + max_cd;
02568
02569 flags.maxAtomMovement = max_tol;
02570
02571
02572
02573 if ( max_tol > ( (1. - simParams->pairlistTrigger) *
02574 doPairlistCheck_newTolerance ) ) {
02575 doPairlistCheck_newTolerance *= (1. + simParams->pairlistGrow);
02576 }
02577
02578 if ( max_tol > doPairlistCheck_newTolerance ) {
02579 doPairlistCheck_newTolerance = max_tol / (1. - simParams->pairlistTrigger);
02580 }
02581
02582 }
02583
02584 void HomePatch::doGroupSizeCheck()
02585 {
02586 if ( ! flags.doNonbonded ) return;
02587
02588 SimParameters *simParams = Node::Object()->simParameters;
02589 BigReal hgcut = 0.5 * simParams->hgroupCutoff; hgcut *= hgcut;
02590 BigReal maxrad2 = 0.;
02591
02592 FullAtomList::iterator p_i = atom.begin();
02593 FullAtomList::iterator p_e = atom.end();
02594
02595 while ( p_i != p_e ) {
02596 const int hgs = p_i->hydrogenGroupSize;
02597 int ngs = hgs;
02598 if ( ngs > 5 ) ngs = 5;
02599 BigReal x = p_i->position.x;
02600 BigReal y = p_i->position.y;
02601 BigReal z = p_i->position.z;
02602 int i;
02603 for ( i = 1; i < ngs; ++i ) {
02604 p_i[i].nonbondedGroupSize = 0;
02605 BigReal dx = p_i[i].position.x - x;
02606 BigReal dy = p_i[i].position.y - y;
02607 BigReal dz = p_i[i].position.z - z;
02608 BigReal r2 = dx * dx + dy * dy + dz * dz;
02609 if ( r2 > hgcut ) break;
02610 else if ( r2 > maxrad2 ) maxrad2 = r2;
02611 }
02612 p_i->nonbondedGroupSize = i;
02613 for ( ; i < hgs; ++i ) {
02614 p_i[i].nonbondedGroupSize = 1;
02615 }
02616 p_i += hgs;
02617 }
02618
02619 flags.maxGroupRadius = sqrt(maxrad2);
02620
02621 }
02622
02623 void HomePatch::doMarginCheck()
02624 {
02625 SimParameters *simParams = Node::Object()->simParameters;
02626
02627 BigReal sysdima = lattice.a_r().unit() * lattice.a();
02628 BigReal sysdimb = lattice.b_r().unit() * lattice.b();
02629 BigReal sysdimc = lattice.c_r().unit() * lattice.c();
02630
02631 BigReal minSize = simParams->patchDimension - simParams->margin;
02632
02633 if ( ( aAwayDist*sysdima < minSize*0.9999 ) ||
02634 ( bAwayDist*sysdimb < minSize*0.9999 ) ||
02635 ( cAwayDist*sysdimc < minSize*0.9999 ) ) {
02636
02637 NAMD_die("Periodic cell has become too small for original patch grid!\n"
02638 "Possible solutions are to restart from a recent checkpoint,\n"
02639 "increase margin, or disable useFlexibleCell for liquid simulation.");
02640 }
02641
02642 BigReal cutoff = simParams->cutoff;
02643
02644 BigReal margina = 0.5 * ( aAwayDist - cutoff / sysdima );
02645 BigReal marginb = 0.5 * ( bAwayDist - cutoff / sysdimb );
02646 BigReal marginc = 0.5 * ( cAwayDist - cutoff / sysdimc );
02647
02648 if ( (margina < -0.0001) || (marginb < -0.0001) || (marginc < -0.0001) ) {
02649 NAMD_die("Periodic cell has become too small for original patch grid!\n"
02650 "There are probably many margin violations already reported.\n"
02651 "Possible solutions are to restart from a recent checkpoint,\n"
02652 "increase margin, or disable useFlexibleCell for liquid simulation.");
02653 }
02654
02655 BigReal minx = min.x - margina;
02656 BigReal miny = min.y - marginb;
02657 BigReal minz = min.z - marginc;
02658 BigReal maxx = max.x + margina;
02659 BigReal maxy = max.y + marginb;
02660 BigReal maxz = max.z + marginc;
02661
02662 int xdev, ydev, zdev;
02663 int problemCount = 0;
02664
02665 FullAtomList::iterator p_i = atom.begin();
02666 FullAtomList::iterator p_e = atom.end();
02667 for ( ; p_i != p_e; ++p_i ) {
02668
02669 ScaledPosition s = lattice.scale(p_i->position);
02670
02671
02672 if (s.x < minx) xdev = 0;
02673 else if (maxx <= s.x) xdev = 2;
02674 else xdev = 1;
02675
02676 if (s.y < miny) ydev = 0;
02677 else if (maxy <= s.y) ydev = 2;
02678 else ydev = 1;
02679
02680 if (s.z < minz) zdev = 0;
02681 else if (maxz <= s.z) zdev = 2;
02682 else zdev = 1;
02683
02684 if (mInfo[xdev][ydev][zdev]) {
02685 ++problemCount;
02686 }
02687
02688 }
02689
02690 marginViolations = problemCount;
02691
02692
02693
02694
02695
02696 }
02697
02698
02699 void
02700 HomePatch::doAtomMigration()
02701 {
02702 int i;
02703
02704 for (i=0; i<numNeighbors; i++) {
02705 realInfo[i].mList.resize(0);
02706 }
02707
02708
02709 atomMapper->unregisterIDsFullAtom(atom.begin(),atom.end());
02710
02711
02712
02713 BigReal minx = min.x;
02714 BigReal miny = min.y;
02715 BigReal minz = min.z;
02716 BigReal maxx = max.x;
02717 BigReal maxy = max.y;
02718 BigReal maxz = max.z;
02719
02720 int xdev, ydev, zdev;
02721 int delnum = 0;
02722
02723 FullAtomList::iterator atom_i = atom.begin();
02724 FullAtomList::iterator atom_e = atom.end();
02725
02726
02727 #if NAMD_SeparateWaters != 0
02728 FullAtomList::iterator atom_first = atom_i;
02729 int numLostWaterAtoms = 0;
02730 #endif
02731
02732 while ( atom_i != atom_e ) {
02733 if ( atom_i->migrationGroupSize ) {
02734 Position pos = atom_i->position;
02735 if ( atom_i->migrationGroupSize != atom_i->hydrogenGroupSize ) {
02736 int mgs = atom_i->migrationGroupSize;
02737 int c = 1;
02738 for ( int j=atom_i->hydrogenGroupSize; j<mgs;
02739 j+=(atom_i+j)->hydrogenGroupSize ) {
02740 pos += (atom_i+j)->position;
02741 ++c;
02742 }
02743 pos *= 1./c;
02744
02745 }
02746
02747 ScaledPosition s = lattice.scale(pos);
02748
02749
02750 if (s.x < minx) xdev = 0;
02751 else if (maxx <= s.x) xdev = 2;
02752 else xdev = 1;
02753
02754 if (s.y < miny) ydev = 0;
02755 else if (maxy <= s.y) ydev = 2;
02756 else ydev = 1;
02757
02758 if (s.z < minz) zdev = 0;
02759 else if (maxz <= s.z) zdev = 2;
02760 else zdev = 1;
02761
02762 }
02763
02764 if (mInfo[xdev][ydev][zdev]) {
02765
02766
02767
02768 MigrationList &mCur = mInfo[xdev][ydev][zdev]->mList;
02769 DebugM(3,"Migrating atom " << atom_i->id << " from patch "
02770 << patchID << " with position " << atom_i->position << "\n");
02771 mCur.add(*atom_i);
02772
02773 ++delnum;
02774
02775
02776
02777 #if NAMD_SeparateWaters != 0
02778
02779
02780
02781
02782
02783
02784
02785 int atomIndex = atom_i - atom_first;
02786 if (atomIndex < numWaterAtoms)
02787 numLostWaterAtoms++;
02788 #endif
02789
02790
02791 } else {
02792
02793 if ( delnum ) { *(atom_i-delnum) = *atom_i; }
02794
02795 }
02796
02797 ++atom_i;
02798 }
02799
02800
02801 #if NAMD_SeparateWaters != 0
02802 numWaterAtoms -= numLostWaterAtoms;
02803 #endif
02804
02805
02806 int delpos = numAtoms - delnum;
02807 DebugM(4,"numAtoms " << numAtoms << " deleted " << delnum << "\n");
02808 atom.del(delpos,delnum);
02809
02810 numAtoms = atom.size();
02811
02812 PatchMgr::Object()->sendMigrationMsgs(patchID, realInfo, numNeighbors);
02813
02814
02815 inMigration = true;
02816
02817
02818 for (i=0; i<numMlBuf; i++) {
02819 DebugM(1, "Draining migration buffer ("<<i<<","<<numMlBuf<<")\n");
02820 depositMigration(msgbuf[i]);
02821 }
02822 numMlBuf = 0;
02823
02824 if (!allMigrationIn) {
02825 DebugM(3,"All Migrations NOT in, we are suspending patch "<<patchID<<"\n");
02826 migrationSuspended = true;
02827 sequencer->suspend();
02828 migrationSuspended = false;
02829 }
02830 allMigrationIn = false;
02831
02832 inMigration = false;
02833 marginViolations = 0;
02834 }
02835
02836 void
02837 HomePatch::depositMigration(MigrateAtomsMsg *msg)
02838 {
02839
02840 if (!inMigration) {
02841
02842 msgbuf[numMlBuf++] = msg;
02843 return;
02844 }
02845
02846
02847
02848 #if NAMD_SeparateWaters != 0
02849
02850
02851
02852
02853
02854
02855 mergeAtomList(msg->migrationList);
02856
02857
02858 #else
02859
02860
02861
02862
02863 {
02864 MigrationList &migrationList = msg->migrationList;
02865 MigrationList::iterator mi;
02866 Transform mother_transform;
02867 for (mi = migrationList.begin(); mi != migrationList.end(); mi++) {
02868 DebugM(1,"Migrating atom " << mi->id << " to patch "
02869 << patchID << " with position " << mi->position << "\n");
02870 if ( mi->migrationGroupSize ) {
02871 if ( mi->migrationGroupSize != mi->hydrogenGroupSize ) {
02872 Position pos = mi->position;
02873 int mgs = mi->migrationGroupSize;
02874 int c = 1;
02875 for ( int j=mi->hydrogenGroupSize; j<mgs;
02876 j+=(mi+j)->hydrogenGroupSize ) {
02877 pos += (mi+j)->position;
02878 ++c;
02879 }
02880 pos *= 1./c;
02881
02882 mother_transform = mi->transform;
02883 pos = lattice.nearest(pos,center,&mother_transform);
02884 mi->position = lattice.reverse_transform(mi->position,mi->transform);
02885 mi->position = lattice.apply_transform(mi->position,mother_transform);
02886 mi->transform = mother_transform;
02887 } else {
02888 mi->position = lattice.nearest(mi->position,center,&(mi->transform));
02889 mother_transform = mi->transform;
02890 }
02891 } else {
02892 mi->position = lattice.reverse_transform(mi->position,mi->transform);
02893 mi->position = lattice.apply_transform(mi->position,mother_transform);
02894 mi->transform = mother_transform;
02895 }
02896 atom.add(*mi);
02897 }
02898 }
02899
02900
02901 #endif // if (NAMD_SeparateWaters != 0)
02902
02903
02904 numAtoms = atom.size();
02905 delete msg;
02906
02907 DebugM(3,"Counter on " << patchID << " = " << patchMigrationCounter << "\n");
02908 if (!--patchMigrationCounter) {
02909 DebugM(3,"All Migrations are in for patch "<<patchID<<"\n");
02910 allMigrationIn = true;
02911 patchMigrationCounter = numNeighbors;
02912 if (migrationSuspended) {
02913 DebugM(3,"patch "<<patchID<<" is being awakened\n");
02914 migrationSuspended = false;
02915 sequencer->awaken();
02916 return;
02917 }
02918 else {
02919 DebugM(3,"patch "<<patchID<<" was not suspended\n");
02920 }
02921 }
02922 }
02923
02924
02925
02926
02927 #if NAMD_SeparateWaters != 0
02928
02929
02930
02931
02932 void HomePatch::separateAtoms() {
02933
02934
02935
02936
02937
02938
02939
02940
02941
02942
02943 if (atom.size() < 0) return;
02944
02945
02946 tempAtom.resize(numAtoms);
02947
02948
02949 int wathgsize = 3;
02950 if (simParams->watmodel == WAT_TIP4) wathgsize = 4;
02951 else if (simParams->watmodel == WAT_SWM4) wathgsize = 5;
02952
02953
02954 int i = 0;
02955 int waterIndex = 0;
02956 int nonWaterIndex = 0;
02957 while (i < numAtoms) {
02958
02959 FullAtom &atom_i = atom[i];
02960 Mass mass = atom_i.mass;
02961 int hgs = atom_i.hydrogenGroupSize;
02962
02963 if (IS_HYDROGEN_GROUP_WATER(hgs, mass)) {
02964
02965
02966 if (waterIndex != i) {
02967 atom[waterIndex ] = atom[i ];
02968 atom[waterIndex + 1] = atom[i + 1];
02969 atom[waterIndex + 2] = atom[i + 2];
02970 if (wathgsize > 3) atom[waterIndex + 3] = atom[i + 3];
02971 if (wathgsize > 4) atom[waterIndex + 4] = atom[i + 4];
02972
02973 }
02974
02975 waterIndex += wathgsize;
02976 i += wathgsize;
02977
02978 } else {
02979
02980
02981 for (int j = 0; j < hgs; j++)
02982 tempAtom[nonWaterIndex + j] = atom[i + j];
02983
02984 nonWaterIndex += hgs;
02985 i += hgs;
02986 }
02987
02988 }
02989
02990
02991
02992
02993
02994
02995
02996 for (i = 0; i < nonWaterIndex; i++)
02997 atom[waterIndex + i] = tempAtom[i];
02998
02999
03000 numWaterAtoms = waterIndex;
03001 }
03002
03003
03004
03005
03006
03007
03008
03009 void HomePatch::mergeAtomList(FullAtomList &al) {
03010
03011
03012 if (al.size() <= 0) return;
03013
03014 const int orig_atomSize = atom.size();
03015 const int orig_alSize = al.size();
03016
03017
03018 atom.resize(orig_atomSize + orig_alSize);
03019
03020
03021 #if 0 // version where non-waters are moved to scratch first
03022
03023
03024
03025
03026
03027
03028
03029
03030
03031
03032
03033
03034
03035
03036
03037
03038
03039
03040
03041
03042
03043 const int orig_atom_numNonWaters = orig_atomSize - numWaterAtoms;
03044 tempAtom.resize(orig_atom_numNonWaters + al.size());
03045 for (int i = 0; i < orig_atom_numNonWaters; i++)
03046 tempAtom[i] = atom[numWaterAtoms + i];
03047
03048
03049
03050 int atom_waterIndex = numWaterAtoms;
03051 int atom_nonWaterIndex = orig_atom_numNonWaters;
03052 int i = 0;
03053 while (i < orig_alSize) {
03054
03055 FullAtom &atom_i = al[i];
03056 int hgs = atom_i.hydrogenGroupSize;
03057 if ( hgs != atom_i.migrationGroupSize ) {
03058 NAMD_bug("HomePatch::mergeAtomList() not updated for migration groups!");
03059 }
03060 Mass mass = atom_i.mass;
03061
03062 if (IS_HYDROGEN_GROUP_WATER(hgs, mass)) {
03063
03064
03065
03066
03067 al[i].position = lattice.nearest(al[i].position, center, &(al[i].transform));
03068 Transform mother_transform = al[i].transform;
03069
03070
03071 al[i+1].position = lattice.reverse_transform(al[i+1].position, al[i+1].transform);
03072 al[i+1].position = lattice.apply_transform(al[i+1].position, mother_transform);
03073 al[i+1].transform = mother_transform;
03074
03075
03076 al[i+2].position = lattice.reverse_transform(al[i+2].position, al[i+2].transform);
03077 al[i+2].position = lattice.apply_transform(al[i+2].position, mother_transform);
03078 al[i+2].transform = mother_transform;
03079
03080
03081 atom[atom_waterIndex ] = al[i ];
03082 atom[atom_waterIndex + 1] = al[i + 1];
03083 atom[atom_waterIndex + 2] = al[i + 2];
03084
03085 atom_waterIndex += 3;
03086 i += 3;
03087
03088 } else {
03089
03090
03091
03092
03093 al[i].position = lattice.nearest(al[i].position, center, &(al[i].transform));
03094 Transform mother_transform = al[i].transform;
03095
03096
03097 for (int j = 1; j < hgs; j++) {
03098 al[i+j].position = lattice.reverse_transform(al[i+j].position, al[i+j].transform);
03099 al[i+j].position = lattice.apply_transform(al[i+j].position, mother_transform);
03100 al[i+j].transform = mother_transform;
03101 }
03102
03103
03104 for (int j = 0; j < hgs; j++)
03105 tempAtom[atom_nonWaterIndex + j] = al[i + j];
03106
03107 atom_nonWaterIndex += hgs;
03108 i += hgs;
03109 }
03110
03111 }
03112
03113
03114 for (int i = 0; i < atom_nonWaterIndex; i++)
03115 atom[atom_waterIndex + i] = tempAtom[i];
03116
03117
03118 numWaterAtoms = atom_waterIndex;
03119 numAtoms = atom.size();
03120
03121
03122 #else
03123
03124
03125
03126
03127
03128
03129
03130
03131
03132
03133
03134
03135
03136
03137 int wathgsize = 3;
03138 if (simParams->watmodel == WAT_TIP4) wathgsize = 4;
03139 else if (simParams->watmodel == WAT_SWM4) wathgsize = 5;
03140
03141
03142 int al_numWaterAtoms = 0;
03143 int i = 0;
03144 while (i < orig_alSize) {
03145
03146 FullAtom &atom_i = al[i];
03147 int hgs = atom_i.hydrogenGroupSize;
03148 Mass mass = atom_i.mass;
03149
03150 if (IS_HYDROGEN_GROUP_WATER(hgs, mass)) {
03151 al_numWaterAtoms += wathgsize;
03152 }
03153
03154 i += hgs;
03155 }
03156
03157
03158
03159 if (al_numWaterAtoms > 0) {
03160 for (i = orig_atomSize - 1; i >= numWaterAtoms; i--) {
03161 atom[i + al_numWaterAtoms] = atom[i];
03162 }
03163 }
03164
03165
03166
03167
03168 int atom_waterIndex = numWaterAtoms;
03169 int atom_nonWaterIndex = orig_atomSize + al_numWaterAtoms;
03170 i = 0;
03171 while (i < orig_alSize) {
03172
03173 FullAtom &atom_i = al[i];
03174 int hgs = atom_i.hydrogenGroupSize;
03175 if ( hgs != atom_i.migrationGroupSize ) {
03176 NAMD_bug("HomePatch::mergeAtomList() not updated for migration groups!");
03177 }
03178 Mass mass = atom_i.mass;
03179
03180 if (IS_HYDROGEN_GROUP_WATER(hgs, mass)) {
03181
03182
03183
03184
03185 al[i].position = lattice.nearest(al[i].position, center, &(al[i].transform));
03186 Transform mother_transform = al[i].transform;
03187
03188
03189 al[i+1].position = lattice.reverse_transform(al[i+1].position, al[i+1].transform);
03190 al[i+1].position = lattice.apply_transform(al[i+1].position, mother_transform);
03191 al[i+1].transform = mother_transform;
03192
03193
03194 al[i+2].position = lattice.reverse_transform(al[i+2].position, al[i+2].transform);
03195 al[i+2].position = lattice.apply_transform(al[i+2].position, mother_transform);
03196 al[i+2].transform = mother_transform;
03197
03198
03199 atom[atom_waterIndex ] = al[i ];
03200 atom[atom_waterIndex + 1] = al[i + 1];
03201 atom[atom_waterIndex + 2] = al[i + 2];
03202
03203 if (wathgsize > 3) atom[atom_waterIndex + 3] = al[i + 3];
03204
03205 atom_waterIndex += wathgsize;
03206 i += wathgsize;
03207
03208 } else {
03209
03210
03211
03212
03213 al[i].position = lattice.nearest(al[i].position, center, &(al[i].transform));
03214 Transform mother_transform = al[i].transform;
03215
03216
03217 for (int j = 1; j < hgs; j++) {
03218 al[i+j].position = lattice.reverse_transform(al[i+j].position, al[i+j].transform);
03219 al[i+j].position = lattice.apply_transform(al[i+j].position, mother_transform);
03220 al[i+j].transform = mother_transform;
03221 }
03222
03223
03224 for (int j = 0; j < hgs; j++)
03225 atom[atom_nonWaterIndex + j] = al[i + j];
03226
03227 atom_nonWaterIndex += hgs;
03228 i += hgs;
03229 }
03230
03231 }
03232
03233
03234 numWaterAtoms = atom_waterIndex;
03235 numAtoms = atom_nonWaterIndex;
03236
03237 #endif
03238 }
03239
03240 #endif
03241
03242
03243
03244 inline void lubksb(HGMatrixBigReal &a, int n, HGArrayInt &indx,
03245 HGArrayBigReal &b)
03246 {
03247 int i,ii=-1,ip,j;
03248 double sum;
03249
03250 for (i=0;i<n;i++) {
03251 ip=indx[i];
03252 sum=b[ip];
03253 b[ip]=b[i];
03254 if (ii >= 0)
03255 for (j=ii;j<i;j++) sum -= a[i][j]*b[j];
03256 else if (sum) ii=i;
03257 b[i]=sum;
03258 }
03259 for (i=n-1;i>=0;i--) {
03260 sum=b[i];
03261 for (j=i+1;j<n;j++) sum -= a[i][j]*b[j];
03262 b[i]=sum/a[i][i];
03263 }
03264 }
03265
03266
03267 inline void ludcmp(HGMatrixBigReal &a, int n, HGArrayInt &indx, BigReal *d)
03268 {
03269
03270 int i,imax,j,k;
03271 double big,dum,sum,temp;
03272 HGArrayBigReal vv;
03273 *d=1.0;
03274 for (i=0;i<n;i++) {
03275 big=0.0;
03276 for (j=0;j<n;j++)
03277 if ((temp=fabs(a[i][j])) > big) big=temp;
03278 if (big == 0.0) NAMD_die("Singular matrix in routine ludcmp\n");
03279 vv[i]=1.0/big;
03280 }
03281 for (j=0;j<n;j++) {
03282 for (i=0;i<j;i++) {
03283 sum=a[i][j];
03284 for (k=0;k<i;k++) sum -= a[i][k]*a[k][j];
03285 a[i][j]=sum;
03286 }
03287 big=0.0;
03288 for (i=j;i<n;i++) {
03289 sum=a[i][j];
03290 for (k=0;k<j;k++)
03291 sum -= a[i][k]*a[k][j];
03292 a[i][j]=sum;
03293 if ( (dum=vv[i]*fabs(sum)) >= big) {
03294 big=dum;
03295 imax=i;
03296 }
03297 }
03298 if (j != imax) {
03299 for (k=0;k<n;k++) {
03300 dum=a[imax][k];
03301 a[imax][k]=a[j][k];
03302 a[j][k]=dum;
03303 }
03304 *d = -(*d);
03305 vv[imax]=vv[j];
03306 }
03307 indx[j]=imax;
03308 if (a[j][j] == 0.0) a[j][j]=TINY;
03309 if (j != n-1) {
03310 dum=1.0/(a[j][j]);
03311 for (i=j+1;i<n;i++) a[i][j] *= dum;
03312 }
03313 }
03314 }
03315
03316
03317 inline void G_q(const HGArrayVector &refab, HGMatrixVector &gqij,
03318 const int n, const int m, const HGArrayInt &ial, const HGArrayInt &ibl) {
03319 int i;
03320
03321 for(i=0;i<m;i++) {
03322 gqij[i][ial[i]]=2.0*refab[i];
03323 gqij[i][ibl[i]]=-gqij[i][ial[i]];
03324 }
03325 };
03326
03327
03328
03329 int average(CompAtom *qtilde,const HGArrayVector &q,BigReal *lambda,const int n,const int m, const HGArrayBigReal &imass, const HGArrayBigReal &length2, const HGArrayInt &ial, const HGArrayInt &ibl, const HGArrayVector &refab, const BigReal tolf, const int ntrial) {
03330
03331
03332
03333
03334
03335
03336
03337
03338
03339
03340
03341
03342
03343 int k,k1,i,j;
03344 BigReal errx,errf,d,tolx;
03345
03346 HGArrayInt indx;
03347 HGArrayBigReal p;
03348 HGArrayBigReal fvec;
03349 HGMatrixBigReal fjac;
03350 HGArrayVector avgab;
03351 HGMatrixVector grhs;
03352 HGMatrixVector auxrhs;
03353 HGMatrixVector glhs;
03354
03355
03356 tolx=tolf;
03357
03358
03359
03360 for (i=0;i<m;i++) {
03361 lambda[i]=0.0;
03362 }
03363
03364
03365
03366
03367 G_q(refab,grhs,n,m,ial,ibl);
03368 for (k=1;k<=ntrial;k++) {
03369
03370 HGArrayBigReal gij;
03371
03372
03373 {
03374 BigReal multiplier;
03375 HGArrayVector tmp;
03376 for (i=0;i<m;i++) {
03377 multiplier = lambda[i];
03378
03379 for (j=0;j<n;j++) {
03380 auxrhs[i][j]=multiplier*imass[j]*grhs[i][j];
03381 }
03382 }
03383 for (j=0;j<n;j++) {
03384
03385 for (i=0;i<m;i++) {
03386 tmp[j]+=auxrhs[i][j];
03387 }
03388 }
03389
03390 for (j=0;j<n;j++) {
03391 qtilde[j].position=q[j]+tmp[j];
03392 }
03393
03394 }
03395
03396 for ( i = 0; i < m; i++ ) {
03397 avgab[i] = qtilde[ial[i]].position - qtilde[ibl[i]].position;
03398 }
03399
03400
03401
03402 {
03403
03404
03405 HGMatrixVector grhs2;
03406
03407 G_q(avgab,glhs,n,m,ial,ibl);
03408 #ifdef DEBUG0
03409 iout<<iINFO << "G_q:" << std::endl<<endi;
03410 for (i=0;i<m;i++) {
03411 iout<<iINFO << glhs[i*n+0] << " " << glhs[i*n+1] << " " << glhs[i*n+2] << std::endl<<endi;
03412 }
03413 #endif
03414
03415
03416 for (j=0; j<n; j++) {
03417 for (i=0; i<m;i++) {
03418 grhs2[i][j] = grhs[i][j]*imass[j];
03419 }
03420 }
03421
03422
03423
03424 for (i=0;i<m;i++) {
03425 for (j=0;j<m;j++) {
03426 fjac[i][j] = 0;
03427 for (k1=0;k1<n;k1++) {
03428 fjac[i][j] += glhs[i][k1]*grhs2[j][k1];
03429 }
03430 }
03431 }
03432 #ifdef DEBUG0
03433 iout<<iINFO << "glhs" <<endi;
03434 for(i=0;i<9;i++) {
03435 iout<<iINFO << glhs[i] << ","<<endi;
03436 }
03437 iout<<iINFO << std::endl<<endi;
03438 for(i=0;i<9;i++) {
03439 iout<<iINFO << grhs2[i] << ","<<endi;
03440 }
03441 iout<<iINFO << std::endl<<endi;
03442 #endif
03443
03444 }
03445
03446 #ifdef DEBUG0
03447 iout<<iINFO << "Jac" << std::endl<<endi;
03448 for (i=0;i<m;i++)
03449 for (j=0;j<m;j++)
03450 iout<<iINFO << fjac[i][j] << " "<<endi;
03451 iout<< std::endl<<endi;
03452 #endif
03453
03454
03455 for (i=0;i<m;i++) {
03456 gij[i]=avgab[i]*avgab[i]-length2[i];
03457 }
03458 #ifdef DEBUG0
03459 iout<<iINFO << "G" << std::endl<<endi;
03460 iout<<iINFO << "( "<<endi;
03461 for(i=0;i<m-1;i++) {
03462 iout<<iINFO << gij[i] << ", "<<endi;
03463 }
03464 iout<<iINFO << gij[m-1] << ")" << std::endl<<endi;
03465 #endif
03466
03467 for(i=0;i<m;i++) {
03468 fvec[i] = gij[i];
03469 }
03470
03471
03472
03473 errf=0.0;
03474 for (i=0;i<m;i++) errf += fabs(fvec[i]);
03475 #ifdef DEBUG0
03476 iout<<iINFO << "errf: " << errf << std::endl<<endi;
03477 #endif
03478 if (errf <= tolf) {
03479 break;
03480 }
03481 for (i=0;i<m;i++) p[i] = -fvec[i];
03482
03483 ludcmp(fjac,m,indx,&d);
03484 lubksb(fjac,m,indx,p);
03485
03486 errx=0.0;
03487 for (i=0;i<m;i++) {
03488 errx += fabs(p[i]);
03489 }
03490 for (i=0;i<m;i++)
03491 lambda[i] += p[i];
03492
03493 #ifdef DEBUG0
03494 iout<<iINFO << "lambda:" << lambda[0]
03495 << " " << lambda[1] << " " << lambda[2] << std::endl<<endi;
03496 iout<<iINFO << "errx: " << errx << std::endl<<endi;
03497 #endif
03498 if (errx <= tolx) break;
03499 #ifdef DEBUG0
03500 iout<<iINFO << "Qtilde:" << std::endl<<endi;
03501 iout<<iINFO << qtilde[0].position << " " << qtilde[1].position << " " << qtilde[2].position << std::endl<<endi;
03502 #endif
03503 }
03504 #ifdef DEBUG
03505 iout<<iINFO << "LAMBDA:" << lambda[0] << " " << lambda[1] << " " << lambda[2] << std::endl<<endi;
03506 #endif
03507
03508 return k;
03509 }
03510
03511 void mollify(CompAtom *qtilde,const HGArrayVector &q0,const BigReal *lambda, HGArrayVector &force,const int n, const int m, const HGArrayBigReal &imass,const HGArrayInt &ial,const HGArrayInt &ibl,const HGArrayVector &refab) {
03512 int i,j,k;
03513 BigReal d;
03514 HGMatrixBigReal fjac;
03515 Vector zero(0.0,0.0,0.0);
03516
03517 HGArrayVector tmpforce;
03518 HGArrayVector tmpforce2;
03519 HGArrayVector y;
03520 HGMatrixVector grhs;
03521 HGMatrixVector glhs;
03522 HGArrayBigReal aux;
03523 HGArrayInt indx;
03524
03525 for(i=0;i<n;i++) {
03526 tmpforce[i]=imass[i]*force[i];
03527 }
03528
03529 HGMatrixVector grhs2;
03530 HGArrayVector avgab;
03531
03532 for ( i = 0; i < m; i++ ) {
03533 avgab[i] = qtilde[ial[i]].position - qtilde[ibl[i]].position;
03534 }
03535
03536 G_q(avgab,glhs,n,m,ial,ibl);
03537 G_q(refab,grhs,n,m,ial,ibl);
03538
03539 for (j=0; j<n; j++) {
03540 for (i=0; i<m;i++) {
03541 grhs2[i][j] = grhs[i][j]*imass[j];
03542 }
03543 }
03544
03545
03546
03547 for (i=0;i<m;i++) {
03548 for (j=0;j<m;j++) {
03549 fjac[j][i] = 0;
03550 for (k=0;k<n;k++) {
03551 fjac[j][i] += glhs[i][k]*grhs2[j][k];
03552 }
03553 }
03554 }
03555
03556
03557
03558
03559 for(i=0;i<m;i++) {
03560 aux[i]=0.0;
03561 for(j=0;j<n;j++) {
03562 aux[i]+=grhs[i][j]*tmpforce[j];
03563 }
03564 }
03565
03566 ludcmp(fjac,m,indx,&d);
03567 lubksb(fjac,m,indx,aux);
03568
03569 for(j=0;j<n;j++) {
03570 y[j] = zero;
03571 for(i=0;i<m;i++) {
03572 y[j] += aux[i]*glhs[i][j];
03573 }
03574 }
03575 for(i=0;i<n;i++) {
03576 y[i]=force[i]-y[i];
03577 }
03578
03579
03580 for(i=0;i<n;i++) {
03581 tmpforce2[i]=imass[i]*y[i];
03582 }
03583
03584
03585 for (i=0;i<n;i++) {
03586 tmpforce[i]=zero;
03587 }
03588
03589 for (j=0;j<m;j++) {
03590 Vector tmpf = 2.0*lambda[j]*(tmpforce2[ial[j]]-tmpforce2[ibl[j]]);
03591 tmpforce[ial[j]] += tmpf;
03592 tmpforce[ibl[j]] -= tmpf;
03593 }
03594
03595
03596 for(i=0;i<n;i++) {
03597 force[i]=tmpforce[i]+y[i];
03598 }
03599
03600 }
03601