PhreeqcRM
Public Member Functions | Data Fields
PhreeqcRM Fortran Module Reference

Public Member Functions

integer function RM_Abort (id, irm_result, err_str)
 
integer function RM_CloseFiles (id)
 
integer function RM_Concentrations2Utility (id, c, n, tc, p_atm)
 
integer function RM_Create (nxyz, nthreads)
 
integer function RM_Createmapping (id, grid2chem)
 
integer function RM_DecodeError (id, e)
 
integer function RM_Destroy (id)
 
integer function RM_DumpModule (id, dump_on, append)
 
integer function RM_ErrorMessage (id, errstr)
 
integer function RM_FindComponents (id)
 
integer function RM_GetBackwardMapping (id, n, list, size)
 
integer function RM_GetChemistryCellCount (id)
 
integer function RM_GetComponent (id, num, comp_name)
 
integer function RM_GetComponentcount (id)
 
integer function RM_GetConcentrations (id, c)
 
integer function RM_GetDensity (id, density)
 
integer function RM_GetEndCell (id, ec)
 
integer function RM_GetErrorString (id, errstr)
 
integer function RM_GetErrorStringlength (id)
 
integer function RM_GetFilePrefix (id, prefix)
 
integer function RM_GetGfw (id, gfw)
 
integer function RM_GetGridCellCount (id)
 
integer function RM_GetIPhreeqcId (id, i)
 
integer function RM_GetMpiMyself (id)
 
integer function RM_GetMpiTasks (id)
 
integer function RM_GetNthSelectedOutputUserNumber (id, n)
 
integer function RM_GetSaturation (id, sat_calc)
 
integer function RM_GetSelectedOutput (id, so)
 
integer function RM_GetSelectedOutputcolumncount (id)
 
integer function RM_GetSelectedOutputcount (id)
 
integer function RM_GetSelectedOutputheading (id, icol, heading)
 
integer function RM_GetSelectedOutputrowcount (id)
 
integer function RM_GetSolutionVolume (id, vol)
 
integer function RM_GetSpeciesConcentrations (id, species_conc)
 
integer function RM_GetSpeciesCount (id)
 
integer function RM_GetSpeciesD25 (id, diffc)
 
integer function RM_GetSpeciesName (id, i, name)
 
integer function RM_GetSpeciesSaveOn (id)
 
integer function RM_GetSpeciesZ (id, z)
 
integer function RM_GetStartCell (id, sc)
 
integer function RM_GetThreadCount (id)
 
double precision function RM_GetTime (id)
 
double precision function RM_GetTimeconversion (id)
 
double precision function RM_GetTimestep (id)
 
integer function RM_InitialPhreeqc2Concentrations (id, bc_conc, n_boundary, bc1, bc2, f1)
 
integer function RM_InitialPhreeqc2Module (id, ic1, ic2, f1)
 
integer function RM_InitialPhreeqc2SpeciesConcentrations (id, bc_conc, n_boundary, bc1, bc2, f1)
 
integer function RM_InitialPhreeqcCell2Module (id, n_user, cell_numbers, n_cell)
 
integer function RM_LoadDatabase (id, db_name)
 
integer function RM_LogMessage (id, str)
 
integer function RM_MpiWorker (id)
 
integer function RM_MpiWorkerbreak (id)
 
integer function RM_OpenFiles (id)
 
integer function RM_OutputMessage (id, str)
 
integer function RM_RunCells (id)
 
integer function RM_RunFile (id, workers, initial_phreeqc, utility, chem_name)
 
integer function RM_RunString (id, initial_phreeqc, workers, utility, input_string)
 
integer function RM_ScreenMessage (id, str)
 
integer function RM_SetComponentH2O (id, tf)
 
integer function RM_SetConcentrations (id, c)
 
integer function RM_SetCurrentSelectedOutputUserNumber (id, n_user)
 
integer function RM_SetDensity (id, density)
 
integer function RM_SetDumpFileName (id, dump_name)
 
integer function RM_SetErrorHandlerMode (id, mode)
 
integer function RM_SetFilePrefix (id, prefix)
 
integer function RM_SetMpiWorkerCallback (id, fcn)
 
integer function RM_SetPartitionUZSolids (id, tf)
 
integer function RM_SetPorosity (id, por)
 
integer function RM_SetPressure (id, p)
 
integer function RM_SetPrintChemistryMask (id, cell_mask)
 
integer function RM_SetPrintChemistryOn (id, workers, initial_phreeqc, utility)
 
integer function RM_SetRebalanceByCell (id, method)
 
integer function RM_SetRebalanceFraction (id, f)
 
integer function RM_SetRepresentativeVolume (id, rv)
 
integer function RM_SetSaturation (id, sat)
 
integer function RM_SetScreenOn (id, tf)
 
integer function RM_SetSelectedOutputOn (id, tf)
 
integer function RM_SetSpeciesSaveOn (id, save_on)
 
integer function RM_SetTemperature (id, t)
 
integer function RM_SetTime (id, time)
 
integer function RM_SetTimeconversion (id, conv_factor)
 
integer function RM_SetTimestep (id, time_step)
 
integer function RM_SetUnitsExchange (id, option)
 
integer function RM_SetUnitsGasPhase (id, option)
 
integer function RM_SetUnitsKinetics (id, option)
 
integer function RM_SetUnitsPPassemblage (id, option)
 
integer function RM_SetUnitsSolution (id, option)
 
integer function RM_SetUnitsSSassemblage (id, option)
 
integer function RM_SetUnitsSurface (id, option)
 
integer function RM_SpeciesConcentrations2Module (id, species_conc)
 
integer function RM_UseSolutionDensityVolume (id, tf)
 
integer function RM_WarningMessage (id, warn_str)
 

Data Fields

logical rmf_debug =.true.
 

Detailed Description

Fortran Documentation for the geochemical reaction module PhreeqcRM.

"USE PhreeqcRM" is included in Fortran source code to define the PhreeqcRM functions. For Windows, define the module by including the file RM_interface.F90 in your project. For Linux, configure, compile, and install the PhreeqcRM library and module file. You will need installed include directory (-I) added to the project) to reference the module file. You will need to link to the library to produce the executable for your code.

Member Function/Subroutine Documentation

integer function RM_Abort ( integer, intent(in)  id,
integer, intent(in)  irm_result,
character(len=*), intent(in)  err_str 
)

Abort the program. irm_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.
irm_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.
Fortran Example:
 iphreeqc_id = RM_Concentrations2Utility(id, c_well, 1, tc, p_atm)
 string = "SELECTED_OUTPUT 5; -pH;RUN_CELLS; -cells 1"
 status = RunString(iphreeqc_id, string)
 if (status .ne. 0) status = RM_Abort(id, status, "IPhreeqc RunString failed")
 
MPI:
Called by root or workers.
integer function RM_CloseFiles ( integer, intent(in)  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
Fortran Example:
 status = RM_CloseFiles(id)
 
MPI:
Called only by root.
integer function RM_Concentrations2Utility ( integer, intent(in)  id,
double precision, dimension(:,:), intent(in)  c,
integer, intent(in)  n,
double precision, dimension(:), intent(in)  tc,
double precision, dimension(:), intent(in)  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 size is (n, ncomps) where ncomps is the number of components (RM_GetComponentcount).
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).
Fortran Example:
 allocate (c_well(1,ncomps))
 do i = 1, ncomps
   c_well(1,i) = 0.5 * c(1,i) + 0.5 * c(10,i)
 enddo
 allocate(tc(1), p_atm(1))
 tc(1) = 15.0
 p_atm(1) = 3.0
 iphreeqc_id = RM_Concentrations2Utility(id, c_well, 1, tc, p_atm)
 string = "SELECTED_OUTPUT 5; -pH; RUN_CELLS; -cells 1"
 status = RunString(iphreeqc_id, string)
 status = SetCurrentSelectedOutputUserNumber(iphreeqc_id, 5)
 status = GetSelectedOutputValue(iphreeqc_id, 1, 1, vtype, pH, svalue)
 
MPI:
Called only by root.
integer function RM_Create ( integer, intent(in)  nxyz,
integer, intent(in)  nthreads 
)

Creates a reaction module. 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
Fortran Example:
 nxyz = 40
 #ifdef USE_MPI
   id = RM_Create(nxyz, MPI_COMM_WORLD)
   call MPI_Comm_rank(MPI_COMM_WORLD, mpi_myself, status)
   if (status .ne. MPI_SUCCESS) then
     stop "Failed to get mpi_myself"
   endif
   if (mpi_myself > 0) then
     status = RM_MpiWorker(id)
     status = RM_Destroy(id)
     return
   endif
 #else
   nthreads = 3
   id = RM_Create(nxyz, nthreads)
 #endif
 
MPI:
Called by root and workers.
integer function RM_Createmapping ( integer, intent(in)  id,
integer, dimension(:), intent(in)  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 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).
Fortran Example:
 ! For demonstation, two equivalent rows by symmetry
 allocate(grid2chem(nxyz))
 do i = 1, nxyz/2
   grid2chem(i) = i - 1
   grid2chem(i+nxyz/2) = i - 1
 enddo
 status = RM_Createmapping(id, grid2chem)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_DecodeError ( integer, intent(in)  id,
integer, intent(in)  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;
 
Fortran Example:
 status = RM_Createmapping(id, grid2chem)
 if (status < 0) status = RM_DecodeError(id, status)
 
MPI:
Can be called by root and (or) workers.
integer function RM_Destroy ( integer, intent(in)  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
Fortran Example:
 status = RM_Destroy(id)
 
MPI:
Called by root and workers.
integer function RM_DumpModule ( integer, intent(in)  id,
integer, intent(in)  dump_on,
integer, intent(in)  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
Fortran Example:
      dump_on = 1
 append = 0
 status = RM_SetDumpFileName(id, "advection_f90.dmp")
 status = RM_DumpModule(id, dump_on, append)
 
MPI:
Called by root; workers must be in the loop of RM_MpiWorker.
integer function RM_ErrorMessage ( integer, intent(in)  id,
character(len=*), intent(in)  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_OpenFiles, RM_LogMessage, RM_OutputMessage, RM_ScreenMessage, RM_WarningMessage.
Fortran Example:
 status = RM_ErrorMessage(id, "Goodbye world")
 
MPI:
Called by root and (or) workers; root writes to output and log files.
integer function RM_FindComponents ( integer, intent(in)  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_SetSpeciesSaveOn, RM_GetSpeciesConcentrations, RM_SpeciesConcentrations2Module, RM_GetSpeciesCount, RM_GetSpeciesName, RM_GetSpeciesD25, RM_GetSpeciesZ, RM_SetComponentH2O.
Fortran Example:
 ! Get list of components
 ncomps = RM_FindComponents(id)
 allocate(components(ncomps))
 do i = 1, ncomps
   status = RM_GetComponent(id, i, components(i))
 enddo
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetBackwardMapping ( integer, intent(in)  id,
integer, intent(in)  n,
integer, dimension(*), intent(in)  list,
integer, intent(inout)  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) .eq. 0) then
   if (fstr(1:l) .eq. "HYDRAULIC_K") then
     my_basic_fortran_callback = K_ptr(list(1)+1) 
   endif
 endif
 
MPI:
Called by root and (or) workers.
integer function RM_GetChemistryCellCount ( integer, intent(in)  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.
Fortran Example:
 status = RM_Createmapping(id, grid2chem)
 nchem = RM_GetChemistryCellCount(id)
 
MPI:
Called by root and (or) workers.
integer function RM_GetComponent ( integer, intent(in)  id,
integer, intent(in)  num,
character(len=*), intent(inout)  comp_name 
)

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. Fortran, 1 based.
comp_nameThe string value associated with component num.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetComponentcount
Fortran Example:
 ! Get list of components
 ncomps = RM_FindComponents(id)
 allocate(components(ncomps))
 do i = 1, ncomps
   status = RM_GetComponent(id, i, components(i))
 enddo
 
MPI:
Called by root and (or) workers.
integer function RM_GetComponentcount ( integer, intent(in)  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.
Fortran Example:
 ncomps1 = RM_GetComponentcount(id)
 
MPI:
Called by root.
integer function RM_GetConcentrations ( integer, intent(in)  id,
double precision, dimension(:,:), intent(out)  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_SetSaturation), porosity (RM_SetPorosity), and representative volume (RM_SetRepresentativeVolume), and the mass of solution is volume times density as defined by RM_SetDensity. 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_GetSaturation, RM_SetConcentrations, RM_SetDensity, RM_SetRepresentativeVolume, RM_SetSaturation, RM_SetUnitsSolution, RM_UseSolutionDensityVolume.
Fortran Example:
 allocate(c(nxyz, ncomps))
 status = RM_RunCells(id)
 status = RM_GetConcentrations(id, c)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetDensity ( integer, intent(in)  id,
double precision, dimension(:), intent(out)  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. 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).
Fortran Example:
 allocate(density(nxyz))
 status = RM_RunCells(id)
 status = RM_GetDensity(id, density)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetEndCell ( integer, intent(in)  id,
integer, dimension(:), intent(out)  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.
Fortran Example:
 n = RM_GetThreadCount(id) * RM_GetMpiTasks(id)
 allocate(ec(n))
 status = RM_GetEndCell(id, ec)
 
MPI:
Called by root and (or) workers.
integer function RM_GetErrorString ( integer, intent(in)  id,
character(len=*), intent(out)  errstr 
)

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.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
Fortran Example:
 character(len=:), allocatable :: errstr
 integer l
 if (status .ne. 0) then
   l = RM_GetErrorStringLength(id)
   allocate (character(len=l) :: errstr)
   status = RM_GetErrorString(id, errstr)
   write(*,"(A)") errstr
   deallocate(errstr)
   status = RM_Destroy(id)
   stop
 endif 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetErrorStringlength ( integer, intent(in)  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.
Fortran Example:
 character(len=:), allocatable :: errstr
 integer l
 if (status .ne. 0) then
   l = RM_GetErrorStringLength(id)
   allocate (character(len=l) :: errstr)
   status = RM_GetErrorString(id, errstr)
   write(*,"(A)") errstr
   deallocate(errstr)
   status = RM_Destroy(id)
   stop
 endif 
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetFilePrefix ( integer, intent(in)  id,
character(len=*), intent(inout)  prefix 
)

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.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_SetFilePrefix.
Fortran Example:
 character(100) :: string
 character(200) :: string1
 status = RM_GetFilePrefix(id, string)
 string1 = "File prefix: "//string
 status = RM_OutputMessage(id, string1)
 
MPI:
Called by root and (or) workers.
integer function RM_GetGfw ( integer, intent(in)  id,
double precision, dimension(:), intent(out)  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_GetComponentcount, RM_GetComponent.
Fortran Example:
 character(100),   dimension(:), allocatable   :: components
 double precision, dimension(:), allocatable   :: gfw
 ncomps = RM_FindComponents(id)
 allocate(components(ncomps))
 allocate(gfw(ncomps))
 status = RM_GetGfw(id, gfw)
 do i = 1, ncomps
   status = RM_GetComponent(id, i, components(i))
   write(string,"(A10, F15.4)") components(i), gfw(i)
   status = RM_OutputMessage(id, string)
 enddo
 
MPI:
Called by root.
integer function RM_GetGridCellCount ( integer, intent(in)  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.
Fortran Example:
 nxyz = RM_GetGridCellCount(id)
 write(string1, "(A,I)") "Number of grid cells in the user's model: ", nxyz
 status = RM_OutputMessage(id, trim(string1))
 
MPI:
Called by root and (or) workers.
integer function RM_GetIPhreeqcId ( integer, intent(in)  id,
integer, intent(in)  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.
Fortran Example:
 ! Utility pointer is worker number nthreads + 1
 iphreeqc_id1 = RM_GetIPhreeqcId(id, RM_GetThreadCount(id) + 1)
 
MPI:
Called by root and (or) workers.
integer function RM_GetMpiMyself ( integer, intent(in)  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.
Fortran Example:
 write(string1, "(A,I)") "MPI task number: ", RM_GetMpiMyself(id)
 status = RM_OutputMessage(id, string1)
 
MPI:
Called by root and (or) workers.
integer function RM_GetMpiTasks ( integer, intent(in)  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_GetMpiMyself, RM_Create.
Fortran Example:
 mpi_tasks = RM_GetMpiTasks(id)
 write(string1, "(A,I)") "Number of MPI processes: ", mpi_tasks
 status = RM_OutputMessage(id, string1)
 
MPI:
Called by root and (or) workers.
integer function RM_GetNthSelectedOutputUserNumber ( integer, intent(in)  id,
integer, intent(in)  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. Fortran, 1 based.
Return values
Theuser number of the nth selected-output definition, negative is failure (See RM_DecodeError).
See also
RM_GetSelectedOutput, RM_GetSelectedOutputcolumncount, RM_GetSelectedOutputcount, RM_GetSelectedOutputheading, RM_GetSelectedOutputrowcount, RM_SetCurrentSelectedOutputUserNumber, RM_SetSelectedOutputOn.
Fortran Example:
 do isel = 1, RM_GetSelectedOutputcount(id)
   n_user = RM_GetNthSelectedOutputUserNumber(id, isel)
   status = RM_SetCurrentSelectedOutputUserNumber(id, n_user)
   write(*,*) "Selected output sequence number: ", isel)
   write(*,*) "Selected output user number:     ", n_user)
   col = RM_GetSelectedOutputcolumncount(id)
   allocate(selected_out(nxyz,col))
   status = RM_GetSelectedOutput(id, selected_out)
   ! Process results here
   deallocate(selected_out)
 enddo
 
MPI:
Called by root.
integer function RM_GetSaturation ( integer, intent(in)  id,
double precision, dimension(:), intent(out)  sat_calc 
)

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

Parameters
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_SetSaturation.
Fortran Example:
 allocate(sat_calc(nxyz))
 status = RM_RunCells(id)
 status = RM_GetSaturation(id, sat_calc)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetSelectedOutput ( integer, intent(in)  id,
double precision, dimension(:,:), intent(out)  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, 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_GetNthSelectedOutputUserNumber, RM_GetSelectedOutputcolumncount, RM_GetSelectedOutputcount, RM_GetSelectedOutputheading, RM_GetSelectedOutputrowcount, RM_SetCurrentSelectedOutputUserNumber, RM_SetSelectedOutputOn.
Fortran Example:
 do isel = 1, RM_GetSelectedOutputcount(id)
   n_user = RM_GetNthSelectedOutputUserNumber(id, isel)
   status = RM_SetCurrentSelectedOutputUserNumber(id, n_user)
   col = RM_GetSelectedOutputcolumncount(id)
   allocate(selected_out(nxyz,col))
   status = RM_GetSelectedOutput(id, selected_out)
   ! Process results here
   deallocate(selected_out)
 enddo
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetSelectedOutputcolumncount ( integer, intent(in)  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_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputcount, RM_GetSelectedOutputheading, RM_GetSelectedOutputrowcount, RM_SetCurrentSelectedOutputUserNumber, RM_SetSelectedOutputOn.
Fortran Example:
 do isel = 1, RM_GetSelectedOutputcount(id)
   n_user = RM_GetNthSelectedOutputUserNumber(id, isel)
   status = RM_SetCurrentSelectedOutputUserNumber(id, n_user)
   col = RM_GetSelectedOutputcolumncount(id)
   allocate(selected_out(nxyz,col))
   status = RM_GetSelectedOutput(id, selected_out)
   ! Process results here
   deallocate(selected_out)
 enddo
 
MPI:
Called by root.
integer function RM_GetSelectedOutputcount ( integer, intent(in)  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_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputcolumncount, RM_GetSelectedOutputheading, RM_GetSelectedOutputrowcount, RM_SetCurrentSelectedOutputUserNumber, RM_SetSelectedOutputOn.
Fortran Example:
 do isel = 1, RM_GetSelectedOutputcount(id)
   n_user = RM_GetNthSelectedOutputUserNumber(id, isel)
   status = RM_SetCurrentSelectedOutputUserNumber(id, n_user)
   col = RM_GetSelectedOutputcolumncount(id)
   allocate(selected_out(nxyz,col))
   status = RM_GetSelectedOutput(id, selected_out)
   ! Process results here
   deallocate(selected_out)
 enddo
 
MPI:
Called by root.
integer function RM_GetSelectedOutputheading ( integer, intent(in)  id,
integer, intent(in)  icol,
character(len=*), intent(out)  heading 
)

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. Fortran, 1 based.
headingA string buffer to receive the heading.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputcolumncount, RM_GetSelectedOutputcount, RM_GetSelectedOutputrowcount, RM_SetCurrentSelectedOutputUserNumber, RM_SetSelectedOutputOn.
Fortran Example:
 do isel = 1, RM_GetSelectedOutputcount(id)
   n_user = RM_GetNthSelectedOutputUserNumber(id, isel)
   status = RM_SetCurrentSelectedOutputUserNumber(id, n_user)
   col = RM_GetSelectedOutputcolumncount(id)
   do j = 1, col
     status = RM_GetSelectedOutputheading(id, j, heading)
     write(*,'(10x,i2,A2,A10,A2,f10.4)') j, " ", trim(heading)
   enddo
 enddo
 
MPI:
Called by root.
integer function RM_GetSelectedOutputrowcount ( integer, intent(in)  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_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputcolumncount, RM_GetSelectedOutputcount, RM_GetSelectedOutputheading, RM_SetCurrentSelectedOutputUserNumber, RM_SetSelectedOutputOn.
Fortran Example:
 do isel = 1, RM_GetSelectedOutputcount(id)
   n_user = RM_GetNthSelectedOutputUserNumber(id, isel)
   status = RM_SetCurrentSelectedOutputUserNumber(id, n_user)
   col = RM_GetSelectedOutputcolumncount(id)
   allocate(selected_out(nxyz,col))
   status = RM_GetSelectedOutput(id, selected_out)
   ! Print results
   do i = 1, RM_GetSelectedOutputrowcount(id)
     write(*,*) "Cell number ", i
     write(*,*) "     Selected output: "
     do j = 1, col
       status = RM_GetSelectedOutputheading(id, j, heading)
       write(*,'(10x,i2,A2,A10,A2,f10.4)') j, " ", trim(heading),": ", selected_out(i,j)
     enddo
   enddo
   deallocate(selected_out)
 enddo
 
MPI:
Called by root.
integer function RM_GetSolutionVolume ( integer, intent(in)  id,
double precision, dimension(:), intent(out)  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_GetSaturation.
Fortran Example:
 allocate(volume(nxyz))
 status = RM_RunCells(id)
 status = RM_GetSolutionVolume(id, volume)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetSpeciesConcentrations ( integer, intent(in)  id,
double precision, dimension(:,:), intent(out)  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_GetSpeciesZ, RM_GetSpeciesName, RM_SpeciesConcentrations2Module, RM_GetSpeciesSaveOn, RM_SetSpeciesSaveOn.
Fortran Example:
 status = RM_SetSpeciesSaveOn(id, 1)
 ncomps = RM_FindComponents(id)
 nspecies = RM_GetSpeciesCount(id)
 nxyz = RM_GetGridCellCount(id)
 allocate(species_c(nxyz, nspecies))
 status = RM_RunCells(id)
 status = RM_GetSpeciesConcentrations(id, species_c)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_GetSpeciesCount ( integer, intent(in)  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_GetSpeciesZ, RM_GetSpeciesName, RM_SpeciesConcentrations2Module, RM_GetSpeciesSaveOn, RM_SetSpeciesSaveOn.
Fortran Example:
 status = RM_SetSpeciesSaveOn(id, 1)
 ncomps = RM_FindComponents(id)
 nspecies = RM_GetSpeciesCount(id)
 nxyz = RM_GetGridCellCount(id)
 allocate(species_c(nxyz, nspecies))
 status = RM_RunCells(id)
 status = RM_GetSpeciesConcentrations(id, species_c)
 
MPI:
Called by root and (or) workers.
integer function RM_GetSpeciesD25 ( integer, intent(in)  id,
double precision, dimension(:), intent(out)  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_GetSpeciesZ, RM_GetSpeciesName, RM_SpeciesConcentrations2Module, RM_GetSpeciesSaveOn, RM_SetSpeciesSaveOn.
Fortran Example:
 status = RM_SetSpeciesSaveOn(id, 1)
 ncomps = RM_FindComponents(id)
 nspecies = RM_GetSpeciesCount(id)
 allocate(diffc(nspecies))
 status = RM_GetSpeciesD25(id, diffc)
 
MPI:
Called by root and (or) workers.
integer function RM_GetSpeciesName ( integer, intent(in)  id,
integer, intent(in)  i,
character(len=*), intent(out)  name 
)

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. Fortran, 1 based.
nameCharacter array to receive the species name.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSpeciesConcentrations, RM_GetSpeciesCount, RM_GetSpeciesD25, RM_GetSpeciesZ, RM_SpeciesConcentrations2Module, RM_GetSpeciesSaveOn, RM_SetSpeciesSaveOn.
Fortran Example:
 char*100 name
 ...
 status = RM_SetSpeciesSaveOn(id, 1)
 ncomps = RM_FindComponents(id)
 nspecies = RM_GetSpeciesCount(id)
 do i = 1, nspecies
   status = RM_GetSpeciesName(id, i, name)
   write(*,*) name
 enddo
 
MPI:
Called by root and (or) workers.
integer function RM_GetSpeciesSaveOn ( integer, intent(in)  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_GetSpeciesZ, RM_GetSpeciesName, RM_SpeciesConcentrations2Module, RM_SetSpeciesSaveOn.
Fortran Example:
 save_on = RM_GetSpeciesSaveOn(id)
 if (save_on .ne. 0) then
   write(*,*) "Reaction module is saving species concentrations"
 else
   write(*,*) "Reaction module is not saving species concentrations"
 end
 
MPI:
Called by root and (or) workers.
integer function RM_GetSpeciesZ ( integer, intent(in)  id,
double precision, dimension(:), intent(out)  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_GetSpeciesName, RM_SpeciesConcentrations2Module, RM_GetSpeciesSaveOn, RM_SetSpeciesSaveOn.
Fortran Example:
 status = RM_SetSpeciesSaveOn(id, 1)
 ncomps = RM_FindComponents(id)
 nspecies = RM_GetSpeciesCount(id)
 allocate(z(nspecies))
 status = RM_GetSpeciesZ(id, z)
 
MPI:
Called by root and (or) workers.
integer function RM_GetStartCell ( integer, intent(in)  id,
integer, dimension(:), intent(out)  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.
Fortran Example:
 n = RM_GetThreadCount(id) * RM_GetMpiTasks(id)
 allocate(sc(n))
 status = RM_GetStartCell(id, sc)
 
MPI:
Called by root and (or) workers.
integer function RM_GetThreadCount ( integer, intent(in)  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.
Fortran Example:
 write(string1, "(A,I)") "Number of threads: ", RM_GetThreadCount(id)
 status = RM_OutputMessage(id, string1)
 
MPI:
Called by root and (or) workers; result is always 1.
double precision function RM_GetTime ( integer, intent(in)  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.
Fortran Example:
 write(string, "(A32,F15.1,A)") "Beginning transport calculation ", &
       RM_GetTime(id) * RM_GetTimeconversion(id), " days"
 status = RM_LogMessage(id, string)
 
MPI:
Called by root and (or) workers.
double precision function RM_GetTimeconversion ( integer, intent(in)  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.
Fortran Example:
 write(string, "(A32,F15.1,A)") "Beginning transport calculation ", &
       RM_GetTime(id) * RM_GetTimeconversion(id), " days"
 status = RM_LogMessage(id, string)
 
MPI:
Called by root and (or) workers.
double precision function RM_GetTimestep ( integer, intent(in)  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.
Fortran Example:
 write(string, "(A32,F15.1,A)") "          Time step             ", &
       RM_GetTimestep(id) * RM_GetTimeconversion(id), " days"
 status = RM_LogMessage(id, string)
 
MPI:
Called by root and (or) workers.
integer function RM_InitialPhreeqc2Concentrations ( integer, intent(in)  id,
double precision, dimension(:,:), intent(out)  bc_conc,
integer, intent(in)  n_boundary,
integer, dimension(:), intent(in)  bc1,
integer, dimension(:), intent(in), optional  bc2,
double precision, dimension(:), intent(in), optional  f1 
)

Fills an array (bc_conc) 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 bc1, then the highest numbered solution in the InitialPhreeqc instance will be used for that cell. Concentrations may be a mixture of two solutions, bc1 and bc2, with a mixing fraction for bc1 1 of f1 and mixing fraction for bc2 of (1 - f1). A negative value for bc2 implies no mixing, and the associated value for f1 is ignored. If bc2 and f1 are omitted, no mixing is used; concentrations are derived from bc1 only.

Parameters
idThe instance id returned from RM_Create.
bc_concArray of concentrations extracted from the InitialPhreeqc instance. The dimension of bc_conc is (n_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.
bc1Array of solution index numbers that refer to solutions in the InitialPhreeqc instance. Size is n_boundary.
bc2Array of solution index numbers that that refer to solutions in the InitialPhreeqc instance and are defined to mix with bc1. Size is n_boundary. Optional in Fortran.
f1Fraction of bc1 that mixes with (1-f1) of bc2. Size is (n_boundary). Optional in Fortran.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetComponentcount.
Fortran Example:
 nbound = 1
 allocate(bc1(nbound), bc2(nbound), f1(nbound))
 allocate(bc_conc(nbound, ncomps))
 bc1 = 0           ! solution 0 from InitialPhreeqc instance
 bc2 = -1          ! no bc2 solution for mixing
 f1 = 1.0          ! mixing fraction for bc1
 status = RM_InitialPhreeqc2Concentrations(id, bc_conc, nbound, bc1, bc2, f1)
 
MPI:
Called by root.
integer function RM_InitialPhreeqc2Module ( integer, intent(in)  id,
integer, dimension(:,:), intent(in)  ic1,
integer, dimension(:,:), intent(in), optional  ic2,
double precision, dimension(:,:), intent(in), optional  f1 
)

Transfer solutions and reactants from the InitialPhreeqc instance to the reaction-module workers, possibly with mixing. In its simplest form, ic1 is used to select initial conditions, including solutions and reactants, for each cell of the model, without mixing. ic1 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 Fortran, ic1(100, 4) = 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, ic2 contains numbers for a second entity that mixes with the entity defined in ic1. F1 contains the mixing fraction for ic1, whereas (1 - f1) is the mixing fraction for ic2. In Fortran, ic1(100, 4) = 2, initial_conditions2(100, 4) = 3, f1(100, 4) = 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 ic2 is negative, no mixing occurs. If ic2 and f1 are omitted, no mixing is used, and initial conditions are derived solely from ic1.

Parameters
idThe instance id returned from RM_Create.
ic1Array 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.
ic2Array of solution and reactant index numbers that refer to definitions in the InitialPhreeqc instance. Nonnegative values of ic2 result in mixing with the entities defined in ic1. Negative values result in no mixing. Size is (nxyz,7). The order of definitions is given above. Optional in Fortran; omitting results in no mixing.
f1Fraction of ic1 that mixes with (1-f1) of ic2. Size is (nxyz,7). The order of definitions is given above. Optional in Fortran; omitting results in no mixing.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialPhreeqcCell2Module.
Fortran Example:
 allocate(ic1(nxyz,7), ic2(nxyz,7), f1(nxyz,7))
 ic1 = -1
 ic2 = -1
 f1 = 1.0
 do i = 1, nxyz
   ic1(i,1) = 1       ! Solution 1
   ic1(i,2) = -1      ! Equilibrium phases none
   ic1(i,3) = 1       ! Exchange 1
   ic1(i,4) = -1      ! Surface none
   ic1(i,5) = -1      ! Gas phase none
   ic1(i,6) = -1      ! Solid solutions none
   ic1(i,7) = -1      ! Kinetics none
 enddo
 status = RM_InitialPhreeqc2Module(id, ic1, ic2, f1)1))
 ! No mixing is defined, so the following is equivalent
 status = RM_InitialPhreeqc2Module(id, ic1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_InitialPhreeqc2SpeciesConcentrations ( integer, intent(in)  id,
double precision, dimension(:,:), intent(out)  bc_conc,
integer, intent(in)  n_boundary,
integer, dimension(:), intent(in)  bc1,
integer, dimension(:), intent(in), optional  bc2,
double precision, dimension(:), intent(in), optional  f1 
)

Fills an array (bc_conc) 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 bc1, then the highest numbered solution in the InitialPhreeqc instance will be used for that cell. Concentrations may be a mixture of two solutions, bc1 and bc2, with a mixing fraction for bc1 of f1 and mixing fraction for bc2 of (1 - f1). A negative value for bc2 implies no mixing, and the associated value for f1 is ignored. If bc2 and f1 are omitted, no mixing is used; concentrations are derived from bc1 only.

Parameters
idThe instance id returned from RM_Create.
bc_concArray 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.
bc1Array of solution index numbers that refer to solutions in the InitialPhreeqc instance. Size is n_boundary.
bc2Array of solution index numbers that that refer to solutions in the InitialPhreeqc instance and are defined to mix with bc1. Size is n_boundary. Optional in Fortran.
f1Fraction of bc1 that mixes with (1-f1) of bc2. Size is n_boundary. Optional in Fortran.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_FindComponents, RM_GetSpeciesCount, RM_SetSpeciesSaveOn.
Fortran Example:
 nbound = 1
 allocate(bc1(nbound), bc2(nbound), f1(nbound))
 allocate(bc_conc(nbound, ncomps))
 bc1 = 0           ! solution 0 from InitialPhreeqc instance
 bc2 = -1          ! no bc2 solution for mixing
 f1 = 1.0          ! mixing fraction for bc1
 status = RM_InitialPhreeqc2SpeciesConcentrations(id, bc_conc, nbound, bc1, bc2, f1)
 
MPI:
Called by root.
integer function RM_InitialPhreeqcCell2Module ( integer, intent(in)  id,
integer, intent(in)  n_user,
integer, dimension(:), intent(in)  cell_numbers,
integer, intent(in)  n_cell 
)

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

Parameters
idThe instance id returned from RM_Create.
n_userCell 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.
cell_numbersA list of cell numbers in the user's grid-cell numbering system that will be populated with cell n_user from the InitialPhreeqc instance.
n_cellThe number of cell numbers in the cell_numbers list.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_InitialPhreeqc2Module.
Fortran Example:
 allocate (cell_numbers(2))
 cell_numbers(1) = 18
 cell_numbers(2) = 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, cell_numbers, 2)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_LoadDatabase ( integer, intent(in)  id,
character(len=*), intent(in)  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.
Fortran Example:
 status = RM_LoadDatabase(id, "phreeqc.dat")
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_LogMessage ( integer, intent(in)  id,
character(len=*), intent(in)  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_OpenFiles, RM_ErrorMessage, RM_OutputMessage, RM_ScreenMessage, RM_WarningMessage.
Fortran Example:
 write(string, "(A32,F15.1,A)") "Beginning transport calculation ", &
       RM_GetTime(id) * RM_GetTimeconversion(id), " days"
 status = RM_LogMessage(id, string)
 
MPI:
Called by root.
integer function RM_MpiWorker ( integer, intent(in)  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.
Fortran Example:
 status = RM_MpiWorker(id)
 
MPI:
Called by all workers.
integer function RM_MpiWorkerbreak ( integer, intent(in)  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.
Fortran Example:
 status = RM_MpiWorkerbreak(id)
 
MPI:
Called by root.
integer function RM_OpenFiles ( integer, intent(in)  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_SetFilePrefix, RM_GetFilePrefix, RM_CloseFiles, RM_ErrorMessage, RM_LogMessage, RM_OutputMessage, RM_WarningMessage.
Fortran Example:
 status = RM_SetFilePrefix(id, "Advect_f90")
 status = RM_OpenFiles(id)
 
MPI:
Called by root.
integer function RM_OutputMessage ( integer, intent(in)  id,
character(len=*), intent(in)  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.
Fortran Example:
 write(string1, "(A,I10)") "Number of threads:                                ", RM_GetThreadCount(id)
 status = RM_OutputMessage(id, string1)
 write(string1, "(A,I10)") "Number of MPI processes:                          ", RM_GetMpiTasks(id)
 status = RM_OutputMessage(id, string1)
 write(string1, "(A,I10)") "MPI task number:                                  ", RM_GetMpiMyself(id)
 status = RM_OutputMessage(id, string1)
 status = RM_GetFilePrefix(id, string)
 write(string1, "(A,A)") "File prefix:                                      ", string
 status = RM_OutputMessage(id, trim(string1))
 write(string1, "(A,I10)") "Number of grid cells in the user's model:         ", RM_GetGridCellCount(id)
 status = RM_OutputMessage(id, trim(string1))
 write(string1, "(A,I10)") "Number of chemistry cells in the reaction module: ", RM_GetChemistryCellCount(id)
 status = RM_OutputMessage(id, trim(string1))
 write(string1, "(A,I10)") "Number of components for transport:               ", RM_GetComponentcount(id)
 status = RM_OutputMessage(id, trim(string1))
 
MPI:
Called by root.
integer function RM_RunCells ( integer, intent(in)  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_SetSaturation), 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_SetSaturation, RM_SetTemperature, RM_SetTimestep.
Fortran Example:
 status = RM_SetPorosity(id, por)                ! If pore volume changes
 status = RM_SetSaturation(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_GetDensity(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.
integer function RM_RunFile ( integer, intent(in)  id,
integer, intent(in)  workers,
integer, intent(in)  initial_phreeqc,
integer, intent(in)  utility,
character(len=*), intent(in)  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.
Fortran Example:
 status = RM_RunFile(id, 1, 1, 1, "advect.pqi")
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_RunString ( integer, intent(in)  id,
integer, intent(in)  initial_phreeqc,
integer, intent(in)  workers,
integer, intent(in)  utility,
character(len=*), intent(in)  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.
Fortran Example:
 string = "DELETE; -all"
 status = RM_RunString(id, 1, 0, 1, string)  ! workers, initial_phreeqc, utility
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_ScreenMessage ( integer, intent(in)  id,
character(len=*), intent(in)  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.
Fortran Example:
 write(string, "(A32,F15.1,A)") "Beginning reaction calculation  ", &
       time * RM_GetTimeconversion(id), " days"
 status = RM_ScreenMessage(id, string)
 
MPI:
Called by root and (or) workers.
integer function RM_SetComponentH2O ( integer, intent(in)  id,
integer, intent(in)  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.
Fortran Example:
 status = RM_SetComponentH2O(id, 0)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetConcentrations ( integer, intent(in)  id,
double precision, dimension(:,:), intent(in)  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_SetSaturation), 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_SetDensity) 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_SetDensity, RM_SetPorosity, RM_SetRepresentativeVolume, RM_SetSaturation, RM_SetUnitsSolution.
Fortran Example:
 allocate(c(nxyz, ncomps))
 ...
 call advect_f90(c, bc_conc, ncomps, nxyz)
 status = RM_SetPorosity(id, por)               ! If porosity changes
 status = RM_SetSaturation(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_GetDensity(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.
integer function RM_SetCurrentSelectedOutputUserNumber ( integer, intent(in)  id,
integer, intent(in)  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_GetNthSelectedOutputUserNumber, RM_GetSelectedOutput, RM_GetSelectedOutputcolumncount, RM_GetSelectedOutputcount, RM_GetSelectedOutputrowcount, RM_GetSelectedOutputheading, RM_SetSelectedOutputOn.
Fortran Example:
 do isel = 1, RM_GetSelectedOutputcount(id)
   n_user = RM_GetNthSelectedOutputUserNumber(id, isel)
   status = RM_SetCurrentSelectedOutputUserNumber(id, n_user)
   col = RM_GetSelectedOutputcolumncount(id)
   allocate(selected_out(nxyz,col))
   status = RM_GetSelectedOutput(id, selected_out)
   ! Process results here
   deallocate(selected_out)
 enddo
 
MPI:
Called by root. */
integer function RM_SetDensity ( integer, intent(in)  id,
double precision, dimension(:), intent(in)  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.
Fortran Example:
 allocate(density(nxyz))
 density = 1.0
 status = RM_SetDensity(id, density(1))
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetDumpFileName ( integer, intent(in)  id,
character(len=*), intent(in)  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.
Fortran Example:
 status = RM_SetDumpFileName(id, "advection_f90.dmp")
 dump_on = 1
 append = 0
 status = RM_DumpModule(id, dump_on, append)
 
MPI:
Called by root.
integer function RM_SetErrorHandlerMode ( integer, intent(in)  id,
integer, intent(in)  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).
Fortran Example:
 id = RM_create(nxyz, nthreads)
 status = RM_SetErrorHandlerMode(id, 2)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetFilePrefix ( integer, intent(in)  id,
character(len=*), intent(in)  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_OpenFiles, RM_CloseFiles.
Fortran Example:
 status = RM_SetFilePrefix(id, "Advect_f90")
 status = RM_OpenFiles(id)
 
MPI:
Called by root.
integer function RM_SetMpiWorkerCallback ( integer, intent(in)  id,
  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.

For Fortran, the functions that are called from the callback function can use USE statements to find the data necessary to perform the tasks, and the only argument to the callback function is an integer message argument. RM_SetMpiWorkerCallback 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.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_MpiWorker, RM_MpiWorkerbreak.
Fortran Example:
 Code executed by root:
 status = do_something()
 
 Code executed by workers:
 status = RM_SetMpiWorkerCallback(id, worker_tasks_f)
 status = RM_MpiWorker(id)
 
 Code executed by root and workers:    
 integer function do_something
   implicit none
   INCLUDE 'mpif.h'
   integer status
   integer i, method_number, mpi_myself, mpi_task, mpi_tasks, worker_number;
   method_number = 1000
   call MPI_Comm_size(MPI_COMM_WORLD, mpi_tasks, status)
   call MPI_Comm_rank(MPI_COMM_WORLD, mpi_myself, status)
   if (mpi_myself .eq. 0) then     
     CALL MPI_Bcast(method_number, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, status)
     write(*,*) "I am root."
     do i = 1, mpi_tasks-1
       CALL MPI_Recv(worker_number, 1, MPI_INTEGER, i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE, status)
       write(*,*) "Recieved data from worker number ", worker_number, "."
     enddo
   else
     CALL MPI_Send(mpi_myself, 1, MPI_INTEGER, 0, 0, MPI_COMM_WORLD, status)
   endif
   do_something = 0
 end function do_something
 
 Code called by workers from method MpiWorker:
 integer(kind=C_INT) function worker_tasks_f(method_number) BIND(C, NAME='worker_tasks_f')
   USE ISO_C_BINDING
   implicit none
   interface
     integer function do_something
     end function do_something
   end interface
   integer(kind=c_int), intent(in) :: method_number
   integer :: status
   if (method_number .eq. 1000) then
     status = do_something()
   endif
   worker_tasks_f = 0
 end function worker_tasks_f
 
MPI:
Called by workers, before call to RM_MpiWorker.
Parameters
[in]id
integer function RM_SetPartitionUZSolids ( integer, intent(in)  id,
integer, intent(in)  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).
Fortran Example:
 status = RM_SetPartitionUZSolids(id, 0)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetPorosity ( integer, intent(in)  id,
double precision, dimension(:), intent(in)  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_SetSaturation), 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_GetSaturation, RM_SetRepresentativeVolume, RM_SetSaturation.
Fortran Example:
 allocate(por(nxyz))
 por = 0.2
 status = RM_SetPorosity(id, por(1))
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetPressure ( integer, intent(in)  id,
double precision, dimension(:), intent(in)  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.
Fortran Example:
 allocate(pressure(nxyz))
 pressure = 2.0
 status = RM_SetPressure(id, pressure(1))
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetPrintChemistryMask ( integer, intent(in)  id,
integer, dimension(:), intent(in)  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.
Fortran Example:
 allocate(print_chemistry_mask(nxyz))
   do i = 1, nxyz/2
   print_chemistry_mask(i) = 1
   print_chemistry_mask(i+nxyz/2) = 0
 enddo
 status = RM_SetPrintChemistryMask(id, print_chemistry_mask(1))
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetPrintChemistryOn ( integer, intent(in)  id,
integer, intent(in)  workers,
integer, intent(in)  initial_phreeqc,
integer, intent(in)  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.
Fortran 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.
integer function RM_SetRebalanceByCell ( integer, intent(in)  id,
integer, intent(in)  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.
Fortran Example:
 status = RM_SetRebalanceByCell(id, 1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetRebalanceFraction ( integer, intent(in)  id,
double precision, intent(in)  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.
Fortran Example:
 status = RM_SetRebalanceFraction(id, 0.5d0)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetRepresentativeVolume ( integer, intent(in)  id,
double precision, dimension(:), intent(in)  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_SetSaturation). The numerical method of PHREEQC is more robust if the water volume for a reaction cell is within a couple orders of magnitude of 1.0. Small water volumes caused by small porosities and (or) small saturations (and (or) small representative volumes) may cause non-convergence of the numerical method. In these cases, a larger representative volume may help. Note that increasing the representative volume also increases the number of moles of the reactants in the reaction cell (minerals, surfaces, exchangers, and others), which are defined as moles per representative volume.

Parameters
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_SetSaturation.
Fortran Example:
 double precision, dimension(:), allocatable   :: rv
 allocate(rv(nxyz))
 rv = 1.0
 status = RM_SetRepresentativeVolume(id, rv(1))
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetSaturation ( integer, intent(in)  id,
double precision, dimension(:), intent(in)  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_SetSaturation), 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_GetDensity, RM_GetSolutionVolume, and RM_GetSaturation 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_GetDensity, RM_GetSaturation, RM_GetSolutionVolume, RM_SetPorosity, RM_SetRepresentativeVolume.
Fortran Example:
 allocate(sat(nxyz))
 sat = 1.0
 status = RM_SetSaturation(id, sat(1))
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetScreenOn ( integer, intent(in)  id,
integer, intent(in)  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.
Fortran Example:
 status = RM_SetScreenOn(rm_id, 1)
 
MPI:
Called by root.
integer function RM_SetSelectedOutputOn ( integer, intent(in)  id,
integer, intent(in)  tf 
)

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.
tf0, disable selected output; 1, enable selected output.
Return values
IRM_RESULT0 is success, negative is failure (See RM_DecodeError).
See also
RM_GetSelectedOutput, RM_SetPrintChemistryOn.
Fortran Example:
 status = RM_SetSelectedOutputOn(id, 1)        ! enable selected output
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetSpeciesSaveOn ( integer, intent(in)  id,
integer, intent(in)  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_GetSpeciesSaveOn, RM_GetSpeciesZ, RM_GetSpeciesName, RM_SpeciesConcentrations2Module.
Fortran Example:
 save_on = RM_SetSpeciesSaveOn(id, 1)
 
MPI:
Called by root and (or) workers.
integer function RM_SetTemperature ( integer, intent(in)  id,
double precision, dimension(:), intent(in)  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.
Fortran Example:
 allocate(temperature(nxyz))
 temperature = 20.0
 status = RM_SetTemperature(id, temperature(1))
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetTime ( integer, intent(in)  id,
double precision, intent(in)  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_SetTimestep, RM_SetTimeconversion.
Fortran Example:
 status = RM_SetTime(id, time)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetTimeconversion ( integer, intent(in)  id,
double precision, intent(in)  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.
Fortran Example:
 status = RM_SetTimeconversion(id, dble(1.0 / 86400.0)) ! days
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetTimestep ( integer, intent(in)  id,
double precision, intent(in)  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.
Fortran Example:
 status = RM_SetTimestep(id, time_step)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetUnitsExchange ( integer, intent(in)  id,
integer, intent(in)  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.
Fortran Example:
 status = RM_SetUnitsExchange(id, 1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetUnitsGasPhase ( integer, intent(in)  id,
integer, intent(in)  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.
Fortran Example:
 status = RM_SetUnitsGasPhase(id, 1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetUnitsKinetics ( integer, intent(in)  id,
integer, intent(in)  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_SetSaturation), 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_SetSaturation.
Fortran Example:
 status = RM_SetUnitsKinetics(id, 1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetUnitsPPassemblage ( integer, intent(in)  id,
integer, intent(in)  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.
Fortran Example:
 status = RM_SetUnitsPPassemblage(id, 1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetUnitsSolution ( integer, intent(in)  id,
integer, intent(in)  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_SetSaturation), 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_SetDensity) and multiplied by the solution volume.

To convert from moles of element in the representative volume of a reaction cell to mg/L, the number of moles of an element is divided by the solution volume resulting in mol/L, and then converted to mg/L. To convert from moles of element in a cell to mol/L, the number of moles of an element is divided by the solution volume resulting in mol/L. To convert from moles of element in a cell to mass fraction, the number of moles of an element is converted to kg and divided by the total mass of the solution. Two options are available for the volume and mass of solution that are used in converting to transport concentrations: (1) the volume and mass of solution are calculated by PHREEQC, or (2) the volume of solution is the product of porosity (RM_SetPorosity), saturation (RM_SetSaturation), and representative volume (RM_SetRepresentativeVolume), and the mass of solution is volume times density as defined by RM_SetDensity. 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_SetDensity, RM_SetPorosity, RM_SetRepresentativeVolume, RM_SetSaturation, RM_UseSolutionDensityVolume.
Fortran Example:
 status = RM_SetUnitsSolution(id, 1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetUnitsSSassemblage ( integer, intent(in)  id,
integer, intent(in)  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.

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.

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.

Fortran Example:
 status = RM_SetUnitsSSassemblage(id, 1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SetUnitsSurface ( integer, intent(in)  id,
integer, intent(in)  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.
Fortran Example:
 status = RM_SetUnitsSurface(id, 1)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_SpeciesConcentrations2Module ( integer, intent(in)  id,
double precision, dimension(:,:), intent(in)  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_GetSpeciesZ, RM_GetSpeciesName, RM_GetSpeciesSaveOn, RM_SetSpeciesSaveOn.
Fortran Example:
 status = RM_SetSpeciesSaveOn(id, 1)
 ncomps = RM_FindComponents(id)
 nspecies = RM_GetSpeciesCount(id)
 nxyz = RM_GetGridCellCount(id)
 allocate(species_c(nxyz, nspecies))
 ...
 status = RM_SpeciesConcentrations2Module(id, species_c(1,1))
 status = RM_RunCells(id)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_UseSolutionDensityVolume ( integer, intent(in)  id,
integer, intent(in)  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_SetDensity) and solution volume are defined by the product of saturation (RM_SetSaturation), 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_SetDensity and the volume determined by the product of RM_SetSaturation, RM_SetPorosity, and RM_SetRepresentativeVolume, will be used to calculate concentrations retrieved by RM_GetConcentrations.
See also
RM_GetConcentrations, RM_SetDensity, RM_SetPorosity, RM_SetRepresentativeVolume, RM_SetSaturation.
Fortran Example:
 status = RM_UseSolutionDensityVolume(id, 0)
 
MPI:
Called by root, workers must be in the loop of RM_MpiWorker.
integer function RM_WarningMessage ( integer, intent(in)  id,
character(len=*), intent(in)  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_OpenFiles, RM_LogMessage, RM_OutputMessage, RM_ScreenMessage, RM_ErrorMessage.
Fortran 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.

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