SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
Classes | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
SSAGES::ForwardFlux Class Referenceabstract

ForwardFlux sampling method. More...

#include <ForwardFlux.h>

Inheritance diagram for SSAGES::ForwardFlux:
Inheritance graph
[legend]

Classes

class  FFSConfigID
 Nested class to store different FFS Config IDs. More...
 

Public Member Functions

 ForwardFlux (const MPI_Comm &world, const MPI_Comm &comm, double ninterfaces, std::vector< double > interfaces, unsigned int N0Target, std::vector< unsigned int > M, bool initialFluxFlag, bool saveTrajectories, unsigned int currentInterface, std::string output_directory, unsigned int frequency)
 Constructor. More...
 
void PreSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override
 Pre-simulation hook. More...
 
virtual void PostIntegration (Snapshot *snapshot, const class CVManager &cvmanager)=0
 Post-integration hook. More...
 
void PostSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override
 Post-simulation hook. More...
 
- Public Member Functions inherited from SSAGES::Method
 Method (uint frequency, const MPI_Comm &world, const MPI_Comm &comm)
 Constructor. More...
 
void SetCVMask (const std::vector< uint > &mask)
 Sets the collective variable mask.
 
virtual ~Method ()
 Destructor.
 
- Public Member Functions inherited from SSAGES::EventListener
 EventListener (uint frequency)
 Constructor. More...
 
uint GetFrequency () const
 Get frequency of event listener. More...
 
virtual ~EventListener ()
 Destructor.
 

Static Public Member Functions

static ForwardFluxBuild (const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path)
 
- Static Public Member Functions inherited from SSAGES::Method
static MethodBuildMethod (const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path)
 Build a derived method from JSON node. More...
 

Protected Member Functions

void CheckInitialStructure (const CVList &)
 Function that checks the initial structure that user provides.
 
void WriteInitialFlux ()
 Function to compute and write the initial flux.
 
void AddNewIDsToQueue ()
 
bool HasReturnedToA (double)
 Function checks if configuration has returned to A.
 
int HasCrossedInterface (double, double, unsigned int interface)
 Function checks if configuration has crossed interface specified since the last check. More...
 
void WriteFFSConfiguration (Snapshot *snapshot, FFSConfigID &ffsconfig, bool wassuccess)
 Write a file corresponding to FFSConfigID from current snapshot.
 
void ReadFFSConfiguration (Snapshot *, FFSConfigID &, bool)
 Read a file corresponding to a FFSConfigID into current snapshot.
 
void ComputeInitialFlux (Snapshot *, const CVList &)
 Compute Initial Flux.
 
virtual void CheckForInterfaceCrossings (Snapshot *, const class CVManager &)=0
 Function that checks if interfaces have been crossed (different for each FFS flavor)
 
virtual void InitializeQueue (Snapshot *, const CVList &)=0
 Initialize the Queue.
 
void ComputeTransitionProbabilities ()
 Compute the probability of going from each lambda_i to lambda_{i+1}. More...
 
void PrintQueue ()
 Print the queue, useful for debugging.
 
void PopQueueMPI (Snapshot *, const CVList &, unsigned)
 Pop the queue, do MPI so that all procs maintain the same queue.
 
void FluxBruteForce (Snapshot *, const CVList &)
 Compute the flux via brute force. More...
 
void ReconstructTrajectories (Snapshot *)
 When simulation is finished, parse through the trajectories that reached B, and reconstruct the complete trajectory from where it started at A (lambda0)
 
void ComputeCommittorProbability (Snapshot *)
 When simulation is finished, recursively parse through the trajectories that reached B or failed back to A and calculate the Commitor Probability of that state going to B (_pB)
 
void AppendTrajectoryFile (Snapshot *, std::ofstream &)
 
void OpenTrajectoryFile (std::ofstream &)
 Take the current config in snapshot and append it to the provided ofstream.
 

Protected Attributes

double _ninterfaces
 
std::vector< double > _interfaces
 FFS Interfaces.
 
bool _interfaces_increase
 Interfaces must monotonically increase (or decrease), this determines whether going to the 'next' interface will be higher values of CV, or lower ones.
 
double _cvvalue_previous
 Previous cv position, used to determine if you've crossed an interface since last time.
 
double _cvvalue
 current cv position
 
