OpenAtom  Version1.5a
Atoms

Handles coordinate and force integration for Atoms, keeping a distributed coordinate cache in AtomsCache and computing updates in AtomsCompute. More...

Classes

class  AtomsCache
 
class  AtomsCompute
 

Functions

 AtomsCache::AtomsCache (int, int, Atom *, UberCollection thisInstance)
 
void AtomsCache::contributeforces ()
 Contributeforces Every proc assigned an electronic structure chare of the bead=i temper=j (i,j) uber instance that generates an atom force MUST be in (i,j) atom group and so contribute to the reduction below. More...
 
void AtomsCache::acceptAtoms (AtomMsg *)
 acceptAtoms More...
 
void AtomsCache::atomsDone (CkReductionMsg *)
 AtomsDone. More...
 
void AtomsCache::atomsDone ()
 = Needs to have each proc invoke directly doneMovingAtoms method of the GSpaceDrivers which are mapped to it. More...
 
void AtomsCache::releaseGSP ()
 = Use eesCache's registered GSPs to More...
 
 AtomsCache::~AtomsCache ()
 Destructor. More...
 
 AtomsCompute::AtomsCompute (int, int, int, int, int, int, int, double, Atom *, AtomNHC *, int nChareAtoms, UberCollection thisInstance)
 
 AtomsCompute::~AtomsCompute ()
 Destructor. More...
 
void AtomsCompute::init ()
 = Bead root is the 0th element. More...
 
void AtomsCompute::recvContribute (CkReductionMsg *)
 recvContribute Every Array member is sent all the forces at present. More...
 
void AtomsCompute::recvContributeForces (CkReductionMsg *)
 recvContributeForces Every Array member has all the forces at present. More...
 
void AtomsCompute::handleForces ()
 
void AtomsCompute::integrateAtoms ()
 Integrate atoms. More...
 
void AtomsCompute::startRealSpaceForces ()
 trigger force computation based on real forces available in each processor's chare then contribute to global group reduction of all atom forces (within this bead) -> recvContribute More...
 
void AtomsCompute::outputAtmEnergy ()
 = Atom energy output
 
void AtomsCompute::send_PIMD_Fx ()
 
void AtomsCompute::send_PIMD_Fx_and_x ()
 
void AtomsCompute::sendAtoms (double, double, double, double, int, int, int)
 This thing is sending its updated atoms to all elements. More...
 
void AtomsCompute::bcastAtomsToAtomCache ()
 For dynamics we have to update all the caches with the new coordinates. More...
 
void AtomsCompute::acceptAtoms (AtomMsg *)
 Take packed message of a chunk of the atoms with updated positions. More...
 
void AtomsCompute::atomsDone (CkReductionMsg *)
 
void AtomsCompute::accept_PIMD_CM (AtomXYZMsg *m)
 
void AtomsCompute::send_PIMD_u ()
 
void AtomsCompute::accept_PIMD_Fu (double _fxu, double _fyu, double _fzu, int atomI)
 = is broadcast to us
 
void AtomsCompute::accept_PIMD_Fu_and_u (double _fxu, double _fyu, double _fzu, double _xu, double _yu, double _zu, int atomI)
 = is broadcast to us
 
void AtomsCompute::accept_PIMD_x (double _x, double _y, double _z, int atomI)
 
void AtomsCompute::send_PIMD_x ()
 = done during initialization in 1st iteration
 
void AtomsCompute::acceptNewTemperature (double temp)
 
void AtomsCompute::accept_PIMD_u (double _ux, double _uy, double _uz, int atomI)
 = done during initialization in 1st iteration
 

Detailed Description

Handles coordinate and force integration for Atoms, keeping a distributed coordinate cache in AtomsCache and computing updates in AtomsCompute.

AtomsCache is a mostly passive structure for handling data replication of the atom coordinates and write-back collation of force contributions.

various chares use CkLocal to access the atoms: Rho_RHart RhoGHart StructureFactor RealParticlePlane ParticlePlane -> fastAtoms read coords write forces RhoRHart -> natm read -> fastAtoms read coords write forces -> iteration read RhoGHart -> natm read -> iteration read coords write forces -> fastAtoms read RealParticlePlane -> temperScreenFile* read -> iteration read -> fastAtoms read coords write forces StructureFactor -> iteration read -> fastAtoms read coords write forces InstanceController -> iteration read -> temperScreenFile* read GSpacePlane -> natm read -> fastAtoms read coords for genwave -> temperScreenFile* read

EJB: this could perhaps be better implemented in an MSA, once they shake performance bugs out of MSA, and pay attention to partitioning considerations, we'll revisit it.

CkCache could be considered here, but our use case is so simple that CkCache is wildly over engineered for our purposes. The number of atoms is always relatively small, so we have no eviction requirements. Similarly the number of clients is static and all clients need access to all the atom coordinates every iteration, so there is no segmentation and no window of opportunity for clever just in time delivery.

Nodegroup? AtomsCache is an excellent candidate for Nodegroup. All force updates are a purely additive operation, each computation will use the coordinates and its own chare local data to produce its force contribution. Force updates are handled via by reference passage into the PINY compute routines. Updates are order independent, so as long as the memory subsystem doesn't lose its mind maintaining associative additive consistency (almost the weakest consistency one could ask for in a shared memory context) then correctness is guaranteed without locks or explicit fencing. False sharing shouldn't be an issue given that each force component is in its own vector as are the coordinates. Having only one of these per node gives us a slight memory footprint reduction. More importantly it requires numprocs/node fewer update messages. If you think this gives us grief just replace nodegroup with Group and recompile, the nodegroup advantages are all implicit.

However, AtomsCache directly uses the registration in the eesCache to execute releaseGSP. So, AtomsCache's nodegroupness needs to be 1:1 with eesCache's nodegroupness, or an alternate launch scheme put in place for releaseGSP to break that dependency.

Function Documentation

void AtomsCompute::accept_PIMD_x ( double  _x,
double  _y,
double  _z,
int  atomI 
)

= CkPrintf("{%d}[%d] AtomsCompute::accept_PIMD_x iteration %d\n ", thisInstance.proxyOffset, thisIndex, *iteration);

Definition at line 1155 of file AtomsCompute.C.

void AtomsCache::acceptAtoms ( AtomMsg msg)

acceptAtoms

AtomCompute gives us a coordinate update after force integration

Definition at line 137 of file AtomsCache.C.

void AtomsCompute::acceptAtoms ( AtomMsg msg)

Take packed message of a chunk of the atoms with updated positions.

Update local copy of atoms. Update local energyGroup members. Print atom energies when we have all of them. Do file output of atoms if desired.

= CkPrintf("{%d}[%d] AtomsCompute::acceptAtoms.\n ", thisInstance.proxyOffset, CkMyPe());

= unpack atom position and velocity : For PIMD x is reconstructed from xu

= unpack energy

= Delete the message

= Copy to the fast vectors and phone home : If PIMD Transform U to X

Definition at line 896 of file AtomsCompute.C.

References AtomsCompute::outputAtmEnergy(), and AtomsCompute::send_PIMD_u().

AtomsCompute::AtomsCompute ( int  n,
int  n_nl,
int  len_nhc_,
int  iextended_on_,
int  cp_min_opt_,
int  cp_wave_opt_,
int  isokin_opt_,
double  kT_,
Atom a,
AtomNHC aNHC,
int  nChareAtoms,
UberCollection  thisInstance 
)

= parameters, options and energies

= Initial positions, forces, velocities

= A copy of the atoms for fast routines : only need to copy charge once

= Number of messages to be received when atoms are moved : simple parallel

= PIMD set up : Even if classical, this can run. It is harmless and prevents careless bugs

Definition at line 66 of file AtomsCompute.C.

void AtomsCache::atomsDone ( )

= Needs to have each proc invoke directly doneMovingAtoms method of the GSpaceDrivers which are mapped to it.

Without migration, we have that map at startup. With migration, one must write an enroll/dropout routine. All

= Increment iteration counters CkPrintf("{%d}[%d] AtomsCache::atomsDone()\n ", thisInstance.proxyOffset, CkMyPe());

Definition at line 188 of file AtomsCache.C.

References AtomsCache::releaseGSP().

Referenced by AtomsCache::atomsDone().

void AtomsCache::atomsDone ( CkReductionMsg *  msg)

AtomsDone.

everyone has the new coordinates in dynamics

Definition at line 171 of file AtomsCache.C.

References AtomsCache::atomsDone().

void AtomsCompute::bcastAtomsToAtomCache ( )

For dynamics we have to update all the caches with the new coordinates.

= Malloc the message CkPrintf("{%d}[%d] AtomsCompute::bcastAtomsToAtomCache.\n ", thisInstance.proxyOffset, thisIndex);

= pack atom positions : new for use : old for output

= Send the message to everyone in the group

FIX THIS: should use delegated section to minimize stupid anytime migration overhead.

Definition at line 838 of file AtomsCompute.C.

void AtomsCache::contributeforces ( )

Contributeforces Every proc assigned an electronic structure chare of the bead=i temper=j (i,j) uber instance that generates an atom force MUST be in (i,j) atom group and so contribute to the reduction below.

Otherwise, all the pieces of the atom forces of (i,j) will not be collected and the user will be very sad!

Definition at line 105 of file AtomsCache.C.

void AtomsCompute::handleForces ( )

Model forces : This is fine for path integral checking too

Check the forces

= if classical go on to integration, otherwise Fx -> Fu

Definition at line 317 of file AtomsCompute.C.

References AtomsCompute::integrateAtoms(), AtomsCompute::send_PIMD_Fx(), and AtomsCompute::send_PIMD_Fx_and_x().

Referenced by AtomsCompute::recvContribute(), and AtomsCompute::recvContributeForces().

void AtomsCompute::init ( )

= Bead root is the 0th element.

The AtomsCompute map will lookup the GSP map to ensure that the 0th AtomCompute is co-mapped with GSP(0,0) in each bead instance.

