PhreeqcRM
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes
PhreeqcRM C++ Class Reference

Public Member Functions

 PhreeqcRM (int nxyz, int thread_count_or_communicator, PHRQ_io *io=NULL)
 
IRM_RESULT CloseFiles (void)
 
IPhreeqc * Concentrations2Utility (std::vector< double > &c, std::vector< double > tc, std::vector< double > p_atm)
 
IRM_RESULT CreateMapping (std::vector< int > &grid2chem)
 
void DecodeError (int result)
 
IRM_RESULT DumpModule (bool dump_on, bool append=false)
 
void ErrorHandler (int result, const std::string &e_string)
 
void ErrorMessage (const std::string &error_string, bool prepend=true)
 
int FindComponents ()
 
const std::vector< std::vector
< int > > & 
GetBackwardMapping (void)
 
int GetChemistryCellCount (void) const
 
int GetComponentCount (void) const
 
const std::vector< std::string > & GetComponents (void) const
 
IRM_RESULT GetConcentrations (std::vector< double > &c)
 
std::string GetDatabaseFileName (void)
 
IRM_RESULT GetDensity (std::vector< double > &density)
 
const std::vector< int > & GetEndCell (void) const
 
const std::vector< std::string > & GetEquilibriumPhases (void) const
 
int GetEquilibriumPhasesCount (void) const
 
int GetErrorHandlerMode (void)
 
std::string GetErrorString (void)
 
const std::vector< std::string > & GetExchangeNames (void) const
 
const std::vector< std::string > & GetExchangeSpecies (void) const
 
int GetExchangeSpeciesCount (void) const
 
std::string GetFilePrefix (void)
 
const std::vector< int > & GetForwardMapping (void)
 
const std::vector< std::string > & GetGasComponents (void) const
 
int GetGasComponentsCount (void) const
 
IRM_RESULT GetGasCompMoles (std::vector< double > &gas_moles)
 
IRM_RESULT GetGasCompPressures (std::vector< double > &gas_pressure)
 
IRM_RESULT GetGasCompPhi (std::vector< double > &gas_phi)
 
IRM_RESULT GetGasPhaseVolume (std::vector< double > &gas_volume)
 
const std::vector< double > & GetGfw (void)
 
int GetGridCellCount (void)
 
IPhreeqc * GetIPhreeqcPointer (int i)
 
const std::vector< std::string > & GetKineticReactions (void) const
 
int GetKineticReactionsCount (void) const
 
int GetMpiMyself (void) const
 
int GetMpiTasks (void) const
 
int GetNthSelectedOutputUserNumber (int n)
 
bool GetPartitionUZSolids (void) const
 
const std::vector< double > & GetPressure (void)
 
const std::vector< int > & GetPrintChemistryMask (void)
 
const std::vector< bool > & GetPrintChemistryOn (void) const
 
bool GetRebalanceByCell (void) const
 
double GetRebalanceFraction (void) const
 
IRM_RESULT GetSaturation (std::vector< double > &sat)
 
IRM_RESULT GetSelectedOutput (std::vector< double > &so)
 
int GetSelectedOutputColumnCount (void)
 
int GetSelectedOutputCount (void)
 
IRM_RESULT GetSelectedOutputHeading (int icol, std::string &heading)
 
bool GetSelectedOutputOn (void)
 
int GetSelectedOutputRowCount (void)
 
int GetSICount (void) const
 
const std::vector< std::string > & GetSINames (void) const
 
const std::vector< std::string > & GetSolidSolutionComponents (void) const
 
int GetSolidSolutionComponentsCount (void) const
 
const std::vector< std::string > & GetSolidSolutionNames (void) const
 
const std::vector< double > & GetSolutionVolume (void)
 
IRM_RESULT GetSpeciesConcentrations (std::vector< double > &species_conc)
 
int GetSpeciesCount (void)
 
const std::vector< double > & GetSpeciesD25 (void)
 
IRM_RESULT GetSpeciesLog10Gammas (std::vector< double > &species_log10gammas)
 
IRM_RESULT GetSpeciesLog10Molalities (std::vector< double > &species_log10molalities)
 
const std::vector< std::string > & GetSpeciesNames (void)
 
bool GetSpeciesSaveOn (void)
 
const std::vector
< cxxNameDouble > & 
GetSpeciesStoichiometry (void)
 
const std::vector< double > & GetSpeciesZ (void)
 
const std::vector< int > & GetStartCell (void) const
 
const std::vector< std::string > & GetSurfaceNames (void) const
 
const std::vector< std::string > & GetSurfaceSpecies (void) const
 
int GetSurfaceSpeciesCount (void) const
 
const std::vector< std::string > & GetSurfaceTypes (void) const
 
const std::vector< double > & GetTemperature (void)
 
int GetThreadCount ()
 
double GetTime (void) const
 
double GetTimeConversion (void)
 
double GetTimeStep (void)
 
int GetUnitsExchange (void)
 
int GetUnitsGasPhase (void)
 
int GetUnitsKinetics (void)
 
int GetUnitsPPassemblage (void)
 
int GetUnitsSolution (void)
 
int GetUnitsSSassemblage (void)
 
int GetUnitsSurface (void)
 
const std::vector
< IPhreeqcPhast * > & 
GetWorkers ()
 
IRM_RESULT InitialPhreeqc2Concentrations (std::vector< double > &destination_c, std::vector< int > &boundary_solution1)
 
IRM_RESULT InitialPhreeqc2Concentrations (std::vector< double > &destination_c, std::vector< int > &boundary_solution1, std::vector< int > &boundary_solution2, std::vector< double > &fraction1)
 
IRM_RESULT InitialPhreeqc2Module (std::vector< int > &initial_conditions1)
 
IRM_RESULT InitialPhreeqc2Module (std::vector< int > &initial_conditions1, std::vector< int > &initial_conditions2, std::vector< double > &fraction1)
 
IRM_RESULT InitialPhreeqc2SpeciesConcentrations (std::vector< double > &destination_c, std::vector< int > &boundary_solution1)
 
IRM_RESULT InitialPhreeqc2SpeciesConcentrations (std::vector< double > &destination_c, std::vector< int > &boundary_solution1, std::vector< int > &boundary_solution2, std::vector< double > &fraction1)
 
IRM_RESULT InitialPhreeqcCell2Module (int n, const std::vector< int > &cell_numbers)
 
IRM_RESULT LoadDatabase (const std::string &database)
 
void LogMessage (const std::string &str)
 
int MpiAbort ()
 
IRM_RESULT MpiWorker ()
 
IRM_RESULT MpiWorkerBreak ()
 
IRM_RESULT OpenFiles (void)
 
void OutputMessage (const std::string &str)
 
IRM_RESULT RunCells (void)
 
IRM_RESULT ReturnHandler (IRM_RESULT result, const std::string &e_string)
 
IRM_RESULT RunFile (bool workers, bool initial_phreeqc, bool utility, const std::string &chemistry_name)
 
IRM_RESULT RunString (bool workers, bool initial_phreeqc, bool utility, const std::string &input_string)
 
void ScreenMessage (const std::string &str)
 
IRM_RESULT SetComponentH2O (bool tf)
 
IRM_RESULT SetConcentrations (const std::vector< double > &c)
 
IRM_RESULT SetCurrentSelectedOutputUserNumber (int n_user)
 
IRM_RESULT SetDensity (const std::vector< double > &density)
 
IRM_RESULT SetDumpFileName (const std::string &dump_name)
 
IRM_RESULT SetErrorHandlerMode (int mode)
 
IRM_RESULT SetErrorOn (bool tf)
 
IRM_RESULT SetFilePrefix (const std::string &prefix)
 
IRM_RESULT SetGasCompMoles (const std::vector< double > &gas_moles)
 
IRM_RESULT SetGasPhaseVolume (const std::vector< double > &gas_volume)
 
IRM_RESULT SetMpiWorkerCallbackC (int(*fcn)(int *method, void *cookie))
 
IRM_RESULT SetMpiWorkerCallbackCookie (void *cookie)
 
IRM_RESULT SetMpiWorkerCallbackFortran (int(*fcn)(int *method))
 
IRM_RESULT SetPartitionUZSolids (bool tf)
 
IRM_RESULT SetPorosity (const std::vector< double > &por)
 
IRM_RESULT SetPressure (const std::vector< double > &p)
 
IRM_RESULT SetPrintChemistryMask (std::vector< int > &cell_mask)
 
IRM_RESULT SetPrintChemistryOn (bool workers, bool initial_phreeqc, bool utility)
 
IRM_RESULT SetRebalanceByCell (bool tf)
 
IRM_RESULT SetRebalanceFraction (double f)
 
IRM_RESULT SetRepresentativeVolume (const std::vector< double > &rv)
 
IRM_RESULT SetSaturation (const std::vector< double > &sat)
 
IRM_RESULT SetScreenOn (bool tf)
 
IRM_RESULT SetSelectedOutputOn (bool tf)
 
IRM_RESULT SetSpeciesSaveOn (bool save_on)
 
IRM_RESULT SetTemperature (const std::vector< double > &t)
 
IRM_RESULT SetTime (double time)
 
IRM_RESULT SetTimeConversion (double conv_factor)
 
IRM_RESULT SetTimeStep (double time_step)
 
IRM_RESULT SetUnitsExchange (int option)
 
IRM_RESULT SetUnitsGasPhase (int option)
 
IRM_RESULT SetUnitsKinetics (int option)
 
IRM_RESULT SetUnitsPPassemblage (int option)
 
IRM_RESULT SetUnitsSolution (int option)
 
IRM_RESULT SetUnitsSSassemblage (int option)
 
IRM_RESULT SetUnitsSurface (int option)
 
IRM_RESULT SpeciesConcentrations2Module (std::vector< double > &species_conc)
 
IRM_RESULT StateSave (int istate)
 
IRM_RESULT StateApply (int istate)
 
IRM_RESULT StateDelete (int istate)
 
void UseSolutionDensityVolume (bool tf)
 
void WarningMessage (const std::string &warnstr)
 

Static Public Member Functions

static void CleanupReactionModuleInstances (void)
 
static int CreateReactionModule (int nxyz, int nthreads)
 
static IRM_RESULT DestroyReactionModule (int n)
 
static PhreeqcRMGetInstance (int n)
 
static std::string Char2TrimString (const char *str, size_t l=0)
 
static bool FileExists (const std::string &name)
 
static void FileRename (const std::string &temp_name, const std::string &name, const std::string &backup_name)
 
static IRM_RESULT Int2IrmResult (int r, bool positive_ok)
 

Protected Member Functions

IRM_RESULT CellInitialize (int i, int n_user_new, int *initial_conditions1, int *initial_conditions2, double *fraction1, std::set< std::string > &error_set)
 
IRM_RESULT CheckCells ()
 
int CheckSelectedOutput ()
 
IPhreeqc * Concentrations2UtilityH2O (std::vector< double > &c_in, std::vector< double > &t_in, std::vector< double > &p_in)
 
IPhreeqc * Concentrations2UtilityNoH2O (std::vector< double > &c_in, std::vector< double > &t_in, std::vector< double > &p_in)
 
void Concentrations2Solutions (int n, std::vector< double > &c)
 
void Concentrations2SolutionsH2O (int n, std::vector< double > &c)
 
void Concentrations2SolutionsNoH2O (int n, std::vector< double > &c)
 
void cxxSolution2concentration (cxxSolution *cxxsoln_ptr, std::vector< double > &d, double v, double dens)
 
void cxxSolution2concentrationH2O (cxxSolution *cxxsoln_ptr, std::vector< double > &d, double v, double dens)
 
void cxxSolution2concentrationNoH2O (cxxSolution *cxxsoln_ptr, std::vector< double > &d, double v, double dens)
 
void GatherNchem (std::vector< double > &source, std::vector< double > &destination)
 
cxxStorageBin & Get_phreeqc_bin (void)
 
IRM_RESULT HandleErrorsInternal (std::vector< int > &r)
 
void PartitionUZ (int n, int iphrq, int ihst, double new_frac)
 
void RebalanceLoad (void)
 
void RebalanceLoadPerCell (void)
 
IRM_RESULT RunCellsThread (int i)
 
IRM_RESULT RunFileThread (int n)
 
IRM_RESULT RunStringThread (int n, std::string &input)
 
IRM_RESULT RunCellsThreadNoPrint (int n)
 
void Scale_solids (int n, int iphrq, double frac)
 
void ScatterNchem (double *d_array)
 
void ScatterNchem (int *i_array)
 
void ScatterNchem (std::vector< double > &source, std::vector< double > &destination)
 
void ScatterNchem (std::vector< int > &source, std::vector< int > &destination)
 
IRM_RESULT SetChemistryFileName (const char *prefix=NULL)
 
IRM_RESULT SetDatabaseFileName (const char *db=NULL)
 
void SetEndCells (void)
 
void SetEndCellsHeterogeneous (void)
 
double TimeStandardTask (void)
 
IRM_RESULT TransferCells (cxxStorageBin &t_bin, int old, int nnew)
 
IRM_RESULT TransferCellsUZ (std::ostringstream &raw_stream, int old, int nnew)
 

Protected Attributes

bool component_h2o
 
std::string database_file_name
 
std::string chemistry_file_name
 
std::string dump_file_name
 
std::string file_prefix
 
cxxStorageBin * phreeqc_bin
 
int mpi_myself
 
int mpi_tasks
 
std::vector< std::string > components
 
std::vector< double > gfw
 
double gfw_water
 
bool partition_uz_solids
 
int nxyz
 
int count_chemistry
 
double time
 
double time_step
 
double time_conversion
 
std::vector< double > old_saturation_root
 
std::vector< double > old_saturation_worker
 
std::vector< double > saturation_root
 
std::vector< double > saturation_worker
 
std::vector< double > pressure_root
 
std::vector< double > pressure_worker
 
std::vector< double > rv_root
 
std::vector< double > rv_worker
 
std::vector< double > porosity_root
 
std::vector< double > porosity_worker
 
std::vector< double > tempc_root
 
std::vector< double > tempc_worker
 
std::vector< double > density_root
 
std::vector< double > density_worker
 
std::vector< double > solution_volume_root
 
std::vector< double > solution_volume_worker
 
std::vector< int > print_chem_mask_root
 
std::vector< int > print_chem_mask_worker
 
bool rebalance_by_cell
 
double rebalance_fraction
 
int units_Solution
 
int units_PPassemblage
 
int units_Exchange
 
int units_Surface
 
int units_GasPhase
 
int units_SSassemblage
 
int units_Kinetics
 
std::vector< int > forward_mapping_root
 
std::vector< std::vector< int > > backward_mapping
 
bool use_solution_density_volume
 
std::vector< bool > print_chemistry_on
 
bool selected_output_on
 
int error_count
 
int error_handler_mode
 
bool need_error_check
 
std::string phreeqcrm_error_string
 
int nthreads
 
std::vector< IPhreeqcPhast * > workers
 
std::vector< int > start_cell
 
std::vector< int > end_cell
 
PHRQ_io * phreeqcrm_io
 
bool delete_phreeqcrm_io
 
int(* mpi_worker_callback_fortran )(int *method)
 
int(* mpi_worker_callback_c )(int *method, void *cookie)
 
void * mpi_worker_callback_cookie
 
bool species_save_on
 
std::vector< std::string > species_names
 
std::vector< double > species_z
 
std::vector< double > species_d_25
 
std::vector< cxxNameDouble > species_stoichiometry
 
std::map< int, int > s_num2rm_species_num
 
std::vector< double > standard_task_vector
 
std::vector< std::string > ExchangeSpeciesNamesList
 
std::vector< std::string > ExchangeNamesList
 
std::vector< std::string > SurfaceSpeciesNamesList
 
std::vector< std::string > SurfaceTypesList
 
std::vector< std::string > SurfaceNamesList
 
std::vector< std::string > EquilibriumPhasesList
 
std::vector< std::string > GasComponentsList
 
std::vector< std::string > KineticReactionsList
 
std::vector< std::string > SolidSolutionComponentsList
 
std::vector< std::string > SolidSolutionNamesList
 
std::vector< std::string > SINamesList
 

Detailed Description

Geochemical reaction module.

Constructor & Destructor Documentation

PhreeqcRM ( int  nxyz,
int  thread_count_or_communicator,
PHRQ_io *  io = NULL 
)

Constructor for the PhreeqcRM reaction module. If the code is compiled with the preprocessor directive USE_OPENMP, the reaction module use OPENMP and multiple threads. If the code is compiled with the preprocessor directive USE_MPI, the reaction module will use MPI and multiple processes. If neither preprocessor directive is used, the reaction module will be serial (unparallelized).

Parameters
nxyzThe number of grid cells in the users model.
thread_count_or_communicatorIf multithreaded, the number of threads to use in parallel segments of the code. If thread_count_or_communicator is <= 0, the number of threads is set equal to the number of processors in the computer. If multiprocessor, the MPI communicator to use within the reaction module.
ioOptionally, a PHRQ_io input/output object can be provided to the constructor. By default a PHRQ_io object is constructed to handle reading and writing files.
C++ Example:
int nxyz = 40;
#ifdef USE_MPI
  PhreeqcRM phreeqc_rm(nxyz, MPI_COMM_WORLD);
  int mpi_myself;
  if (MPI_Comm_rank(MPI_COMM_WORLD, &mpi_myself) != MPI_SUCCESS)
  {
    exit(4);
  }
  if (mpi_myself > 0)
  {
    phreeqc_rm.MpiWorker();
    return EXIT_SUCCESS;
  }
#else
  int nthreads = 3;
  PhreeqcRM phreeqc_rm(nxyz, nthreads);
#endif
MPI:
Called by root and all workers.

Member Function Documentation

IRM_RESULT CloseFiles ( void  )

Close the output and log files.

Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
OpenFiles, SetFilePrefix
C++ Example:
status = phreeqc_rm.CloseFiles();
MPI:
Called only by root.
IPhreeqc* Concentrations2Utility ( std::vector< double > &  c,
std::vector< double >  tc,
std::vector< double >  p_atm 
)

N sets of component concentrations are converted to SOLUTIONs numbered 1-n in the Utility IPhreeqc. The solutions can be reacted and manipulated with the methods of IPhreeqc. If solution concentration units (SetUnitsSolution) are per liter, one liter of solution is created in the Utility instance; if solution concentration units are mass fraction, one kilogram of solution is created in the Utility instance. The motivation for this method is the mixing of solutions in wells, where it may be necessary to calculate solution properties (pH for example) or react the mixture to form scale minerals. The code fragment below makes a mixture of concentrations and then calculates the pH of the mixture.

