PhreeqcRM
PhreeqcRM C Reference
Include dependency graph for RM_interface_C.h:

Go to the source code of this file.

Macros

#define RM_INTERFACE_C_H
 

Functions

IRM_DLL_EXPORT int RM_BmiCreate (int nxyz, int nthreads)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiDestroy (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiAddOutputVars (int id, char *option, char *def)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiFinalize (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetComponentName (int id, char *component_name, int l)
 
IRM_DLL_EXPORT double RM_BmiGetCurrentTime (int id)
 
IRM_DLL_EXPORT double RM_BmiGetEndTime (int id)
 
IRM_DLL_EXPORT int RM_BmiGetGridRank (int id, int grid)
 
IRM_DLL_EXPORT int RM_BmiGetGridSize (int id, int grid)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetGridType (int id, int grid, char *str, int l)
 
IRM_DLL_EXPORT int RM_BmiGetInputItemCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetInputVarName (int id, int i, char *name, int l)
 
IRM_DLL_EXPORT int RM_BmiGetOutputItemCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetOutputVarName (int id, int i, char *name, int l)
 
IRM_DLL_EXPORT int RM_BmiGetPointableItemCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetPointableVarName (int id, int i, char *name, int l)
 
IRM_DLL_EXPORT double RM_BmiGetStartTime (int id)
 
IRM_DLL_EXPORT double RM_BmiGetTime (int id)
 
IRM_DLL_EXPORT double RM_BmiGetTimeStep (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetTimeUnits (int id, char *units, int l)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetValueInt (int id, char *var, int *dest)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetValueDouble (int id, char *var, double *dest)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetValueChar (int id, char *var, char *dest, int l)
 
IRM_DLL_EXPORT void * RM_BmiGetValuePtr (int id, char *var)
 
IRM_DLL_EXPORT int RM_BmiGetVarGrid (int id, char *var)
 
IRM_DLL_EXPORT int RM_BmiGetVarItemsize (int id, char *name)
 
IRM_DLL_EXPORT int RM_BmiGetVarNbytes (int id, char *name)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetVarType (int id, char *name, char *vtype, int l)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiGetVarUnits (int id, char *name, char *units, int l)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiInitialize (int id, char *config_file)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiSetValueChar (int id, char *name, const char *src)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiSetValueDouble (int id, char *name, double src)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiSetValueDoubleArray (int id, char *name, double *src)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiSetValueInt (int id, char *name, int src)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiUpdate (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_BmiUpdateUntil (int id, double end_time)
 
IRM_DLL_EXPORT void RM_BmiGetValueAtIndices (int id, char *name, void *dest, int *inds, int count)
 
IRM_DLL_EXPORT void RM_BmiSetValueAtIndices (int id, char *name, int *inds, int count, void *src)
 
IRM_DLL_EXPORT void RM_BmiGetGridShape (int id, const int grid, int *shape)
 
IRM_DLL_EXPORT void RM_BmiGetGridSpacing (int id, const int grid, double *spacing)
 
IRM_DLL_EXPORT void RM_BmiGetGridOrigin (int id, const int grid, double *origin)
 
IRM_DLL_EXPORT void RM_BmiGetGridX (int id, const int grid, double *x)
 
IRM_DLL_EXPORT void RM_BmiGetGridY (int id, const int grid, double *y)
 
IRM_DLL_EXPORT void RM_BmiGetGridZ (int id, const int grid, double *z)
 
IRM_DLL_EXPORT int RM_BmiGetGridNodeCount (int id, const int grid)
 
IRM_DLL_EXPORT int RM_BmiGetGridEdgeCount (int id, const int grid)
 
IRM_DLL_EXPORT int RM_BmiGetGridFaceCount (int id, const int grid)
 
IRM_DLL_EXPORT void RM_BmiGetGridEdgeNodes (int id, const int grid, int *edge_nodes)
 
IRM_DLL_EXPORT void RM_BmiGetGridFaceEdges (int id, const int grid, int *face_edges)
 
IRM_DLL_EXPORT void RM_BmiGetGridFaceNodes (int id, const int grid, int *face_nodes)
 
IRM_DLL_EXPORT void RM_BmiGetGridNodesPerFace (int id, const int grid, int *nodes_per_face)
 
IRM_DLL_EXPORT IRM_RESULT RM_Abort (int id, int result, const char *err_str)
 
IRM_DLL_EXPORT IRM_RESULT RM_CloseFiles (int id)
 
IRM_DLL_EXPORT int RM_Concentrations2Utility (int id, double *c, int n, double *tc, double *p_atm)
 
IRM_DLL_EXPORT int RM_Create (int nxyz, int nthreads)
 
IRM_DLL_EXPORT IRM_RESULT RM_CreateMapping (int id, int *grid2chem)
 
IRM_DLL_EXPORT IRM_RESULT RM_DecodeError (int id, int e)
 
IRM_DLL_EXPORT IRM_RESULT RM_Destroy (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_DumpModule (int id, int dump_on, int append)
 
IRM_DLL_EXPORT IRM_RESULT RM_ErrorMessage (int id, const char *errstr)
 
IRM_DLL_EXPORT int RM_FindComponents (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetBackwardMapping (int id, int n, int *list, int *size)
 
IRM_DLL_EXPORT int RM_GetChemistryCellCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetComponent (int id, int num, char *chem_name, int l)
 
IRM_DLL_EXPORT int RM_GetComponentCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetConcentrations (int id, double *c)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetIthConcentration (int id, int i, double *c)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetIthSpeciesConcentration (int id, int i, double *c)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetIthConcentration (int id, int i, double *c)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetIthSpeciesConcentration (int id, int i, double *c)
 
IRM_DLL_EXPORT int RM_GetCurrentSelectedOutputUserNumber (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetDensityCalculated (int id, double *density)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetDensity (int id, double *density)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetEndCell (int id, int *ec)
 
IRM_DLL_EXPORT int RM_GetEquilibriumPhasesCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetEquilibriumPhasesName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetErrorString (int id, char *errstr, int l)
 
IRM_DLL_EXPORT int RM_GetErrorStringLength (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetExchangeName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT int RM_GetExchangeSpeciesCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetExchangeSpeciesName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetFilePrefix (int id, char *prefix, int l)
 
IRM_DLL_EXPORT int RM_GetGasComponentsCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetGasComponentsName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetGasCompMoles (int id, double *gas_moles)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetGasCompPressures (int id, double *gas_pressure)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetGasCompPhi (int id, double *gas_phi)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetGasPhaseVolume (int id, double *gas_volume)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetGfw (int id, double *gfw)
 
IRM_DLL_EXPORT int RM_GetGridCellCount (int id)
 
IRM_DLL_EXPORT int RM_GetGridCellCountYAML (const char *config_file)
 
IRM_DLL_EXPORT int RM_GetIPhreeqcId (int id, int i)
 
IRM_DLL_EXPORT int RM_GetKineticReactionsCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetKineticReactionsName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT int RM_GetMpiMyself (int id)
 
IRM_DLL_EXPORT int RM_GetMpiTasks (int id)
 
IRM_DLL_EXPORT int RM_GetNthSelectedOutputUserNumber (int id, int n)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetPorosity (int id, double *porosity)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetPressure (int id, double *pressure)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSaturationCalculated (int id, double *sat_calc)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSaturation (int id, double *sat_calc)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSelectedOutput (int id, double *so)
 
IRM_DLL_EXPORT int RM_GetSelectedOutputColumnCount (int id)
 
IRM_DLL_EXPORT int RM_GetSelectedOutputCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSelectedOutputHeading (int id, int icol, char *heading, int length)
 
IRM_DLL_EXPORT int RM_GetSelectedOutputRowCount (int id)
 
IRM_DLL_EXPORT int RM_GetSICount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSIName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT int RM_GetSolidSolutionComponentsCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSolidSolutionComponentsName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSolidSolutionName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSolutionVolume (int id, double *vol)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesConcentrations (int id, double *species_conc)
 
IRM_DLL_EXPORT int RM_GetSpeciesCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesD25 (int id, double *diffc)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesLog10Gammas (int id, double *species_log10gammas)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesLog10Molalities (int id, double *species_log10molalities)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesName (int id, int i, char *name, int length)
 
IRM_DLL_EXPORT int RM_GetSpeciesSaveOn (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesZ (int id, double *z)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetStartCell (int id, int *sc)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetTemperature (int id, double *temperature)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSurfaceName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT int RM_GetSurfaceSpeciesCount (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSurfaceSpeciesName (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetSurfaceType (int id, int num, char *name, int l1)
 
IRM_DLL_EXPORT int RM_GetThreadCount (int id)
 
IRM_DLL_EXPORT double RM_GetTime (int id)
 
IRM_DLL_EXPORT double RM_GetTimeConversion (int id)
 
IRM_DLL_EXPORT double RM_GetTimeStep (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_GetViscosity (int id, double *viscosity)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitializeYAML (int id, const char *yamlfile)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialPhreeqc2Concentrations (int id, double *c, int n_boundary, int *boundary_solution1, int *boundary_solution2, double *fraction1)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialSolutions2Module (int id, int *solutions)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialEquilibriumPhases2Module (int id, int *equilibrium_phases)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialExchanges2Module (int id, int *exchanges)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialSurfaces2Module (int id, int *surfaces)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialGasPhases2Module (int id, int *gas_phases)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialSolidSolutions2Module (int id, int *solid_solutions)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialKinetics2Module (int id, int *kinetics)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialPhreeqc2Module (int id, int *initial_conditions1, int *initial_conditions2, double *fraction1)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialPhreeqc2SpeciesConcentrations (int id, double *species_c, int n_boundary, int *boundary_solution1, int *boundary_solution2, double *fraction1)
 
IRM_DLL_EXPORT IRM_RESULT RM_InitialPhreeqcCell2Module (int id, int n, int *module_numbers, int dim_module_numbers)
 
IRM_DLL_EXPORT IRM_RESULT RM_LoadDatabase (int id, const char *db_name)
 
IRM_DLL_EXPORT IRM_RESULT RM_LogMessage (int id, const char *str)
 
IRM_DLL_EXPORT IRM_RESULT RM_MpiWorker (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_MpiWorkerBreak (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_OpenFiles (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_OutputMessage (int id, const char *str)
 
IRM_DLL_EXPORT IRM_RESULT RM_RunCells (int id)
 
IRM_DLL_EXPORT IRM_RESULT RM_RunFile (int id, int workers, int initial_phreeqc, int utility, const char *chem_name)
 
IRM_DLL_EXPORT IRM_RESULT RM_RunString (int id, int workers, int initial_phreeqc, int utility, const char *input_string)
 
IRM_DLL_EXPORT IRM_RESULT RM_ScreenMessage (int id, const char *str)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetComponentH2O (int id, int tf)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetConcentrations (int id, double *c)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetCurrentSelectedOutputUserNumber (int id, int n_user)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetDensityUser (int id, double *density)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetDensity (int id, double *density)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetDumpFileName (int id, const char *dump_name)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetErrorHandlerMode (int id, int mode)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetErrorOn (int id, int tf)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetFilePrefix (int id, const char *prefix)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetGasCompMoles (int id, double *gas_moles)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetGasPhaseVolume (int id, double *gas_volume)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetMpiWorkerCallback (int id, int(*fcn)(int *x1, void *cookie))
 
IRM_DLL_EXPORT IRM_RESULT RM_SetMpiWorkerCallbackCookie (int id, void *cookie)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetNthSelectedOutput (int id, int n)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetPartitionUZSolids (int id, int tf)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetPorosity (int id, double *por)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetPressure (int id, double *p)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetPrintChemistryMask (int id, int *cell_mask)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetPrintChemistryOn (int id, int workers, int initial_phreeqc, int utility)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetRebalanceByCell (int id, int method)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetRebalanceFraction (int id, double f)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetRepresentativeVolume (int id, double *rv)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetSaturationUser (int id, double *sat)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetSaturation (int id, double *sat)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetScreenOn (int id, int tf)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetSelectedOutputOn (int id, int selected_output)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetSpeciesSaveOn (int id, int save_on)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetTemperature (int id, double *t)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetTime (int id, double time)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetTimeConversion (int id, double conv_factor)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetTimeStep (int id, double time_step)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsExchange (int id, int option)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsGasPhase (int id, int option)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsKinetics (int id, int option)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsPPassemblage (int id, int option)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsSolution (int id, int option)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsSSassemblage (int id, int option)
 
IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsSurface (int id, int option)
 
IRM_DLL_EXPORT IRM_RESULT RM_SpeciesConcentrations2Module (int id, double *species_conc)
 
IRM_DLL_EXPORT IRM_RESULT RM_StateSave (int id, int istate)
 
IRM_DLL_EXPORT IRM_RESULT RM_StateApply (int id, int istate)
 
IRM_DLL_EXPORT IRM_RESULT RM_StateDelete (int id, int istate)
 
IRM_DLL_EXPORT IRM_RESULT RM_UseSolutionDensityVolume (int id, int tf)
 
IRM_DLL_EXPORT IRM_RESULT RM_WarningMessage (int id, const char *warn_str)
 

Detailed Description

C header file for module BMIPhreeqcRM. RM_BmiCreate creates a module with access to all BMI methods. For backward compatibility, the deprecated method RM_Create creates an old-style PhreeqcRM instance, which does not have access to BMI methods.

Macro Definition Documentation

◆ RM_INTERFACE_C_H

#define RM_INTERFACE_C_H

Function Documentation

◆ RM_Abort()

IRM_DLL_EXPORT IRM_RESULT RM_Abort ( int  id,
int  result,
const char *  err_str 
)

Abort the program. Result will be interpreted as an IRM_RESULT value and decoded; err_str will be printed; and the reaction module will be destroyed. If using MPI, an MPI_Abort message will be sent before the reaction module is destroyed. If the id is an invalid instance, RM_Abort will return a value of IRM_BADINSTANCE, otherwise the program will exit with a return code of 4.

Parameters
idThe instance id returned from RM_Create.
resultInteger treated as an IRM_RESULT return code.
err_strString to be printed as an error message.
Return values
IRM_RESULTProgram will exit before returning unless id is an invalid reaction module id.
See also

RM_Destroy, RM_ErrorMessage.
C Example:
iphreeqc_id = RM_Concentrations2Utility(id, c_well, 1, tc, p_atm);
Utilities::strcpy_safe(str, MAX_LENGTH, "SELECTED_OUTPUT 5; -pH; RUN_CELLS; -cells 1");
status = RunString(iphreeqc_id, str);
if (status != 0) status = RM_Abort(id, status, "IPhreeqc RunString failed");
MPI:
Called by root or workers.

◆ RM_BmiAddOutputVars()

IRM_DLL_EXPORT IRM_RESULT RM_BmiAddOutputVars ( int  id,
char *  option,
char *  def 
)

RM_BmiAddOutputVars allows selection of sets of variables that can be retieved by the RM_BmiGetValue methods. Sets of variables can be included or excluded with multiple calls to this method. All calls must precede the final call to the PhreeqcRM method FindComponents. FindComponents generates SELECTED_OUTPUT 333 and USER_PUNCH 333 data blocks that make the variables accessible. Variables will only be accessible if the system includes the given reactant; for example, no gas variables will be Created if there are no GAS_PHASEs in the model.

Parameters
idId number returned by RM_BmiCreate.
optionA string value, among those listed below, that includes or excludes variables from RM_BmiGetOutputVarName, RM_BmiGetValue methods, and other BMI methods.
defA string value that can be "false", "true", or a list of items to be included as accessible variables. A value of "false", excludes all variables of the given type; a value of "true" includes all variables of the given type for the current system; a list specifies a subset of items of the given type.
Return values
0is success, 0 is failure.

Values for the the parameter option:


AddOutputVars: False excludes all variables; True causes the settings for each variable group to determine the variables that will be defined. Default True;
SolutionProperties: False excludes all solution property variables; True includes variables pH, pe, alkalinity, ionic strength, water mass, charge balance, percent error, and specific conductance. Default True.
SolutionTotalMolalities: False excludes all total element and element redox state variables; True includes all elements and element redox state variables for the system defined for the calculation; list restricts variables to the specified elements and redox states. Default True.
ExchangeMolalities: False excludes all variables related to exchange; True includes all variables related to exchange; list includes variables for the specified exchange species. Default True.
SurfaceMolalities: False excludes all variables related to surfaces; True includes all variables related to surfaces; list includes variables for the specified surface species. Default True.
EquilibriumPhases: False excludes all variables related to equilibrium phases; True includes all variables related to equilibrium phases; list includes variables for the specified equilibiurm phases. Default True.
Gases: False excludes all variables related to gases; True includes all variables related to gases; list includes variables for the specified gas components. Default True.
KineticReactants: False excludes all variables related to kinetic reactants; True includes all variables related to kinetic reactants; list includes variables for the specified kinetic reactants. Default True.
SolidSolutions: False excludes all variables related to solid solutions; True includes all variables related to solid solutions; list includes variables for the specified solid solutions components. Default True.
CalculateValues: False excludes all calculate values; True includes all calculate values; list includes the specified calculate values. CALCLUATE_VALUES can be used to calculate geochemical quantities not available in the other sets of variables. Default True.
SolutionActivities: False excludes all aqueous species; True includes all aqueous species; list includes only the specified aqueous species. Default False.
SolutionMolalities: False excludes all aqueous species; True includes all aqueous species; list includes only the specified aqueous species. Default False.
SaturationIndices: False excludes all saturation indices; True includes all saturation indices; list includes only the specified saturation indices. Default False.

C example:
status = RM_BmiAddOutputVars(id, "SolutionMolalities", "True");
status = RM_BmiAddOutputVars(id, "SaturationIndices", "Calcite Dolomite");

◆ RM_BmiCreate()

IRM_DLL_EXPORT int RM_BmiCreate ( int  nxyz,
int  nthreads 
)

RM_BmiCreate Creates a BMI reaction module, which allows use of all of the RM_Bmi methods. If the code is compiled with the preprocessor directive USE_OPENMP, the reaction module is multithreaded. 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 user's model.
nthreads(or comm, MPI) When using OPENMP, the argument (nthreads) is the number of worker threads to be used. If nthreads <= 0, the number of threads is set equal to the number of processors of the computer. When using MPI, the argument (comm) is the MPI communicator to use within the reaction module.
Return values
Idof the BMIPhreeqcRM instance, negative is failure.
See also
RM_BmiFinalize.
C example:
nxyz = 40;
nthreads = 3;
id = RM_BmiCreate(nxyz, nthreads);
MPI:
Called by root and workers.

◆ RM_BmiDestroy()

IRM_DLL_EXPORT IRM_RESULT RM_BmiDestroy ( int  id)

RM_BmiDestroy Destroys a BMI reaction module; same as RM_BmiFinalize.

Parameters
idId number returned by RM_BmiCreate.
Return values
0is success, 0 is failure.
See also
RM_BmiCreate.
C example:
status = RM_BmiDestroy(id);
MPI:
Called by root and workers.

◆ RM_BmiFinalize()

IRM_DLL_EXPORT IRM_RESULT RM_BmiFinalize ( int  id)

RM_BmiFinalize Destroys a reaction module.

Parameters
idId number returned by RM_BmiCreate.
Return values
0is success, 0 is failure.
See also
RM_BmiCreate.
C example:
status = RM_BmiFinalize(id);
MPI:
Called by root and workers.

◆ RM_BmiGetComponentName()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetComponentName ( int  id,
char *  component_name,
int  l 
)

RM_BmiGetComponentName returns the component name–"BMIPhreeqcRM".

Parameters
idId number returned by RM_BmiCreate.
component_nameReturns "BMIPhreeqcRM", the name of the component.
lLength of string buffer component_name.
Return values
0is success, 1 is failure; negative indicates buffer is too small.
C example:
status = RM_BmiGetComponentName(id, component_name);
MPI:
Called by root.

◆ RM_BmiGetCurrentTime()

IRM_DLL_EXPORT double RM_BmiGetCurrentTime ( int  id)

RM_BmiGetCurrentTime returns the current simulation time, in seconds.

Parameters
idId number returned by RM_BmiCreate.
Return values
Thecurrent simulation time, in seconds.
See also
RM_BmiGetEndTime, RM_BmiGetTimeStep, RM_BmiGetTime.
C example:
now = RM_BmiGetCurrentTime(id);
MPI:
Called by root.

◆ RM_BmiGetEndTime()

IRM_DLL_EXPORT double RM_BmiGetEndTime ( int  id)

RM_BmiGetEndTime returns RM_BmiGetCurrentTime plus RM_BmiGetTimeStep, in seconds.

Parameters
idId number returned by RM_BmiCreate.
Return values
Theend of the time step, in seconds.
See also
RM_BmiGetCurrentTime, RM_BmiGetTimeStep.
C example:
status = RM_BmiGetEndTime(id);
MPI:
Called by root.

◆ RM_BmiGetGridEdgeCount()

IRM_DLL_EXPORT int RM_BmiGetGridEdgeCount ( int  id,
const int  grid 
)

RM_BmiGetGridEdgeCount is not implemented.

◆ RM_BmiGetGridEdgeNodes()

IRM_DLL_EXPORT void RM_BmiGetGridEdgeNodes ( int  id,
const int  grid,
int *  edge_nodes 
)

RM_BmiGetGridEdgeNodes is not implemented.

◆ RM_BmiGetGridFaceCount()

IRM_DLL_EXPORT int RM_BmiGetGridFaceCount ( int  id,
const int  grid 
)

RM_BmiGetGridFaceCount is not implemented.

◆ RM_BmiGetGridFaceEdges()

IRM_DLL_EXPORT void RM_BmiGetGridFaceEdges ( int  id,
const int  grid,
int *  face_edges 
)

RM_BmiGetGridFaceEdges is not implemented.

◆ RM_BmiGetGridFaceNodes()

IRM_DLL_EXPORT void RM_BmiGetGridFaceNodes ( int  id,
const int  grid,
int *  face_nodes 
)

RM_BmiGetGridFaceNodes is not implemented.

◆ RM_BmiGetGridNodeCount()

IRM_DLL_EXPORT int RM_BmiGetGridNodeCount ( int  id,
const int  grid 
)

RM_BmiGetGridNodeCount is not implemented.

◆ RM_BmiGetGridNodesPerFace()

IRM_DLL_EXPORT void RM_BmiGetGridNodesPerFace ( int  id,
const int  grid,
int *  nodes_per_face 
)

RM_BmiGetGridNodesPerFace is not implemented.

◆ RM_BmiGetGridOrigin()

IRM_DLL_EXPORT void RM_BmiGetGridOrigin ( int  id,
const int  grid,
double *  origin 
)

RM_BmiGetGridOrigin is not implemented.

◆ RM_BmiGetGridRank()

IRM_DLL_EXPORT int RM_BmiGetGridRank ( int  id,
int  grid 
)

RM_BmiGetGridRank returns a rank of 1 for grid 0. BMIPhreeqcRM has a 1D series of cells; any grid or spatial information must be found in the user's model.

Parameters
idId number returned by RM_BmiCreate.
gridGrid number, only grid 0 is considered.
Return values
Rankof 1 is returned for grid 0; 0 for all other values of grid.
C example:
rank = RM_BmiGetGridRank(id, grid)
MPI:
Called by root.

◆ RM_BmiGetGridShape()

IRM_DLL_EXPORT void RM_BmiGetGridShape ( int  id,
const int  grid,
int *  shape 
)

RM_BmiGetGridShape is not implemented.

◆ RM_BmiGetGridSize()

IRM_DLL_EXPORT int RM_BmiGetGridSize ( int  id,
int  grid 
)

RM_BmiGetGridSize returns the number of cells specified at creation of the BMIPhreeqcRM instance.

Parameters
idId number returned by RM_BmiCreate.
gridGrid number, only grid 0 is considered.
Return values
Numberof cells. Same value as RM_BmiGetValueInt(id, "GridCellCount") is returned for grid 0; 0 for all other values of grid.
C example:
nxyz = RM_BmiGetGridSize(id, grid);
MPI:
Called by root.

◆ RM_BmiGetGridSpacing()

IRM_DLL_EXPORT void RM_BmiGetGridSpacing ( int  id,
const int  grid,
double *  spacing 
)

RM_BmiGetGridSpacing is not implemented.

◆ RM_BmiGetGridType()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetGridType ( int  id,
int  grid,
char *  str,
int  l 
)

RM_BmiGetGridType defines the grid to be points. No grid information is available in BMIPhreeqcRM; all grid information must be found in the user's model.

Parameters
idId number returned by RM_BmiCreate.
gridGrid number, only grid 0 is considered.
str"Points" is returned for grid 0; "Undefined grid identifier" is returned for all other values of grid.
lLength of string buffer str.
Return values
0is success, 1 is failure, negative indicates the buffer is too small.
C example:
status = RM_BmiGetGridType(id, grid, str, l)
MPI:
Called by root.

◆ RM_BmiGetGridX()

IRM_DLL_EXPORT void RM_BmiGetGridX ( int  id,
const int  grid,
double *  x 
)

RM_BmiGetGridX is not implemented.

◆ RM_BmiGetGridY()

IRM_DLL_EXPORT void RM_BmiGetGridY ( int  id,
const int  grid,
double *  y 
)

RM_BmiGetGridY is not implemented.

◆ RM_BmiGetGridZ()

IRM_DLL_EXPORT void RM_BmiGetGridZ ( int  id,
const int  grid,
double *  z 
)

RM_BmiGetGridZ is not implemented.

◆ RM_BmiGetInputItemCount()

IRM_DLL_EXPORT int RM_BmiGetInputItemCount ( int  id)

RM_BmiGetInputItemCount returns count of variables that can be set with RM_BmiSetValue methods.

Parameters
idId number returned by RM_BmiCreate.
Return values
Numberof input variables that can be set with RM_BmiSetValue methods.
See also
RM_BmiGetInputVarName, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
count = RM_BmiGetInputItemCount(id);
MPI:
Called by root.

◆ RM_BmiGetInputVarName()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetInputVarName ( int  id,
int  i,
char *  name,
int  l 
)

RM_BmiGetInputVarName returns the ith variable name that can be set with RM_BmiSetValue methods.

Parameters
idId number returned by RM_BmiCreate.
i0-based index of variable name to retrieve.
nameRetrieved variable name.
lLength of buffer for name.
Return values
0is success, 1 is failure; negative indicates buffer is too small.
See also
RM_BmiGetInputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
char name[256];
status = RM_BmiGetInputVarName(id, 0, name, 256);
MPI:
Called by root.

◆ RM_BmiGetOutputItemCount()

IRM_DLL_EXPORT int RM_BmiGetOutputItemCount ( int  id)

RM_BmiGetOutputItemCount returns count of output variables that can be retrieved with RM_BmiGetValue methods.

Parameters
idId number returned by RM_BmiCreate.
Return values
Numberof output variables that can be retrieved with RM_BmiGetValue methods.
See also
RM_BmiGetOutputVarName, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
count = RM_BmiGetOutputItemCount(id);
MPI:
Called by root.

◆ RM_BmiGetOutputVarName()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetOutputVarName ( int  id,
int  i,
char *  name,
int  l 
)

RM_BmiGetOutputVarName returns ith variable name for which a pointer can be retrieved with RM_BmiGetValue methods.

Parameters
idId number returned by RM_BmiCreate.
i0-based index of variable name to retrieve.
nameRetrieved variable name.
lLength of buffer for name.
Return values
0is success, 1 is failure; negative indicates buffer is too small.
See also
RM_BmiGetOutputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
char name[256]
status = RM_BmiGetOutputVarName(id, 0, name, 256);
MPI:
Called by root.

◆ RM_BmiGetPointableItemCount()

IRM_DLL_EXPORT int RM_BmiGetPointableItemCount ( int  id)

RM_BmiGetPointableItemCount returns count of pointable variables that can be retrieved with RM_BmiGetValuePtr.

Parameters
idId number returned by RM_BmiCreate.
Return values
Numberof output variables that can be retrieved with RM_BmiGetValuePtr.
See also
RM_BmiGetPointableVarName, RM_BmiGetValuePtr, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
count = RM_BmiGetPointableItemCount(id);
MPI:
Called by root.

◆ RM_BmiGetPointableVarName()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetPointableVarName ( int  id,
int  i,
char *  name,
int  l 
)

RM_BmiGetPointableVarName returns ith variable name for which a pointer can be retrieved with RM_BmiGetValuePtr.

Parameters
idId number returned by RM_BmiCreate.
i0-based index of variable name to retrieve.
nameRetrieved variable name.
lLength of buffer for name.
Return values
0is success, 1 is failure; negative indicates buffer is too small.
See also
RM_BmiGetPointableItemCount, RM_BmiGetValuePtr, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
char name[256];
status = RM_BmiGetPointableVarName(id, 0, name, 256);
MPI:
Called by root.

◆ RM_BmiGetStartTime()

IRM_DLL_EXPORT double RM_BmiGetStartTime ( int  id)

RM_BmiGetStartTime returns the current simulation time, in seconds. (Same as RM_BmiGetCurrentTime.)

Parameters
idId number returned by RM_BmiCreate.
Return values
Thecurrent simulation time, in seconds.

◆ RM_BmiGetTime()

IRM_DLL_EXPORT double RM_BmiGetTime ( int  id)

RM_BmiGetTime returns the current simulation time, in seconds. (Same as RM_BmiGetCurrentTime.)

Parameters
idId number returned by RM_BmiCreate.
Return values
Thecurrent simulation time, in seconds.

◆ RM_BmiGetTimeStep()

IRM_DLL_EXPORT double RM_BmiGetTimeStep ( int  id)

RM_BmiGetTimeStep returns the current simulation time step, in seconds.

Parameters
idId number returned by RM_BmiCreate.
Return values
Thecurrent simulation time step, in seconds.
See also
RM_BmiGetCurrentTime, RM_BmiGetEndTime.
C example:
time_step = RM_BmiGetTimeStep(id)
MPI:
Called by root.

◆ RM_BmiGetTimeUnits()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetTimeUnits ( int  id,
char *  units,
int  l 
)

RM_BmiGetTimeUnits returns the time units of PhreeqcRM. All time units are seconds for PhreeqcRM.

Parameters
idId number returned by RM_BmiCreate.
unitsReturns the string "seconds".
lLength of the string buffer units.
Return values
0is success, 1 failure; negative indicates buffer is too small.
See also
RM_BmiGetCurrentTime, RM_BmiGetEndTime, RM_BmiGetTimeStep.
C example:
char time_units[256]
status = RM_BmiGetTimeUnits(id, time_units, 256)
MPI:
Called by root.

◆ RM_BmiGetValueAtIndices()

IRM_DLL_EXPORT void RM_BmiGetValueAtIndices ( int  id,
char *  name,
void *  dest,
int *  inds,
int  count 
)

RM_BmiGetValueAtIndices is not implemented

◆ RM_BmiGetValueChar()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetValueChar ( int  id,
char *  var,
char *  dest,
int  l 
)

RM_BmiGetValueChar retrieves char model variables. Only variables in the list provided by RM_BmiGetOutputVarName can be retrieved.

Parameters
idId number returned by RM_BmiCreate.
varName of the variable to retrieve.
destVariable in which to place results.
lLength of the string buffer dest.
Return values
0is success, 1 is failure; negative indicates buffer is too small.

The buffer length must be at least one character greater than the value returned by RM_BmiGetVarNbytes to allow for null termination. "ErrorString" and "FilePrefix" return single strings. "Components" and "SelectedOutputHeadings" retrieve a string that is a concatenated list of components or selected-output headings. The length of each item in a list is given by RM_BmiGetVarItemsize. The concatenated list must be processed to extract each component or heading and a null termination must be appended. Alternatively, the components can be retrieved one at a time with RM_GetComponent or RM_GetSelectedOutputHeading.

Variable names for the second argument (var).


"Components"
"ErrorString"
"FilePrefix"
"SelectedOutputHeadings".

See also
RM_BmiGetOutputVarName, RM_BmiGetOutputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
char string[256];
status = RM_BmiGetValueChar(id, "FilePrefix", string);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiGetValueDouble()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetValueDouble ( int  id,
char *  var,
double *  dest 
)

RM_BmiGetValueDouble retrieves model variables. Only variables in the list provided by RM_BmiGetOutputVarName can be retrieved.

Parameters
idId number returned by RM_BmiCreate.
varName of the variable to retrieve.
destVariable in which to place results.
Return values
0is success, 1 is failure.

Variables in addition to the ones listed below may be retrieved by this method, depending on variables selected by RM_BmiAddOutputVars. All variables added by RM_BmiAddOutputVars will be double arrays of size equal to the number of model cells [RM_BmiGetValueInt(id, "GridCellCount")].

Variable names for the second argument (var).


"Concentrations"
"DensityCalculated"
"Gfw"
"Porosity"
"Pressure"
"SaturationCalculated"
"SelectedOutput"
"SolutionVolume"
"Temperature"
"Time"
"TimeStep"
"Viscosity".

See also
RM_BmiGetOutputVarName, RM_BmiGetOutputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
density = (double *)malloc(nxyz*sizeof(double));
status = RM_BmiGetValueDouble(id, "DensityCalculated", density);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiGetValueInt()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetValueInt ( int  id,
char *  var,
int *  dest 
)

RM_BmiGetValueInt retrieves int model variables. Only variables in the list provided by RM_BmiGetOutputVarName can be retrieved.

Parameters
idId number returned by RM_BmiCreate.
varName of the variable to retrieve.
destVariable in which to place results.
Return values
0is success, 1 is failure.

Variable names for the second argument (var).


"ComponentCount"
"CurrentSelectedOutputUserNumber"
"GridCellCount"
"SelectedOutputColumnCount"
"SelectedOutputCount"
"SelectedOutputOn"
"SelectedOutputRowCount".

See also
RM_BmiGetOutputVarName, RM_BmiGetOutputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
status = RM_BmiGetValueInt(id, "ComponentCount", &count);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiGetValuePtr()

IRM_DLL_EXPORT void * RM_BmiGetValuePtr ( int  id,
char *  var 
)

RM_BmiGetValuePtr retrieves pointers to model variables. Only variables in the list provided by RM_BmiGetPointableVarName can be pointed to.

Parameters
idId number returned by RM_BmiCreate.
varName of the variable to retrieve.
Return values
Pointerto an up-to-date copy of the variable's data.

The following list gives the name in the second argument (var) and the data type the pointer:


"ComponentCount"
"Concentrations"
"DensityCalculated"
"Gfw"
"GridCellCount"
"Porosity"
"Pressure"
"SaturationCalculated"
"SelectedOutputOn"
"SolutionVolume"
"Temperature"
"Time"
"TimeStep"
"Viscosity"

MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiGetVarGrid()

IRM_DLL_EXPORT int RM_BmiGetVarGrid ( int  id,
char *  var 
)

RM_BmiGetVarGrid returns a value of 1, indicating points. BMIPhreeqcRM does not have a grid of its own. The cells of BMIPhreeqcRM are associated with the user's model grid, and all spatial characterists are assigned by the user's model.

Parameters
idId number returned by RM_BmiCreate.
varVaraiable name. (Return value is the same regardless of value of @ var.)
Return values
1(points). BMIPhreeqcRM cells derive meaning from the user's model.

◆ RM_BmiGetVarItemsize()

IRM_DLL_EXPORT int RM_BmiGetVarItemsize ( int  id,
char *  name 
)

RM_BmiGetVarItemsize retrieves the size, in bytes, of a variable that can be set with RM_BmiSetValue methods, retrieved with RM_BmiGetValue methods, or pointed to with RM_BmiGetValuePtr. Sizes may be the size of an integer, double, or a character length for string variables.

Parameters
idId number returned by RM_BmiCreate.
nameName of the variable to retrieve the item size.
Return values
Size,inbytes, of one element of the variable.
See also
RM_BmiGetInputVarName, RM_BmiGetInputItemCount, RM_BmiGetOutputVarName, RM_BmiGetOutputItemCount, RM_BmiGetPointableVarName, RM_BmiGetPointableItemCount, RM_BmiGetValuePtr, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
for(i = 0; i < GetInputVarCount(id); i++)
{
    itemsize = GetVarItemsize(id, name);
}
MPI:
Called by root.

◆ RM_BmiGetVarNbytes()

IRM_DLL_EXPORT int RM_BmiGetVarNbytes ( int  id,
char *  name 
)

RM_BmiGetVarNbytes retrieves the total number of bytes needed for a variable that can be set with RM_BmiSetValue methods, retrieved with RM_BmiGetValue methods, or pointed to with RM_BmiGetValuePtr.

Parameters
idId number returned by RM_BmiCreate.
nameName of the variable to retrieve the number of bytes needed to retrieve or store the variable.
Return values
Totalnumber of bytes needed for the variable.
See also
RM_BmiGetInputVarName, RM_BmiGetInputItemCount, RM_BmiGetOutputVarName, RM_BmiGetOutputItemCount, RM_BmiGetPointableVarName, RM_BmiGetPointableItemCount, RM_BmiGetValuePtr, RM_BmiGetVarItemsize, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
for(i = 0; i < GetInputVarCount(id); i++)
{
    nbytes = GetVarNbytes(id, name);
}
MPI:
Called by root.

◆ RM_BmiGetVarType()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetVarType ( int  id,
char *  name,
char *  vtype,
int  l 
)

RM_BmiGetVarType retrieves the type of a variable that can be set with RM_BmiSetValue methods, retrieved with RM_BmiGetValue methods, or pointed to with RM_BmiGetValuePtr. Types are "char", "double", or "int", or an array of these types.

Parameters
idId number returned by RM_BmiCreate.
nameName of the variable to retrieve the type.
vtypeType of the variable.
lLength of string buffer vtype.
Return values
0is success, 1 is failure; negative indicates buffer is too small.
See also
RM_BmiGetInputVarName, RM_BmiGetInputItemCount, RM_BmiGetOutputVarName, RM_BmiGetOutputItemCount, RM_BmiGetPointableVarName, RM_BmiGetPointableItemCount, RM_BmiGetValuePtr, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarUnits.
C example:
char string[256];
for(i = 0; i < GetInputVarCount(id); i++)
{
    status = GetVarType(id, i, string, 256);
}
MPI:
Called by root.

◆ RM_BmiGetVarUnits()

IRM_DLL_EXPORT IRM_RESULT RM_BmiGetVarUnits ( int  id,
char *  name,
char *  units,
int  l 
)

RM_BmiGetVarType retrieves the units of a variable that can be set with RM_BmiSetValue methods, retrieved with RM_BmiGetValue methods, or pointed to with RM_BmiGetValuePtr.

Parameters
idId number returned by RM_BmiCreate.
nameName of the variable to retrieve the type.
unitsUnits of the variable.
lLength of string buffer units.
Return values
0is success, 1 is failure; negative indicates buffer is too small.
See also
RM_BmiGetInputVarName, RM_BmiGetInputItemCount, RM_BmiGetOutputVarName, RM_BmiGetOutputItemCount, RM_BmiGetPointableVarName, RM_BmiGetPointableItemCount, RM_BmiGetValuePtr, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType.
C example:
char string[256];
for(i = 0; i < GetInputVarCount(id); i++)
{
    status = GetVarUnits(id, i, string, 256);
}
MPI:
Called by root.

◆ RM_BmiInitialize()

IRM_DLL_EXPORT IRM_RESULT RM_BmiInitialize ( int  id,
char *  config_file 
)

RM_BmiInitialize uses a YAML file to initialize an instance of BMIPhreeqcRM.

Parameters
idId number returned by RM_BmiCreate.
config_fileString containing the YAML file name.
Return values
0is success, 1 is failure.

The file contains a YAML map of PhreeqcRM methods and the arguments corresponding to the methods. For example,

- key: LoadDatabase
  database: phreeqc.dat
- key: RunFile
  workers: true
  initial_phreeqc: true
  utility: true
  chemistry_name: advect.pqi

RM_BmiInitialize will read the YAML file and execute the specified methods with the specified arguments. Using YAML terminology, the argument(s) for a method may be a scalar, a sequence, or a map, depending if the argument is a single item, a single vector, or there are multiple arguments. In the case of a map, the name associated with each argument (for example "chemistry_name" above) is arbitrary. The names of the map keys for map arguments are not used in parsing the YAML file; only the order of the arguments is important.

The following list gives the PhreeqcRM methods that can be specified in a YAML file and the arguments that are required. The arguments are described with C++ formats, which are sufficient to identify which arguments are YAML scalars (single bool/logical, int, double, string/character argument), sequences (single vector argument), or maps (multiple arguments).

CloseFiles();
CreateMapping(std::vector< int >& grid2chem);
DumpModule();
FindComponents();
InitialEquilibriumPhases2Module(std::vector< int > equilibrium_phases);
InitialExchanges2Module(std::vector< int > exchanges);
InitialGasPhases2Module(std::vector< int > gas_phases);
InitialKineticss2Module(std::vector< int > kinetics);
InitialSolidSolutions2Module(std::vector< int > solid_solutions);
InitialSolutions2Module(std::vector< int > solutions);
InitialSurfaces2Module(std::vector< int > surfaces);
InitialPhreeqc2Module(std::vector< int > initial_conditions1);
InitialPhreeqc2Module(std::vector< int > initial_conditions1,
std::vector< int > initial_conditions2, std::vector< double > fraction1);
InitialPhreeqcCell2Module(int n, std::vector< int > cell_numbers);
LoadDatabase(std::string database);
OpenFiles();
OutputMessage(std::string str);
RunCells();
RunFile(bool workers, bool initial_phreeqc, bool utility, std::string chemistry_name);
RunString(bool workers, bool initial_phreeqc, bool utility, std::string input_string);
ScreenMessage(std::string str);
SetComponentH2O(bool tf);
SetConcentrations(std::vector< double > c);
SetCurrentSelectedOutputUserNumber(int n_user);
SetDensityUser(std::vector< double > density);
SetDumpFileName(std::string dump_name);
SetErrorHandlerMode(int mode);
SetErrorOn(bool tf);
SetFilePrefix(std::string prefix);
SetGasCompMoles(std::vector< double > gas_moles);
SetGasPhaseVolume(std::vector< double > gas_volume);
SetPartitionUZSolids(bool tf);
SetPorosity(std::vector< double > por);
SetPressure(std::vector< double > p);
SetPrintChemistryMask(std::vector< int > cell_mask);
SetPrintChemistryOn(bool workers, bool initial_phreeqc, bool utility);
SetRebalanceByCell(bool tf);
SetRebalanceFraction(double f);
SetRepresentativeVolume(std::vector< double > rv);
SetSaturationUser(std::vector< double > sat);
SetScreenOn(bool tf);
SetSelectedOutputOn(bool tf);
SetSpeciesSaveOn(bool save_on);
SetTemperature(std::vector< double > t);
SetTime(double time);
SetTimeConversion(double conv_factor);
SetTimeStep(double time_step);
SetUnitsExchange(int option);
SetUnitsGasPhase(int option);
SetUnitsKinetics(int option);
SetUnitsPPassemblage(int option);
SetUnitsSolution(int option);
SetUnitsSSassemblage(int option);
SetUnitsSurface(int option);
SpeciesConcentrations2Module(std::vector< double > species_conc);
StateSave(int istate);
StateApply(int istate);
StateDelete(int istate);
UseSolutionDensityVolume(bool tf);
WarningMessage(std::string warnstr);

C example:
id = RM_BmiCreate(nxyz, nthreads);
status = RM_BmiInitializeYAML(id, "myfile.yaml");
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiSetValueAtIndices()

IRM_DLL_EXPORT void RM_BmiSetValueAtIndices ( int  id,
char *  name,
int *  inds,
int  count,
void *  src 
)

RM_BmiSetValueAtIndices is not implemented.

◆ RM_BmiSetValueChar()

IRM_DLL_EXPORT IRM_RESULT RM_BmiSetValueChar ( int  id,
char *  name,
const char *  src 
)

RM_BmiSetValueChar sets model character variables. Only variables in the list provided by RM_BmiGetInputVarName can be set.

Parameters
idId number returned by RM_BmiCreate.
nameName of variable to set.
srcData to use to set the variable.
Return values
0is success, 1 is failure.

Variable names for the second argument (var):


"FilePrefix"

See also
RM_BmiGetInputVarName, RM_BmiGetInputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
status = RM_BmiSetValueChar(id, "FilePrefix", "my_prefix");
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiSetValueDouble()

IRM_DLL_EXPORT IRM_RESULT RM_BmiSetValueDouble ( int  id,
char *  name,
double  src 
)

RM_BmiSetValueDouble sets model double variables. Only variables in the list provided by RM_BmiGetInputVarName can be set.

Parameters
idId number returned by RM_BmiCreate.
nameName of variable to set.
srcData to use to set the variable.
Return values
0is success, 1 is failure.

Variable names for the second argument (var):


"Time"
"TimeStep".

See also
RM_BmiGetInputVarName, RM_BmiGetInputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
status = RM_BmiSetValueDouble(id, "TimeStep", 86400.0);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiSetValueDoubleArray()

IRM_DLL_EXPORT IRM_RESULT RM_BmiSetValueDoubleArray ( int  id,
char *  name,
double *  src 
)

RM_BmiSetValueDoubleArray sets model double array variables. Only variables in the list provided by RM_BmiGetInputVarName can be set.

Parameters
idId number returned by RM_BmiCreate.
nameName of variable to set.
srcData to use to set the variable.
Return values
0is success, 1 is failure.

Variable names for the second argument (var):


"Concentrations"
"DensityUser"
"Porosity"
"Pressure"
"SaturationUser"
"Temperature".

See also
RM_BmiGetInputVarName, RM_BmiGetInputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
tc = (double *)malloc(nxyz*sizeof(double));
for(i=0; i < nxyz; i++) tc[i] = 28.0e0;
status = RM_BmiSetValueDoubleArray(id, "Temperature", tc);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiSetValueInt()

IRM_DLL_EXPORT IRM_RESULT RM_BmiSetValueInt ( int  id,
char *  name,
int  src 
)

RM_BmiSetValueInt sets model int variables. Only variables in the list provided by RM_BmiGetInputVarName can be set.

Parameters
idId number returned by RM_BmiCreate.
nameName of variable to set.
srcData to use to set the variable.
Return values
0is success, 1 is failure.

Variable names for the second argument (var):


"NthSelectedOutput"
"SelectedOutputOn".

See also
RM_BmiGetInputVarName, RM_BmiGetInputItemCount, RM_BmiGetVarItemsize, RM_BmiGetVarNbytes, RM_BmiGetVarType, RM_BmiGetVarUnits.
C example:
status = RM_BmiSetValueInt(id, "SelectedOutputOn", 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiUpdate()

IRM_DLL_EXPORT IRM_RESULT RM_BmiUpdate ( int  id)

RM_BmiUpdate runs a reaction step for all of the cells in the reaction module.

Parameters
idId number returned by RM_BmiCreate.
Return values
0is success, 1 is failure.

Tranported concentrations are transferred to the reaction cells (Concentrations) before reaction calculations are run. The length of time over which kinetic reactions are integrated is set by TimeStep. Other properties that may need to be updated as a result of the transport calculations include porosity, pressure, saturation, temperature.

C example:
status = RM_BmiSetValue(id, "Porosity", por);                ! If pore volume changes
status = RM_BmiSetValue(id, "SaturationUser", sat);          ! If saturation changes
status = RM_BmiSetValue(id, "Temperature", temperature);     ! If temperature changes
status = RM_BmiSetValue(id, "Pressure", pressure);           ! If pressure changes
status = RM_BmiSetValue(id, "Concentrations", c);            ! Transported concentrations
status = RM_BmiSetValue(id, "TimeStep", time_step);          ! Time step for kinetic reactions
status = RM_BmiUpdate(id);
status = RM_BmiGetValue(id, "Concentrations", c);            ! Concentrations after reaction
status = RM_BmiGetValue(id, "DensityCalculated", density);   ! Density after reaction
status = RM_BmiGetValue(id, "SolutionVolume", volume);       ! Solution volume after reaction
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_BmiUpdateUntil()

IRM_DLL_EXPORT IRM_RESULT RM_BmiUpdateUntil ( int  id,
double  end_time 
)

RM_BmiUpdateUntil is the same as RM_BmiUpdate, except the time step is calculated from the argument end_time. The time step is calculated to be end_time minus the current time (RM_BmiGetCurrentTime).

Parameters
idId number returned by RM_BmiCreate..
end_timeTime at the end of the time step.
See also
RM_BmiInitialize, RM_BmiUpdate.
C example:
status = RM_BmiSetValue(id, "Time", time);
status = RM_BmiSetValue(id, "Concentrations", c);
status = RM_BmiUpdateUntil(id, time + 86400.0);
status = RM_BmiGetValue(id, "Concentrations", c);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_CloseFiles()

IRM_DLL_EXPORT IRM_RESULT RM_CloseFiles ( int  id)

Close the output and log files.

Parameters
idThe instance id returned from RM_Create.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_OpenFiles, RM_SetFilePrefix.
C Example:
status = RM_CloseFiles(id);
MPI:
Called only by root.

◆ RM_Concentrations2Utility()

IRM_DLL_EXPORT int RM_Concentrations2Utility ( int  id,
double *  c,
int  n,
double *  tc,
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 (RM_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 fragments below make a mixture of concentrations and then calculate the pH of the mixture.

Parameters
idThe instance id returned from RM_Create.
cArray of concentrations to be made SOLUTIONs in Utility IPhreeqc. Array storage is n * ncomps.
nThe number of sets of concentrations.
tcArray of temperatures to apply to the SOLUTIONs, in degree C. Array of size n.
p_atmArray of pressures to apply to the SOLUTIONs, in atm. Array of size n.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
C Example:
c_well = (double *) malloc((size_t) ((size_t) (1 * ncomps * sizeof(double))));
for (i = 0; i < ncomps; i++)
{
  c_well[i] = 0.5 * c[0 + nxyz*i] + 0.5 * c[9 + nxyz*i];
}
tc = (double *) malloc((size_t) (1 * sizeof(double)));
p_atm = (double *) malloc((size_t) (1 * sizeof(double)));
tc[0] = 15.0;
p_atm[0] = 3.0;
iphreeqc_id = RM_Concentrations2Utility(id, c_well, 1, tc, p_atm);
Utilities::strcpy_safe(str, MAX_LENGTH, "SELECTED_OUTPUT 5; -pH; RUN_CELLS; -cells 1");
status = RunString(iphreeqc_id, str);
status = SetCurrentSelectedOutputUserNumber(iphreeqc_id, 5);
status = GetSelectedOutputValue2(iphreeqc_id, 1, 0, &vtype, &pH, svalue, 100);
MPI:
Called only by root.

◆ RM_Create()

IRM_DLL_EXPORT int RM_Create ( int  nxyz,
int  nthreads 
)

Creates a reaction module without BMI methods. This method is deprecated and included only for backward compatibility. Use RM_BmiCreate to create a reaction module with access to all RM_Bmi methods.

If the code is compiled with the preprocessor directive USE_OPENMP, the reaction module is multithreaded. 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 user's model.
nthreads(or comm, MPI) When using OPENMP, the argument (nthreads) is the number of worker threads to be used. If nthreads <= 0, the number of threads is set equal to the number of processors of the computer. When using MPI, the argument (comm) is the MPI communicator to use within the reaction module.
Return values
Idof the PhreeqcRM instance, negative is failure (See RM_DecodeError).
See also

RM_Destroy.
C Example:
nxyz = 40;
#ifdef USE_MPI
  id = RM_Create(nxyz, MPI_COMM_WORLD);
  if (MPI_Comm_rank(MPI_COMM_WORLD, &mpi_myself) != MPI_SUCCESS)
  {
    exit(4);
  }
  if (mpi_myself > 0)
  {
    status = RM_MpiWorker(id);
    status = RM_Destroy(id);
    return;
  }
#else
  nthreads = 3;
  id = RM_Create(nxyz, nthreads);
#endif
MPI:
Called by root and workers.

◆ RM_CreateMapping()

IRM_DLL_EXPORT IRM_RESULT RM_CreateMapping ( int  id,
int *  grid2chem 
)

Provides a mapping from grid cells in the user's model to reaction cells in PhreeqcRM. 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, RM_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
idThe instance id returned from RM_Create.
grid2chemAn array of integers: Nonnegative is a reaction cell number (0 based), negative is an inactive cell. Array of size nxyz (number of grid cells).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
C Example:
// For demonstation, two equivalent rows by symmetry
grid2chem = (int *) malloc((size_t) (nxyz * sizeof(int)));
for (i = 0; i < nxyz/2; i++)
{
  grid2chem[i] = i;
  grid2chem[i+nxyz/2] = i;
}
status = RM_CreateMapping(id, grid2chem);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_DecodeError()

IRM_DLL_EXPORT IRM_RESULT RM_DecodeError ( int  id,
int  e 
)

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

Parameters
idThe instance id returned from RM_Create.
eAn IRM_RESULT value returned by one of the reaction-module methods.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
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 = RM_CreateMapping(id, grid2chem);
if (status < 0) status = RM_DecodeError(id, status);
MPI:
Can be called by root and (or) workers.

◆ RM_Destroy()

IRM_DLL_EXPORT IRM_RESULT RM_Destroy ( int  id)

Destroys a reaction module.

Parameters
idThe instance id returned from RM_Create.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_Create.
C Example:
status = RM_Destroy(id);
MPI:
Called by root and workers.

◆ RM_DumpModule()

IRM_DLL_EXPORT IRM_RESULT RM_DumpModule ( int  id,
int  dump_on,
int  append 
)

Writes the contents of all workers to file in _RAW formats, including SOLUTIONs and all reactants.

Parameters
idThe instance id returned from RM_Create.
dump_onSignal for writing the dump file: 1 true, 0 false.
appendSignal to append to the contents of the dump file: 1 true, 0 false.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetDumpFileName.
C Example:
dump_on = 1;
append = 0;
status = RM_SetDumpFileName(id, "advection_c.dmp");
status = RM_DumpModule(id, dump_on, append);
MPI:
Called by root; workers must be in the loop of RM_MpiWorker.

◆ RM_ErrorMessage()

IRM_DLL_EXPORT IRM_RESULT RM_ErrorMessage ( int  id,
const char *  errstr 
)

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

Parameters
idThe instance id returned from RM_Create.
errstrString to be printed.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_LogMessage,
RM_OpenFiles,
RM_OutputMessage, RM_ScreenMessage, RM_WarningMessage.
C Example:
status = RM_ErrorMessage(id, "Goodbye world");
MPI:
Called by root and (or) workers; root writes to output and log files.

◆ RM_FindComponents()

IRM_DLL_EXPORT int RM_FindComponents ( int  id)

Returns the number of items in the list of all elements in the InitialPhreeqc instance. Elements are those that have been defined in a solution or any other reactant (EQUILIBRIUM_PHASE, KINETICS, and others). The 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 (RM_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 (RM_GetSpeciesConcentrations) and to set new solution concentrations after MCD by using individual species concentrations (RM_SpeciesConcentrations2Module). To use these methods the save-species property needs to be turned on (RM_SetSpeciesSaveOn). If the save-species property is on, RM_FindComponents will generate a list of aqueous species (RM_GetSpeciesCount, RM_GetSpeciesName), their diffusion coefficients at 25 C (RM_GetSpeciesD25), their charge (RM_GetSpeciesZ).

Parameters
idThe instance id returned from RM_Create.
Return values
Numberof components currently in the list, or IRM_RESULT error code (see RM_DecodeError).
See also

RM_GetComponent, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName, RM_GetSpeciesZ, RM_SetComponentH2O, RM_SetSpeciesSaveOn, RM_SpeciesConcentrations2Module.
The RM_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
RM_GetEquilibriumPhasesName, RM_GetEquilibriumPhasesCount, RM_GetExchangeName, RM_GetExchangeSpeciesName, RM_GetExchangeSpeciesCount, RM_GetGasComponentsName, RM_GetGasComponentsCount, RM_GetKineticReactionsName, RM_GetKineticReactionsCount, RM_GetSICount, RM_GetSIName, RM_GetSolidSolutionComponentsName, RM_GetSolidSolutionComponentsCount, RM_GetSolidSolutionName, RM_GetSurfaceName, RM_GetSurfaceSpeciesName, RM_GetSurfaceSpeciesCount, RM_GetSurfaceType.
C Example:
// Get list of components
ncomps = RM_FindComponents(id);
components = (char **) malloc((size_t) (ncomps * sizeof(char *)));
for (i = 0; i < ncomps; i++)
{
  components[i] = (char *) malloc((size_t) (100 * sizeof(char *)));
  status = RM_GetComponent(id, i, components[i], 100);
}
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetBackwardMapping()

IRM_DLL_EXPORT IRM_RESULT RM_GetBackwardMapping ( int  id,
int  n,
int *  list,
int *  size 
)

Fills an array with the cell numbers in the user's numbering sytstem that map to a cell in the PhreeqcRM numbering system. The mapping is defined by RM_CreateMapping.

Parameters
idThe instance id returned from RM_Create.
nA cell number in the PhreeqcRM numbering system (0 <= n < RM_GetChemistryCellCount).
listArray to store the user cell numbers mapped to PhreeqcRM cell n.
sizeInput, the allocated size of list; it is an error if the array is too small. Output, the number of cells mapped to cell n.
Return values
IRM_RESULTerror code (see RM_DecodeError).
See also

RM_CreateMapping, RM_GetChemistryCellCount, RM_GetGridCellCount.
C Example:
if (RM_GetBackwardMapping(rm_id, rm_cell_number, list, &size) == 0)
{
  if (strcmp(str, "HYDRAULIC_K") == 0)
  {
    return data->K_ptr[list[0]];
  }
}
MPI:
Called by root and (or) workers.

◆ RM_GetChemistryCellCount()

IRM_DLL_EXPORT int RM_GetChemistryCellCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Numberof chemistry cells, or IRM_RESULT error code (see RM_DecodeError).
See also

RM_CreateMapping, RM_GetGridCellCount.
C Example:
status = RM_CreateMapping(id, grid2chem);
nchem = RM_GetChemistryCellCount(id);
MPI:
Called by root and (or) workers.

◆ RM_GetComponent()

IRM_DLL_EXPORT IRM_RESULT RM_GetComponent ( int  id,
int  num,
char *  chem_name,
int  l 
)

Retrieves an item from the reaction-module component list that was generated by calls to RM_FindComponents.

Parameters
idThe instance id returned from RM_Create.
numThe number of the component to be retrieved. C, 0 based.
chem_nameThe string value associated with component num.
lThe length of the maximum number of characters for chem_name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetComponentCount
C Example:
// Get list of components
ncomps = RM_FindComponents(id);
components = (char **) malloc((size_t) (ncomps * sizeof(char *)));
for (i = 0; i < ncomps; i++)
{
  components[i] = (char *) malloc((size_t) (100 * sizeof(char *)));
  status = RM_GetComponent(id, i, components[i], 100);
}
MPI:
Called by root and (or) workers.

◆ RM_GetComponentCount()

IRM_DLL_EXPORT int RM_GetComponentCount ( int  id)

Returns the number of components in the reaction-module component list. The component list is generated by calls to RM_FindComponents. The return value from the last call to RM_FindComponents is equal to the return value from RM_GetComponentCount.

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of components in the reaction-module component list, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetComponent.
C Example:
ncomps1 = RM_GetComponentCount(id);
MPI:
Called by root.

◆ RM_GetConcentrations()

IRM_DLL_EXPORT IRM_RESULT RM_GetConcentrations ( int  id,
double *  c 
)

Transfer solution concentrations from each reaction cell to the concentration array given in the argument list (c). Units of concentration for c are defined by RM_SetUnitsSolution. For concentration units of per liter, the 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 (RM_SetSaturationUser), porosity (RM_SetPorosity), and representative volume (RM_SetRepresentativeVolume), and the mass of solution is volume times density as defined by RM_SetDensityUser. RM_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
idThe instance id returned from RM_Create.
cArray to receive the concentrations. Dimension of the array is nxyz * ncomps, where nxyz is the number of user grid cells and ncomps is the result of RM_FindComponents or RM_GetComponentCount. Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetComponentCount, RM_GetSaturationCalculated, RM_SetConcentrations, RM_SetDensityUser, RM_SetRepresentativeVolume, RM_SetSaturationUser, RM_SetUnitsSolution, RM_UseSolutionDensityVolume.
C Example:
c = (double *) malloc((size_t) (ncomps * nxyz * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetConcentrations(id, c);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetCurrentSelectedOutputUserNumber()

IRM_DLL_EXPORT int RM_GetCurrentSelectedOutputUserNumber ( int  id)

Returns the user number of the current selected-output definition. RM_SetCurrentSelectedOutputUserNumber or RM_SetNthSelectedOutput specifies which of the selected-output definitions is used.

Return values
Usernumber of the the current selected-output definition, negative is failure (See RM_DecodeError).
See also
RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputCount, RM_GetSelectedOutputHeading, RM_GetSelectedOutputRowCount, RM_SetCurrentSelectedOutputUserNumber, RM_SetNthSelectedOutput, RM_SetSelectedOutputOn.
C Example:
    for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
    {
      status = RM_SetNthSelectedOutputUser(id, isel);
      n_user = RM_GetCurrentSelectedOutputUserNumber(id);
      col = RM_GetSelectedOutputColumnCount(id);
      selected_out = (double *) malloc((size_t) (col * nxyz * sizeof(double)));
      status = RM_GetSelectedOutput(id, selected_out);
      // Process results here
      free(selected_out);
    }
    
MPI:
Called by root.

◆ RM_GetDensity()

IRM_DLL_EXPORT IRM_RESULT RM_GetDensity ( int  id,
double *  density 
)

Deprecated equivalent of RM_GetDensityCalculated.

◆ RM_GetDensityCalculated()

IRM_DLL_EXPORT IRM_RESULT RM_GetDensityCalculated ( int  id,
double *  density 
)

Transfer solution densities from the reaction cells to the array given in the argument list (density). Densities are those calculated by the reaction module. This method always returns the calculated densities; RM_SetDensityUser does not affect the result. Only the following databases distributed with PhreeqcRM have molar volume information needed to accurately calculate density: phreeqc.dat, Amm.dat, and pitzer.dat.

Parameters
idThe instance id returned from RM_Create.
densityArray to receive the densities. Dimension of the array is nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
C Example:
density = (double *) malloc((size_t) (nxyz * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetDensityCalculated(id, density);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetEndCell()

IRM_DLL_EXPORT IRM_RESULT RM_GetEndCell ( int  id,
int *  ec 
)

Returns an array with the ending cell numbers from the range of cell numbers assigned to each worker.

Parameters
idThe instance id returned from RM_Create.
ecArray to receive the ending cell numbers. Dimension of the array is the number of threads (OpenMP) or the number of processes (MPI).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_Create, RM_GetMpiTasks, RM_GetStartCell, RM_GetThreadCount.
C Example:
n = RM_GetThreadCount(id) * RM_GetMpiTasks(id);
ec = (int *) malloc((size_t) (n * sizeof(int)));
status = RM_GetEndCell(id, ec);
MPI:
Called by root and (or) workers.

◆ RM_GetEquilibriumPhasesCount()

IRM_DLL_EXPORT int RM_GetEquilibriumPhasesCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of equilibrium phases in the initial-phreeqc module.
See also
RM_FindComponents, RM_GetEquilibriumPhasesName.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -equilibrium_phases\n");
for (i = 0; i < RM_GetEquilibriumPhasesCount(id); i++)
{
status = RM_GetEquilibriumPhasesName(id, i, line1, 100);
sprintf(line, "%4s%20s\n", "    ", line1);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetEquilibriumPhasesName()

IRM_DLL_EXPORT IRM_RESULT RM_GetEquilibriumPhasesName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves an item from the equilibrium phase list. The list includes all phases included in any EQUILIBRIUM_PHASES definitions in the initial-phreeqc module. RM_FindComponents must be called before RM_GetEquilibriumPhasesName. This method may be useful when generating selected output definitions related to equilibrium phases.

Parameters
idThe instance id returned from RM_Create.
numThe number of the equilibrium phase name to be retrieved. (0 basedindex.)
nameThe equilibrium phase name at number num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetEquilibriumPhasesCount.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -equilibrium_phases\n");
for (i = 0; i < RM_GetEquilibriumPhasesCount(id); i++)
{
status = RM_GetEquilibriumPhasesName(id, i, line1, 100);
sprintf(line, "%4s%20s\n", "    ", line1);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetErrorString()

IRM_DLL_EXPORT IRM_RESULT RM_GetErrorString ( int  id,
char *  errstr,
int  l 
)

Returns a string containing error messages related to the last call to a PhreeqcRM method to the character argument (errstr).

Parameters
idThe instance id returned from RM_Create.
errstrThe error string related to the last call to a PhreeqcRM method.
lMaximum number of characters that can be written to the argument (errstr).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
C Example:
if (status != IRM_OK)
{
  l = RM_GetErrorStringLength(id);
  errstr = (char *) malloc((size_t) (l * sizeof(char) + 1));
  RM_GetErrorString(id, errstr, l+1);
  fprintf(stderr,"%s", errstr);
  free(errstr);
  RM_Destroy(id);
  exit(1);
}
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetErrorStringLength()

IRM_DLL_EXPORT int RM_GetErrorStringLength ( int  id)

Returns the length of the string that contains error messages related to the last call to a PhreeqcRM method.

Parameters
idThe instance id returned from RM_Create.
Return values
intLength of the error message string (for C, equivalent to strlen, does not include terminating \0).
C Example:
if (status != IRM_OK)
{
  l = RM_GetErrorStringLength(id);
  errstr = (char *) malloc((size_t) (l * sizeof(char) + 1));
  RM_GetErrorString(id, errstr, l+1);
  fprintf(stderr,"%s", errstr);
  free(errstr);
  RM_Destroy(id);
  exit(1);
}
MPI:
Called by root, workers must be in the loop of RM_MpiWorker..

◆ RM_GetExchangeName()

IRM_DLL_EXPORT IRM_RESULT RM_GetExchangeName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves an item from the exchange name list. RM_FindComponents must be called before RM_GetExchangeName. The exchange names vector is the same length as the exchange species names vector and provides the corresponding exchange site (for example, X corresponing to NaX). This method may be useful when generating selected output definitions related to exchangers.

Parameters
idThe instance id returned from RM_Create.
numThe number of the exchange name to be retrieved. (0 based index.)
nameThe exchange name associated with exchange species num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetExchangeSpeciesCount, RM_GetExchangeSpeciesName.
C Example:
for (i = 0; i < RM_GetExchangeSpeciesCount(id); i++)
{
Utilities::strcpy_safe(line, MAX_LENGTH, "");
status = RM_GetExchangeSpeciesName(id, i, line1, 100);
status = RM_GetExchangeName(id, i, line2, 100);
sprintf(line, "%4s%20s%3s%20s\n", "    ", line1, " # ", line2);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetExchangeSpeciesCount()

IRM_DLL_EXPORT int RM_GetExchangeSpeciesCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of exchange species in the initial-phreeqc module.
See also
RM_FindComponents, RM_GetExchangeSpeciesName, RM_GetExchangeName.
C Example:
for (i = 0; i < RM_GetExchangeSpeciesCount(id); i++)
{
Utilities::strcpy_safe(line, MAX_LENGTH), "");
status = RM_GetExchangeSpeciesName(id, i, line1, 100);
status = RM_GetExchangeName(id, i, line2, 100);
sprintf(line, "%4s%20s%3s%20s\n", "    ", line1, " # ", line2);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetExchangeSpeciesName()

IRM_DLL_EXPORT IRM_RESULT RM_GetExchangeSpeciesName ( int  id,
int  num,
char *  name,
int  l1 
)

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

Parameters
idThe instance id returned from RM_Create.
numThe number of the exchange species to be retrieved. (0 based index.)
nameThe exchange species name at number num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetExchangeSpeciesCount, RM_GetExchangeName.
C Example:
for (i = 0; i < RM_GetExchangeSpeciesCount(id); i++)
{
Utilities::strcpy_safe(line, MAX_LENGTH, "");
status = RM_GetExchangeSpeciesName(id, i, line1, 100);
status = RM_GetExchangeName(id, i, line2, 100);
sprintf(line, "%4s%20s%3s%20s\n", "    ", line1, " # ", line2);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetFilePrefix()

IRM_DLL_EXPORT IRM_RESULT RM_GetFilePrefix ( int  id,
char *  prefix,
int  l 
)

Returns the reaction-module file prefix to the character argument (prefix).

Parameters
idThe instance id returned from RM_Create.
prefixCharacter string where the prefix is written.
lMaximum number of characters that can be written to the argument (prefix).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetFilePrefix.
C Example:
char str[100], str1[200];
status = RM_GetFilePrefix(id, str, 100);
Utilities::strcpy_safe(str1, MAX_LENGTH, "File prefix: ");
Utilities::strcat_safe(str1, MAX_LENGTH, str);
Utilities::strcat_safe(str1, MAX_LENGTH, "\n");
status = RM_OutputMessage(id, str1);
MPI:
Called by root and (or) workers.

◆ RM_GetGasCompMoles()

IRM_DLL_EXPORT IRM_RESULT RM_GetGasCompMoles ( int  id,
double *  gas_moles 
)

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

Parameters
idThe instance id returned from RM_Create.
gas_molesVector to receive the moles of gas components. Dimension of the vector must be ngas_comps times nxyz, where, ngas_comps is the result of RM_GetGasComponentsCount, and nxyz is the number of user grid cells (RM_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 RM_DecodeError).
See also
RM_FindComponents, RM_GetGasComponentsCount, RM_GetGasCompPressures, RM_GetGasCompPhi, RM_GetGasPhaseVolume, RM_SetGasCompMoles, RM_SetGasPhaseVolume.
C Example:
ngas_comps = RM_GetGasComponentsCount();
gas_moles = (double *) malloc((size_t) (ngas_comps * nxyz * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetGasCompMoles(id, gas_moles);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetGasComponentsCount()

IRM_DLL_EXPORT int RM_GetGasComponentsCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of gas phase components in the initial-phreeqc module.
See also
RM_FindComponents, RM_GetGasComponentsName.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -gases\n");
for (i = 0; i < RM_GetGasComponentsCount(id); i++)
{
status = RM_GetGasComponentsName(id, i, line1, 100);
sprintf(line, "%4s%20s\n", "    ", line1);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetGasComponentsName()

IRM_DLL_EXPORT IRM_RESULT RM_GetGasComponentsName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves an item from the gas components list. The list includes all gas components included in any GAS_PHASE definitions in the initial-phreeqc module. RM_FindComponents must be called before RM_GetGasComponentsName. This method may be useful when generating selected output definitions related to gas phases.

Parameters
idThe instance id returned from RM_Create.
numThe number of the gas component name to be retrieved. (0 based index.)
nameThe gas component name at number num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetGasComponentsCount.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -gases\n");
for (i = 0; i < RM_GetGasComponentsCount(id); i++)
{
status = RM_GetGasComponentsName(id, i, line1, 100);
sprintf(line, "%4s%20s\n", "    ", line1);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetGasCompPhi()

IRM_DLL_EXPORT IRM_RESULT RM_GetGasCompPhi ( int  id,
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 of a gas component is equal to its pressure times the fugacity coefficient.

Parameters
idThe instance id returned from RM_Create.
gas_phiVector to receive the fugacity coefficients of gas components. Dimension of the vector must be ngas_comps times nxyz, where, ngas_comps is the result of RM_GetGasComponentsCount, and nxyz is the number of user grid cells (RM_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 RM_DecodeError).
See also
RM_FindComponents, RM_GetGasComponentsCount, RM_GetGasCompMoles, RM_GetGasCompPressures, RM_GetGasPhaseVolume, RM_SetGasCompMoles, RM_SetGasPhaseVolume.
C Example:
ngas_comps = RM_GetGasComponentsCount();
gas_phi = (double *) malloc((size_t) (ngas_comps * nxyz * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetGasCompPhi(id, gas_phi);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetGasCompPressures()

IRM_DLL_EXPORT IRM_RESULT RM_GetGasCompPressures ( int  id,
double *  gas_pressure 
)

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

Parameters
idThe instance id returned from RM_Create.
gas_pressureVector to receive the pressures of gas components. Dimension of the vector must be ngas_comps times nxyz, where, ngas_comps is the result of RM_GetGasComponentsCount, and nxyz is the number of user grid cells (RM_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 RM_DecodeError).
See also
RM_FindComponents, RM_GetGasComponentsCount, RM_GetGasCompMoles, RM_GetGasCompPhi, RM_GetGasPhaseVolume, RM_SetGasCompMoles, RM_SetGasPhaseVolume.
C Example:
ngas_comps = RM_GetGasComponentsCount();
gas_pressure = (double *) malloc((size_t) (ngas_comps * nxyz * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetGasCompPressures(id, gas_pressure);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetGasPhaseVolume()

IRM_DLL_EXPORT IRM_RESULT RM_GetGasPhaseVolume ( int  id,
double *  gas_volume 
)

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

Parameters
idThe instance id returned from RM_Create.
gas_volumeArray to receive the gas phase volumes. Dimension of the vector must be nxyz, where, nxyz is the number of user grid cells (RM_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 RM_DecodeError).
See also
RM_FindComponents, RM_GetGasComponentsCount, RM_GetGasCompMoles, RM_GetGasCompPressures, RM_GetGasCompPhi, RM_SetGasCompMoles, RM_SetGasPhaseVolume.
C Example:
gas_volume = (double *) malloc((size_t) (nxyz * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetGasPhaseVolume(id, gas_volume);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetGfw()

IRM_DLL_EXPORT IRM_RESULT RM_GetGfw ( int  id,
double *  gfw 
)

Returns the gram formula weights (g/mol) for the components in the reaction-module component list.

Parameters
idThe instance id returned from RM_Create.
gfwArray to receive the gram formula weights. Dimension of the array is ncomps, where ncomps is the number of components in the component list.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents,
RM_GetComponent, RM_GetComponentCount.
C Example:
ncomps = RM_FindComponents(id);
components = (char **) malloc((size_t) (ncomps * sizeof(char *)));
gfw = (double *) malloc((size_t) (ncomps * sizeof(double)));
status = RM_GetGfw(id, gfw);
for (i = 0; i < ncomps; i++)
{
  components[i] = (char *) malloc((size_t) (100 * sizeof(char *)));
  status = RM_GetComponent(id, i, components[i], 100);
  sprintf(str,"%10s    %10.3f\n", components[i], gfw[i]);
  status = RM_OutputMessage(id, str);
}
MPI:
Called by root.

◆ RM_GetGridCellCount()

IRM_DLL_EXPORT int RM_GetGridCellCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Numberof grid cells in the user's model, negative is failure (See RM_DecodeError).
See also

RM_Create,
RM_CreateMapping.
C Example:
nxyz = RM_GetGridCellCount(id);
sprintf(str1, "Number of grid cells in the user's model: %d\n", nxyz);
status = RM_OutputMessage(id, str1);
MPI:
Called by root and (or) workers.

◆ RM_GetGridCellCountYAML()

IRM_DLL_EXPORT int RM_GetGridCellCountYAML ( const char *  config_file)

GetGridCellCountYAML will read a YAML file and extract the value of GridCellCount, which can be used to construct a PhreeqcRM instance. RM_Create requires a value for the number of cells. If a GUI or preprocessor is used to write a YAML file to initialize PhreeqcRM, the number of cells can be written to the YAML file and extracted with this method.

Parameters
config_fileString containing the YAML file name.
Return values
Numberof grid cells specified in the YAML file; returns zero if GridCellCount is not defined.
See also
RM_Create, and RM_InitializeYAML.
C Example:
int nthreads = 0;
int nxyz;
nxyz = RM_GetGridCellCountYAML("myfile.yaml");
status = RM_Create(nxyz, nthreads);
status = RM_InitializeYAML("myfile.yaml");
Sequence:
Called before RM_Create.

◆ RM_GetIPhreeqcId()

IRM_DLL_EXPORT int RM_GetIPhreeqcId ( int  id,
int  i 
)

Returns an IPhreeqc id for 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 (RM_Create). The number of threads can be determined by RM_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
idThe instance id returned from RM_Create.
iThe number of the IPhreeqc instance to be retrieved (0 based).
Return values
IPhreeqcid for the ith IPhreeqc instance, negative is failure (See RM_DecodeError).
See also

RM_Create, RM_GetThreadCount. See IPhreeqc documentation for descriptions of IPhreeqc methods.
C Example:
// Utility pointer is worker number nthreads + 1
iphreeqc_id1 = RM_GetIPhreeqcId(id, RM_GetThreadCount(id) + 1);
MPI:
Called by root and (or) workers.

◆ RM_GetIthConcentration()

IRM_DLL_EXPORT IRM_RESULT RM_GetIthConcentration ( int  id,
int  i,
double *  c 
)

Transfer the concentration from each cell for one component to the array given in the argument list (c). The concentrations are those resulting from the last call to RM_RunCells. Units of concentration for c are defined by RM_SetUnitsSolution.

Parameters
idThe instance id returned from RM_Create.
iZero-based index for the component to retrieve. Indices refer to the order produced by RM_GetComponent. The total number of components is given by RM_GetComponentCount.
cAllocated array to receive the component concentrations. Dimension of the array must be at least nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetComponent, RM_GetComponentCount, RM_GetConcentrations.
C Example:
c = (double*)malloc(nxyz * sizeof(double));
status = RM_RunCells(id);
status = RM_GetIthConcentration(id, 0, c)
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetIthSpeciesConcentration()

IRM_DLL_EXPORT IRM_RESULT RM_GetIthSpeciesConcentration ( int  id,
int  i,
double *  c 
)

Transfer the concentrations for one species from each cell to the array given in the argument list (c). The concentrations are those resulting from the last call to RM_RunCells. Units of concentration for c are mol/L. To retrieve species concentrations, RM_SetSpeciesSaveOn must be set to 1. This method is for use with multicomponent diffusion calculations.

Parameters
idThe instance id returned from RM_Create.
iZero-based index for the species to retrieve. Indices refer to the order given by RM_GetSpeciesName. The total number of species is given by RM_GetSpeciesCount.
cAllocated array to receive the species concentrations. Dimension of the array must be at least nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSpeciesCount, RM_GetSpeciesName, RM_GetSpeciesConcentrations, RM_SetSpeciesSaveOn.
C Example:
c = (double*)malloc(nxyz*sizeof(double));
status = RM_RunCells(id);
status = RM_GetIthSpeciesConcentration(id, 0, c);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetKineticReactionsCount()

IRM_DLL_EXPORT int RM_GetKineticReactionsCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of kinetic reactions in the initial-phreeqc module.
See also
RM_FindComponents, RM_GetKineticReactionsName.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -kinetics\n");
for (i = 0; i < RM_GetKineticReactionsCount(id); i++)
{
status = RM_GetKineticReactionsName(id, i, line1, 100);
sprintf(line, "%4s%20s\n", "    ", line1);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetKineticReactionsName()

IRM_DLL_EXPORT IRM_RESULT RM_GetKineticReactionsName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves an item from the kinetic reactions list. The list includes all kinetic reactions included in any KINETICS definitions in the initial-phreeqc module. RM_FindComponents must be called before RM_GetKineticReactionsName. This method may be useful when generating selected output definitions related to kinetic reactions.

Parameters
idThe instance id returned from RM_Create.
numThe number of the kinetic reaction name to be retrieved. (0 based index.)
nameThe kinetic reaction name at number num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetKineticReactionsCount.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -kinetics\n");
for (i = 0; i < RM_GetKineticReactionsCount(id); i++)
{
status = RM_GetKineticReactionsName(id, i, line1, 100);
sprintf(line, "%4s%20s\n", "    ", line1);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetMpiMyself()

IRM_DLL_EXPORT int RM_GetMpiMyself ( int  id)

Returns the MPI task number. For the OPENMP version, the task number is always zero and the result of RM_GetMpiTasks is one. For the MPI version, the root task number is zero, and all workers have a task number greater than zero. The number of tasks can be obtained with RM_GetMpiTasks. The number of tasks and computer hosts are 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 (RM_Create).

Parameters
idThe instance id returned from RM_Create.
Return values
TheMPI task number for a process, negative is failure (See RM_DecodeError).
See also

RM_GetMpiTasks.
C Example:
sprintf(str1, "MPI task number:  %d\n", RM_GetMpiMyself(id));
status = RM_OutputMessage(id, str1);
MPI:
Called by root and (or) workers.

◆ RM_GetMpiTasks()

IRM_DLL_EXPORT int RM_GetMpiTasks ( int  id)

Returns the number of MPI processes (tasks) assigned to the reaction module. For the OPENMP version, the number of tasks is always one (although there may be multiple threads, RM_GetThreadCount), and the task number returned by RM_GetMpiMyself is zero. For the MPI version, the number of tasks and computer hosts are determined at run time by the mpiexec command. An MPI communicator is used in constructing reaction modules for MPI. The communicator may define a subset of the total number of MPI processes. The root task number is zero, and all workers have a task number greater than zero.

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of MPI processes assigned to the reaction module, negative is failure (See RM_DecodeError).
See also

RM_Create,
RM_GetMpiMyself.
C Example:
mpi_tasks = RM_GetMpiTasks(id);
sprintf(str1, "Number of MPI processes: %d\n", mpi_tasks);
status = RM_OutputMessage(id, str1);
MPI:
Called by root and (or) workers.

◆ RM_GetNthSelectedOutputUserNumber()

IRM_DLL_EXPORT int RM_GetNthSelectedOutputUserNumber ( int  id,
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 RM_GetSelectedOutputCount. To cycle through all of the definitions, RM_GetNthSelectedOutputUserNumber can be used to identify the user number for each selected-output definition in sequence. RM_SetCurrentSelectedOutputUserNumber is then used to select that user number for selected-output processing.

Parameters
idThe instance id returned from RM_Create.
nThe sequence number of the selected-output definition for which the user number will be returned. C, 0 based.
Return values
Theuser number of the nth selected-output definition, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputCount, RM_GetSelectedOutputHeading, RM_GetSelectedOutputRowCount, RM_SetCurrentSelectedOutputUserNumber, RM_SetNthSelectedOutput, RM_SetSelectedOutputOn.
C Example:
for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
{
  n_user = RM_GetNthSelectedOutputUserNumber(id, isel);
  status = RM_SetCurrentSelectedOutputUserNumber(id, n_user);
  fprintf(stderr, "Selected output sequence number: %d\n", isel);
  fprintf(stderr, "Selected output user number:     %d\n", n_user);
  col = RM_GetSelectedOutputColumnCount(id);
  selected_out = (double *) malloc((size_t) (col * nxyz * sizeof(double)));
  status = RM_GetSelectedOutput(id, selected_out);
  // Process results here
  free(selected_out);
}
MPI:
Called by root.

◆ RM_GetPorosity()

IRM_DLL_EXPORT IRM_RESULT RM_GetPorosity ( int  id,
double *  porosity 
)

Transfer current porosities to the array given in the argument list (porosity). Porosity is not changed by PhreeqcRM; the values are either the default values or the values set by the last call to RM_SetPorosity.

Parameters
idThe instance id returned from RM_Create.
porosityArray to receive the porosities. Dimension of the array must be nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
C Example:
porosity = (double*)malloc(nxyz*sizeof(double));
status = RM_GetPorosity(id, porosity);
MPI:
Called by root.

◆ RM_GetPressure()

IRM_DLL_EXPORT IRM_RESULT RM_GetPressure ( int  id,
double *  pressure 
)

Transfer current pressures to the array given in the argument list (pressure). Pressure is not usually calculated by PhreeqcRM; the values are either the default values or the values set by the last call to RM_SetPressure. Pressures can be calculated by PhreeqcRM if a fixed-volume GAS_PHASE is used.

Parameters
idThe instance id returned from RM_Create.
pressureArray to receive the porosities. Dimension of the array must be nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
C Example:
pressure = (double*)malloc(nxyz*sizeof(double));
status = RM_GetPressure(id, pressure);
MPI:
Called by root.

◆ RM_GetSaturation()

IRM_DLL_EXPORT IRM_RESULT RM_GetSaturation ( int  id,
double *  sat_calc 
)

Deprecated equivalent of RM_GetSaturationCalculated.

◆ RM_GetSaturationCalculated()

IRM_DLL_EXPORT IRM_RESULT RM_GetSaturationCalculated ( int  id,
double *  sat_calc 
)

Returns a vector of saturations (sat) as calculated by the reaction module. This method always returns solution_volume/(rv * porosity); the method RM_SetSaturationUser has no effect on the values returned. 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 (RM_GetSolutionVolume) divided by the product of representative volume (RM_SetRepresentativeVolume) and the porosity (RM_SetPorosity). The cell saturation returned by RM_GetSaturationCalculated may be less than or greater than the saturation set by the transport code (RM_SetSaturationUser 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
idThe instance id returned from RM_Create.
sat_calcVector to receive the saturations. Dimension of the array is set to nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_GetSolutionVolume, RM_SetPorosity, RM_SetRepresentativeVolume, RM_SetSaturationUser.
C Example:
sat_calc = (double *) malloc((size_t) (nxyz * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetSaturationCalculated(id, sat_calc);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetSelectedOutput()

IRM_DLL_EXPORT IRM_RESULT RM_GetSelectedOutput ( int  id,
double *  so 
)

Populates an array with values from the current selected-output definition. RM_SetCurrentSelectedOutputUserNumber determines which of the selected-output definitions is used to populate the array.

Parameters
idThe instance id returned from RM_Create.
soAn array to contain the selected-output values. Size of the array is nxyz x col, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount), and col is the number of columns in the selected-output definition (RM_GetSelectedOutputColumnCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputCount, RM_GetSelectedOutputHeading, RM_GetSelectedOutputRowCount, RM_SetCurrentSelectedOutputUserNumber, RM_SetNthSelectedOutput, RM_SetSelectedOutputOn.
C Example:
for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
{
  n_user = RM_GetNthSelectedOutputUserNumber(id, isel);
  status = RM_SetCurrentSelectedOutputUserNumber(id, n_user);
  col = RM_GetSelectedOutputColumnCount(id);
  selected_out = (double *) malloc((size_t) (col * nxyz * sizeof(double)));
  status = RM_GetSelectedOutput(id, selected_out);
  // Process results here
  free(selected_out);
}
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetSelectedOutputColumnCount()

IRM_DLL_EXPORT int RM_GetSelectedOutputColumnCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Numberof columns in the current selected-output definition, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputCount, RM_GetSelectedOutputHeading, RM_GetSelectedOutputRowCount, RM_SetCurrentSelectedOutputUserNumber, RM_SetNthSelectedOutput, RM_SetSelectedOutputOn.
C Example:
for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
{
  n_user = RM_GetNthSelectedOutputUserNumber(id, isel);
  status = RM_SetCurrentSelectedOutputUserNumber(id, n_user);
  col = RM_GetSelectedOutputColumnCount(id);
  selected_out = (double *) malloc((size_t) (col * nxyz * sizeof(double)));
  status = RM_GetSelectedOutput(id, selected_out);
  // Process results here
  free(selected_out);
}
MPI:
Called by root.

◆ RM_GetSelectedOutputCount()

IRM_DLL_EXPORT int RM_GetSelectedOutputCount ( int  id)

Returns the number of selected-output definitions. RM_SetCurrentSelectedOutputUserNumber determines which of the selected-output definitions is used.

Parameters
idThe instance id returned from RM_Create.
Return values
Numberof selected-output definitions, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputHeading, RM_GetSelectedOutputRowCount, RM_SetCurrentSelectedOutputUserNumber, RM_SetNthSelectedOutput, RM_SetSelectedOutputOn.
C Example:
for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
{
  n_user = RM_GetNthSelectedOutputUserNumber(id, isel);
  status = RM_SetCurrentSelectedOutputUserNumber(id, n_user);
  col = RM_GetSelectedOutputColumnCount(id);
  selected_out = (double *) malloc((size_t) (col * nxyz * sizeof(double)));
  status = RM_GetSelectedOutput(id, selected_out);
  // Process results here
  free(selected_out);
}
MPI:
Called by root.

◆ RM_GetSelectedOutputHeading()

IRM_DLL_EXPORT IRM_RESULT RM_GetSelectedOutputHeading ( int  id,
int  icol,
char *  heading,
int  length 
)

Returns a selected output heading. The number of headings is determined by RM_GetSelectedOutputColumnCount. RM_SetCurrentSelectedOutputUserNumber determines which of the selected-output definitions is used.

Parameters
idThe instance id returned from RM_Create.
icolThe sequence number of the heading to be retrieved. C, 0 based.
headingA string buffer to receive the heading.
lengthThe maximum number of characters that can be written to the string buffer.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputCount, RM_GetSelectedOutputRowCount, RM_SetCurrentSelectedOutputUserNumber, RM_SetNthSelectedOutput, RM_SetSelectedOutputOn.
C Example:
char heading[100];
for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
{
  n_user = RM_GetNthSelectedOutputUserNumber(id, isel);
  status = RM_SetCurrentSelectedOutputUserNumber(id, n_user);
  col = RM_GetSelectedOutputColumnCount(id);
  for (j = 0; j < col; j++)
  {
    status = RM_GetSelectedOutputHeading(id, j, heading, 100);
    fprintf(stderr, "          %2d %10s\n", j, heading);
  }
}
MPI:
Called by root.

◆ RM_GetSelectedOutputRowCount()

IRM_DLL_EXPORT int RM_GetSelectedOutputRowCount ( int  id)

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, and is equal to RM_GetGridCellCount.

Parameters
idThe instance id returned from RM_Create.
Return values
Numberof rows in the current selected-output definition, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputCount, RM_GetSelectedOutputHeading, RM_SetCurrentSelectedOutputUserNumber, RM_SetNthSelectedOutput, RM_SetSelectedOutputOn.
C Example:
for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
{
  n_user = RM_GetNthSelectedOutputUserNumber(id, isel);
  status = RM_SetCurrentSelectedOutputUserNumber(id, n_user);
  col = RM_GetSelectedOutputColumnCount(id);
  selected_out = (double *) malloc((size_t) (col * nxyz * sizeof(double)));
  status = RM_GetSelectedOutput(id, selected_out);
  // Print results
  for (i = 0; i < RM_GetSelectedOutputRowCount(id)/2; i++)
  {
    fprintf(stderr, "Cell number %d\n", i);
    fprintf(stderr, "     Selected output: \n");
    for (j = 0; j < col; j++)
    {
      status = RM_GetSelectedOutputHeading(id, j, heading, 100);
      fprintf(stderr, "          %2d %10s: %10.4f\n", j, heading, selected_out[j*nxyz + i]);
    }
  }
  free(selected_out);
}
MPI:
Called by root.

◆ RM_GetSICount()

IRM_DLL_EXPORT int RM_GetSICount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of phases in the initial-phreeqc module for which saturation indices could be calculated.
See also
RM_FindComponents, RM_GetSIName.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -saturation_indices\n");
for (i = 0; i < RM_GetSICount(id); i++)
{
status = RM_GetSIName(id, i, line1, 100);
sprintf(line, "%4s%20s\n", "    ", line1);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetSIName()

IRM_DLL_EXPORT IRM_RESULT RM_GetSIName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves an item from the list of all phases for which saturation indices can 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. RM_FindComponents must be called before RM_GetSIName. This method may be useful when generating selected output definitions related to saturation indices.

Parameters
idThe instance id returned from RM_Create.
numThe number of the saturation-index-phase name to be retrieved. (0 based index.)
nameThe saturation-index-phase name at number num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSICount.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -saturation_indices\n");
for (i = 0; i < RM_GetSICount(id); i++)
{
status = RM_GetSIName(id, i, line1, 100);
sprintf(line, "%4s%20s\n", "    ", line1);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetSolidSolutionComponentsCount()

IRM_DLL_EXPORT int RM_GetSolidSolutionComponentsCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of solid solution components in the initial-phreeqc module.
See also
RM_FindComponents, RM_GetSolidSolutionComponentsName, RM_GetSolidSolutionName.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -solid_solutions\n");
for (i = 0; i < RM_GetSolidSolutionComponentsCount(id); i++)
{
status = RM_GetSolidSolutionComponentsName(id, i, line1, 100);
status = RM_GetSolidSolutionName(id, i, line2, 100);
sprintf(line, "%4s%20s%3s%20s\n", "    ", line1, " # ", line2);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetSolidSolutionComponentsName()

IRM_DLL_EXPORT IRM_RESULT RM_GetSolidSolutionComponentsName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves an item from the solid solution components list. The list includes all solid solution components included in any SOLID_SOLUTIONS definitions in the initial-phreeqc module. RM_FindComponents must be called before RM_GetSolidSolutionComponentsName. This method may be useful when generating selected output definitions related to solid solutions.

Parameters
idThe instance id returned from RM_Create.
numThe number of the solid solution components name to be retrieved. (0 based index.)
nameThe solid solution compnent name at number num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSolidSolutionComponentsCount, RM_GetSolidSolutionName.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -solid_solutions\n");
for (i = 0; i < RM_GetSolidSolutionComponentsCount(id); i++)
{
status = RM_GetSolidSolutionComponentsName(id, i, line1, 100);
status = RM_GetSolidSolutionName(id, i, line2, 100);
sprintf(line, "%4s%20s%3s%20s\n", "    ", line1, " # ", line2);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetSolidSolutionName()

IRM_DLL_EXPORT IRM_RESULT RM_GetSolidSolutionName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves an item from the solid solution names list. The list includes solid solution names included in SOLID_SOLUTIONS definitions in the initial-phreeqc module. 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. RM_FindComponents must be called before RM_GetSolidSolutionName. This method may be useful when generating selected output definitions related to solid solutions.

Parameters
idThe instance id returned from RM_Create.
numThe number of the solid solution name to be retrieved. (0 based index.)
nameThe solid solution name at number num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSolidSolutionComponentsCount, RM_GetSolidSolutionComponentsName.
C Example:
Utilities::strcat_safe(input, MAX_LENGTH, "  -solid_solutions\n");
for (i = 0; i < RM_GetSolidSolutionComponentsCount(id); i++)
{
status = RM_GetSolidSolutionComponentsName(id, i, line1, 100);
status = RM_GetSolidSolutionName(id, i, line2, 100);
sprintf(line, "%4s%20s%3s%20s\n", "    ", line1, " # ", line2);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetSolutionVolume()

IRM_DLL_EXPORT IRM_RESULT RM_GetSolutionVolume ( int  id,
double *  vol 
)

Transfer solution volumes from the reaction cells to the array given in the argument list (vol). Solution volumes are those 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
idThe instance id returned from RM_Create.
volArray to receive the solution volumes. Dimension of the array is (nxyz), where nxyz is the number of user grid cells. Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_GetSaturationCalculated.
C Example:
volume = (double *) malloc((size_t) (nxyz * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetSolutionVolume(id, volume);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetSpeciesConcentrations()

IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesConcentrations ( int  id,
double *  species_conc 
)

Transfer concentrations of aqueous species to the array argument (species_conc) This method is intended for use with multicomponent-diffusion transport calculations, and RM_SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by RM_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
idThe instance id returned from RM_Create.
species_concArray to receive the aqueous species concentrations. Dimension of the array is (nxyz, nspecies), where nxyz is the number of user grid cells (RM_GetGridCellCount), and nspecies is the number of aqueous species (RM_GetSpeciesCount). Concentrations are moles per liter. Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName, RM_GetSpeciesSaveOn, RM_GetSpeciesZ,
RM_SetSpeciesSaveOn, RM_SpeciesConcentrations2Module.
C Example:
status = RM_SetSpeciesSaveOn(id, 1);
ncomps = RM_FindComponents(id);
nspecies = RM_GetSpeciesCount(id);
nxyz = RM_GetGridCellCount(id);
species_c = (double *) malloc((size_t) (nxyz * nspecies * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetSpeciesConcentrations(id, species_c);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetSpeciesCount()

IRM_DLL_EXPORT int RM_GetSpeciesCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
IRM_RESULTThe number of aqueous species, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName, RM_GetSpeciesSaveOn, RM_GetSpeciesZ, RM_SpeciesConcentrations2Module, RM_SetSpeciesSaveOn.
C Example:
status = RM_SetSpeciesSaveOn(id, 1);
ncomps = RM_FindComponents(id);
nspecies = RM_GetSpeciesCount(id);
nxyz = RM_GetGridCellCount(id);
species_c = (double *) malloc((size_t) (nxyz * nspecies * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetSpeciesConcentrations(id, species_c);
MPI:
Called by root and (or) workers.

◆ RM_GetSpeciesD25()

IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesD25 ( int  id,
double *  diffc 
)

Transfers diffusion coefficients at 25C to the array argument (diffc). This method is intended for use with multicomponent-diffusion transport calculations, and RM_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.

Parameters
idThe instance id returned from RM_Create.
diffcArray to receive the diffusion coefficients at 25 C, m^2/s. Dimension of the array is nspecies, where nspecies is is the number of aqueous species (RM_GetSpeciesCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName,
RM_GetSpeciesSaveOn, RM_GetSpeciesZ, RM_SetSpeciesSaveOn, RM_SpeciesConcentrations2Module.
C Example:
status = RM_SetSpeciesSaveOn(id, 1);
ncomps = RM_FindComponents(id);
nspecies = RM_GetSpeciesCount(id);
diffc = (double *) malloc((size_t) (nspecies * sizeof(double)));
status = RM_GetSpeciesD25(id, diffc);
MPI:
Called by root and (or) workers.

◆ RM_GetSpeciesLog10Gammas()

IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesLog10Gammas ( int  id,
double *  species_log10gammas 
)

Transfer aqueous-species log10 activity coefficients to the array argument (species_log10gammas) This method is intended for use with multicomponent-diffusion transport calculations, and RM_SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by RM_FindComponents and includes all aqueous species that can be made from the set of components.

Parameters
idThe instance id returned from RM_Create.
species_log10gammasArray to receive the aqueous species log10 activity coefficients. Dimension of the array is (nxyz, nspecies), where nxyz is the number of user grid cells (RM_GetGridCellCount), and nspecies is the number of aqueous species (RM_GetSpeciesCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName, RM_GetSpeciesSaveOn, RM_GetSpeciesZ, RM_SetSpeciesSaveOn, RM_SpeciesConcentrations2Module.
C Example:
status = RM_SetSpeciesSaveOn(id, 1);
ncomps = RM_FindComponents(id);
nspecies = RM_GetSpeciesCount(id);
nxyz = RM_GetGridCellCount(id);
species_log10gammas = (double *) malloc((size_t) (nxyz * nspecies * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetSpeciesLog10Gammas(id, species_log10gammas);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetSpeciesLog10Molalities()

IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesLog10Molalities ( int  id,
double *  species_log10molalities 
)

Transfer aqueous-species log10 molalities to the array argument (species_log10molalities) To use this method RM_SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by RM_FindComponents and includes all aqueous species that can be made from the set of components.

Parameters
idThe instance id returned from RM_Create.
species_log10molalitiesArray to receive the aqueous species log10 molalities. Dimension of the array is (nxyz, nspecies), where nxyz is the number of user grid cells (RM_GetGridCellCount), and nspecies is the number of aqueous species (RM_GetSpeciesCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesName, RM_GetSpeciesSaveOn, RM_GetSpeciesZ, RM_SetSpeciesSaveOn, RM_SpeciesConcentrations2Module.
C Example:
status = RM_SetSpeciesSaveOn(id, 1);
ncomps = RM_FindComponents(id);
nspecies = RM_GetSpeciesCount(id);
nxyz = RM_GetGridCellCount(id);
species_log10molalities = (double *) malloc((size_t) (nxyz * nspecies * sizeof(double)));
status = RM_RunCells(id);
status = RM_GetSpeciesLog10Molalities(id, species_log10molalities);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetSpeciesName()

IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesName ( int  id,
int  i,
char *  name,
int  length 
)

Transfers the name of the ith aqueous species to the character argument (name). This method is intended for use with multicomponent-diffusion transport calculations, and RM_SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by RM_FindComponents and includes all aqueous species that can be made from the set of components.

Parameters
idThe instance id returned from RM_Create.
iSequence number of the species in the species list. C, 0 based.
nameCharacter array to receive the species name.
lengthMaximum length of string that can be stored in the character array.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesSaveOn, RM_GetSpeciesZ, RM_SetSpeciesSaveOn, RM_SpeciesConcentrations2Module.
C Example:
char name[100];
...
status = RM_SetSpeciesSaveOn(id, 1);
ncomps = RM_FindComponents(id);
nspecies = RM_GetSpeciesCount(id);
for (i = 0; i < nspecies; i++)
{
  status = RM_GetSpeciesName(id, i, name, 100);
  fprintf(stderr, "%s\n", name);
}
MPI:
Called by root and (or) workers.

◆ RM_GetSpeciesSaveOn()

IRM_DLL_EXPORT int RM_GetSpeciesSaveOn ( int  id)

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 RM_GetSpeciesConcentrations, and solution compositions to be set with RM_SpeciesConcentrations2Module.

Parameters
idThe instance id returned from RM_Create.
Return values
IRM_RESULT0, species are not saved; 1, species are saved; negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName, RM_GetSpeciesZ, RM_SetSpeciesSaveOn, RM_SpeciesConcentrations2Module.
C Example:
save_on = RM_GetSpeciesSaveOn(id);
if (save_on .ne. 0)
{
  fprintf(stderr, "Reaction module is saving species concentrations\n");
}
else
{
  fprintf(stderr, "Reaction module is not saving species concentrations\n");
}
MPI:
Called by root and (or) workers.

◆ RM_GetSpeciesZ()

IRM_DLL_EXPORT IRM_RESULT RM_GetSpeciesZ ( int  id,
double *  z 
)

Transfers the charge of each aqueous species to the array argument (z). This method is intended for use with multicomponent-diffusion transport calculations, and RM_SetSpeciesSaveOn must be set to true.

Parameters
idThe instance id returned from RM_Create.
zArray that receives the charge for each aqueous species. Dimension of the array is nspecies, where nspecies is is the number of aqueous species (RM_GetSpeciesCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName, RM_GetSpeciesSaveOn, RM_SetSpeciesSaveOn, RM_SpeciesConcentrations2Module.
C Example:
status = RM_SetSpeciesSaveOn(id, 1);
ncomps = RM_FindComponents(id);
nspecies = RM_GetSpeciesCount(id);
z = (double *) malloc((size_t) (nspecies * sizeof(double)));
status = RM_GetSpeciesZ(id, z);
MPI:
Called by root and (or) workers.

◆ RM_GetStartCell()

IRM_DLL_EXPORT IRM_RESULT RM_GetStartCell ( int  id,
int *  sc 
)

Returns an array with the starting cell numbers from the range of cell numbers assigned to each worker.

Parameters
idThe instance id returned from RM_Create.
scArray to receive the starting cell numbers. Dimension of the array is the number of threads (OpenMP) or the number of processes (MPI).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_Create, RM_GetEndCell, RM_GetMpiTasks, RM_GetThreadCount.
C Example:
n = RM_GetThreadCount(id) * RM_GetMpiTasks(id);
sc = (int *) malloc((size_t) (n * sizeof(int)));
status = RM_GetStartCell(id, sc);
MPI:
Called by root and (or) workers.

◆ RM_GetSurfaceName()

IRM_DLL_EXPORT IRM_RESULT RM_GetSurfaceName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves the surface name (such as "Hfo") that corresponds with the surface species name. The lists of surface species names and surface names are the same length. RM_FindComponents must be called before RM_GetSurfaceName. This method may be useful when generating selected output definitions related to surfaces.

Parameters
idThe instance id returned from RM_Create.
numThe number of the surface name to be retrieved. (0 based index.)
nameThe surface name associated with surface species num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSurfaceSpeciesCount, RM_GetSurfaceSpeciesName, RM_GetSurfaceType.
C Example:
for (i = 0; i < RM_GetSurfaceSpeciesCount(id); i++)
{
status = RM_GetSurfaceSpeciesName(id, i, line1, 100);
status = RM_GetSurfaceType(id, i, line2, 100);
status = RM_GetSurfaceName(id, i, line3, 100);
sprintf(line, "%4s%20s%3s%20s%20s\n", "    ", line1, " # ", line2, line3);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetSurfaceSpeciesCount()

IRM_DLL_EXPORT int RM_GetSurfaceSpeciesCount ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of surface species in the initial-phreeqc module.
See also
RM_FindComponents, RM_GetSurfaceSpeciesName, RM_GetSurfaceType, RM_GetSurfaceName.
C Example:
for (i = 0; i < RM_GetSurfaceSpeciesCount(id); i++)
{
status = RM_GetSurfaceSpeciesName(id, i, line1, 100);
status = RM_GetSurfaceType(id, i, line2, 100);
status = RM_GetSurfaceName(id, i, line3, 100);
sprintf(line, "%4s%20s%3s%20s%20s\n", "    ", line1, " # ", line2, line3);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetSurfaceSpeciesName()

IRM_DLL_EXPORT IRM_RESULT RM_GetSurfaceSpeciesName ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves an item from the surface species list. The list of surface species (for example, "Hfo_wOH") is derived from the list of components (RM_FindComponents) and the list of all surface types (such as "Hfo_w") that are included in SURFACE definitions in the initial-phreeqc module. RM_FindComponents must be called before RM_GetSurfaceSpeciesName. This method may be useful when generating selected output definitions related to surfaces.

Parameters
idThe instance id returned from RM_Create.
numThe number of the surface type to be retrieved. (0 based index.)
nameThe surface species name at number num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSurfaceSpeciesCount, RM_GetSurfaceType, RM_GetSurfaceName.
C Example:
for (i = 0; i < RM_GetSurfaceSpeciesCount(id); i++)
{
status = RM_GetSurfaceSpeciesName(id, i, line1, 100);
status = RM_GetSurfaceType(id, i, line2, 100);
status = RM_GetSurfaceName(id, i, line3, 100);
sprintf(line, "%4s%20s%3s%20s%20s\n", "    ", line1, " # ", line2, line3);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetSurfaceType()

IRM_DLL_EXPORT IRM_RESULT RM_GetSurfaceType ( int  id,
int  num,
char *  name,
int  l1 
)

Retrieves the surface site type (such as "Hfo_w") that corresponds with the surface species name. The lists of surface species names and surface species types are the same length. RM_FindComponents must be called before RM_GetSurfaceType. This method may be useful when generating selected output definitions related to surfaces.

Parameters
idThe instance id returned from RM_Create.
numThe number of the surface type to be retrieved. (0 based index.)
nameThe surface type associated with surface species num.
l1The length of the maximum number of characters for name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSurfaceSpeciesCount, RM_GetSurfaceSpeciesName, RM_GetSurfaceName.
C Example:
for (i = 0; i < RM_GetSurfaceSpeciesCount(id); i++)
{
status = RM_GetSurfaceSpeciesName(id, i, line1, 100);
status = RM_GetSurfaceType(id, i, line2, 100);
status = RM_GetSurfaceName(id, i, line3, 100);
sprintf(line, "%4s%20s%3s%20s%20s\n", "    ", line1, " # ", line2, line3);
Utilities::strcat_safe(input, MAX_LENGTH, line);
}
MPI:
Called by root.

◆ RM_GetTemperature()

IRM_DLL_EXPORT IRM_RESULT RM_GetTemperature ( int  id,
double *  temperature 
)

Returns an array of temperatures (temperature) from the reaction module. Reactions do not change the temperature, so the temperatures are either the temperatures at initialization, or the values set with the last call to RM_SetTemperature.

Parameters
idThe instance id returned from RM_Create.
temperatureAllocatable array to receive the temperatures. Dimension of the array must be nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_SetTemperature.
C Example:
temperature = (double*)malloc(nxyz*sizeof(double));
status = RM_GetTemperature(id, temperature);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_GetThreadCount()

IRM_DLL_EXPORT int RM_GetThreadCount ( int  id)

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 RM_Create. For the MPI version, the number of threads is always one for each process.

Parameters
idThe instance id returned from RM_Create.
Return values
Thenumber of threads, negative is failure (See RM_DecodeError).
See also

RM_GetMpiTasks.
C Example:
sprintf(str1, "Number of threads: %d\n", RM_GetThreadCount(id));
status = RM_OutputMessage(id, str1);
MPI:
Called by root and (or) workers; result is always 1.

◆ RM_GetTime()

IRM_DLL_EXPORT double RM_GetTime ( int  id)

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 RM_SetTime.

Parameters
idThe instance id returned from RM_Create.
Return values
Thecurrent simulation time in seconds.
See also

RM_GetTimeConversion, RM_GetTimeStep, RM_SetTime, RM_SetTimeConversion, RM_SetTimeStep.
C Example:
sprintf(str, "%s%10.1f%s", "Beginning reaction calculation ",
        RM_GetTime(id) * RM_GetTimeConversion(id), " days\n");
status = RM_LogMessage(id, str);
MPI:
Called by root and (or) workers.

◆ RM_GetTimeConversion()

IRM_DLL_EXPORT double RM_GetTimeConversion ( int  id)

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 (RM_SetTimeConversion) and retrieve it with RM_GetTimeConversion. The reaction module only uses the conversion factor when printing the long version of cell chemistry (RM_SetPrintChemistryOn), which is rare. Default conversion factor is 1.0.

Parameters
idThe instance id returned from RM_Create.
Return values
Multiplierto convert seconds to another time unit.
See also

RM_GetTime, RM_GetTimeStep, RM_SetTime, RM_SetTimeConversion, RM_SetTimeStep.
C Example:
sprintf(str, "%s%10.1f%s", "Beginning reaction calculation ",
        RM_GetTime(id) * RM_GetTimeConversion(id), " days\n");
status = RM_LogMessage(id, str);
MPI:
Called by root and (or) workers.

◆ RM_GetTimeStep()

IRM_DLL_EXPORT double RM_GetTimeStep ( int  id)

Returns the current simulation time step in seconds. This is the time over which kinetic reactions are integrated in a call to RM_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 RM_SetTimeStep.

Parameters
idThe instance id returned from RM_Create.
Return values
Thecurrent simulation time step in seconds.
See also

RM_GetTime, RM_GetTimeConversion, RM_SetTime, RM_SetTimeConversion, RM_SetTimeStep.
C Example:
sprintf(str, "%s%10.1f%s", "          Time step                  ",
        RM_GetTimeStep(id) * RM_GetTimeConversion(id), " days\n");
status = RM_LogMessage(id, str);
MPI:
Called by root and (or) workers.

◆ RM_GetViscosity()

IRM_DLL_EXPORT IRM_RESULT RM_GetViscosity ( int  id,
double *  viscosity 
)

Transfer current viscosities to the array given in the argument list (viscosity).

Parameters
idThe instance id returned from RM_Create.
viscosityAllocated array to receive the viscosities. Dimension of the array must be nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are set to 1e30.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
C Example:
viscosity = (double*)malloc(nxyz*sizeof(double));
status = RM_GetViscosity(id, viscosity);
MPI:
Called by root.

◆ RM_InitialEquilibriumPhases2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialEquilibriumPhases2Module ( int  id,
int *  equilibrium_phases 
)

Transfer EQUILIBRIUM_PHASES definitions from the InitialPhreeqc instance to the reaction-module workers. equilibrium_phases is used to select EQUILIBRIUM_PHASES definitions for each cell of the model. equilibrium_phases is dimensioned nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).

Parameters
idThe instance id returned from RM_Create.
equilibrium_phasesArray of EQUILIBRIUM_PHASES index numbers that refer to definitions in the InitialPhreeqc instance. Size is nxyz. Negative values are ignored, resulting in no transfer of an EQUILIBRIUM_PHASES definition for that cell. (Note that an EQUILIBRIUM_PHASES definition for a cell could be defined by other calls to RM_InitialEquilibriumPhases2Module, RM_InitialPhreeqc2Module, or RM_InitialPhreeqcCell2Module.)
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialSolutions2Module, RM_InitialExchanges2Module, RM_InitialGasPhases2Module, RM_InitialKinetics2Module, RM_InitialSolidSolutions2Module, RM_InitialSurfaces2Module, RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module.
C Example:
equilibrium_phases = (double*)malloc(nxyz*sizeof(double));
for (i=0; i < nxyz; i++) equilibrium_phases[i] = 1;
status = RM_InitialEquilibriumPhases2Module(id, equilibrium_phases);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitialExchanges2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialExchanges2Module ( int  id,
int *  exchanges 
)

Transfer EXCHANGE definitions from the InitialPhreeqc instance to the reaction-module workers. exchanges is used to select EXCHANGE definitions for each cell of the model. exchanges is dimensioned nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).

Parameters
idThe instance id returned from RM_Create.
exchangesVector of EXCHANGE index numbers that refer to definitions in the InitialPhreeqc instance. Size is nxyz. Negative values are ignored, resulting in no transfer of an EXCHANGE definition for that cell. (Note that an EXCHANGE definition for a cell could be defined by other calls to RM_InitialExchanges2Module, RM_InitialPhreeqc2Module, or RM_InitialPhreeqcCell2Module.)
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialSolutions2Module, RM_InitialEquilibriumPhases2Module, RM_InitialGasPhases2Module, RM_InitialKinetics2Module, RM_InitialSolidSolutions2Module, RM_InitialSurfaces2Module, RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module.
C Example:
exchanges = (double*)malloc(nxyz*sizeof(double));
for (i=0; i < nxyz; i++) exchanges[i] = 1;
status = RM_InitialExchanges2Module(id, exchanges);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitialGasPhases2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialGasPhases2Module ( int  id,
int *  gas_phases 
)

Transfer GAS_PHASE definitions from the InitialPhreeqc instance to the reaction-module workers. gas_phases is used to select GAS_PHASE definitions for each cell of the model. gas_phases is dimensioned nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).

Parameters
idThe instance id returned from RM_Create.
gas_phasesVector of GAS_PHASE index numbers that refer to definitions in the InitialPhreeqc instance.Size is nxyz. Negative values are ignored, resulting in no transfer of a GAS_PHASE definition for that cell. (Note that an GAS_PHASE definition for a cell could be defined by other calls to RM_InitialGasPhases2Module, RM_InitialPhreeqc2Module, or RM_InitialPhreeqcCell2Module.)
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialSolutions2Module, RM_InitialEquilibriumPhases2Module, RM_InitialExchanges2Module, RM_InitialKinetics2Module, RM_InitialSolidSolutions2Module, RM_InitialSurfaces2Module, RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module.
C Example:
gas_phases = (double*)malloc(nxyz*sizeof(double));
for (i=0; i < nxyz; i++) gas_phases[i] = 1;
status = RM_InitialGasPhases2Module(id, gas_phases);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitializeYAML()

IRM_DLL_EXPORT IRM_RESULT RM_InitializeYAML ( int  id,
const char *  yamlfile 
)

A YAML file can be used to initialize an instance of PhreeqcRM.

Parameters
idThe instance id returned from RM_Create.
yamlfileString containing the YAML file name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).

The file contains a YAML map of PhreeqcRM methods and the arguments corresponding to the methods. Note that the PhreeqcRM methods do not have the "RM_" prefix and the id argument is not included. For example,

- key: LoadDatabase
  database: phreeqc.dat
- key: RunFile
  workers: true
  initial_phreeqc: true
  utility: true
  chemistry_name: advect.pqi

RM_InitializeYAML will read the YAML file and execute the specified methods with the specified arguments. Using YAML terminology, the argument(s) for a method may be a scalar, a sequence, or a map, depending if the argument is a single item, a single vector, or there are multiple arguments. In the case of a map, the name associated with each argument (for example "chemistry_name" above) is arbitrary. The names of the map keys for map arguments are not used in parsing the YAML file; only the order of the arguments is important.

The following list gives the PhreeqcRM methods that can be specified in a YAML file and the arguments that are required. The arguments are described with C++ formats, which are sufficient to identify which arguments are YAML scalars (single bool, int, double, string argument), sequences (single vector argument), or maps (multiple arguments).

CloseFiles(void);
CreateMapping(std::vector< int >& grid2chem);
DumpModule();
FindComponents();
InitialPhreeqc2Module(std::vector< int > initial_conditions1);
InitialPhreeqc2Module(std::vector< int > initial_conditions1, std::vector< int > initial_conditions2, std::vector< double > fraction1);
InitialPhreeqcCell2Module(int n, std::vector< int > cell_numbers);
LoadDatabase(std::string database);
OpenFiles(void);
OutputMessage(std::string str);
RunCells(void);
RunFile(bool workers, bool initial_phreeqc, bool utility, std::string chemistry_name);
RunString(bool workers, bool initial_phreeqc, bool utility, std::string input_string);
ScreenMessage(std::string str);
SetComponentH2O(bool tf);
SetConcentrations(std::vector< double > c);
SetCurrentSelectedOutputUserNumber(int n_user);
SetDensityUser(std::vector< double > density);
SetDumpFileName(std::string dump_name);
SetErrorHandlerMode(int mode);
SetErrorOn(bool tf);
SetFilePrefix(std::string prefix);
SetGasCompMoles(std::vector< double > gas_moles);
SetGasPhaseVolume(std::vector< double > gas_volume);
SetPartitionUZSolids(bool tf);
SetPorosity(std::vector< double > por);
SetPressure(std::vector< double > p);
SetPrintChemistryMask(std::vector< int > cell_mask);
SetPrintChemistryOn(bool workers, bool initial_phreeqc, bool utility);
SetRebalanceByCell(bool tf);
SetRebalanceFraction(double f);
SetRepresentativeVolume(std::vector< double > rv);
SetSaturationUser(std::vector< double > sat);
SetScreenOn(bool tf);
SetSelectedOutputOn(bool tf);
SetSpeciesSaveOn(bool save_on);
SetTemperature(std::vector< double > t);
SetTime(double time);
SetTimeConversion(double conv_factor);
SetTimeStep(double time_step);
SetUnitsExchange(int option);
SetUnitsGasPhase(int option);
SetUnitsKinetics(int option);
SetUnitsPPassemblage(int option);
SetUnitsSolution(int option);
SetUnitsSSassemblage(int option);
SetUnitsSurface(int option);
SpeciesConcentrations2Module(std::vector< double > species_conc);
StateSave(int istate);
StateApply(int istate);
StateDelete(int istate);
UseSolutionDensityVolume(bool tf);
WarningMessage(std::string warnstr);

C Example:
        id = RM_Create(nxyz, MPI_COMM_WORLD)
        status = RM_InitializeYAML(id, "myfile.yaml")
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitialKinetics2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialKinetics2Module ( int  id,
int *  kinetics 
)

Transfer KINETICS definitions from the InitialPhreeqc instance to the reaction-module workers. kinetics is used to select KINETICS definitions for each cell of the model. kinetics is dimensioned nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).

Parameters
idThe instance id returned from RM_Create.
kineticsArray of KINETICS index numbers that refer to definitions in the InitialPhreeqc instance. Size is nxyz. Negative values are ignored, resulting in no transfer of a KINETICS definition for that cell. (Note that an KINETICS definition for a cell could be defined by other calls to RM_InitialKinetics2Module, RM_InitialPhreeqc2Module, or RM_InitialPhreeqcCell2Module.)
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialSolutions2Module, RM_InitialEquilibriumPhases2Module, RM_InitialExchanges2Module, RM_InitialGasPhases2Module, RM_InitialSolidSolutions2Module, RM_InitialSurfaces2Module, RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module.
C Example:
kinetics = (double*)malloc(nxyz*sizeof(double));
for (i=0; i < nxyz; i++) kinetics[i] = 1;
status = RM_InitialKinetics2Module(id, kinetics);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitialPhreeqc2Concentrations()

IRM_DLL_EXPORT IRM_RESULT RM_InitialPhreeqc2Concentrations ( int  id,
double *  c,
int  n_boundary,
int *  boundary_solution1,
int *  boundary_solution2,
double *  fraction1 
)

Fills an array (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. Concentrations may be a mixture of two solutions, boundary_solution1 and boundary_solution2, with a mixing fraction for boundary_solution1 1 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. If boundary_solution2 and fraction1 are NULL, no mixing is used; concentrations are derived from boundary_solution1 only.

Parameters
idThe instance id returned from RM_Create.
cArray of concentrations extracted from the InitialPhreeqc instance. The dimension of c is en_boundary * ncomp, where ncomp is the number of components returned from RM_FindComponents or RM_GetComponentCount.
n_boundaryThe number of boundary condition solutions that need to be filled.
boundary_solution1Array of solution index numbers that refer to solutions in the InitialPhreeqc instance. Size is n_boundary.
boundary_solution2Array 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. May be NULL in C.
fraction1Fraction of boundary_solution1 that mixes with (1-fraction1) of boundary_solution2. Size is (n_boundary). May be NULL in C.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetComponentCount.
C Example:
nbound = 1;
bc1 = (int *) malloc((size_t) (nbound * sizeof(int)));
bc2 = (int *) malloc((size_t) (nbound * sizeof(int)));
bc_f1 = (double *) malloc((size_t) (nbound * sizeof(double)));
bc_conc = (double *) malloc((size_t) (ncomps * nbound * sizeof(double)));
for (i = 0; i < nbound; i++)
{
  bc1[i]          = 0;       // Solution 0 from InitialPhreeqc instance
  bc2[i]          = -1;      // no bc2 solution for mixing
  bc_f1[i]        = 1.0;     // mixing fraction for bc1
}
status = RM_InitialPhreeqc2Concentrations(id, bc_conc, nbound, bc1, bc2, bc_f1);
MPI:
Called by root.

◆ RM_InitialPhreeqc2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialPhreeqc2Module ( int  id,
int *  initial_conditions1,
int *  initial_conditions2,
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 (nxyz, 7), where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount). The dimension of 7 refers to solutions and reactants in the following order: (1) SOLUTIONS, (2) EQUILIBRIUM_PHASES, (3) EXCHANGE, (4) SURFACE, (5) GAS_PHASE, (6) SOLID_SOLUTIONS, and (7) KINETICS. In C, initial_solution1[3*nxyz + 99] = 2, indicates that cell 99 (0 based) contains the SURFACE definition with user number 2 that has been defined in the InitialPhreeqc instance (either by RM_RunFile or RM_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. In C, 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. If initials_conditions2 and fraction1 are NULL, no mixing is used, and initial conditions are derived solely from initials_conditions1.

Parameters
idThe instance id returned from RM_Create.
initial_conditions1Array of solution and reactant index numbers that refer to definitions in the InitialPhreeqc instance. Size is (nxyz,7). The order of definitions is given above. Negative values are ignored, resulting in no definition of that entity for that cell.
initial_conditions2Array 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 (nxyz,7). The order of definitions is given above. May be NULL in C; setting to NULL results in no mixing.
fraction1Fraction of initial_conditions1 that mixes with (1-fraction1) of initial_conditions2. Size is (nxyz,7). The order of definitions is given above. May be NULL in C; setting to NULL results in no mixing.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_InitialPhreeqcCell2Module.
C Example:
ic1 = (int *) malloc((size_t) (7 * nxyz * sizeof(int)));
ic2 = (int *) malloc((size_t) (7 * nxyz * sizeof(int)));
f1 = (double *) malloc((size_t) (7 * nxyz * sizeof(double)));
for (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

  ic2[i]          = -1;      // Solution none
  ic2[nxyz + i]   = -1;      // Equilibrium phases none
  ic2[2*nxyz + i] = -1;      // Exchange none
  ic2[3*nxyz + i] = -1;      // Surface none
  ic2[4*nxyz + i] = -1;      // Gas phase none
  ic2[5*nxyz + i] = -1;      // Solid solutions none
  ic2[6*nxyz + i] = -1;      // Kinetics none

  f1[i]          = 1.0;      // Mixing fraction ic1 Solution
  f1[nxyz + i]   = 1.0;      // Mixing fraction ic1 Equilibrium phases
  f1[2*nxyz + i] = 1.0;      // Mixing fraction ic1 Exchange 1
  f1[3*nxyz + i] = 1.0;      // Mixing fraction ic1 Surface
  f1[4*nxyz + i] = 1.0;      // Mixing fraction ic1 Gas phase
  f1[5*nxyz + i] = 1.0;      // Mixing fraction ic1 Solid solutions
  f1[6*nxyz + i] = 1.0;      // Mixing fraction ic1 Kinetics
}
status = RM_InitialPhreeqc2Module(id, ic1, ic2, f1);
// No mixing is defined, so the following is equivalent
status = RM_InitialPhreeqc2Module(id, ic1, NULL, NULL);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitialPhreeqc2SpeciesConcentrations()

IRM_DLL_EXPORT IRM_RESULT RM_InitialPhreeqc2SpeciesConcentrations ( int  id,
double *  species_c,
int  n_boundary,
int *  boundary_solution1,
int *  boundary_solution2,
double *  fraction1 
)

Fills an array (species_c) with aqueous species concentrations from solutions in the InitialPhreeqc instance. This method is intended for use with multicomponent-diffusion transport calculations, and RM_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 1 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. If boundary_solution2 and fraction1 are NULL, no mixing is used; concentrations are derived from boundary_solution1 only.

Parameters
idThe instance id returned from RM_Create.
species_cArray of aqueous concentrations extracted from the InitialPhreeqc instance. The dimension of species_c is n_boundary * nspecies, where nspecies is the number of aqueous species returned from RM_GetSpeciesCount.
n_boundaryThe number of boundary condition solutions that need to be filled.
boundary_solution1Array of solution index numbers that refer to solutions in the InitialPhreeqc instance. Size is n_boundary.
boundary_solution2Array 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. May be NULL in C.
fraction1Fraction of boundary_solution1 that mixes with (1-fraction1) of boundary_solution2. Size is n_boundary. May be NULL in C.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesCount, RM_SetSpeciesSaveOn.
C Example:
nbound = 1;
nspecies = RM_GetSpeciesCount(id);
bc1 = (int *) malloc((size_t) (nbound * sizeof(int)));
bc2 = (int *) malloc((size_t) (nbound * sizeof(int)));
bc_f1 = (double *) malloc((size_t) (nbound * sizeof(double)));
bc_conc = (double *) malloc((size_t) (nspecies * nbound * sizeof(double)));
for (i = 0; i < nbound; i++)
{
  bc1[i]          = 0;       // Solution 0 from InitialPhreeqc instance
  bc2[i]          = -1;      // no bc2 solution for mixing
  bc_f1[i]        = 1.0;     // mixing fraction for bc1
}
status = RM_InitialPhreeqc2SpeciesConcentrations(id, bc_conc, nbound, bc1, bc2, bc_f1);
MPI:
Called by root.

◆ RM_InitialPhreeqcCell2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialPhreeqcCell2Module ( int  id,
int  n,
int *  module_numbers,
int  dim_module_numbers 
)

A cell numbered n in the InitialPhreeqc instance is selected to populate a series of 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 module_numbers are removed before the cell definition is copied from the InitialPhreeqc instance to the workers.

Parameters
idThe instance id returned from RM_Create.
nCell number refers to a solution or MIX and associated reactants in the InitialPhreeqc instance. A negative number indicates the largest solution or MIX number in the InitialPhreeqc instance will be used.
module_numbersA list of cell numbers in the user's grid-cell numbering system that will be populated with cell n from the InitialPhreeqc instance.
dim_module_numbersThe number of cell numbers in the module_numbers list.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialPhreeqc2Module.
C Example:
module_cells = (int *) malloc((size_t) (2 * sizeof(int)));
module_cells[0] = 18;
module_cells[1] = 19;
// n will be the largest SOLUTION number in InitialPhreeqc instance
// copies solution and reactants to cells 18 and 19
status = RM_InitialPhreeqcCell2Module(id, -1, module_cells, 2);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitialSolidSolutions2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialSolidSolutions2Module ( int  id,
int *  solid_solutions 
)

Transfer SOLID_SOLUTIONS definitions from the InitialPhreeqc instance to the reaction-module workers. solid_solutions is used to select SOLID_SOLUTIONS definitions for each cell of the model. solid_solutions is dimensioned nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).

Parameters
idThe instance id returned from RM_Create.
solid_solutionsArray of SOLID_SOLUTIONS index numbers that refer to definitions in the InitialPhreeqc instance. Size is nxyz. Negative values are ignored, resulting in no transfer of a SOLID_SOLUTIONS definition for that cell. (Note that an SOLID_SOLUTIONS definition for a cell could be defined by other calls to RM_InitialSolidSolutions2Module, RM_InitialPhreeqc2Module, or RM_InitialPhreeqcCell2Module.)
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialSolutions2Module, RM_InitialEquilibriumPhases2Module, RM_InitialExchanges2Module, RM_InitialGasPhases2Module, RM_InitialKinetics2Module, RM_InitialSurfaces2Module, RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module.
C Example:
solid_solutions = (double*)malloc(nxyz*sizeof(double));
for (i=0; i < nxyz; i++) solid_solutions[i] = 1;
status = RM_InitialSolidSolutions2Module(id, solid_solutions);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitialSolutions2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialSolutions2Module ( int  id,
int *  solutions 
)

Transfer SOLUTION definitions from the InitialPhreeqc instance to the reaction- module workers. solutions is used to select SOLUTION definitions for each cell of the model. solutions is dimensioned nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).

Parameters
idThe instance id returned from RM_Create.
solutionsArray of SOLUTION index numbers that refer to definitions in the InitialPhreeqc instance. Size is nxyz. Negative values are ignored, resulting in no transfer of a SOLUTION definition for that cell. (Note that all cells must have a SOLUTION definition, which could be defined by other calls to RM_InitialSolutions2Module, RM_InitialPhreeqc2Module, or RM_InitialPhreeqcCell2Module.)
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialEquilibriumPhases2Module, RM_InitialExchanges2Module, RM_InitialGasPhases2Module, RM_InitialKinetics2Module, RM_InitialSolidSolutions2Module, RM_InitialSurfaces2Module, RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module.
C Example:
solutions = (double*)malloc(nxyz*sizeof(double));
for (i=0; i < nxyz; i++) solutions[i] = 1;
status = RM_InitialSolutions2Module(id, solutions);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_InitialSurfaces2Module()

IRM_DLL_EXPORT IRM_RESULT RM_InitialSurfaces2Module ( int  id,
int *  surfaces 
)

Transfer SURFACE definitions from the InitialPhreeqc instance to the reaction-module workers. surfaces is used to select SURFACE definitions for each cell of the model. surfaces is dimensioned nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).

Parameters
idThe instance id returned from RM_Create.
surfacesArray of SURFACE index numbers that refer to definitions in the InitialPhreeqc instance. Size is nxyz. Negative values are ignored, resulting in no transfer of a SURFACE definition for that cell. (Note that an SURFACE definition for a cell could be defined by other calls to RM_InitialSurfaces2Module, RM_InitialPhreeqc2Module, or RM_InitialPhreeqcCell2Module.)
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialSolutions2Module, RM_InitialEquilibriumPhases2Module, RM_InitialExchanges2Module, RM_InitialGasPhases2Module, RM_InitialKinetics2Module, RM_InitialSolidSolutions2Module, RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module.
C Example:
surfaces = (double*)malloc(nxyz*sizeof(double));
for (i=0; i < nxyz; i++) surfaces[i] = 1;
status = RM_InitialSurfaces2Module(id, surfaces);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_LoadDatabase()

IRM_DLL_EXPORT IRM_RESULT RM_LoadDatabase ( int  id,
const char *  db_name 
)

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
idThe instance id returned from RM_Create.
db_nameString containing the database name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_Create.
C Example:
status = RM_LoadDatabase(id, "phreeqc.dat");
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_LogMessage()

IRM_DLL_EXPORT IRM_RESULT RM_LogMessage ( int  id,
const char *  str 
)

Print a message to the log file.

Parameters
idThe instance id returned from RM_Create.
strString to be printed.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_ErrorMessage,
RM_OpenFiles, RM_OutputMessage, RM_ScreenMessage, RM_WarningMessage.
C Example:
sprintf(str, "%s%10.1f%s", "Beginning transport calculation      ",
        RM_GetTime(id) * RM_GetTimeConversion(id), " days\n");
status = RM_LogMessage(id, str);
MPI:
Called by root.

◆ RM_MpiWorker()

IRM_DLL_EXPORT IRM_RESULT RM_MpiWorker ( int  id)

MPI only. Workers (processes with RM_GetMpiMyself > 0) must call RM_MpiWorker to be able to respond to messages from the root to accept data, perform calculations, and (or) return data. RM_MpiWorker contains a loop that reads a message from root, performs a task, and waits for another message from root. RM_SetConcentrations, RM_RunCells, and RM_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 RM_MpiWorker" in the MPI section of the method documentation. The workers will continue to respond to messages from root until root calls RM_MpiWorkerBreak.

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

Parameters
idThe instance id returned from RM_Create.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError). RM_MpiWorker returns a value only when RM_MpiWorkerBreak is called by root.
See also

RM_MpiWorkerBreak, RM_SetMpiWorkerCallback, RM_SetMpiWorkerCallbackCookie.
C Example:
status = RM_MpiWorker(id);
MPI:
Called by all workers.

◆ RM_MpiWorkerBreak()

IRM_DLL_EXPORT IRM_RESULT RM_MpiWorkerBreak ( int  id)

MPI only. This method is called by root to force workers (processes with RM_GetMpiMyself > 0) to return from a call to RM_MpiWorker. RM_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 RM_MpiWorker" in the MPI section of the method documentation. The workers will continue to respond to messages from root until root calls RM_MpiWorkerBreak.

Parameters
idThe instance id returned from RM_Create.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_MpiWorker, RM_SetMpiWorkerCallback, RM_SetMpiWorkerCallbackCookie.
C Example:
status = RM_MpiWorkerBreak(id);
MPI:
Called by root.

◆ RM_OpenFiles()

IRM_DLL_EXPORT IRM_RESULT RM_OpenFiles ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_CloseFiles,
RM_ErrorMessage, RM_GetFilePrefix, RM_LogMessage, RM_OutputMessage,
RM_SetFilePrefix, RM_WarningMessage.
C Example:
status = RM_SetFilePrefix(id, "Advect_c");
status = RM_OpenFiles(id);
MPI:
Called by root.

◆ RM_OutputMessage()

IRM_DLL_EXPORT IRM_RESULT RM_OutputMessage ( int  id,
const char *  str 
)

Print a message to the output file.

Parameters
idThe instance id returned from RM_Create.
strString to be printed.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_ErrorMessage, RM_LogMessage, RM_ScreenMessage, RM_WarningMessage.
C Example:
sprintf(str1, "Number of threads:                                %d\n", RM_GetThreadCount(id));
status = RM_OutputMessage(id, str1);
sprintf(str1, "Number of MPI processes:                          %d\n", RM_GetMpiTasks(id));
status = RM_OutputMessage(id, str1);
sprintf(str1, "MPI task number:                                  %d\n", RM_GetMpiMyself(id));
status = RM_OutputMessage(id, str1);
status = RM_GetFilePrefix(id, str, 100);
sprintf(str1, "File prefix:                                      %s\n", str);
status = RM_OutputMessage(id, str1);
sprintf(str1, "Number of grid cells in the user's model:         %d\n", RM_GetGridCellCount(id));
status = RM_OutputMessage(id, str1);
sprintf(str1, "Number of chemistry cells in the reaction module: %d\n", RM_GetChemistryCellCount(id));
status = RM_OutputMessage(id, str1);
sprintf(str1, "Number of components for transport:               %d\n", RM_GetComponentCount(id));
status = RM_OutputMessage(id, str1);
MPI:
Called by root.

◆ RM_RunCells()

IRM_DLL_EXPORT IRM_RESULT RM_RunCells ( int  id)

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

Parameters
idThe instance id returned from RM_Create.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetConcentrations,
RM_SetPorosity, RM_SetPressure, RM_SetSaturationUser, RM_SetTemperature, RM_SetTimeStep.
C Example:
status = RM_SetPorosity(id, por);              // If porosity changes
status = RM_SetSaturationUser(id, sat);            // If saturation changes
status = RM_SetTemperature(id, temperature);   // If temperature changes
status = RM_SetPressure(id, pressure);         // If pressure changes
status = RM_SetConcentrations(id, c);          // Transported concentrations
status = RM_SetTimeStep(id, time_step);        // Time step for kinetic reactions
status = RM_RunCells(id);
status = RM_GetConcentrations(id, c);          // Concentrations after reaction
status = RM_GetDensityCalculated(id, density);           // Density after reaction
status = RM_GetSolutionVolume(id, volume);     // Solution volume after reaction
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_RunFile()

IRM_DLL_EXPORT IRM_RESULT RM_RunFile ( int  id,
int  workers,
int  initial_phreeqc,
int  utility,
const char *  chem_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
idThe instance id returned from RM_Create.
workers1, the workers will run the file; 0, the workers will not run the file.
initial_phreeqc1, the InitialPhreeqc instance will run the file; 0, the InitialPhreeqc will not run the file.
utility1, the Utility instance will run the file; 0, the Utility instance will not run the file.
chem_nameName of the file to run.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_RunString.
C Example:
status = RM_RunFile(id, 1, 1, 1, "advect.pqi");
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_RunString()

IRM_DLL_EXPORT IRM_RESULT RM_RunString ( int  id,
int  workers,
int  initial_phreeqc,
int  utility,
const char *  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
idThe instance id returned from RM_Create.
workers1, the workers will run the string; 0, the workers will not run the string.
initial_phreeqc1, the InitialPhreeqc instance will run the string; 0, the InitialPhreeqc will not run the string.
utility1, the Utility instance will run the string; 0, the Utility instance will not run the string.
input_stringString containing PHREEQC input.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_RunFile.
C Example:
Utilities::strcpy_safe(str, MAX_LENGTH, "DELETE; -all");
status = RM_RunString(id, 1, 0, 1, str);    // workers, initial_phreeqc, utility
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_ScreenMessage()

IRM_DLL_EXPORT IRM_RESULT RM_ScreenMessage ( int  id,
const char *  str 
)

Print message to the screen.

Parameters
idThe instance id returned from RM_Create.
strString to be printed.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_ErrorMessage,
RM_LogMessage, RM_OutputMessage, RM_WarningMessage.
C Example:
sprintf(str, "%s%10.1f%s", "Beginning transport calculation      ",
        time * RM_GetTimeConversion(id), " days\n");
status = RM_ScreenMessage(id, str);
MPI:
Called by root and (or) workers.

◆ RM_SetComponentH2O()

IRM_DLL_EXPORT IRM_RESULT RM_SetComponentH2O ( int  id,
int  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. RM_SetComponentH2O must be called before RM_FindComponents.

Parameters
idThe instance id returned from RM_Create.
tf0, total H and O are included in the component list; 1, excess H, excess O, and water are included in the component list.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents.
C Example:
status = RM_SetComponentH2O(id, 0);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetConcentrations()

IRM_DLL_EXPORT IRM_RESULT RM_SetConcentrations ( int  id,
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 (RM_SetPorosity), saturation (RM_SetSaturationUser), and reference volume (RM_SetRepresentativeVolume). The moles of each component are determined by the volume of water and per liter concentrations. If concentration units (RM_SetUnitsSolution) are mass fraction, the density (as specified by RM_SetDensityUser) is used to convert from mass fraction to per mass per liter.

Parameters
idThe instance id returned from RM_Create.
cArray of component concentrations. Size of array is nxyz * ncomps, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount), and ncomps is the number of components as determined by RM_FindComponents or RM_GetComponentCount.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetDensityUser, RM_SetPorosity, RM_SetRepresentativeVolume, RM_SetSaturationUser, RM_SetUnitsSolution.
C Example:
c = (double *) malloc((size_t) (ncomps * nxyz * sizeof(double)));
...
advect_c(c, bc_conc, ncomps, nxyz, nbound);
status = RM_SetPorsity(id, por);               // If porosity changes
status = RM_SetSaturationUser(id, sat);            // If saturation changes
status = RM_SetTemperature(id, temperature);   // If temperature changes
status = RM_SetPressure(id, pressure);         // If pressure changes
status = RM_SetConcentrations(id, c);          // Transported concentrations
status = RM_SetTimeStep(id, time_step);        // Time step for kinetic reactions
status = RM_SetTime(id, time);                 // Current time
status = RM_RunCells(id);
status = RM_GetConcentrations(id, c);          // Concentrations after reaction
status = RM_GetDensityCalculated(id, density);           // Density after reaction
status = RM_GetSolutionVolume(id, volume);     // Solution volume after reaction
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetCurrentSelectedOutputUserNumber()

IRM_DLL_EXPORT IRM_RESULT RM_SetCurrentSelectedOutputUserNumber ( int  id,
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
idThe instance id returned from RM_Create.
n_userUser number of the SELECTED_OUTPUT data block that is to be used.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputCount, RM_GetSelectedOutputHeading, RM_GetSelectedOutputRowCount, RM_SetNthSelectedOutput, RM_SetSelectedOutputOn.
C Example:
for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
{
  n_user = RM_GetNthSelectedOutputUserNumber(id, isel);
  status = RM_SetCurrentSelectedOutputUserNumber(id, n_user);
  col = RM_GetSelectedOutputColumnCount(id);
  selected_out = (double *) malloc((size_t) (col * nxyz * sizeof(double)));
  status = RM_GetSelectedOutput(id, selected_out);
  // Process results here
  free(selected_out);
}
MPI:
Called by root.

◆ RM_SetDensity()

IRM_DLL_EXPORT IRM_RESULT RM_SetDensity ( int  id,
double *  density 
)

Deprecated equivalent of RM_SetDensityUser.

◆ RM_SetDensityUser()

IRM_DLL_EXPORT IRM_RESULT RM_SetDensityUser ( int  id,
double *  density 
)

Set the density for each reaction cell. These density values are used when converting from transported mass fraction concentrations (RM_SetUnitsSolution) to produce per liter concentrations during a call to RM_SetConcentrations. They are also used when converting from module concentrations to transport concentrations of mass fraction (RM_GetConcentrations), if RM_UseSolutionDensityVolume is set to false.

Parameters
idThe instance id returned from RM_Create.
densityArray of densities. Size of array is nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_GetConcentrations, RM_SetConcentrations, RM_SetUnitsSolution, RM_UseSolutionDensityVolume.
C Example:
density = (double *) malloc((size_t) (nxyz * sizeof(double)));
for (i = 0; i < nxyz; i++)
{
    density[i] = 1.0;
}
status = RM_SetDensityUser(id, density);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetDumpFileName()

IRM_DLL_EXPORT IRM_RESULT RM_SetDumpFileName ( int  id,
const char *  dump_name 
)

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

Parameters
idThe instance id returned from RM_Create.
dump_nameName of dump file.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_DumpModule.
C Example:
status = RM_SetDumpFileName(id, "advection_c.dmp");
dump_on = 1;
append = 0;
status = RM_DumpModule(id, dump_on, append);
MPI:
Called by root.

◆ RM_SetErrorHandlerMode()

IRM_DLL_EXPORT IRM_RESULT RM_SetErrorHandlerMode ( int  id,
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
idThe instance id returned from RM_Create.
modeError handling mode: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
C Example:
id = RM_Create(nxyz, nthreads);
status = RM_SetErrorHandlerMode(id, 2);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetErrorOn()

IRM_DLL_EXPORT IRM_RESULT RM_SetErrorOn ( int  id,
int  tf 
)

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

Parameters
idThe instance id returned from RM_Create.
tf1, enable error messages; 0, disable error messages. Default is 1.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_ErrorMessage, RM_ScreenMessage.
C Example:
status = RM_SetErrorOn(id, 1);
MPI:
Called by root.

◆ RM_SetFilePrefix()

IRM_DLL_EXPORT IRM_RESULT RM_SetFilePrefix ( int  id,
const char *  prefix 
)

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

Parameters
idThe instance id returned from RM_Create.
prefixPrefix used when opening the output and log files.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_CloseFiles,
RM_OpenFiles.
C Example:
status = RM_SetFilePrefix(id, "Advect_c");
status = RM_OpenFiles(id);
MPI:
Called by root.

◆ RM_SetGasCompMoles()

IRM_DLL_EXPORT IRM_RESULT RM_SetGasCompMoles ( int  id,
double *  gas_moles 
)

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

Parameters
idThe instance id returned from RM_Create.
gas_molesVector of moles of gas components. Dimension of the vector must be ngas_comps times nxyz, where, ngas_comps is the result of RM_GetGasComponentsCount, and nxyz is the number of user grid cells (RM_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 RM_DecodeError).
See also

RM_FindComponents, RM_GetGasComponentsCount, RM_GetGasCompMoles, RM_GetGasCompPressures, RM_GetGasPhaseVolume, RM_GetGasCompPhi, RM_SetGasPhaseVolume.
C Example:
ngas_comps = RM_GetGasComponentsCount();
gas_moles = (double *) malloc((size_t) (ngas_comps * nxyz * sizeof(double)));
...
status = RM_SetGasCompMoles(id, gas_moles);
status = RM_RunCells(id)
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetGasPhaseVolume()

IRM_DLL_EXPORT IRM_RESULT RM_SetGasPhaseVolume ( int  id,
double *  gas_volume 
)

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

Parameters
idThe instance id returned from RM_Create.
gas_volumeVector of volumes for each gas phase. Dimension of the vector must be nxyz, where, nxyz is the number of user grid cells (RM_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 RM_DecodeError).
See also
RM_FindComponents, RM_GetGasComponentsCount, RM_GetGasCompMoles, RM_GetGasCompPressures, RM_GetGasPhaseVolume, RM_GetGasCompPhi, RM_SetGasCompMoles.
C Example:
gas_volume = (double *) malloc((size_t) (nxyz * sizeof(double)));
...
status = RM_SetGasPhaseVolume(id, gas_moles);
status = RM_RunCells(id)
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetIthConcentration()

IRM_DLL_EXPORT IRM_RESULT RM_SetIthConcentration ( int  id,
int  i,
double *  c 
)

Transfer the concentrations for one component given by the vector c to each reaction cell. Units of concentration for c are defined by RM_SetUnitsSolution. It is required that RM_SetIthConcentration be called for each component in the system before RM_RunCells is called.

Parameters
idThe instance id returned from RM_Create.
iZero-based index for the component to transfer. Indices refer to the order produced by RM_GetComponent. The total number of components is given by RM_GetComponentCount.
cArray of concentrations to transfer to the reaction cells. Dimension of the vector is nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are ignored.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetComponentCount, RM_GetComponent, RM_SetConcentrations.
C Example:
status = RM_SetIthConcentration(id, i, c); ! repeat for all components
...
status = phreeqc_rm.RunCells(id);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetIthSpeciesConcentration()

IRM_DLL_EXPORT IRM_RESULT RM_SetIthSpeciesConcentration ( int  id,
int  i,
double *  c 
)

Transfer the concentrations for one aqueous species given by the vector c to each reaction cell. Units of concentration for c are mol/L. To set species concentrations, RM_SetSpeciesSaveOn must be set to true. It is required that RM_SetIthSpeciesConcentration be called for each aqueous species in the system before RM_RunCells is called. This method is for use with multicomponent diffusion calculations.

Parameters
idThe instance id returned from RM_Create.
iZero-based index for the species to transfer. Indices refer to the order produced by RM_GetSpeciesName. The total number of species is given by RM_GetSpeciesCount.
cArray of concentrations to transfer to the reaction cells. Dimension of the array is nxyz, where nxyz is the number of user grid cells (RM_GetGridCellCount). Values for inactive cells are ignored.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSpeciesCount, RM_GetSpeciesName, RM_SpeciesConcentrations2Module, RM_SetSpeciesSaveOn.
C Example:
status = RM_SetIthSpeciesConcentration(id, i, c); ! repeat for all species
...
status = RM_RunCells(id);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetMpiWorkerCallback()

IRM_DLL_EXPORT IRM_RESULT RM_SetMpiWorkerCallback ( int  id,
int(*)(int *x1, void *cookie)  fcn 
)

MPI only. Defines a callback function that allows additional tasks to be done by the workers. The method RM_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. RM_SetMpiWorkerCallback allows the developer to add another function that responds to additional integer messages by calling developer-defined functions corresponding to those integers. RM_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 RM_MpiWorker to wait for another message.

In C, an additional pointer can be supplied to find the data necessary to do the task. A void pointer may be set with RM_SetMpiWorkerCallbackCookie. This pointer is passed to the callback function through a void pointer argument in addition to the integer message argument. The void pointer may be a pointer to a struct that contains pointers to additional data. RM_SetMpiWorkerCallbackCookie must be called by each worker before RM_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 RM_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 RM_MpiWorker. Alternatively, it is possible for the workers to return from RM_MpiWorker by a call to RM_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 RM_MpiWorker.

Parameters
idThe instance id returned from RM_Create.
fcnA function that returns an integer and has an integer argument. C has an additional void * argument.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_MpiWorker, RM_MpiWorkerBreak, RM_SetMpiWorkerCallbackCookie.
C Example:
Code executed by root:
// root calls a function that will involve the workers
status = do_something(&comm);

Code executed by workers:
status = RM_SetMpiWorkerCallback(id, worker_tasks_c);
status = RM_SetMpiWorkerCallbackCookie(id, &comm);
status = RM_MpiWorker(id);

Code executed by root and workers:
int do_something(void *cookie)
{
    MPI_Status status;
    MPI_Comm *comm = (MPI_Comm *) cookie;
    int i, method_number, mpi_myself, mpi_tasks, worker_number;
    method_number = 1000;
    MPI_Comm_size(MPI_COMM_WORLD, &mpi_tasks);
    MPI_Comm_rank(MPI_COMM_WORLD, &mpi_myself);
    if (mpi_myself == 0)
    {
        MPI_Bcast(&method_number, 1, MPI_INT, 0, *comm);
        fprintf(stderr, "I am root.\n");
        for (i = 1; i < mpi_tasks; i++)
        {
            MPI_Recv(&worker_number, 1, MPI_INT, i, 0, *comm, &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_c(int *method_number, void * cookie)
{
    if (*method_number == 1000)
    {
        do_something(cookie);
    }
    return 0;
}
MPI:
Called by workers, before call to RM_MpiWorker.

◆ RM_SetMpiWorkerCallbackCookie()

IRM_DLL_EXPORT IRM_RESULT RM_SetMpiWorkerCallbackCookie ( int  id,
void *  cookie 
)

MPI and C only. Defines a void pointer that can be used by C functions called from the callback function (RM_SetMpiWorkerCallback) to locate data for a task. The C callback function that is registered with RM_SetMpiWorkerCallback has two arguments, an integer message to identify a task, and a void pointer. RM_SetMpiWorkerCallbackCookie sets the value of the void pointer that is passed to the callback function.

Parameters
idThe instance id returned from RM_Create.
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 RM_DecodeError).
See also

RM_MpiWorker, RM_MpiWorkerBreak, RM_SetMpiWorkerCallback.
C Example:
Code executed by root:
// root calls a function that will involve the workers
status = do_something(&comm);

Code executed by workers:
status = RM_SetMpiWorkerCallback(id, worker_tasks_c);
status = RM_SetMpiWorkerCallbackCookie(id, &comm);
status = RM_MpiWorker(id);

Code executed by root and workers:
int do_something(void *cookie)
{
    MPI_Status status;
    MPI_Comm *comm = (MPI_Comm *) cookie;
    int i, method_number, mpi_myself, mpi_tasks, worker_number;
    method_number = 1000;
    MPI_Comm_size(MPI_COMM_WORLD, &mpi_tasks);
    MPI_Comm_rank(MPI_COMM_WORLD, &mpi_myself);
    if (mpi_myself == 0)
    {
        MPI_Bcast(&method_number, 1, MPI_INT, 0, *comm);
        fprintf(stderr, "I am root.\n");
        for (i = 1; i < mpi_tasks; i++)
        {
            MPI_Recv(&worker_number, 1, MPI_INT, i, 0, *comm, &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_c(int *method_number, void * cookie)
{
    if (*method_number == 1000)
    {
        do_something(cookie);
    }
    return 0;
}
MPI:
Called by workers, before call to RM_MpiWorker.

◆ RM_SetNthSelectedOutput()

IRM_DLL_EXPORT IRM_RESULT RM_SetNthSelectedOutput ( int  id,
int  n 
)

Specify the current selected output by sequence number. The user may define multiple SELECTED_OUTPUT data blocks for the workers. A user number is specified for each data block, and the blocks are stored in user-number order. The value of the argument n selects the sequence number of the SELECTED_OUTPUT definition that will be used for selected-output operations.

Parameters
idThe instance id returned from RM_Create.
nSequence number of the SELECTED_OUTPUT data block that is to be used.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputCount, RM_GetSelectedOutputHeading, RM_GetSelectedOutputRowCount, RM_SetCurrentSelectedOutputUserNumber, RM_SetSelectedOutputOn.
C Example:
    for (isel = 0; isel < RM_GetSelectedOutputCount(id); isel++)
    {
      status = RM_SetNthSelectedOutputUser(id, isel);
      n_user = RM_GetCurrentSelectedOutputUserNumber(id);
      col = RM_GetSelectedOutputColumnCount(id);
      selected_out = (double *) malloc((size_t) (col * nxyz * sizeof(double)));
      status = RM_GetSelectedOutput(id, selected_out);
      // Process results here
      free(selected_out);
    }
    
MPI:
Called by root.

◆ RM_SetPartitionUZSolids()

IRM_DLL_EXPORT IRM_RESULT RM_SetPartitionUZSolids ( int  id,
int  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 RM_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
idThe instance id returned from RM_Create.
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 RM_DecodeError).
C Example:
status = RM_SetPartitionUZSolids(id, 0);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetPorosity()

IRM_DLL_EXPORT IRM_RESULT RM_SetPorosity ( int  id,
double *  por 
)

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

Parameters
idThe instance id returned from RM_Create.
porArray of porosities, unitless. Default is 0.1. Size of array is nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_GetSaturationCalculated, RM_SetRepresentativeVolume, RM_SetSaturationUser.
C Example:
por = (double *) malloc((size_t) (nxyz * sizeof(double)));
for (i = 0; i < nxyz; i++) por[i] = 0.2;
status = RM_SetPorosity(id, por);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetPressure()

IRM_DLL_EXPORT IRM_RESULT RM_SetPressure ( int  id,
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
idThe instance id returned from RM_Create.
pArray of pressures, in atm. Size of array is nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetTemperature.
C Example:
pressure = (double *) malloc((size_t) (nxyz * sizeof(double)));
for (i = 0; i < nxyz; i++) pressure[i] = 2.0;
status = RM_SetPressure(id, pressure);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetPrintChemistryMask()

IRM_DLL_EXPORT IRM_RESULT RM_SetPrintChemistryMask ( int  id,
int *  cell_mask 
)

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

Parameters
idThe instance id returned from RM_Create.
cell_maskArray of integers. Size of array is nxyz, where nxyz is the number of grid cells in the user's model (RM_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 RM_DecodeError).
See also

RM_SetPrintChemistryOn.
C Example:
print_chemistry_mask = (int *) malloc((size_t) (nxyz * sizeof(int)));
for (i = 0; i < nxyz/2; i++)
{
  print_chemistry_mask[i] = 1;
  print_chemistry_mask[i + nxyz/2] = 0;
}
status = RM_SetPrintChemistryMask(id, print_chemistry_mask);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetPrintChemistryOn()

IRM_DLL_EXPORT IRM_RESULT RM_SetPrintChemistryOn ( int  id,
int  workers,
int  initial_phreeqc,
int  utility 
)

Setting to enable or disable printing detailed output from reaction calculations to the output file for a set of cells defined by RM_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, RM_OpenFiles). For the worker instances, the output can be limited to a set of cells (RM_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 RM_RunFile or RM_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
idThe instance id returned from RM_Create.
workers0, disable detailed printing in the worker instances, 1, enable detailed printing in the worker instances.
initial_phreeqc0, disable detailed printing in the InitialPhreeqc instance, 1, enable detailed printing in the InitialPhreeqc instances.
utility0, disable detailed printing in the Utility instance, 1, enable detailed printing in the Utility instance.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetPrintChemistryMask.
C Example:
status = RM_SetPrintChemistryOn(id, 0, 1, 0); // workers, initial_phreeqc, utility
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetRebalanceByCell()

IRM_DLL_EXPORT IRM_RESULT RM_SetRebalanceByCell ( int  id,
int  method 
)

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 RM_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
idThe instance id returned from RM_Create.
method0, indicates average times are used in rebalancing; 1 indicates individual cell times are used in rebalancing (default).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetRebalanceFraction.
C Example:
status = RM_SetRebalanceByCell(id, 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetRebalanceFraction()

IRM_DLL_EXPORT IRM_RESULT RM_SetRebalanceFraction ( int  id,
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 RM_RunCells calculation. The rebalancing transfers cell calculations among threads or processes to try to achieve an optimum balance. RM_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
idThe instance id returned from RM_Create.
fFraction from 0.0 to 1.0.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetRebalanceByCell.
C Example:
status = RM_SetRebalanceFraction(id, 0.5);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetRepresentativeVolume()

IRM_DLL_EXPORT IRM_RESULT RM_SetRepresentativeVolume ( int  id,
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 procuct of the representative volume, the porosity (RM_SetPorosity), and the saturation (RM_SetSaturationUser). 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
idThe instance id returned from RM_Create.
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 (RM_GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetPorosity, RM_SetSaturationUser.
C Example:
double * rv;
rv = (double *) malloc((size_t) (nxyz * sizeof(double)));
for (i = 0; i < nxyz; i++) rv[i] = 1.0;
status = RM_SetRepresentativeVolume(id, rv);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetSaturation()

IRM_DLL_EXPORT IRM_RESULT RM_SetSaturation ( int  id,
double *  sat 
)

Deprecated equivalent of RM_SetSaturationUser.

◆ RM_SetSaturationUser()

IRM_DLL_EXPORT IRM_RESULT RM_SetSaturationUser ( int  id,
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 (RM_SetPorosity), saturation (RM_SetSaturationUser), and representative volume (RM_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 RM_GetDensityCalculated, RM_GetSolutionVolume, and RM_GetSaturationCalculated can be used to account for these changes in the succeeding transport calculation. RM_SetRepresentativeVolume should be called before initial conditions are defined for the reaction cells.

Parameters
idThe instance id returned from RM_Create.
satArray of saturations, unitless. Size of array is nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_GetDensityCalculated, RM_GetSaturationCalculated, RM_GetSolutionVolume, RM_SetPorosity, RM_SetRepresentativeVolume.
C Example:
sat = (double *) malloc((size_t) (nxyz * sizeof(double)));
for (i = 0; i < nxyz; i++) sat[i] = 1.0;
status = RM_SetSaturationUser(id, sat);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetScreenOn()

IRM_DLL_EXPORT IRM_RESULT RM_SetScreenOn ( int  id,
int  tf 
)

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

Parameters
idThe instance id returned from RM_Create.
tf1, enable screen messages; 0, disable screen messages. Default is 1.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_RunCells, RM_ScreenMessage.
C Example:
status = RM_SetScreenOn(id, 1);
MPI:
Called by root.

◆ RM_SetSelectedOutputOn()

IRM_DLL_EXPORT IRM_RESULT RM_SetSelectedOutputOn ( int  id,
int  selected_output 
)

Setting determines whether selected-output results are available to be retrieved with RM_GetSelectedOutput. 1 indicates that selected-output results will be accumulated during RM_RunCells and can be retrieved with RM_GetSelectedOutput; 0 indicates that selected-output results will not be accumulated during RM_RunCells.

Parameters
idThe instance id returned from RM_Create.
selected_output0, disable selected output; 1, enable selected output.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_GetCurrentSelectedOutputUserNumber, RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputColumnCount, RM_GetSelectedOutputCount, RM_GetSelectedOutputHeading, RM_GetSelectedOutputRowCount, RM_SetCurrentSelectedOutputUserNumber, RM_SetNthSelectedOutput.
C Example:
status = RM_SetSelectedOutputOn(id, 1);       // enable selected output
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetSpeciesSaveOn()

IRM_DLL_EXPORT IRM_RESULT RM_SetSpeciesSaveOn ( int  id,
int  save_on 
)

Sets the value of the species-save property. This method enables use of PhreeqcRM with multicomponent-diffusion transport calculations. By default, concentrations of aqueous species are not saved. Setting the species-save property to 1 allows aqueous species concentrations to be retrieved with RM_GetSpeciesConcentrations, and solution compositions to be set with RM_SpeciesConcentrations2Module. RM_SetSpeciesSaveOn must be called before calls to RM_FindComponents.

Parameters
idThe instance id returned from RM_Create.
save_on0, indicates species concentrations are not saved; 1, indicates species concentrations are saved.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName, RM_GetSpeciesSaveOn, RM_GetSpeciesZ, RM_SpeciesConcentrations2Module.
C Example:
status = RM_SetSpeciesSaveOn(id, 1);
MPI:
Called by root and (or) workers.

◆ RM_SetTemperature()

IRM_DLL_EXPORT IRM_RESULT RM_SetTemperature ( int  id,
double *  t 
)

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

Parameters
idThe instance id returned from RM_Create.
tArray of temperatures, in degrees C. Size of array is nxyz, where nxyz is the number of grid cells in the user's model (RM_GetGridCellCount).
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module, RM_SetPressure.
C Example:
temperature = (double *) malloc((size_t) (nxyz * sizeof(double)));
for (i = 0; i < nxyz; i++)
{
  temperature[i] = 20.0;
}
status = RM_SetTemperature(id, temperature);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetTime()

IRM_DLL_EXPORT IRM_RESULT RM_SetTime ( int  id,
double  time 
)

Set current simulation time for the reaction module.

Parameters
idThe instance id returned from RM_Create.
timeCurrent simulation time, in seconds.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetTimeConversion,
RM_SetTimeStep.
C Example:
status = RM_SetTime(id, time);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetTimeConversion()

IRM_DLL_EXPORT IRM_RESULT RM_SetTimeConversion ( int  id,
double  conv_factor 
)

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

Parameters
idThe instance id returned from RM_Create.
conv_factorFactor to convert seconds to user time units.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetTime, RM_SetTimeStep.
C Example:
status = RM_SetTimeConversion(id, 1.0 / 86400.0); // days
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetTimeStep()

IRM_DLL_EXPORT IRM_RESULT RM_SetTimeStep ( int  id,
double  time_step 
)

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

Parameters
idThe instance id returned from RM_Create.
time_stepCurrent time step, in seconds.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetTime, RM_SetTimeConversion.
C Example:
time_step = 86400.0;
status = RM_SetTimeStep(id, time_step);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetUnitsExchange()

IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsExchange ( int  id,
int  option 
)

Sets input units for exchangers. In PHREEQC input, exchangers are defined by moles of exchange sites (Mp). RM_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 (RM_SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (RM_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
idThe instance id returned from RM_Create.
optionUnits option for exchangers: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module, RM_SetPorosity, RM_SetRepresentativeVolume.
C Example:
status = RM_SetUnitsExchange(id, 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetUnitsGasPhase()

IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsGasPhase ( int  id,
int  option 
)

Set input units for gas phases. In PHREEQC input, gas phases are defined by moles of component gases (Mp). RM_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 (RM_SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (RM_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
idThe instance id returned from RM_Create.
optionUnits option for gas phases: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module, RM_SetPorosity, RM_SetRepresentativeVolume.
C Example:
status = RM_SetUnitsGasPhase(id, 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetUnitsKinetics()

IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsKinetics ( int  id,
int  option 
)

Set input units for kinetic reactants.

In PHREEQC input, kinetics are defined by moles of kinetic reactants (Mp). RM_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 (RM_SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (RM_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 (RM_SetPorosity), the saturation (RM_SetSaturationUser), and representative volume (RM_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
idThe instance id returned from RM_Create.
optionUnits option for kinetic reactants: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module, RM_SetPorosity, RM_SetRepresentativeVolume, RM_SetSaturationUser.
C Example:
status = RM_SetUnitsKinetics(id, 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetUnitsPPassemblage()

IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsPPassemblage ( int  id,
int  option 
)

Set input units for pure phase assemblages (equilibrium phases). In PHREEQC input, equilibrium phases are defined by moles of each phase (Mp). RM_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 (RM_SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (RM_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
idThe instance id returned from RM_Create.
optionUnits option for equilibrium phases: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module, RM_SetPorosity, RM_SetRepresentativeVolume.
C Example:
status = RM_SetUnitsPPassemblage(id, 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetUnitsSolution()

IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsSolution ( int  id,
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 (RM_SetPorosity), saturation (RM_SetSaturationUser), and representative volume (RM_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 (RM_SetDensityUser) 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 (RM_SetPorosity), saturation (RM_SetSaturationUser), and representative volume (RM_SetRepresentativeVolume), and the mass of solution is volume times density as defined by RM_SetDensityUser. Which option is used is determined by RM_UseSolutionDensityVolume.

Parameters
idThe instance id returned from RM_Create.
optionUnits option for solutions: 1, 2, or 3, default is 1, mg/L.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_SetDensityUser, RM_SetPorosity, RM_SetRepresentativeVolume, RM_SetSaturationUser, RM_UseSolutionDensityVolume.
C Example:
status = RM_SetUnitsSolution(id, 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetUnitsSSassemblage()

IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsSSassemblage ( int  id,
int  option 
)

Set input units for solid-solution assemblages. In PHREEQC, solid solutions are defined by moles of each component (Mp). RM_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 (RM_SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (RM_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
idThe instance id returned from RM_Create.
optionUnits option for solid solutions: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module, RM_SetPorosity, RM_SetRepresentativeVolume.
C Example:
status = RM_SetUnitsSSassemblage(id, 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SetUnitsSurface()

IRM_DLL_EXPORT IRM_RESULT RM_SetUnitsSurface ( int  id,
int  option 
)

Set input units for surfaces. In PHREEQC input, surfaces are defined by moles of surface sites (Mp). RM_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 (RM_SetRepresentativeVolume); 1, Mp is mol/L of water in the RV, Mc = Mp*P*RV, where P is porosity (RM_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
idThe instance id returned from RM_Create.
optionUnits option for surfaces: 0, 1, or 2.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_InitialPhreeqc2Module, RM_InitialPhreeqcCell2Module, RM_SetPorosity, RM_SetRepresentativeVolume.
C Example:
status = RM_SetUnitsSurface(id, 1);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_SpeciesConcentrations2Module()

IRM_DLL_EXPORT IRM_RESULT RM_SpeciesConcentrations2Module ( int  id,
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 RM_SetSpeciesSaveOn must be set to true. The list of aqueous species is determined by RM_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
idThe instance id returned from RM_Create.
species_concArray of aqueous species concentrations. Dimension of the array is (nxyz, nspecies), where nxyz is the number of user grid cells (RM_GetGridCellCount), and nspecies is the number of aqueous species (RM_GetSpeciesCount). Concentrations are moles per liter.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesLog10Gammas, RM_GetSpeciesLog10Molalities, RM_GetSpeciesName, RM_GetSpeciesSaveOn, RM_GetSpeciesZ, RM_SetSpeciesSaveOn.
C Example:
status = RM_SetSpeciesSaveOn(id, 1);
ncomps = RM_FindComponents(id);
nspecies = RM_GetSpeciesCount(id);
nxyz = RM_GetGridCellCount(id);
species_c = (double *) malloc((size_t) (nxyz * nspecies * sizeof(double)));
...
status = RM_SpeciesConcentrations2Module(id, species_c);
status = RM_RunCells(id);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_StateApply()

IRM_DLL_EXPORT IRM_RESULT RM_StateApply ( int  id,
int  istate 
)

Reset the state of the module to a state previously saved with RM_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 workersand 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
idThe instance id returned from RM_Create.
istateInteger identifying the state that is to be applied.
Return values
IRM_RESULT0 is success, negative is failure(See RM_DecodeError).
See also
RM_StateSave and RM_StateDelete.
C Example:
status = RM_StateSave(id, 1);
...
status = RM_StateApply(id, 1);
status = RM_StateDelete(id, 1);
MPI :
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_StateDelete()

IRM_DLL_EXPORT IRM_RESULT RM_StateDelete ( int  id,
int  istate 
)

Delete a state previously saved with RM_StateSave.

Parameters
idThe instance id returned from RM_Create.
istateInteger identifying the state that is to be deleted.
Return values
IRM_RESULT0 is success, negative is failure(See RM_DecodeError).
See also
RM_StateSave and ref RM_StateApply.
C Example:
status = RM_StateSave(id, 1);
...
status = RM_StateApply(id, 1);
status = RM_StateDelete(id, 1);
MPI :
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_StateSave()

IRM_DLL_EXPORT IRM_RESULT RM_StateSave ( int  id,
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 workersand the chemistry of fully or partially unsaturated cells are also saved.The state is saved in memory; use RM_DumpModule to save the state to file.PhreeqcRM can be reset to this state by using RM_StateApply. A state is identified by an integer, and multiple states can be saved.

Parameters
idThe instance id returned from RM_Create.
istateInteger identifying the state that is saved.
Return values
IRM_RESULT0 is success, negative is failure(See RM_DecodeError).
See also
RM_DumpModule, RM_StateApply, and RM_StateDelete.
C Example:
status = RM_StateSave(id, 1);
...
status = RM_StateApply(id, 1);
status = RM_StateDelete(id, 1);
MPI :
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_UseSolutionDensityVolume()

IRM_DLL_EXPORT IRM_RESULT RM_UseSolutionDensityVolume ( int  id,
int  tf 
)

Determines the volume and density to use when converting from the reaction-module concentrations to transport concentrations (RM_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 (RM_SetDensityUser) and solution volume are defined by the product of saturation (RM_SetSaturationUser), porosity (RM_SetPorosity), and representative volume (RM_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 transport units of mass fraction.

Parameters
idThe instance id returned from RM_Create.
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 RM_SetDensityUser and the volume determined by the product of RM_SetSaturationUser, RM_SetPorosity, and RM_SetRepresentativeVolume, will be used to calculate concentrations retrieved by RM_GetConcentrations.
See also

RM_GetConcentrations, RM_SetDensityUser, RM_SetPorosity, RM_SetRepresentativeVolume, RM_SetSaturationUser.
C Example:
status = RM_UseSolutionDensityVolume(id, 0);
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.

◆ RM_WarningMessage()

IRM_DLL_EXPORT IRM_RESULT RM_WarningMessage ( int  id,
const char *  warn_str 
)

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

Parameters
idThe instance id returned from RM_Create.
warn_strString to be printed.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also

RM_ErrorMessage,
RM_LogMessage,
RM_OpenFiles, RM_OutputMessage, RM_ScreenMessage.
C Example:
status = RM_WarningMessage(id, "Parameter is out of range, using default");
MPI:
Called by root and (or) workers; only root writes to the log file.