21 #ifndef GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSSOLVER_HPP_
22 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSSOLVER_HPP_
28 #include "constitutive/solid/CoupledSolidBase.hpp"
29 #include "constitutive/solid/PorousSolid.hpp"
30 #include "constitutive/contact/HydraulicApertureBase.hpp"
33 #include "codingUtilities/Utilities.hpp"
38 namespace stabilization
40 enum class StabilizationType :
integer
54 template<
typename FLOW_SOLVER,
typename MECHANICS_SOLVER = Sol
idMechanicsLagrangianFEM >
82 :
Base( name, parent ),
86 setApplyDefaultValue( 0 ).
88 setDescription(
"Flag indicating whether the problem is thermal or not. Set isThermal=\"1\" to enable the thermal coupling" );
91 setApplyDefaultValue(
false ).
93 setDescription(
"Flag to indicate that the solver is going to perform stress initialization" );
97 setDescription(
"StabilizationType. Options are:\n" +
98 toString( stabilization::StabilizationType::None ) +
"- Add no stabilization to mass equation \n" +
99 toString( stabilization::StabilizationType::Global ) +
"- Add jump stabilization to all faces \n" +
100 toString( stabilization::StabilizationType::Local ) +
"- Add jump stabilization on interior of macro elements" );
103 setRTTypeName( rtTypes::CustomTypes::groupNameRefArray ).
105 setDescription(
"Regions where stabilization is applied." );
108 setApplyDefaultValue( 1.0 ).
110 setDescription(
"Constant multiplier of stabilization strength" );
118 GEOS_FMT(
"{} {}: The attribute `{}` of the flow solver `{}` must be set to 1 since the poromechanics solver is thermal",
133 hydraulicApertureModelName = PhysicsSolverBase::getConstitutiveName< constitutive::HydraulicApertureBase >( subRegion );
134 GEOS_ERROR_IF( hydraulicApertureModelName.empty(), GEOS_FMT(
"{}: HydraulicApertureBase model not found on subregion {}",
146 ": Local stabilization has been temporarily disabled",
149 DomainPartition & domain = this->
template getGroupByPath< DomainPartition >(
"/Problem/domain" );
155 ElementRegionManager & elementRegionManager = mesh.getElemManager();
156 elementRegionManager.forElementSubRegions< CellElementSubRegion >( regionNames,
157 [&]( localIndex const,
158 ElementSubRegionBase & subRegion )
160 string & porousName = subRegion.getReference< string >( viewKeyStruct::porousMaterialNamesString() );
161 porousName = this->template getConstitutiveName< constitutive::CoupledSolidBase >( subRegion );
162 GEOS_THROW_IF( porousName.empty(),
163 GEOS_FMT(
"{} {} : Solid model not found on subregion {}",
164 this->getCatalogName(), this->getDataContext().toString(), subRegion.getName() ),
167 string & porosityModelName = subRegion.getReference< string >( constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString() );
168 porosityModelName = this->template getConstitutiveName< constitutive::PorosityBase >( subRegion );
169 GEOS_THROW_IF( porosityModelName.empty(),
170 GEOS_FMT(
"{} {} : Porosity model not found on subregion {}",
171 this->getCatalogName(), this->getDataContext().toString(), subRegion.getName() ),
174 if( subRegion.hasField< fields::poromechanics::bulkDensity >() )
177 constitutive::CoupledSolidBase const & solid = this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, porousName );
178 subRegion.getField< fields::poromechanics::bulkDensity >().resizeDimension< 1 >( solid.getDensity().size( 1 ) );
191 solidMechanicsSolver()->enableFixedStressPoromechanicsUpdate();
193 flowSolver()->enableFixedStressPoromechanicsUpdate();
196 if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
198 flowSolver()->enableJumpStabilization();
211 subRegion.registerWrapper<
string >( viewKeyStruct::porousMaterialNamesString() ).
214 setSizedFromParent( 0 );
217 subRegion.registerWrapper<
string >( constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString() ).
220 setSizedFromParent( 0 );
226 subRegion.registerField< fields::poromechanics::bulkDensity >( this->getName() );
229 if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
231 subRegion.registerField< fields::flow::macroElementIndex >( this->getName());
232 subRegion.registerField< fields::flow::elementStabConstant >( this->getName());
242 flowSolver()->setKeepVariablesConstantDuringInitStep( m_performStressInitialization );
244 if( this->m_stabilizationType == stabilization::StabilizationType::Global || this->m_stabilizationType == stabilization::StabilizationType::Local )
246 this->updateStabilizationParameters( domain );
249 Base::implicitStepSetup( time_n, dt, domain );
258 solidMechanicsSolver()->setupDofs( domain, dofManager );
259 flowSolver()->setupDofs( domain, dofManager );
260 this->setupCoupling( domain, dofManager );
263 virtual bool checkSequentialConvergence(
int const & iter,
269 auto & subcycling = this->getNonlinearSolverParameters().m_subcyclingOption;
270 auto const subcycling_orig = subcycling;
271 if( m_performStressInitialization )
274 bool isConverged = Base::checkSequentialConvergence( iter, time_n, dt, domain );
277 subcycling = subcycling_orig;
288 return std::get< toUnderlying( SolverType::SolidMechanics ) >( m_solvers );
297 return std::get< toUnderlying( SolverType::Flow ) >( m_solvers );
304 void setStressInitialization(
integer const performStressInitialization )
306 m_performStressInitialization = performStressInitialization;
338 for(
string const & regionName : m_stabilizationRegionNames )
340 regionFilter.insert( regionName );
344 this->
template forDiscretizationOnMeshTargets( domain.
getMeshBodies(), [&] (
string const &,
346 arrayView1d< string const >
const & targetRegionNames )
349 array1d< string > filteredTargetRegionNames;
350 filteredTargetRegionNames.reserve( targetRegionNames.size() );
352 for( string const & targetRegionName : targetRegionNames )
354 if( regionFilter.count( targetRegionName ) )
356 filteredTargetRegionNames.emplace_back( targetRegionName );
361 mesh.getElemManager().forElementSubRegions( filteredTargetRegionNames.toViewConst(), [&](
localIndex const,
362 ElementSubRegionBase & subRegion )
364 arrayView1d< integer > const macroElementIndex = subRegion.getField< fields::flow::macroElementIndex >();
365 arrayView1d< real64 > const elementStabConstant = subRegion.getField< fields::flow::elementStabConstant >();
367 geos::constitutive::CoupledSolidBase const & porousSolid =
368 this->template getConstitutiveModel< geos::constitutive::CoupledSolidBase >( subRegion, subRegion.getReference< string >( viewKeyStruct::porousMaterialNamesString() ) );
370 arrayView1d< real64 const > const bulkModulus = porousSolid.getBulkModulus();
371 arrayView1d< real64 const > const shearModulus = porousSolid.getShearModulus();
372 arrayView1d< real64 const > const biotCoefficient = porousSolid.getBiotCoefficient();
374 real64 const stabilizationMultiplier = m_stabilizationMultiplier;
376 forAll< parallelDevicePolicy<> >( subRegion.size(), [bulkModulus,
379 stabilizationMultiplier,
381 elementStabConstant] GEOS_HOST_DEVICE ( localIndex const ei )
384 real64 const bM = bulkModulus[ei];
385 real64 const sM = shearModulus[ei];
386 real64 const bC = biotCoefficient[ei];
388 macroElementIndex[ei] = 1;
389 elementStabConstant[ei] = stabilizationMultiplier * 9.0 * (bC * bC) / (32.0 * (10.0 * sM / 3.0 + bM));
398 template<
typename CONSTITUTIVE_BASE,
399 typename KERNEL_WRAPPER,
400 typename ... PARAMS >
401 real64 assemblyLaunch( MeshLevel & mesh,
402 DofManager
const & dofManager,
403 arrayView1d< string const >
const & regionNames,
404 string const & materialNamesString,
405 CRSMatrixView< real64, globalIndex const >
const & localMatrix,
406 arrayView1d< real64 >
const & localRhs,
408 PARAMS && ... params )
412 NodeManager
const & nodeManager = mesh.getNodeManager();
414 string const dofKey = dofManager.getKey( fields::solidMechanics::totalDisplacement::key() );
415 arrayView1d< globalIndex const >
const & dofNumber = nodeManager.getReference<
globalIndex_array >( dofKey );
417 real64 const gravityVectorData[3] = LVARRAY_TENSOROPS_INIT_LOCAL_3( PhysicsSolverBase::gravityVector() );
419 KERNEL_WRAPPER kernelWrapper( dofNumber,
420 dofManager.rankOffset(),
425 std::forward< PARAMS >( params )... );
427 return finiteElement::
428 regionBasedKernelApplication< parallelDevicePolicy< >,
430 CellElementSubRegion >( mesh,
432 this->solidMechanicsSolver()->getDiscretizationName(),
439 void recordAverageMeanTotalStressIncrement( DomainPartition & domain,
440 array1d< real64 > & averageMeanTotalStressIncrement )
442 averageMeanTotalStressIncrement.resize( 0 );
443 PhysicsSolverBase::forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&](
string const &,
445 arrayView1d< string const >
const & regionNames ) {
446 mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
449 string const solidName = subRegion.template getReference< string >(
"porousMaterialNames" );
450 constitutive::CoupledSolidBase & solid = PhysicsSolverBase::getConstitutiveModel< constitutive::CoupledSolidBase >(
451 subRegion, solidName );
453 arrayView1d< const real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
454 for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
456 averageMeanTotalStressIncrement.emplace_back( averageMeanTotalStressIncrement_k[k] );
462 void applyAcceleratedAverageMeanTotalStressIncrement( DomainPartition & domain,
463 array1d< real64 > & averageMeanTotalStressIncrement )
466 PhysicsSolverBase::forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&](
string const &,
468 arrayView1d< string const >
const & regionNames ) {
469 mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
472 string const solidName = subRegion.template getReference< string >(
"porousMaterialNames" );
473 constitutive::CoupledSolidBase & solid = PhysicsSolverBase::getConstitutiveModel< constitutive::CoupledSolidBase >(
474 subRegion, solidName );
475 auto & porosityModel = dynamic_cast< constitutive::BiotPorosity const & >( solid.getBasePorosityModel());
476 arrayView1d< real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
477 for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
479 porosityModel.updateAverageMeanTotalStressIncrement( k, averageMeanTotalStressIncrement[i] );
486 real64 computeAitkenRelaxationFactor( array1d< real64 >
const & s0,
487 array1d< real64 >
const & s1,
488 array1d< real64 >
const & s1_tilde,
489 array1d< real64 >
const & s2_tilde,
492 array1d< real64 > r1 = axpy( s1_tilde, s0, -1.0 );
493 array1d< real64 > r2 = axpy( s2_tilde, s1, -1.0 );
496 array1d< real64 > diff = axpy( r2, r1, -1.0 );
498 real64 const denom = dot( diff, diff );
499 real64 const numer = dot( r1, diff );
502 if( !isZero( denom ))
504 omega1 = -1.0 * omega0 * numer / denom;
509 array1d< real64 > computeUpdate( array1d< real64 >
const & s1,
510 array1d< real64 >
const & s2_tilde,
513 return axpy( scale( s1, 1.0 - omega1 ),
514 scale( s2_tilde, omega1 ),
518 void startSequentialIteration(
integer const & iter,
519 DomainPartition & domain )
override
521 if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
525 recordAverageMeanTotalStressIncrement( domain, m_s1 );
531 m_s1_tilde = m_s2_tilde;
537 void finishSequentialIteration(
integer const & iter,
538 DomainPartition & domain )
override
540 if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
549 m_omega1 = computeAitkenRelaxationFactor( m_s0, m_s1, m_s1_tilde, m_s2_tilde, m_omega0 );
550 m_s2 = computeUpdate( m_s1, m_s2_tilde, m_omega1 );
551 applyAcceleratedAverageMeanTotalStressIncrement( domain, m_s2 );
561 if( solverType ==
static_cast< integer >( SolverType::Flow ) )
563 this->
template forDiscretizationOnMeshTargets( domain.
getMeshBodies(), [&](
string const &,
568 mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
574 updateBulkDensity( subRegion );
580 if( solverType ==
static_cast< integer >( SolverType::SolidMechanics )
581 && !m_performStressInitialization )
584 averageMeanTotalStressIncrement( domain );
586 this->
template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&](
string const &,
591 mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
596 flowSolver()->updatePorosityAndPermeability( subRegion );
598 updateBulkDensity( subRegion );
604 if( solverType ==
static_cast< integer >( SolverType::SolidMechanics ) &&
605 this->getNonlinearSolverParameters().m_nonlinearAccelerationType== NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
607 recordAverageMeanTotalStressIncrement( domain, m_s2_tilde );
617 this->
template forDiscretizationOnMeshTargets( domain.
getMeshBodies(), [&](
string const &,
621 mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
625 string const solidName = subRegion.template getReference< string >( viewKeyStruct::porousMaterialNamesString() );
626 constitutive::CoupledSolidBase & solid = this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
628 arrayView2d< real64 const > const meanTotalStressIncrement_k = solid.getMeanTotalStressIncrement_k();
629 arrayView1d< real64 > const averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
631 finiteElement::FiniteElementBase & subRegionFE =
632 subRegion.template getReference< finiteElement::FiniteElementBase >( solidMechanicsSolver()->getDiscretizationName() );
635 finiteElement::FiniteElementDispatchHandler< BASE_FE_TYPES >::
636 dispatch3D( subRegionFE, [&] ( auto const finiteElement )
638 using FE_TYPE = decltype( finiteElement );
641 AverageOverQuadraturePoints1DKernelFactory::
642 createAndLaunch< CellElementSubRegion,
644 parallelDevicePolicy<> >( mesh.getNodeManager(),
645 mesh.getEdgeManager(),
646 mesh.getFaceManager(),
649 meanTotalStressIncrement_k,
650 averageMeanTotalStressIncrement_k );
656 virtual void updateBulkDensity( ElementSubRegionBase & subRegion ) = 0;
658 virtual void validateNonlinearAcceleration()
override
662 GEOS_ERROR(
"Nonlinear acceleration is not implemented for MPI runs" );
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
#define GEOS_ERROR_IF(EXP, msg)
Conditionally raise a hard error and terminate the program.
#define GEOS_THROW_IF(EXP, msg, TYPE)
Conditionally throw an exception.
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
std::tuple< SOLVERS *... > m_solvers
Pointers of the single-physics solvers.
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
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.
The ElementRegionManager class provides an interface to ObjectManagerBase in order to manage ElementR...
void forElementSubRegions(LAMBDA &&lambda)
This function is used to launch kernel function over the element subregions of all the subregion type...
Class facilitating the representation of a multi-level discretization of a MeshBody.
ElementRegionManager const & getElemManager() const
Get the element region manager.
@ Sequential
Sequential coupling.
virtual void registerDataOnMesh(Group &MeshBodies) override
Register wrappers that contain data on the mesh objects.
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.
array1d< string > m_stabilizationRegionNames
Names of regions where stabilization applied.
array1d< real64 > m_s0
Member variables needed for Nonlinear Acceleration ( Aitken ). Naming convention follows ( Jiang & Tc...
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override
function to perform setup for implicit timestep
real64 m_stabilizationMultiplier
Stabilization Multiplier.
virtual void initializePreSubGroups() override
Called by Initialize() prior to initializing sub-Groups.
stabilization::StabilizationType m_stabilizationType
Type of stabilization used.
FLOW_SOLVER * flowSolver() const
accessor for the pointer to the flow solver
PoromechanicsSolver(const string &name, dataRepository::Group *const parent)
main constructor for CoupledSolver Objects
virtual void initializePostInitialConditionsPreSubGroups() override
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
virtual void setConstitutiveNamesCallSuper(ElementSubRegionBase &subRegion) const override final
This function sets constitutive name fields on an ElementSubRegionBase, and calls the base function i...
MECHANICS_SOLVER * solidMechanicsSolver() const
accessor for the pointer to the solid mechanics solver
virtual void mapSolutionBetweenSolvers(DomainPartition &domain, integer const solverType) override
Maps the solution obtained from one solver to the fields used by the other solver(s)
integer m_isThermal
Flag to determine whether or not this is a thermal simulation.
integer m_performStressInitialization
Flag to indicate that the solver is going to perform stress initialization.
virtual void setupDofs(DomainPartition const &domain, DofManager &dofManager) const override
Populate degree-of-freedom manager with fields relevant to this solver.
virtual void registerDataOnMesh(dataRepository::Group &meshBodies) override
Register wrappers that contain data on the mesh objects.
void averageMeanTotalStressIncrement(DomainPartition &domain)
Helper function to average the mean total stress increment over quadrature points.
static string coupledSolverAttributePrefix()
String used to form the solverName used to register solvers in CoupledSolver.
virtual void initializePreSubGroups()
Called by Initialize() prior to initializing sub-Groups.
Wrapper< TBASE > & registerWrapper(string const &name, wrapperMap::KeyIndex::index_type *const rkey=nullptr)
Create and register a Wrapper around a new object.
DataContext const & getDataContext() const
void setRestartFlags(RestartFlags flags)
Set flags that control restart output of this group.
string const & getName() const
Get group name.
Group & setSizedFromParent(int val)
Set whether this wrapper is resized when its parent is resized.
DataContext const & getWrapperDataContext(KEY key) const
virtual void initializePostInitialConditionsPreSubGroups()
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
@ NOPLOT
Do not ever write to plot file.
@ OPTIONAL
Optional in input.
@ FALSE
Not read from input.
@ NO_WRITE
Do not write into restart.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
int MPI_COMM_GEOS
Global MPI communicator used by GEOSX.
array1d< globalIndex > globalIndex_array
A 1-dimensional array of geos::globalIndex types.
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "preconditioner")
Declare strings associated with enumeration values.
std::set< T > set
A set of local indices.
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.
@ None
No Schur complement - just block-GS/block-Jacobi preconditioner.
Array< T, 1 > array1d
Alias for 1D array.
constexpr static char const * porousMaterialNamesString()
Names of the porous materials.
constexpr static const char * stabilizationRegionNamesString()
Name of regions on which stabilization applied.
constexpr static char const * performStressInitializationString()
Flag to indicate that the solver is going to perform stress initialization.
static constexpr char const * hydraulicApertureRelationNameString()
Name of the hydraulicApertureRelationName.
constexpr static char const * isThermalString()
Flag to indicate that the simulation is thermal.
constexpr static char const * stabilizationTypeString()
Type of pressure stabilization.
constexpr static const char * stabilizationMultiplierString()
Multiplier on stabilization strength.
Structure to hold scoped key names.