Parameters
cVector of concentrations to be made SOLUTIONs in Utility IPhreeqc. Vector contains n values for each component (GetComponentCount) in sequence.
tcVector of temperatures to apply to the SOLUTIONs, in degrees C. Vector of size n.
p_atmVector of pressures to apply to the SOLUTIONs, in atm. Vector of size n.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
C++ Example:
std::vector  c_well;
c_well.resize(1*ncomps, 0.0);
for (int i = 0; i < ncomps; i++)
{
  c_well[i] = 0.5 * c[0 + nxyz*i] + 0.5 * c[9 + nxyz*i];
}
std::vector tc, p_atm;
tc.resize(1, 15.0);
p_atm.resize(1, 3.0);
IPhreeqc * util_ptr = phreeqc_rm.Concentrations2Utility(c_well, tc, p_atm);
input = "SELECTED_OUTPUT 5; -pH;RUN_CELLS; -cells 1";
int iphreeqc_result;
util_ptr->SetOutputFileName("utility_cpp.txt");
util_ptr->SetOutputFileOn(true);
iphreeqc_result = util_ptr->RunString(input.c_str());
phreeqc_rm.ErrorHandler(iphreeqc_result, "IPhreeqc RunString failed");
int vtype;
double pH;
char svalue[100];
util_ptr->SetCurrentSelectedOutputUserNumber(5);
iphreeqc_result = util_ptr->GetSelectedOutputValue2(1, 0, &vtype, &pH, svalue, 100);
MPI:
Called only by root.
IRM_RESULT CreateMapping ( std::vector< int > &  grid2chem)

Provides a mapping from grid cells in the user's model to reaction cells for which chemistry needs to be run. The mapping is used to eliminate inactive cells and to use symmetry to decrease the number of cells for which chemistry must be run. The array grid2chem of size nxyz (the number of grid cells, GetGridCellCount) must contain the set of all integers 0 <= i < count_chemistry, where count_chemistry is a number less than or equal to nxyz. Inactive cells are assigned a negative integer. The mapping may be many-to-one to account for symmetry. Default is a one-to-one mapping–all user grid cells are reaction cells (equivalent to grid2chem values of 0,1,2,3,...,nxyz-1).

