changa  3.5
 All Classes Files Functions Variables Typedefs Enumerations Friends Macros Groups Pages
Public Member Functions | Public Attributes | Friends | List of all members
TreePiece Class Reference

Fundamental structure that holds particle and tree data. More...

#include <ParallelGravity.h>

Inheritance diagram for TreePiece:

Public Member Functions

void memCacheStats (const CkCallback &cb)
 
void addActiveWalk (int iAwi, TreeWalk *tw, Compute *c, Opt *o, State *s)
 
void markWalkDone ()
 Called when walk on the current TreePiece is done.
 
void finishWalk ()
 Called when walk on all TreePieces is done.
 
void markSmoothWalkDone ()
 Called when smooth walk on the current TreePiece is done.
 
void finishSmoothWalk ()
 Called when smooth walk on all TreePieces is done.
 
int getIndex ()
 
void addToNodeInterRemote (int chunk, int howmany)
 accumulate node interaction count for statistics
 
void addToParticleInterRemote (int chunk, int howmany)
 accumulate particle interaction count for statistics
 
void addToNodeInterLocal (int howmany)
 accumulate node interaction count for statistics
 
void addToParticleInterLocal (int howmany)
 accumulate particle interaction count for statistics
 
void initiatePrefetch (int chunk)
 
void startRemoteChunk ()
 Start a new remote computation upon prefetch finished.
 
int getNumParticles ()
 Return the number of particles on this TreePiece.
 
GravityParticlegetParticles ()
 Return the pointer to the particles on this TreePiece.
 
void continueStartRemoteChunk (int chunk)
 Main work of StartRemoteChunk() Schedule a TreePiece::calculateGravityRemote() then start prefetching for the next chunk.
 
void continueWrapUp ()
 
void getBucketsBeneathBounds (GenericTreeNode *&node, int &start, int &end)
 get range of bucket numbers beneath a given TreeNode. More...
 
void updateBucketState (int start, int end, int n, int chunk, State *state)
 
void updateUnfinishedBucketState (int start, int end, int n, int chunk, State *state)
 
GenericTreeNodegetStartAncestor (int current, int previous, GenericTreeNode *dflt)
 return the largest node which contains current bucket, but which does not contain previous bucket. If previous is -1 then return the root.
 
GenericTreeNodekeyToNode (const Tree::NodeKey k)
 convert a key to a node using the nodeLookupTable
 
GenericTreeNodegetBucket (int i)
 
void finishedChunk (int chunk)
 
void EwaldGPU ()
 
void EwaldGPUComplete ()
 
void quiescence ()
 
NodeLookupType & getNodeLookupTable ()
 
int getNumNodes ()
 
void deleteTree ()
 delete treenodes if allocated
 
void calculateRemoteMoments (GenericTreeNode *node)
 
void checkTree (GenericTreeNode *node)
 
bool nodeOwnership (const Tree::NodeKey nkey, int &firstOwner, int &lastOwner)
 
void initBuckets ()
 
template<class Tsmooth >
void initBucketsSmooth (Tsmooth tSmooth)
 
void smoothNextBucket ()
 
void reSmoothNextBucket ()
 
void markSmoothNextBucket ()
 
void smoothBucketComputation ()
 
void startNextBucket ()
 Start the treewalk for the next bucket among those belonging to me. The buckets are simply ordered in a vector.
 
void doAllBuckets ()
 Start a full step of bucket computation, it sends a message to trigger nextBucket() which will loop over all the buckets.
 
void reconstructNodeLookup (GenericTreeNode *node)
 
 TreePiece (CkMigrateMessage *m)
 
void setPeriodic (int nReplicas, Vector3D< cosmoType > fPeriod, int bEwald, double fEwCut, double fEwhCut, int bPeriod, int bComove, double dRhoFac)
 
void BucketEwald (GenericTreeNode *req, int nReps, double fEwCut)
 
void EwaldInit ()
 
void calculateEwald (dummyMsg *m)
 
void calculateEwaldUsingCkLoop (dummyMsg *msg, int yield_num)
 
void callBucketEwald (int id)
 
void doParallelNextBucketWork (int id, LoopParData *lpdata)
 
void initCoolingData (const CkCallback &cb)
 
void velScale (double dScale, const CkCallback &cb)
 
void loadNChilada (const std::string &filename, const double dTuFac, const CkCallback &cb)
 Load I.C. from NChilada file. More...
 
void readFloatBinary (OutputParams &params, int bParaRead, const CkCallback &cb)
 
void loadTipsy (const std::string &filename, const double dTuFac, const bool bDoublePos, const bool bDoubleVel, const CkCallback &cb)
 Load I.C. from Tipsy file. More...
 
void readTipsyArray (OutputParams &params, const CkCallback &cb)
 read a tipsy array file (binary or ascii)
 
