PAPRECA hybrid off-lattice kMC/MD simulator
2.0.0 (17 September 2024)
|
Classes | |
class | Bond |
Custom bond class (not to be confused with a LAMMPS bond). More... | |
class | Event |
Parent class for PAPRECA events. More... | |
class | Reaction |
Child of PAPRECA:Event and parent of PAPRECA::BondBreak and PAPRECA::BondForm classes. More... | |
class | BondBreak |
child of PAPRECA::Reaction dedicated to bond-breaking events. More... | |
class | BondForm |
child of PAPRECA::Reaction dedicated to bond formation events. More... | |
class | Deposition |
Child of PAPRECA::Event dedicated to monoatomic or molecular adsorption. More... | |
class | Diffusion |
child of PAPRECA::Event class dedicated to diffusion events. More... | |
class | MonoatomicDesorption |
child of PAPRECA::Event dedicated to monoatomic desorption events (i.e., where a single atom is ejected from the system). More... | |
class | PredefinedReaction |
class associated with reaction events. This is a Predefined template for PAPRECA::BondBreak and PAPRECA::BondForm events. More... | |
class | PredefinedBondForm |
child class of PAPRECA::PredefinedReaction associated with bond formation events. This is a Predefined template for PAPRECA::BondForm events. More... | |
class | PredefinedDiffusionHop |
class associated with diffusion events. This is a Predefined template for PAPRECA::DiffusionHop events. More... | |
class | PredefinedDeposition |
class associated with deposition events. This is a Predefined template for PAPRECA::Deposition events. More... | |
class | PredefinedMonoatomicDesorption |
associated with monoatomic desorption events. This is a Predefined template for PAPRECA::MonoatomicDesorption events. More... | |
class | PredefinedEventsCatalog |
General class that stores ALL the predefined events in the system. You can consider this as the Predefined Catalog of Events. More... | |
class | File |
General parent file class. Any PAPRECA export file should be a child to this class. More... | |
class | Log |
Child class of File, manages papreca.log files. More... | |
class | HeightVtime |
Child class of File, manages heightVtime.log files. More... | |
class | SurfaceCoverage |
Child class of File, manages surface_coverage.log files. More... | |
class | ElementalDistribution |
Child class of File, manages distribution.log files. More... | |
class | ExecTime |
Child class of File, manages execTimes.log files. More... | |
class | PaprecaConfig |
Class storing settings and global variables for the PAPRECA run. More... | |
class | PairHash |
Utility hash function (struct) for custom typedefs. More... | |
Typedefs | |
typedef std::vector< Bond > | BOND_VECTOR |
typedef std::unordered_map< LAMMPS_NS::tagint, BOND_VECTOR > | ATOM2BONDS_MAP |
maps atom IDs to their associated PAPRECA::BOND_VECTOR to allow easy access of bonds and bond types. More... | |
typedef std::unordered_map< int, PredefinedReaction * > | TYPE2REACTION_MAP |
typedef std::unordered_map< INT_PAIR, PredefinedBondForm *, PairHash > | PAIR2BONDFORM_MAP |
typedef std::unordered_map< int, PredefinedDiffusionHop * > | TYPE2DIFFUSION_MAP |
typedef std::unordered_map< int, PredefinedDeposition * > | TYPE2DEPOSITION_MAP |
typedef std::unordered_map< int, PredefinedMonoatomicDesorption * > | TYPE2MONODES_MAP |
typedef std::pair< int, int > | INT_PAIR |
typedef std::unordered_set< INT_PAIR, PairHash > | PAIR_SET |
typedef std::unordered_map< INT_PAIR, double, PairHash > | INTPAIR2DOUBLE_MAP |
typedef std::unordered_map< INT_PAIR, int, PairHash > | INTPAIR2INT_MAP |
typedef std::unordered_map< int, int > | INT2INT_MAP |
typedef std::unordered_map< int, INT2INT_MAP > | INT2INTSMAP_MAP |
typedef std::unordered_set< int > | INT_SET |
typedef std::unordered_set< LAMMPS_NS::tagint > | TAGINT_SET |
typedef std::array< double, 3 > | ARRAY3D |
typedef std::pair< double, int > | DOUBLE2INTPAIR |
typedef std::vector< DOUBLE2INTPAIR > | DOUBLE2INTPAIR_VEC |
Functions | |
void | debugPrintBondMapPairs (ATOM2BONDS_MAP const &bonds_map, const int &proc_id) |
void | debugPrintBasicAtomInfo (LAMMPS_NS::LAMMPS *lmp, const int &proc_id) |
void | debugPrintNeighborLists (LAMMPS_NS::LAMMPS *lmp, const int &proc_id) |
void | debugCheckBondsInNeibLists (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, ATOM2BONDS_MAP &atomID2bonds) |
void | debugPrintBondsList (LAMMPS_NS::tagint *bonds_list, LAMMPS_NS::bigint &bonds_num, const int &proc_id) |
void | debugPrintType2SigmaMap (INTPAIR2DOUBLE_MAP &types2sigma) |
void | debugPrintEventInfo (Event *event, const int &proc_id) |
void | debugCheckDeposition () |
void | fillDelidsLocalVec (LAMMPS_NS::LAMMPS *lmp, const double &desorb_cut, std::vector< LAMMPS_NS::tagint > &delids_local, ATOM2BONDS_MAP &atomID2bonds) |
bool | delidsLocalVectorsAreEmpty (std::vector< LAMMPS_NS::tagint > &delids_local) |
void | gatherAndTrimDelIdsOnDriverProc (const int &proc_id, const int &nprocs, std::vector< LAMMPS_NS::tagint > &delids_local, std::vector< LAMMPS_NS::tagint > &delids_global) |
int | fillDelidsVec (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, const double &desorb_cut, std::vector< LAMMPS_NS::tagint > &delids, ATOM2BONDS_MAP &atomID2bonds) |
void | broadcastDelidsFromMasterProc (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, int &delids_num, std::vector< LAMMPS_NS::tagint > &delids) |
void | deleteDesorbedAtoms (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, int &proc_id, const int &nprocs, double &film_height, ATOM2BONDS_MAP &atomID2bonds) |
void | equilibrateFluidAtoms (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, double &time) |
void | equilibrate (LAMMPS_NS::LAMMPS *lmp, int &proc_id, const int &nprocs, double &time, PaprecaConfig &papreca_config, double &film_height, int &zero_rate, const int &KMC_loopid, ATOM2BONDS_MAP &atomID2bonds) |
const bool | feCandidateHas4PO4Neibs (PaprecaConfig &papreca_config, PredefinedDiffusionHop *diff_template, LAMMPS_NS::tagint *atom_ids, int *atom_types, int *neighbors, int &neighbors_num, ATOM2BONDS_MAP &atomID2bonds) |
void | getDiffPointCandidateCoords (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, double *iatom_xyz, double *candidate_xyz, PredefinedDiffusionHop *diff_template) |
const bool | candidateDiffHasCollisions (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, int *neighbors, int &neighbors_num, double *candidate_xyz, const int &diffused_type, double *iatom_xyz, const int &iatom_type) |
void | getDiffEventsFromAtom (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, const int &iatom, int *neighbors, int &neighbors_num, std::vector< Event * > &events_local, ATOM2BONDS_MAP &atomID2bonds) |
const bool | atomIsInDepoScanRange (PaprecaConfig &papreca_config, double *iatom_xyz, double &film_height) |
void | getDepoPointCandidateCoords (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, double *iatom_xyz, double *candidate_xyz, PredefinedDeposition *depo_template) |
const bool | depoCandidateIsBelowRejectionHeight (PaprecaConfig &papreca_config, double *candidate_xyz, const double &film_height) |
void | getMolCoords (LAMMPS_NS::LAMMPS *lmp, double **mol_xyz, double **mol_dx, const int &mol_natoms, double *candidate_center) |
void | initMolCoordsArr (double ***mol_xyz, const int &mol_natoms) |
void | deleteMolCoordsArr (double **mol_xyz, const int &mol_natoms) |
bool | atomHasCollisionWithMolAtoms (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, double *atom_xyz, const int &atom_type, const int &mol_natoms, double **mol_xyz, int *mol_atomtype) |
bool | candidateDepoHasCollisions (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, const int &nprocs, PaprecaConfig &papreca_config, int *neighbors, int neighbors_num, double *candidate_center, double *iatom_xyz, const int &iatom_type, PredefinedDeposition *depo_template) |
void | getDepoEventsFromAtom (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, const int &proc_id, const int &nprocs, const int &iatom, int *neighbors, int &neighbors_num, double &film_height, std::vector< Event * > &events_local) |
const bool | headAtomIsCatalyzed (PredefinedReaction *reaction_template, int *atom_types, int *neighbors, int &neighbors_num) |
const bool | bondLengthIsWithinBreakLimits (LAMMPS_NS::LAMMPS *lmp, PredefinedReaction *break_template, const int &iatom, const LAMMPS_NS::tagint &jatom_id) |
void | getBondBreakingEventsFromAtom (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, const int &iatom, int *neighbors, int &neighbors_num, std::vector< Event * > &events_local, ATOM2BONDS_MAP &atomID2bonds) |
const bool | atomsBelong2TheSameMol (const LAMMPS_NS::tagint &iatom_mol, const LAMMPS_NS::tagint &jneib_mol) |
const bool | atomHasMaxBonds (PaprecaConfig &papreca_config, ATOM2BONDS_MAP &atomID2bonds, const LAMMPS_NS::tagint &atom_id, const int atom_type) |
bool | bondBetweenAtomsExists (ATOM2BONDS_MAP &atomID2bonds, const LAMMPS_NS::tagint &atom1_id, const LAMMPS_NS::tagint &atom2_id) |
const bool | atomCandidatesAreLone (const LAMMPS_NS::tagint atom1_id, const LAMMPS_NS::tagint atom2_id, ATOM2BONDS_MAP &atomID2bonds) |
const bool | atomHasMaxBondTypes (PaprecaConfig &papreca_config, ATOM2BONDS_MAP &atomID2bonds, const LAMMPS_NS::tagint &atom_id, const int &atom_type, const int &bond_type) |
void | getBondFormEventsFromAtom (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, const int &iatom, int *neighbors, int &neighbors_num, std::vector< Event * > &events_local, ATOM2BONDS_MAP &atomID2bonds) |
void | getMonoDesEventsFromAtom (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, const int &iatom, std::vector< Event * > &events_local, ATOM2BONDS_MAP &atomID2bonds) |
void | loopAtomsAndIdentifyEvents (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, int &nprocs, const int &KMC_loopid, PaprecaConfig &papreca_config, std::vector< Event * > &events_local, ATOM2BONDS_MAP &atomID2bonds, double &film_height) |
const bool | bondLengthIsWithinBreakLimits (LAMMPS_NS::LAMMPS *lmp, PredefinedBondForm *bonding_template, const int &iatom, const LAMMPS_NS::tagint &jatom_id) |
void | fillFormTransferDataArr (BondForm *bond_form, int *form_data) |
void | deserializeFormTransferDataArr (int *form_data, int &bond_type, int &delete_atoms) |
void | executeBondForm (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, int &KMC_loopid, double &time, const int &proc_id, const int &nprocs, const int &event_proc, Event *selected_event) |
void | executeBondBreak (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, int &KMC_loopid, double &time, const int &proc_id, const int &nprocs, const int &event_proc, Event *selected_event, ATOM2BONDS_MAP &atomID2bonds) |
void | fillDepoDataTransfArr (double *depo_data, Deposition *depo) |
void | deserializeDepoTransfData (double *depo_data, double *site_pos, double *rot_pos, double &rot_theta, double &insertion_vel) |
void | executeDeposition (LAMMPS_NS::LAMMPS *lmp, int &KMC_loopid, double &time, PaprecaConfig &papreca_config, const int &proc_id, const int &nprocs, const int &event_proc, Event *selected_event) |
void | fillIntegerDiffDataTransfArray (int *diff_intdata, Diffusion *diff) |
void | fillDoubleDiffDataTransfArray (double *diff_doubledata, Diffusion *diff) |
void | deserializeIntegerDiffDataArr (int *diff_intdata, int &parent_type, int &is_displacive, int &diffused_type) |
void | deserializeDoubleDiffDataArr (double *diff_doubledata, double *vac_pos, double &insertion_vel) |
void | executeDiffusion (LAMMPS_NS::LAMMPS *lmp, int &KMC_loopid, double &time, PaprecaConfig &papreca_config, const int &proc_id, const int &nprocs, const int &event_proc, Event *selected_event) |
void | executeMonoatomicDesorption (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, int &KMC_loopid, double &time, const int &proc_id, const int &nprocs, const int &event_proc, Event *selected_event) |
void | printStepInfo (PaprecaConfig &papreca_config, const int &KMC_loopid, const double &time, const double &film_height, const double &proc_rates_sum) |
void | executeEvent (LAMMPS_NS::LAMMPS *lmp, int &KMC_loopid, double &time, PaprecaConfig &papreca_config, const int &proc_id, const int &nprocs, const int &event_proc, const int &event_num, char *event_type, std::vector< Event * > &events_local, ATOM2BONDS_MAP &atomID2bonds) |
int | selectAndExecuteEvent (LAMMPS_NS::LAMMPS *lmp, int &KMC_loopid, double &time, char *event_type, int &proc_id, int &nprocs, PaprecaConfig &papreca_config, std::vector< Event * > &events_local, ATOM2BONDS_MAP &atomID2bonds, double &film_height) |
double | getLocalRate (std::vector< Event * > &events_local, PaprecaConfig &papreca_config) |
void | fillAndSortIndexedRatesVec (double *arr, const int &arr_size, DOUBLE2INTPAIR_VEC &rates_indexed) |
int | selectProcessStochastically (double *arr, const int &arr_size, double &rnum, double &rates_sum) |
void | calcLocalMassAndFillMassProfile (LAMMPS_NS::LAMMPS *lmp, double **mass_profiles, double &local_mass, const int &atom_type, double *atom_xyz, const double &atom_mass, const double &bin_width, const int &bins_num) |
double ** | initMassProfilesArr (const int &types_num, const int &bins_num) |
void | deleteMassProfilesArr (double **mass_profiles, const int &bins_num) |
void | fillMassProfilesTotalArrFromMassProfilesLocal (const int &bins_num, const int &types_num, double **mass_profiles_total, double **mass_profiles_local) |
void | getFilmHeightFromMassBinsMethod (PaprecaConfig &papreca_config, LAMMPS_NS::LAMMPS *lmp, const int &proc_id, double &film_height, double **mass_profiles_total, const double &local_mass, double *atom_mass, const int &bins_num, const int &types_num, const double &bin_width) |
void | calcFilmHeight (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, const int &KMC_loopid, PaprecaConfig &papreca_config, double &film_height) |
const bool | atomsCollide (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config, double *atom1_xyz, const int &atom1_type, double *atom2_xyz, const int &atom2_type) |
void | setTimeUnitsConversionConstant (LAMMPS_NS::LAMMPS *lmp, PaprecaConfig &papreca_config) |
void | checkForAcceptableKeywordsUsedMultipleTimes (std::vector< std::string > &commands, const std::string &keyword) |
void | check4AcceptableKeywords (std::vector< std::string > &commands, const int &start, std::unordered_set< std::string > &acceptable_keywords, const bool &accept_bool) |
void | processCatalyzedOption (std::vector< std::string > &commands, int ¤t_pos, std::vector< int > &catalyzing_types) |
void | processBondLimitOption (std::vector< std::string > &commands, int ¤t_pos, double &length_equil, double &length_perc) |
void | processSigmaMixOptions (std::vector< std::string > &commands, PaprecaConfig &papreca_config, int ¤t_pos) |
void | processBinWidthOptionForElementalDistributions (std::vector< std::string > &commands, PaprecaConfig &papreca_config, int ¤t_pos) |
double | getStickingCoeffFromDepositionEventOptions (std::vector< std::string > &commands, int ¤t_pos) |
double | getRateFromInputRateOptions (std::vector< std::string > &commands, int ¤t_pos) |
void | processCustomDiffEventOptions (std::vector< std::string > &commands, int ¤t_pos, std::string &custom_style, std::vector< int > &style_atomtypes) |
void | executeKMCstepsCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeKMCperMDCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeTimeEndCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeRandomSeedCommand (LAMMPS_NS::LAMMPS *lmp, std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeFluidAtomTypesCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeFrozenAtomTypesCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeDesorptionCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeHeightCalculationCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeSpeciesMaxBondsCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeSpeciesMaxBondTypesCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeMinimizePriorCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeMinimizeAfterCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeTrajectoryDurationCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeDepoheightsCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeRandomDepovecsCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeRandomDiffvecsCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeCreateBondBreakCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeCreateBondFormCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeCreateDiffusionHopCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeCreateDepositionCommand (LAMMPS_NS::LAMMPS *lmp, std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeCreateMonoatomicDesorptionCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeExportHeightVtimeCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeExportSurfaceCoverageCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeExportElementalDistributionsCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeExportExecutionTimesCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeRestartFreqCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeSigmasOptionsCommand (LAMMPS_NS::LAMMPS *lmp, std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executeInitSigmaCommand (LAMMPS_NS::LAMMPS *lmp, std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | executePaprecaCommand (LAMMPS_NS::LAMMPS *lmp, std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
std::vector< std::string > | processLine (char *line) |
void | warn4IllegalRuns (const int &proc_id, PaprecaConfig &papreca_config) |
void | abortIllegalRun (const int &proc_id, PaprecaConfig &papreca_config) |
void | readInputAndInitPaprecaConfig (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, const char *file_name, PaprecaConfig &papreca_config) |
void | processSigmaMixOptions (std::vector< std::string > &commands, int ¤t_pos) |
double | getBinWidthFromElementalDistributions (std::vector< std::string > &commands, int ¤t_pos) |
void | executeCreateDepositionCommand (std::vector< std::string > &commands, PaprecaConfig &papreca_config) |
void | initializeLMP (LAMMPS_NS::LAMMPS **lmp) |
Initialize LAMMPS. More... | |
void | readLMPinput (const std::string &lmp_input, LAMMPS_NS::LAMMPS *lmp) |
void | runLammps (LAMMPS_NS::LAMMPS *lmp, const int ×teps_num) |
void | remap3DArrayInPeriodicBox (LAMMPS_NS::LAMMPS *lmp, double *arr) |
void | resetMobileAtomsGroups (LAMMPS_NS::LAMMPS *lmp, const std::vector< int > &fluid_atomtypes) |
void | deleteAtoms (LAMMPS_NS::LAMMPS *lmp, LAMMPS_NS::tagint *atom_ids, const int &num_atoms, const std::string &delete_bonds, const std::string &delete_molecule) |
void | deleteAtoms (LAMMPS_NS::LAMMPS *lmp, std::vector< LAMMPS_NS::tagint > &atom_ids, const std::string &delete_bonds, const std::string &delete_molecule) |
void | deleteAtomsInBoxRegion (LAMMPS_NS::LAMMPS *lmp, double &boxxlo, double &boxxhi, double &boxylo, double &boxyhi, double &boxzlo, double &boxzhi, const std::string &delete_bonds, const std::string &delete_molecule) |
void | createAtom (LAMMPS_NS::LAMMPS *lmp, const double atom_pos[3], const int &atom_type) |
void | deleteBond (LAMMPS_NS::LAMMPS *lmp, const LAMMPS_NS::tagint &atom1id, const LAMMPS_NS::tagint &atom2id, const bool special) |
void | formBond (LAMMPS_NS::LAMMPS *lmp, const LAMMPS_NS::tagint &atom1id, const LAMMPS_NS::tagint &atom2id, const int &bond_type) |
void | insertMolecule (LAMMPS_NS::LAMMPS *lmp, const double site_pos[3], const double rot_pos[3], const double &rot_theta, const int &type_offset, const char *mol_name) |
void | diffuseAtom (LAMMPS_NS::LAMMPS *lmp, const double vac_pos[3], const LAMMPS_NS::tagint &parent_id, const int &parent_type, const int &is_displacive, const int &diffused_type) |
void | initType2SigmaFromLammpsPairCoeffs (LAMMPS_NS::LAMMPS *lmp, INTPAIR2DOUBLE_MAP &type2sigma) |
int | getMaskedNeibIndex (int *neighbors, int &j) |
void | initAndGatherBondsList (LAMMPS_NS::LAMMPS *lmp, LAMMPS_NS::tagint **bonds_list, LAMMPS_NS::bigint &bonds_num) |
const int | getMolIndexFromMolName (LAMMPS_NS::LAMMPS *lmp, std::string mol_name) |
void | computeMolCenter (LAMMPS_NS::LAMMPS *lmp, std::string mol_name) |
void | dumpRestart (LAMMPS_NS::LAMMPS *lmp, const int &KMC_loopid, const int &dump_freq) |
double | get3DSqrDistWithPBC (LAMMPS_NS::LAMMPS *lmp, const double *x1, const double *x2) |
void | MPIBcastAndExecuteCommand (LAMMPS_NS::LAMMPS *lmp, std::string &command) |
void | setupMPI (int *narg, char ***arg, int *nprocs, int *proc_id) |
const int | getMPIRank (MPI_Comm communicator) |
void | warnOne (MPI_Comm communicator, const std::string &message) |
void | warnAll (MPI_Comm communicator, const std::string &message) |
void | allAbort (MPI_Comm communicator) |
void | allAbortWithMessage (MPI_Comm communicator, const std::string &message) |
double | getRateFromArrhenius (const double &activation_energy, const double &attempt_freq, const double &temperature) |
double | getDepoRateFromHertzKnudsen (const double &pressure_in, const double &ads_area_in, const double &ads_mass_in, const double &temperature_in) |
double | getDesorptionRate (const double &activation_energy, const double &temperature) |
void | advanceSimClockFromKMC (PaprecaConfig &papreca_config, const double &proc_rates_sum, double &time) |
void | advanceSimClockFromLAMMPS (PaprecaConfig &papreca_config, double &time) |
double | doubleArrSum (double *arr, const int &size) |
void | copyDoubleArray3D (double *copy, const double *source, const int start, const int end) |
const std::string | getConcatenatedStringFromStringsVector (const std::vector< std::string > &strings, size_t start=-1, size_t end=-1) |
const std::string | getConcatenatedStringWithSpacesFromStringsVector (const std::vector< std::string > &strings, size_t start=-1, size_t end=-1) |
const int | getStringPosInStringVec (const std::string &string2check, const std::vector< std::string > &strings) |
bool | stringIsNumber (const std::string &string) |
bool | stringIsIntNumber (const std::string &string) |
const bool | stringIsBool (const std::string &string) |
const unsigned long int | string2UnsignedLongInt (std::string &string) |
const int | string2Int (std::string &string) |
const double | string2Double (std::string &string) |
const bool | string2Bool (std::string &string) |
const int | boolString2Int (std::string &string) |
template<typename T > | |
T | getSumOfVecElements (const std::vector< T > &vec) |
template<typename T , size_t N> | |
void | pushArrayToVector (const T *arr, std::vector< T > &vec) |
template<typename T > | |
std::vector< T > | getSubVectorFromVector (const std::vector< T > &vec, size_t start, size_t end) |
template<typename T , typename Hash > | |
bool | elementIsInUnorderedSet (const std::unordered_set< T, Hash > &set, const T &element) |
template<typename KeyType , typename ValueType , typename Hash > | |
bool | mappingExists (const std::unordered_map< KeyType, ValueType, Hash > &map, const KeyType &key) |
template<typename T > | |
bool | elementIsInVector (const std::vector< T > &vec, const T &element) |
typedef std::array< double , 3 > PAPRECA::ARRAY3D |
typedef std::unordered_map< LAMMPS_NS::tagint , BOND_VECTOR > PAPRECA::ATOM2BONDS_MAP |
maps atom IDs to their associated PAPRECA::BOND_VECTOR to allow easy access of bonds and bond types.
typedef std::vector< Bond > PAPRECA::BOND_VECTOR |
typedef std::pair< double , int > PAPRECA::DOUBLE2INTPAIR |
typedef std::vector< DOUBLE2INTPAIR > PAPRECA::DOUBLE2INTPAIR_VEC |
typedef std::unordered_map< int , int > PAPRECA::INT2INT_MAP |
typedef std::unordered_map< int , INT2INT_MAP > PAPRECA::INT2INTSMAP_MAP |
typedef std::pair< int , int > PAPRECA::INT_PAIR |
typedef std::unordered_set< int > PAPRECA::INT_SET |
typedef std::unordered_map< INT_PAIR , double , PairHash > PAPRECA::INTPAIR2DOUBLE_MAP |
typedef std::unordered_map< INT_PAIR , int , PairHash > PAPRECA::INTPAIR2INT_MAP |
typedef std::unordered_map< INT_PAIR , PredefinedBondForm* , PairHash > PAPRECA::PAIR2BONDFORM_MAP |
typedef std::unordered_set< INT_PAIR , PairHash > PAPRECA::PAIR_SET |
typedef std::unordered_set< LAMMPS_NS::tagint > PAPRECA::TAGINT_SET |
typedef std::unordered_map< int , PredefinedDeposition* > PAPRECA::TYPE2DEPOSITION_MAP |
typedef std::unordered_map< int , PredefinedDiffusionHop* > PAPRECA::TYPE2DIFFUSION_MAP |
typedef std::unordered_map< int , PredefinedMonoatomicDesorption* > PAPRECA::TYPE2MONODES_MAP |
typedef std::unordered_map< int , PredefinedReaction* > PAPRECA::TYPE2REACTION_MAP |
void PAPRECA::abortIllegalRun | ( | const int & | proc_id, |
PaprecaConfig & | papreca_config | ||
) |
After processing all commands, this function is called abort illegal PAPRECA runs.
[in] | proc_id | ID of current MPI process. |
[in] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::advanceSimClockFromKMC | ( | PaprecaConfig & | papreca_config, |
const double & | proc_rates_sum, | ||
double & | time | ||
) |
Stochastically advances the simulation clock by a KMC step.
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | proc_rates_sum | total (gathered/sum from all MPI processes) event rate. |
[in,out] | time | time on current PAPRECA step. |
void PAPRECA::advanceSimClockFromLAMMPS | ( | PaprecaConfig & | papreca_config, |
double & | time | ||
) |
Advances the simulation clock by timestep*dt (as defined in the LAMMPS and PAPRECA inputs).
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in,out] | time | current time. |
void PAPRECA::allAbort | ( | MPI_Comm | communicator | ) |
Aborts all MPI processes.
[in] | communicator | communicator |
void PAPRECA::allAbortWithMessage | ( | MPI_Comm | communicator, |
const std::string & | message | ||
) |
Aborts all MPI processes and throws a message on the screen.
[in] | communicator | MPI communicator. |
[in] | message | std::string of the error message. |
const bool PAPRECA::atomCandidatesAreLone | ( | const LAMMPS_NS::tagint | atom1_id, |
const LAMMPS_NS::tagint | atom2_id, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Checks if two atoms are lone (i.e., if they do not have any bonds). This is done by retrieving the PAPRECA:Bond objects vector for each atom from the atomID2bonds container.
[in] | atom1_id | ID of the first atom. |
[in] | atom2_id | ID of the second atom. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
bool PAPRECA::atomHasCollisionWithMolAtoms | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
double * | atom_xyz, | ||
const int & | atom_type, | ||
const int & | mol_natoms, | ||
double ** | mol_xyz, | ||
int * | mol_atomtype | ||
) |
Checks if the parent event atom has collisions with any of the inserted molecule atoms.
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | atom_xyz | coordinates of the parent atom. |
[in] | atom_type | atom type of parent atom. |
[in] | mol_natoms | total number of atoms in the inserted molecule. |
[in] | mol_xyz | a mol_natomsx3 array storing the coordinates of all inserted molecule atoms. |
[in] | mol_atomtype | array storing the atom types of all inserted molecule atoms. |
const bool PAPRECA::atomHasMaxBonds | ( | PaprecaConfig & | papreca_config, |
ATOM2BONDS_MAP & | atomID2bonds, | ||
const LAMMPS_NS::tagint & | atom_id, | ||
const int | atom_type | ||
) |
Determines the total number of bonds for an atom by retrieving its PAPRECA::Bond objects vector. Then compares it to the user defined species max bonds.
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
[in] | atom_id | ID of atom. |
[in] | atom_type | type of current atom. |
const bool PAPRECA::atomHasMaxBondTypes | ( | PaprecaConfig & | papreca_config, |
ATOM2BONDS_MAP & | atomID2bonds, | ||
const LAMMPS_NS::tagint & | atom_id, | ||
const int & | atom_type, | ||
const int & | bond_type | ||
) |
Checks if the number of bonds of a specific bond type for the current atom is greater or equal than the user defined max bond types of species. This is done by retrieving the PAPRECA::Bond objects vector of the current atom from the atomID2bonds map.
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
[in] | atom_id | id of current atom. |
[in] | atom_type | of current atom. |
[in] | bond_type | to be checked |
const bool PAPRECA::atomIsInDepoScanRange | ( | PaprecaConfig & | papreca_config, |
double * | iatom_xyz, | ||
double & | film_height | ||
) |
Checks if the parent atom is within the acceptable scan range for deposition events. The acceptable scan ranges can be modified through the PAPRECA input file.
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | iatom_xyz | coordinates of candidate atom |
[in] | film_height | height in current PAPRECA step. |
const bool PAPRECA::atomsBelong2TheSameMol | ( | const LAMMPS_NS::tagint & | iatom_mol, |
const LAMMPS_NS::tagint & | jneib_mol | ||
) |
Uses LAMMPS molecule IDs to decide if two atoms pbelong to the same mol.
[in] | iatom_mol | ID of the first atom. |
[in] | jneib_mol | ID of the second atom (which is typically one of the neighbors of iatom). |
const bool PAPRECA::atomsCollide | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
double * | atom1_xyz, | ||
const int & | atom1_type, | ||
double * | atom2_xyz, | ||
const int & | atom2_type | ||
) |
Checks for collisions between two atoms. It is assumed that two atoms of types i and j "collide" if the distance between them is smaller than sigma_ij.
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | atom1_xyz | coordinates of the first atom (x,y, and z). |
[in] | atom1_type | atom type of the first atom. |
[in] | atom2_xyz | coordinates of the second atom (x,y, and z). |
[in] | atom2_type | atom type of the second atom. |
bool PAPRECA::bondBetweenAtomsExists | ( | ATOM2BONDS_MAP & | atomID2bonds, |
const LAMMPS_NS::tagint & | atom1_id, | ||
const LAMMPS_NS::tagint & | atom2_id | ||
) |
The PAPRECA::Bond objects vector is retrieved from the atomID2bonds container for two atoms to determine if the two atoms are already bonded.
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
[in] | atom1_id | ID of the first atom. |
[in] | atom2_id | ID of the second atom. |
const bool PAPRECA::bondLengthIsWithinBreakLimits | ( | LAMMPS_NS::LAMMPS * | lmp, |
PredefinedBondForm * | bonding_template, | ||
const int & | iatom, | ||
const LAMMPS_NS::tagint & | jatom_id | ||
) |
const bool PAPRECA::bondLengthIsWithinBreakLimits | ( | LAMMPS_NS::LAMMPS * | lmp, |
PredefinedReaction * | break_template, | ||
const int & | iatom, | ||
const LAMMPS_NS::tagint & | jatom_id | ||
) |
Checks if the current length of the bond is within the user-defined limits for bond-breaking.
[in] | lmp | pointer to LAMMPS object. |
[in] | break_template | pointer to predefined bond-breaking object |
[in] | iatom | local index of atom i |
[in] | jatom_id | id of bonded atom to atom i |
const int PAPRECA::boolString2Int | ( | std::string & | string | ) |
This function receives an std::string which is yes/no and converts it to an int (1/0, respectively).
[in] | string | std::string to be converted. |
void PAPRECA::broadcastDelidsFromMasterProc | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id, | ||
int & | delids_num, | ||
std::vector< LAMMPS_NS::tagint > & | delids | ||
) |
Broadcasts the delids vector (containing the IDs of atoms marked for deletion) from the master MPI process to all other processes. This function is only called when the delete_desorbed option is active and the gather_all algorithm is selected.
[in] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | delids_num | number of atom IDs marked for deletion. |
[in,out] | delids | vector of collected atom IDs initialized by the master MPI process and broadcasted to all other MPI processes. |
void PAPRECA::calcFilmHeight | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id, | ||
const int & | KMC_loopid, | ||
PaprecaConfig & | papreca_config, | ||
double & | film_height | ||
) |
Calculates the film height based on the user defined film calculation method.
[in] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | KMC_loopid | PAPRECA (kMC) step number. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in,out] | film_height | height at current PAPRECA (kMC) step. |
void PAPRECA::calcLocalMassAndFillMassProfile | ( | LAMMPS_NS::LAMMPS * | lmp, |
double ** | mass_profiles, | ||
double & | local_mass, | ||
const int & | atom_type, | ||
double * | atom_xyz, | ||
const double & | atom_mass, | ||
const double & | bin_width, | ||
const int & | bins_num | ||
) |
Appends atom_mass to the relevant bin of mass_profiles array.
[in] | lmp | pointer to LAMMPS object. |
[in,out] | mass_profiles | 2D (i.e., mass_profiles[bins_num][atomtypes_num]) array containing the total mass of a specific atom type for each x-y slice (bin) in the system. |
[in,out] | local_mass | total mass in current MPI process. |
[in] | atom_type | type of the current atom. |
[in] | atom_xyz | array containing the coordinates of the current atom. |
[in] | atom_mass | mass of current atom. |
[in] | bin_width | width of current x-y slice (bin). |
[in] | bins_num | total number of x-y slices (bins). |
bool PAPRECA::candidateDepoHasCollisions | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id, | ||
const int & | nprocs, | ||
PaprecaConfig & | papreca_config, | ||
int * | neighbors, | ||
int | neighbors_num, | ||
double * | candidate_center, | ||
double * | iatom_xyz, | ||
const int & | iatom_type, | ||
PredefinedDeposition * | depo_template | ||
) |
Checks if the inserted molecule atoms have collisions with 1) the parent atom, 2) all neighbors of the parent event atom.
[in] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | neighbors | IDs of all neighbors of the parent event atom. |
[in] | neighbors_num | number of atoms. |
[in] | candidate_center | coordinates of the center of mass of the candidate inserted molecule. |
[in] | iatom_xyz | coordinates of parent atom. |
[in] | iatom_type | atom type of parent atom. |
[in] | depo_template | Deposition template (PAPRECA::PredefinedDeposition) as initialized by the user (in the PAPRECA input file). |
const bool PAPRECA::candidateDiffHasCollisions | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
int * | neighbors, | ||
int & | neighbors_num, | ||
double * | candidate_xyz, | ||
const int & | diffused_type, | ||
double * | iatom_xyz, | ||
const int & | iatom_type | ||
) |
Checks for collisions between the diffusion atom (i.e., atom moving to the vacancy) and existing atoms in the simulation.
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | neighbors | IDs of neighbors of the parent atom. |
[in] | neighbors_num | number of neighbors of the parent atom. |
[in] | candidate_xyz | coordinates (x,y, and z) of vacancy (point towards which the atom diffuses). |
[in] | diffused_type | atom type of diffused atom. |
[in] | iatom_xyz | coordinates of the parent atom. |
[in] | iatom_type | atom type of parent atom. |
void PAPRECA::check4AcceptableKeywords | ( | std::vector< std::string > & | commands, |
const int & | start, | ||
std::unordered_set< std::string > & | acceptable_keywords, | ||
const bool & | accept_bool | ||
) |
Checks if a specific command line includes acceptable keyword. The "acceptable_keyword" unordered set has to be initialized correctly by the user, before the call to this function. The present function aborts if it detected a non-acceptable keyword (i.e, a keyword not included in the acceptable_keywords set).
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in] | start | current index position in the commands vector. Only check keywords for indices that are equal or greater than start. |
[in] | acceptable_keywords | unordered_set containing the acceptable (for this command) keywords. Has to be initialized by the user before the call to this function. |
[in] | accept_bool | if true, then yes and no keywords are also accepted in the command line. |
void PAPRECA::checkForAcceptableKeywordsUsedMultipleTimes | ( | std::vector< std::string > & | commands, |
const std::string & | keyword | ||
) |
Receives the commands vector of strings and aborts if the provided "keyword" is used multiple times.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in] | keyword | specific keyword to be checked. |
void PAPRECA::computeMolCenter | ( | LAMMPS_NS::LAMMPS * | lmp, |
std::string | mol_name | ||
) |
Invokes the LAMMPS function compute_center for the mol_name.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | mol_name | valid molecule names (as defined in the LAMMPS input file. e.g., if you used this command: "molecule mmmTCP ./TCP.xyz", then mol_name should be "mmmTCP"). |
void PAPRECA::copyDoubleArray3D | ( | double * | copy, |
const double * | source, | ||
const int | start, | ||
const int | end | ||
) |
Receives an array of doubles and copies it to a different array of doubles.
[in,out] | copy | copied array. |
[in] | source | original array. |
[in] | start | copy starting from that index. |
[in] | end | copy all elements whose index is less than this variable (i.e., copy[i] with i < end ). |
void PAPRECA::createAtom | ( | LAMMPS_NS::LAMMPS * | lmp, |
const double | atom_pos[3], | ||
const int & | atom_type | ||
) |
Inserts an atom in the simulation.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | atom_pos | coordinates (x,y, and z) of insertion point. |
[in] | atom_type | type of insert atom |
void PAPRECA::debugCheckBondsInNeibLists | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
void PAPRECA::debugCheckDeposition | ( | ) |
void PAPRECA::debugPrintBasicAtomInfo | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id | ||
) |
void PAPRECA::debugPrintBondMapPairs | ( | ATOM2BONDS_MAP const & | bonds_map, |
const int & | proc_id | ||
) |
void PAPRECA::debugPrintBondsList | ( | LAMMPS_NS::tagint * | bonds_list, |
LAMMPS_NS::bigint & | bonds_num, | ||
const int & | proc_id | ||
) |
void PAPRECA::debugPrintEventInfo | ( | Event * | event, |
const int & | proc_id | ||
) |
void PAPRECA::debugPrintNeighborLists | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id | ||
) |
void PAPRECA::debugPrintType2SigmaMap | ( | INTPAIR2DOUBLE_MAP & | types2sigma | ) |
void PAPRECA::deleteAtoms | ( | LAMMPS_NS::LAMMPS * | lmp, |
LAMMPS_NS::tagint * | atom_ids, | ||
const int & | num_atoms, | ||
const std::string & | delete_bonds, | ||
const std::string & | delete_molecule | ||
) |
Receives an array of atom IDs (atom_ids) and deletes all atoms in the simulation with the corresponding IDs.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | atom_ids | array of IDs of atoms to be deleted. |
[in] | num_atoms | number of atoms to be deleted (length of atom_ids vector). |
[in] | delete_bonds | if "yes", then the LAMMPS delete_atoms command is called with the "bond yes" keyword. |
[in] | delete_molecule | if "yes", then the LAMMPS delete atoms command is called with the "mol yes" option. |
void PAPRECA::deleteAtoms | ( | LAMMPS_NS::LAMMPS * | lmp, |
std::vector< LAMMPS_NS::tagint > & | atom_ids, | ||
const std::string & | delete_bonds, | ||
const std::string & | delete_molecule | ||
) |
Receives a vector of atom IDs and deletes all atoms in the simulation with the corresponding IDs.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | atom_ids | std::vector of IDs of atoms to be deleted. |
[in] | delete_bonds | if "yes", then the LAMMPS delete_atoms command is called with the "bond yes" keyword. |
[in] | delete_molecule | if "yes", then the LAMMPS delete atoms command is called with the "mol yes" option. |
void PAPRECA::deleteAtomsInBoxRegion | ( | LAMMPS_NS::LAMMPS * | lmp, |
double & | boxxlo, | ||
double & | boxxhi, | ||
double & | boxylo, | ||
double & | boxyhi, | ||
double & | boxzlo, | ||
double & | boxzhi, | ||
const std::string & | delete_bonds, | ||
const std::string & | delete_molecule | ||
) |
Deletes all atoms in a block region defined by the input region bounds.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | boxxlo | low deletion region x-boundary. |
[in] | boxxhi | high deletion region x-boundary. |
[in] | boxylo | low deletion region y-boundary. |
[in] | boxyhi | high deletion region y-boundary. |
[in] | boxzlo | low deletion region z-boundary. |
[in] | boxzhi | high deletion region z-boundary. |
[in] | delete_bonds | if "yes" then the delete_atoms LAMMPS command is called with the "bond yes" keyword. The delete_atoms command is called to delete atoms in the deletion region. |
[in] | delete_molecule | if "yes" then the delete_atoms LAMMPS command is called with the "mol yes" keyword. The delete_atoms command is called to delete atoms in the deletion region. |
void PAPRECA::deleteBond | ( | LAMMPS_NS::LAMMPS * | lmp, |
const LAMMPS_NS::tagint & | atom1id, | ||
const LAMMPS_NS::tagint & | atom2id, | ||
const bool | special | ||
) |
Deletes an existing bond between atom1 and atom2.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | atom1id | ID of the first atom. |
[in] | atom2id | ID of the second atom. |
[in] | special | if true, then the delete_bonds LAMMPS command is called with the "special" keyword. |
void PAPRECA::deleteDesorbedAtoms | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
int & | proc_id, | ||
const int & | nprocs, | ||
double & | film_height, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
This function is always called but performs computations only if the user has set a desorption height cutoff in the PAPRECA input file. The present function compares the z-coordinate of each atom with the desorption height cutoff. If any z-coordinate value is greater than or equal to the desorption height cutoff, the associated atom (along with its bonded atoms) is deleted. Currently, the user can select between two different algorithms: 1) gather_local (see fillDelidsLocalVec() function description/notes), and 2) gather_all (see fillDelidsVec() function description/notes). A comparison of the performance between the two algorithms is not currently available. However, as a quick note, it can be mentioned that gather_all is expected to be more memory intensive, since it calls the LAMMPS function lammps_gather_atoms_concat() to gather the coordinates of non-consecutive IDs on the master proc.
[in,out] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | film_height | height at current PAPRECA step. |
[in,out] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
void PAPRECA::deleteMassProfilesArr | ( | double ** | mass_profiles, |
const int & | bins_num | ||
) |
Deletes a previously-initialized 2D array of mass profiles.
[in,out] | mass_profiles | 2D array of mass profiles |
[in] | bins_num | number of x-y slices (bins). |
void PAPRECA::deleteMolCoordsArr | ( | double ** | mol_xyz, |
const int & | mol_natoms | ||
) |
Deletes the molecule coordinates array used to check for collisions between the deposited molecule and system atoms.
[in,out] | mol_xyz | array storing the x,y, and z distances of each molecule atom from the molecule center. See molecule.h and molecule.cpp files in the LAMMPS source directory for more information. |
[in] | mol_natoms | total number of molecule atoms. |
bool PAPRECA::delidsLocalVectorsAreEmpty | ( | std::vector< LAMMPS_NS::tagint > & | delids_local | ) |
Communicates information between MPI processes to check if every local (i.e., on every MPI process) delids_local container is empty. Effectively, this means that the z-coordinates of all atoms in the simulation are lower than desorb_cut and that no deletions have to be performed.
[in] | delids_local | vector of collected atom IDs. |
const bool PAPRECA::depoCandidateIsBelowRejectionHeight | ( | PaprecaConfig & | papreca_config, |
double * | candidate_xyz, | ||
const double & | film_height | ||
) |
Checks if the deposition candidate point is below the rejection height (as set in the PAPRECA input).
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | candidate_xyz | array storing the coordinates of the candidate deposition point. |
[in] | film_height | film height at the current PAPRECA step. |
void PAPRECA::deserializeDepoTransfData | ( | double * | depo_data, |
double * | site_pos, | ||
double * | rot_pos, | ||
double & | rot_theta, | ||
double & | insertion_vel | ||
) |
Deserializes (post-processes) data after calls to MPI function related to the execution of PAPRECA::Deposition events.
[in] | depo_data | 8-element array of serialized (in fillDepoDataTransfArr() function) data. |
[in,out] | site_pos | array containing the coordinates of the center-of-mass of the inserted molecule/atom. |
[in,out] | rot_pos | array containing the coordinates of the center-of-rotation of the inserted molecule. |
[in,out] | rot_theta | angle of rotation of inserted molecule |
[in,out] | insertion_vel | velocity of inserted molecule. |
void PAPRECA::deserializeDoubleDiffDataArr | ( | double * | diff_doubledata, |
double * | vac_pos, | ||
double & | insertion_vel | ||
) |
Deserializes (post-processes) double data after calls to MPI function related to the execution of PAPRECA::Diffusion events.
[in] | diff_doubledata | 4-element array of serialized (in fillDoubleDiffDataTransfArray() function) data. |
[in,out] | vac_pos | coordinates of vacancy. |
[in,out] | insertion_vel | velocity of diffused atom. |
void PAPRECA::deserializeFormTransferDataArr | ( | int * | form_data, |
int & | bond_type, | ||
int & | delete_atoms | ||
) |
Deserializes (post-processes) data after calls to MPI function related to the execution of PAPRECA::BondForm events.
[in] | form_data | 2-element array of serialized (in fillFormTransferDataArr() function) data. |
[in,out] | bond_type | type of bond. |
[in,out] | delete_atoms | 0 if atoms are not deleted after the execution of the bond formation event, or 1 if atoms are deleted. |
void PAPRECA::deserializeIntegerDiffDataArr | ( | int * | diff_intdata, |
int & | parent_type, | ||
int & | is_displacive, | ||
int & | diffused_type | ||
) |
Deserializes (post-processes) integer data after calls to MPI function related to the execution of PAPRECA::Diffusion events.
[in] | diff_intdata | 3-element array of serialized (in fillIntegerDiffDataTransfArray() function) data. |
[in,out] | parent_type | atom type of parent atom. |
[in,out] | is_displacive | 0 or 1 depending on whether the diffusion event is diplacive or not. |
[in,out] | diffused_type | atom type of diffused atom. |
void PAPRECA::diffuseAtom | ( | LAMMPS_NS::LAMMPS * | lmp, |
const double | vac_pos[3], | ||
const LAMMPS_NS::tagint & | parent_id, | ||
const int & | parent_type, | ||
const int & | is_displacive, | ||
const int & | diffused_type | ||
) |
Executes a diffusion operation (i.e., moves if the diffusion is displacive, or creates an atom if the diffusion is non-displacive).
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | vac_pos | coordinates (x,y, and z) of vacancy. |
[in] | parent_id | ID of parent candidate atom. |
[in] | parent_type | atom type of parent candidate atom. |
[in] | is_displacive | 1 if the diffusion is displacive (i.e., if the parent atom moves), or 0 if the diffusion is non-displacive (i.e., if a new atom is created on the vacancy position). |
[in] | diffused_type | type of diffused atom. Can be the same as parent type or can be set to a different type if you wish to change the atom type after performing a diffusion event. |
double PAPRECA::doubleArrSum | ( | double * | arr, |
const int & | size | ||
) |
Receives a C-style array of doubles and returns the sum of its elements.
[in] | arr | pointer to 1-D array of doubles. |
[in] | size | size (length) of arr. |
void PAPRECA::dumpRestart | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | KMC_loopid, | ||
const int & | dump_freq | ||
) |
Dumps a LAMMPS simulation restart file (can be used to restart a PAPRECA simulation).
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | KMC_loopid | current PAPRECA step number. |
[in] | dump_freq | dump a restart file every dump_freq PAPRECA step. |
bool PAPRECA::elementIsInUnorderedSet | ( | const std::unordered_set< T, Hash > & | set, |
const T & | element | ||
) |
Searches if an element is included in an unordered set.
[in] | set | unordered set to be searched. |
[in] | element | element to be searched in the unordered set. |
bool PAPRECA::elementIsInVector | ( | const std::vector< T > & | vec, |
const T & | element | ||
) |
Searches if an element is included in an std::vector.
[in] | vec | vector to be searched. |
[in] | element | element to be searched in the vector. |
void PAPRECA::equilibrate | ( | LAMMPS_NS::LAMMPS * | lmp, |
int & | proc_id, | ||
const int & | nprocs, | ||
double & | time, | ||
PaprecaConfig & | papreca_config, | ||
double & | film_height, | ||
int & | zero_rate, | ||
const int & | KMC_loopid, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
This function performs a LAMMPS run on the current system configuration, every KMC_per_MD (as set by the user in the PAPRECA input file). Then, it deletes atoms whose z-coordinate is equal to or greater than the desorption height cutoff (defined in the PAPRECA input file).
[in,out] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in,out] | time | current PAPRECA simulation time. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | film_height | current PAPRECA simulation height. |
[in] | zero_rate | 0 if the total event rate at the current step is zero, or 1 otherwise. |
[in] | KMC_loopid | current PAPRECA simulation step. |
[in,out] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
void PAPRECA::equilibrateFluidAtoms | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
double & | time | ||
) |
Performs a LAMMPS simulation on the fluid atom types (as defined in the PAPRECA input). Then, updates the simulation clock by timestep*trajectory_duration (as defined by the user in the LAMMPS and PAPRECA inputs). Additionally, might perform minimizations before/after the LAMMPS trajectory (if an appropriate LAMMPS minimization command is defined by the user).
[in,out] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in,out] | time | current time. |
void PAPRECA::executeBondBreak | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
int & | KMC_loopid, | ||
double & | time, | ||
const int & | proc_id, | ||
const int & | nprocs, | ||
const int & | event_proc, | ||
Event * | selected_event, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Executes PAPRECA::Bondbreak event. This is done by 1) communicating information from the MPI process that detected the event to all other MPI processes, and 2) calling deleteBond() after the appropriate data have been communicated.
[in,out] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | event_proc | MPI process calling the function (i.e., proc that detected this event). |
[in] | selected_event | selected (for execution) PAPRECA::Event. Casted to the correct type (i.e., PAPRECA::BondBreak) before execution. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
void PAPRECA::executeBondForm | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
int & | KMC_loopid, | ||
double & | time, | ||
const int & | proc_id, | ||
const int & | nprocs, | ||
const int & | event_proc, | ||
Event * | selected_event | ||
) |
Executes PAPRECA::BondForm event. This is done by 1) communicating information from the MPI process that detected the event to all other MPI processes, and 2) calling formBond() and potentially deleteAtoms() after the appropriate data have been communicated.
[in,out] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | event_proc | MPI process calling the function (i.e., proc that detected this event). |
[in] | selected_event | PAPRECA::Event selected (to be executed). Casted to the correct type (i.e., PAPRECA::BondForm) before execution. |
void PAPRECA::executeCreateBondBreakCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Initializes a PAPRECA::PredefinedBondBreak event template and inserts it in the PAPRECA::PredefinedEventsCatalog of the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeCreateBondFormCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Initializes a PAPRECA::PredefinedBondForm event template and inserts it in the PAPRECA::PredefinedEventsCatalog of the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeCreateDepositionCommand | ( | LAMMPS_NS::LAMMPS * | lmp, |
std::vector< std::string > & | commands, | ||
PaprecaConfig & | papreca_config | ||
) |
Initializes a PAPRECA::PredefinedDeposition event template and inserts it in the PAPRECA::PredefinedEventsCatalog of the PAPRECA::PaprecaConfig object.
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeCreateDepositionCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
void PAPRECA::executeCreateDiffusionHopCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Initializes a PAPRECA::PredefinedDiffusionHop event template and inserts it in the PAPRECA::PredefinedEventsCatalog of the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeCreateMonoatomicDesorptionCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Initializes a PAPRECA::PredefinedMonoatomicDesorption event template and inserts it in the PAPRECA::PredefinedEventsCatalog of the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeDepoheightsCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the c_height_scan and c_height_reject parameters in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeDeposition | ( | LAMMPS_NS::LAMMPS * | lmp, |
int & | KMC_loopid, | ||
double & | time, | ||
PaprecaConfig & | papreca_config, | ||
const int & | proc_id, | ||
const int & | nprocs, | ||
const int & | event_proc, | ||
Event * | selected_event | ||
) |
Executes PAPRECA::Deposition event. This is done by 1) communicating information from the MPI process that detected the event to all other MPI processes, and 2) calling insertMolecule() after the appropriate data have been communicated.
[in,out] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | event_proc | MPI process calling the function (i.e., proc that detected this event). |
[in] | selected_event | PAPRECA::Event selected (to be executed). Casted to the correct type (i.e., PAPRECA::Deposition) before execution. |
void PAPRECA::executeDesorptionCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the necessary parameters for the desorption command in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeDiffusion | ( | LAMMPS_NS::LAMMPS * | lmp, |
int & | KMC_loopid, | ||
double & | time, | ||
PaprecaConfig & | papreca_config, | ||
const int & | proc_id, | ||
const int & | nprocs, | ||
const int & | event_proc, | ||
Event * | selected_event | ||
) |
Executes PAPRECA::Diffusion event. This is done by 1) communicating information from the MPI process that detected the event to all other MPI processes, and 2) calling diffuseAtom() after the appropriate data have been communicated.
[in,out] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | event_proc | MPI process calling the function (i.e., proc that detected this event). |
[in] | selected_event | PAPRECA::Event selected (to be executed). Casted to the correct type (i.e., PAPRECA::Diffusion) before execution. |
void PAPRECA::executeEvent | ( | LAMMPS_NS::LAMMPS * | lmp, |
int & | KMC_loopid, | ||
double & | time, | ||
PaprecaConfig & | papreca_config, | ||
const int & | proc_id, | ||
const int & | nprocs, | ||
const int & | event_proc, | ||
const int & | event_num, | ||
char * | event_type, | ||
std::vector< Event * > & | events_local, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Casts event type from parent MPI process to all other processes and then executes event. Currently, only PAPRECA::BondForm, PAPRECA::BondBreak, PAPRECA::Deposition, PAPRECA::Diffusion, PAPRECA::MonoatomicDesorption events are supported.
[in,out] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | event_proc | ID of proc that discovered the event. |
[in] | event_num | local (on the MPI process that discovered the event) event index (in the PAPRECA::Event objects vector). |
[in,out] | event_type | type of selected event. |
[in] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
void PAPRECA::executeExportElementalDistributionsCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Enables the generation of PAPRECA::ElementalDistribution files with a give frequency in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeExportExecutionTimesCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Enables the generation of PAPRECA::ExecTime files with a give frequency in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeExportHeightVtimeCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Enables the generation of PAPRECA::HeightVtime files with a give frequency in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeExportSurfaceCoverageCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Enables the generation of PAPRECA::SurfaceCoverage files with a give frequency in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeFluidAtomTypesCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Initializes the fluid atom types in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeFrozenAtomTypesCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Initializes the frozen atom types in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeHeightCalculationCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets all parameters relevant to height calculation for the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeInitSigmaCommand | ( | LAMMPS_NS::LAMMPS * | lmp, |
std::vector< std::string > & | commands, | ||
PaprecaConfig & | papreca_config | ||
) |
Sets the sigma value for a pair of atom types in the PAPRECA::PaprecaConfig object.
[in] | lmp | pointer to previously initialized LAMMPS object. |
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeKMCperMDCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the frequency of MD (LAMMPS) steps (i.e., how many KMC stages are performed for every MD stages) in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeKMCstepsCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the total number of KMC steps in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeMinimizeAfterCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the after minimization command in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeMinimizePriorCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets a prior minimization command in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeMonoatomicDesorption | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
int & | KMC_loopid, | ||
double & | time, | ||
const int & | proc_id, | ||
const int & | nprocs, | ||
const int & | event_proc, | ||
Event * | selected_event | ||
) |
Executes PAPRECA::MonoatomicDesorption event. This is done by 1) communicating information from the MPI process that detected the event to all other MPI processes, and 2) calling deleteAtoms() after the appropriate data have been communicated.
[in,out] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | event_proc | MPI process calling the function (i.e., proc that detected this event). |
[in] | selected_event | PAPRECA::Event selected (to be executed). Casted to the correct type (i.e., PAPRECA::MonoatomicDesorption) before execution. |
void PAPRECA::executePaprecaCommand | ( | LAMMPS_NS::LAMMPS * | lmp, |
std::vector< std::string > & | commands, | ||
PaprecaConfig & | papreca_config | ||
) |
General PAPRECA input function that redirects commands to the relevant command function.
[in] | lmp | pointer to previously initialized LAMMPS object. |
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeRandomDepovecsCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the type of deposition vector (random or above) in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeRandomDiffvecsCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the type of diffusion vector in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeRandomSeedCommand | ( | LAMMPS_NS::LAMMPS * | lmp, |
std::vector< std::string > & | commands, | ||
PaprecaConfig & | papreca_config | ||
) |
Initializes the LAMMPS random number generator (RanMars) of the PAPRECA::PaprecaConfig object using the input random seed.
[in] | lmp | pointer to LAMMPS object. |
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeRestartFreqCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Enables the generation of LAMMPS restart files with a give frequency in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeSigmasOptionsCommand | ( | LAMMPS_NS::LAMMPS * | lmp, |
std::vector< std::string > & | commands, | ||
PaprecaConfig & | papreca_config | ||
) |
Sets the sigma options in the PAPRECA::PaprecaConfig object.
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeSpeciesMaxBondsCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets a maximum number of bonds for a species in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeSpeciesMaxBondTypesCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets a maximum number of bonds of specific bond type for an atom type in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeTimeEndCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the ending target time in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::executeTrajectoryDurationCommand | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config | ||
) |
Sets the LAMMPS trajectory duration in the PAPRECA::PaprecaConfig object.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
const bool PAPRECA::feCandidateHas4PO4Neibs | ( | PaprecaConfig & | papreca_config, |
PredefinedDiffusionHop * | diff_template, | ||
LAMMPS_NS::tagint * | atom_ids, | ||
int * | atom_types, | ||
int * | neighbors, | ||
int & | neighbors_num, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Scans neighbors of a diffusion parent atom to check whether 4 distinct PO4 (phosphate) structures exist in the neighbors list. This function is only called when the custom PAPRECA::PredefinedDiffusionHop style: Fe_4PO4neib is active.
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | diff_template | Diffusion Hop template (PAPRECA::PredefinedDiffusionHop) as initialized by the user (in the PAPRECA input file). |
[in] | atom_ids | LAMMPS atom ids. |
[in] | atom_types | LAMMPS atom types. |
[in] | neighbors | LAMMPS array storing the IDs of neighbors of the current (parent) atom type. |
[in] | neighbors_num | total number of neighbors. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
void PAPRECA::fillAndSortIndexedRatesVec | ( | double * | arr, |
const int & | arr_size, | ||
DOUBLE2INTPAIR_VEC & | rates_indexed | ||
) |
Uses std::sort to sort a double array. Typically used to sort rates vector in ascending order. Sorting helps with faster event selection. It also ensures repeatability (if the same random number seed is used), because the event enumeration is maintained even if the domain decomposition numbering changes (i.e., if MPI decides to decompose the domain differently). The sorting is not performed on the array itself but on an auxiliary std::vector< std::pair< double , int > > containers that helps us sort while maintaining the index of events. This is necessary, since each time we stochastically select an event using the rates array, we have to link it (by index) to the correct PAPRECA::Event.
[in] | arr | double array container. |
[in] | arr_size | size of arr. |
[in,out] | rates_indexed | std::vector< std::pair< double , int > > container used for ascending sorting. |
void PAPRECA::fillDelidsLocalVec | ( | LAMMPS_NS::LAMMPS * | lmp, |
const double & | desorb_cut, | ||
std::vector< LAMMPS_NS::tagint > & | delids_local, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Called by deleteDesorbedAtoms() and only when the delete_desorbed algorithm is set to gather_local. The function compares the z-coordinate of an atom. If the atom z-value is higher than desorb_cut, the atom ID is marked for deletion (i.e., inserted in the delids_local container) and it is later deleted along with its bonded atoms (retrieved from atomID2bonds map).
[in] | lmp | pointer to LAMMPS object. |
[in] | desorb_cut | cutoff distance for atom deletion. Atoms whose z-coordinate is equal to or greater than desorb_cut are marked for deletion. |
[in,out] | delids_local | vector of collected atom IDs (delids_local will be different on each MPI process). |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
int PAPRECA::fillDelidsVec | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id, | ||
const double & | desorb_cut, | ||
std::vector< LAMMPS_NS::tagint > & | delids, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Called by deleteDesorbedAtoms() and only when the delete_desorbed algorithm is set to gather_all. The function compares the z-coordinate of an atom. If the atom z-coordinate is higher than desorb_cut, the atom ID is marked for deletion (i.e., inserted in the delids container) and it is deleted along with its bonded atoms (retrieved from atomID2bonds map). Here, a gather operation collects data from all atoms on the master MPI process (i.e., proc_id==0) before comparing the z-coordinates of atoms with desorb_cut.
[in] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | desorb_cut | cutoff distance for atom deletion. Atoms whose z-coordinate is equal to or greater than desorb_cut are marked for deletion. |
[in,out] | delids | vector of collected atom IDs on the master MPI process. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
void PAPRECA::fillDepoDataTransfArr | ( | double * | depo_data, |
Deposition * | depo | ||
) |
Serializes (prepares) data for calls to MPI function related to the execution of PAPRECA::Deposition events.
[in,out] | depo_data | 8-element array of double data for transfer. |
[in] | depo | PAPRECA::Deposition event to be executed. |
void PAPRECA::fillDoubleDiffDataTransfArray | ( | double * | diff_doubledata, |
Diffusion * | diff | ||
) |
Serializes (prepares) double data for calls to MPI function related to the execution of PAPRECA::Diffusion events.
[in,out] | diff_doubledata | 4-element array of double data for transfer. |
[in] | diff | PAPRECA::Diffusion event to be executed. |
void PAPRECA::fillFormTransferDataArr | ( | BondForm * | bond_form, |
int * | form_data | ||
) |
Serializes (prepares) data for calls to MPI function related to the execution of PAPRECA::BondForm events.
[in] | bond_form | PAPRECA::BondForm event to be executed. |
[in,out] | form_data | 2-element array of int data. |
void PAPRECA::fillIntegerDiffDataTransfArray | ( | int * | diff_intdata, |
Diffusion * | diff | ||
) |
Serializes (prepares) integer data for calls to MPI function related to the execution of PAPRECA::Diffusion events.
[in,out] | diff_intdata | 3-element array of integer data for transfer. |
[in] | diff | PAPRECA::Diffusion event to be executed. |
void PAPRECA::fillMassProfilesTotalArrFromMassProfilesLocal | ( | const int & | bins_num, |
const int & | types_num, | ||
double ** | mass_profiles_total, | ||
double ** | mass_profiles_local | ||
) |
Fills global array of mass profiles on master process (i.e., proc_id==0). This is done by invoking MPI_Reduce on the local mass profile arrays.
[in] | bins_num | total number of x-y slices (bins). |
[in] | types_num | total number of atom types. |
[in,out] | mass_profiles_total | global 2D mass profiles array initialized and filled ONLY for the master process (i.e., proc_id==0). |
[in] | mass_profiles_local | local (per MPI process) 2D array of mass profiles. |
void PAPRECA::formBond | ( | LAMMPS_NS::LAMMPS * | lmp, |
const LAMMPS_NS::tagint & | atom1id, | ||
const LAMMPS_NS::tagint & | atom2id, | ||
const int & | bond_type | ||
) |
Bonds atom1 and atom2 with a specific bond type (as defined in the LAMMPS input file).
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | atom1id | ID of the first atom. |
[in] | atom2id | ID of the second atom. |
[in] | bond_type | type of bond to be formed. |
void PAPRECA::gatherAndTrimDelIdsOnDriverProc | ( | const int & | proc_id, |
const int & | nprocs, | ||
std::vector< LAMMPS_NS::tagint > & | delids_local, | ||
std::vector< LAMMPS_NS::tagint > & | delids_global | ||
) |
Called by deleteDesorbedAtoms() and only when the delete_desorbed algorithm is set to gather_local. The gather_local algorithm compares the z-coordinates of atoms with desorb_cut locally (i.e., individually for each MPI process). After all comparisons are performed, this function (i.e., gatherAndTrimDelIdsOnDriverProc() ) communicates information among the MPI processes. Then, the function fills delids_global with the IDs of all atoms marked for deletion and makes sure that no duplicate deletions are performed.
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | number of MPI processes. |
[in] | delids_local | contains the IDs of atoms marked for deletion on a specific MPI process (i.e., the delids_local container is different on each MPI process). |
[in,out] | delids_global | vector containing the IDs of all atoms marked for deletion (contains the same data on each MPI process as casted from the master process). |
Different MPI processes can return a different number of delids so we need to calculate the counts/displacements first and then use Gatherv
double PAPRECA::get3DSqrDistWithPBC | ( | LAMMPS_NS::LAMMPS * | lmp, |
const double * | x1, | ||
const double * | x2 | ||
) |
Calculates and returns the distances between 2 points, while accounting for any Periodic Boundary Conditions (PBC) in the system.
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | x1 | array of coordinates of atom1. |
[in] | x2 | array of coordinates of atom2. |
double PAPRECA::getBinWidthFromElementalDistributions | ( | std::vector< std::string > & | commands, |
int & | current_pos | ||
) |
void PAPRECA::getBondBreakingEventsFromAtom | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
const int & | iatom, | ||
int * | neighbors, | ||
int & | neighbors_num, | ||
std::vector< Event * > & | events_local, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Scans all atoms in the current MPI process for PAPRECA::BondBreak (a.k.a. PAPRECA::PredefinedReaction) events. Discovered events are inserted in the events_local vector of PAPRECA::Event objects (storing all events of the current MPI process).
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | iatom | local (LAMMPS, per MPI process) index of candidate atom. |
[in] | neighbors | array storing the IDs of LAMMPS neighbors of the candidate atom. |
[in] | neighbors_num | total number of neighbors of the candidate atom. |
[in,out] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
void PAPRECA::getBondFormEventsFromAtom | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
const int & | iatom, | ||
int * | neighbors, | ||
int & | neighbors_num, | ||
std::vector< Event * > & | events_local, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Scans a half neighbor list and detect bond formation (a.k.a. PAPRECA::PredefinedBondForm). Detected events are inserted in events_local vector, which is a container of PAPRECA::Event object and stores the events detected on the current MPI process.
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | iatom | local index of current atom. |
[in] | neighbors | IDs of neighbors of the current atom. |
[in] | neighbors_num | number of neighbors of the current atom. |
[in,out] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
const std::string PAPRECA::getConcatenatedStringFromStringsVector | ( | const std::vector< std::string > & | strings, |
size_t | start = -1 , |
||
size_t | end = -1 |
||
) |
Receives an std::vector of strings and returns a concatenated std::string (i.e., a string combining all elements of the std::vector< std::string > container)
[in] | strings | std::vector< std::string > container |
[in] | start | concatenate strings starting for this index in the strings vector. |
[in] | end | concatenate strings whose index in the strings vector is less than end. |
const std::string PAPRECA::getConcatenatedStringWithSpacesFromStringsVector | ( | const std::vector< std::string > & | strings, |
size_t | start = -1 , |
||
size_t | end = -1 |
||
) |
Same as getConcatenatedStringFromStringsVector(). However, the present function also adds a space character between each concatenated string (received from std::vector< string > ).
[in] | strings | std::vector< std::string > container |
[in] | start | concatenate strings starting for this index in the strings vector. |
[in] | end | concatenate strings whose index in the strings vector is less than end. |
void PAPRECA::getDepoEventsFromAtom | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
const int & | proc_id, | ||
const int & | nprocs, | ||
const int & | iatom, | ||
int * | neighbors, | ||
int & | neighbors_num, | ||
double & | film_height, | ||
std::vector< Event * > & | events_local | ||
) |
Checks if a system atom is parent to a deposition event (a.k.a. PAPRECA::PredefinedDeposition). Detected deposition events are inserted in the events_local vector of PAPRECA::Event objects (storing all events of the current MPI process).
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | iatom | local index of current atom. |
[in] | neighbors | array containing the IDs of the neighbors to the current atom. |
[in] | neighbors_num | total number of neighbors |
[in] | film_height | film height at the current PAPRECA step. |
[in,out] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
void PAPRECA::getDepoPointCandidateCoords | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
double * | iatom_xyz, | ||
double * | candidate_xyz, | ||
PredefinedDeposition * | depo_template | ||
) |
Fills candidate_xyz array with the coordinates of the deposition candidate. The deposition candidate coordinates DO NOT coincide with the parent atom coordinates. Depending on the user input settings in PAPRECA input, the deposition candidate coordinates can either be directly above the parent atom or on the surface of the upper hemisphere centered at the parent atom coordinates.
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | iatom_xyz | coordinates of parent atom. |
[in,out] | candidate_xyz | coordinates of candidate point for the deposition event. |
[in] | depo_template | Deposition template (PAPRECA::PredefinedDeposition) as initialized by the user (in the PAPRECA input file). |
double PAPRECA::getDepoRateFromHertzKnudsen | ( | const double & | pressure_in, |
const double & | ads_area_in, | ||
const double & | ads_mass_in, | ||
const double & | temperature_in | ||
) |
This function outputs the deposition base rate using the Hertz-Knudsen equation (kinetic theory of gases).
[in] | pressure_in | Partial Pressure in [Bar] units. |
[in] | ads_area_in | Adsorption Site Area in [Angstrom^2] units. |
[in] | ads_mass_in | Molecule mass in [g/mol] units. |
[in] | temperature_in | temperature in [K]. |
double PAPRECA::getDesorptionRate | ( | const double & | activation_energy, |
const double & | temperature | ||
) |
Calculates a desorption rate based on an Arrhenius-type equation.
[in] | activation_energy | Activation Energy in [kcal/mol] units. |
[in] | temperature | Temperature in [K] units. |
void PAPRECA::getDiffEventsFromAtom | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
const int & | iatom, | ||
int * | neighbors, | ||
int & | neighbors_num, | ||
std::vector< Event * > & | events_local, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Checks if the current atom is candidate to diffusion events (a.k.a. PAPRECA::PredefinedDiffusionHop). Detected deposition events are inserted in the events_local vector of PAPRECA::Event objects (storing all events of the current MPI process).
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | iatom | local index of current atom. |
[in] | neighbors | IDs of neighbors of the current atom. |
[in] | neighbors_num | number of neighbors of the current atom. |
[in,out] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
void PAPRECA::getDiffPointCandidateCoords | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
double * | iatom_xyz, | ||
double * | candidate_xyz, | ||
PredefinedDiffusionHop * | diff_template | ||
) |
Calculates the diffusion point (vacancy) coordinates for a given parent atom. Depending on the user settings, the diffusion point can be directly above the parent atom or at the surface of a sphere centered at the coordinates of the parent atom.
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | iatom_xyz | coordinates of parent atom. |
[in,out] | candidate_xyz | coordinates of diffusion point. |
[in] | diff_template | Diffusion template (PAPRECA::PredefinedDiffusionHop) as initialized by the user (in the PAPRECA input file). |
void PAPRECA::getFilmHeightFromMassBinsMethod | ( | PaprecaConfig & | papreca_config, |
LAMMPS_NS::LAMMPS * | lmp, | ||
const int & | proc_id, | ||
double & | film_height, | ||
double ** | mass_profiles_total, | ||
const double & | local_mass, | ||
double * | atom_mass, | ||
const int & | bins_num, | ||
const int & | types_num, | ||
const double & | bin_width | ||
) |
Uses the global 2D mass profiles array (initialized/filled in the master process) to calculate the film height in the current step.
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | lmp | pointer to LAMMPS object. |
[in] | proc_id | id of current MPI process. |
[in,out] | film_height | height of current step. |
[in] | mass_profiles_total | global 2D array of mass_profiles (initialized/filled in master MPI process). |
[in] | local_mass | total mass in the current MPI process. |
[in] | atom_mass | LAMMPS array storing the masses of atom_types. |
[in] | bins_num | total number of x-y slices (bins). |
[in] | types_num | total number of atom types. |
[in] | bin_width | width of each x-y slice (bin). |
double PAPRECA::getLocalRate | ( | std::vector< Event * > & | events_local, |
PaprecaConfig & | papreca_config | ||
) |
Calculates the total local (on the current MPI process) rate by summing the rates of the detected local events (contained in the PAPRECA::Event objects vector, events_local).
[in] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
int PAPRECA::getMaskedNeibIndex | ( | int * | neighbors, |
int & | j | ||
) |
Neighbor indexes in neighbor lists are masked with extra bits. This function unmasks the index and returns a valid local index to avoid segmentation faults.
[in] | neighbors | array of IDs of neighbors. |
[in] | j | neighbor index. |
void PAPRECA::getMolCoords | ( | LAMMPS_NS::LAMMPS * | lmp, |
double ** | mol_xyz, | ||
double ** | mol_dx, | ||
const int & | mol_natoms, | ||
double * | candidate_center | ||
) |
Fills mol_xyz with the coordinates of the candidate molecule for deposition.
[in] | lmp | pointer to LAMMPS object. |
[in,out] | mol_xyz | Temporary array storing the candidate (for deposition) molecule coordinates. Container mol_xyz is only used to check for possible collisions between the inserted molecule and system molecules. |
[in] | mol_dx | array storing the x,y, and z distances of each molecule atom from the molecule center. See molecule.h and molecule.cpp files in the LAMMPS source directory for more information. |
[in] | mol_natoms | total number of molecule atoms. |
[in] | candidate_center | coordinates of the candidate deposition point. |
|
inline |
Receives a molecule name and returns its molecule index. This allows one to access a LAMMPS molecule by using: lmp->atom->molecules[imol] (where imol is the molecule index).
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | mol_name | valid molecule names (as defined in the LAMMPS input file. e.g., if you used this command: "molecule mmmTCP ./TCP.xyz", then mol_name should be "mmmTCP"). |
void PAPRECA::getMonoDesEventsFromAtom | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config, | ||
const int & | iatom, | ||
std::vector< Event * > & | events_local, | ||
ATOM2BONDS_MAP & | atomID2bonds | ||
) |
Checks if the parent atom is candidate to monoatomic desorption events (a.k.a. PAPRECA::PredefinedMonoatomicDesorption). Detected events are inserted in events_local vector, which is a container of PAPRECA::Event object and stores the events detected on the current MPI process.
[in] | lmp | pointer to LAMMPS object. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | iatom | local index of current atom. |
[in,out] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
const int PAPRECA::getMPIRank | ( | MPI_Comm | communicator | ) |
Wrapper around MPI_Comm_rank. Returns the MPI process ID of the specific process calling the function.
[in] | communicator | MPI communicator. |
double PAPRECA::getRateFromArrhenius | ( | const double & | activation_energy, |
const double & | attempt_freq, | ||
const double & | temperature | ||
) |
This function uses the classic Hertzian equation to output the rate of a rare event.
[in] | activation_energy | Activation energy in [kcal/mol] units. |
[in] | attempt_freq | Attempt frequency in [Hz (1/s)] units. |
[in] | temperature | in [K] units. |
double PAPRECA::getRateFromInputRateOptions | ( | std::vector< std::string > & | commands, |
int & | current_pos | ||
) |
Calculates the predefined event rate by calling the chosen (in the command-line) function in rates_calc.h and returns it to the caller function.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | current_pos | position on the commands vector. This variable gets updated after the call to the function to denote the exiting position in the commands vector. |
double PAPRECA::getStickingCoeffFromDepositionEventOptions | ( | std::vector< std::string > & | commands, |
int & | current_pos | ||
) |
Extracts the sticking coefficient from the commands vector returns it to the caller function.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | current_pos | position on the commands vector. This variable gets updated after the call to the function to denote the exiting position in the commands vector. |
const int PAPRECA::getStringPosInStringVec | ( | const std::string & | string2check, |
const std::vector< std::string > & | strings | ||
) |
Returns the position of a string in a vector of strings (std::vector< std::string > container).
[in] | string2check | string to identify. |
[in] | strings | std::vector< std::string > container. |
std::vector< T > PAPRECA::getSubVectorFromVector | ( | const std::vector< T > & | vec, |
size_t | start, | ||
size_t | end | ||
) |
Receives an std::vector and returns a new std::vector containing all the elements of the original vector whose indexes are equal to or greater than "start" and smaller than "end".
[in] | vec | input std::vector. |
[in] | start | starting index (i.e., include all elements whose indexes are equal to or greater than start). |
[in] | end | ending index (i.e., include elements whose indexes are smaller than end). |
T PAPRECA::getSumOfVecElements | ( | const std::vector< T > & | vec | ) |
Receives an std::vector and returns the sum of its elements.
[in] | vec | std::vector |
const bool PAPRECA::headAtomIsCatalyzed | ( | PredefinedReaction * | reaction_template, |
int * | atom_types, | ||
int * | neighbors, | ||
int & | neighbors_num | ||
) |
Scans neighbors of parent predefined reaction candidate to check if catalyzing types exist in the neighborhood. The catalyzing types are provided in the PAPRECA input file.
[in] | reaction_template | Reaction template (PAPRECA::PredefinedReaction or PAPRECA::PredefinedBondForm) as initialized by the user (in the PAPRECA input file). |
[in] | atom_types | number of LAMMPS atom types. |
[in] | neighbors | array storing the neighbor IDs of the parent PAPRECA::PredefinedReaction candidate atom. |
[in] | neighbors_num | number of neighbors of the parent candidate atom. |
void PAPRECA::initAndGatherBondsList | ( | LAMMPS_NS::LAMMPS * | lmp, |
LAMMPS_NS::tagint ** | bonds_list, | ||
LAMMPS_NS::bigint & | bonds_num | ||
) |
Calls lammps_gather_bonds to fill array bonds_list. The container bonds_list contains the IDs of bonds of the LAMMPS object (lmp).
[in] | lmp | pointer to LAMMPS object. |
[in,out] | bonds_list | array containing the IDs of bonds. |
[in,out] | bonds_num | number of bonds in the bonds_list array. |
void PAPRECA::initializeLMP | ( | LAMMPS_NS::LAMMPS ** | lmp | ) |
Initialize LAMMPS.
Instantiates a LAMMPS object
[in,out] | lmp | pointer to LAMMPS object. |
double ** PAPRECA::initMassProfilesArr | ( | const int & | types_num, |
const int & | bins_num | ||
) |
Initializes a 2D array of mass profiles (i.e., arr[bins_num][types_num]).
[in] | types_num | total number of atom types. |
[in] | bins_num | total number of x-y slices (bins). |
void PAPRECA::initMolCoordsArr | ( | double *** | mol_xyz, |
const int & | mol_natoms | ||
) |
Initializes the molecule coordinates array used to check for collisions between the deposited molecule and system atoms.
[in,out] | mol_xyz | array storing the x,y, and z distances of each molecule atom from the molecule center. See molecule.h and molecule.cpp files in the LAMMPS source directory for more information. |
[in] | mol_natoms | total number of molecule atoms. |
void PAPRECA::initType2SigmaFromLammpsPairCoeffs | ( | LAMMPS_NS::LAMMPS * | lmp, |
INTPAIR2DOUBLE_MAP & | type2sigma | ||
) |
Retrieves the pairstyle sigmas from LAMMPS to initialize the type sigmas in PAPRECA.
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in,out] | type2sigma | PAPRECA::INTPAIR2DOUBLE_MAP mapping a pair of atom types to their corresponding sigma. |
void PAPRECA::insertMolecule | ( | LAMMPS_NS::LAMMPS * | lmp, |
const double | site_pos[3], | ||
const double | rot_pos[3], | ||
const double & | rot_theta, | ||
const int & | type_offset, | ||
const char * | mol_name | ||
) |
Inserts a molecule in the current LAMMPS system.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | site_pos | array of insertion site coordinates. LAMMPS will insert a new molecule such that its center-of-mass coincides with site_pos. |
[in] | rot_pos | array of coordinates of the center of rotation. |
[in] | rot_theta | angle of molecule rotation. |
[in] | type_offset | add type_offset to the atom types of the molecule (as defined in the LAMMPS molecule template file). WARNING: For almost all cases you should use a value of 0 to ensure that the inserted atom types are identical to the atom types as defined in the LAMMPS molecule input file. |
[in] | mol_name | name of molecule (as defined in the LAMMPS input command. e.g., if you used this command: "molecule mmmTCP ./TCP.xyz", your mol_name should be "mmmTCP" ). |
void PAPRECA::loopAtomsAndIdentifyEvents | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id, | ||
int & | nprocs, | ||
const int & | KMC_loopid, | ||
PaprecaConfig & | papreca_config, | ||
std::vector< Event * > & | events_local, | ||
ATOM2BONDS_MAP & | atomID2bonds, | ||
double & | film_height | ||
) |
1) Calculates film height (if that is requested by the user). 2) Discovers events and inserts them in the PAPRECA::Event objects vector (storing all events detected on the current MPI process). PAPRECA::BondBreak, PAPRECA::Deposition, PAPRECA::MonoatomicDesorption, and PAPRECA::Diffusion are discovered through a full neighbors list. PAPRECA::BondForm are detected using a half-neighbors list. This is done because bondbreak, deposition, monoatomic desorption, and diffusion events perform collision checks which require full neighbor lists (because we need to know all the neighbors of each atom). On the other hand, it is more efficient to check for PAPRECA:BondForm events using a half neighbors list that include each neighbor pair once.
[in] | lmp | pointer to LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | KMC_loopid | current PAPRECA step. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in,out] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
[in,out] | film_height | film height at current PAPRECA step. |
bool PAPRECA::mappingExists | ( | const std::unordered_map< KeyType, ValueType, Hash > & | map, |
const KeyType & | key | ||
) |
Searches if an element is included in an unordered map.
[in] | map | unordered_map to be searched. |
[in] | key | element to be searched in the unordered map. |
void PAPRECA::MPIBcastAndExecuteCommand | ( | LAMMPS_NS::LAMMPS * | lmp, |
std::string & | command | ||
) |
void PAPRECA::printStepInfo | ( | PaprecaConfig & | papreca_config, |
const int & | KMC_loopid, | ||
const double & | time, | ||
const double & | film_height, | ||
const double & | proc_rates_sum | ||
) |
Prints essential PAPRECA information to the screen/terminal.
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | KMC_loopid | number of current PAPRECA step. |
[in] | time | current time. |
[in] | film_height | current height. |
[in] | proc_rates_sum | total (gathered/sum from all MPI processes) event rate. |
void PAPRECA::processBinWidthOptionForElementalDistributions | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config, | ||
int & | current_pos | ||
) |
Sets the bind width variable of the PAPRECA::PaprecaConfig object (based on the input command-line).
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
[in,out] | current_pos | position on the commands vector. This variable gets updated after the call to the function to denote the exiting position in the commands vector. |
void PAPRECA::processBondLimitOption | ( | std::vector< std::string > & | commands, |
int & | current_pos, | ||
double & | length_equil, | ||
double & | length_perc | ||
) |
void PAPRECA::processCatalyzedOption | ( | std::vector< std::string > & | commands, |
int & | current_pos, | ||
std::vector< int > & | catalyzing_types | ||
) |
Checks if the optional "catalyzed" keyword is present in a command. If yes then the catalyzing types are inserted in the catalyzing_types std::vector< int > container.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | current_pos | position on the commands vector. This variable gets updated after the call to the function to denote the exiting position in the commands vector. |
[in,out] | catalyzing_types | vector containing the catalyzing types for a specific command (typically create_BondBreak or create_BondForm). |
void PAPRECA::processCustomDiffEventOptions | ( | std::vector< std::string > & | commands, |
int & | current_pos, | ||
std::string & | custom_style, | ||
std::vector< int > & | style_atomtypes | ||
) |
Initializes custom_style and style_atomtypes for custom diffusion events.
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | current_pos | position on the commands vector. This variable gets updated after the call to the function to denote the exiting position in the commands vector. |
[in,out] | custom_style | type of custom style. |
[in,out] | style_atomtypes | potentially useful auxiliary vector storing atom types and passing them to the main function of papreca.cpp. |
std::vector< std::string > PAPRECA::processLine | ( | char * | line | ) |
Receives a line from the PAPRECA input file, removes blank characters and comments (i.e., text to the right of a '#' character). Then, converts the char* to an std::vector< std::string > container that stores each word/command/keyword on a separate std::string.
[in] | line | raw (unprocessed command-line sent from the PAPRECA input file. |
void PAPRECA::processSigmaMixOptions | ( | std::vector< std::string > & | commands, |
int & | current_pos | ||
) |
void PAPRECA::processSigmaMixOptions | ( | std::vector< std::string > & | commands, |
PaprecaConfig & | papreca_config, | ||
int & | current_pos | ||
) |
Sets the sigma mix variable of the PAPRECA::PaprecaConfig object (based on the input command-line).
[in] | commands | trimmed/processed vector of strings. This is effectively the entire command line with each vector element (i.e., std::string) being a single word/number. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
[in,out] | current_pos | position on the commands vector. This variable gets updated after the call to the function to denote the exiting position in the commands vector. |
void PAPRECA::pushArrayToVector | ( | const T * | arr, |
std::vector< T > & | vec | ||
) |
Receives an std::vector and an array and pushes all the array elements to the vector.
[in] | arr | input array. |
[in,out] | vec | std::vector container. |
void PAPRECA::readInputAndInitPaprecaConfig | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | proc_id, | ||
const char * | file_name, | ||
PaprecaConfig & | papreca_config | ||
) |
Reads PAPRECA input line-by-line on the master MPI process and casts each line to all other processes. Then, all MPI processes execute the relevant PAPRECA command in executePaprecaCommand(). In the end, all MPI processes initialize an identical copy of PAPRECA::PaprecaConfig (containing the settings and global variables for the PAPRECA run).
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | proc_id | ID of current MPI process. |
[in] | file_name | name of PAPRECA input file. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::readLMPinput | ( | const std::string & | lmp_input, |
LAMMPS_NS::LAMMPS * | lmp | ||
) |
Reads LAMMPS input and initializes LAMMPS parameters.
[in] | lmp_input | name of LAMMPS input file (command-line argument). |
[in,out] | lmp | pointer to LAMMPS object. |
void PAPRECA::remap3DArrayInPeriodicBox | ( | LAMMPS_NS::LAMMPS * | lmp, |
double * | arr | ||
) |
Receives an 1D array of 3 elements (i.e., x,y, and z coordinates) and remaps it inside the existing periodic box (as set in the provided LAMMPS object).
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | arr | array of coordinates. |
void PAPRECA::resetMobileAtomsGroups | ( | LAMMPS_NS::LAMMPS * | lmp, |
const std::vector< int > & | fluid_atomtypes | ||
) |
Clears the fluid group (containing the fluid atom types, as defined in the PAPRECA input file) and redefines it to. This ensures that all atoms of fluid types are included in the fluid group.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | fluid_atomtypes | vector containing the fluid atom types. |
void PAPRECA::runLammps | ( | LAMMPS_NS::LAMMPS * | lmp, |
const int & | timesteps_num | ||
) |
Runs LAMMPS trajectory for given (previously initialized) LAMMPS object.
[in,out] | lmp | pointer to previously instantiated LAMMPS object. |
[in] | timesteps_num | number of timesteps to execute. |
int PAPRECA::selectAndExecuteEvent | ( | LAMMPS_NS::LAMMPS * | lmp, |
int & | KMC_loopid, | ||
double & | time, | ||
char * | event_type, | ||
int & | proc_id, | ||
int & | nprocs, | ||
PaprecaConfig & | papreca_config, | ||
std::vector< Event * > & | events_local, | ||
ATOM2BONDS_MAP & | atomID2bonds, | ||
double & | film_height | ||
) |
Gets total event rate by gathering all local rates (i.e., sum of rates on a single MPI process). Then, selects an event MPI process using the N-FOLD way. Afterwards, an event is chosen from the selected MPI process and executed on all procs.
[in,out] | lmp | pointer to LAMMPS object. |
[in] | KMC_loopid | current PAPRECA step. |
[in,out] | time | current time. |
[in,out] | event_type | type of selected event. |
[in] | proc_id | ID of current MPI process. |
[in] | nprocs | total number of MPI processes. |
[in] | papreca_config | object of the PAPRECA::PaprecaConfig class that stores global variables and settings for the current PAPRECA run. |
[in] | events_local | vector containing all the PAPRECA::Event objects for a specific MPI process. |
[in] | atomID2bonds | PAPRECA::ATOM2BONDS_MAP container (i.e., std::unordered_map< LAMMPS_NS::tagint parent_atomID , std::vector< PAPRECA::Bond > >). The atomID2bonds container provides direct access to all the bonds of the parent atom. |
[in] | film_height | current height. |
int PAPRECA::selectProcessStochastically | ( | double * | arr, |
const int & | arr_size, | ||
double & | rnum, | ||
double & | rates_sum | ||
) |
Selects a process stochastically (through the classic N-FOLD way selection process). Can be used to select an MPI proc that fires an event or can be used to select the next kMC event inside the proc. The event selection is performed on a std::vector< std::pair< double , int > > container that holds the rate and the associated index of events.
[in] | arr | array of rates. |
[in] | arr_size | size of arr. |
[in] | rnum | uniformly distributed (between 0 and 1) pseudorandom number (usually drawn on master proc). |
[in] | rates_sum | total sum of rates |
void PAPRECA::setTimeUnitsConversionConstant | ( | LAMMPS_NS::LAMMPS * | lmp, |
PaprecaConfig & | papreca_config | ||
) |
Retrieves the LAMMPS units style from the provided LAMMPS object (lmp) and sets the conversion time variable of the PAPRECA::PaprecaConfig object.
[in] | lmp | pointer to previously instantiated LAMMPS object. |
[in,out] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::setupMPI | ( | int * | narg, |
char *** | arg, | ||
int * | nprocs, | ||
int * | proc_id | ||
) |
Initializes MPI environment using command-line arguments passed during the program invocation from the terminal.
[in] | narg | number of command-line arguments. |
[in] | arg | array containing command-line arguments. |
[in,out] | nprocs | number of MPI processes. |
[in,out] | proc_id | identifier of current MPI process. |
const bool PAPRECA::string2Bool | ( | std::string & | string | ) |
Converts an input std::string to a boolean.
[in] | string | std::string to be converted. |
const double PAPRECA::string2Double | ( | std::string & | string | ) |
Converts an input std::string to a double number.
[in] | string | std::string to be converted. |
const int PAPRECA::string2Int | ( | std::string & | string | ) |
Converts an std::string to an integer.
[in] | string | std::string to be converted. |
const unsigned long int PAPRECA::string2UnsignedLongInt | ( | std::string & | string | ) |
Converts input string to unsigned long int.
[in] | string | std::string to be converted. |
const bool PAPRECA::stringIsBool | ( | const std::string & | string | ) |
Checks if the string is a bool (i.e., yes or no).
[in] | string | std::string to be checked. |
bool PAPRECA::stringIsIntNumber | ( | const std::string & | string | ) |
This function checks if the string is a number.
[in] | string | std::string to be checked. |
bool PAPRECA::stringIsNumber | ( | const std::string & | string | ) |
This function checks if the string is a number. We allow negative/positive signs and decimal points in there so we can potentially indetify double numbers or negative numbers or exponentials as well.
[in] | string | input std::string. |
void PAPRECA::warn4IllegalRuns | ( | const int & | proc_id, |
PaprecaConfig & | papreca_config | ||
) |
After processing all commands, this function is called to warn all MPI processes for potentially illegal PAPRECA runs.
[in] | proc_id | ID of current MPI process. |
[in] | papreca_config | previously instantiated PAPRECA::PaprecaConfig object storing the settings and global variables for the PAPRECA simulation. |
void PAPRECA::warnAll | ( | MPI_Comm | communicator, |
const std::string & | message | ||
) |
Throws a warning to all MPI processes.
[in] | communicator | MPI communicator. |
[in] | message | std::string of the warning message. |
void PAPRECA::warnOne | ( | MPI_Comm | communicator, |
const std::string & | message | ||
) |
Throws a warning on a specific MPI process.
[in] | communicator | MPI communicator. |
[in] | message | std::string of the warning message. |