PAPRECA hybrid off-lattice kMC/MD simulator  2.0.0 (17 September 2024)
PAPRECA::PaprecaConfig Class Reference

Class storing settings and global variables for the PAPRECA run. More...

#include <papreca_config.h>

Collaboration diagram for PAPRECA::PaprecaConfig:

Public Member Functions

 PaprecaConfig ()
 
 ~PaprecaConfig ()
 
void setKMCsteps (const unsigned long int &KMC_steps_in)
 
const unsigned long int & getKMCsteps () const
 
void setKMCperMD (const unsigned long int &KMC_per_MD_in)
 
const unsigned long int & getKMCperMD () const
 
void setTimeEnd (const double &time_end_in)
 
const double & getTimeEnd ()
 
void initRanNumGenerator (LAMMPS_NS::LAMMPS *lmp, const int &seed)
 
const double getUniformRanNum () const
 
const bool ranNumGeneratorIsInitialized () const
 
void setFluidAtomTypes (const std::vector< int > &fluid_atomtypes_in)
 
const std::vector< int > & getFluidAtomTypes () const
 
void setFrozenAtomTypes (const std::vector< int > &frozen_atomtypes_in)
 
const std::vector< int > & getFrozenAtomTypes () const
 
PredefinedReactiongetReactionFromBondType (const int &bond_type)
 
PredefinedBondFormgetBondFormFromAtomTypesPair (const INT_PAIR &types_pair)
 
int getMaxBondsFromSpecies (const int &atom_type)
 
int getMaxBondTypesOfSpecies (const int &atom_type, const int &bond_type)
 
PredefinedDiffusionHopgetDiffusionHopFromAtomType (const int &atom_type)
 
PredefinedDepositiongetDepositionFromParentAtomType (const int &atom_type)
 
PredefinedMonoatomicDesorptiongetMonoatomicDesorptionFromAtomType (const int &atom_type)
 
void initPredefinedReaction (const int &atom1_type, const int &atom2_type, const int &bond_type, const double &rate, const std::vector< int > &catalyzing_types)
 
void initPredefinedBondForm (const int &atom1_type, const int &atom2_type, const int &bond_type, const double &bond_dist, const int &delete_atoms, const int &lone_candidates, const bool &same_mol, const double &rate, const std::vector< int > &catalyzing_types)
 
void initPredefinedDiffusionHop (const int &parent_type, const double &insertion_vel, const double &diff_dist, const bool &is_displacive, const int &diffused_type, const double &rate, const std::string &custom_style, const std::vector< int > &custom_atomtypes)
 
void initPredefinedDeposition (LAMMPS_NS::LAMMPS *lmp, const int &parent_type, const double &depo_offset, const double &insertion_vel, const std::string &adsorbate_name, const double &rate, const bool &variable_sticking, const double &sticking_coeff)
 
void initPredefinedMonoatomicDesorption (const int &parent_type, const double &rate)
 
void setSpeciesMaxBonds (const int &species, const int &bonds_max)
 
void setSpeciesMaxBondTypes (const int &species, const int &bond_type, const int &bonds_max)
 
void calcStickingCoeffs ()
 
const bool predefinedCatalogHasBondBreakEvents () const
 
const bool predefinedCatalogHasBondFormEvents () const
 
const bool predefinedCatalogHasDiffusionHopEvents () const
 
const bool predefinedCatalogHasDepositionEvents () const
 
const bool predefinedCatalogHasMonoDesEvents () const
 
const bool predefinedCatalogIsEmpty () const
 
void setRandomDepoVecs (const bool &random_depovecs_in)
 
const bool & depoVecsAreRandom () const
 
void setRandomDiffVecs (const bool &random_diffvecs_in)
 
const bool & diffVecsAreRandom () const
 
void setRandomDiffVecsStyle (const std::string &diffvecs_style_in)
 
const std::string & getRandomDiffVecsStyle () const
 
void setDepoHeights (const double &height_deposcan_in, const double &height_deporeject_in)
 
const double & getHeightDepoScan () const
 
const double & getHeightDepoReject () const
 
void setDesorptionHeight (const double &desorb_cut_in)
 
const double & getDesorptionHeight () const
 
void setDesorbDelMax (const int &desorb_delmax_in)
 
const int & getDesorbDelMax () const
 
void setDesorptionStyle (const std::string &desorb_style_in)
 
const std::string & getDesorptionStyle () const
 
void setHeightMethod (const std::string &height_method_in)
 
const std::string & getHeightMethod () const
 
void setHeightPercentage (const double &height_percentage_in)
 
const double & getHeightPercentage () const
 
void setBinWidth (const double &bin_width_in)
 
const double & getBinWidth () const
 
void initSigmasFromLammps (LAMMPS_NS::LAMMPS *lmp)
 
void setSpeciesPair2Sigma (const int &species1, const int &species2, const double &sigma)
 
void setSigmaStyle (const std::string &sigmastyle_in)
 
const std::string & getSigmaStyle () const
 
void setSigmaMix (const std::string &mixstyle_in)
 
const std::string & getSigmaMixStyle () const
 
void mixSigmas (LAMMPS_NS::LAMMPS *lmp)
 
const double getSigmaFromAtomTypes (const int &atom1_type, const int &atom2_type)
 
const bool type2SigmaMapIsEmpty () const
 
void setMinimize1 (const std::string &minimize1_in)
 
const std::string & getMinimize1 () const
 
void setMinimize2 (const std::string &minimize2_in)
 
const std::string & getMinimize2 () const
 
void setTrajDuration (const int &traj_duration_in)
 
const int & getTrajDuration () const
 
void setCtimeConvert (const double &c_convert_in)
 
const double & getCtimeConvert ()
 
void setNeibLists (const std::string &neiblist_half_in, const std::string &neiblist_full_in)
 
const std::string & getHalfNeibListName () const
 
const std::string & getFullNeibListName () const
 
LoggetLogFile ()
 
HeightVtimegetHeightVtimeFile ()
 
SurfaceCoveragegetSurfaceCoverageFile ()
 
void calcSurfaceCoverage ()
 
ElementalDistributiongetElementalDistributionsFile ()
 
ExecTimegetExecTimeFile ()
 
void setupExportFiles (const int &proc_id)
 
void setHybridStartTimeStamp4ExecTimeFile (const int &KMC_loopid)
 
void calcHybridAndKMCTimes4ExecTimeFile (const int &nprocs, const int &KMC_loopid)
 
void setMDTimeStamp4ExecTimeFile (const int &KMC_loopid)
 
void calcMDTime4ExecTimeFile (const int &nprocs, const int &KMC_loopid)
 
void appendExportFiles (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, const double &time, const char *event_type, const double &film_height, const int &KMC_loopid)
 
void dumpElementalDistributionFile (LAMMPS_NS::LAMMPS *lmp, const int &proc_id, const int &KMC_loopid, double **mass_profiles_total, double *atom_mass, const int &bins_num, const int &types_num)
 
void closeExportFiles (const int &proc_id)
 
void setRestartDumpFreq (const int &restart_dumpfreq_in)
 
const int & getRestartDumpFreq () const
 
void dumpLAMMPSRestart (LAMMPS_NS::LAMMPS *lmp, const int KMC_loopid)
 

Protected Attributes

LAMMPS_NS::RanMars * rnum_gen = NULL
 Points to RanMars object as implemented in LAMMPS (see random_mars.h). The LAMMPS RanMars object is initialized from a user defined random seed in the PAPRECA input file (seed > 0 && seed < 900000000). More...
 
unsigned long int KMC_steps = 0
 perform that many PAPRECA steps. More...
 
unsigned long int KMC_per_MD = 0
 Perform that many KMC(PAPRECA) steps for every MD(LAMMPS) step. More...
 
double time_end = std::numeric_limits< double >::max( )
 Optional parameter. The user can define an ending simulation time. If not set it stays at limits of max, so you effectively never go beyond time_end;. More...
 
std::vector< int > fluid_atomtypes
 stores fluid atom types (as defined in the PAPRECA input file). More...
 
std::vector< int > frozen_atomtypes
 stores frozen atom types (as defined in the PAPRECA input file). More...
 
PredefinedEventsCatalog predefined_catalog
 stores a PAPRECA::PredefinedEventsCatalog. More...
 
bool random_depovecs = false
 Controls deposition sites. If true, the deposition sites are not directly above the parent atom, but on the surface of a sphere of radius depo_offset. More...
 
bool random_diffvecs = false
 Controls diffusion sites. If true, the diffusion sites are not directly above the parent atom, but on the surface of a sphere centered at the parent atom. More...
 
std::string diffvecs_style = "3D"
 Type 2D means that the random diffusion vector is always ABOVE the parent type. Type 3D means that the random diffusion vector can be anywhere in 3D space. Default is 3D as it is the more general case. Note that: we do not use 2D/3D random deposition vectors as it does not really make sense to get a random vector below the parent atom type. More...
 
double height_deposcan = -1
 Scan for deposition events only +- above/below the current film height. Default at -1 which means scan everywhere. More...
 
double height_deporeject = -1
 Reject deposition event above height_current + height_deporeject. Default at -1 which means do not reject anything. More...
 
double desorb_cut = -1
 Atoms above film_height + desorb_cut are deleted. The default value is -1 which means that desorption is disabled. More...
 
int desorb_delmax = std::numeric_limits< int >::max( )
 Maximum number of atoms that can be deleted at once. Initialized at max limits of int so if the user does not set that, the maximum number of deleted atoms will be unlimited. More...
 
std::string desorb_style = ""
 User defined desorption algorithm. Currently, can be gather_local or gather_all (defined in the PAPRECA input file). More...
 
std::string height_method = ""
 Algorithm to calculate height. Currently, only the mass_bins method is supported. More...
 
double height_percentage = 0.0
 Only used for height calculation method: mass_bins. Defines the mass percentage cutoff. The film height is calculated as the first bin whose cumulative mass is above the percentage cutoff. More...
 
double bin_width = 1.0
 Bin width for height calculation or for ElementaDistribution files. Initialized at 1.0, so even if the user does not use the height_caulcation or heightbin_width command, you would be able to calculate mass bins with a width of 1.0. More...
 
INTPAIR2DOUBLE_MAP type2sigma
 maps types of atom types to their corresponding sigma. More...
 
std::string sigma_style = ""
 method for initialization of sigma values (can be manual/LAMMPS). More...
 
std::string sigma_mix =""
 
std::string minimize1 = ""
 stores a valid LAMMPS minimization command to be executed before the LAMMPS trajectory. See: https://docs.lammps.org/minimize.html. More...
 
std::string minimize2 = ""
 stores a valid LAMMPS minimization command to be executed after the LAMMPS trajectory. See: https://docs.lammps.org/minimize.html. More...
 
int traj_duration = 0
 This control the duration of the LAMMPS minimization. The specific LAMMPS fixes have to be defined in the LAMMPS input file. More...
 
double c_time_convert = -1
 This variable is initialized from lmp->update->unit_style to allow conversion of time units to second. The constant is also pre-multiplied with the LAMMPS timestep to allow easy conversion of timesteps (i.e., traj_duration) to time interval (in seconds). More...
 
std::string neiblist_half = ""
 Name of LAMMPS half neighbors list. More...
 
std::string neiblist_full = ""
 Name of LAMMPS full neighbors list. More...
 
Log log_file
 stores a PAPRECA::Log file The log file is always active (no need to be set to active from the input file). More...
 
HeightVtime heightVtime_file
 stores a PAPRECA::HeightVtime file. More...
 
SurfaceCoverage surfcoverage_file
 stores a PAPRECA::SurfaceCoverage file. More...
 
double surface_coverage = 0.0
 stores a surface coverage for easier printing (if necessary). More...
 
ElementalDistribution elementalDistribution_files
 stores the PAPRECA::ElementalDistribution files generated in the simulation. More...
 
ExecTime execTime_file
 stores a PAPRECA::ExecTime file. More...
 
int restart_dumpfreq = std::numeric_limits< int >::max( )
 dump a restart every restart_dumpfreq PAPRECA steps. Initialized at int limits, so if it is not set you virtually never dump restarts (see how restarts are dumped in lammps_wrappers.h of papreca lib). More...
 

Detailed Description

Class storing settings and global variables for the PAPRECA run.

Constructor & Destructor Documentation

◆ PaprecaConfig()

PAPRECA::PaprecaConfig::PaprecaConfig ( )

◆ ~PaprecaConfig()

PAPRECA::PaprecaConfig::~PaprecaConfig ( )

Member Function Documentation

◆ appendExportFiles()

void PAPRECA::PaprecaConfig::appendExportFiles ( LAMMPS_NS::LAMMPS *  lmp,
const int &  proc_id,
const double &  time,
const char *  event_type,
const double &  film_height,
const int &  KMC_loopid 
)

◆ calcHybridAndKMCTimes4ExecTimeFile()

void PAPRECA::PaprecaConfig::calcHybridAndKMCTimes4ExecTimeFile ( const int &  nprocs,
const int &  KMC_loopid 
)

◆ calcMDTime4ExecTimeFile()

void PAPRECA::PaprecaConfig::calcMDTime4ExecTimeFile ( const int &  nprocs,
const int &  KMC_loopid 
)

◆ calcStickingCoeffs()

void PAPRECA::PaprecaConfig::calcStickingCoeffs ( )

◆ calcSurfaceCoverage()

void PAPRECA::PaprecaConfig::calcSurfaceCoverage ( )

◆ closeExportFiles()

void PAPRECA::PaprecaConfig::closeExportFiles ( const int &  proc_id)

◆ depoVecsAreRandom()

const bool & PAPRECA::PaprecaConfig::depoVecsAreRandom ( ) const

◆ diffVecsAreRandom()

const bool & PAPRECA::PaprecaConfig::diffVecsAreRandom ( ) const

◆ dumpElementalDistributionFile()

void PAPRECA::PaprecaConfig::dumpElementalDistributionFile ( LAMMPS_NS::LAMMPS *  lmp,
const int &  proc_id,
const int &  KMC_loopid,
double **  mass_profiles_total,
double *  atom_mass,
const int &  bins_num,
const int &  types_num 
)

◆ dumpLAMMPSRestart()

void PAPRECA::PaprecaConfig::dumpLAMMPSRestart ( LAMMPS_NS::LAMMPS *  lmp,
const int  KMC_loopid 
)

◆ getBinWidth()

const double & PAPRECA::PaprecaConfig::getBinWidth ( ) const

◆ getBondFormFromAtomTypesPair()

PredefinedBondForm * PAPRECA::PaprecaConfig::getBondFormFromAtomTypesPair ( const INT_PAIR types_pair)

Returns a PAPRECA::PredefinedBondForm object as defined in the PAPRECA::PredefinedEventsCatalog.

Parameters
[in]types_pairstd::pair( int , int ) of atom types.
Returns
mapped PAPRECA::PredefinedBondForm object.
Note
each pair of types can only be mapped to one PAPRECA::PredefinedBondForm object.

◆ getCtimeConvert()

const double & PAPRECA::PaprecaConfig::getCtimeConvert ( )

◆ getDepositionFromParentAtomType()

PredefinedDeposition * PAPRECA::PaprecaConfig::getDepositionFromParentAtomType ( const int &  atom_type)

Returns a PAPRECA::PredefinedDeposition object as defined in the PAPRECA::PredefinedEventsCatalog.

Parameters
[in]atom_typetype of atom.
Returns
mapped PAPRECA::PredefinedDeposition as defined in the PAPRECA::PredefinedEventsCatalog.
Note
each atom type can only be mapped to one PAPRECA::PredefinedDeposition.

◆ getDesorbDelMax()

const int & PAPRECA::PaprecaConfig::getDesorbDelMax ( ) const

◆ getDesorptionHeight()

const double & PAPRECA::PaprecaConfig::getDesorptionHeight ( ) const

◆ getDesorptionStyle()

const std::string & PAPRECA::PaprecaConfig::getDesorptionStyle ( ) const

◆ getDiffusionHopFromAtomType()

PredefinedDiffusionHop * PAPRECA::PaprecaConfig::getDiffusionHopFromAtomType ( const int &  atom_type)

Returns a PAPRECA::PredefinedDiffusionHop object as defined in the PAPRECA::PredefinedEventsCatalog.

Parameters
[in]atom_typetype of atom.
Returns
mapped PAPRECA::PredefinedDiffusionHop object.
Note
each atom type can only be mapped to one PAPRECA::PredefinedDiffusionHop.

◆ getElementalDistributionsFile()

ElementalDistribution & PAPRECA::PaprecaConfig::getElementalDistributionsFile ( )

◆ getExecTimeFile()

ExecTime & PAPRECA::PaprecaConfig::getExecTimeFile ( )

◆ getFluidAtomTypes()

const std::vector< int > & PAPRECA::PaprecaConfig::getFluidAtomTypes ( ) const

◆ getFrozenAtomTypes()

const std::vector< int > & PAPRECA::PaprecaConfig::getFrozenAtomTypes ( ) const

◆ getFullNeibListName()

const std::string & PAPRECA::PaprecaConfig::getFullNeibListName ( ) const

◆ getHalfNeibListName()

const std::string & PAPRECA::PaprecaConfig::getHalfNeibListName ( ) const

◆ getHeightDepoReject()

const double & PAPRECA::PaprecaConfig::getHeightDepoReject ( ) const

◆ getHeightDepoScan()

const double & PAPRECA::PaprecaConfig::getHeightDepoScan ( ) const

◆ getHeightMethod()

const std::string & PAPRECA::PaprecaConfig::getHeightMethod ( ) const

◆ getHeightPercentage()

const double & PAPRECA::PaprecaConfig::getHeightPercentage ( ) const

◆ getHeightVtimeFile()

HeightVtime & PAPRECA::PaprecaConfig::getHeightVtimeFile ( )

◆ getKMCperMD()

const unsigned long int & PAPRECA::PaprecaConfig::getKMCperMD ( ) const