void resetMetals (const CkCallback &cb)
 Set total metals based on Ox and Fe mass fractions.
 
void getMaxIOrds (const CkCallback &cb)
 return maximum iOrders. Like the above, but the file has already been read
 
void RestartEnergy (double dTuFac, const CkCallback &cb)
 
void findTotalMass (const CkCallback &cb)
 
void recvTotalMass (CkReductionMsg *msg)
 
void writeTipsy (Tipsy::TipsyWriter &w, const double dvFac, const double duTFac, const bool bDoublePos, const bool bDoubleVel, const int bCool)
 
void setupWrite (int iStage, u_int64_t iPrevOffset, const std::string &filename, const double dTime, const double dvFac, const double duTFac, const bool bDoublePos, const bool bDoubleVel, const int bCool, const CkCallback &cb)
 Find starting offsets and begin parallel write. More...
 
void parallelWrite (int iPass, const CkCallback &cb, const std::string &filename, const double dTime, const double dvFac, const double duTFac, const bool bDoublePos, const bool bDoubleVel, const int bCool)
 Control the parallelism in the tipsy output by breaking it up into nIOProcessor pieces. More...
 
void serialWrite (u_int64_t iPrevOffset, const std::string &filename, const double dTime, const double dvFac, const double duTFac, const bool bDoublePos, const bool bDoubleVel, const int bCool, const CkCallback &cb)
 Serial output of tipsy file. More...
 
void oneNodeWrite (int iIndex, int iOutParticles, int iOutSPH, int iOutStar, GravityParticle *particles, extraSPHData *pGas, extraStarData *pStar, int *piSPH, int *piStar, const u_int64_t iPrevOffset, const std::string &filename, const double dTime, const double dvFac, const double duTFac, const bool bDoublePos, const bool bDoubleVel, const int bCool, const CkCallback &cb)
 write out the particles I have been sent
 
void reOrder (int64_t nMaxOrder, const CkCallback &cb)
 Reorder particles for output. More...
 
void ioShuffle (CkReductionMsg *msg)
 Perform the shuffle for reOrder.
 
void ioAcceptSortedParticles (ParticleShuffleMsg *)
 Accept particles from other TreePieces once the sorting has finished.
 
void resetObjectLoad (const CkCallback &cb)
 Set the load balancing data after a restart from checkpoint.
 
void assignKeys (CkReductionMsg *m)
 After the bounding box has been found, we can assign keys to the particles.
 
void evaluateBoundaries (SFC::Key *keys, const int n, int isRefine, const CkCallback &cb)
 
void unshuffleParticles (CkReductionMsg *m)
 
void acceptSortedParticles (ParticleShuffleMsg *)
 Accept particles from other TreePieces once the sorting has finished.
 
void shuffleAfterQD ()
 
void unshuffleParticlesWoDD (const CkCallback &cb)
 
void acceptSortedParticlesFromOther (ParticleShuffleMsg *)
 
void setNumExpectedNeighborMsgs ()
 
void initORBPieces (const CkCallback &cb)
 Initialize stuff before doing ORB decomposition.
 
void initBeforeORBSend (unsigned int myCount, unsigned int myCountGas, unsigned int myCountStar, const CkCallback &cb, const CkCallback &cback)
 
void sendORBParticles ()
 
void acceptORBParticles (const GravityParticle *particles, const int n)
 
void acceptORBParticles (const GravityParticle *particles, const int n, const extraSPHData *pGas, const int nGasIn, const extraStarData *pStar, const int nStarIn)
 Accept particles from other TreePieces once the sorting has finished.
 
void finalizeBoundaries (ORBSplittersMsg *splittersMsg)
 Determine my boundaries at the end of ORB decomposition.
 
void evaluateParticleCounts (ORBSplittersMsg *splittersMsg)
 Evaluate particle counts for ORB decomposition. More...
 
void kick (int iKickRung, double dDelta[MAXRUNG+1], int bClosing, int bNeedVPred, int bGasIsothermal, double dMaxEnergy, double duDelta[MAXRUNG+1], double gammam1, double dThermalCondSatCoeff, double dMultiPhaseMaxTime, double dMultiPhaseMinTemp, double dEvapCoeff, const CkCallback &cb)
 
void drift (double dDelta, int bNeedVPred, int bGasIsothermal, double dvDelta, double duDelta, int nGrowMass, bool buildTree, double dMaxEnergy, const CkCallback &cb)
 
void initAccel (int iKickRung, const CkCallback &cb)
 
void applyFrameAcc (int iKickRung, Vector3D< double > frameAcc, const CkCallback &cb)
 
void externalGravity (int activeRung, const ExternalGravity exGrav, const CkCallback &cb)
 Apply an external gravitational force. More...
 
