ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
(Generate patch)

Comparing branches/development/src/brains/ForceManager.cpp (file contents):
Revision 1579 by gezelter, Thu Jun 9 20:26:29 2011 UTC vs.
Revision 1581 by gezelter, Mon Jun 13 22:13:12 2011 UTC

# Line 475 | Line 475 | namespace OpenMD {
475    }
476    
477    void ForceManager::longRangeInteractions() {
478 <    // some of this initial stuff will go away:
478 >
479      Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
480      DataStorage* config = &(curSnapshot->atomData);
481      DataStorage* cgConfig = &(curSnapshot->cgData);
482    RealType* frc = config->getArrayPointer(DataStorage::dslForce);
483    RealType* pos = config->getArrayPointer(DataStorage::dslPosition);
484    RealType* trq = config->getArrayPointer(DataStorage::dslTorque);
485    RealType* A = config->getArrayPointer(DataStorage::dslAmat);
486    RealType* electroFrame = config->getArrayPointer(DataStorage::dslElectroFrame);
487    RealType* particlePot = config->getArrayPointer(DataStorage::dslParticlePot);
482  
483 <    // new stuff starts here:
484 <    
483 >    //calculate the center of mass of cutoff group
484 >
485 >    SimInfo::MoleculeIterator mi;
486 >    Molecule* mol;
487 >    Molecule::CutoffGroupIterator ci;
488 >    CutoffGroup* cg;
489 >
490 >    if(info_->getNCutoffGroups() > 0){      
491 >      for (mol = info_->beginMolecule(mi); mol != NULL;
492 >           mol = info_->nextMolecule(mi)) {
493 >        for(cg = mol->beginCutoffGroup(ci); cg != NULL;
494 >            cg = mol->nextCutoffGroup(ci)) {
495 >          cg->updateCOM();
496 >        }
497 >      }      
498 >    } else {
499 >      // center of mass of the group is the same as position of the atom  
500 >      // if cutoff group does not exist
501 >      cgConfig->position = config->position;
502 >    }
503 >
504      fDecomp_->zeroWorkArrays();
505      fDecomp_->distributeData();
506      
# Line 496 | Line 509 | namespace OpenMD {
509      RealType rgrpsq, rgrp, r2, r;
510      RealType electroMult, vdwMult;
511      RealType vij;
512 <    Vector3d fij, fg;
512 >    Vector3d fij, fg, f1;
513      tuple3<RealType, RealType, RealType> cuts;
514      RealType rCutSq;
515      bool in_switching_region;
# Line 512 | Line 525 | namespace OpenMD {
525  
526      int loopStart, loopEnd;
527  
528 +    idat.vdwMult = &vdwMult;
529 +    idat.electroMult = &electroMult;
530 +    idat.pot = &pot;
531 +    idat.vpair = &vpair;
532 +    idat.f1 = &f1;
533 +    idat.sw = &sw;
534 +
535      loopEnd = PAIR_LOOP;
536      if (info_->requiresPrepair() ) {
537        loopStart = PREPAIR_LOOP;
# Line 519 | Line 539 | namespace OpenMD {
539        loopStart = PAIR_LOOP;
540      }
541      
522
542      for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
543      
544        if (iLoop == loopStart) {
# Line 551 | Line 570 | namespace OpenMD {
570            
571            in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
572                                                       rgrp);
554
555          idat.sw = &sw;
573                
574            atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
575            atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
# Line 565 | Line 582 | namespace OpenMD {
582                   jb != atomListColumn.end(); ++jb) {              
583                atom2 = (*jb);
584                
568              cerr << "doing atoms " << atom1 << " " << atom2 << "\n";
585                if (!fDecomp_->skipAtomPair(atom1, atom2)) {
586                  
587                  vpair = 0.0;
588 +                f1 = V3Zero;
589  
590 <                cerr << "filling idat atoms " << atom1 << " " << atom2 << "\n";
574 <                idat = fDecomp_->fillInteractionData(atom1, atom2);
575 <                cerr << "done with idat\n";
590 >                fDecomp_->fillInteractionData(idat, atom1, atom2);
591                  
592                  topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
593                  vdwMult = vdwScale_[topoDist];
594                  electroMult = electrostaticScale_[topoDist];
595  
581                idat.vdwMult = &vdwMult;
582                idat.electroMult = &electroMult;
583                idat.pot = &pot;
584                idat.vpair = &vpair;
585
596                  if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
597                    idat.d = &d_grp;
598                    idat.r2 = &rgrpsq;
# Line 594 | Line 604 | namespace OpenMD {
604                    idat.r2 = &r2;
605                  }
606                  
607 <                cerr << "d = " << d << "\n";
598 <                cerr << "r2 = " << r2 << "\n";
599 <                r = sqrt( r2 );
607 >                r = sqrt( *(idat.r2) );
608                  idat.rij = &r;
609                
610                  if (iLoop == PREPAIR_LOOP) {
611                    interactionMan_->doPrePair(idat);
612                  } else {
605                  cerr << "doing doPair " << atom1 << " " << atom2 << " " << r << "\n";
613                    interactionMan_->doPair(idat);
614                    fDecomp_->unpackInteractionData(idat, atom1, atom2);
615 <                  vij += *(idat.vpair);
616 <                  fij += *(idat.f1);
617 <                  tau -= outProduct( *(idat.d), *(idat.f1));
615 >                  vij += vpair;
616 >                  fij += f1;
617 >                  tau -= outProduct( *(idat.d), f1);
618                  }
619                }
620              }
# Line 673 | Line 680 | namespace OpenMD {
680            fDecomp_->collectIntermediateData();
681  
682            for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
683 <            sdat = fDecomp_->fillSelfData(atom1);
683 >            fDecomp_->fillSelfData(sdat, atom1);
684              interactionMan_->doPreForce(sdat);
685            }
686  
# Line 695 | Line 702 | namespace OpenMD {
702               jb != skipList.end(); ++jb) {        
703      
704            atom2 = (*jb);
705 <          idat = fDecomp_->fillSkipData(atom1, atom2);
705 >          fDecomp_->fillSkipData(idat, atom1, atom2);
706            interactionMan_->doSkipCorrection(idat);
707  
708          }
# Line 705 | Line 712 | namespace OpenMD {
712      if (info_->requiresSelfCorrection()) {
713  
714        for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {          
715 <        sdat = fDecomp_->fillSelfData(atom1);
715 >        fDecomp_->fillSelfData(sdat, atom1);
716          interactionMan_->doSelfCorrection(sdat);
717        }
718  

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines