Solvers

These classes were developed to help the user control the execution of the solver used in their simulation and to handle certain pygeos wrappers easily. Depending on the specifics of the solver targeted, methods have been added to acces the values of GEOS fields and create outputs of different formats.

The base class for every other solver class is called Solver.

Solver

Solver class which is the base class for every other OtherSolver classes. The driving methods for pygeosx such as initialize and execute, and get/set methods for pygeosx wrappers are defined.

Warning

This does not handle coupled solvers simulations.

Methods ending with the name “For1RegionWith1CellBlock” are designed to only work when dealing with mesh containing only hexahedral elements and only 1 CellElementRegion. In every other cases, do not use these methods.

Todo

If possible, add the capabilities to handle coupled solvers.

class geos.pygeos_tools.solvers.Solver.Solver(solverType, **kwargs)[source]

Bases: object

Solver class containing the main methods of a GEOS solver

alreadyInitialized

Tells if the Solver has been initialized

Type:

bool

collections

geosx group

Type:

List[ pygeosx.Group ]

collectionsTargets

Name of TimeHistory

Type:

List[ str ]

dt

Time step for simulation

Type:

float

geosx

Problem group

Type:

pygeosx.Group

geosxArgs

Object containing GEOSX launching options

Type:

GeosxArgs

hdf5Outputs

geosx group

Type:

List[ pygeosx.Group ]

hdf5Targets

Outputs of TimeHistory

Type:

List[ str ]

maxTime

End time of the simulation

Type:

float

meshName

Name of the mesh in the GEOS XML deck

Type:

str

minTime

Begin time of the simulation

Type:

float

name

Name of the solver defined in the GEOS XML deck

Type:

str

timeVariables

Possible time variables found in the GEOS XML deck that are link to the solver

Type:

Dict[ str, float ]

type

The solverType targeted in GEOS XML deck

Type:

str

vtkOutputs

geosx group

Type:

List[ pygeosx.Group ]

vtkTargets

Outputs of vtk

Type:

List[ str ]

xml

XML object

Type:

XML

applyInitialConditions()[source]

Apply the initial conditions after GEOS (re)initialization

bcastFieldFor1RegionWith1CellBlock(fullField, comm, root=0, **kwargs)[source]

Broadcast a field to local ranks with GEOS local to global map WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:
  • fullField (numpy array) – Full field

  • comm (MPI.COMM_WORLD) – MPI communicator

  • root (int) – MPI rank used for the gather Default is rank 0

Returns:

field – Local field

Return type:

numpy array

cleanup(*args, **kwargs)
execute(*args, **kwargs)
filterGhostRankFor1RegionWith1CellBlock(field, **kwargs)[source]

Filter the ghost rank from a GEOS field WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:

field (numpy array) – Field to filter

Returns:

field – Filtered field

Return type:

numpy array

finalize()[source]

Terminate GEOSX

gatherFieldFor1RegionWith1CellBlock(field, comm, root=0, **kwargs)[source]

Gather a full GEOS field from all local ranks WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:
  • field (numpy array) – Local field

  • comm (MPI.COMM_WORLD) – MPI communicator

  • root (int) – MPI rank used for the gather Default is rank 0

getAllGeosWrapperByName(*args, **kwargs)
getCollections(*args, **kwargs)
getDiscretization(*args, **kwargs)
getDt(*args, **kwargs)
getElementCenterFor1RegionWith1CellBlock(filterGhost=False, **kwargs)[source]

Get element center position as numpy array WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:

filterGhost (bool)

Returns:

elementCenter – Element center coordinates

Return type:

array-like

getElementCenterZFor1RegionWith1CellBlock(filterGhost=False, **kwargs)[source]

Get the z coordinate of the element center WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:

filterGhost (str)

Returns:

elementCenterZ – Element center z coordinates

Return type:

array-like

getGeosWrapperByName(*args, **kwargs)
getGeosWrapperByTargetKey(*args, **kwargs)
getGeosx(*args, **kwargs)
getGhostRankFor1RegionWith1CellBlock(**kwargs)[source]

Get the local ghost ranks WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:

filters (list(str)) – Keywords that can be used to restrict more the research of ‘ghostRank’ field in GEOS.

Returns:

ghostRank – Local ghost ranks

Return type:

npt.NDArray

getHdf5Outputs(*args, **kwargs)
getLocalToGlobalMapFor1RegionWith1CellBlock(filterGhost=False, **kwargs)[source]

Get the local rank element id list WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:

filterGhost (str)

Returns:

Numpy Array

Return type:

Array containing the element id list for the local rank

getMaxTime(*args, **kwargs)
getMeshName(*args, **kwargs)
getMinTime(*args, **kwargs)
getName(*args, **kwargs)
getSolverFieldWithPrefix(fieldName, **kwargs)[source]

Get the requested field as numpy array. WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:

fieldName (str) – Name of the field in GEOSX

Returns:

field – Field requested

Return type:

npt.NDArray

getTargetRegions(*args, **kwargs)
getTimeVariables(*args, **kwargs)
getType(*args, **kwargs)
getVtkOutputs(*args, **kwargs)
initialize(rank=0, xml=None)[source]

Initialization or reinitialization of GEOS that will update these parameters: - the solver name stored in self.name - the solver pygeosx.Group stored in self.solver - the discretization used by the solver - the name of the mesh - the regions targeted by the solver - the different possible outputs which are self.collections, self.hdf5Outputs, self.vtkOutputs - the available time variables defined in the XML and that are relevant to the current solver

Parameters:
  • rank (int) – Process rank

  • xml (XML) – XML object containing parameters for GEOS initialization. Only required if not set in the __init__ OR if different from it

outputVtk(*args, **kwargs)
reinitSolver(*args, **kwargs)
setDt(value)[source]
setDtFromTimeVariable(*args, **kwargs)
setGeosWrapperValueByName(*args, **kwargs)
setGeosWrapperValueByTargetKey(*args, **kwargs)
setHdf5OutputsName(*args, **kwargs)
setMaxTime(new_maxTime)[source]
setMinTime(new_minTime)[source]
setTimeVariable(*args, **kwargs)
setVtkOutputsName(*args, **kwargs)
setXml(xml)[source]

Sets the new XML object.

Parameters:

xml (XML) – XML object corresponding to GEOSX input

updateDiscretization(*args, **kwargs)
updateMeshName(*args, **kwargs)
updateOutputs(*args, **kwargs)
updateSolverGroup(*args, **kwargs)
updateSolverName(*args, **kwargs)
updateTargetRegions(*args, **kwargs)
updateTimeVariables(*args, **kwargs)

WaveSolver

WaveSolver class which is the base class for every AcousticSolver and ElasticSolver classes, and inherits the Solver capabilities.

This adds more methods to the Solver class to handle seismic sources and receivers.

class geos.pygeos_tools.solvers.WaveSolver.WaveSolver(solverType, dt=None, minTime=0.0, maxTime=None, dtSeismo=None, dtWaveField=None, sourceType=None, sourceFreq=None, **kwargs)[source]

Bases: Solver

WaveSolver Object containing methods useful for simulation using wave solvers in GEOS

The ones herited from Solver class and
dtSeismo

Time step to save pressure for seismic trace

Type:

float

dtWaveField

Time step to save fields

Type:

float

minTime

Min time to consider

Type:

float

maxTime

End time to consider

Type:

float

minTimeSim

Starting time of simulation

Type:

float

maxTimeSim

End Time of simulation

Type:

float

sourceType

Type of source

Type:

str

sourceFreq

Frequency of the source

Type:

float

Parameters:
  • solverType (str) – The solverType targeted in GEOS XML deck

  • dt (float) – Time step for simulation

  • minTime (float) – Starting time of simulation Default is 0

  • maxTime (float) – End Time of simulation

  • dtSeismo (float) – Time step to save pressure for seismic trace

  • dtWaveField (float) – Time step to save fields

  • sourceType (str) – Type of source Default is None

  • sourceFreq (float) – Frequency of the source Default is None

  • kwargs (keyword args) –

    geosx_argvlist

    GEOSX arguments or command line as a splitted line