Parameters
grid2chemA vector of integers: Nonnegative is a reaction-cell number (0 based), negative is an inactive cell. Vector is of size nxyz (number of grid cells, GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
C++ Example:
// For demonstation, two equivalent rows by symmetry
std::vector grid2chem;
grid2chem.resize(nxyz, -1);
for (int i = 0; i < nxyz/2; i++)
{
  grid2chem[i] = i;
  grid2chem[i + nxyz/2] = i;
}
status = phreeqc_rm.CreateMapping(grid2chem);
if (status < 0) phreeqc_rm.DecodeError(status);
int nchem = phreeqc_rm.GetChemistryCellCount();
MPI:
Called by root, workers must be in the loop of MpiWorker.
void DecodeError ( int  result)

If result is negative, this method prints an error message corresponding to IRM_RESULT result. If result is non-negative, no action is taken.

Parameters
resultAn IRM_RESULT value returned by one of the reaction-module methods.
IRM_RESULT definition:
typedef enum {
  IRM_OK            =  0,  //Success
  IRM_OUTOFMEMORY   = -1,  //Failure, Out of memory
  IRM_BADVARTYPE    = -2,  //Failure, Invalid VAR type
  IRM_INVALIDARG    = -3,  //Failure, Invalid argument
  IRM_INVALIDROW    = -4,  //Failure, Invalid row
  IRM_INVALIDCOL    = -5,  //Failure, Invalid column
  IRM_BADINSTANCE   = -6,  //Failure, Invalid rm instance id
  IRM_FAIL          = -7,  //Failure, Unspecified
} IRM_RESULT;
C++ Example:
status = phreeqc_rm.CreateMapping(grid2chem);
phreeqc_rm.DecodeError(status);
MPI:
Can be called by root and (or) workers.
IRM_RESULT DumpModule ( bool  dump_on,
bool  append = false 
)

Writes the contents of all workers to file in _RAW formats (see appendix of PHREEQC version 3 manual), including SOLUTIONs and all reactants.

Parameters
dump_onSignal for writing the dump file, true or false.
appendSignal to append to the contents of the dump file, true or false.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetDumpFileName.
C++ Example:
bool dump_on = true;
bool append = false;
status = phreeqc_rm.SetDumpFileName("advection_cpp.dmp");
status = phreeqc_rm.DumpModule(dump_on, append);
MPI:
Called by root; workers must be in the loop of MpiWorker.
void ErrorHandler ( int  result,
const std::string &  e_string 
)

Checks result for an error code. If result is negative, the result is decoded (DecodeError), and printed as an error message along with the e_string, and an exception is thrown. If the result is nonnegative, no action is taken.

Parameters
resultIRM_RESULT to be checked for an error.
e_stringString to be printed if an error is found.
See also
DecodeError, ErrorMessage.
C++ Example:
iphreeqc_result = util_ptr->RunString(input.c_str());
if (iphreeqc_result != 0)
{
  phreeqc_rm.ErrorHandler(IRM_FAIL, "IPhreeqc RunString failed");
}
MPI:
Called by root and (or) workers.
void ErrorMessage ( const std::string &  error_string,
bool  prepend = true 
)

Send an error message to the screen, the output file, and the log file.

Parameters
error_stringString to be printed.
prependTrue, prepends error_string with "Error: "; false, error_string is used with no prepended text.
See also
OpenFiles, LogMessage, ScreenMessage, WarningMessage.
C++ Example:
phreeqc_rm.ErrorMessage("Goodbye world");
MPI:
Called by root and (or) workers; root writes to output and log files.
int FindComponents ( )

This method accumulates a list of elements. Elements are those that have been defined in a solution or any other reactant (EQUILIBRIUM_PHASE, KINETICS, and others), including charge imbalance. This method can be called multiple times and the list that is created is cummulative. The list is the set of components that needs to be transported. By default the list includes water, excess H and excess O (the H and O not contained in water); alternatively, the list may be set to contain total H and total O (SetComponentH2O), which requires transport results to be accurate to eight or nine significant digits. If multicomponent diffusion (MCD) is to be modeled, there is a capability to retrieve aqueous species concentrations (GetSpeciesConcentrations) and to set new solution concentrations after MCD by using individual species concentrations (SpeciesConcentrations2Module). To use these methods the save-species property needs to be turned on (SetSpeciesSaveOn). If the save-species property is on, FindComponents will generate a list of aqueous species (GetSpeciesCount, GetSpeciesNames), their diffusion coefficients at 25 C (GetSpeciesD25), and their charge (GetSpeciesZ).

Return values
Numberof components currently in the list, or IRM_RESULT error code (negative value, see DecodeError).
See also
GetComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SetComponentH2O, SetSpeciesSaveOn, SpeciesConcentrations2Module.
The FindComponents method also generates lists of reactants–equilibrium phases,
exchangers, gas components, kinetic reactants, solid solution components, and surfaces. The lists are cumulative, including all reactants that were defined in the initial phreeqc instance at any time FindComponents was called. In addition, a list of phases is generated for which saturation indices may be calculated from the cumulative list of components.
See also
also GetEquilibriumPhases, GetEquilibriumPhasesCount, GetExchangeNames, GetExchangeSpecies, GetExchangeSpeciesCount, GetGasComponents, GetGasComponentsCount, GetKineticReactions, GetKineticReactionsCount, GetSICount, GetSINames, GetSolidSolutionComponents, GetSolidSolutionComponentsCount, GetSolidSolutionNames, GetSurfaceNames, GetSurfaceSpecies, GetSurfaceSpeciesCount, GetSurfaceTypes.
C++ Example:
int ncomps = phreeqc_rm.FindComponents();
const std::vector &components = phreeqc_rm.GetComponents();
const std::vector < double > & gfw = phreeqc_rm.GetGfw();
for (int i = 0; i < ncomps; i++)
{
  std::ostringstream strm;
  strm.width(10);
  strm << components[i] << "    " << gfw[i] << "\n";
  phreeqc_rm.OutputMessage(strm.str());
}
MPI:
Called by root, workers must be in the loop of MpiWorker.
const std::vector< std::vector <int> >& GetBackwardMapping ( void  )
inline

Returns a vector of vectors, where the nth vector is a vector of grid-cell numbers that are mapped to reaction-cell number n. Each reaction-cell number has a vector of one or more grid-cell numbers.

Return values
Vectorof vectors of ints. For each reaction cell n, the nth vector in the vector of vectors contains the grid-cell numbers that map to the reaction cell.
See also
CreateMapping, GetForwardMapping, GetChemistryCellCount, GetGridCellCount.
C++ Example:
const std::vector < std::vector  > & back = phreeqcrm_ptr->GetBackwardMapping();
if (option == "HYDRAULIC_K")
{
  return (*data_ptr->hydraulic_K)[back[rm_cell_number][0]];
}
MPI:
Called by root or workers.
int GetChemistryCellCount ( void  ) const
inline

Returns the number of reaction cells in the reaction module. The number of reaction cells is defined by the set of non-negative integers in the mapping from grid cells (CreateMapping), or, by default, the number of grid cells (GetGridCellCount). The number of reaction cells is less than or equal to the number of grid cells in the user's model.

Return values
Numberof reaction cells.
See also
CreateMapping, GetGridCellCount.
C++ Example:
status = phreeqc_rm.CreateMapping(grid2chem);
std::ostringstream oss;
oss << "Number of reaction cells in the reaction module: "
    << phreeqc_rm.GetChemistryCellCount() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
int GetComponentCount ( void  ) const
inline

Returns the number of components in the reaction-module component list.

Return values
Thenumber of components in the reaction-module component list. The component list is generated by calls to FindComponents. The return value from the last call to FindComponents is equal to the return value from GetComponentCount.
See also
FindComponents, GetComponents.
C++ Example:
std::ostringstream oss;
oss << "Number of components for transport: " << phreeqc_rm.GetComponentCount() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root.
const std::vector<std::string>& GetComponents ( void  ) const
inline

Returns a reference to the reaction-module component list that was generated by calls to FindComponents.

Return values
conststd::vector<std::string>& A vector of strings; each string is a component name.
See also
FindComponents, GetComponentCount
C++ Example:
const std::vector &components = phreeqc_rm.GetComponents();
const std::vector < double > & gfw = phreeqc_rm.GetGfw();
int ncomps = phreeqc_rm.GetComponentCount();
for (int i = 0; i < ncomps; i++)
{
  std::ostringstream strm;
  strm.width(10);
  strm << components[i] << "    " << gfw[i] << "\n";
  phreeqc_rm.OutputMessage(strm.str());
}
MPI:
Called by root and (or) workers.
IRM_RESULT GetConcentrations ( std::vector< double > &  c)

Transfer solution concentrations from each reaction cell to the concentration vector given in the argument list (c). Units of concentration for c are defined by SetUnitsSolution. For per liter concentration units, solution volume is used to calculate the concentrations for c. For mass-fraction concentration units, the solution mass is used to calculate concentrations for c. Two options are available for the volume and mass of solution that are used in converting to transport concentrations: (1) the volume and mass of solution are calculated by PHREEQC, or (2) the volume of solution is the product of saturation (SetSaturation), porosity (SetPorosity), and representative volume (SetRepresentativeVolume), and the mass of solution is volume times density as defined by SetDensity. UseSolutionDensityVolume determines which option is used. For option 1, the databases that have partial molar volume definitions needed to accurately calculate solution volume are phreeqc.dat, Amm.dat, and pitzer.dat.

Parameters
cVector to receive the concentrations. Dimension of the vector is set to ncomps times nxyz, where, ncomps is the result of FindComponents or GetComponentCount, and nxyz is the number of user grid cells (GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetComponentCount, GetSaturation, SetConcentrations, SetDensity, SetRepresentativeVolume, SetSaturation, SetUnitsSolution, UseSolutionDensityVolume.
C++ Example:
std::vector c;
status = phreeqc_rm.RunCells();
status = phreeqc_rm.GetConcentrations(c);
MPI:
Called by root, workers must be in the loop of MpiWorker.
std::string GetDatabaseFileName ( void  )
inline

Returns the file name of the database. Should be called after LoadDatabase.

Return values
std::stringThe file name defined in LoadDatabase.
See also
LoadDatabase.
C++ Example:
std::ostringstream oss;
oss << "Database: " << phreeqc_rm.GetDatabaseFileName() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
IRM_RESULT GetDensity ( std::vector< double > &  density)

Transfer solution densities from the reaction-module workers to the vector given in the argument list (density).

Parameters
densityVector to receive the densities. Dimension of the array is set to nxyz, where nxyz is the number of user grid cells (GetGridCellCount). Values for inactive cells are set to 1e30. Densities are those calculated by the reaction module. Only the following databases distributed with PhreeqcRM have molar volume information needed to accurately calculate density: phreeqc.dat, Amm.dat, and pitzer.dat.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetSolutionVolume, SetDensity.
C++ Example:
status = phreeqc_rm.RunCells();
status = phreeqc_rm.GetConcentrations(c);
std::vector density;
status = phreeqc_rm.GetDensity(density);
MPI:
Called by root, workers must be in the loop of MpiWorker.
const std::vector< int>& GetEndCell ( void  ) const
inline

Returns a vector of integers that contains the largest reaction-cell number assigned to each worker. Each worker is assigned a range of reaction-cell numbers that are run during a call to RunCells. The range of reaction cells for a worker may vary as load rebalancing occurs. At any point in the calculations, the first cell and last cell to be run by a worker can be found in the vectors returned by GetStartCell and GetEndCell. Each method returns a vector of integers that has length of the number of threads (GetThreadCount), if using OPENMP, or the number of processes (GetMpiTasks), if using MPI.

Return values
IRM_RESULTVector of integers, one for each worker, that gives the last reaction cell to be run by each worker.
See also
GetStartCell, GetThreadCount, GetMpiTasks.
C++ Example:
std::ostringstream oss;
oss << "Current distribution of cells for workers\n";
oss << "Worker First Cell   Last Cell\n";
int n;
n = phreeqc_rm.GetThreadCount() * phreeqc_rm.GetMpiTasks();
for (int i = 0; i < n; i++)
{
    oss << i << "      "
        << phreeqc_rm.GetStartCell()[i]
        << "            "
        << phreeqc_rm.GetEndCell()[i] << "\n";
}
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
const std::vector<std::string>& GetEquilibriumPhases ( void  ) const
inline

Returns a reference to the vector of all equilibrium phases. The list includes all phases included in any EQUILIBRIUM_PHASES definitions in the initial-phreeqc module. FindComponents must be called before GetEquilibriumPhases. This method may be useful when generating selected output definitions related to equilibrium phases.

Return values
conststd::vector<std::string>& A vector of strings; each string is a unique equilibrium phases name.
See also
FindComponents, GetEquilibriumPhasesCount.
C++ Example:
oss << "  -equilibrium_phases " << "\n";
// equilibrium phases
const std::vector &eq_phases = phreeqc_rm.GetEquilibriumPhases();
for (size_t i = 0; i < phreeqc_rm.GetEquilibriumPhasesCount(); i++)
{
oss << "    " << eq_phases[i] << "\n";
}
MPI:
Called by root.
int GetEquilibriumPhasesCount ( void  ) const
inline

Returns the number of equilibrium phases in the initial-phreeqc module. FindComponents must be called before GetEquilibriumPhasesCount. This method may be useful when generating selected output definitions related to equilibrium phases.

Return values
Thenumber of equilibrium phases in the initial-phreeqc module.
See also
FindComponents, GetEquilibriumPhases.
C++ Example:
oss << "  -equilibrium_phases " << "\n";
// equilibrium phases
const std::vector &eq_phases = phreeqc_rm.GetEquilibriumPhases();
for (size_t i = 0; i < phreeqc_rm.GetEquilibriumPhasesCount(); i++)
{
oss << "    " << eq_phases[i] << "\n";
}
MPI:
Called by root.
int GetErrorHandlerMode ( void  )
inline

Get the setting for the action to be taken when the reaction module encounters an error. Options are 0, return to calling program with an error return code (default); 1, throw an exception, which can be caught in C++ (for C and Fortran, the program will exit); 2, attempt to exit gracefully.

Return values
IRM_RESULTCurrent setting for the error handling mode: 0, 1, or 2.
See also
SetErrorHandlerMode.
C++ Example:
std::ostringstream oss;
oss << "Error handler mode: " << phreeqc_rm.GetErrorHandlerMode() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
std::string GetErrorString ( void  )

Returns a standard string containing error messages related to the last call to a PhreeqcRM method.

Return values
Errormessages related to the last call to a PhreeqcRM method.
See also
GetErrorHandlerMode.
C++ Example:
if (status != IRM_OK)
{
  std::cerr << phreeqc_rm.GetErrorString();
}
MPI:
Called by root, workers must be in the loop of MpiWorker.
const std::vector<std::string>& GetExchangeNames ( void  ) const
inline

Returns a reference to the vector of exchange names (such as "X") that correspond with the exchange species names. FindComponents must be called before GetExchangeNames. The exchange names vector is the same length as the exchange species names vector and provides the corresponding exchange site. This method may be useful when generating selected output definitions related to exchangers.

Return values
conststd::vector<std::string>& A vector of strings; each string is an exchange name corresponding to the exchange species vector; an exchange name may occur multiple times.
See also
FindComponents, GetExchangeSpeciesCount, GetExchangeSpecies.
C++ Example:
// molalities of exchange species
const std::vector &ex_species = phreeqc_rm.GetExchangeSpecies();
const std::vector &ex_names = phreeqc_rm.GetExchangeNames();
for (size_t i = 0; i < phreeqc_rm.GetExchangeSpeciesCount(); i++)
{

oss << "    ";
oss.width(15);
oss << std::left << ex_species[i];
oss << " # " << ex_names[i] << "\n";
}
MPI:
Called by root.
const std::vector<std::string>& GetExchangeSpecies ( void  ) const
inline

Returns a reference to the vector of exchange species names (such as "NaX"). The list of exchange species (such as "NaX") is derived from the list of components (FindComponents) and the list of all exchange names (such as "X") that are included in EXCHANGE definitions in the initial-phreeqc module. FindComponents must be called before GetExchangeSpecies. This method may be useful when generating selected output definitions related to exchangers.

Return values
conststd::vector<std::string>& A vector of strings; each string is a unique exchange species name.
See also
FindComponents, GetExchangeSpeciesCount, GetExchangeNames.
C++ Example:
// molalities of exchange species
const std::vector &ex_species = phreeqc_rm.GetExchangeSpecies();
const std::vector &ex_names = phreeqc_rm.GetExchangeNames();
for (size_t i = 0; i < phreeqc_rm.GetExchangeSpeciesCount(); i++)
{

oss << "    ";
oss.width(15);
oss << std::left << ex_species[i];
oss << " # " << ex_names[i] << "\n";
}
MPI:
Called by root.
int GetExchangeSpeciesCount ( void  ) const
inline

Returns the number of exchange species in the initial-phreeqc module. FindComponents must be called before GetExchangeSpeciesCount. This method may be useful when generating selected output definitions related to exchangers.

Return values
Thenumber of exchange species in the initial-phreeqc module.
See also
FindComponents, GetExchangeSpecies, GetExchangeNames.
C++ Example:
// molalities of exchange species
const std::vector &ex_species = phreeqc_rm.GetExchangeSpecies();
const std::vector &ex_names = phreeqc_rm.GetExchangeNames();
for (size_t i = 0; i < phreeqc_rm.GetExchangeSpeciesCount(); i++)
{

oss << "    ";
oss.width(15);
oss << std::left << ex_species[i];
oss << " # " << ex_names[i] << "\n";
}
MPI:
Called by root.
std::string GetFilePrefix ( void  )
inline

Returns the file prefix for the output (.chem.txt) and log files (.log.txt).

Return values
std::stringThe file prefix as set by SetFilePrefix, or "myrun", by default.
See also
SetFilePrefix.
C++ Example:
std::ostringstream oss;
oss << "File prefix: " << phreeqc_rm.GetFilePrefix() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
const std::vector< int >& GetForwardMapping ( void  )
inline

Returns a reference to a vector of ints that is a mapping from grid cells to reaction cells. The mapping is used to eliminate cells that are inactive and cells that are unnecessary because of symmetry from the list of cells for which reactions must be run. The mapping may be many-to-one to account for symmetry. The mapping is set by CreateMapping, or, by default, is a one-to-one mapping–all grid cells are reaction cells (vector contains 0,1,2,3,...,nxyz-1).

Return values
conststd::vector < int >& A vector of integers of size nxyz (number of grid cells, GetGridCellCount). Nonnegative is a reaction-cell number (0 based), negative is an inactive cell.
C++ Example:
const std::vector &f_map = phreeqc_rm.GetForwardMapping();
MPI:
Called by root.
IRM_RESULT GetGasCompMoles ( std::vector< double > &  gas_moles)

Transfer moles of gas components from each reaction cell to the vector given in the argument list (gas_moles).

Parameters
gas_molesVector to receive the moles of gas components. Dimension of the vector is set to ngas_comps times nxyz, where, ngas_comps is the result of GetGasComponentsCount, and nxyz is the number of user grid cells (GetGridCellCount). If a gas component is not defined for a cell, the number of moles is set to -1. Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetGasComponentsCount, GetGasCompPressures, GetGasCompPhi, GetGasPhaseVolume, SetGasCompMoles, SetGasPhaseVolume.
C++ Example:
std::vector gas_moles;
status = phreeqc_rm.RunCells();
status = phreeqc_rm.GetGasCompMoles(gas_moles);
MPI:
Called by root, workers must be in the loop of MpiWorker.
const std::vector<std::string>& GetGasComponents ( void  ) const
inline

Returns a reference to the vector of all gas components in the initial-phreeqc module. The list includes all gas components included in any GAS_PHASE definitions in the initial-phreeqc module. FindComponents must be called before GetGasComponents. This method may be useful when generating selected output definitions related to gas phases.

Return values
conststd::vector<std::string>& A vector of strings; each string is a unique gas component name.
See also
FindComponents, GetGasComponentsCount.
C++ Example:
oss << "  -gases " << "\n";
// gas components
const std::vector &gas_phases = phreeqc_rm.GetGasComponents();
for (size_t i = 0; i < phreeqc_rm.GetGasComponentsCount(); i++)
{
oss << "    " << gas_phases[i] << "\n";
}
MPI:
Called by root.
int GetGasComponentsCount ( void  ) const
inline

Returns the number of gas phase components in the initial-phreeqc module. FindComponents must be called before GetGasComponentsCount. This method may be useful when generating selected output definitions related to gas phases.

Return values
Thenumber of gas phase components in the initial-phreeqc module.
See also
FindComponents, GetGasComponents.
C++ Example:
oss << "  -gases " << "\n";
// gas components
const std::vector &gas_phases = phreeqc_rm.GetGasComponents();
for (size_t i = 0; i < phreeqc_rm.GetGasComponentsCount(); i++)
{
oss << "    " << gas_phases[i] << "\n";
}
MPI:
Called by root.
IRM_RESULT GetGasCompPhi ( std::vector< double > &  gas_phi)

Transfer fugacity coefficients (phi) of gas components from each reaction cell to the vector given in the argument list (gas_phi). Fugacity is equal to the gas component pressure times the fugacity coefficient.

Parameters
gas_phiVector to receive the fugacity coefficients of gas components. Dimension of the vector is set to ngas_comps times nxyz, where, ngas_comps is the result of GetGasComponentsCount, and nxyz is the number of user grid cells (GetGridCellCount). If a gas component is not defined for a cell, the fugacity coefficient is set to -1. Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetGasComponentsCount, GetGasCompMoles, GetGasCompPressures, GetGasPhaseVolume, SetGasCompMoles, SetGasPhaseVolume.
C++ Example:
std::vector gas_phi;
status = phreeqc_rm.RunCells();
status = phreeqc_rm.GetGasCompPhi(gas_phi);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT GetGasCompPressures ( std::vector< double > &  gas_pressure)

Transfer pressures of gas components from each reaction cell to the vector given in the argument list (gas_pressure).

Parameters
gas_pressureVector to receive the pressures of gas components. Dimension of the vector is set to ngas_comps times nxyz, where, ngas_comps is the result of GetGasComponentsCount, and nxyz is the number of user grid cells (GetGridCellCount). If a gas component is not defined for a cell, the pressure is set to -1. Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetGasComponentsCount, GetGasCompMoles, GetGasCompPhi, GetGasPhaseVolume, SetGasCompMoles, SetGasPhaseVolume.
C++ Example:
std::vector gas_pressure;
status = phreeqc_rm.RunCells();
status = phreeqc_rm.GetGasCompPressures(gas_pressure);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT GetGasPhaseVolume ( std::vector< double > &  gas_volume)

Transfer volume of gas phase from each reaction cell to the vector given in the argument list (gas_volume).

Parameters
gas_volumeVector to receive the gas phase volumes. Dimension of the vector is set to nxyz, where, nxyz is the number of user grid cells (GetGridCellCount). If a gas phase is not defined for a cell, the volume is set to -1. Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetGasComponentsCount, GetGasCompMoles, GetGasCompPressures, GetGasCompPhi, SetGasCompMoles, SetGasPhaseVolume.
C++ Example:
std::vector gas_volume;
status = phreeqc_rm.RunCells();
status = phreeqc_rm.GetGasPhaseVolume(gas_volume);
MPI:
Called by root, workers must be in the loop of MpiWorker.
const std::vector< double >& GetGfw ( void  )
inline

Returns a reference to a vector of doubles that contains the gram-formula weight of each component. Called after FindComponents. Order of weights corresponds to the list of components from GetComponents.

Return values
conststd::vector<double>& A vector of doubles; each value is a component gram-formula weight, g/mol.
See also
FindComponents, GetComponents.
C++ Example:
const std::vector &components = phreeqc_rm.GetComponents();
const std::vector < double > & gfw = phreeqc_rm.GetGfw();
int ncomps = phreeqc_rm.GetComponentCount();
for (int i = 0; i < ncomps; i++)
{
  std::ostringstream strm;
  strm.width(10);
  strm << components[i] << "    " << gfw[i] << "\n";
  phreeqc_rm.OutputMessage(strm.str());
}
MPI:
Called by root and (or) workers.
int GetGridCellCount ( void  )
inline

Returns the number of grid cells in the user's model, which is defined in the call to the constructor for the reaction module. The mapping from grid cells to reaction cells is defined by CreateMapping. The number of reaction cells may be less than the number of grid cells if there are inactive regions or there is symmetry in the model definition.

Return values
Numberof grid cells in the user's model.
See also
PhreeqcRM::PhreeqcRM , CreateMapping.
C++ Example:
std::ostringstream oss;
oss << "Number of grid cells in the user's model: " << phreeqc_rm.GetGridCellCount() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
IPhreeqc* GetIPhreeqcPointer ( int  i)

Returns an IPhreeqc pointer to the ith IPhreeqc instance in the reaction module. For the threaded version, there are nthreads + 2 IPhreeqc instances, where nthreads is defined in the constructor (PhreeqcRM::PhreeqcRM). The number of threads can be determined by GetThreadCount. The first nthreads (0 based) instances will be the workers, the next (nthreads) is the InitialPhreeqc instance, and the next (nthreads + 1) is the Utility instance. Getting the IPhreeqc pointer for one of these instances allows the user to use any of the IPhreeqc methods on that instance. For MPI, each process has exactly three IPhreeqc instances, one worker (number 0), one InitialPhreeqc instance (number 1), and one Utility instance (number 2).

Parameters
iThe number of the IPhreeqc instance (0 based) to be retrieved.
Return values
IPhreeqcpointer to the ith IPhreeqc instance (0 based) in the reaction module.
See also
PhreeqcRM::PhreeqcRM, GetThreadCount. See IPhreeqc documentation for descriptions of IPhreeqc methods.
C++ Example:
// Utility pointer is worker nthreads + 1
IPhreeqc * util_ptr = phreeqc_rm.GetIPhreeqcPointer(phreeqc_rm.GetThreadCount() + 1);
MPI:
Called by root and (or) workers.
const std::vector<std::string>& GetKineticReactions ( void  ) const
inline

Returns a reference to the vector of all kinetic reactions in the initial-phreeqc module. The list includes all kinetic reactions included in any KINETICS definitions in the reaction model. FindComponents must be called before GetKineticReactions. This method may be useful when generating selected output definitions related to kinetic reactions.

Return values
conststd::vector<std::string>& A vector of strings; each string is a unique kinetic reaction name.
See also
FindComponents, GetKineticReactionsCount.
C++ Example:
oss << "  -kinetics " << "\n";
// kinetic reactions
const std::vector &kin_reactions = phreeqc_rm.GetKineticReactions();
for (size_t i = 0; i < phreeqc_rm.GetKineticReactionsCount(); i++)
{
oss << "    " << kin_reactions[i] << "\n";
}
MPI:
Called by root.
int GetKineticReactionsCount ( void  ) const
inline

Returns the number of kinetic reactions in the initial-phreeqc module. FindComponents must be called before GetKineticReactionsCount. This method may be useful when generating selected output definitions related to kinetic reactions.

Return values
Thenumber of kinetic reactions in the initial-phreeqc module.
See also
FindComponents, GetKineticReactions.
C++ Example:
oss << "  -kinetics " << "\n";
// kinetic reactions
const std::vector &kin_reactions = phreeqc_rm.GetKineticReactions();
for (size_t i = 0; i < phreeqc_rm.GetKineticReactionsCount(); i++)
{
oss << "    " << kin_reactions[i] << "\n";
}
MPI:
Called by root.
int GetMpiMyself ( void  ) const
inline

Returns the MPI process (task) number. For the MPI version, the root task number is zero, and all MPI tasks have unique task numbers greater than zero. The number of tasks can be obtained with GetMpiTasks. The number of tasks and computer hosts is determined at run time by the mpiexec command, and the number of reaction-module processes is defined by the communicator used in constructing the reaction modules (PhreeqcRM::PhreeqcRM). For the OPENMP version, the task number is always zero, and the result of GetMpiTasks is one.

Return values
TheMPI task number for a process.
See also
GetMpiTasks, PhreeqcRM::PhreeqcRM.
C++ Example:
std::ostringstream oss;
oss << "MPI task number: " << phreeqc_rm.GetMpiMyself() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
int GetMpiTasks ( void  ) const
inline

Returns the number of MPI processes (tasks) assigned to the reaction module. For the MPI version, the number of tasks and computer hosts is specified at run time by the mpiexec command. The number of MPI processes used for reaction calculations is determined by the MPI communicator used in constructing the reaction modules. The communicator may define a subset of the total number of MPI processes. The root task number is zero, and all other MPI tasks have unique task numbers greater than zero. For the OPENMP version, the number of tasks is one, and the task number returned by GetMpiMyself is zero.

Return values
Thenumber of MPI processes assigned to the reaction module.
See also
GetMpiMyself, PhreeqcRM::PhreeqcRM.
C++ Example:
std::ostringstream oss;
oss << "Number of MPI processes: " << phreeqc_rm.GetMpiTasks() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
int GetNthSelectedOutputUserNumber ( int  n)

Returns the user number for the nth selected-output definition. Definitions are sorted by user number. Phreeqc allows multiple selected-output definitions, each of which is assigned a nonnegative integer identifier by the user. The number of definitions can be obtained by GetSelectedOutputCount. To cycle through all of the definitions, GetNthSelectedOutputUserNumber can be used to identify the user number for each selected-output definition in sequence. SetCurrentSelectedOutputUserNumber is then used to select that user number for selected-output processing.

Parameters
nThe sequence number of the selected-output definition for which the user number will be returned. Fortran, 1 based; C, 0 based.
Return values
Theuser number of the nth selected-output definition, negative is failure (See DecodeError).
See also
GetSelectedOutput, GetSelectedOutputColumnCount, GetSelectedOutputCount, GetSelectedOutputHeading, GetSelectedOutputRowCount, SetCurrentSelectedOutputUserNumber, SetSelectedOutputOn.
C++ Example:
for (int isel = 0; isel < phreeqc_rm.GetSelectedOutputCount(); isel++)
{
  int n_user = phreeqc_rm.GetNthSelectedOutputUserNumber(isel);
  status = phreeqc_rm.SetCurrentSelectedOutputUserNumber(n_user);
  std::cerr << "Selected output sequence number: " << isel << "\n";
  std::cerr << "Selected output user number:     " << n_user << "\n";
  std::vector so;
  int col = phreeqc_rm.GetSelectedOutputColumnCount();
  status = phreeqc_rm.GetSelectedOutput(so);
  // Process results here
}
MPI:
Called by root.
bool GetPartitionUZSolids ( void  ) const
inline

Returns the setting for partitioning solids between the saturated and unsaturated parts of a partially saturated cell. The option is intended to be used by saturated-only flow codes that allow a variable water table. The value has meaning only when saturations less than 1.0 are encountered. The partially saturated cells may have a small water-to-rock ratio that causes reactions to proceed slowly relative to fully saturated cells. By setting SetPartitionUZSolids to true, the amounts of solids and gases are partioned according to the saturation. If a cell has a saturation of 0.5, then the water interacts with only half of the solids and gases; the other half is unreactive until the water table rises. As the saturation in a cell varies, solids and gases are transferred between the saturated and unsaturated (unreactive) reservoirs of the cell. Unsaturated-zone flow and transport codes will probably use the default (false), which assumes all gases and solids are reactive regardless of saturation.

Return values
boolTrue, the fraction of solids and gases available for reaction is equal to the saturation; False (default), all solids and gases are reactive regardless of saturation.
See also
SetPartitionUZSolids.
C++ Example:
oss << "Partioning of UZ solids: " << phreeqc_rm.GetPartitionUZSolids();
MPI:
Called by root and (or) workers.
const std::vector<double>& GetPressure ( void  )

Returns the pressure for each cell. By default, the pressure vector is initialized with 1 atm; if SetPressure has not been called, worker solutions will have pressures as defined in input files (RunFile) or input strings (RunString); if SetPressure has been called, worker solutions will have the pressures as defined by SetPressure. Pressure effects are considered by three PHREEQC databases: phreeqc.dat, Amm.dat, and pitzer.dat.

Return values
conststd::vector<double>& A vector reference to the pressures in each cell, in atm. Size of vector is nxyz, the number of grid cells in the user's model (GetGridCellCount).
See also
SetPressure, GetTemperature, SetTemperature.
C++ Example:
const std::vector & p_atm = phreeqc_rm.GetPressure();
MPI:
Called by root and (or) workers.
const std::vector<int>& GetPrintChemistryMask ( void  )
inline

Return a reference to the vector of print flags that enable or disable detailed output for each cell. Printing for a cell will occur only when the printing is enabled with SetPrintChemistryOn, and the value in the vector for the cell is 1.

Return values
std::vector<int>& Vector of integers. Size of vector is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount). A value of 0 for a cell indicates printing is disabled; a value of 1 for a cell indicates printing is enabled.
See also
SetPrintChemistryOn.
C++ Example:
const std::vector & print_chemistry_mask1 = phreeqc_rm.GetPrintChemistryMask();
MPI:
Called by root.
const std::vector<bool>& GetPrintChemistryOn ( void  ) const
inline

Returns a vector reference to the current print flags for detailed output for the three sets of IPhreeqc instances: the workers, the InitialPhreeqc instance, and the Utility instance. Dimension of the vector is 3. Printing of detailed output from reaction calculations to the output file is enabled when the vector value is true, disabled when false. The detailed output prints all of the output typical of a PHREEQC reaction calculation, which includes solution descriptions and the compositions of all other reactants. The output can be several hundred lines per cell, which can lead to a very large output file (prefix.chem.txt, OpenFiles). For the worker instances, the output can be limited to a set of cells (SetPrintChemistryMask) and, in general, the amount of information printed can be limited by use of options in the PRINT data block of PHREEQC (applied by using RunFile or RunString). Printing the detailed output for the workers is generally used only for debugging, and PhreeqcRM will run faster when printing detailed output for the workers is disabled (SetPrintChemistryOn).

Return values
conststd::vector<bool> & Print flag for the workers, InitialPhreeqc, and Utility IPhreeqc instances, respectively.
See also
SetPrintChemistryOn, SetPrintChemistryMask.
C++ Example:
const std::vector & print_on = phreeqc_rm.GetPrintChemistryOn();
MPI:
Called by root and (or) workers.
bool GetRebalanceByCell ( void  ) const
inline

Get the load-rebalancing method used for parallel processing. PhreeqcRM attempts to rebalance the load of each thread or process such that each thread or process takes the same amount of time to run its part of a RunCells calculation. Two algorithms are available: one accounts for cells that were not run because saturation was zero (true), and the other uses the average time to run all of the cells assigned to a process or thread (false), . The methods are similar, but preliminary results indicate the default is better in most cases.

Return values
boolTrue indicates individual cell run times are used in rebalancing (default); False, indicates average run times are used in rebalancing.
See also
GetRebalanceFraction, SetRebalanceByCell, SetRebalanceFraction.
C++ Example:
bool rebalance = phreeqc_rm.GetRebalanceByCell();
MPI:
Called by root and (or) workers.
double GetRebalanceFraction ( void  ) const
inline

Get the fraction used to determine the number of cells to transfer among threads or processes. PhreeqcRM attempts to rebalance the load of each thread or process such that each thread or process takes the same amount of time to run its part of a RunCells calculation. The rebalancing transfers cell calculations among threads or processes to try to achieve an optimum balance. SetRebalanceFraction adjusts the calculated optimum number of cell transfers by a fraction from 0 to 1.0 to determine the number of cell transfers that actually are made. A value of zero eliminates load rebalancing. A value less than 1.0 is suggested to avoid possible oscillations, where too many cells are transferred at one iteration, requiring reverse transfers at the next iteration. Default is 0.5.

Return values
intFraction used in rebalance, 0.0 to 1.0.
See also
GetRebalanceByCell, SetRebalanceByCell, SetRebalanceFraction.
C++ Example:
double f_rebalance = phreeqc_rm.GetRebalanceFraction();
MPI:
Called by root.
IRM_RESULT GetSaturation ( std::vector< double > &  sat)

Returns a vector of saturations (sat) as calculated by the reaction module. Reactions will change the volume of solution in a cell. The transport code must decide whether to ignore or account for this change in solution volume due to reactions. Following reactions, the cell saturation is calculated as solution volume (GetSolutionVolume) divided by the product of representative volume (SetRepresentativeVolume) and the porosity (SetPorosity). The cell saturation returned by GetSaturation may be less than or greater than the saturation set by the transport code (SetSaturation), and may be greater than or less than 1.0, even in fully saturated simulations. Only the following databases distributed with PhreeqcRM have molar volume information needed to accurately calculate solution volume and saturation: phreeqc.dat, Amm.dat, and pitzer.dat.

Parameters
satVector to receive the saturations. Dimension of the array is set to nxyz, where nxyz is the number of user grid cells (GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetSolutionVolume, SetPorosity, SetRepresentativeVolume, SetSaturation.
C++ Example:
std::vector sat;
status = phreeqc_rm.GetSaturation(sat);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT GetSelectedOutput ( std::vector< double > &  so)

Returns the array of selected-output values for the current selected-output definition. SetCurrentSelectedOutputUserNumber specifies which of the selected-output definitions is returned to the vector (so).

Parameters
soA vector to contain the selected-output values. Size of the vector is set to col times nxyz, where col is the number of columns in the selected-output definition (GetSelectedOutputColumnCount), and nxyz is the number of grid cells in the user's model (GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetNthSelectedOutputUserNumber, GetSelectedOutputColumnCount, GetSelectedOutputCount, GetSelectedOutputHeading, GetSelectedOutputRowCount, SetCurrentSelectedOutputUserNumber, SetSelectedOutputOn.
C++ Example:
for (int isel = 0; isel < phreeqc_rm.GetSelectedOutputCount(); isel++)
{
  int n_user = phreeqc_rm.GetNthSelectedOutputUserNumber(isel);
  status = phreeqc_rm.SetCurrentSelectedOutputUserNumber(n_user);
  std::vector so;
  int col = phreeqc_rm.GetSelectedOutputColumnCount();
  status = phreeqc_rm.GetSelectedOutput(so);
  // Process results here
}
MPI:
Called by root, workers must be in the loop of MpiWorker.
int GetSelectedOutputColumnCount ( void  )

Returns the number of columns in the current selected-output definition. SetCurrentSelectedOutputUserNumber specifies which of the selected-output definitions is used.

Return values
Numberof columns in the current selected-output definition, negative is failure (See DecodeError).
See also
GetNthSelectedOutputUserNumber, GetSelectedOutput, GetSelectedOutputCount, GetSelectedOutputHeading, GetSelectedOutputRowCount, SetCurrentSelectedOutputUserNumber, SetSelectedOutputOn.
C++ Example:
for (int isel = 0; isel < phreeqc_rm.GetSelectedOutputCount(); isel++)
{
  int n_user = phreeqc_rm.GetNthSelectedOutputUserNumber(isel);
  status = phreeqc_rm.SetCurrentSelectedOutputUserNumber(n_user);
  std::vector so;
  int col = phreeqc_rm.GetSelectedOutputColumnCount();
  status = phreeqc_rm.GetSelectedOutput(so);
  // Print results
  for (int i = 0; i < phreeqc_rm.GetSelectedOutputRowCount()/2; i++)
  {
    std::vector headings;
    headings.resize(col);
    std::cerr << "     Selected output: " << "\n";
    for (int j = 0; j < col; j++)
    {
      status = phreeqc_rm.GetSelectedOutputHeading(j, headings[j]);
      std::cerr << "          " << j << " " << headings[j] << ": " << so[j*nxyz + i] << "\n";
    }
  }
}
MPI:
Called by root.
int GetSelectedOutputCount ( void  )

Returns the number of selected-output definitions. SetCurrentSelectedOutputUserNumber specifies which of the selected-output definitions is used.

Return values
Numberof selected-output definitions, negative is failure (See DecodeError).
See also
GetNthSelectedOutputUserNumber, GetSelectedOutput, GetSelectedOutputColumnCount, GetSelectedOutputHeading, GetSelectedOutputRowCount, SetCurrentSelectedOutputUserNumber, SetSelectedOutputOn.
C++ Example:
for (int isel = 0; isel < phreeqc_rm.GetSelectedOutputCount(); isel++)
{
  int n_user = phreeqc_rm.GetNthSelectedOutputUserNumber(isel);
  status = phreeqc_rm.SetCurrentSelectedOutputUserNumber(n_user);
  std::vector so;
  int col = phreeqc_rm.GetSelectedOutputColumnCount();
  status = phreeqc_rm.GetSelectedOutput(so);
  // Process results here
}
MPI:
Called by root.
IRM_RESULT GetSelectedOutputHeading ( int  icol,
std::string &  heading 
)

Returns a selected-output heading. The number of headings is determined by GetSelectedOutputColumnCount. SetCurrentSelectedOutputUserNumber specifies which of the selected-output definitions is used.

Parameters
icolThe sequence number of the heading to be retrieved, 0 based.
headingA string to receive the heading.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetNthSelectedOutputUserNumber, GetSelectedOutput, GetSelectedOutputColumnCount, GetSelectedOutputCount, GetSelectedOutputRowCount, SetCurrentSelectedOutputUserNumber, SetSelectedOutputOn.
C++ Example:
for (int isel = 0; isel < phreeqc_rm.GetSelectedOutputCount(); isel++)
{
  int n_user = phreeqc_rm.GetNthSelectedOutputUserNumber(isel);
  status = phreeqc_rm.SetCurrentSelectedOutputUserNumber(n_user);
  std::vector so;
  int col = phreeqc_rm.GetSelectedOutputColumnCount();
  status = phreeqc_rm.GetSelectedOutput(so);
  std::vector headings;
  headings.resize(col);
  std::cerr << "     Selected output: " << "\n";
  for (int j = 0; j < col; j++)
  {
    status = phreeqc_rm.GetSelectedOutputHeading(j, headings[j]);
    std::cerr << "          " << j << " " << headings[j] << "\n";
  }
}
MPI:
Called by root.
bool GetSelectedOutputOn ( void  )
inline

Returns the current value of the selected-output property. A value of true for this property indicates that selected output data will be requested this time step. A value of false indicates that selected output will not be retrieved for this time step; processing the selected output is avoided with some time savings.

Return values
boolTrue, selected output will be requested; false, selected output will not be retrieved.
See also
SetSelectedOutputOn.
C++ Example:
bool so_on = phreeqc_rm.GetSelectedOutputOn();
MPI:
Called by root and (or) workers.
int GetSelectedOutputRowCount ( void  )

Returns the number of rows in the current selected-output definition. However, the method is included only for convenience; the number of rows is always equal to the number of grid cells in the user's model (GetGridCellCount).

Return values
Numberof rows in the current selected-output definition, negative is failure (See DecodeError).
See also
GetNthSelectedOutputUserNumber, GetSelectedOutput, GetSelectedOutputColumnCount, GetSelectedOutputCount, GetSelectedOutputHeading, SetCurrentSelectedOutputUserNumber, SetSelectedOutputOn.
C++ Example:
for (int isel = 0; isel < phreeqc_rm.GetSelectedOutputCount(); isel++)
{
  int n_user = phreeqc_rm.GetNthSelectedOutputUserNumber(isel);
  status = phreeqc_rm.SetCurrentSelectedOutputUserNumber(n_user);
  std::vector so;
  int col = phreeqc_rm.GetSelectedOutputColumnCount();
  status = phreeqc_rm.GetSelectedOutput(so);
  // Print results
  for (int i = 0; i < phreeqc_rm.GetSelectedOutputRowCount()/2; i++)
  {
    std::vector headings;
    headings.resize(col);
    std::cerr << "     Selected output: " << "\n";
    for (int j = 0; j < col; j++)
    {
      status = phreeqc_rm.GetSelectedOutputHeading(j, headings[j]);
      std::cerr << "          " << j << " " << headings[j] << ": " << so[j*nxyz + i] << "\n";
    }
  }
}
MPI:
Called by root.
int GetSICount ( void  ) const
inline

Returns the number of phases in the initial-phreeqc module for which saturation indices could be calculated. FindComponents must be called before GetSICount. This method may be useful when generating selected output definitions related to saturation indices.

Return values
Thenumber of phases in the initial-phreeqc module for which saturation indices could be calculated.
See also
FindComponents, GetSINames.
C++ Example:
const std::vector &si = phreeqc_rm.GetSINames();
for (size_t i = 0; i < phreeqc_rm.GetSICount(); i++)
{
oss << "    " << si[i] << "\n";
}
MPI:
Called by root.
const std::vector<std::string>& GetSINames ( void  ) const
inline

Returns a reference to the vector of the names of all phases for which saturation indices (SIs) could be calculated. The list includes all phases that contain only elements included in the components in the initial-phreeqc module. The list assumes that all components are present to be able to calculate the entire list of SIs; it may be that one or more components are missing in any specific cell. FindComponents must be called before GetSINames. This method may be useful when generating selected output definitions related to saturation indices.

Return values
conststd::vector<std::string>& A vector of strings; each string is a unique phase name.
See also
FindComponents, GetSICount.
C++ Example:
oss << "  -saturation_indices " << "\n";
// molalities of aqueous species
const std::vector &si = phreeqc_rm.GetSINames();
for (size_t i = 0; i < phreeqc_rm.GetSICount(); i++)
{
oss << "    " << si[i] << "\n";
}
MPI:
Called by root.
const std::vector<std::string>& GetSolidSolutionComponents ( void  ) const
inline

Returns a reference to the vector of solid solution components. The list of solid solution components includes all components in any SOLID_SOLUTION definitions in the initial-phreeqc module. FindComponents must be called before GetSolidSolutionComponents. This method may be useful when generating selected output definitions related to solid solutions.

Return values
conststd::vector<std::string>& A vector of strings; each string is a unique solid solution component.
See also
FindComponents, GetSolidSolutionComponentsCount, GetSolidSolutionNames.
C++ Example:
oss << "  -solid_solutions " << "\n";
// solid solutions
const std::vector &ss_comps = phreeqc_rm.GetSolidSolutionComponents();
const std::vector &ss_names = phreeqc_rm.GetSolidSolutionNames();
for (size_t i = 0; i < phreeqc_rm.GetSolidSolutionComponentsCount(); i++)
{

oss << "    ";
oss.width(15);
oss  << std::left << ss_comps[i];
oss << " # " << ss_names[i] << "\n";
}
MPI:
Called by root.
int GetSolidSolutionComponentsCount ( void  ) const
inline

Returns the number of solid solution components in the initial-phreeqc module. FindComponents must be called before GetSolidSolutionComponentsCount. This method may be useful when generating selected output definitions related to solid solutions.

Return values
Thenumber of solid solution components in the initial-phreeqc module.
See also
FindComponents, GetSolidSolutionComponents, GetSolidSolutionNames.
C++ Example:
oss << "  -solid_solutions " << "\n";
// solid solutions
const std::vector &ss_comps = phreeqc_rm.GetSolidSolutionComponents();
const std::vector &ss_names = phreeqc_rm.GetSolidSolutionNames();
for (size_t i = 0; i < phreeqc_rm.GetSolidSolutionComponentsCount(); i++)
{

oss << "    ";
oss.width(15);
oss  << std::left << ss_comps[i];
oss << " # " << ss_names[i] << "\n";
}
MPI:
Called by root.
const std::vector<std::string>& GetSolidSolutionNames ( void  ) const
inline

Returns a reference to the vector of solid solution names that correspond with the solid solution components. FindComponents must be called before GetSolidSolutionNames. The solid solution names vector is the same length as the solid solution components vector and provides the corresponding name of solid solution containing the component. This method may be useful when generating selected output definitions related to solid solutions.

Return values
conststd::vector<std::string>& A vector of strings; each string is a solid solution name corresponding to the solid solution components vector; a solid solution name may occur multiple times.
See also
FindComponents, GetSolidSolutionComponentsCount, GetSolidSolutionComponents.
C++ Example:
oss << "  -solid_solutions " << "\n";
// solid solutions
const std::vector &ss_comps = phreeqc_rm.GetSolidSolutionComponents();
const std::vector &ss_names = phreeqc_rm.GetSolidSolutionNames();
for (size_t i = 0; i < phreeqc_rm.GetSolidSolutionComponentsCount(); i++)
{

oss << "    ";
oss.width(15);
oss  << std::left << ss_comps[i];
oss << " # " << ss_names[i] << "\n";
}
MPI:
Called by root.
const std::vector<double>& GetSolutionVolume ( void  )

Return a vector reference to the current solution volumes as calculated by the reaction module. Dimension of the vector will be nxyz, where nxyz is the number of user grid cells. Values for inactive cells are set to 1e30. Only the following databases distributed with PhreeqcRM have molar volume information needed to accurately calculate solution volume: phreeqc.dat, Amm.dat, and pitzer.dat.

Return values
Vectorreference to current solution volumes.
See also
GetSaturation.
C++ Example:
status = phreeqc_rm.RunCells();
const std::vector &volume = phreeqc_rm.GetSolutionVolume();
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT GetSpeciesConcentrations ( std::vector< double > &  species_conc)

Returns a vector reference to aqueous species concentrations (species_conc). This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by FindComponents and includes all aqueous species that can be made from the set of components. Solution volumes used to calculate mol/L are calculated by the reaction module. Only the following databases distributed with PhreeqcRM have molar volume information needed to accurately calculate solution volume: phreeqc.dat, Amm.dat, and pitzer.dat.

Parameters
species_concVector to receive the aqueous species concentrations. Dimension of the vector is set to nspecies times nxyz, where nspecies is the number of aqueous species (GetSpeciesCount), and nxyz is the number of grid cells (GetGridCellCount). Concentrations are moles per liter. Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int npecies = phreeqc_rm.GetSpeciesCount();
status = phreeqc_rm.RunCells();
std::vector c;
status = phreeqc_rm.GetSpeciesConcentrations(c);
MPI:
Called by root, workers must be in the loop of MpiWorker.
int GetSpeciesCount ( void  )
inline

Returns the number of aqueous species used in the reaction module. This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by FindComponents and includes all aqueous species that can be made from the set of components.

Return values
intThe number of aqueous species.
See also
FindComponents, GetSpeciesConcentrations, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int npecies = phreeqc_rm.GetSpeciesCount();
MPI:
Called by root and (or) workers.
const std::vector<double>& GetSpeciesD25 ( void  )
inline

Returns a vector reference to diffusion coefficients at 25C for the set of aqueous species. This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true. Diffusion coefficients are defined in SOLUTION_SPECIES data blocks, normally in the database file. Databases distributed with the reaction module that have diffusion coefficients defined are phreeqc.dat, Amm.dat, and pitzer.dat.

Return values
Vectorreference to the diffusion coefficients at 25 C, m^2/s. Dimension of the vector is nspecies, where nspecies is the number of aqueous species (GetSpeciesCount).
See also
FindComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int npecies = phreeqc_rm.GetSpeciesCount();
const std::vector < double > & species_d = phreeqc_rm.GetSpeciesD25();
MPI:
Called by root and (or) workers.
IRM_RESULT GetSpeciesLog10Gammas ( std::vector< double > &  species_log10gammas)

Returns a vector reference to log10 aqueous species activity coefficients (species_log10gammas). This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by FindComponents and includes all aqueous species that can be made from the set of components.

Parameters
species_log10gammasVector to receive the log10 aqueous species activity coefficients. Dimension of the vector is set to nspecies times nxyz, where nspecies is the number of aqueous species (GetSpeciesCount), and nxyz is the number of grid cells (GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int npecies = phreeqc_rm.GetSpeciesCount();
status = phreeqc_rm.RunCells();
std::vector species_gammas;
status = phreeqc_rm.GetSpeciesLog10Gammas(species_gammas);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT GetSpeciesLog10Molalities ( std::vector< double > &  species_log10molalities)

Returns a vector reference to log10 aqueous species molalities (species_log10molalities). To use this method SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by FindComponents and includes all aqueous species that can be made from the set of components.

Parameters
species_log10molalitiesVector to receive the log10 aqueous species molalites. Dimension of the vector is set to nspecies times nxyz, where nspecies is the number of aqueous species (GetSpeciesCount), and nxyz is the number of grid cells (GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int npecies = phreeqc_rm.GetSpeciesCount();
status = phreeqc_rm.RunCells();
std::vector species_molalities;
status = phreeqc_rm.GetSpeciesLog10Molalities(species_molalities);
MPI:
Called by root, workers must be in the loop of MpiWorker.
const std::vector<std::string>& GetSpeciesNames ( void  )
inline

Returns a vector reference to the names of the aqueous species. This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by FindComponents and includes all aqueous species that can be made from the set of components.

Return values
namesVector of strings containing the names of the aqueous species. Dimension of the vector is nspecies, where nspecies is the number of aqueous species (GetSpeciesCount).
See also
FindComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int npecies = phreeqc_rm.GetSpeciesCount();
const std::vector &species = phreeqc_rm.GetSpeciesNames();
MPI:
Called by root and (or) workers.
bool GetSpeciesSaveOn ( void  )
inline

Returns the value of the species-save property. By default, concentrations of aqueous species are not saved. Setting the species-save property to true allows aqueous species concentrations to be retrieved with GetSpeciesConcentrations, and solution compositions to be set with SpeciesConcentrations2Module.

Return values
Trueindicates solution species concentrations are saved and can be used for multicomponent-diffusion calculations; False indicates that solution species concentrations are not saved.
See also
FindComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesSaveOn, GetSpeciesZ, GetSpeciesNames, SpeciesConcentrations2Module.
FindComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesStoichiometry, GetSpeciesZ, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int npecies = phreeqc_rm.GetSpeciesCount();
bool species_on = phreeqc_rm.GetSpeciesSaveOn();
MPI:
Called by root and (or) workers.
const std::vector<cxxNameDouble>& GetSpeciesStoichiometry ( void  )
inline

Returns a vector reference to the stoichiometry of each aqueous species. This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true.

Return values
Vectorof cxxNameDouble instances (maps) that contain the component names and associated stoichiometric coefficients for each aqueous species. Dimension of the vector is nspecies, where nspecies is the number of aqueous species (GetSpeciesCount).
See also
FindComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesZ, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
const std::vector &species = phreeqc_rm.GetSpeciesNames();
const std::vector < double > & species_z = phreeqc_rm.GetSpeciesZ();
const std::vector < double > & species_d = phreeqc_rm.GetSpeciesD25();
bool species_on = phreeqc_rm.GetSpeciesSaveOn();
int nspecies = phreeqc_rm.GetSpeciesCount();
for (int i = 0; i < nspecies; i++)
{
  std::ostringstream strm;
  strm << species[i] << "\n";
  strm << "    Charge: " << species_z[i] << std::endl;
  strm << "    Dw:     " << species_d[i] << std::endl;
  cxxNameDouble::const_iterator it = phreeqc_rm.GetSpeciesStoichiometry()[i].begin();
  for (; it != phreeqc_rm.GetSpeciesStoichiometry()[i].end(); it++)
  {
    strm << "          " << it->first << "   " << it->second << "\n";
  }
  phreeqc_rm.OutputMessage(strm.str());
}
MPI:
Called by root and (or) workers.
const std::vector<double>& GetSpeciesZ ( void  )
inline

Returns a vector reference to the charge on each aqueous species. This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true.

Return values
Vectorcontaining the charge on each aqueous species. Dimension of the vector is nspecies, where nspecies is the number of aqueous species (GetSpeciesCount).
See also
FindComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, SetSpeciesSaveOn, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int npecies = phreeqc_rm.GetSpeciesCount();
const std::vector < double > & species_z = phreeqc_rm.GetSpeciesZ();
MPI:
Called by root and (or) workers.
const std::vector< int>& GetStartCell ( void  ) const
inline

Returns a vector of integers that contains the smallest reaction-cell number assigned to each worker. Each worker is assigned a range of reaction-cell numbers that are run during a call to RunCells. The range of reaction cell numbers for a worker may vary as load rebalancing occurs. At any point in the calculations, the first cell and last cell to be run by a worker can be found in the vectors returned by GetStartCell and GetEndCell. Each method returns a vector of integers that has size of the number of threads (GetThreadCount), if using OPENMP, or the number of processes (GetMpiTasks), if using MPI.

Return values
IRM_RESULTVector of integers, one for each worker, that gives the first reaction cell to be run by each worker.
See also
GetEndCell, GetThreadCount, GetMpiTasks, RunCells.
C++ Example:
std::ostringstream oss;
oss << "Current distribution of cells for workers\n";
oss << "Worker First Cell   Last Cell\n";
int n;
n = phreeqc_rm.GetThreadCount() * phreeqc_rm.GetMpiTasks();
for (int i = 0; i < n; i++)
{
    oss << i << "      "
        << phreeqc_rm.GetStartCell()[i]
        << "            "
        << phreeqc_rm.GetEndCell()[i] << "\n";
}
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers.
const std::vector<std::string>& GetSurfaceNames ( void  ) const
inline

Returns a reference to the vector of surface names (such as "Hfo") that correspond with the surface species names. The vectors referenced by GetSurfaceSpecies and GetSurfaceNames are the same length. FindComponents must be called before GetSurfaceNames. This method may be useful when generating selected output definitions related to surfaces.

Return values
conststd::vector<std::string>& A vector of strings; each string is a surface name corresponding to the surface species vector; a surface name may occur multiple times.
See also
FindComponents, GetSurfaceSpeciesCount, GetSurfaceSpecies, GetSurfaceTypes.
C++ Example:
// molalities of surface species
const std::vector &surf_species = phreeqc_rm.GetSurfaceSpecies();
const std::vector &surf_types = phreeqc_rm.GetSurfaceTypes();
const std::vector &surf_names = phreeqc_rm.GetSurfaceNames();
for (size_t i = 0; i < phreeqc_rm.GetSurfaceSpeciesCount(); i++)
{
oss << "    ";
oss.width(15);
oss << std::left << surf_species[i];
oss << " # ";
oss.width(15);
oss << surf_types[i] << "   " << surf_names[i] << "\n";
}
MPI:
Called by root.
const std::vector<std::string>& GetSurfaceSpecies ( void  ) const
inline

Returns a reference to the vector of surface species names (such as "Hfo_wOH"). The list of surface species is derived from the list of components (FindComponents) and the list of all surface site types (such as "Hfo_w") that are included in SURFACE definitions in the initial-phreeqc module. FindComponents must be called before GetSurfaceSpecies. This method may be useful when generating selected output definitions related to surfaces.

Return values
conststd::vector<std::string>& A vector of strings; each string is a unique surface species name.
See also
FindComponents, GetSurfaceSpeciesCount, GetSurfaceTypes, GetSurfaceNames.
C++ Example:
// molalities of surface species
const std::vector &surf_species = phreeqc_rm.GetSurfaceSpecies();
const std::vector &surf_types = phreeqc_rm.GetSurfaceTypes();
const std::vector &surf_names = phreeqc_rm.GetSurfaceNames();
for (size_t i = 0; i < phreeqc_rm.GetSurfaceSpeciesCount(); i++)
{
oss << "    ";
oss.width(15);
oss << std::left << surf_species[i];
oss << " # ";
oss.width(15);
oss << surf_types[i] << "   " << surf_names[i] << "\n";
}
MPI:
Called by root.
int GetSurfaceSpeciesCount ( void  ) const
inline

Returns the number of surface species (such as "Hfo_wOH") in the initial-phreeqc module. FindComponents must be called before GetSurfaceSpeciesCount. This method may be useful when generating selected output definitions related to surfaces.

Return values
Thenumber of surface species in the initial-phreeqc module.
See also
FindComponents, GetSurfaceSpecies, GetSurfaceTypes, GetSurfaceNames.
C++ Example:
// molalities of surface species
const std::vector &surf_species = phreeqc_rm.GetSurfaceSpecies();
const std::vector &surf_types = phreeqc_rm.GetSurfaceTypes();
const std::vector &surf_names = phreeqc_rm.GetSurfaceNames();
for (size_t i = 0; i < phreeqc_rm.GetSurfaceSpeciesCount(); i++)
{
oss << "    ";
oss.width(15);
oss << std::left << surf_species[i];
oss << " # ";
oss.width(15);
oss << surf_types[i] << "   " << surf_names[i] << "\n";
}
MPI:
Called by root.
const std::vector<std::string>& GetSurfaceTypes ( void  ) const
inline

Returns a reference to the vector of surface site types (such as "Hfo_w") that correspond with the surface species names. The vectors referenced by GetSurfaceSpecies and GetSurfaceTypes are the same length. FindComponents must be called before GetSurfaceTypes. This method may be useful when generating selected output definitions related to surfaces.

Return values
conststd::vector<std::string>& A vector of strings; each string is a surface site type for the corresponding species in the surface species vector; a surface site type may occur multiple times.
See also
FindComponents, GetSurfaceSpeciesCount, GetSurfaceSpecies, GetSurfaceNames.
C++ Example:
// molalities of surface species
const std::vector &surf_species = phreeqc_rm.GetSurfaceSpecies();
const std::vector &surf_types = phreeqc_rm.GetSurfaceTypes();
const std::vector &surf_names = phreeqc_rm.GetSurfaceNames();
for (size_t i = 0; i < phreeqc_rm.GetSurfaceSpeciesCount(); i++)
{
oss << "    ";
oss.width(15);
oss << std::left << surf_species[i];
oss << " # ";
oss.width(15);
oss << surf_types[i] << "   " << surf_names[i] << "\n";
}
MPI:
Called by root.
const std::vector<double>& GetTemperature ( void  )

Vector reference to the current temperatures of the cells. By default, the temperature vector is initialized to 25 C; if SetTemperature has not been called, worker solutions will have temperatures as defined in input files (RunFile) or input strings (RunString); if SetTemperature has been called, worker solutions will have the temperatures as defined by SetTemperature.

Return values
Vectorof temperatures, in degrees C. Size of vector is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount).
See also
SetTemperature, GetPressure, SetPressure.
C++ Example:
const std::vector &  tempc = phreeqc_rm.GetTemperature();
MPI:
Called by root.
int GetThreadCount ( )
inline

Returns the number of threads, which is equal to the number of workers used to run in parallel with OPENMP. For the OPENMP version, the number of threads is set implicitly or explicitly with the constructor (PhreeqcRM::PhreeqcRM). For the MPI version, the number of threads is always one for each process.

Return values
Thenumber of threads used for OPENMP parallel processing.
See also
GetMpiTasks.
C++ Example:
std::ostringstream oss;
oss << "Number of threads: " << phreeqc_rm.GetThreadCount() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root and (or) workers; result is always 1.
double GetTime ( void  ) const
inline

Returns the current simulation time in seconds. The reaction module does not change the time value, so the returned value is equal to the default (0.0) or the last time set by SetTime.

Return values
Thecurrent simulation time, in seconds.
See also
GetTimeConversion, GetTimeStep, SetTime, SetTimeConversion, SetTimeStep.
C++ Example:
std::ostringstream strm;
strm << "Beginning transport calculation "
     << phreeqc_rm.GetTime() * phreeqc_rm.GetTimeConversion()
     << " days\n";
phreeqc_rm.LogMessage(strm.str());
MPI:
Called by root and (or) workers.
double GetTimeConversion ( void  )
inline

Returns a multiplier to convert time from seconds to another unit, as specified by the user. The reaction module uses seconds as the time unit. The user can set a conversion factor (SetTimeConversion) and retrieve it with GetTimeConversion. The reaction module only uses the conversion factor when printing the long version of cell chemistry (SetPrintChemistryOn), which is rare. Default conversion factor is 1.0.

Return values
Multiplierto convert seconds to another time unit.
See also
GetTime, GetTimeStep, SetTime, SetTimeConversion, SetTimeStep.
C++ Example:
std::ostringstream strm;
strm << "Beginning transport calculation "
     <<   phreeqc_rm.GetTime() * phreeqc_rm.GetTimeConversion()
     << " days\n";
phreeqc_rm.LogMessage(strm.str());
MPI:
Called by root and (or) workers.
double GetTimeStep ( void  )
inline

Returns the current simulation time step in seconds. This is the time over which kinetic reactions are integrated in a call to RunCells. The reaction module does not change the time-step value, so the returned value is equal to the default (0.0) or the last time step set by SetTimeStep.

Return values
Thecurrent simulation time step, in seconds.
See also
GetTime, GetTimeConversion, SetTime, SetTimeConversion, SetTimeStep.
C++ Example:
std::ostringstream strm;
strm << "Time step "
     << phreeqc_rm.GetTimeStep() * phreeqc_rm.GetTimeConversion()
     << " days\n";
phreeqc_rm.LogMessage(strm.str());
MPI:
Called by root and (or) workers.
int GetUnitsExchange ( void  )
inline

Returns the input units for exchangers. In PHREEQC input, exchangers are defined by moles of exchange sites (Mp). SetUnitsExchange specifies how the number of moles of exchange sites in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

Return values
Inputunits for exchangers.
See also
SetUnitsExchange, SetPorosity, SetRepresentativeVolume.
C++ Example:
int units_exchange = phreeqc_rm.GetUnitsExchange();
MPI:
Called by root and (or) workers.
int GetUnitsGasPhase ( void  )
inline

Returns the input units for gas phases. In PHREEQC input, gas phases are defined by moles of component gases (Mp). SetUnitsGasPhase specifies how the number of moles of component gases in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

Return values
Inputunits for gas phases (0, 1, or 2).
See also
SetUnitsGasPhase, SetPorosity, SetRepresentativeVolume.
C++ Example:
int units_gas_phase = phreeqc_rm.GetUnitsGasPhase();
MPI:
Called by root and (or) workers.
int GetUnitsKinetics ( void  )
inline

Returns the input units for kinetic reactants. In PHREEQC input, kinetics are defined by moles of kinetic reactants (Mp). SetUnitsKinetics specifies how the number of moles of kinetic reactants in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

Return values
Inputunits for kinetic reactants (0, 1, or 2).
See also
SetUnitsKinetics, SetPorosity, SetRepresentativeVolume.
C++ Example:
int units_kinetics = phreeqc_rm.GetUnitsKinetics();
MPI:
Called by root and (or) workers.
int GetUnitsPPassemblage ( void  )
inline

Returns the input units for pure phase assemblages (equilibrium phases). In PHREEQC input, equilibrium phases are defined by moles of each phase (Mp). SetUnitsPPassemblage specifies how the number of moles of phases in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

Return values
Inputunits for equilibrium phases (0, 1, or 2).
See also
SetUnitsPPassemblage, SetPorosity, SetRepresentativeVolume.
C++ Example:
int units_pp_assemblage = phreeqc_rm.GetUnitsPPassemblage();
MPI:
Called by root and (or) workers.
int GetUnitsSolution ( void  )
inline

Returns the units of concentration used by the transport model. Options are 1, mg/L; 2 mol/L; or 3, mass fraction, kg/kgs. In PHREEQC, solutions are defined by the number of moles of each element in the solution. The units of transport concentration are used when transport concentrations are converted to solution moles by SetConcentrations and Concentrations2Utility. The units of solution concentration also are used when solution moles are converted to transport concentrations by GetConcentrations.

To convert from mg/L to moles of element in the representative volume of a reaction cell, mg/L is converted to mol/L and multiplied by the solution volume, which is the product of porosity (SetPorosity), saturation (SetSaturation), and representative volume (SetRepresentativeVolume). To convert from mol/L to moles of element in a cell, mol/L is multiplied by the solution volume. To convert from mass fraction to moles of element in a cell, kg/kgs is converted to mol/kgs, multiplied by density (SetDensity) and multiplied by the solution volume.

To convert from moles of element in the representative volume of a reaction cell to mg/L, the number of moles of an element is divided by the solution volume resulting in mol/L, and then converted to mg/L. To convert from moles of element in the representative volume of a reaction cell to mol/L, the number of moles of an element is divided by the solution volume resulting in mol/L. To convert from moles of element in the representative volume of a reaction cell to mass fraction, the number of moles of an element is converted to kg and divided by the total mass of the solution. Two options are available for the volume and mass of solution that are used in converting to transport concentrations: (1) the volume and mass of solution are calculated by PHREEQC, or (2) the volume of solution is the product of porosity, saturation, and representative volume, and the mass of solution is volume times density as defined by SetDensity. Which option is used is determined by UseSolutionDensityVolume.

Return values
Unitsfor concentrations in transport.
See also
Concentrations2Utility, GetConcentrations, SetConcentrations, SetDensity, SetPorosity, SetRepresentativeVolume, SetUnitsSolution, UseSolutionDensityVolume.
C++ Example:
int units_solution = phreeqc_rm.GetUnitsSolution();
MPI:
Called by root and (or) workers.
int GetUnitsSSassemblage ( void  )
inline

Returns the input units for solid-solution assemblages. In PHREEQC input, solid solutions are defined by moles of each component (Mp). SetUnitsSSassemblage specifies how the number of moles in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

Return values
Inputunits for solid solutions (0, 1, or 2).
See also
SetUnitsSSassemblage, SetPorosity, SetRepresentativeVolume.
C++ Example:
int units_ss_exchange = phreeqc_rm.GetUnitsSSassemblage();
MPI:
Called by root and (or) workers.
int GetUnitsSurface ( void  )
inline

Returns the input units for surfaces. In PHREEQC input, surfaces are defined by moles of surface sites (Mp). SetUnitsSurface specifies how the number of moles of surface sites in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

Return values
Inputunits for solid surfaces (0, 1, or 2).
See also
SetUnitsSurface, SetPorosity, SetRepresentativeVolume.
C++ Example:
int units_surface = phreeqc_rm.GetUnitsSurface();
MPI:
Called by root and (or) workers.
const std::vector<IPhreeqcPhast *>& GetWorkers ( )
inline

Returns a reference to the vector of IPhreeqcPhast instances. IPhreeqcPhast inherits from IPhreeqc, and the vector can be interpreted as a vector of pointers to the worker, InitialPhreeqc, and Utility IPhreeqc instances. For OPENMP, there are nthreads workers, where nthreads is defined in the constructor (PhreeqcRM::PhreeqcRM). For MPI, there is a single worker. For OPENMP and MPI, there is one InitialPhreeqc and one Utility instance.

Return values
Vectorof IPhreeqcPhast instances.
See also
PhreeqcRM::PhreeqcRM.
C++ Example:
const std::vector < IPhreeqcPhast *> & w = phreeqc_rm.GetWorkers();
w[0]->AccumulateLine("Delete; -all");
int iphreeqc_result = w[0]->RunAccumulated();
MPI:
Called by root and (or) workers.
IRM_RESULT InitialPhreeqc2Concentrations ( std::vector< double > &  destination_c,
std::vector< int > &  boundary_solution1 
)

Fills a vector (destination_c) with concentrations from solutions in the InitialPhreeqc instance. The method is used to obtain concentrations for boundary conditions. If a negative value is used for a cell in boundary_solution1, then the highest numbered solution in the InitialPhreeqc instance will be used for that cell.

Parameters
destination_cVector to receive the concentrations.The dimension of destination_c is set to ncomps times n_boundary, where ncomps is the number of components returned from FindComponents or GetComponentCount, and n_boundary is the dimension of the vector boundary_solution1.
boundary_solution1Vector of solution index numbers that refer to solutions in the InitialPhreeqc instance. Size is n_boundary.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetComponentCount.
C++ Example:
std::vector bc_conc;
std::vector bc1;
int nbound = 1;
bc1.resize(nbound, 0);                      // solution 0 from InitialIPhreeqc instance
status = phreeqc_rm.InitialPhreeqc2Concentrations(bc_conc, bc1);
MPI:
Called by root.
IRM_RESULT InitialPhreeqc2Concentrations ( std::vector< double > &  destination_c,
std::vector< int > &  boundary_solution1,
std::vector< int > &  boundary_solution2,
std::vector< double > &  fraction1 
)

Fills a vector (destination_c) with concentrations from solutions in the InitialPhreeqc instance. The method is used to obtain concentrations for boundary conditions that are mixtures of solutions. If a negative value is used for a cell in boundary_solution1, then the highest numbered solution in the InitialPhreeqc instance will be used for that cell. Concentrations may be a mixture of two solutions, boundary_solution1 and boundary_solution2, with a mixing fraction for boundary_solution1 of fraction1 and mixing fraction for boundary_solution2 of (1 - fraction1). A negative value for boundary_solution2 implies no mixing, and the associated value for fraction1 is ignored.

Parameters
destination_cVector of concentrations extracted from the InitialPhreeqc instance. The dimension of destination_c is set to ncomps times n_boundary, where ncomps is the number of components returned from FindComponents or GetComponentCount, and n_boundary is the dimension of the vectors boundary_solution1, boundary_solution2, and fraction1.
boundary_solution1Vector of solution index numbers that refer to solutions in the InitialPhreeqc instance. Size is n_boundary.
boundary_solution2Vector of solution index numbers that that refer to solutions in the InitialPhreeqc instance and are defined to mix with boundary_solution1. Size is n_boundary.
fraction1Fraction of boundary_solution1 that mixes with (1 - fraction1) of boundary_solution2. Size is n_boundary.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetComponentCount.
C++ Example:
std::vector bc_conc, bc_f1;
std::vector bc1, bc2;
int nbound = 1;
bc1.resize(nbound, 0);                      // solution 0 from InitialIPhreeqc instance
bc2.resize(nbound, 1);                      // solution 1 from InitialIPhreeqc instance
bc_f1.resize(nbound, 0.4);                  // mixing fraction for bc1, result is 0.4/0.6 mix
status = phreeqc_rm.InitialPhreeqc2Concentrations(bc_conc, bc1, bc2, bc_f1);
MPI:
Called by root.
IRM_RESULT InitialPhreeqc2Module ( std::vector< int > &  initial_conditions1)

Transfer solutions and reactants from the InitialPhreeqc instance to the reaction-module workers. Initial_conditions1 is used to select initial conditions, including solutions and reactants, for each cell of the model, without mixing. Initial_conditions1 is dimensioned 7 times nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount). The dimension of 7 refers to solutions and reactants in the following order: (0) SOLUTIONS, (1) EQUILIBRIUM_PHASES, (2) EXCHANGE, (3) SURFACE, (4) GAS_PHASE, (5) SOLID_SOLUTIONS, and (6) KINETICS. The definition initial_solution1[3*nxyz + 99] = 2, indicates that cell 99 (0 based) contains the SURFACE definition (index 3) defined by SURFACE 2 in the InitialPhreeqc instance (created in the InitialPhreeqc instance either by RunFile or RunString).

Parameters
initial_conditions1Vector of solution and reactant index numbers that refer to definitions in the InitialPhreeqc instance. Size is 7 times nxyz. The order of definitions is given above. Negative values are ignored, resulting in no definition of that entity for that cell.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
InitialPhreeqcCell2Module.
C++ Example:
std::vector ic1;
ic1.resize(nxyz*7, -1);
for (int i = 0; i < nxyz; i++)
{
  ic1[i] = 1;              // Solution 1
  ic1[nxyz + i] = -1;      // Equilibrium phases none
  ic1[2*nxyz + i] = 1;     // Exchange 1
  ic1[3*nxyz + i] = -1;    // Surface none
  ic1[4*nxyz + i] = -1;    // Gas phase none
  ic1[5*nxyz + i] = -1;    // Solid solutions none
  ic1[6*nxyz + i] = -1;    // Kinetics none
}
status = phreeqc_rm.InitialPhreeqc2Module(ic1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT InitialPhreeqc2Module ( std::vector< int > &  initial_conditions1,
std::vector< int > &  initial_conditions2,
std::vector< double > &  fraction1 
)

Transfer solutions and reactants from the InitialPhreeqc instance to the reaction-module workers, possibly with mixing. In its simplest form, initial_conditions1 is used to select initial conditions, including solutions and reactants, for each cell of the model, without mixing. Initial_conditions1 is dimensioned 7 times nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount). The dimension of 7 refers to solutions and reactants in the following order: (0) SOLUTIONS, (1) EQUILIBRIUM_PHASES, (2) EXCHANGE, (3) SURFACE, (4) GAS_PHASE, (5) SOLID_SOLUTIONS, and (6) KINETICS. The definition initial_solution1[3*nxyz + 99] = 2, indicates that cell 99 (0 based) contains the SURFACE definition (index 3) defined by SURFACE 2 in the InitialPhreeqc instance (either by RunFile or RunString).

It is also possible to mix solutions and reactants to obtain the initial conditions for cells. For mixing, initials_conditions2 contains numbers for a second entity that mixes with the entity defined in initial_conditions1. Fraction1 contains the mixing fraction for initial_conditions1, whereas (1 - fraction1) is the mixing fraction for initial_conditions2. The definitions initial_solution1[3*nxyz + 99] = 2, initial_solution2[3*nxyz + 99] = 3, fraction1[3*nxyz + 99] = 0.25 indicates that cell 99 (0 based) contains a mixture of 0.25 SURFACE 2 and 0.75 SURFACE 3, where the surface compositions have been defined in the InitialPhreeqc instance. If the user number in initial_conditions2 is negative, no mixing occurs.

Parameters
initial_conditions1Vector of solution and reactant index numbers that refer to definitions in the InitialPhreeqc instance. Size is 7 times nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount). The order of definitions is given above. Negative values are ignored, resulting in no definition of that entity for that cell.
initial_conditions2Vector of solution and reactant index numbers that refer to definitions in the InitialPhreeqc instance. Nonnegative values of initial_conditions2 result in mixing with the entities defined in initial_conditions1. Negative values result in no mixing. Size is 7 times nxyz. The order of definitions is given above.
fraction1Fraction of initial_conditions1 that mixes with (1 - fraction1) of initial_conditions2. Size is 7 times nxyz. The order of definitions is given above.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
InitialPhreeqcCell2Module.
C++ Example:
std::vector ic1, ic2;
ic1.resize(nxyz*7, -1);
ic2.resize(nxyz*7, -1);
std::vector f1;
f1.resize(nxyz*7, 1.0);
for (int i = 0; i < nxyz; i++)
{
  ic1[i] = 1;              // Solution 1
  ic1[nxyz + i] = -1;      // Equilibrium phases none
  ic1[2*nxyz + i] = 1;     // Exchange 1
  ic1[3*nxyz + i] = -1;    // Surface none
  ic1[4*nxyz + i] = -1;    // Gas phase none
  ic1[5*nxyz + i] = -1;    // Solid solutions none
  ic1[6*nxyz + i] = -1;    // Kinetics none
}
status = phreeqc_rm.InitialPhreeqc2Module(ic1, ic2, f1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT InitialPhreeqc2SpeciesConcentrations ( std::vector< double > &  destination_c,
std::vector< int > &  boundary_solution1 
)

Fills a vector destination_c with aqueous species concentrations from solutions in the InitialPhreeqc instance. This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true. The method is used to obtain aqueous species concentrations for boundary conditions. If a negative value is used for a cell in boundary_solution1, then the highest numbered solution in the InitialPhreeqc instance will be used for that cell.

Parameters
destination_cVector of aqueous concentrations extracted from the InitialPhreeqc instance. The dimension of species_c is nspecies times n_boundary, where nspecies is the number of aqueous species returned from GetSpeciesCount, and n_boundary is the dimension of boundary_solution1.
boundary_solution1Vector of solution index numbers that refer to solutions in the InitialPhreeqc instance.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetSpeciesCount, SetSpeciesSaveOn.
C++ Example:
std::vector bc_conc, bc_f1;
std::vector bc1, bc2;
int nbound = 1;
bc1.resize(nbound, 0);                      // solution 0 from Initial IPhreeqc instance
status = phreeqc_rm.InitialPhreeqc2SpeciesConcentrations(bc_conc, bc1);
MPI:
Called by root.
IRM_RESULT InitialPhreeqc2SpeciesConcentrations ( std::vector< double > &  destination_c,
std::vector< int > &  boundary_solution1,
std::vector< int > &  boundary_solution2,
std::vector< double > &  fraction1 
)

Fills a vector destination_c with aqueous species concentrations from solutions in the InitialPhreeqc instance. This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true. The method is used to obtain aqueous species concentrations for boundary conditions. If a negative value is used for a cell in boundary_solution1, then the highest numbered solution in the InitialPhreeqc instance will be used for that cell. Concentrations may be a mixture of two solutions, boundary_solution1 and boundary_solution2, with a mixing fraction for boundary_solution1 of fraction1 and mixing fraction for boundary_solution2 of (1 - fraction1). A negative value for boundary_solution2 implies no mixing, and the associated value for fraction1 is ignored.

Parameters
destination_cVector of aqueous concentrations extracted from the InitialPhreeqc instance. The dimension of species_c is nspecies times n_boundary, where nspecies is the number of aqueous species returned from GetSpeciesCount, and n_boundary is the dimension of boundary_solution1.
boundary_solution1Vector of solution index numbers that refer to solutions in the InitialPhreeqc instance.
boundary_solution2Vector of solution index numbers that refer to solutions in the InitialPhreeqc instance and are defined to mix with boundary_solution1. Size is same as boundary_solution1.
fraction1Vector of fractions of boundary_solution1 that mix with (1 - fraction1) of boundary_solution2. Size is same as boundary_solution1.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetSpeciesCount, SetSpeciesSaveOn.
C++ Example:
std::vector bc_conc, bc_f1;
std::vector bc1, bc2;
int nbound = 1;
bc1.resize(nbound, 0);                      // solution 0 from Initial IPhreeqc instance
bc2.resize(nbound, -1);                     // no bc2 solution for mixing
bc_f1.resize(nbound, 1.0);                  // mixing fraction for bc1
status = phreeqc_rm.InitialPhreeqc2SpeciesConcentrations(bc_conc, bc1, bc2, bc_f1);
MPI:
Called by root.
IRM_RESULT InitialPhreeqcCell2Module ( int  n,
const std::vector< int > &  cell_numbers 
)

A cell numbered n in the InitialPhreeqc instance is selected to populate a series of transport cells. All reactants with the number n are transferred along with the solution. If MIX n exists, it is used for the definition of the solution. If n is negative, n is redefined to be the largest solution or MIX number in the InitialPhreeqc instance. All reactants for each cell in the list cell_numbers are removed before the cell definition is copied from the InitialPhreeqc instance to the workers.

Parameters
nNumber that refers to a solution or MIX and associated reactants in the InitialPhreeqc instance.
cell_numbersA vector of grid-cell numbers (user's grid-cell numbering system) that will be populated with cell n from the InitialPhreeqc instance.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
InitialPhreeqc2Module.
C++ Example:
std::vector module_cells;
module_cells.push_back(18);
module_cells.push_back(19);
status = phreeqc_rm.InitialPhreeqcCell2Module(-1, module_cells);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT LoadDatabase ( const std::string &  database)

Load a database for all IPhreeqc instances–workers, InitialPhreeqc, and Utility. All definitions of the reaction module are cleared (SOLUTION_SPECIES, PHASES, SOLUTIONs, etc.), and the database is read.

Parameters
databaseString containing the database name.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
C++ Example:
status = phreeqc_rm.LoadDatabase("phreeqc.dat");
MPI:
Called by root, workers must be in the loop of MpiWorker.
void LogMessage ( const std::string &  str)

Print a message to the log file.

Parameters
strString to be printed.
See also
OpenFiles, ErrorMessage, OutputMessage, ScreenMessage, WarningMessage.
C++ Example:
std::ostringstream strm;
strm << "Beginning transport calculation " <<   phreeqc_rm.GetTime() * phreeqc_rm.GetTimeConversion() << " days\n";
strm << "          Time step             " <<   phreeqc_rm.GetTimeStep() * phreeqc_rm.GetTimeConversion() << " days\n";
phreeqc_rm.LogMessage(strm.str());
MPI:
Called by root.
int MpiAbort ( )

MPI only. Calls MPI_Abort, which aborts MPI, and makes the reaction module unusable. Should be used only on encountering an unrecoverable error.

C++ Example:
int status = phreeqc_rm.MPI_Abort();
MPI:
Called by root or workers.
IRM_RESULT MpiWorker ( )

MPI only. Nonroot processes (processes with GetMpiMyself > 0) must call MpiWorker to be able to respond to messages from the root to accept data, perform calculations, and (or) return data within the reaction module. MpiWorker contains a loop that reads a message from root, performs a task, and waits for another message from root. SetConcentrations, RunCells, and GetConcentrations are examples of methods that send a message from root to get the workers to perform a task. The workers will respond to all methods that are designated "workers must be in the loop of MpiWorker" in the MPI section of the method documentation. The workers will continue to respond to messages from root until root calls MpiWorkerBreak.

(Advanced) The list of tasks that the workers perform can be extended by using SetMpiWorkerCallbackC. It is then possible to use the MPI processes to perform other developer-defined tasks, such as transport calculations, without exiting from the MpiWorker loop. Alternatively, root calls MpiWorkerBreak to allow the workers to continue past a call to MpiWorker. The workers perform developer-defined calculations, and then MpiWorker is called again to respond to requests from root to perform reaction-module tasks.

Return values
IRM_RESULT0 is success, negative is failure (See DecodeError). MpiWorker returns a value only when MpiWorkerBreak is called by root.
See also
MpiWorkerBreak, SetMpiWorkerCallbackC, SetMpiWorkerCallbackCookie.
C++ Example:
PhreeqcRM phreeqc_rm(nxyz, MPI_COMM_WORLD);
int mpi_myself;
if (MPI_Comm_rank(MPI_COMM_WORLD, &mpi_myself) != MPI_SUCCESS)
{
  exit(4);
}
if (mpi_myself > 0)
{
  phreeqc_rm.MpiWorker();
  return EXIT_SUCCESS;
}
MPI:
Called by all workers.
IRM_RESULT MpiWorkerBreak ( )

MPI only. This method is called by root to force nonroot processes (processes with GetMpiMyself > 0) to return from a call to MpiWorker. MpiWorker contains a loop that reads a message from root, performs a task, and waits for another message from root. The workers respond to all methods that are designated "workers must be in the loop of MpiWorker" in the MPI section of the method documentation. The workers will continue to respond to messages from root until root calls MpiWorkerBreak.

Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
MpiWorker, SetMpiWorkerCallbackC, SetMpiWorkerCallbackCookie (C only).
C++ Example:
status = phreeqc_rm.MpiWorkerBreak();
MPI:
Called by root.
IRM_RESULT OpenFiles ( void  )

Opens the output and log files. Files are named prefix.chem.txt and prefix.log.txt based on the prefix defined by SetFilePrefix.

Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetFilePrefix, GetFilePrefix, CloseFiles, ErrorMessage, LogMessage, OutputMessage, WarningMessage.
C++ Example:
status = phreeqc_rm.SetFilePrefix("Advect_cpp");
phreeqc_rm.OpenFiles();
MPI:
Called by root.
void OutputMessage ( const std::string &  str)

Print a message to the output file.

Parameters
strString to be printed.
See also
ErrorMessage, LogMessage, ScreenMessage, WarningMessage.
C++ Example:
std::ostringstream oss;
oss << "Database:                                         " << phreeqc_rm.GetDatabaseFileName().c_str() << "\n";
oss << "Number of threads:                                " << phreeqc_rm.GetThreadCount() << "\n";
oss << "Number of MPI processes:                          " << phreeqc_rm.GetMpiTasks() << "\n";
oss << "MPI task number:                                  " << phreeqc_rm.GetMpiMyself() << "\n";
oss << "File prefix:                                      " << phreeqc_rm.GetFilePrefix() << "\n";
oss << "Number of grid cells in the user's model:         " << phreeqc_rm.GetGridCellCount() << "\n";
oss << "Number of chemistry cells in the reaction module: " << phreeqc_rm.GetChemistryCellCount() << "\n";
oss << "Number of components for transport:               " << phreeqc_rm.GetComponentCount() << "\n";
oss << "Error handler mode:                               " << phreeqc_rm.GetErrorHandlerMode() << "\n";
phreeqc_rm.OutputMessage(oss.str());
MPI:
Called by root.
IRM_RESULT ReturnHandler ( IRM_RESULT  result,
const std::string &  e_string 
)

Process an IRM_RESULT return code. If the return code is nonnegative, no action is taken. If the return code is negative, the return code is decoded and printed as an error message along with the second argument (std::string). On an error, the method will return the same return code, throw an exception, or exit the program depending on the setting for SetErrorHandlerMode.

Parameters
resultReturn code to be processed.
e_stringError message to be printed in case of an error.
Return values
IRM_RESULTThe first argument to the method is returned.
See also
SetErrorHandlerMode.
C++ Example:
status = phreeqc_rm.ReturnHandler(irm_result, "Previous method failed.");
MPI:
Called by root or workers.
IRM_RESULT RunCells ( void  )

Runs a reaction step for all reaction cells in the reaction module. Normally, tranport concentrations are transferred to the reaction cells (SetConcentrations) before reaction calculations are run. The length of time over which kinetic reactions are integrated is set by SetTimeStep. Other properties that may need to be updated as a result of the transport calculations include porosity (SetPorosity), saturation (SetSaturation), temperature (SetTemperature), and pressure (SetPressure).

Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetConcentrations, SetPorosity, SetTemperature, SetPressure, SetSaturation, SetTimeStep.
C++ Example:
status = phreeqc_rm.SetSelectedOutputOn(print_selected_output_on);
status = phreeqc_rm.SetPrintChemistryOn(print_chemistry_on, false, false);
status = phreeqc_rm.SetPorosity(por);             // If porosity changes
status = phreeqc_rm.SetSaturation(sat);           // If saturation changes
status = phreeqc_rm.SetTemperature(temperature);  // If temperature changes
status = phreeqc_rm.SetPressure(pressure);        // If pressure changes
status = phreeqc_rm.SetConcentrations(c);         // Transported concentrations
status = phreeqc_rm.SetTimeStep(time_step);       // Time step for kinetic reactions
time = time + time_step;
status = phreeqc_rm.SetTime(time);
status = phreeqc_rm.RunCells();
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT RunFile ( bool  workers,
bool  initial_phreeqc,
bool  utility,
const std::string &  chemistry_name 
)

Run a PHREEQC input file. The first three arguments determine which IPhreeqc instances will run the file–the workers, the InitialPhreeqc instance, and (or) the Utility instance. Input files that modify the thermodynamic database should be run by all three sets of instances. Files with SELECTED_OUTPUT definitions that will be used during the time-stepping loop need to be run by the workers. Files that contain initial conditions or boundary conditions should be run by the InitialPhreeqc instance.

Parameters
workersTrue, the workers will run the file; False, the workers will not run the file.
initial_phreeqcTrue, the InitialPhreeqc instance will run the file; False, the InitialPhreeqc will not run the file.
utilityTrue, the Utility instance will run the file; False, the Utility instance will not run the file.
chemistry_nameName of the file to run.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
RunString.
C++ Example:
status = phreeqc_rm.RunFile(true, true, true, "advect.pqi");
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT RunString ( bool  workers,
bool  initial_phreeqc,
bool  utility,
const std::string &  input_string 
)

Run a PHREEQC input string. The first three arguments determine which IPhreeqc instances will run the string–the workers, the InitialPhreeqc instance, and (or) the Utility instance. Input strings that modify the thermodynamic database should be run by all three sets of instances. Strings with SELECTED_OUTPUT definitions that will be used during the time-stepping loop need to be run by the workers. Strings that contain initial conditions or boundary conditions should be run by the InitialPhreeqc instance.

Parameters
workersTrue, the workers will run the string; False, the workers will not run the string.
initial_phreeqcTrue, the InitialPhreeqc instance will run the string; False, the InitialPhreeqc will not run the string.
utilityTrue, the Utility instance will run the string; False, the Utility instance will not run the string.
input_stringString containing PHREEQC input.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
RunFile.
C++ Example:
std::string input = "DELETE; -all";
status = phreeqc_rm.RunString(true, false, true, input.c_str());
MPI:
Called by root, workers must be in the loop of MpiWorker.
void ScreenMessage ( const std::string &  str)

Print message to the screen.

Parameters
strString to be printed.
See also
ErrorMessage, LogMessage, OutputMessage, WarningMessage.
C++ Example:
std::ostringstream strm;
strm << "Beginning transport calculation "
     <<   phreeqc_rm.GetTime() * phreeqc_rm.GetTimeConversion()
     << " days\n";
strm << "          Time step             "
     <<   phreeqc_rm.GetTimeStep() * phreeqc_rm.GetTimeConversion()
     << " days\n";
phreeqc_rm.ScreenMessage(strm.str());
MPI:
Called by root and (or) workers.
IRM_RESULT SetComponentH2O ( bool  tf)

Select whether to include H2O in the component list. The concentrations of H and O must be known accurately (8 to 10 significant digits) for the numerical method of PHREEQC to produce accurate pH and pe values. Because most of the H and O are in the water species, it may be more robust (require less accuracy in transport) to transport the excess H and O (the H and O not in water) and water. The default setting (true) is to include water, excess H, and excess O as components. A setting of false will include total H and total O as components. SetComponentH2O must be called before FindComponents.

Parameters
tfTrue (default), excess H, excess O, and water are included in the component list; False, total H and O are included in the component list.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents.
C++ Example:
status = phreeqc_rm.SetComponentH2O(true);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetConcentrations ( const std::vector< double > &  c)

Use the vector of concentrations (c) to set the moles of components in each reaction cell. The volume of water in a cell is the product of porosity (SetPorosity), saturation (SetSaturation), and reference volume (SetRepresentativeVolume). The moles of each component are determined by the volume of water and per liter concentrations. If concentration units (SetUnitsSolution) are mass fraction, the density (as specified by SetDensity) is used to convert from mass fraction to per mass per liter.

Parameters
cVector of component concentrations. Size of vector is ncomps times nxyz, where ncomps is the number of components as determined by FindComponents or GetComponentCount and nxyz is the number of grid cells in the user's model (GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetDensity, SetPorosity, SetRepresentativeVolume, SetSaturation, SetUnitsSolution.
C++ Example:
std::vector c;
c.resize(nxyz * components.size());
...
AdvectCpp(c, bc_conc, ncomps, nxyz, nbound);
status = phreeqc_rm.SetPorosity(por);             // If porosity changes
status = phreeqc_rm.SetSaturation(sat);           // If saturation changes
status = phreeqc_rm.SetTemperature(temperature);  // If temperature changes
status = phreeqc_rm.SetPressure(pressure);        // If pressure changes
status = phreeqc_rm.SetConcentrations(c);         // Transported concentrations
status = phreeqc_rm.SetTimeStep(time_step);       // Time step for kinetic reactions
time = time + time_step;
status = phreeqc_rm.SetTime(time);
status = phreeqc_rm.RunCells();
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetCurrentSelectedOutputUserNumber ( int  n_user)

Select the current selected output by user number. The user may define multiple SELECTED_OUTPUT data blocks for the workers. A user number is specified for each data block. The value of the argument n_user selects which of the SELECTED_OUTPUT definitions will be used for selected-output operations.

Parameters
n_userUser number of the SELECTED_OUTPUT data block that is to be used.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetNthSelectedOutputUserNumber, GetSelectedOutput, GetSelectedOutputColumnCount, GetSelectedOutputCount, GetSelectedOutputRowCount, GetSelectedOutputHeading, SetSelectedOutputOn.
C++ Example:
for (int isel = 0; isel < phreeqc_rm.GetSelectedOutputCount(); isel++)
{
  int n_user = phreeqc_rm.GetNthSelectedOutputUserNumber(isel);
  status = phreeqc_rm.SetCurrentSelectedOutputUserNumber(n_user);
  std::cerr << "Selected output sequence number: " << isel << "\n";
  std::cerr << "Selected output user number:     " << n_user << "\n";
  std::vector so;
  status = phreeqc_rm.GetSelectedOutput(so);
  // Process results here
}
MPI:
Called by root.
IRM_RESULT SetDensity ( const std::vector< double > &  density)

Set the density for each reaction cell. These density values are used when converting from transported mass-fraction concentrations (SetUnitsSolution) to produce per liter concentrations during a call to SetConcentrations. They are also used when converting from reaction-cell concentrations to transport concentrations (GetConcentrations), if UseSolutionDensityVolume is set to false.

Parameters
densityVector of densities. Size of vector is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetConcentrations, GetGridCellCount, SetConcentrations, SetUnitsSolution, UseSolutionDensityVolume.
C++ Example:
std::vector initial_density;
initial_density.resize(nxyz, 1.0);
phreeqc_rm.SetDensity(initial_density);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetDumpFileName ( const std::string &  dump_name)

Set the name of the dump file. It is the name used by DumpModule.

Parameters
dump_nameName of dump file.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
DumpModule.
C++ Example:
status = phreeqc_rm.SetDumpFileName("advection_cpp.dmp");
bool dump_on = true;
bool append = false;
status = phreeqc_rm.DumpModule(dump_on, append);
MPI:
Called by root.
IRM_RESULT SetErrorHandlerMode ( int  mode)

Set the action to be taken when the reaction module encounters an error. Options are 0, return to calling program with an error return code (default); 1, throw an exception, in C++, the exception can be caught, for C and Fortran, the program will exit; or 2, attempt to exit gracefully.

Parameters
modeError handling mode: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
C++ Example:
PhreeqcRM phreeqc_rm(nxyz, nthreads);
IRM_RESULT status;
status = phreeqc_rm.SetErrorHandlerMode(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetErrorOn ( bool  tf)

Set the property that controls whether error messages are generated and displayed. Messages include PHREEQC "ERROR" messages, and any messages written with ErrorMessage.

Parameters
tfTrue, enable error messages; False, disable error messages. Default is true.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
ErrorMessage, ScreenMessage.
C++ Example:
status = phreeqc_rm.SetErrorOn(true);
MPI:
Called by root.
IRM_RESULT SetFilePrefix ( const std::string &  prefix)

Set the prefix for the output (prefix.chem.txt) and log (prefix.log.txt) files. These files are opened by OpenFiles.

Parameters
prefixPrefix used when opening the output and log files.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
OpenFiles, CloseFiles.
C++ Example:
status = phreeqc_rm.SetFilePrefix("Advect_cpp");
phreeqc_rm.OpenFiles();
MPI:
Called by root.
IRM_RESULT SetGasCompMoles ( const std::vector< double > &  gas_moles)

Transfer moles of gas components from the vector given in the argument list (gas_moles) to each reaction cell.

Parameters
gas_molesVector of moles of gas components. Dimension of the vector is set to ngas_comps times nxyz, where, ngas_comps is the result of GetGasComponentsCount, and nxyz is the number of user grid cells (GetGridCellCount). If the number of moles is set to a negative number, the gas component will not be defined for the GAS_PHASE of the reaction cell.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetGasComponentsCount, GetGasCompMoles, GetGasCompPressures, GetGasPhaseVolume, GetGasCompPhi, SetGasPhaseVolume.
C++ Example:
std::vector gas_moles;
gas_moles.resize(nxyz*ngas);
...
status = phreeqc_rm.SetGasCompMoles(gas_moles);
status = phreeqc_rm.RunCells();
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetGasPhaseVolume ( const std::vector< double > &  gas_volume)

Transfer volumes of gas phases from the vector given in the argument list (gas_volume) to each reaction cell. The gas-phase volume affects the gas-component pressures calculated for fixed-volume gas phases. If a gas-phase volume is defined with this methood for a GAS_PHASE in a cell, the gas phase is forced to be a fixed-volume gas phase.

Parameters
gas_volumeVector of volumes for each gas phase. Dimension of the vector is nxyz, where nxyz is the number of user grid cells (GetGridCellCount). If the volume is set to a negative number for a cell, the gas-phase volume for that cell is not changed.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetGasComponentsCount, GetGasCompMoles, GetGasCompPressures, GetGasCompPhi, GetGasPhaseVolume, SetGasCompMoles.
C++ Example:
std::vector gas_volume;
gas_volume.resize(nxyz, -1.0);
...
status = phreeqc_rm.SetGasPhaseVolume(gas_volume);
status = phreeqc_rm.RunCells();
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetMpiWorkerCallbackC ( int(*)(int *method, void *cookie)  fcn)

MPI and C/C++ only. Defines a callback function that allows additional tasks to be done by the workers. The method MpiWorker contains a loop, where the workers receive a message (an integer), run a function corresponding to that integer, and then wait for another message. SetMpiWorkerCallbackC allows C or C++ developers to add another function that responds to additional integer messages by calling developer-defined functions corresponding to those integers. MpiWorker calls the callback function when the message number is not one of the PhreeqcRM message numbers. Messages are unique integer numbers. PhreeqcRM uses integers in a range beginning at 0. It is suggested that developers use message numbers starting at 1000 or higher for their tasks. The callback function calls a developer-defined function specified by the message number and then returns to MpiWorker to wait for another message.

In C and C++, an additional pointer can be supplied to find the data necessary to do the task. A void pointer may be set with SetMpiWorkerCallbackCookie. This pointer is passed to the callback function through a void pointer argument in addition to the integer message argument. The pointer may be to a struct or class instance that provides a number of additional pointers to data. SetMpiWorkerCallbackCookie must be called by each worker before MpiWorker is called.

The motivation for this method is to allow the workers to perform other tasks, for instance, parallel transport calculations, within the structure of MpiWorker. The callback function can be used to allow the workers to receive data, perform transport calculations, and (or) send results, without leaving the loop of MpiWorker. Alternatively, it is possible for the workers to return from MpiWorker by a call to MpiWorkerBreak by root. The workers could then call subroutines to receive data, calculate transport, and send data, and then resume processing PhreeqcRM messages from root with another call to MpiWorker.

Parameters
fcnA function that returns an integer and has an integer argument and a void * argument.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
MpiWorker, MpiWorkerBreak, SetMpiWorkerCallbackCookie.
C++ Example:
Code executed by root:
// root calls a function that will involve the workers
int istatus = do_something(&comm);

Code executed by workers:
phreeqc_rm.SetMpiWorkerCallbackC(worker_tasks_cc);
phreeqc_rm.SetMpiWorkerCallbackCookie(&comm);
phreeqc_rm.MpiWorker();

Code executed by root and workers:
int do_something(void *cookie)
{
    int method_number = 1000;
    MP_TYPE *comm = (MP_TYPE *) cookie;
    int mpi_tasks, mpi_myself, worker_number;
    MPI_Comm_size(MPI_COMM_WORLD, &mpi_tasks);
    MPI_Comm_rank(MPI_COMM_WORLD, &mpi_myself);
    std::stringstream msg;
    if (mpi_myself == 0)
    {
        MPI_Bcast(&method_number, 1, MPI_INT, 0, *comm);
        fprintf(stderr, "I am root.\n");
        for (int i = 1; i < mpi_tasks; i++)
        {
            MPI_Status status;
            MPI_Recv(&worker_number, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &status);
            fprintf(stderr, "Recieved data from worker number %d.\n", worker_number);
        }
    }
    else
    {
        MPI_Send(&mpi_myself, 1, MPI_INT, 0, 0, *comm);
    }
    return 0;
}

Code called by workers from method MpiWorker:
int worker_tasks_cc(int *task_number, void * cookie)
{
    if (*task_number == 1000)
    {
        do_something(cookie);
    }
    return 0;
}
MPI:
Called by workers, before call to MpiWorker.
IRM_RESULT SetMpiWorkerCallbackCookie ( void *  cookie)

MPI and C/C++ only. Defines a void pointer that can be used by C and C++ functions called from the callback function (SetMpiWorkerCallbackC) to locate data for a task. The C callback function that is registered with SetMpiWorkerCallbackC has two arguments, an integer message to identify a task, and a void pointer. SetMpiWorkerCallbackCookie sets the value of the void pointer that is passed to the callback function. The void pointer may be a pointer to a struct of class instance that contains additonal pointers to data.

Parameters
cookieVoid pointer that can be used by subroutines called from the callback function to locate data needed to perform a task.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
MpiWorker, MpiWorkerBreak, SetMpiWorkerCallbackC.
C++ Example:
Code executed by root:
// root calls a function that will involve the workers
int istatus = do_something(&comm);

Code executed by workers:
phreeqc_rm.SetMpiWorkerCallbackC(worker_tasks_cc);
phreeqc_rm.SetMpiWorkerCallbackCookie(&comm);
phreeqc_rm.MpiWorker();

Code executed by root and workers:
int do_something(void *cookie)
{
    int method_number = 1000;
    MP_TYPE *comm = (MP_TYPE *) cookie;
    int mpi_tasks, mpi_myself, worker_number;
    MPI_Comm_size(MPI_COMM_WORLD, &mpi_tasks);
    MPI_Comm_rank(MPI_COMM_WORLD, &mpi_myself);
    std::stringstream msg;
    if (mpi_myself == 0)
    {
        MPI_Bcast(&method_number, 1, MPI_INT, 0, *comm);
        fprintf(stderr, "I am root.\n");
        for (int i = 1; i < mpi_tasks; i++)
        {
            MPI_Status status;
            MPI_Recv(&worker_number, 1, MPI_INT, i, 0, MPI_COMM_WORLD, &status);
            fprintf(stderr, "Recieved data from worker number %d.\n", worker_number);
        }
    }
    else
    {
        MPI_Send(&mpi_myself, 1, MPI_INT, 0, 0, *comm);
    }
    return 0;
}

Code called by workers from method MpiWorker:
int worker_tasks_cc(int *task_number, void * cookie)
{
    if (*task_number == 1000)
    {
        do_something(cookie);
    }
    return 0;
}
MPI:
Called by workers, before call to MpiWorker.
IRM_RESULT SetMpiWorkerCallbackFortran ( int(*)(int *method)  fcn)

MPI and Fortran only. Defines a callback function that allows additional tasks to be done by the workers. See documentation of PhreeqcRM for Fortran, method SetMpiWorkerCallback.

IRM_RESULT SetPartitionUZSolids ( bool  tf)

Sets the property for partitioning solids between the saturated and unsaturated parts of a partially saturated cell.

The option is intended to be used by saturated-only flow codes that allow a variable water table. The value has meaning only when saturations less than 1.0 are encountered. The partially saturated cells may have a small water-to-rock ratio that causes reactions to proceed differently relative to fully saturated cells. By setting SetPartitionUZSolids to true, the amounts of solids and gases are partioned according to the saturation. If a cell has a saturation of 0.5, then the water interacts with only half of the solids and gases; the other half is unreactive until the water table rises. As the saturation in a cell varies, solids and gases are transferred between the saturated and unsaturated (unreactive) reservoirs of the cell. Unsaturated-zone flow and transport codes will probably use the default (false), which assumes all gases and solids are reactive regardless of saturation.

Parameters
tfTrue, the fraction of solids and gases available for reaction is equal to the saturation; False (default), all solids and gases are reactive regardless of saturation.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetPartitionUZSolids.
C++ Example:
phreeqc_rm.SetPartitionUZSolids(false);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetPorosity ( const std::vector< double > &  por)

Set the porosity for each reaction cell. The volume of water in a reaction cell is the product of porosity, saturation (SetSaturation), and representative volume (SetRepresentativeVolume).

Parameters
porVector of porosities, unitless. Default is 0.1. Size of vector is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetSaturation, SetRepresentativeVolume, SetSaturation.
C++ Example:
std::vector por;
por.resize(nxyz, 0.2);
status = phreeqc_rm.SetPorosity(por);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetPressure ( const std::vector< double > &  p)

Set the pressure for each reaction cell. Pressure effects are considered only in three of the databases distributed with PhreeqcRM: phreeqc.dat, Amm.dat, and pitzer.dat.

Parameters
pVector of pressures, in atm. Size of vector is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetTemperature.
C++ Example:
std::vector pressure;
pressure.resize(nxyz, 2.0);
phreeqc_rm.SetPressure(pressure);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetPrintChemistryMask ( std::vector< int > &  cell_mask)

Enable or disable detailed output for each reaction cell. Printing for a reaction cell will occur only when the printing is enabled with SetPrintChemistryOn and the cell_mask value is 1.

Parameters
cell_maskVector of integers. Size of vector is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount). A value of 0 will disable printing detailed output for the cell; a value of 1 will enable printing detailed output for a cell.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetPrintChemistryOn.
C++ Example:
std::vector print_chemistry_mask;
print_chemistry_mask.resize(nxyz, 0);
for (int i = 0; i < nxyz/2; i++)
{
  print_chemistry_mask[i] = 1;
}
status = phreeqc_rm.SetPrintChemistryMask(print_chemistry_mask);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetPrintChemistryOn ( bool  workers,
bool  initial_phreeqc,
bool  utility 
)

Set property that enables or disables printing detailed output from reaction calculations to the output file for a set of cells defined by SetPrintChemistryMask. The detailed output prints all of the output typical of a PHREEQC reaction calculation, which includes solution descriptions and the compositions of all other reactants. The output can be several hundred lines per cell, which can lead to a very large output file (prefix.chem.txt, OpenFiles). For the worker instances, the output can be limited to a set of cells (SetPrintChemistryMask) and, in general, the amount of information printed can be limited by use of options in the PRINT data block of PHREEQC (applied by using RunFile or RunString). Printing the detailed output for the workers is generally used only for debugging, and PhreeqcRM will run significantly faster when printing detailed output for the workers is disabled.

Parameters
workersTrue, enable detailed printing in the worker instances; False, disable detailed printing in the worker instances.
initial_phreeqcTrue, enable detailed printing in the InitialPhreeqc instance; False, disable detailed printing in the InitialPhreeqc instance.
utilityTrue, enable detailed printing in the Utility instance; False, disable detailed printing in the Utility instance.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
OpenFiles, RunFile, RunString, SetPrintChemistryMask.
C++ Example:
status = phreeqc_rm.SetPrintChemistryOn(false, true, false);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetRebalanceByCell ( bool  tf)

Set the load-balancing algorithm. PhreeqcRM attempts to rebalance the load of each thread or process such that each thread or process takes the same amount of time to run its part of a RunCells calculation. Two algorithms are available; one uses individual times for each cell and accounts for cells that were not run because saturation was zero (default), and the other assigns an average time to all cells. The methods are similar, but limited testing indicates the default method performs better.

Parameters
tfTrue, indicates individual cell times are used in rebalancing (default); False, indicates average times are used in rebalancing.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetRebalanceFraction.
C++ Example:
status = phreeqc_rm.SetRebalanceByCell(true);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetRebalanceFraction ( double  f)

Sets the fraction of cells that are transferred among threads or processes when rebalancing. PhreeqcRM attempts to rebalance the load of each thread or process such that each thread or process takes the same amount of time to run its part of a RunCells calculation. The rebalancing transfers cell calculations among threads or processes to try to achieve an optimum balance. SetRebalanceFraction adjusts the calculated optimum number of cell transfers by a fraction from 0 to 1.0 to determine the actual number of cell transfers. A value of zero eliminates load rebalancing. A value less than 1.0 is suggested to slow the approach to the optimum cell distribution and avoid possible oscillations when too many cells are transferred at one iteration, requiring reverse transfers at the next iteration. Default is 0.5.

Parameters
fFraction from 0.0 to 1.0.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetRebalanceByCell.
C++ Example:
status = phreeqc_rm.SetRebalanceFraction(0.5);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetRepresentativeVolume ( const std::vector< double > &  rv)

Set the representative volume of each reaction cell. By default the representative volume of each reaction cell is 1 liter. The volume of water in a reaction cell is determined by the product of the representative volume, the porosity (SetPorosity), and the saturation (SetSaturation). The numerical method of PHREEQC is more robust if the water volume for a reaction cell is within a couple orders of magnitude of 1.0. Small water volumes caused by small porosities and (or) small saturations (and (or) small representative volumes) may cause non-convergence of the numerical method. In these cases, a larger representative volume may help. Note that increasing the representative volume also increases the number of moles of the reactants in the reaction cell (minerals, surfaces, exchangers, and others), which are defined as moles per representative volume.

Parameters
rvVector of representative volumes, in liters. Default is 1.0 liter. Size of array is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetPorosity, SetSaturation.
C++ Example:
std::vector rv;
rv.resize(nxyz, 2.0);
status = phreeqc_rm.SetRepresentativeVolume(rv);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetSaturation ( const std::vector< double > &  sat)

Set the saturation of each reaction cell. Saturation is a fraction ranging from 0 to 1. The volume of water in a cell is the product of porosity (SetPorosity), saturation (SetSaturation), and representative volume (SetRepresentativeVolume). As a result of a reaction calculation, solution properties (density and volume) will change; the databases phreeqc.dat, Amm.dat, and pitzer.dat have the molar volume data to calculate these changes. The methods GetDensity, GetSolutionVolume, and GetSaturation can be used to account for these changes in the succeeding transport calculation. SetRepresentativeVolume should be called before initial conditions are defined for the reaction cells.

Parameters
satVector of saturations, unitless. Default 1.0. Size of vector is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetDensity, GetSaturation, GetSolutionVolume, SetPorosity, SetRepresentativeVolume.
C++ Example:
std::vector sat;
sat.resize(nxyz, 1.0);
status = phreeqc_rm.SetSaturation(sat);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetScreenOn ( bool  tf)

Set the property that controls whether messages are written to the screen. Messages include information about rebalancing during RunCells, and any messages written with ScreenMessage.

Parameters
tfTrue, enable screen messages; False, disable screen messages. Default is true.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
RunCells, ScreenMessage.
C++ Example:
status = phreeqc_rm.SetScreenOn(true);
MPI:
Called by root.
IRM_RESULT SetSelectedOutputOn ( bool  tf)

Set the property that controls whether selected-output results are available to be retrieved with GetSelectedOutput. True indicates that selected-output results will be accumulated during RunCells and can be retrieved with GetSelectedOutput; False indicates that selected-output results will not be accumulated during RunCells.

Parameters
tfTrue, enable selected output; False, disable selected output.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetSelectedOutput, SetPrintChemistryOn.
C++ Example:
status = phreeqc_rm.SetSelectedOutputOn(true);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetSpeciesSaveOn ( bool  save_on)

Sets the value of the species-save property. This method enables or disables use of PhreeqcRM with multicomponent-diffusion transport calculations. By default, concentrations of aqueous species are not saved. Setting the species-save property to true allows aqueous species concentrations to be retrieved with GetSpeciesConcentrations, and solution compositions to be set with SpeciesConcentrations2Module. SetSpeciesSaveOn must be called before calls to FindComponents.

Parameters
save_onTrue indicates species concentrations are saved; False indicates species concentrations are not saved.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SpeciesConcentrations2Module.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
MPI:
Called by root and (or) workers.
IRM_RESULT SetTemperature ( const std::vector< double > &  t)

Set the temperature for each reaction cell. If SetTemperature is not called, worker solutions will have temperatures as defined by initial conditions (InitialPhreeqc2Module and InitialPhreeqcCell2Module).

Parameters
tVector of temperatures, in degrees C. Size of vector is nxyz, where nxyz is the number of grid cells in the user's model (GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetPressure, InitialPhreeqc2Module, InitialPhreeqcCell2Module, SetPressure, GetTemperature.
C++ Example:
std::vector temperature;
temperature.resize(nxyz, 20.0);
phreeqc_rm.SetTemperature(temperature);
MPI:
Called by root and (or) workers.
IRM_RESULT SetTime ( double  time)

Set current simulation time for the reaction module.

Parameters
timeCurrent simulation time, in seconds.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetTimeStep, SetTimeConversion.
C++ Example:
time += time_step;
status = phreeqc_rm.SetTime(time);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetTimeConversion ( double  conv_factor)

Set a factor to convert from seconds to user time units. Factor times seconds produces user time units.

Parameters
conv_factorFactor to convert seconds to user time units.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetTime, SetTimeStep.
C++ Example:
double time_conversion = 1.0 / 86400;
status = phreeqc_rm.SetTimeConversion(time_conversion);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetTimeStep ( double  time_step)

Set current time step for the reaction module. This is the length of time over which kinetic reactions are integrated.

Parameters
time_stepTime step, in seconds.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetTime, SetTimeConversion.
C++ Example:
time_step = 86400.;
status = phreeqc_rm.SetTimeStep(time_step);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetUnitsExchange ( int  option)

Sets input units for exchangers. In PHREEQC input, exchangers are defined by moles of exchange sites (Mp). SetUnitsExchange specifies how the number of moles of exchange sites in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

If a single EXCHANGE definition is used for cells with different initial porosity, the three options scale quite differently. For option 0, the number of moles of exchangers will be the same regardless of porosity. For option 1, the number of moles of exchangers will be vary directly with porosity and inversely with rock volume. For option 2, the number of moles of exchangers will vary directly with rock volume and inversely with porosity.

Parameters
optionUnits option for exchangers: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetUnitsExchange, InitialPhreeqc2Module, InitialPhreeqcCell2Module, SetPorosity, SetRepresentativeVolume.
C++ Example:
status = phreeqc_rm.SetUnitsExchange(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetUnitsGasPhase ( int  option)

Set input units for gas phases. In PHREEQC input, gas phases are defined by moles of component gases (Mp). SetUnitsGasPhase specifies how the number of moles of component gases in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

If a single GAS_PHASE definition is used for cells with different initial porosity, the three options scale quite differently. For option 0, the number of moles of a gas component will be the same regardless of porosity. For option 1, the number of moles of a gas component will be vary directly with porosity and inversely with rock volume. For option 2, the number of moles of a gas component will vary directly with rock volume and inversely with porosity.

Parameters
optionUnits option for gas phases: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetUnitsGasPhase, InitialPhreeqc2Module, InitialPhreeqcCell2Module, SetPorosity, SetRepresentativeVolume.
C++ Example:
status = phreeqc_rm.SetUnitsGasPhase(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetUnitsKinetics ( int  option)

Set input units for kinetic reactants.

In PHREEQC input, kinetics are defined by moles of kinetic reactants (Mp). SetUnitsKinetics specifies how the number of moles of kinetic reactants in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

If a single KINETICS definition is used for cells with different initial porosity, the three options scale quite differently. For option 0, the number of moles of kinetic reactants will be the same regardless of porosity. For option 1, the number of moles of kinetic reactants will be vary directly with porosity and inversely with rock volume. For option 2, the number of moles of kinetic reactants will vary directly with rock volume and inversely with porosity.

Note that the volume of water in a cell in the reaction module is equal to the product of porosity (SetPorosity), the saturation (SetSaturation), and representative volume (SetRepresentativeVolume), which is usually less than 1 liter. It is important to write the RATES definitions for homogeneous (aqueous) kinetic reactions to account for the current volume of water, often by calculating the rate of reaction per liter of water and multiplying by the volume of water (Basic function SOLN_VOL).

Rates that depend on surface area of solids, are not dependent on the volume of water. However, it is important to get the correct surface area for the kinetic reaction. To scale the surface area with the number of moles, the specific area (m^2 per mole of reactant) can be defined as a parameter (KINETICS; -parm), which is multiplied by the number of moles of reactant (Basic function M) in RATES to obtain the surface area.

Parameters
optionUnits option for kinetic reactants: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetUnitsKinetics, InitialPhreeqc2Module, InitialPhreeqcCell2Module, SetPorosity, SetRepresentativeVolume, SetSaturation.
C++ Example:
status = phreeqc_rm.SetUnitsKinetics(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetUnitsPPassemblage ( int  option)

Set input units for pure phase assemblages (equilibrium phases). In PHREEQC input, equilibrium phases are defined by moles of each phase (Mp). SetUnitsPPassemblage specifies how the number of moles of phases in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

If a single EQUILIBRIUM_PHASES definition is used for cells with different initial porosity, the three options scale quite differently. For option 0, the number of moles of a mineral will be the same regardless of porosity. For option 1, the number of moles of a mineral will be vary directly with porosity and inversely with rock volume. For option 2, the number of moles of a mineral will vary directly with rock volume and inversely with porosity.

Parameters
optionUnits option for equilibrium phases: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetUnitsPPassemblage, InitialPhreeqc2Module, InitialPhreeqcCell2Module, SetPorosity, SetRepresentativeVolume.
C++ Example:
status = phreeqc_rm.SetUnitsPPassemblage(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetUnitsSolution ( int  option)

Solution concentration units used by the transport model. Options are 1, mg/L; 2 mol/L; or 3, mass fraction, kg/kgs. PHREEQC defines solutions by the number of moles of each element in the solution.

To convert from mg/L to moles of element in the representative volume of a reaction cell, mg/L is converted to mol/L and multiplied by the solution volume, which is the product of porosity (SetPorosity), saturation (SetSaturation), and representative volume (SetRepresentativeVolume). To convert from mol/L to moles of element in the representative volume of a reaction cell, mol/L is multiplied by the solution volume. To convert from mass fraction to moles of element in the representative volume of a reaction cell, kg/kgs is converted to mol/kgs, multiplied by density (SetDensity) and multiplied by the solution volume.

To convert from moles of element in the representative volume of a reaction cell to mg/L, the number of moles of an element is divided by the solution volume resulting in mol/L, and then converted to mg/L. To convert from moles of element in a cell to mol/L, the number of moles of an element is divided by the solution volume resulting in mol/L. To convert from moles of element in a cell to mass fraction, the number of moles of an element is converted to kg and divided by the total mass of the solution. Two options are available for the volume and mass of solution that are used in converting to transport concentrations: (1) the volume and mass of solution are calculated by PHREEQC, or (2) the volume of solution is the product of porosity (SetPorosity), saturation (SetSaturation), and representative volume (SetRepresentativeVolume), and the mass of solution is volume times density as defined by SetDensity. Which option is used is determined by UseSolutionDensityVolume.

Parameters
optionUnits option for solutions: 1, 2, or 3, default is 1, mg/L.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
SetDensity, SetPorosity, SetRepresentativeVolume, SetSaturation, UseSolutionDensityVolume.
C++ Example:
status = phreeqc_rm.SetUnitsSolution(2);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetUnitsSSassemblage ( int  option)

Set input units for solid-solution assemblages. In PHREEQC, solid solutions are defined by moles of each component (Mp). SetUnitsSSassemblage specifies how the number of moles of solid-solution components in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-@ P)*RV.

If a single SOLID_SOLUTION definition is used for cells with different initial porosity, the three options scale quite differently. For option 0, the number of moles of a solid-solution component will be the same regardless of porosity. For option 1, the number of moles of a solid-solution component will be vary directly with porosity and inversely with rock volume. For option 2, the number of moles of a solid-solution component will vary directly with rock volume and inversely with porosity.

Parameters
optionUnits option for solid solutions: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetUnitsSSassemblage, InitialPhreeqc2Module, InitialPhreeqcCell2Module, SetPorosity, SetRepresentativeVolume.
C++ Example:
status = phreeqc_rm.SetUnitsSSassemblage(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SetUnitsSurface ( int  option)

Set input units for surfaces. In PHREEQC input, surfaces are defined by moles of surface sites (Mp). SetUnitsSurface specifies how the number of moles of surface sites in a reaction cell (Mc) is calculated from the input value (Mp).

Options are 0, Mp is mol/L of RV (default), Mc = Mp*RV, where RV is the representative volume (SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (SetPorosity); or 2, Mp is mol/L of rock in the RV, Mc = Mp*(1-P)*RV.

If a single SURFACE definition is used for cells with different initial porosity, the three options scale quite differently. For option 0, the number of moles of surface sites will be the same regardless of porosity. For option 1, the number of moles of surface sites will be vary directly with porosity and inversely with rock volume. For option 2, the number of moles of surface sites will vary directly with rock volume and inversely with porosity.

Parameters
optionUnits option for surfaces: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
GetUnitsSurface, InitialPhreeqc2Module, InitialPhreeqcCell2Module, SetPorosity, SetRepresentativeVolume.
C++ Example:
status = phreeqc_rm.SetUnitsSurface(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT SpeciesConcentrations2Module ( std::vector< double > &  species_conc)

Set solution concentrations in the reaction cells based on the vector of aqueous species concentrations (species_conc). This method is intended for use with multicomponent-diffusion transport calculations, and SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by FindComponents and includes all aqueous species that can be made from the set of components. The method determines the total concentration of a component by summing the molarities of the individual species times the stoichiometric coefficient of the element in each species. Solution compositions in the reaction cells are updated with these component concentrations.

Parameters
species_concVector of aqueous species concentrations. Dimension of the array is nspecies times nxyz, where nspecies is the number of aqueous species (GetSpeciesCount), and nxyz is the number of user grid cells (GetGridCellCount). Concentrations are moles per liter.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
FindComponents, GetSpeciesConcentrations, GetSpeciesCount, GetSpeciesD25, GetSpeciesLog10Gammas, GetSpeciesLog10Molalities, GetSpeciesNames, GetSpeciesSaveOn, GetSpeciesStoichiometry, GetSpeciesZ, SetSpeciesSaveOn.
C++ Example:
status = phreeqc_rm.SetSpeciesSaveOn(true);
int ncomps = phreeqc_rm.FindComponents();
int nspecies = phreeqc_rm.GetSpeciesCount();
std::vector c;
status = phreeqc_rm.GetSpeciesConcentrations(c);
...
SpeciesAdvectCpp(c, bc_conc, nspecies, nxyz, nbound);
status = phreeqc_rm.SpeciesConcentrations2Module(c);
status = phreeqc_rm.RunCells();
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT StateApply ( int  istate)

Reset the state of the module to a state previously saved with StateSave. The chemistry of all model cells are reset, including SOLUTIONs, EQUILIBRIUM_PHASES, EXCHANGEs, GAS_PHASEs, KINETICS, SOLID_SOLUTIONs, and SURFACEs. MIXes, REACTIONs, REACTION_PRESSUREs, and REACTION_TEMPERATUREs will be reset for each cell, if they were defined in the worker IPhreeqc instances at the time the state was saved. The distribution of cells among the workers and the chemistry of fully or partially unsaturated cells are also reset to the saved state. The state to be applied is identified by an integer.

Parameters
istateInteger identifying the state that is to be applied.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
StateSave and StateDelete.
C++ Example:
status = phreeqc_rm.StateSave(1);
...
status = phreeqc_rm.StateApply(1);
status = phreeqc_rm.StateDelete(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT StateDelete ( int  istate)

Delete a state previously saved with StateSave.

Parameters
istateInteger identifying the state that is to be deleted.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
StateSave and StateApply.
C++ Example:
status = phreeqc_rm.StateSave(1);
...
status = phreeqc_rm.StateApply(1);
status = phreeqc_rm.StateDelete(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
IRM_RESULT StateSave ( int  istate)

Save the state of the chemistry in all model cells, including SOLUTIONs, EQUILIBRIUM_PHASES, EXCHANGEs, GAS_PHASEs, KINETICS, SOLID_SOLUTIONs, and SURFACEs. Although not generally used, MIXes, REACTIONs, REACTION_PRESSUREs, and REACTION_TEMPERATUREs will be saved for each cell, if they have been defined in the worker IPhreeqc instances. The distribution of cells among the workers and the chemistry of fully or partially unsaturated cells are also saved. The state is saved in memory; use DumpModule to save the state to file. PhreeqcRM can be reset to this state by using StateApply. A state is identified by an integer, and multiple states can be saved.

Parameters
istateInteger identifying the state that is saved.
Return values
IRM_RESULT0 is success, negative is failure (See DecodeError).
See also
DumpModule, StateApply, and StateDelete.
C++ Example:
status = phreeqc_rm.StateSave(1);
...
status = phreeqc_rm.StateApply(1);
status = phreeqc_rm.StateDelete(1);
MPI:
Called by root, workers must be in the loop of MpiWorker.
void UseSolutionDensityVolume ( bool  tf)

Determines the volume and density to use when converting from the reaction-cell concentrations to transport concentrations (GetConcentrations). Two options are available to convert concentration units: (1) the density and solution volume calculated by PHREEQC are used, or (2) the specified density (SetDensity) and solution volume are determined by the product of saturation (SetSaturation), porosity (SetPorosity), and representative volume (SetRepresentativeVolume). Transport models that consider density-dependent flow will probably use the PHREEQC-calculated density and solution volume (default), whereas transport models that assume constant-density flow will probably use specified values of density and solution volume. Only the following databases distributed with PhreeqcRM have molar-volume information needed to accurately calculate density and solution volume: phreeqc.dat, Amm.dat, and pitzer.dat. Density is only used when converting to or from transport units of mass fraction.

Parameters
tfTrue indicates that the solution density and volume as calculated by PHREEQC will be used to calculate concentrations. False indicates that the solution density set by SetDensity and the volume determined by the product of SetSaturation, SetPorosity, and SetRepresentativeVolume, will be used to calculate concentrations retrieved by GetConcentrations.
See also
GetConcentrations, SetDensity, SetPorosity, SetRepresentativeVolume, SetSaturation.
C++ Example:
phreeqc_rm.UseSolutionDensityVolume(false);
MPI:
Called by root, workers must be in the loop of MpiWorker.
void WarningMessage ( const std::string &  warnstr)

Print a warning message to the screen and the log file.

Parameters
warnstrString to be printed.
See also
OpenFiles, ErrorMessage, LogMessage, OutputMessage, ScreenMessage.
C++ Example:
phreeqc_rm.WarningMessage("Parameter is out of range, using default");
MPI:
Called by root and (or) workers; only root writes to the log file.

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