◆ getKMCsteps()

const unsigned long int & PAPRECA::PaprecaConfig::getKMCsteps ( ) const

◆ getLogFile()

Log & PAPRECA::PaprecaConfig::getLogFile ( )

◆ getMaxBondsFromSpecies()

int PAPRECA::PaprecaConfig::getMaxBondsFromSpecies ( const int &  atom_type)

Returns the maximum number of bonds of a specific atom type as defined in the PAPRECA::PredefinedEventsCatalog.

Parameters
[in]atom_typetype of atom.
Returns
maximum number of bonds for a specific atom type.

◆ getMaxBondTypesOfSpecies()

int PAPRECA::PaprecaConfig::getMaxBondTypesOfSpecies ( const int &  atom_type,
const int &  bond_type 
)

Returns the maximum number of bonds of a specific bond type, for a specific atom type.

Parameters
[in]atom_typetype of atom.
[in]bond_typetype of bond to be checked.
Returns
maximum number of bonds of a specific bond type for an atom type.

◆ getMinimize1()

const std::string & PAPRECA::PaprecaConfig::getMinimize1 ( ) const

◆ getMinimize2()

const std::string & PAPRECA::PaprecaConfig::getMinimize2 ( ) const

◆ getMonoatomicDesorptionFromAtomType()

PredefinedMonoatomicDesorption * PAPRECA::PaprecaConfig::getMonoatomicDesorptionFromAtomType ( const int &  atom_type)

Returns a PAPRECA::PredefinedMonoatomicDesorption object as defined in the PAPRECA::PredefinedEventsCatalog.

Parameters
[in]atom_typetype of atom.
Returns
mapped PAPRECA::PredefinedMonoatomicDesorption as defined in the PAPRECA::PredefinedEventsCatalog.
Note
each atom type can only be mapped to one PAPRECA::PredefinedMonoatomicDesorption.

◆ getRandomDiffVecsStyle()

const std::string & PAPRECA::PaprecaConfig::getRandomDiffVecsStyle ( ) const

◆ getReactionFromBondType()

PredefinedReaction * PAPRECA::PaprecaConfig::getReactionFromBondType ( const int &  bond_type)

Returns a PAPRECA::PredefinedReaction object as defined in the PAPRECA:PredefinedEventsCatalog. param[in] bond_type type of bond.

Returns
mapped PAPRECA::PredefinedReaction object.
Note
each bond can only be mapped to one PAPRECA::PredefinedReaction event.

◆ getRestartDumpFreq()

const int & PAPRECA::PaprecaConfig::getRestartDumpFreq ( ) const

◆ getSigmaFromAtomTypes()

const double PAPRECA::PaprecaConfig::getSigmaFromAtomTypes ( const int &  atom1_type,
const int &  atom2_type 
)

◆ getSigmaMixStyle()

const std::string & PAPRECA::PaprecaConfig::getSigmaMixStyle ( ) const

◆ getSigmaStyle()

const std::string & PAPRECA::PaprecaConfig::getSigmaStyle ( ) const

◆ getSurfaceCoverageFile()

SurfaceCoverage & PAPRECA::PaprecaConfig::getSurfaceCoverageFile ( )

◆ getTimeEnd()

const double & PAPRECA::PaprecaConfig::getTimeEnd ( )

◆ getTrajDuration()

const int & PAPRECA::PaprecaConfig::getTrajDuration ( ) const

◆ getUniformRanNum()

const double PAPRECA::PaprecaConfig::getUniformRanNum ( ) const

◆ initPredefinedBondForm()

void PAPRECA::PaprecaConfig::initPredefinedBondForm ( const int &  atom1_type,
const int &  atom2_type,
const int &  bond_type,
const double &  bond_dist,
const int &  delete_atoms,
const int &  lone_candidates,
const bool &  same_mol,
const double &  rate,
const std::vector< int > &  catalyzing_types 
)

◆ initPredefinedDeposition()

void PAPRECA::PaprecaConfig::initPredefinedDeposition ( LAMMPS_NS::LAMMPS *  lmp,
const int &  parent_type,
const double &  depo_offset,
const double &  insertion_vel,
const std::string &  adsorbate_name,
const double &  rate,
const bool &  variable_sticking,
const double &  sticking_coeff 
)

◆ initPredefinedDiffusionHop()

void PAPRECA::PaprecaConfig::initPredefinedDiffusionHop ( const int &  parent_type,
const double &  insertion_vel,
const double &  diff_dist,
const bool &  is_displacive,
const int &  diffused_type,
const double &  rate,
const std::string &  custom_style,
const std::vector< int > &  custom_atomtypes 
)