evaluateSource()[source]

Evaluate source and update on GEOS Only ricker order {0 - 4} accepted

filterSource(fmax)[source]

Filter the source value and give the value to GEOSX. Note that is can also modify the start and end time of simulation to avoid discontinuity.

Parameters:

fmax (float/string) – Max frequency of the source wanted. The source then have frequencies in the interval [0,fmax+1]

getVelocityModel(velocityName, filterGhost=False, **kwargs)[source]

Get the velocity values WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:
  • velocityName (str) – Name of velocity array in GEOS

  • filterGhost (bool) – Filter the ghost ranks

Returns:

Numpy Array

Return type:

Array containing the velocity values

initialize(rank=0, xml=None)[source]

Initialization or reinitialization of GEOSX

Parameters:
  • rank (int) – Process rank

  • xml (XML) – XML object containing parameters for GEOSX initialization. Only required if not set in the __init__ OR if different from it

outputWaveField(time, cycleNumber)[source]

Trigger the wavefield output

Parameters:
  • time (float) – Current time of simulation

  • cycleNumber (int) – Current cycle number

setSourceAndReceivers(sourcesCoords=[], receiversCoords=[])[source]

Update sources and receivers positions in GEOS

Parameters:
  • sourcesCoords (list) – List of coordinates for the sources

  • receiversCoords (list) – List of coordinates for the receivers

setSourceFrequency(freq)[source]

Overwrite GEOSX source frequency and set self.sourceFreq

Parameters:

freq (float) – Frequency of the source in Hz

setSourceValue(value)[source]

Set the value of the source in GEOS

Parameters:

value (array/list) – List/array containing the value of the source at each time step

updateSourceProperties()[source]

Updates the frequency and type of source to match the XML

AcousticSolver

AcousticSolver class inherits from WaveSolver class.

This adds method to read / set pressures at receiver, compute partial gradients and accessors for Density and Velocity models.

class geos.pygeos_tools.solvers.AcousticSolver.AcousticSolver(solverType='AcousticSEM', dt=None, minTime=0.0, maxTime=None, dtSeismo=None, dtWaveField=None, sourceType=None, sourceFreq=None, **kwargs)[source]

Bases: WaveSolver

AcousticSolver Object containing all methods to run AcousticSEM simulation with GEOSX

The ones inherited from WaveSolver class
modelForGradient

Gradient model used

Type:

str

Parameters:
  • solverType (str) – The solverType targeted in GEOS XML deck. Defaults to “AcousticSEM”

  • dt (float) – Time step for simulation

  • minTime (float) – Starting time of simulation Default is 0

  • maxTime (float) – End Time of simulation

  • dtSeismo (float) – Time step to save pressure for seismic trace

  • dtWaveField (float) – Time step to save fields

  • sourceType (str) – Type of source Default is None

  • sourceFreq (float) – Frequency of the source Default is None

  • kwargs (keyword args) –

    geosx_argvlist

    GEOSX arguments or command line as a splitted line

computePartialGradientFor1RegionWith1CellBlock(shotId, minDepth, comm, velocityName, gradDirectory='partialGradient', filterGhost=False, **kwargs)[source]

Compute the partial Gradient WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:
  • shotId (string) – Number of the shot as string

  • minDepth (float) – Depth at which gradient values are kept, otherwise it is set to 0. NOTE : this is done in this routine to avoid storage of elementCenter coordinates in the .hdf5 but might be problem for WolfeConditions later on if minDepth is too large

  • comm (MPI_COMM) – MPI communicators

  • velocity (str) – Name of the velocity model in GEOS

  • gradDirectory (str, optional) – Partial gradient directory Default is partialGradient

getFullPressureAtReceivers(comm)[source]

Return all pressures at receivers values on all ranks Note that for a too large 2d array this may not work.

Parameters:

commMPI_COMM

MPI communicators

getFullWaveFieldAtReceivers(comm)[source]
getModelForGradient(*args, **kwargs)
getPartialGradientFor1RegionWith1CellBlock(filterGhost=False, **kwargs)[source]

Get the local rank gradient value WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Returns:

partialGradient – Array containing the element id list for the local rank