Bead root should be the 0th element. The AtomsCompute map will use the GSP map to ensure that the 0th AtomCompute is co-mapped with GSP(0,0) in each bead instance.

Definition at line 197 of file AtomsCompute.C.

void AtomsCompute::integrateAtoms ( )

Integrate atoms.

This is parallelized so that a subset of the atoms are computed on each processor and their results sent to AtomCompute->acceptAtoms().

Local pointers

Atom parallelization: This decompostion is done dynamically because its a pretty trivial choice and it allows us to support atom migration into and out of the QM region as long that atom migration mechanism updates natm and the atom array accordingly.

Zero some local copies energies and flags

DEBUGGING : Compute the distribution function for model

Path integral : Overwrite variables to keep integrator clean: Add the chain force to transformed forces. Scale the masses to fict bead masses This is PIMD Respa implementation friendly

Integrate the atoms : Path Integral Ready

Path integral : Rescale to the physical masses

Debug output : Not correct for PIMD which needs chain PE

Get ready for the next iteration : zero forces, outputAtmEnergy, atomsDone

Definition at line 420 of file AtomsCompute.C.

References AtomsCompute::outputAtmEnergy(), and AtomsCompute::sendAtoms().

Referenced by AtomsCompute::accept_PIMD_Fu(), AtomsCompute::accept_PIMD_Fu_and_u(), AtomsCompute::accept_PIMD_u(), and AtomsCompute::handleForces().

void AtomsCompute::recvContribute ( CkReductionMsg *  msg)

recvContribute Every Array member is sent all the forces at present.

= Local pointers

Copy out the reduction of energy and forces and nuke msg

Tuck things that can be tucked.

Definition at line 250 of file AtomsCompute.C.

References AtomsCompute::handleForces().

void AtomsCompute::recvContributeForces ( CkReductionMsg *  msg)

recvContributeForces Every Array member has all the forces at present.

= Local pointers

Definition at line 294 of file AtomsCompute.C.

References AtomsCompute::handleForces().

void AtomsCache::releaseGSP ( )

= Use eesCache's registered GSPs to

= Use the cool new data caching system to say we're done.

Definition at line 219 of file AtomsCache.C.

Referenced by AtomsCache::atomsDone().

void AtomsCompute::send_PIMD_Fx ( )

= CkPrintf("{%d}[%d] AtomsCompute::send_PIMD_fx iteration %d\n ", thisInstance.proxyOffset, CkMyPe(), *iteration); every BOC has all Fx might as well just bcast from beadroot

Definition at line 711 of file AtomsCompute.C.

Referenced by AtomsCompute::handleForces().

void AtomsCompute::send_PIMD_Fx_and_x ( )

= CkPrintf("{%d}[%d] AtomsCompute::send_PIMD_fx iteration %d\n ", thisInstance.proxyOffset, CkMyPe(), *iteration); every BOC has all Fx might as well just bcast from beadroot

Definition at line 736 of file AtomsCompute.C.

Referenced by AtomsCompute::handleForces().

void AtomsCompute::send_PIMD_u ( )

= CkPrintf("{%d}[%d] AtomsCompute::send_PIMD_u iteration %d\n ", thisInstance.proxyOffset, CkMyPe(), *iteration);

Definition at line 1079 of file AtomsCompute.C.

Referenced by AtomsCompute::acceptAtoms().

void AtomsCompute::sendAtoms ( double  eKinetic_loc,
double  eKineticNhc_loc,
double  potNhc_loc,
double  potPIMDChain_loc,
int  natmNow,
int  natmStr,
int  natmEnd 
)

This thing is sending its updated atoms to all elements.

Effectively an all-to-all implementation of what is basically an all-gather. Cost of this isn't too bad for the array version if we use a section proxy for the multicast it constrains the number of destinations fairly managably. Formerly, this degraded to nAtoms global broadcasts in the old group implementation. If you've ever wondered why network fabrics quake in fear at the approach of an openatom run, this sort of thing should give you an inkling of the abuse we dish out. The only saving grace is that the number of atoms is relatively small.

= Malloc the message CkPrintf("{%d}[%d] AtomsCompute::sendAtoms.\n ", thisInstance.proxyOffset, CkMyPe());

= pack atom positions : new for use : old for output

= pack the 3 energies

= Send the message to everyone in the group

FIX THIS: should use delegated section to minimize stupid anytime migration overhead. Also, this is really an all-gather.

Definition at line 780 of file AtomsCompute.C.

Referenced by AtomsCompute::integrateAtoms().

void AtomsCompute::startRealSpaceForces ( )

trigger force computation based on real forces available in each processor's chare then contribute to global group reduction of all atom forces (within this bead) -> recvContribute

= Atom parallelization : same as for integrate

= Get the real space atom forces plus these 4 quantities

= Everybody contributes to the reduction (see contribute forces comment)

Definition at line 602 of file AtomsCompute.C.

AtomsCache::~AtomsCache ( )

Destructor.

Definition at line 264 of file AtomsCache.C.

AtomsCompute::~AtomsCompute ( )

Destructor.

Definition at line 167 of file AtomsCompute.C.