void adjust (int iKickRung, int bEpsAccStep, int bGravStep, int bSphStep, int bViscosityLimitdt, double dEta, double dEtaCourant, double dEtauDot, double dDiffCoeff, double dEtaDiffusion, double dDelta, double dAccFac, double dCosmoFac, double dhMinOverSoft, double dResolveJeans, int bDoGas, const CkCallback &cb)
 
void truncateRung (int iCurrMaxRung, const CkCallback &cb)
 Truncate the highest rung. More...
 
void rungStats (const CkCallback &cb)
 
void countActive (int activeRung, const CkCallback &cb)
 
void countType (int iType, const CkCallback &cb)
 count total number of particles of given type
 
void outputBlackHoles (const std::string &pszFileName, double dvFac, long lFPos, const CkCallback &cb)
 processor specific routine to output black hole orbits. More...
 
void SetSink (double dSinkMassMin, const CkCallback &cb)
 set sink type based on formation time. More...
 
void SinkStep (int iCurrSinkRung, int iKickRung, const CkCallback &cb)
 set sink timesteps.
 
void formSinks (int bJeans, double dJConst2, int bDensity, double dDensityCut, double dTime, int iKickRung, int bSimple, const CkCallback &cb)
 Form sink particles: per processor routine.
 
void emergencyAdjust (int iRung, double dDelta, double dDeltaThresh, const CkCallback &cb)
 Look for gas particles reporting small new timesteps and move their timesteps down, appropriately adjusting predicted quantities. More...
 
void assignDomain (const CkCallback &cb)
 assign domain number to each particle for diagnostic
 
void starCenterOfMass (const CkCallback &cb)
 
void calcEnergy (const CkCallback &cb)
 
void newParticle (GravityParticle *p)
 add new particle
 
void adjustTreePointers (GenericTreeNode *node, GravityParticle *newParts)
 Adjust particlePointer attribute of all the nodes in the tree. More...
 
void colNParts (const CkCallback &cb)
 Count add/deleted particles, and compact main particle storage. More...
 
void newOrder (const NewMaxOrder *nStarts, const int n, const CkCallback &cb)
 Assign iOrders to recently added particles. More...
 
void setNParts (int64_t _nTotalSPH, int64_t _nTotalDark, int64_t _nTotalStar, const CkCallback &cb)
 Update total particle numbers. More...
 
void setSoft (const double dSoft, const CkCallback &cb)
 Set the gravitational softening on all particles. More...
 
void physicalSoft (const double dSoftMax, const double dFac, const int bSoftMaxMul, const CkCallback &cb)
 Adjust comoving softening to maintain constant physical softening.
 
void growMass (int nGrowMass, double dDeltaM, const CkCallback &cb)
 
void InitEnergy (double dTuFac, double z, double dTime, double gammam1, const CkCallback &cb)
 
void updateuDot (int activeRung, double duDelta[MAXRUNG+1], double dStartTime[MAXRUNG+1], int bCool, int bAll, int bUpdateState, double gammam1, double dResolveJeans, const CkCallback &cb)
 Update the cooling rate (uDot) More...
 
void ballMax (int activeRung, double dFac, const CkCallback &cb)
 
void sphViscosityLimiter (int bOn, int activeRung, const CkCallback &cb)
 
void getAdiabaticGasPressure (double gamma, double gammam1, double dTuFac, double dThermalCondCoeff, double dThermalCond2Coeff, double dThermalCondSatCoeff, double dThermalCond2SatCoeff, double dEvapMinTemp, double dDtCourantFac, double dResolveJeans, const CkCallback &cb)
 
void getCoolingGasPressure (double gamma, double gammam1, double dThermalCondCoeff, double dThermalCond2Coeff, double dThermalCondSatCoeff, double dThermalCond2SatCoeff, double dEvapMinTemp, double dDtCourantFac, double dResolveJeans, const CkCallback &cb)
 
COOLCool ()
 
void initRand (int iRand, const CkCallback &cb)
 initialize random seed for star formation
 
void FormStars (Stfm param, double dTime, double dDelta, double dCosmoFac, const CkCallback &cb)
 
void flushStarLog (const CkCallback &cb)
 flush starlog table to disk.
 
void Feedback (const Fdbk &fb, double dTime, double dDelta, const CkCallback &cb)
 
void massMetalsEnergyCheck (int bPreDist, const CkCallback &cb)
 total feedback quantities for conservation check More...
 
void SetTypeFromFileSweep (int iSetMask, char *file, struct SortStruct *ss, int nss, int *pniOrder, int *pnSet)
 Read file of iOrders to set particle type. More...
 
void setTypeFromFile (int iSetMask, char *file, const CkCallback &cb)
 
void getCOM (const CkCallback &cb, int bLiveViz)
 
void getCOMByType (int iType, const CkCallback &cb, int bLiveViz)
 