double _rate
 rate constant
 
std::vector< FFSConfigIDLambda0ConfigLibrary
 Data structure that holds a Library N0 configurations at lambda0.
 
double _N0TotalSimTime
 Total Simulation Time spent in accumulating \ _N0.
 
unsigned int _N0Target
 Number of configurations to store at lambda0, target.
 
double _fluxA0
 Flux of trajectories out of state A. Denoted PhiA0 over h_A in Allen2009.
 
std::vector< unsigned int > _M
 
std::vector< unsigned int > _A
 Number of attempts from interface i.
 
std::vector< double > _P
 Flag to determine wheter fluxA0 should be calculated, seems not using this. More...
 
std::vector< unsigned int > _S
 
std::vector< unsigned int > _N
 
bool _pop_tried_but_empty_queue
 
bool _initialFluxFlag
 if 1 compute initial flux
 
bool initializeQueueFlag
 
FFSConfigID myFFSConfigID
 The current FFSConfigID of this MPI process.
 
bool _saveTrajectories
 should the FFS trajectories be saved
 
unsigned int _nfailure_total
 
std::vector< std::vector
< double > > 
_pB
 
unsigned int _current_interface
 Current Interface.
 
std::default_random_engine _generator
 random number generator
 
std::deque< FFSConfigIDFFSConfigIDQueue
 
std::string _output_directory
 Directory of FFS output.
 
std::ofstream _trajectory_file
 file to which the current trajectory is written to
 
uint iteration_
 Method iteration counter/.
 
- Protected Attributes inherited from SSAGES::Method
mxx::comm world_
 Global MPI communicator.
 
mxx::comm comm_
 Local MPI communicator.
 
std::vector< uint > cvmask_
 Mask which identifies which CVs to act on.
 

Detailed Description

ForwardFlux sampling method.

The notation used here is drawn largely from Allen, Valeriani and Rein ten Wolde. J. Phys.: Condens. Matter (2009) 21:463102. We recommend referring to this review if the reader is unfamiliar with the method, or our variable naming conventions.

Definition at line 38 of file ForwardFlux.h.

Constructor & Destructor Documentation

SSAGES::ForwardFlux::ForwardFlux ( const MPI_Comm &  world,
const MPI_Comm &  comm,
double  ninterfaces,
std::vector< double >  interfaces,
unsigned int  N0Target,
std::vector< unsigned int >  M,
bool  initialFluxFlag,
bool  saveTrajectories,
unsigned int  currentInterface,
std::string  output_directory,
unsigned int  frequency 
)
inline

Constructor.

Parameters
worldMPI global communicator.
commMPI local communicator.
frequencyFrequency with which this method is invoked.

Create instance of Forward Flux

Definition at line 246 of file ForwardFlux.h.

References _A, _initialFluxFlag, _interfaces, _interfaces_increase, _N, _N0Target, _N0TotalSimTime, _nfailure_total, _ninterfaces, _output_directory, _P, _pop_tried_but_empty_queue, _S, Lambda0ConfigLibrary, and SSAGES::Method::world_.

251  :
252  Method(frequency, world, comm), _ninterfaces(ninterfaces), _interfaces(interfaces), _N0Target(N0Target),
253  _M(M), _initialFluxFlag(initialFluxFlag), _saveTrajectories(saveTrajectories), _current_interface(currentInterface),
254  _output_directory(output_directory), _generator(1), iteration_(0)
255  {
256  //_output_directory = "FFSoutput";
257  mkdir(_output_directory.c_str(),S_IRWXU); //how to make directory?
258 
259  //_current_interface = 0;
260  _N0TotalSimTime = 0;
261 
262  if (!_initialFluxFlag)
263  initializeQueueFlag = true;
264 
265  _A.resize(_ninterfaces);
266  _P.resize(_ninterfaces);
267  _S.resize(_ninterfaces);
268  _N.resize(_ninterfaces);
269 
270  //_N[_current_interface] = _N0Target;
271  /*_M.resize(_ninterfaces);
272 
273  //code for setting up simple simulation and debugging
274  _ninterfaces = 5;
275  int i;
276  _interfaces.resize(_ninterfaces);
277  _interfaces[0]=-1.0;
278  _interfaces[1]=-0.95;
279  _interfaces[2]=-0.8;
280  _interfaces[3]= 0;
281  _interfaces[4]= 1.0;
282 
283  _saveTrajectories = true;
284 
285  _initialFluxFlag = true;
286 
287  for(i=0;i<_ninterfaces;i++) _M[i] = 50;
288 
289  _N0Target = 100;*/
290  //_N0Target = 100;
291  _nfailure_total = 0;
292 
293 
294  //check if interfaces monotonically increase or decrease
295  bool errflag = false;
296  if (_interfaces[0]< _interfaces[1]) _interfaces_increase = true;
297  else if (_interfaces[0] > _interfaces[1]) _interfaces_increase = false;
298  else errflag = true;
299  for (unsigned int i=0;i<_ninterfaces-1;i++)
300  {
301  if((_interfaces_increase) && (_interfaces[i] >= _interfaces[i+1]))
302  errflag = true;
303  else if ((!_interfaces_increase) && (_interfaces[i] <= _interfaces[i+1]))
304  errflag = true;
305  }
306  if (errflag)
307  {
308  std::cerr << "Error! The interfaces are poorly defined. They must be monotonically increasing or decreasing and cannot equal one another! Please fix this.\n";
309  for (auto interface : _interfaces)
310  std::cerr << interface << " ";
311  std::cerr << "\n";
312  MPI_Abort(world_, EXIT_FAILURE);
313  }
314 
315  // This is to generate an artificial Lambda0ConfigLibrary, Hadi's code does this for real
316  // THIS SHOULD BE SOMEWHERE ELSE!!!
318  std::normal_distribution<double> distribution(0,1);
319  for (unsigned int i = 0; i < _N0Target ; i++)
320  {
321  Lambda0ConfigLibrary[i].l = 0;
322  Lambda0ConfigLibrary[i].n = i;
323  Lambda0ConfigLibrary[i].a = 0;
324  Lambda0ConfigLibrary[i].lprev = 0;
325  Lambda0ConfigLibrary[i].nprev = i;
326  Lambda0ConfigLibrary[i].aprev = 0;
327  //FFSConfigID ffsconfig = Lambda0ConfigLibrary[i];
328  }
329  /*// Write the dump file out
330  std::ofstream file;
331  std::string filename = _output_directory + "/l" + std::to_string(ffsconfig.l) + "-n" + std::to_string(ffsconfig.n) + ".dat";
332  file.open(filename.c_str());
333 
334  //first line gives ID of where it came from
335  file << ffsconfig.lprev << " " << ffsconfig.nprev << " " << ffsconfig.aprev << "\n";
336  //write position and velocity
337  file << "1 -1 0 0 " << distribution(_generator) << " " << distribution(_generator) << " 0\n";
338 
339  }
340  */
342  }
std::vector< unsigned int > _M
Definition: ForwardFlux.h:101
std::vector< unsigned int > _S
Definition: ForwardFlux.h:114
std::vector< unsigned int > _A
Number of attempts from interface i.
Definition: ForwardFlux.h:104
unsigned int _current_interface
Current Interface.
Definition: ForwardFlux.h:147
bool _pop_tried_but_empty_queue
Definition: ForwardFlux.h:124
std::vector< unsigned int > _N
Definition: ForwardFlux.h:119
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
uint iteration_
Method iteration counter/.
Definition: ForwardFlux.h:166
bool _initialFluxFlag
if 1 compute initial flux
Definition: ForwardFlux.h:127
Method(uint frequency, const MPI_Comm &world, const MPI_Comm &comm)
Constructor.
Definition: Method.h:61
std::string _output_directory
Directory of FFS output.
Definition: ForwardFlux.h:160
std::vector< double > _P
Flag to determine wheter fluxA0 should be calculated, seems not using this.
Definition: ForwardFlux.h:110
unsigned int _nfailure_total
Definition: ForwardFlux.h:140
unsigned int _N0Target
Number of configurations to store at lambda0, target.
Definition: ForwardFlux.h:94
bool _saveTrajectories
should the FFS trajectories be saved
Definition: ForwardFlux.h:135
std::default_random_engine _generator
random number generator
Definition: ForwardFlux.h:150
double _N0TotalSimTime
Total Simulation Time spent in accumulating \ _N0.
Definition: ForwardFlux.h:91
std::vector< double > _interfaces
FFS Interfaces.
Definition: ForwardFlux.h:73
bool _interfaces_increase
Interfaces must monotonically increase (or decrease), this determines whether going to the 'next' int...
Definition: ForwardFlux.h:76
std::vector< FFSConfigID > Lambda0ConfigLibrary
Data structure that holds a Library N0 configurations at lambda0.
Definition: ForwardFlux.h:88