Return type:

Numpy Array-

getPressureAtReceivers()[source]

Get the pressure values at receivers coordinates

Returns:

numpy Array

Return type:

Array containing the pressure values at all time step at all receivers coordinates

getWaveField()[source]
resetPressureAtReceivers()[source]

Reinitialize pressure values at receivers to 0

resetWaveField(**kwargs)[source]

Reinitialize all pressure values on the Wavefield to zero in GEOSX

setModelForGradient(modelForGradient)[source]
updateDensityModel(density)[source]

Update density values in GEOS

Parameters:

density (array) – New values for the density

updateVelocityModel(vel)[source]

Update velocity value in GEOS

Parameters:

vel (float/array) – Value(s) for velocity field

updateVelocityModelFromHDF5(filename, low, high, comm, velocityModelName, **kwargs)[source]

Update velocity model WARNING: this function aims to work in the specific case of having only 1 CellElementRegion in your XML file and that this CellElementRegion contains only one cellBlock.

Parameters:
  • filename (str) – .hdf5 file where to get the new model

  • low (float) – Min value threshold. All new values < low are set to low

  • high (float) – Max value threshold. All new values > high are set to high

  • comm (MPI_COMM) – MPI communicators

ElasticSolver

AcousticSolver class inherits from WaveSolver class.

This adds method to read / set displacements at receiver.

class geos.pygeos_tools.solvers.ElasticSolver.ElasticSolver(solverType='ElasticSEM', dt=None, minTime=0.0, maxTime=None, dtSeismo=None, dtWaveField=None, sourceType=None, sourceFreq=None, **kwargs)[source]

Bases: WaveSolver

ElasticSolver Object containing all methods to run ElasticSEM simulation with GEOSX

The ones inherited from WaveSolver class
Parameters:
  • solverType (str) – The solverType targeted in GEOS XML deck. Defaults to “ElasticSEM”

  • dt (float) – Time step for simulation

  • minTime (float) – Starting time of simulation Default is 0

  • maxTime (float) – End Time of simulation

  • dtSeismo (float) – Time step to save pressure for seismic trace

  • dtWaveField (float) – Time step to save fields

  • sourceType (str) – Type of source Default is None

  • sourceFreq (float) – Frequency of the source Default is None

  • kwargs (keyword args) –

    geosx_argvlist

    GEOSX arguments or command line as a splitted line

getAllDisplacementAtReceivers()[source]

Get the displacement for the x, y and z directions at all time step and all receivers coordinates

Returns:

  • displacementX (numpy array) – Component X of the displacement

  • displacementY (numpy array) – Component Y of the displacement

  • displacementZ (numpy array) – Component Z of the displacement

getDASSignalAtReceivers()[source]

Get the DAS signal values at receivers coordinates

Returns:

dassignal – Array containing the DAS signal values at all time step at all receivers coordinates

Return type:

numpy array

getDisplacementAtReceivers(component='X')[source]

Get the displacement values at receivers coordinates for a given direction

Returns:

displacement – Array containing the displacements values at all time step at all receivers coordinates

Return type:

numpy array

getWaveField()[source]
initialize(rank=0, xml=None)[source]

Initialization or reinitialization of GEOSX

Parameters:
  • rank (int) – Process rank

  • xml (XML) – XML object containing parameters for GEOSX initialization. Only required if not set in the __init__ OR if different from it

resetDisplacementAtReceivers()[source]

Reinitialize displacement values at receivers to 0

resetWaveField(**kwargs)[source]

Reinitialize all displacement values on the Wavefield to zero in GEOSX

updateDensityModel(density)[source]

Update density values in GEOS

Parameters:

density (array) – New values for the density

updateVelocityModel(vel, component)[source]

Update velocity value in GEOS

Parameters:
  • vel (float/array) – Value(s) for velocity field

  • component (str) – Vs or Vp

GeomechanicsSolver

AcousticSolver class inherits from Solver class.

This adds accessor methods for stresses, displacements and geomechanics parameters.

class geos.pygeos_tools.solvers.GeomechanicsSolver.GeomechanicsSolver(solverType, dt=None, maxTime=None, **kwargs)[source]

Bases: Solver

Geomechanics solver object containing methods to run geomechanics simulations in GEOS

The ones inherited from Solver class
Parameters:
  • solverType (str) – The solver used in the XML file

  • initialDt (float) – Initial time step Default is None

  • maxTime (float) – End time of simulation Default is None

getBulkModulus()[source]

Get the local bulk modulus for each CellElementRegion and its cellBlocks of the mesh.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each. The first 2 regions are using “shale”,

  • the last region uses “sand”, the result is

  • { “region1/block1/shale” (npt.NDArray, “region1/block1/shale”: npt.NDArray,) – “region1/block2/shale”: npt.NDArray, “region1/block1/shale”: npt.NDArray, “region2/block3/shale”: npt.NDArray, “region2/block1/shale”: npt.NDArray, “region2/block4/shale”: npt.NDArray, “region2/block1/shale”: npt.NDArray, “region3/block5/sand”: npt.NDArray, “region3/block1/sand”: npt.NDArray, “region3/block6/sand”: npt.NDArray, “region3/block1/sand”: npt.NDArray }

getConstitutiveModelData(modelName)[source]

Get the local constitutive model data for each CellElementRegion and its cellBlocks of the mesh.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each. The first 2 regions are using the constituve

  • ”model1”, the last region uses “model2”, the result is

  • { “region1/block1/model1” (npt.NDArray, “region1/block1/model1”: npt.NDArray,) – “region1/block2/model1”: npt.NDArray, “region1/block1/model1”: npt.NDArray, “region2/block3/model1”: npt.NDArray, “region2/block1/model1”: npt.NDArray, “region2/block4/model1”: npt.NDArray, “region2/block1/model1”: npt.NDArray, “region3/block5/model2”: npt.NDArray, “region3/block1/model2”: npt.NDArray, “region3/block6/model2”: npt.NDArray, “region3/block1/model2”: npt.NDArray }

getDensities()[source]

Get the local density for each CellElementRegion and its cellBlocks of the mesh.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each. The first 2 regions are using “shale”,

  • the last region uses “sand”, the result is

  • { “region1/block1/shale” (npt.NDArray, “region1/block1/shale”: npt.NDArray,) – “region1/block2/shale”: npt.NDArray, “region1/block1/shale”: npt.NDArray, “region2/block3/shale”: npt.NDArray, “region2/block1/shale”: npt.NDArray, “region2/block4/shale”: npt.NDArray, “region2/block1/shale”: npt.NDArray, “region3/block5/sand”: npt.NDArray, “region3/block1/sand”: npt.NDArray, “region3/block6/sand”: npt.NDArray, “region3/block1/sand”: npt.NDArray }

getShearModulus()[source]

Get the local shear modulus for each CellElementRegion and its cellBlocks of the mesh.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each. The first 2 regions are using “shale”,

  • the last region uses “sand”, the result is

  • { “region1/block1/shale” (npt.NDArray, “region1/block1/shale”: npt.NDArray,) – “region1/block2/shale”: npt.NDArray, “region1/block1/shale”: npt.NDArray, “region2/block3/shale”: npt.NDArray, “region2/block1/shale”: npt.NDArray, “region2/block4/shale”: npt.NDArray, “region2/block1/shale”: npt.NDArray, “region3/block5/sand”: npt.NDArray, “region3/block1/sand”: npt.NDArray, “region3/block6/sand”: npt.NDArray, “region3/block1/sand”: npt.NDArray }

getStresses()[source]

Get the local stresses for each CellElementRegion and its cellBlocks of the mesh.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each. The first 2 regions are using “shale”,

  • the last region uses “sand”, the result is

  • { “region1/block1/shale” (npt.NDArray, “region1/block1/shale”: npt.NDArray,) – “region1/block2/shale”: npt.NDArray, “region1/block1/shale”: npt.NDArray, “region2/block3/shale”: npt.NDArray, “region2/block1/shale”: npt.NDArray, “region2/block4/shale”: npt.NDArray, “region2/block1/shale”: npt.NDArray, “region3/block5/sand”: npt.NDArray, “region3/block1/sand”: npt.NDArray, “region3/block6/sand”: npt.NDArray, “region3/block1/sand”: npt.NDArray }