void DumpFrame (InDumpFrame in, const CkCallback &cb, int liveVizDump)
 
void liveVizDumpFrameInit (liveVizRequestMsg *msg)
 
void setProjections (int bOn)
 
void buildTree (int bucketSize, const CkCallback &cb)
 Charm entry point to build the tree (called by Main). More...
 
void startOctTreeBuild (CkReductionMsg *m)
 Real tree build, independent of other TreePieces. More...
 
void recvBoundary (SFC::Key key, NborDir dir)
 
void recvdBoundaries (CkReductionMsg *m)
 
void startORBTreeBuild (CkReductionMsg *m)
 
OrientedBox< float > constructBoundingBox (GenericTreeNode *node, int level, int numChild)
 
void buildORBTree (GenericTreeNode *node, int level)
 
bool sendFillReqNodeWhenNull (CkCacheRequestMsg< KeyType > *msg)
 When the node is found to be NULL, forward the request. More...
 
void requestRemoteMoments (const Tree::NodeKey key, int sender)
 Request the moments for this node. More...
 
void receiveRemoteMoments (const Tree::NodeKey key, Tree::NodeType type, int firstParticle, int numParticles, int remIdx, const MultipoleMoments &moments, const OrientedBox< double > &box, const OrientedBox< double > &boxBall, const unsigned int iParticleTypes, const int64_t nSPH)
 response from requestRemoteMoments
 
void calculateGravityLocal ()
 This function could be replaced by the doAllBuckets() call. More...
 
void commenceCalculateGravityLocal ()
 
void calculateGravityRemote (ComputeChunkMsg *msg)
 
void executeCkLoopParallelization (LoopParData *lpdata, int startbucket, int yield_num, int chunkNum, State *gravityState)
 
int doBookKeepingForTargetActive (int curbucket, int end, int chunkNum, bool updatestate, State *gravityState)
 
int doBookKeepingForTargetInactive (int chunkNum, bool updatestate, State *gravityState)
 
void calculateSmoothLocal ()
 As above but for the Smooth operation.
 
void calculateReSmoothLocal ()
 
void calculateMarkSmoothLocal ()
 
void nextBucketSmooth (dummyMsg *msg)
 
void nextBucketReSmooth (dummyMsg *msg)
 
void nextBucketMarkSmooth (dummyMsg *msg)
 
void startGravity (int am, double myTheta, const CkCallback &cb)
 Start a tree based gravity computation. More...
 
void setupSmooth ()
 Setup utility function for all the smooths. Initializes caches.
 
void startSmooth (SmoothParams *p, int iLowhFix, int nSmooth, double dfBall2OverSoft2, const CkCallback &cb)
 
void startReSmooth (SmoothParams *p, const CkCallback &cb)
 
void startMarkSmooth (SmoothParams *p, const CkCallback &cb)
 
void finishNodeCache (const CkCallback &cb)
 We are done with the node Cache.
 
GenericTreeNoderequestNode (int remoteIndex, Tree::NodeKey lookupKey, int chunk, int reqID, int awi, void *source)
 Retrieve the remote node, goes through the cache if present.
 
void fillRequestNode (CkCacheRequestMsg< KeyType > *msg)
 Receive a request for Nodes from a remote processor, copy the data into it, and send back a message.
 
const GenericTreeNodelookupNode (Tree::NodeKey key)
 Receive the node from the cache as following a previous request which returned NULL, and continue the treewalk of the bucket which requested it with this new node. More...
 
const GravityParticlelookupParticles (int begin)
 Find the particles starting at "begin", and return a pointer to it.
 
void finishBucket (int iBucket)
 Check if we have done with the treewalk on a specific bucket, and if we have, check also if we are done with all buckets.
 
void cachedWalkBucketTree (GenericTreeNode *node, int chunk, int reqID)
 Routine which does the tree walk on non-local nodes. It is called back for every incoming node (which are those requested to the cache during previous treewalks), and continue the treewalk from where it had been interrupted. It will possibly made other remote requests.
 
void calculateForcesNode (OffsetNode node, GenericTreeNode *myNode, int level, int chunk)
 
void calculateForces (OffsetNode node, GenericTreeNode *myNode, int level, int chunk)
 
ExternalGravityParticlerequestParticles (Tree::NodeKey key, int chunk, int remoteIndex, int begin, int end, int reqID, int awi, void *source)
 
GravityParticlerequestSmoothParticles (Tree::NodeKey key, int chunk, int remoteIndex, int begin, int end, int reqID, int awi, void *source)
 
void fillRequestParticles (CkCacheRequestMsg< KeyType > *msg)
 
void fillRequestSmoothParticles (CkCacheRequestMsg< KeyType > *msg)
 
void flushSmoothParticles (CkCacheFillMsg< KeyType > *msg)
 
void processReqSmoothParticles ()
 
void getParticleInfoForLB (int64_t active_part, int64_t total_part)
 
void startlb (const CkCallback &cb, int activeRung)
 
void setTreePieceLoad (int activeRung)
 Sets the load of the TreePiece object based on the rung. More...
 
void populateSavedPhaseData (int phase, double tpload, unsigned int activeparts)
 
bool havePhaseData (int phase)
 
void savePhaseData (std::vector< double > &loads, std::vector< unsigned int > &parts_per_phase, double *shuffleloads, unsigned int *shuffleparts, int shufflelen)
 
void ResumeFromSync ()
 
void outputASCII (OutputParams &params, int bParaWrite, const CkCallback &cb)
 
void oneNodeOutVec (OutputParams &params, Vector3D< double > *avOut, int nPart, int iIndex, int bDone, const CkCallback &cb)
 
void oneNodeOutArr (OutputParams &params, double *adOut, int nPart, int iIndex, int bDone, const CkCallback &cb)
 
void oneNodeOutIntArr (OutputParams &params, int *aiOut, int nPart, int iIndex, const CkCallback &cb)
 
void outputBinaryStart (OutputParams &params, int64_t nStart, const CkCallback &cb)
 Determine start offsets for this piece based on previous pieces.
 
void outputBinary (Ck::IO::Session, OutputParams &params)
 Output a Tipsy XDR binary float array file.
 
void minmaxNCOut (OutputParams &params, const CkCallback &cb)
 
void outputStatistics (const CkCallback &cb)
 
void collectStatistics (const CkCallback &cb)
 Collect the total statistics from the various chares.
 
