21 #ifndef GEOS_PHYSICSSOLVERS_MULTIPHYSICS_COUPLEDRESERVOIRANDWELLSBASE_HPP_
22 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_COUPLEDRESERVOIRANDWELLSBASE_HPP_
27 #include "constitutive/permeability/PermeabilityFields.hpp"
28 #include "constitutive/permeability/PermeabilityBase.hpp"
31 #include "physicsSolvers/fluidFlow/wells/WellControls.hpp"
37 namespace coupledReservoirAndWellsInternal
57 string const & resElemDofName,
58 string const & wellElemDofName );
73 template<
typename RESERVOIR_SOLVER,
typename WELL_SOLVER >
102 :
Base( name, parent ),
105 this->
template getWrapper< string >( Base::viewKeyStruct::discretizationString() ).
127 bool const setSparsity =
true )
override
132 reservoirSolver()->setupSystem( domain, dofManager, localMatrix, rhs, solution, setSparsity );
145 for(
localIndex localRow = 0; localRow < patternDiag.numRows(); ++localRow )
147 rowLengths[localRow] = patternDiag.numNonZeros( localRow );
155 pattern.resizeFromRowCapacities< parallelHostPolicy >( patternDiag.numRows(), patternDiag.numColumns(), rowLengths.data() );
158 for(
localIndex localRow = 0; localRow < patternDiag.numRows(); ++localRow )
160 globalIndex const * cols = patternDiag.getColumns( localRow ).dataIfContiguous();
161 pattern.insertNonZeros( localRow, cols, cols + patternDiag.numNonZeros( localRow ) );
168 localMatrix.assimilate< parallelDevicePolicy<> >( std::move( pattern ) );
169 localMatrix.setName( this->
getName() +
"/localMatrix" );
171 rhs.setName( this->
getName() +
"/rhs" );
174 solution.setName( this->
getName() +
"/solution" );
199 DomainPartition & domain = this->
template getGroupByPath< DomainPartition >(
"/Problem/domain" );
202 if( !validateWellPerforations( domain ))
204 GEOS_ERROR( GEOS_FMT(
"{}: well perforations validation failed, bad perforations found", this->
getName()));
229 computeWellTransmissibility( domain );
235 assembleFluxTerms(
real64 const dt,
240 {
reservoirSolver()->assembleFluxTerms( dt, domain, dofManager, localMatrix, localRhs ); }
243 assembleStabilizedFluxTerms(
real64 const dt,
244 DomainPartition
const & domain,
245 DofManager
const & dofManager,
246 CRSMatrixView< real64, globalIndex const >
const & localMatrix,
247 arrayView1d< real64 >
const & localRhs )
const
248 {
reservoirSolver()->assembleStabilizedFluxTerms( dt, domain, dofManager, localMatrix, localRhs ); }
250 real64 updateFluidState( ElementSubRegionBase & subRegion )
const
252 void updatePorosityAndPermeability( CellElementSubRegion & subRegion )
const
254 void updateSolidInternalEnergyModel( ObjectManagerBase & dataGroup )
const
259 void enableJumpStabilization()
262 void enableFixedStressPoromechanicsUpdate()
265 void setKeepVariablesConstantDuringInitStep(
bool const keepVariablesConstantDuringInitStep )
267 reservoirSolver()->setKeepVariablesConstantDuringInitStep( keepVariablesConstantDuringInitStep );
268 wellSolver()->setKeepVariablesConstantDuringInitStep( keepVariablesConstantDuringInitStep );
287 coupledReservoirAndWellsInternal::
288 addCouplingNumNonzeros(
this,
309 virtual void setMGRStrategy()
340 ElementRegionManager & elemManager = meshLevel.getElemManager();
342 ElementRegionManager::ElementViewAccessor< arrayView2d< real64 > > const elemCenter =
343 elemManager.constructViewAccessor< array2d< real64 >, arrayView2d< real64 > >( ElementSubRegionBase::viewKeyStruct::elementCenterString() );
346 elemManager.forElementSubRegions< WellElementSubRegion >( regionNames, [&]( localIndex const,
347 WellElementSubRegion & subRegion )
349 array1d< array1d< arrayView3d< real64 const > > > const permeability =
350 elemManager.constructMaterialFieldAccessor< constitutive::PermeabilityBase,
351 fields::permeability::permeability >();
353 PerforationData & perforationData = *subRegion.getPerforationData();
354 WellControls const & wellControls = wellSolver()->getWellControls( subRegion );
357 perforationData.computeWellTransmissibility( meshLevel, subRegion, permeability );
360 if( wellControls.getLogLevel() >= 2 )
362 arrayView2d< real64 const > const perfLocation =
363 perforationData.getField< fields::perforation::location >();
364 arrayView1d< real64 const > const perfTrans =
365 perforationData.getField< fields::perforation::wellTransmissibility >();
368 arrayView1d< localIndex const > const resElemRegion =
369 perforationData.getField< fields::perforation::reservoirElementRegion >();
370 arrayView1d< localIndex const > const resElemSubRegion =
371 perforationData.getField< fields::perforation::reservoirElementSubRegion >();
372 arrayView1d< localIndex const > const resElemIndex =
373 perforationData.getField< fields::perforation::reservoirElementIndex >();
375 GEOS_UNUSED_VAR( perfLocation );
376 GEOS_UNUSED_VAR( perfTrans );
377 GEOS_UNUSED_VAR( resElemRegion );
378 GEOS_UNUSED_VAR( resElemSubRegion );
379 GEOS_UNUSED_VAR( resElemIndex );
381 forAll< serialPolicy >( perforationData.size(), [=] ( localIndex const iperf )
383 GEOS_UNUSED_VAR( iperf );
384 GEOS_LOG_RANK( GEOS_FMT(
"{}: perforation at ({},{},{}), perforated element center = ({},{},{}), transmissibility = {} [{}]",
385 this->getName(), perfLocation[iperf][0], perfLocation[iperf][1], perfLocation[iperf][2],
386 elemCenter[resElemRegion[iperf]][resElemSubRegion[iperf]][resElemIndex[iperf]][0],
387 elemCenter[resElemRegion[iperf]][resElemSubRegion[iperf]][resElemIndex[iperf]][1],
388 elemCenter[resElemRegion[iperf]][resElemSubRegion[iperf]][resElemIndex[iperf]][2],
389 perfTrans[iperf], getSymbol( units::Transmissibility ) ) );
bool validateWellPerforations(PhysicsSolverBase const *const reservoirSolver, WellSolverBase const *const wellSolver, DomainPartition const &domain)
Validate the well perforations ensuring that each perforation is located in a reservoir region that i...
void addCouplingNumNonzeros(PhysicsSolverBase const *const solver, DomainPartition &domain, DofManager &dofManager, arrayView1d< localIndex > const &rowLengths, integer const resNumDof, integer const wellNumDof, string const &resElemDofName, string const &wellElemDofName)
Utility function for the implementation details of addCouplingNumZeros.
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
bool m_isWellTransmissibilityComputed
Flag to determine whether the well transmissibility needs to be computed.
virtual void initializePostInitialConditionsPreSubGroups() override
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
virtual void addCouplingSparsityPattern(DomainPartition const &domain, DofManager const &dofManager, SparsityPatternView< globalIndex > const &pattern) const =0
static string coupledSolverAttributePrefix()
String used to form the solverName used to register solvers in CoupledSolver.
CoupledReservoirAndWellsBase(const string &name, dataRepository::Group *const parent)
main constructor for ManagedGroup Objects
virtual void postInputInitialization() override
virtual ~CoupledReservoirAndWellsBase() override
default destructor
WELL_SOLVER * wellSolver() const
accessor for the pointer to the well solver
void addCouplingNumNonzeros(DomainPartition &domain, DofManager &dofManager, arrayView1d< localIndex > const &rowLengths) const
virtual void saveSequentialIterationState(DomainPartition &domain) override
Save the state of the solver for sequential iteration.
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override
function to perform setup for implicit timestep
RESERVOIR_SOLVER * reservoirSolver() const
accessor for the pointer to the reservoir solver
virtual void postInputInitialization() override
std::array< string, sizeof...(SOLVERS) > m_names
Names of the single-physics solvers.
std::tuple< SOLVERS *... > m_solvers
Pointers of the single-physics solvers.
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
void setSparsityPattern(SparsityPattern< globalIndex > &pattern) const
Populate sparsity pattern of the entire system matrix.
void reorderByRank()
Finish populating fields and apply appropriate dof renumbering.
localIndex numLocalDofs(string const &fieldName) const
void setDomain(DomainPartition &domain)
Assign a domain.
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
Group const & getMeshBodies() const
Get the mesh bodies, const version.
Class facilitating the representation of a multi-level discretization of a MeshBody.
Base class for all physics solvers.
virtual string getCatalogName() const =0
void forDiscretizationOnMeshTargets(Group const &meshBodies, LAMBDA &&lambda) const
Loop over the target discretization on all mesh targets and apply callback.
CRSMatrix< real64, globalIndex > m_localMatrix
Local system matrix and rhs.
DofManager m_dofManager
Data structure to handle degrees of freedom.
ParallelVector m_solution
System solution vector.
ParallelVector m_rhs
System right-hand side vector.
LinearSolverParametersInput m_linearSolverParameters
Linear solver parameters.
string const & getName() const
Get group name.
virtual void initializePostInitialConditionsPreSubGroups()
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override
virtual void setupSystem(DomainPartition &domain, DofManager &dofManager, CRSMatrix< real64, globalIndex > &localMatrix, ParallelVector &rhs, ParallelVector &solution, bool const setSparsity=true) override
Set up the linear system (DOF indices and sparsity patterns)
void setupDofs(DomainPartition const &domain, DofManager &dofManager) const override
@ FALSE
Not read from input.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
int MPI_COMM_GEOS
Global MPI communicator used by GEOSX.
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
LvArray::SparsityPatternView< COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > SparsityPatternView
Alias for Sparsity pattern View.
LvArray::SparsityPattern< COL_INDEX, INDEX_TYPE, LvArray::ChaiBuffer > SparsityPattern
Alias for Sparsity pattern class.
double real64
64-bit floating point type.
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
std::int32_t integer
Signed integer type.
LvArray::CRSMatrix< T, COL_INDEX, localIndex, LvArray::ChaiBuffer > CRSMatrix
Alias for CRS Matrix class.
Array< T, 1 > array1d
Alias for 1D array.
LAInterface::ParallelVector ParallelVector
Alias for ParallelVector.
@ mgr
Multigrid reduction (Hypre only)
PreconditionerType preconditionerType
Preconditioner type.