Member Function Documentation

void SSAGES::ForwardFlux::AddNewIDsToQueue ( )
protected

Function that adds new FFS configurations to the Queue Different FFS flavors can have differences in this method

void SSAGES::ForwardFlux::AppendTrajectoryFile ( Snapshot snapshot,
std::ofstream &  file 
)
protected

Take the current config in snapshot and append it to the provided ofstream Current format is xyz style (including vx,vy,vz)

Definition at line 279 of file ForwardFlux.cpp.

References SSAGES::Snapshot::GetAtomIDs(), SSAGES::Snapshot::GetPositions(), and SSAGES::Snapshot::GetVelocities().

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings().

280  {
281  const auto& positions = snapshot->GetPositions();
282  const auto& velocities = snapshot->GetVelocities();
283  const auto& atomID = snapshot->GetAtomIDs();
284 
285  //first line gives number of atoms
286  file << atomID.size() << "\n\n";
287 
288  // Then write positions and velocities
289  for(size_t i = 0; i< atomID.size(); i++)
290  {
291  file<<atomID[i]<<" ";
292  file<<positions[i][0]<<" "<<positions[i][1]<<" "<<positions[i][2]<<" ";
293  file<<velocities[i][0]<<" "<<velocities[i][1]<<" "<<velocities[i][2]<<std::endl;
294  }
295  }

Here is the call graph for this function:

Here is the caller graph for this function:

ForwardFlux * SSAGES::ForwardFlux::Build ( const Json::Value &  json,
const MPI_Comm &  world,
const MPI_Comm &  comm,
const std::string &  path 
)
static

Definition at line 629 of file ForwardFlux.cpp.

References Json::Requirement::GetErrors(), Json::Requirement::HasErrors(), Json::ObjectRequirement::Parse(), and Json::ObjectRequirement::Validate().

633  {
634  ObjectRequirement validator;
635  Value schema;
636  Reader reader;
637 
638  reader.parse(JsonSchema::ForwardFluxMethod, schema);
639  validator.Parse(schema, path);
640 
641  // Validate inputs.
642  validator.Validate(json, path);
643  if(validator.HasErrors())
644  throw BuildException(validator.GetErrors());
645 
646  double ninterfaces = json.get("nInterfaces", 2).asDouble();
647  std::vector<double> interfaces;
648  for(auto& s : json["interfaces"])
649  interfaces.push_back(s.asDouble());
650 
651  std::vector<unsigned int> M;
652  for(auto& s : json["trials"])
653  M.push_back(s.asInt());
654 
655  if ((ninterfaces != interfaces.size()) || (ninterfaces != M.size()))
656  throw BuildException({"The size of \"interfaces\" and \"trials\" must be equal to \"nInterfaces\". See documentation for more information"});
657 
658  auto N0Target = json.get("N0Target", 1).asInt();
659  auto initialFluxFlag = json.get("computeInitialFlux", true).asBool();
660  auto saveTrajectories = json.get("saveTrajectories", true).asBool();
661  auto currentInterface = json.get("currentInterface", 0).asInt();
662  auto freq = json.get("frequency", 1).asInt();
663  auto flavor = json.get("flavor", "none").asString();
664  auto output_directory = json.get("outputDirectoryName", "FFSoutput").asString();
665 
666  if(flavor == "DirectForwardFlux")
667  {
668  return new DirectForwardFlux(world, comm, ninterfaces, interfaces, N0Target, M, initialFluxFlag, saveTrajectories, currentInterface, output_directory, freq);
669  }
670  else
671  {
672  throw BuildException({"Unknow flavor of forward flux. The options are \"DirectForwardFlux\""});
673  }
674  }
bool HasErrors()
Check if errors have occured.
Definition: Requirement.h:86
virtual void Parse(Value json, const std::string &path) override
Parse JSON value to generate Requirement(s).
std::vector< std::string > GetErrors()
Get list of error messages.
Definition: Requirement.h:92
Requirements on an object.
virtual void Validate(const Value &json, const std::string &path) override
Validate JSON value.

Here is the call graph for this function:

void SSAGES::ForwardFlux::ComputeTransitionProbabilities ( )
protected

Compute the probability of going from each lambda_i to lambda_{i+1}.

Using number of successes and number of trials This will need to be different for each FFS flavor

Definition at line 214 of file ForwardFlux.cpp.

215  {
216  double Ptotal = 1;
217  for (unsigned int i = 0; i < _ninterfaces -1; i++)
218  Ptotal *= _P[i];
219 
220  _rate = Ptotal*_fluxA0;
221 
222  //write file
223  std::ofstream file;
224  std::string filename = _output_directory + "/rate.dat";
225  file.open(filename.c_str());
226  if (!file) {std::cerr << "Error! Unable to write " << filename << "\n"; exit(1);}
227  file << _rate << "\n";
228  file.close();
229  }
double _rate
rate constant
Definition: ForwardFlux.h:85
std::string _output_directory
Directory of FFS output.
Definition: ForwardFlux.h:160
std::vector< double > _P
Flag to determine wheter fluxA0 should be calculated, seems not using this.
Definition: ForwardFlux.h:110
double _fluxA0
Flux of trajectories out of state A. Denoted PhiA0 over h_A in Allen2009.
Definition: ForwardFlux.h:97
void SSAGES::ForwardFlux::FluxBruteForce ( Snapshot ,
const CVList  
)
protected

Compute the flux via brute force.

Eventually this should be a new class that inherits from ForwardFlux, but for the time being I'll just hard code it This function takes the configurations at lambda0 and run them until they reach B (lambdaN) or return to A

int SSAGES::ForwardFlux::HasCrossedInterface ( double  current,
double  prev,
unsigned int  interface 
)
protected

Function checks if configuration has crossed interface specified since the last check.

Simple function, given current and previous cv position, checks if interface i has been crossed. If crossed in positive direction, return +1, if crossed in negative direction return -1, if nothing crossed return 0

Definition at line 79 of file ForwardFlux.cpp.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings().

80  {
81  double interface_location = _interfaces[i];
83  {
84  if ((prev <= interface_location) && (current >= interface_location))
85  return 1;
86  else if ((prev >= interface_location) && (current <= interface_location))
87  return -1;
88  else
89  return 0;
90  }
91  else
92  {
93  if ((prev >= interface_location) && (current <= interface_location))
94  return 1;
95  else if ((prev <= interface_location) && (current >= interface_location))
96  return -1;
97  else
98  return 0;
99  }
100  }
std::vector< double > _interfaces
FFS Interfaces.
Definition: ForwardFlux.h:73
bool _interfaces_increase
Interfaces must monotonically increase (or decrease), this determines whether going to the 'next' int...
Definition: ForwardFlux.h:76

Here is the caller graph for this function:

virtual void SSAGES::ForwardFlux::PostIntegration ( Snapshot snapshot,
const class CVManager cvmanager 
)
pure virtual

Post-integration hook.

Parameters
snapshotCurrent simulation snapshot.
cvmanagerCollective variable manager.

Implements SSAGES::Method.

Implemented in SSAGES::DirectForwardFlux.

void SSAGES::ForwardFlux::PostSimulation ( Snapshot snapshot,
const class CVManager cvmanager 
)
overridevirtual

Post-simulation hook.

Parameters
snapshotCurrent simulation snapshot.
cvmanagerCollective variable manager.

Implements SSAGES::Method.

Definition at line 48 of file ForwardFlux.cpp.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings().

49  {
50  std::cout << "Post simulation\n";
51 
54 
56  ReconstructTrajectories(snapshot);
57 
58  MPI_Abort(world_, EXIT_FAILURE); //more elegant solution?
59  }
void ReconstructTrajectories(Snapshot *)
When simulation is finished, parse through the trajectories that reached B, and reconstruct the compl...
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
void ComputeTransitionProbabilities()
Compute the probability of going from each lambda_i to lambda_{i+1}.
void ComputeCommittorProbability(Snapshot *)
When simulation is finished, recursively parse through the trajectories that reached B or failed back...
bool _saveTrajectories
should the FFS trajectories be saved
Definition: ForwardFlux.h:135