void nextBucket (dummyMsg *m)
 Entry method used to split the processing of all the buckets in small pieces. It calls startNextBucket() _yieldPeriod number of times, and then it returns to the scheduler after enqueuing a message for itself. After each startNextBucket() call, the state is checked for completed walks (and forces calculated, and finishBucket() is called to clean up.
 
void nextBucketUsingCkLoop (dummyMsg *m)
 
void report ()
 Write a file containing a graphviz dot graph of my tree.
 
void printTreeViz (GenericTreeNode *node, std::ostream &os)
 Print a graphviz version of A tree.
 
void printTree (GenericTreeNode *node, std::ostream &os)
 Print a text version of a tree.
 
void pup (PUP::er &p)
 Fix pup routine to handle correctly the tree
 
GenericTreeNodegetRoot ()
 
Vector3D< cosmoType > decodeOffset (int reqID)
 
void receiveNodeCallback (GenericTreeNode *node, int chunk, int reqID, int awi, void *source)
 
void receiveParticlesCallback (ExternalGravityParticle *egp, int num, int chunk, int reqID, Tree::NodeKey &remoteBucket, int awi, void *source)
 
void receiveParticlesFullCallback (GravityParticle *egp, int num, int chunk, int reqID, Tree::NodeKey &remoteBucket, int awi, void *source)
 
void doAtSync ()
 
void balanceBeforeInitialForces (const CkCallback &cb)
 
void sendRequestForNonLocalMoments (GenericTreeNode *pickedNode)
 
void mergeNonLocalRequestsDone ()
 
std::map< NodeKey,
NonLocalMomentsClientList >
::iterator 
createTreeBuildMomentsEntry (GenericTreeNode *pickedNode)
 

Public Attributes

double dStartTime
 Time read in from input file.
 
int64_t nTotalParticles
 Total Particles in the simulation.
 
int64_t nTotalSPH
 Total Gas Particles.
 
int64_t nTotalDark
 Total Dark Particles.
 
int64_t nTotalStar
 Total Star Particles.
 
bool(* compFuncPtr [3])(GravityParticle, GravityParticle)
 Array of comparison function pointers.
 
double tmpTime
 
double totalTime
 

Friends

class PrefetchCompute
 
class GravityCompute
 
class SmoothCompute
 
class KNearestSmoothCompute
 
class ReSmoothCompute
 
class MarkSmoothCompute
 
class ListCompute
 
class NearNeighborState
 
class ReNearNeighborState
 
class MarkNeighborState
 
class BottomUpTreeWalk
 
class RemoteTreeBuilder
 
class LocalTreeBuilder
 

Detailed Description

Fundamental structure that holds particle and tree data.

Member Function Documentation

void TreePiece::adjust ( int  iKickRung,
int  bEpsAccStep,
int  bGravStep,
int  bSphStep,
int  bViscosityLimitdt,
double  dEta,
double  dEtaCourant,
double  dEtauDot,
double  dDiffCoeff,
double  dEtaDiffusion,
double  dDelta,
double  dAccFac,
double  dCosmoFac,
double  dhMinOverSoft,
double  dResolveJeans,
int  bDoGas,
const CkCallback &  cb 
)

Adjust timesteps of active particles.

Parameters
iKickRungThe rung we are on.
bEpsAccStepUse sqrt(eps/acc) timestepping
bGravStepUse sqrt(r^3/GM) timestepping
bSphStepUse Courant condition
bViscosityLimitdtUse viscosity in Courant condition
dEtaFactor to use in determing timestep
dEtaCourantCourant factor to use in determing timestep
dEtauDotFactor to use in uDot based timestep
dDiffCoeffDiffusion coefficent
dEtaDiffusionFactor to use in diffusion based timestep
dDeltaBase timestep
dAccFacAcceleration scaling for cosmology
dCosmoFacCosmo scaling for Courant
dhMinOverSoftminimum smoothing parameter.
dResolveJeansmultiple of Jeans length to be resolved.
bDoGasWe are calculating gas forces.
cbCallback function reduces currrent maximum rung
void TreePiece::adjustTreePointers ( GenericTreeNode node,
GravityParticle newParts 
)

Adjust particlePointer attribute of all the nodes in the tree.

Call this if "myParticles" is about to change, but you still need the tree.

void TreePiece::buildTree ( int  bucketSize,
const CkCallback &  cb 
)

Charm entry point to build the tree (called by Main).

Overall start of building Tree.

For ORB trees, this continues on to TreePiece::startORBTreeBuild.

void TreePiece::calculateGravityLocal ( )

This function could be replaced by the doAllBuckets() call.

Entry point for the local computation: for each bucket compute the force that its particles see due to the other particles hosted in this TreePiece. The opening angle theta has already been passed through startGravity(). This function just calls doAllBuckets().

void TreePiece::calculateGravityRemote ( ComputeChunkMsg msg)

Entry point for the remote computation: for each bucket compute the force that its particles see due to the other particles NOT hosted by this TreePiece, and belonging to a subset of the global tree (specified by chunkNum).

void TreePiece::calculateRemoteMoments ( GenericTreeNode node)

Compute all the moments for the nodes that are NonLocal, so that during the tree traversal, they contain useful information to decide whether to open or not.

void TreePiece::checkTree ( GenericTreeNode node)

Checks that every particle is indeed included in its bucket node (may not be true due to truncation of the last two bits while generating the 63 bit keys.

Check that all the particles in the tree are really in their boxes. Because the keys are made of only the first 21 out of 23 bits of the floating point representation, there can be particles that are outside their box by tiny amounts. Whether this is bad is not yet known.

void TreePiece::colNParts ( const CkCallback &  cb)

Count add/deleted particles, and compact main particle storage.

Count add/deleted particles, and compact main particle storage. Extradata storage will be reclaimed on the next domain decompose. This contributes to a "set" reduction, which the main chare will iterate through to assign new iOrders.

void TreePiece::commenceCalculateGravityLocal ( )

Do some minor preparation for the local walkk then calculateGravityLocal().

void TreePiece::emergencyAdjust ( int  iRung,
double  dDelta,
double  dDeltaThresh,
const CkCallback &  cb 
)

Look for gas particles reporting small new timesteps and move their timesteps down, appropriately adjusting predicted quantities.

Parameters
iRungCurrent rung
dDeltaBase timestep for rungs.
dDeltaThreshThreshold timestep to adjust timestep.
void TreePiece::evaluateBoundaries ( SFC::Key *  keys,
const int  n,
int  skipEvery,
const CkCallback &  cb 
)

Determine my part of the sorting histograms by counting the number of my particles in each bin. This routine assumes the particles in key order. The parameter skipEvery means that every "skipEvery" bins counted, one must be skipped. When skipEvery is set, the keys are in groups of "skipEvery" size, and only splits within each group need to be evaluated. Hence the counts between the end of one group, and the start of the next group are not evaluated. This feature is used by the Oct decomposition.

find the last place I could put this splitter key in my array of particles

this tells me the number of particles between the last two splitter keys

void TreePiece::evaluateParticleCounts ( ORBSplittersMsg splittersMsg)

Evaluate particle counts for ORB decomposition.

Parameters
mA message containing splitting dimensions and splits, and the callback to contribute Counts the particles of this treepiece on each side of the splits. These counts are summed in a contribution to the specified callback.
void TreePiece::externalGravity ( int  activeRung,
const ExternalGravity  exGrav,
const CkCallback &  cb 
)

Apply an external gravitational force.

Parameters
activeRungThe rung to apply the force.
exGravParamsParameters of the external force
cbCallback function
void TreePiece::Feedback ( const Fdbk fb,
double  dTime,
double  dDelta,
const CkCallback &  cb 
)

processor specific method for stellar feedback

void TreePiece::flushSmoothParticles ( CkCacheFillMsg< KeyType > *  msg)

Combine cached copies with the originals on the treepiece. This function also decrements the count of outstanding cache accesses and does a check to see if the smooth walk is finished.

void TreePiece::FormStars ( Stfm  stfm,
double  dTime,
double  dDelta,
double  dCosmoFac,
const CkCallback &  cb 
)

processor specific method for star formation

void TreePiece::getBucketsBeneathBounds ( GenericTreeNode *&  source,
int &  start,
int &  end 
)

get range of bucket numbers beneath a given TreeNode.

Parameters
sourceGiven TreeNode
startIndex of first bucket (returned)
endIndex of last bucket (returned)
void TreePiece::initBeforeORBSend ( unsigned int  myCount,
unsigned int  myCountGas,
unsigned int  myCountStar,
const CkCallback &  cb,
const CkCallback &  cback 
)

Allocate memory for sorted particles.

Parameters
cbcallback after everything is sorted.
cbackcallback to perform now.
void TreePiece::initBuckets ( )

Initialize all the buckets for the tree walk : Eliminate this redundant copy!

Initialize all particles for gravity force calculation. This includes zeroing out the acceleration and potential.

void TreePiece::initCoolingData ( const CkCallback &  cb)

Per thread initialization

void TreePiece::initiatePrefetch ( int  chunk)

Start prefetching the specfied chunk; prefetch compute calls startRemoteChunk() once chunk prefetch is complete

Starts the prefetching of the specified chunk; once the given chunk has been completely prefetched, the prefetch compute invokes startRemoteChunk(). Chunks in this context is a division of the remote walk so that remote computation can start before all remote prefetches are finished. By default there is only one chunk: the prefetch is done all at once.

void TreePiece::loadNChilada ( const std::string &  filename,
const double  dTuFac,
const CkCallback &  cb 
)

Load I.C. from NChilada file.

Parameters
dTuFacconversion factor from temperature to internal energy
void TreePiece::loadTipsy ( const std::string &  filename,
const double  dTuFac,
const bool  bDoublePos,
const bool  bDoubleVel,
const CkCallback &  cb 
)

Load I.C. from Tipsy file.

Parameters
filenametipsy file
dTuFacconversion factor from temperature to internal energy
bDoublePosPositions are in double precision
bDoubleVelVelocities are in double precision
const GenericTreeNode * TreePiece::lookupNode ( Tree::NodeKey  key)

Receive the node from the cache as following a previous request which returned NULL, and continue the treewalk of the bucket which requested it with this new node.

Find the key in the KeyTable, and copy the node over the passed pointer

void TreePiece::massMetalsEnergyCheck ( int  bPreDist,
const CkCallback &  cb 
)

total feedback quantities for conservation check

Sums are contributed back to main chare.

Parameters
bPreDistIs this before the feedback gets distributed. In this case the "out" quantities need to be summed.
void TreePiece::newOrder ( const NewMaxOrder nStarts,
const int  n,
const CkCallback &  cb 
)

Assign iOrders to recently added particles.

Assign iOrders to recently added particles. Also insure keys are OK

bool TreePiece::nodeOwnership ( const Tree::NodeKey  nkey,
int &  firstOwner,
int &  lastOwner 
)

Given a node, check who is the first owner and the last owner of it. It assumes that there are splitters, and that there is an ordering of them across chares. It works for SFC ordering.

Determine who are all the owners of this node

Returns
true if the caller is part of the owners, false otherwise
void TreePiece::outputBlackHoles ( const std::string &  pszFileName,
double  dvFac,
long  lFPos,
const CkCallback &  cb 
)

processor specific routine to output black hole orbits.

Note this is not parallel.

void TreePiece::parallelWrite ( int  iPass,
const CkCallback &  cb,
const std::string &  filename,
const double  dTime,
const double  dvFac,
const double  duTFac,
const bool  bDoublePos,
const bool  bDoubleVel,
const int  bCool 
)

Control the parallelism in the tipsy output by breaking it up into nIOProcessor pieces.

Parameters
iPassWhat pass we are on in the parallel write. The initial call should be with "0".
void TreePiece::readFloatBinary ( OutputParams params,
int  bParaRead,
const CkCallback &  cb 
)

Generic read of binary (NChilada) array format into floating point or integer particle attribute (NOTE MISNOMER)

void TreePiece::reOrder ( int64_t  _nMaxOrder,
const CkCallback &  cb 
)

Reorder particles for output.

Parameters
_nMaxOrderMaximum iOrder used to calculate TreePiece boundaries for output.
cbCallback for when all the shuffling is done.

After determining iOrder boundaries, the number of particles in each TreePiece is calculated via a reduction.

void TreePiece::requestRemoteMoments ( const Tree::NodeKey  key,
int  sender 
)

Request the moments for this node.

entry method to obtain the moments of a node

bool TreePiece::sendFillReqNodeWhenNull ( CkCacheRequestMsg< KeyType > *  msg)

When the node is found to be NULL, forward the request.

When the node is found to be null, find the owner and forward the request.

void TreePiece::serialWrite ( u_int64_t  iPrevOffset,
const std::string &  filename,
const double  dTime,
const double  dvFac,
const double  duTFac,
const bool  bDoublePos,
const bool  bDoubleVel,
const int  bCool,
const CkCallback &  cb 
)

Serial output of tipsy file.

Send my particles to piece 0 for output.

void TreePiece::setNParts ( int64_t  _nTotalSPH,
int64_t  _nTotalDark,
int64_t  _nTotalStar,
const CkCallback &  cb 
)

Update total particle numbers.

Update total particle numbers.

void TreePiece::SetSink ( double  dSinkMassMin,
const CkCallback &  cb 
)

set sink type based on formation time.

Initial identify sinks; processor specific method.

void TreePiece::setSoft ( const double  dSoft,
const CkCallback &  cb 
)

Set the gravitational softening on all particles.

Parameters
dSoftgravitational softening
cbcallback
void TreePiece::setTreePieceLoad ( int  activeRung)

Sets the load of the TreePiece object based on the rung.

Parameters
activeRungRung to use.
void TreePiece::SetTypeFromFileSweep ( int  iSetMask,
char *  file,
struct SortStruct ss,
int  nss,
int *  pniOrder,
int *  pnSet 
)

Read file of iOrders to set particle type.

This is used by the photogenic code of dumpframe.

void TreePiece::setupWrite ( int  iStage,
u_int64_t  iPrevOffset,
const std::string &  filename,
const double  dTime,
const double  dvFac,
const double  duTFac,
const bool  bDoublePos,
const bool  bDoubleVel,
const int  bCool,
const CkCallback &  cb 
)

Find starting offsets and begin parallel write.

Perform Parallel Scan (a log(p) algorithm) to establish start of parallel writes, then do the writing. Calls writeTipsy() to do the actual writing.

void TreePiece::startGravity ( int  am,
double  myTheta,
const CkCallback &  cb 
)

Start a tree based gravity computation.

Parameters
amthe active rung for the computation
thetathe opening angle
cbthe callback to use after all the computation has finished

This method starts the tree walk and gravity calculation. It first registers with the node and particle caches. It initializes the particle acceleration by calling initBucket(). It initializes treewalk bookkeeping, and starts a prefetch walk which will eventually call a remote walk (a walk on non-local nodes). It then starts the Ewald initialization, and finally starts the local gravity walk.

void TreePiece::startOctTreeBuild ( CkReductionMsg *  m)

Real tree build, independent of other TreePieces.

Actual treebuild for each treepiece. Each treepiece begins its treebuild at the global root. During the treebuild, the nodeLookupTable, a mapping from keys to nodes, is constructed. Also, the bucket list is constructed. Constructing moments for the boundary nodes requires requesting the moments of nodes on other pieces. After all the moments are constructed, the treepiece passes its root to the DataManager via DataManager::notifyPresence.

After the local treebuild is finished DataManager::combineLocalTrees is called, and the treebuild phase is finished.

void TreePiece::startSmooth ( SmoothParams p,
int  iLowhFix,
int  nSmooth,
double  dfBall2OverSoft2,
const CkCallback &  cb 
)

Start a tree based smooth computation.

Parameters
pparameters, including the type, of the smooth.
iLowhFixtrue if a minimum h/softening is used.
nSmoothNumber of particles to smooth over.
dfBall2OverSoft2square of minimum ratio of h to softening.
cbthe callback to use after all the computation has finished SmoothParams
void TreePiece::truncateRung ( int  iCurrMaxRung,
const CkCallback &  cb 
)

Truncate the highest rung.

Parameters
iCurrMaxRungnew maximum rung.
cbcallback.
void TreePiece::unshuffleParticles ( CkReductionMsg *  m)

Once final splitter keys have been decided, I need to give my particles out to the TreePiece responsible for them

void TreePiece::updateuDot ( int  activeRung,
double  duDelta[MAXRUNG+1],
double  dStartTime[MAXRUNG+1],
int  bCool,
int  bUpdateState,
int  bAll,
double  gammam1,
double  dResolveJeans,
const CkCallback &  cb 
)

Update the cooling rate (uDot)

Parameters
activeRung(minimum) rung being updated
duDeltaarray of timesteps of length MAXRUNG+1
dStartTimearray of start times of length MAXRUNG+1
bCoolWhether cooling is on
bUpdateStateWhether the ionization factions need updating
bAllDo all rungs below activeRung
gammam1Isentropic expansion factor/adiabatic index - 1.
dResolveJeansFraction of Pressure to resolve Jeans mass (comoving)
cbCallback.

The documentation for this class was generated from the following files: