GEOS
PoromechanicsSolver.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 TotalEnergies
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
21 #ifndef GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSSOLVER_HPP_
22 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSSOLVER_HPP_
23 
28 #include "constitutive/solid/CoupledSolidBase.hpp"
29 #include "constitutive/solid/PorousSolid.hpp"
30 #include "constitutive/solid/PorousDamageSolid.hpp"
31 #include "constitutive/contact/HydraulicApertureBase.hpp"
32 #include "mesh/DomainPartition.hpp"
34 #include "codingUtilities/Utilities.hpp"
35 
36 namespace geos
37 {
38 
39 namespace stabilization
40 {
41 enum class StabilizationType : integer
42 {
44  Global,
45  Local,
46 };
47 
48 ENUM_STRINGS( StabilizationType,
49  "None",
50  "Global",
51  "Local" );
52 }
53 
54 
55 template< typename FLOW_SOLVER, typename MECHANICS_SOLVER = SolidMechanicsLagrangianFEM >
56 class PoromechanicsSolver : public CoupledSolver< FLOW_SOLVER, MECHANICS_SOLVER >
57 {
58 public:
59 
61  using Base::m_solvers;
62  using Base::m_dofManager;
63  using Base::m_localMatrix;
64  using Base::m_rhs;
65  using Base::m_solution;
66 
67  enum class SolverType : integer
68  {
69  Flow = 0,
70  SolidMechanics = 1
71  };
72 
74  static string coupledSolverAttributePrefix() { return "poromechanics"; }
75 
81  PoromechanicsSolver( const string & name,
82  dataRepository::Group * const parent )
83  : Base( name, parent ),
84  m_isThermal( 0 ),
86  {
88  setApplyDefaultValue( 0 ).
90  setDescription( "Flag indicating whether the problem is thermal or not. Set isThermal=\"1\" to enable the thermal coupling" );
91 
94  setDescription( "StabilizationType. Options are:\n" +
95  toString( stabilization::StabilizationType::None ) + "- Add no stabilization to mass equation \n" +
96  toString( stabilization::StabilizationType::Global ) + "- Add jump stabilization to all faces \n" +
97  toString( stabilization::StabilizationType::Local ) + "- Add jump stabilization on interior of macro elements" );
98 
100  setRTTypeName( rtTypes::CustomTypes::groupNameRefArray ).
101  setInputFlag( dataRepository::InputFlags::OPTIONAL ).
102  setDescription( "Regions where stabilization is applied." );
103 
105  setApplyDefaultValue( 1.0 ).
106  setInputFlag( dataRepository::InputFlags::OPTIONAL ).
107  setDescription( "Constant multiplier of stabilization strength" );
108  }
109 
111  {
113 
114  GEOS_THROW_IF( this->m_isThermal && !this->flowSolver()->isThermal(),
115  GEOS_FMT( "{} {}: The attribute `{}` of the flow solver `{}` must be set to 1 since the poromechanics solver is thermal",
116  this->getCatalogName(), this->getName(), FlowSolverBase::viewKeyStruct::isThermalString(), this->flowSolver()->getName() ),
117  InputError );
118  }
119 
120  virtual void setConstitutiveNamesCallSuper( ElementSubRegionBase & subRegion ) const override final
121  {
122  if( dynamic_cast< SurfaceElementSubRegion * >( &subRegion ) )
123  {
124  subRegion.registerWrapper< string >( viewKeyStruct::hydraulicApertureRelationNameString() ).
125  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
127  setSizedFromParent( 0 );
128 
129  string & hydraulicApertureModelName = subRegion.getReference< string >( viewKeyStruct::hydraulicApertureRelationNameString() );
130  hydraulicApertureModelName = PhysicsSolverBase::getConstitutiveName< constitutive::HydraulicApertureBase >( subRegion );
131  GEOS_ERROR_IF( hydraulicApertureModelName.empty(), GEOS_FMT( "{}: HydraulicApertureBase model not found on subregion {}",
132  this->getDataContext(), subRegion.getDataContext() ) );
133  }
134 
135  }
136 
137  virtual void initializePreSubGroups() override
138  {
140 
141  GEOS_THROW_IF( m_stabilizationType == stabilization::StabilizationType::Local,
143  ": Local stabilization has been temporarily disabled",
144  InputError );
145 
146  DomainPartition & domain = this->template getGroupByPath< DomainPartition >( "/Problem/domain" );
147 
148  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&] ( string const &,
149  MeshLevel & mesh,
150  string_array const & regionNames )
151  {
152  ElementRegionManager & elementRegionManager = mesh.getElemManager();
153  elementRegionManager.forElementSubRegions< CellElementSubRegion >( regionNames,
154  [&]( localIndex const,
155  ElementSubRegionBase & subRegion )
156  {
157  string & porousName = subRegion.getReference< string >( viewKeyStruct::porousMaterialNamesString() );
158  porousName = this->template getConstitutiveName< constitutive::CoupledSolidBase >( subRegion );
159  GEOS_THROW_IF( porousName.empty(),
160  GEOS_FMT( "{} {} : Solid model not found on subregion {}",
161  this->getCatalogName(), this->getDataContext().toString(), subRegion.getName() ),
162  InputError );
163 
164  string & porosityModelName = subRegion.getReference< string >( constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString() );
165  porosityModelName = this->template getConstitutiveName< constitutive::PorosityBase >( subRegion );
166  GEOS_THROW_IF( porosityModelName.empty(),
167  GEOS_FMT( "{} {} : Porosity model not found on subregion {}",
168  this->getCatalogName(), this->getDataContext().toString(), subRegion.getName() ),
169  InputError );
170 
171  if( subRegion.hasField< fields::poromechanics::bulkDensity >() )
172  {
173  // get the solid model to know the number of quadrature points and resize the bulk density
174  constitutive::CoupledSolidBase const & solid = this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, porousName );
175  subRegion.getField< fields::poromechanics::bulkDensity >().resizeDimension< 1 >( solid.getDensity().size( 1 ) );
176  }
177  } );
178  } );
179  }
180 
181  virtual void registerDataOnMesh( dataRepository::Group & meshBodies ) override
182  {
184 
185  if( this->getNonlinearSolverParameters().m_couplingType == NonlinearSolverParameters::CouplingType::Sequential )
186  {
187  // to let the solid mechanics solver that there is a pressure and temperature RHS in the mechanics solve
188  solidMechanicsSolver()->enableFixedStressPoromechanicsUpdate();
189  // to let the flow solver that saving pressure_k and temperature_k is necessary (for the fixed-stress porosity terms)
190  flowSolver()->enableFixedStressPoromechanicsUpdate();
191  }
192 
193  if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
194  {
195  flowSolver()->enableJumpStabilization();
196  }
197 
198  PhysicsSolverBase::forDiscretizationOnMeshTargets( meshBodies, [&] ( string const &,
199  MeshLevel & mesh,
200  string_array const & regionNames )
201  {
202  ElementRegionManager & elemManager = mesh.getElemManager();
203 
204  elemManager.forElementSubRegions< CellElementSubRegion >( regionNames,
205  [&]( localIndex const,
206  ElementSubRegionBase & subRegion )
207  {
208  subRegion.registerWrapper< string >( viewKeyStruct::porousMaterialNamesString() ).
209  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
210  setRestartFlags( dataRepository::RestartFlags::NO_WRITE ).
211  setSizedFromParent( 0 );
212 
213  // This is needed by the way the surface generator currently does things.
214  subRegion.registerWrapper< string >( constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString() ).
215  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
216  setRestartFlags( dataRepository::RestartFlags::NO_WRITE ).
217  setSizedFromParent( 0 );
218 
219  // register the bulk density for use in the solid mechanics solver
220  // ideally we would resize it here as well, but the solid model name is not available yet (see below)
221  subRegion.registerField< fields::poromechanics::bulkDensity >( this->getName() );
222 
223  if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
224  {
225  subRegion.registerField< fields::flow::macroElementIndex >( this->getName());
226  subRegion.registerField< fields::flow::elementStabConstant >( this->getName());
227  }
228  } );
229  } );
230  }
231 
232  virtual void implicitStepSetup( real64 const & time_n,
233  real64 const & dt,
234  DomainPartition & domain ) override
235  {
236  flowSolver()->setKeepVariablesConstantDuringInitStep( m_performStressInitialization );
237 
238  if( this->m_stabilizationType == stabilization::StabilizationType::Global || this->m_stabilizationType == stabilization::StabilizationType::Local )
239  {
240  this->updateStabilizationParameters( domain );
241  }
242 
243  Base::implicitStepSetup( time_n, dt, domain );
244  }
245 
246  virtual void setupDofs( DomainPartition const & domain,
247  DofManager & dofManager ) const override
248  {
249  // note that the order of operations matters a lot here (for instance for the MGR labels)
250  // we must set up dofs for solid mechanics first, and then for flow
251  // that's the reason why this function is here and not in CoupledSolvers.hpp
252  solidMechanicsSolver()->setupDofs( domain, dofManager );
253  flowSolver()->setupDofs( domain, dofManager );
254  this->setupCoupling( domain, dofManager );
255  }
256 
257  virtual bool checkSequentialConvergence( int const & iter,
258  real64 const & time_n,
259  real64 const & dt,
260  DomainPartition & domain ) override
261  {
262  // always force outer loop for initialization
263  auto & subcycling = this->getNonlinearSolverParameters().m_subcyclingOption;
264  auto const subcycling_orig = subcycling;
265  if( m_performStressInitialization )
266  subcycling = 1;
267 
268  bool isConverged = Base::checkSequentialConvergence( iter, time_n, dt, domain );
269 
270  // restore original
271  subcycling = subcycling_orig;
272 
273  return isConverged;
274  }
275 
280  MECHANICS_SOLVER * solidMechanicsSolver() const
281  {
282  return std::get< toUnderlying( SolverType::SolidMechanics ) >( m_solvers );
283  }
284 
289  FLOW_SOLVER * flowSolver() const
290  {
291  return std::get< toUnderlying( SolverType::Flow ) >( m_solvers );
292  }
293 
294  /*
295  * @brief Utility function to set the stress initialization flag
296  * @param[in] performStressInitialization true if the solver has to initialize stress, false otherwise
297  */
298  void setStressInitialization( bool const performStressInitialization )
299  {
300  m_performStressInitialization = performStressInitialization;
301  solidMechanicsSolver()->setStressInitialization( performStressInitialization );
302  }
303 
305  {
307  constexpr static char const * porousMaterialNamesString() { return "porousMaterialNames"; }
308 
310  constexpr static char const * isThermalString() { return "isThermal"; }
311 
313  constexpr static char const * stabilizationTypeString() {return "stabilizationType"; }
314 
316  constexpr static const char * stabilizationRegionNamesString() {return "stabilizationRegionNames"; }
317 
319  constexpr static const char * stabilizationMultiplierString() {return "stabilizationMultiplier"; }
320 
322  static constexpr char const * hydraulicApertureRelationNameString() {return "hydraulicApertureRelationName"; }
323 
324  };
325 
326  void updateStabilizationParameters( DomainPartition & domain ) const
327  {
328  // Step 1: we loop over the regions where stabilization is active and collect their name
329  set< string > regionFilter;
330  for( string const & regionName : m_stabilizationRegionNames )
331  {
332  regionFilter.insert( regionName );
333  }
334 
335  // Step 2: loop over target regions of solver, and tag the elements belonging to the stabilization regions
336  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&] ( string const &,
337  MeshLevel & mesh,
338  string_array const & targetRegionNames )
339  {
340  //keep only target regions in filter
341  string_array filteredTargetRegionNames;
342  filteredTargetRegionNames.reserve( targetRegionNames.size() );
343 
344  for( string const & targetRegionName : targetRegionNames )
345  {
346  if( regionFilter.count( targetRegionName ) )
347  {
348  filteredTargetRegionNames.emplace_back( targetRegionName );
349  }
350  }
351 
352  // Loop over elements and update stabilization constant
353  mesh.getElemManager().forElementSubRegions( filteredTargetRegionNames, [&]( localIndex const,
354  ElementSubRegionBase & subRegion )
355  {
356  arrayView1d< integer > const macroElementIndex = subRegion.getField< fields::flow::macroElementIndex >();
357  arrayView1d< real64 > const elementStabConstant = subRegion.getField< fields::flow::elementStabConstant >();
358 
359  geos::constitutive::CoupledSolidBase const & porousSolid =
360  this->template getConstitutiveModel< geos::constitutive::CoupledSolidBase >( subRegion, subRegion.getReference< string >( viewKeyStruct::porousMaterialNamesString() ) );
361 
362  arrayView1d< real64 const > const bulkModulus = porousSolid.getBulkModulus();
363  arrayView1d< real64 const > const shearModulus = porousSolid.getShearModulus();
364  arrayView1d< real64 const > const biotCoefficient = porousSolid.getBiotCoefficient();
365 
366  real64 const stabilizationMultiplier = m_stabilizationMultiplier;
367 
368  forAll< parallelDevicePolicy<> >( subRegion.size(), [bulkModulus,
369  shearModulus,
370  biotCoefficient,
371  stabilizationMultiplier,
372  macroElementIndex,
373  elementStabConstant] GEOS_HOST_DEVICE ( localIndex const ei )
374 
375  {
376  real64 const bM = bulkModulus[ei];
377  real64 const sM = shearModulus[ei];
378  real64 const bC = biotCoefficient[ei];
379 
380  macroElementIndex[ei] = 1;
381  elementStabConstant[ei] = stabilizationMultiplier * 9.0 * (bC * bC) / (32.0 * (10.0 * sM / 3.0 + bM));
382  } );
383  } );
384  } );
385 
386  }
387 
388  void updateBulkDensity( DomainPartition & domain )
389  {
390  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
391  MeshLevel & mesh,
392  string_array const & regionNames )
393  {
394  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
395  auto & subRegion )
396  {
397  // update the bulk density
398  // TODO: ideally, we would not recompute the bulk density, but a more general "rhs" containing the body force and the
399  // pressure/temperature terms
400  updateBulkDensity( subRegion );
401  } );
402  } );
403  }
404 
405 protected:
406 
407  template< typename CONSTITUTIVE_BASE,
408  typename KERNEL_WRAPPER,
409  typename ... PARAMS >
410  real64 assemblyLaunch( MeshLevel & mesh,
411  DofManager const & dofManager,
412  string_array const & regionNames,
413  string const & materialNamesString,
414  CRSMatrixView< real64, globalIndex const > const & localMatrix,
415  arrayView1d< real64 > const & localRhs,
416  real64 const dt,
417  PARAMS && ... params )
418  {
420 
421  NodeManager const & nodeManager = mesh.getNodeManager();
422 
423  string const dofKey = dofManager.getKey( fields::solidMechanics::totalDisplacement::key() );
424  arrayView1d< globalIndex const > const & dofNumber = nodeManager.getReference< globalIndex_array >( dofKey );
425 
426  real64 const gravityVectorData[3] = LVARRAY_TENSOROPS_INIT_LOCAL_3( PhysicsSolverBase::gravityVector() );
427 
428  KERNEL_WRAPPER kernelWrapper( dofNumber,
429  dofManager.rankOffset(),
430  localMatrix,
431  localRhs,
432  dt,
433  gravityVectorData,
434  std::forward< PARAMS >( params )... );
435 
436  return finiteElement::
437  regionBasedKernelApplication< parallelDevicePolicy< >,
438  CONSTITUTIVE_BASE,
439  CellElementSubRegion >( mesh,
440  regionNames,
441  this->solidMechanicsSolver()->getDiscretizationName(),
442  materialNamesString,
443  kernelWrapper );
444  }
445 
446  /* Implementation of Nonlinear Acceleration (Aitken) of averageMeanTotalStressIncrement */
447 
448  void recordAverageMeanTotalStressIncrement( DomainPartition & domain,
449  array1d< real64 > & averageMeanTotalStressIncrement )
450  {
451  averageMeanTotalStressIncrement.resize( 0 );
452  PhysicsSolverBase::forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
453  MeshLevel & mesh,
454  string_array const & regionNames ) {
455  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
456  auto & subRegion ) {
457  // get the solid model (to access stress increment)
458  string const solidName = subRegion.template getReference< string >( "porousMaterialNames" );
459  constitutive::CoupledSolidBase & solid = PhysicsSolverBase::getConstitutiveModel< constitutive::CoupledSolidBase >(
460  subRegion, solidName );
461 
462  arrayView1d< const real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
463  for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
464  {
465  averageMeanTotalStressIncrement.emplace_back( averageMeanTotalStressIncrement_k[k] );
466  }
467  } );
468  } );
469  }
470 
471  void applyAcceleratedAverageMeanTotalStressIncrement( DomainPartition & domain,
472  array1d< real64 > & averageMeanTotalStressIncrement )
473  {
474  integer i = 0;
475  PhysicsSolverBase::forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
476  MeshLevel & mesh,
477  string_array const & regionNames ) {
478  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
479  auto & subRegion ) {
480  // get the solid model (to access stress increment)
481  string const solidName = subRegion.template getReference< string >( "porousMaterialNames" );
482  constitutive::CoupledSolidBase & solid = PhysicsSolverBase::getConstitutiveModel< constitutive::CoupledSolidBase >(
483  subRegion, solidName );
484  auto & porosityModel = dynamic_cast< constitutive::BiotPorosity const & >( solid.getBasePorosityModel());
485  arrayView1d< real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
486  for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
487  {
488  porosityModel.updateAverageMeanTotalStressIncrement( k, averageMeanTotalStressIncrement[i] );
489  i++;
490  }
491  } );
492  } );
493  }
494 
495  real64 computeAitkenRelaxationFactor( array1d< real64 > const & s0,
496  array1d< real64 > const & s1,
497  array1d< real64 > const & s1_tilde,
498  array1d< real64 > const & s2_tilde,
499  real64 const omega0 )
500  {
501  array1d< real64 > r1 = axpy( s1_tilde, s0, -1.0 );
502  array1d< real64 > r2 = axpy( s2_tilde, s1, -1.0 );
503 
504  // diff = r2 - r1
505  array1d< real64 > diff = axpy( r2, r1, -1.0 );
506 
507  real64 const denom = dot( diff, diff );
508  real64 const numer = dot( r1, diff );
509 
510  real64 omega1 = 1.0;
511  if( !isZero( denom ))
512  {
513  omega1 = -1.0 * omega0 * numer / denom;
514  }
515  return omega1;
516  }
517 
518  array1d< real64 > computeUpdate( array1d< real64 > const & s1,
519  array1d< real64 > const & s2_tilde,
520  real64 const omega1 )
521  {
522  return axpy( scale( s1, 1.0 - omega1 ),
523  scale( s2_tilde, omega1 ),
524  1.0 );
525  }
526 
527  void startSequentialIteration( integer const & iter,
528  DomainPartition & domain ) override
529  {
530  if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
531  {
532  if( iter == 0 )
533  {
534  recordAverageMeanTotalStressIncrement( domain, m_s1 );
535  }
536  else
537  {
538  m_s0 = m_s1;
539  m_s1 = m_s2;
540  m_s1_tilde = m_s2_tilde;
541  m_omega0 = m_omega1;
542  }
543  }
544  }
545 
546  void finishSequentialIteration( integer const & iter,
547  DomainPartition & domain ) override
548  {
549  if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
550  {
551  if( iter == 0 )
552  {
553  m_s2 = m_s2_tilde;
554  m_omega1 = 1.0;
555  }
556  else
557  {
558  m_omega1 = computeAitkenRelaxationFactor( m_s0, m_s1, m_s1_tilde, m_s2_tilde, m_omega0 );
559  m_s2 = computeUpdate( m_s1, m_s2_tilde, m_omega1 );
560  applyAcceleratedAverageMeanTotalStressIncrement( domain, m_s2 );
561  }
562  }
563  }
564 
565  virtual void mapSolutionBetweenSolvers( DomainPartition & domain, integer const solverType ) override
566  {
568 
570  if( solverType == static_cast< integer >( SolverType::Flow ) )
571  {
572  updateBulkDensity( domain );
573  }
574 
576  if( solverType == static_cast< integer >( SolverType::SolidMechanics )
577  && !m_performStressInitialization ) // do not update during poromechanics initialization
578  {
579  // compute the average of the mean total stress increment over quadrature points
580  averageMeanTotalStressIncrement( domain );
581 
582  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
583  MeshLevel & mesh,
584  string_array const & regionNames )
585  {
586 
587  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
588  auto & subRegion )
589  {
590  // update the porosity after a change in displacement (after mechanics solve)
591  // or a change in pressure/temperature (after a flow solve)
592  flowSolver()->updatePorosityAndPermeability( subRegion );
593  // update bulk density to reflect porosity change into mechanics
594  updateBulkDensity( subRegion );
595  } );
596  } );
597  }
598 
599  // needed to perform nonlinear acceleration
600  if( solverType == static_cast< integer >( SolverType::SolidMechanics ) &&
601  this->getNonlinearSolverParameters().m_nonlinearAccelerationType== NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
602  {
603  recordAverageMeanTotalStressIncrement( domain, m_s2_tilde );
604  }
605  }
606 
612  {
613  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
614  MeshLevel & mesh,
615  string_array const & regionNames )
616  {
617  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
618  auto & subRegion )
619  {
620  // get the solid model (to access stress increment)
621  string const solidName = subRegion.template getReference< string >( viewKeyStruct::porousMaterialNamesString() );
622  constitutive::CoupledSolidBase & solid = this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
623 
624  arrayView2d< real64 const > const meanTotalStressIncrement_k = solid.getMeanTotalStressIncrement_k();
625  arrayView1d< real64 > const averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
626 
627  finiteElement::FiniteElementBase & subRegionFE =
628  subRegion.template getReference< finiteElement::FiniteElementBase >( solidMechanicsSolver()->getDiscretizationName() );
629 
630  // determine the finite element type
631  finiteElement::FiniteElementDispatchHandler< BASE_FE_TYPES >::
632  dispatch3D( subRegionFE, [&] ( auto const finiteElement )
633  {
634  using FE_TYPE = decltype( finiteElement );
635 
636  // call the factory and launch the kernel
637  AverageOverQuadraturePoints1DKernelFactory::
638  createAndLaunch< CellElementSubRegion,
639  FE_TYPE,
640  parallelDevicePolicy<> >( mesh.getNodeManager(),
641  mesh.getEdgeManager(),
642  mesh.getFaceManager(),
643  subRegion,
644  finiteElement,
645  meanTotalStressIncrement_k,
646  averageMeanTotalStressIncrement_k );
647  } );
648  } );
649  } );
650  }
651 
652  virtual void updateBulkDensity( ElementSubRegionBase & subRegion ) = 0;
653 
654  virtual void validateNonlinearAcceleration() override
655  {
656  if( MpiWrapper::commSize( MPI_COMM_GEOS ) > 1 )
657  {
658  GEOS_ERROR( "Nonlinear acceleration is not implemented for MPI runs" );
659  }
660  }
661 
664 
667 
669  stabilization::StabilizationType m_stabilizationType;
670 
673 
676 
678  array1d< real64 > m_s0; // Accelerated averageMeanTotalStresIncrement @ outer iteration v ( two iterations ago )
679  array1d< real64 > m_s1; // Accelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( previous iteration )
680  array1d< real64 > m_s1_tilde; // Unaccelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( previous iteration )
681  array1d< real64 > m_s2; // Accelerated averageMeanTotalStresIncrement @ outer iteration v + 2 ( current iteration )
682  array1d< real64 > m_s2_tilde; // Unaccelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( current iteration )
683  real64 m_omega0; // Old Aitken relaxation factor
684  real64 m_omega1; // New Aitken relaxation factor
685 
686 };
687 
688 } /* namespace geos */
689 
690 #endif //GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSSOLVER_HPP_
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
Definition: Logger.hpp:157
#define GEOS_ERROR_IF(EXP, msg)
Conditionally raise a hard error and terminate the program.
Definition: Logger.hpp:142
#define GEOS_THROW_IF(EXP, msg, TYPE)
Conditionally throw an exception.
Definition: Logger.hpp:151
#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,...
Definition: DofManager.hpp:44
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.
Definition: MeshLevel.hpp:42
ElementRegionManager const & getElemManager() const
Get the element region manager.
Definition: MeshLevel.hpp:207
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< 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.
bool m_performStressInitialization
Flag to indicate that the solver is going to perform stress initialization.
string_array m_stabilizationRegionNames
Names of regions where stabilization applied.
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.
Definition: Group.hpp:1542
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
Definition: Group.hpp:1345
void setRestartFlags(RestartFlags flags)
Set flags that control restart output of this group.
Definition: Group.hpp:1426
string const & getName() const
Get group name.
Definition: Group.hpp:1331
Group & setSizedFromParent(int val)
Set whether this wrapper is resized when its parent is resized.
Definition: Group.hpp:1411
DataContext const & getWrapperDataContext(KEY key) const
Definition: Group.hpp:1356
virtual void initializePostInitialConditionsPreSubGroups()
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
Definition: Group.hpp:1554
@ NOPLOT
Do not ever write to plot file.
@ OPTIONAL
Optional in input.
@ NO_WRITE
Do not write into restart.
int MPI_COMM_GEOS
Global MPI communicator used by GEOSX.
array1d< globalIndex > globalIndex_array
A 1-dimensional array of geos::globalIndex types.
Definition: DataTypes.hpp:402
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.
Definition: DataTypes.hpp:263
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:85
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
@ None
No Schur complement - just block-GS/block-Jacobi preconditioner.
std::vector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:393
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:176
Exception class used to report errors in user input.
Definition: Logger.hpp:464
constexpr static char const * porousMaterialNamesString()
Names of the porous materials.
constexpr static const char * stabilizationRegionNamesString()
Name of regions on which stabilization applied.
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.
Definition: Group.hpp:1444