◆ initPredefinedMonoatomicDesorption()

void PAPRECA::PaprecaConfig::initPredefinedMonoatomicDesorption ( const int &  parent_type,
const double &  rate 
)

◆ initPredefinedReaction()

void PAPRECA::PaprecaConfig::initPredefinedReaction ( const int &  atom1_type,
const int &  atom2_type,
const int &  bond_type,
const double &  rate,
const std::vector< int > &  catalyzing_types 
)

◆ initRanNumGenerator()

void PAPRECA::PaprecaConfig::initRanNumGenerator ( LAMMPS_NS::LAMMPS *  lmp,
const int &  seed 
)

◆ initSigmasFromLammps()

void PAPRECA::PaprecaConfig::initSigmasFromLammps ( LAMMPS_NS::LAMMPS *  lmp)

◆ mixSigmas()

void PAPRECA::PaprecaConfig::mixSigmas ( LAMMPS_NS::LAMMPS *  lmp)

◆ predefinedCatalogHasBondBreakEvents()

const bool PAPRECA::PaprecaConfig::predefinedCatalogHasBondBreakEvents ( ) const

◆ predefinedCatalogHasBondFormEvents()

const bool PAPRECA::PaprecaConfig::predefinedCatalogHasBondFormEvents ( ) const

◆ predefinedCatalogHasDepositionEvents()

const bool PAPRECA::PaprecaConfig::predefinedCatalogHasDepositionEvents ( ) const

◆ predefinedCatalogHasDiffusionHopEvents()

const bool PAPRECA::PaprecaConfig::predefinedCatalogHasDiffusionHopEvents ( ) const

◆ predefinedCatalogHasMonoDesEvents()

const bool PAPRECA::PaprecaConfig::predefinedCatalogHasMonoDesEvents ( ) const

◆ predefinedCatalogIsEmpty()

const bool PAPRECA::PaprecaConfig::predefinedCatalogIsEmpty ( ) const

◆ ranNumGeneratorIsInitialized()

const bool PAPRECA::PaprecaConfig::ranNumGeneratorIsInitialized ( ) const

◆ setBinWidth()

void PAPRECA::PaprecaConfig::setBinWidth ( const double &  bin_width_in)

◆ setCtimeConvert()

void PAPRECA::PaprecaConfig::setCtimeConvert ( const double &  c_convert_in)

◆ setDepoHeights()

void PAPRECA::PaprecaConfig::setDepoHeights ( const double &  height_deposcan_in,
const double &  height_deporeject_in 
)

◆ setDesorbDelMax()

void PAPRECA::PaprecaConfig::setDesorbDelMax ( const int &  desorb_delmax_in)

◆ setDesorptionHeight()

void PAPRECA::PaprecaConfig::setDesorptionHeight ( const double &  desorb_cut_in)

◆ setDesorptionStyle()

void PAPRECA::PaprecaConfig::setDesorptionStyle ( const std::string &  desorb_style_in)

◆ setFluidAtomTypes()

void PAPRECA::PaprecaConfig::setFluidAtomTypes ( const std::vector< int > &  fluid_atomtypes_in)

◆ setFrozenAtomTypes()

void PAPRECA::PaprecaConfig::setFrozenAtomTypes ( const std::vector< int > &  frozen_atomtypes_in)

◆ setHeightMethod()

void PAPRECA::PaprecaConfig::setHeightMethod ( const std::string &  height_method_in)

◆ setHeightPercentage()

void PAPRECA::PaprecaConfig::setHeightPercentage ( const double &  height_percentage_in)

◆ setHybridStartTimeStamp4ExecTimeFile()

void PAPRECA::PaprecaConfig::setHybridStartTimeStamp4ExecTimeFile ( const int &  KMC_loopid)

◆ setKMCperMD()

void PAPRECA::PaprecaConfig::setKMCperMD ( const unsigned long int &  KMC_per_MD_in)

◆ setKMCsteps()

void PAPRECA::PaprecaConfig::setKMCsteps ( const unsigned long int &  KMC_steps_in)

◆ setMDTimeStamp4ExecTimeFile()

void PAPRECA::PaprecaConfig::setMDTimeStamp4ExecTimeFile ( const int &  KMC_loopid)

◆ setMinimize1()

void PAPRECA::PaprecaConfig::setMinimize1 ( const std::string &  minimize1_in)

◆ setMinimize2()

void PAPRECA::PaprecaConfig::setMinimize2 ( const std::string &  minimize2_in)

◆ setNeibLists()

void PAPRECA::PaprecaConfig::setNeibLists ( const std::string &  neiblist_half_in,
const std::string &  neiblist_full_in 
)

◆ setRandomDepoVecs()

void PAPRECA::PaprecaConfig::setRandomDepoVecs ( const bool &  random_depovecs_in)

◆ setRandomDiffVecs()

void PAPRECA::PaprecaConfig::setRandomDiffVecs ( const bool &  random_diffvecs_in)

◆ setRandomDiffVecsStyle()

void PAPRECA::PaprecaConfig::setRandomDiffVecsStyle ( const std::string &  diffvecs_style_in)

◆ setRestartDumpFreq()

void PAPRECA::PaprecaConfig::setRestartDumpFreq ( const int &  restart_dumpfreq_in)

◆ setSigmaMix()

void PAPRECA::PaprecaConfig::setSigmaMix ( const std::string &  mixstyle_in)

◆ setSigmaStyle()

void PAPRECA::PaprecaConfig::setSigmaStyle ( const std::string &  sigmastyle_in)

◆ setSpeciesMaxBonds()

void PAPRECA::PaprecaConfig::setSpeciesMaxBonds ( const int &  species,
const int &  bonds_max 
)

◆ setSpeciesMaxBondTypes()

void PAPRECA::PaprecaConfig::setSpeciesMaxBondTypes ( const int &  species,
const int &  bond_type,
const int &  bonds_max 
)

◆ setSpeciesPair2Sigma()

void PAPRECA::PaprecaConfig::setSpeciesPair2Sigma ( const int &  species1,
const int &  species2,
const double &  sigma 
)

◆ setTimeEnd()

void PAPRECA::PaprecaConfig::setTimeEnd ( const double &  time_end_in)

◆ setTrajDuration()

void PAPRECA::PaprecaConfig::setTrajDuration ( const int &  traj_duration_in)

◆ setupExportFiles()

void PAPRECA::PaprecaConfig::setupExportFiles ( const int &  proc_id)

◆ type2SigmaMapIsEmpty()

const bool PAPRECA::PaprecaConfig::type2SigmaMapIsEmpty ( ) const

Member Data Documentation

◆ bin_width

double PAPRECA::PaprecaConfig::bin_width = 1.0
protected

Bin width for height calculation or for ElementaDistribution files. Initialized at 1.0, so even if the user does not use the height_caulcation or heightbin_width command, you would be able to calculate mass bins with a width of 1.0.

◆ c_time_convert

double PAPRECA::PaprecaConfig::c_time_convert = -1
protected

This variable is initialized from lmp->update->unit_style to allow conversion of time units to second. The constant is also pre-multiplied with the LAMMPS timestep to allow easy conversion of timesteps (i.e., traj_duration) to time interval (in seconds).

◆ desorb_cut

double PAPRECA::PaprecaConfig::desorb_cut = -1
protected

Atoms above film_height + desorb_cut are deleted. The default value is -1 which means that desorption is disabled.

◆ desorb_delmax

int PAPRECA::PaprecaConfig::desorb_delmax = std::numeric_limits< int >::max( )
protected

Maximum number of atoms that can be deleted at once. Initialized at max limits of int so if the user does not set that, the maximum number of deleted atoms will be unlimited.

◆ desorb_style

std::string PAPRECA::PaprecaConfig::desorb_style = ""
protected

User defined desorption algorithm. Currently, can be gather_local or gather_all (defined in the PAPRECA input file).

◆ diffvecs_style

std::string PAPRECA::PaprecaConfig::diffvecs_style = "3D"
protected

Type 2D means that the random diffusion vector is always ABOVE the parent type. Type 3D means that the random diffusion vector can be anywhere in 3D space. Default is 3D as it is the more general case. Note that: we do not use 2D/3D random deposition vectors as it does not really make sense to get a random vector below the parent atom type.

◆ elementalDistribution_files

ElementalDistribution PAPRECA::PaprecaConfig::elementalDistribution_files
protected

stores the PAPRECA::ElementalDistribution files generated in the simulation.

◆ execTime_file

ExecTime PAPRECA::PaprecaConfig::execTime_file
protected

stores a PAPRECA::ExecTime file.

◆ fluid_atomtypes

std::vector< int > PAPRECA::PaprecaConfig::fluid_atomtypes
protected

stores fluid atom types (as defined in the PAPRECA input file).

◆ frozen_atomtypes

std::vector< int > PAPRECA::PaprecaConfig::frozen_atomtypes
protected

stores frozen atom types (as defined in the PAPRECA input file).

◆ height_deporeject

double PAPRECA::PaprecaConfig::height_deporeject = -1
protected