Here is the caller graph for this function:

void SSAGES::ForwardFlux::PreSimulation ( Snapshot snapshot,
const class CVManager cvmanager 
)
overridevirtual

Pre-simulation hook.

Parameters
snapshotCurrent simulation snapshot.
cvmanagerCollective variable manager.

Implements SSAGES::Method.

Definition at line 39 of file ForwardFlux.cpp.

References SSAGES::CVManager::GetCVs().

40  {
41  auto cvs = cvmanager.GetCVs(cvmask_);
42  if(cvs.size() > 1)
43  throw BuildException({"Forwardflux currently only works with one cv."});
44 
45  std::cout << "\nWARNING! MAKE SURE LAMMPS GIVES A DIFFERENT RANDOM SEED TO EACH PROCESSOR, OTHERWISE EACH FFS TRAJ WILL BE IDENTICAL!\n";
46  }
std::vector< uint > cvmask_
Mask which identifies which CVs to act on.
Definition: Method.h:50

Here is the call graph for this function:

Member Data Documentation

std::vector<unsigned int> SSAGES::ForwardFlux::_M
protected

Number of trials to attemt from each interface Note _M[0] sets the number of 'branches' for RBFFS and BGFFS?

Definition at line 101 of file ForwardFlux.h.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings(), and SSAGES::DirectForwardFlux::InitializeQueue().

std::vector<unsigned int> SSAGES::ForwardFlux::_N
protected

Current number of configurations currently stored at interface i This is somewhat redundant since _N[i] == _S[i-1], but for clarity N[0] - current number of configurations collected at lambda0 (first interface)

Definition at line 119 of file ForwardFlux.h.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings(), ForwardFlux(), and SSAGES::DirectForwardFlux::InitializeQueue().

unsigned int SSAGES::ForwardFlux::_nfailure_total
protected

Counts the total number of failures eventually if I prune, will want this to be a vector, where it stores the number of failures at each interface however in the absence of pruning, a traj can only fail at lambda0, so this is just a scalar

Definition at line 140 of file ForwardFlux.h.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings(), and ForwardFlux().

double SSAGES::ForwardFlux::_ninterfaces
protected

Number of FFS interfaces note that _ninterfaces = n+1 where n is the interface defining B

Definition at line 70 of file ForwardFlux.h.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings(), and ForwardFlux().

std::vector<double> SSAGES::ForwardFlux::_P
protected

Flag to determine wheter fluxA0 should be calculated, seems not using this.

Probability of going from lambda_{i} to lambda_{i+1}

Definition at line 110 of file ForwardFlux.h.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings(), and ForwardFlux().

std::vector<std::vector<double> > SSAGES::ForwardFlux::_pB
protected

commitor probability. The probability of a given configuration reaching B

Definition at line 144 of file ForwardFlux.h.

bool SSAGES::ForwardFlux::_pop_tried_but_empty_queue
protected

Keep track of jobs that have suceeded or failed but couldn't get reassigned a new task and must wait for the queue to get more jobs This could happen in DFFS once a job has finished but M[i] hasn't been reached (waiting on other jobs) If this is the case I call it a 'zombie job', since the job is running, but isnt doing anything useful. Its just burning cpu cycles waiting for the queue to repopulate

Definition at line 124 of file ForwardFlux.h.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings(), and ForwardFlux().

std::vector<unsigned int> SSAGES::ForwardFlux::_S
protected

Number of successes from lambda_{i} to lambda_{i+1} (might need to be 2d vector if multiple branches are used (with RBFFS)

Definition at line 114 of file ForwardFlux.h.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings(), and ForwardFlux().

std::deque<FFSConfigID> SSAGES::ForwardFlux::FFSConfigIDQueue
protected

Queue When a given processor reaches an interface, it pulls a config from this Queue to figure out what it should do next This object should be syncronized between all FFS walkers (is walker the correct terminology here?) technically this is a double-ended queue, this was mostly for debugging to allow element access of the queue (which std::queue doesn't allow). I use it like a queue though.

Definition at line 157 of file ForwardFlux.h.

Referenced by SSAGES::DirectForwardFlux::CheckForInterfaceCrossings(), and SSAGES::DirectForwardFlux::InitializeQueue().


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