getTotalDisplacement()[source]

Get the local totalDipslacements from the nodes.

Return type:

npt.NDArray of totalDipslacements, shape = ( number of nodes, 3 )

initialize(rank=0, xml=None)[source]

Initialization or reinitialization of GEOS that will update these parameters: - the solver name stored in self.name - the solver pygeosx.Group stored in self.solver - the discretization used by the solver - the name of the mesh - the regions targeted by the solver - the different possible outputs which are self.collections, self.hdf5Outputs, self.vtkOutputs - the available time variables defined in the XML and that are relevant to the current solver

Parameters:
  • rank (int) – Process rank

  • xml (XML) – XML object containing parameters for GEOS initialization. Only required if not set in the __init__ OR if different from it

ReservoirSolver

ReservoirSolver class inherits from Solver class.

This adds accessor methods for pressure, phase volume fractions and porosities.

class geos.pygeos_tools.solvers.ReservoirSolver.ReservoirSolver(solverType, initialDt=None, maxTime=None, **kwargs)[source]

Bases: Solver

Reservoir solver object containing methods to run reservoir simulations in GEOS

The ones inherited from Solver class
Parameters:
  • solverType (str) – The solver used in the XML file

  • initialDt (float) – Initial time step Default is None

  • maxTime (float) – End time of simulation Default is None

getDeltaPressures()[source]

Get the local delta pressure for each CellElementRegion and each cellBlocks of the mesh.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each, the result is

  • { “region1/block1” (npt.NDArray, “region1/block2”: npt.NDArray,) – “region2/block3”: npt.NDArray, “region2/block4”: npt.NDArray, “region3/block5”: npt.NDArray, “region3/block6”: npt.NDArray }

getPhaseVolumeFractions()[source]

Get the local phaseVolumeFraction for each CellElementRegion and each cellBlocks of the mesh and each phase.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each, and two phases ‘oil’ and ‘gas’ give the result

  • { “region1/block1/oil” (npt.NDArray, “region1/block1/gas”: npt.NDArray,) – “region1/block2/oil”: npt.NDArray, “region1/block1/gas”: npt.NDArray, “region2/block3/oil”: npt.NDArray, “region2/block1/gas”: npt.NDArray, “region2/block4/oil”: npt.NDArray, “region2/block1/gas”: npt.NDArray, “region3/block5/oil”: npt.NDArray, “region3/block1/gas”: npt.NDArray, “region3/block6/oil”: npt.NDArray, “region3/block1/gas”: npt.NDArray }

getPorosities()[source]

Get the local porosity for each CellElementRegion and its cellBlocks of the mesh.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each. The first 2 regions are using “burdenPorosity”,

  • the last region uses “sandPorosity”, the result is

  • { “region1/block1/burdenPorosity” (npt.NDArray, “region1/block1/burdenPorosity”: npt.NDArray,) – “region1/block2/burdenPorosity”: npt.NDArray, “region1/block1/burdenPorosity”: npt.NDArray, “region2/block3/burdenPorosity”: npt.NDArray, “region2/block1/burdenPorosity”: npt.NDArray, “region2/block4/burdenPorosity”: npt.NDArray, “region2/block1/burdenPorosity”: npt.NDArray, “region3/block5/sandPorosity”: npt.NDArray, “region3/block1/sandPorosity”: npt.NDArray, “region3/block6/sandPorosity”: npt.NDArray, “region3/block1/sandPorosity”: npt.NDArray }

getPressures()[source]

Get the local pressure for each CellElementRegion and each cellBlocks of the mesh.

Returns:

  • Dict[ str, npt.NDArray ]

  • If your mesh contains 3 regions with 2 cellBlocks in each, the result is

  • { “region1/block1” (npt.NDArray, “region1/block2”: npt.NDArray,) – “region2/block3”: npt.NDArray, “region2/block4”: npt.NDArray, “region3/block5”: npt.NDArray, “region3/block6”: npt.NDArray }