Reject deposition event above height_current + height_deporeject. Default at -1 which means do not reject anything.

◆ height_deposcan

double PAPRECA::PaprecaConfig::height_deposcan = -1
protected

Scan for deposition events only +- above/below the current film height. Default at -1 which means scan everywhere.

◆ height_method

std::string PAPRECA::PaprecaConfig::height_method = ""
protected

Algorithm to calculate height. Currently, only the mass_bins method is supported.

◆ height_percentage

double PAPRECA::PaprecaConfig::height_percentage = 0.0
protected

Only used for height calculation method: mass_bins. Defines the mass percentage cutoff. The film height is calculated as the first bin whose cumulative mass is above the percentage cutoff.

◆ heightVtime_file

HeightVtime PAPRECA::PaprecaConfig::heightVtime_file
protected

stores a PAPRECA::HeightVtime file.

◆ KMC_per_MD

unsigned long int PAPRECA::PaprecaConfig::KMC_per_MD = 0
protected

Perform that many KMC(PAPRECA) steps for every MD(LAMMPS) step.

◆ KMC_steps

unsigned long int PAPRECA::PaprecaConfig::KMC_steps = 0
protected

perform that many PAPRECA steps.

◆ log_file

Log PAPRECA::PaprecaConfig::log_file
protected

stores a PAPRECA::Log file The log file is always active (no need to be set to active from the input file).

◆ minimize1

std::string PAPRECA::PaprecaConfig::minimize1 = ""
protected

stores a valid LAMMPS minimization command to be executed before the LAMMPS trajectory. See: https://docs.lammps.org/minimize.html.

◆ minimize2

std::string PAPRECA::PaprecaConfig::minimize2 = ""
protected

stores a valid LAMMPS minimization command to be executed after the LAMMPS trajectory. See: https://docs.lammps.org/minimize.html.

◆ neiblist_full

std::string PAPRECA::PaprecaConfig::neiblist_full = ""
protected

Name of LAMMPS full neighbors list.

◆ neiblist_half

std::string PAPRECA::PaprecaConfig::neiblist_half = ""
protected

Name of LAMMPS half neighbors list.

◆ predefined_catalog

PredefinedEventsCatalog PAPRECA::PaprecaConfig::predefined_catalog
protected

◆ random_depovecs

bool PAPRECA::PaprecaConfig::random_depovecs = false
protected

Controls deposition sites. If true, the deposition sites are not directly above the parent atom, but on the surface of a sphere of radius depo_offset.

◆ random_diffvecs

bool PAPRECA::PaprecaConfig::random_diffvecs = false
protected

Controls diffusion sites. If true, the diffusion sites are not directly above the parent atom, but on the surface of a sphere centered at the parent atom.

◆ restart_dumpfreq

int PAPRECA::PaprecaConfig::restart_dumpfreq = std::numeric_limits< int >::max( )
protected

dump a restart every restart_dumpfreq PAPRECA steps. Initialized at int limits, so if it is not set you virtually never dump restarts (see how restarts are dumped in lammps_wrappers.h of papreca lib).

◆ rnum_gen

LAMMPS_NS::RanMars* PAPRECA::PaprecaConfig::rnum_gen = NULL
protected

Points to RanMars object as implemented in LAMMPS (see random_mars.h). The LAMMPS RanMars object is initialized from a user defined random seed in the PAPRECA input file (seed > 0 && seed < 900000000).

◆ sigma_mix

std::string PAPRECA::PaprecaConfig::sigma_mix =""
protected

◆ sigma_style

std::string PAPRECA::PaprecaConfig::sigma_style = ""
protected

method for initialization of sigma values (can be manual/LAMMPS).

◆ surface_coverage

double PAPRECA::PaprecaConfig::surface_coverage = 0.0
protected

stores a surface coverage for easier printing (if necessary).

◆ surfcoverage_file

SurfaceCoverage PAPRECA::PaprecaConfig::surfcoverage_file
protected

stores a PAPRECA::SurfaceCoverage file.

◆ time_end

double PAPRECA::PaprecaConfig::time_end = std::numeric_limits< double >::max( )
protected

Optional parameter. The user can define an ending simulation time. If not set it stays at limits of max, so you effectively never go beyond time_end;.

◆ traj_duration

int PAPRECA::PaprecaConfig::traj_duration = 0
protected

This control the duration of the LAMMPS minimization. The specific LAMMPS fixes have to be defined in the LAMMPS input file.

◆ type2sigma

INTPAIR2DOUBLE_MAP PAPRECA::PaprecaConfig::type2sigma
protected

maps types of atom types to their corresponding sigma.


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