API Reference
- madupite.MPI_master()
This helper function allows calling functions only by the master rank, e.g. if you want to write to file by only one process. This function must be called within the PETScContextManager.
- Returns:
True if this is rank 0 on PETSC_COMM_WORLD
- Return type:
bool
- class madupite.PETScContextManager
Bases:
object
Use this context manager to handle initializing and finalizing the PETSc/MPI execution environment.
with madupite.PETScContextManager(): madupite.PyMdp() ...
- class madupite.PyMdp
Bases:
object
This class represents a single MDP.
- generateInitCost(injson='')
Generate an initial cost according to the options given to the PyMdp instance. Options can be passed either one by one (mdp[“discount factor”]=0.9) or by loading an entire .json file.
- Parameters:
injson (String) – provide the filename to a .json file that contains the options for the solver. If the existing options should be used pass the empty string: “”.
- loadG(filename)
Load the stage cost matrix from a PetscBinary file. This file format can be created with madupite.writePETScBinary()
- Parameters:
filename (str) – input filename
- loadP(filename)
Load the transition probability matrix from a PetscBinary file. This file format can be created with madupite.writePETScBinary()
- Parameters:
filename (str) – input filename
- solve(outdir='', injson='', verbose=False)
This calls the solver according to the options given to the PyMdp instance Options can be passed either one by one (mdp[“discount factor”]=0.9) or by loading an entire .json file. An initial cost will be generated if none has been generated before.
- Parameters:
outdir (str, optional) – directory in which the output files should be stored, by default “.”
injson (str, optional) – provide the filename to a .json file that contains the options for the solver. If the existing options should be reused pass the empty string: “”, by default “”
verbose (bool, optional) – verbose option, by default False
- Returns:
Total execution time of solver
- Return type:
float
- madupite.generateArtificialMDP(nstates, mactions, transition_rate, seed=None)
Generate an artificial MDP, whose nonzero transition probabilities are sampled from a uniform distribution. Returns the transition probability in the correct format for PyMdp.loadP().
- Parameters:
nstates (int) – Number of states.
mactions (int) – Number of actions
transition_rate (float in (0,1]) – average rate of transitions to other states. 0.3 means that on average a state has a nonzero probability of transitioning to 30% of the other states.
seed (int, optional) – random seed, by default None
- Returns:
(transition probability matrix, stage cost matrix)
- Return type:
Tuple
- madupite.generateMDP(nstates, mactions, probfunction, costfunction, transprobfilename, stagecostfilename, verbose=False)
Generate the transition probability tensor and stage cost matrix from a probability function and a stagecost function. Output in a PetscBinary file which can be read by the solver. This function must be wrapped in a if madupite.MPI_master() statement. Note: defining probfunction in a numba compatible way provides a significant speedup.
- Parameters:
nstates (int) – Number of states
mactions (int) – Number of actions
probfunction (func) – A function that returns the transition probability to every other state given a state-action pair: \(\mathcal{P}(s\vert s,a)\). Since this is often sparse probfunction(state,action) should return a tuple which consists of an array of the indices of nonzero entries and an array of the values of the nonzero entries.
costfunction (func) – A function that returns the stage cost for a state-action pair. costfunction(state, action) returns a scalar
transprobfilename (str) – output filename for the transition probability file
stagecostfilename (str) – output filename for the stage cost file
verbose (bool, optional) – print potential numba compile error, by default False
- madupite.writePETScBinary(matrix, filename)
Write numpy/scipy matrix as petsc binary sparse format to file https://petsc.org/release/manualpages/Mat/MatLoad/#notes
- Parameters:
matrix (numpy/scipy matrix) – any matrix type that allows calling scipy.sparse.csr_array(matrix)
filename (string) – output filename