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/contact/HydraulicApertureBase.hpp"
31 #include "mesh/DomainPartition.hpp"
33 #include "codingUtilities/Utilities.hpp"
34 
35 namespace geos
36 {
37 
38 namespace stabilization
39 {
40 enum class StabilizationType : integer
41 {
43  Global,
44  Local,
45 };
46 
47 ENUM_STRINGS( StabilizationType,
48  "None",
49  "Global",
50  "Local" );
51 }
52 
53 
54 template< typename FLOW_SOLVER, typename MECHANICS_SOLVER = SolidMechanicsLagrangianFEM >
55 class PoromechanicsSolver : public CoupledSolver< FLOW_SOLVER, MECHANICS_SOLVER >
56 {
57 public:
58 
60  using Base::m_solvers;
61  using Base::m_dofManager;
62  using Base::m_localMatrix;
63  using Base::m_rhs;
64  using Base::m_solution;
65 
66  enum class SolverType : integer
67  {
68  Flow = 0,
69  SolidMechanics = 1
70  };
71 
73  static string coupledSolverAttributePrefix() { return "poromechanics"; }
74 
80  PoromechanicsSolver( const string & name,
81  dataRepository::Group * const parent )
82  : Base( name, parent ),
83  m_isThermal( 0 )
84  {
86  setApplyDefaultValue( 0 ).
88  setDescription( "Flag indicating whether the problem is thermal or not. Set isThermal=\"1\" to enable the thermal coupling" );
89 
91  setApplyDefaultValue( false ).
92  setInputFlag( dataRepository::InputFlags::FALSE ).
93  setDescription( "Flag to indicate that the solver is going to perform stress initialization" );
94 
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" );
101 
103  setRTTypeName( rtTypes::CustomTypes::groupNameRefArray ).
104  setInputFlag( dataRepository::InputFlags::OPTIONAL ).
105  setDescription( "Regions where stabilization is applied." );
106 
108  setApplyDefaultValue( 1.0 ).
109  setInputFlag( dataRepository::InputFlags::OPTIONAL ).
110  setDescription( "Constant multiplier of stabilization strength" );
111  }
112 
114  {
116 
117  GEOS_THROW_IF( this->m_isThermal && !this->flowSolver()->isThermal(),
118  GEOS_FMT( "{} {}: The attribute `{}` of the flow solver `{}` must be set to 1 since the poromechanics solver is thermal",
119  this->getCatalogName(), this->getName(), FlowSolverBase::viewKeyStruct::isThermalString(), this->flowSolver()->getName() ),
120  InputError );
121  }
122 
123  virtual void setConstitutiveNamesCallSuper( ElementSubRegionBase & subRegion ) const override final
124  {
125  if( dynamic_cast< SurfaceElementSubRegion * >( &subRegion ) )
126  {
127  subRegion.registerWrapper< string >( viewKeyStruct::hydraulicApertureRelationNameString() ).
128  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
130  setSizedFromParent( 0 );
131 
132  string & hydraulicApertureModelName = subRegion.getReference< string >( viewKeyStruct::hydraulicApertureRelationNameString() );
133  hydraulicApertureModelName = PhysicsSolverBase::getConstitutiveName< constitutive::HydraulicApertureBase >( subRegion );
134  GEOS_ERROR_IF( hydraulicApertureModelName.empty(), GEOS_FMT( "{}: HydraulicApertureBase model not found on subregion {}",
135  this->getDataContext(), subRegion.getDataContext() ) );
136  }
137 
138  }
139 
140  virtual void initializePreSubGroups() override
141  {
143 
144  GEOS_THROW_IF( m_stabilizationType == stabilization::StabilizationType::Local,
146  ": Local stabilization has been temporarily disabled",
147  InputError );
148 
149  DomainPartition & domain = this->template getGroupByPath< DomainPartition >( "/Problem/domain" );
150 
151  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&] ( string const &,
152  MeshLevel & mesh,
153  arrayView1d< string const > const & regionNames )
154  {
155  ElementRegionManager & elementRegionManager = mesh.getElemManager();
156  elementRegionManager.forElementSubRegions< CellElementSubRegion >( regionNames,
157  [&]( localIndex const,
158  ElementSubRegionBase & subRegion )
159  {
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() ),
165  InputError );
166 
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() ),
172  InputError );
173 
174  if( subRegion.hasField< fields::poromechanics::bulkDensity >() )
175  {
176  // get the solid model to know the number of quadrature points and resize the bulk density
177  constitutive::CoupledSolidBase const & solid = this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, porousName );
178  subRegion.getField< fields::poromechanics::bulkDensity >().resizeDimension< 1 >( solid.getDensity().size( 1 ) );
179  }
180  } );
181  } );
182  }
183 
184  virtual void registerDataOnMesh( dataRepository::Group & meshBodies ) override
185  {
187 
188  if( this->getNonlinearSolverParameters().m_couplingType == NonlinearSolverParameters::CouplingType::Sequential )
189  {
190  // to let the solid mechanics solver that there is a pressure and temperature RHS in the mechanics solve
191  solidMechanicsSolver()->enableFixedStressPoromechanicsUpdate();
192  // to let the flow solver that saving pressure_k and temperature_k is necessary (for the fixed-stress porosity terms)
193  flowSolver()->enableFixedStressPoromechanicsUpdate();
194  }
195 
196  if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
197  {
198  flowSolver()->enableJumpStabilization();
199  }
200 
201  PhysicsSolverBase::forDiscretizationOnMeshTargets( meshBodies, [&] ( string const &,
202  MeshLevel & mesh,
203  arrayView1d< string const > const & regionNames )
204  {
205  ElementRegionManager & elemManager = mesh.getElemManager();
206 
207  elemManager.forElementSubRegions< CellElementSubRegion >( regionNames,
208  [&]( localIndex const,
209  ElementSubRegionBase & subRegion )
210  {
211  subRegion.registerWrapper< string >( viewKeyStruct::porousMaterialNamesString() ).
212  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
213  setRestartFlags( dataRepository::RestartFlags::NO_WRITE ).
214  setSizedFromParent( 0 );
215 
216  // This is needed by the way the surface generator currently does things.
217  subRegion.registerWrapper< string >( constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString() ).
218  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
219  setRestartFlags( dataRepository::RestartFlags::NO_WRITE ).
220  setSizedFromParent( 0 );
221 
222  if( this->getNonlinearSolverParameters().m_couplingType == NonlinearSolverParameters::CouplingType::Sequential )
223  {
224  // register the bulk density for use in the solid mechanics solver
225  // ideally we would resize it here as well, but the solid model name is not available yet (see below)
226  subRegion.registerField< fields::poromechanics::bulkDensity >( this->getName() );
227  }
228 
229  if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
230  {
231  subRegion.registerField< fields::flow::macroElementIndex >( this->getName());
232  subRegion.registerField< fields::flow::elementStabConstant >( this->getName());
233  }
234  } );
235  } );
236  }
237 
238  virtual void implicitStepSetup( real64 const & time_n,
239  real64 const & dt,
240  DomainPartition & domain ) override
241  {
242  flowSolver()->setKeepVariablesConstantDuringInitStep( m_performStressInitialization );
243 
244  if( this->m_stabilizationType == stabilization::StabilizationType::Global || this->m_stabilizationType == stabilization::StabilizationType::Local )
245  {
246  this->updateStabilizationParameters( domain );
247  }
248 
249  Base::implicitStepSetup( time_n, dt, domain );
250  }
251 
252  virtual void setupDofs( DomainPartition const & domain,
253  DofManager & dofManager ) const override
254  {
255  // note that the order of operations matters a lot here (for instance for the MGR labels)
256  // we must set up dofs for solid mechanics first, and then for flow
257  // that's the reason why this function is here and not in CoupledSolvers.hpp
258  solidMechanicsSolver()->setupDofs( domain, dofManager );
259  flowSolver()->setupDofs( domain, dofManager );
260  this->setupCoupling( domain, dofManager );
261  }
262 
263  virtual bool checkSequentialConvergence( int const & iter,
264  real64 const & time_n,
265  real64 const & dt,
266  DomainPartition & domain ) override
267  {
268  // always force outer loop for initialization
269  auto & subcycling = this->getNonlinearSolverParameters().m_subcyclingOption;
270  auto const subcycling_orig = subcycling;
271  if( m_performStressInitialization )
272  subcycling = 1;
273 
274  bool isConverged = Base::checkSequentialConvergence( iter, time_n, dt, domain );
275 
276  // restore original
277  subcycling = subcycling_orig;
278 
279  return isConverged;
280  }
281 
286  MECHANICS_SOLVER * solidMechanicsSolver() const
287  {
288  return std::get< toUnderlying( SolverType::SolidMechanics ) >( m_solvers );
289  }
290 
295  FLOW_SOLVER * flowSolver() const
296  {
297  return std::get< toUnderlying( SolverType::Flow ) >( m_solvers );
298  }
299 
300  /*
301  * @brief Utility function to set the stress initialization flag
302  * @param[in] performStressInitialization true if the solver has to initialize stress, false otherwise
303  */
304  void setStressInitialization( integer const performStressInitialization )
305  {
306  m_performStressInitialization = performStressInitialization;
307  }
308 
310  {
312  constexpr static char const * porousMaterialNamesString() { return "porousMaterialNames"; }
313 
315  constexpr static char const * isThermalString() { return "isThermal"; }
316 
318  constexpr static char const * performStressInitializationString() { return "performStressInitialization"; }
319 
321  constexpr static char const * stabilizationTypeString() {return "stabilizationType"; }
322 
324  constexpr static const char * stabilizationRegionNamesString() {return "stabilizationRegionNames"; }
325 
327  constexpr static const char * stabilizationMultiplierString() {return "stabilizationMultiplier"; }
328 
330  static constexpr char const * hydraulicApertureRelationNameString() {return "hydraulicApertureRelationName"; }
331 
332  };
333 
334  void updateStabilizationParameters( DomainPartition & domain ) const
335  {
336  // Step 1: we loop over the regions where stabilization is active and collect their name
337  set< string > regionFilter;
338  for( string const & regionName : m_stabilizationRegionNames )
339  {
340  regionFilter.insert( regionName );
341  }
342 
343  // Step 2: loop over target regions of solver, and tag the elements belonging to the stabilization regions
344  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&] ( string const &,
345  MeshLevel & mesh,
346  arrayView1d< string const > const & targetRegionNames )
347  {
348  //keep only target regions in filter
349  array1d< string > filteredTargetRegionNames;
350  filteredTargetRegionNames.reserve( targetRegionNames.size() );
351 
352  for( string const & targetRegionName : targetRegionNames )
353  {
354  if( regionFilter.count( targetRegionName ) )
355  {
356  filteredTargetRegionNames.emplace_back( targetRegionName );
357  }
358  }
359 
360  // Loop over elements and update stabilization constant
361  mesh.getElemManager().forElementSubRegions( filteredTargetRegionNames.toViewConst(), [&]( localIndex const,
362  ElementSubRegionBase & subRegion )
363  {
364  arrayView1d< integer > const macroElementIndex = subRegion.getField< fields::flow::macroElementIndex >();
365  arrayView1d< real64 > const elementStabConstant = subRegion.getField< fields::flow::elementStabConstant >();
366 
367  geos::constitutive::CoupledSolidBase const & porousSolid =
368  this->template getConstitutiveModel< geos::constitutive::CoupledSolidBase >( subRegion, subRegion.getReference< string >( viewKeyStruct::porousMaterialNamesString() ) );
369 
370  arrayView1d< real64 const > const bulkModulus = porousSolid.getBulkModulus();
371  arrayView1d< real64 const > const shearModulus = porousSolid.getShearModulus();
372  arrayView1d< real64 const > const biotCoefficient = porousSolid.getBiotCoefficient();
373 
374  real64 const stabilizationMultiplier = m_stabilizationMultiplier;
375 
376  forAll< parallelDevicePolicy<> >( subRegion.size(), [bulkModulus,
377  shearModulus,
378  biotCoefficient,
379  stabilizationMultiplier,
380  macroElementIndex,
381  elementStabConstant] GEOS_HOST_DEVICE ( localIndex const ei )
382 
383  {
384  real64 const bM = bulkModulus[ei];
385  real64 const sM = shearModulus[ei];
386  real64 const bC = biotCoefficient[ei];
387 
388  macroElementIndex[ei] = 1;
389  elementStabConstant[ei] = stabilizationMultiplier * 9.0 * (bC * bC) / (32.0 * (10.0 * sM / 3.0 + bM));
390  } );
391  } );
392  } );
393 
394  }
395 
396 protected:
397 
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,
407  real64 const dt,
408  PARAMS && ... params )
409  {
411 
412  NodeManager const & nodeManager = mesh.getNodeManager();
413 
414  string const dofKey = dofManager.getKey( fields::solidMechanics::totalDisplacement::key() );
415  arrayView1d< globalIndex const > const & dofNumber = nodeManager.getReference< globalIndex_array >( dofKey );
416 
417  real64 const gravityVectorData[3] = LVARRAY_TENSOROPS_INIT_LOCAL_3( PhysicsSolverBase::gravityVector() );
418 
419  KERNEL_WRAPPER kernelWrapper( dofNumber,
420  dofManager.rankOffset(),
421  localMatrix,
422  localRhs,
423  dt,
424  gravityVectorData,
425  std::forward< PARAMS >( params )... );
426 
427  return finiteElement::
428  regionBasedKernelApplication< parallelDevicePolicy< >,
429  CONSTITUTIVE_BASE,
430  CellElementSubRegion >( mesh,
431  regionNames,
432  this->solidMechanicsSolver()->getDiscretizationName(),
433  materialNamesString,
434  kernelWrapper );
435  }
436 
437  /* Implementation of Nonlinear Acceleration (Aitken) of averageMeanTotalStressIncrement */
438 
439  void recordAverageMeanTotalStressIncrement( DomainPartition & domain,
440  array1d< real64 > & averageMeanTotalStressIncrement )
441  {
442  averageMeanTotalStressIncrement.resize( 0 );
443  PhysicsSolverBase::forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
444  MeshLevel & mesh,
445  arrayView1d< string const > const & regionNames ) {
446  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
447  auto & subRegion ) {
448  // get the solid model (to access stress increment)
449  string const solidName = subRegion.template getReference< string >( "porousMaterialNames" );
450  constitutive::CoupledSolidBase & solid = PhysicsSolverBase::getConstitutiveModel< constitutive::CoupledSolidBase >(
451  subRegion, solidName );
452 
453  arrayView1d< const real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
454  for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
455  {
456  averageMeanTotalStressIncrement.emplace_back( averageMeanTotalStressIncrement_k[k] );
457  }
458  } );
459  } );
460  }
461 
462  void applyAcceleratedAverageMeanTotalStressIncrement( DomainPartition & domain,
463  array1d< real64 > & averageMeanTotalStressIncrement )
464  {
465  integer i = 0;
466  PhysicsSolverBase::forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
467  MeshLevel & mesh,
468  arrayView1d< string const > const & regionNames ) {
469  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
470  auto & subRegion ) {
471  // get the solid model (to access stress increment)
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++ )
478  {
479  porosityModel.updateAverageMeanTotalStressIncrement( k, averageMeanTotalStressIncrement[i] );
480  i++;
481  }
482  } );
483  } );
484  }
485 
486  real64 computeAitkenRelaxationFactor( array1d< real64 > const & s0,
487  array1d< real64 > const & s1,
488  array1d< real64 > const & s1_tilde,
489  array1d< real64 > const & s2_tilde,
490  real64 const omega0 )
491  {
492  array1d< real64 > r1 = axpy( s1_tilde, s0, -1.0 );
493  array1d< real64 > r2 = axpy( s2_tilde, s1, -1.0 );
494 
495  // diff = r2 - r1
496  array1d< real64 > diff = axpy( r2, r1, -1.0 );
497 
498  real64 const denom = dot( diff, diff );
499  real64 const numer = dot( r1, diff );
500 
501  real64 omega1 = 1.0;
502  if( !isZero( denom ))
503  {
504  omega1 = -1.0 * omega0 * numer / denom;
505  }
506  return omega1;
507  }
508 
509  array1d< real64 > computeUpdate( array1d< real64 > const & s1,
510  array1d< real64 > const & s2_tilde,
511  real64 const omega1 )
512  {
513  return axpy( scale( s1, 1.0 - omega1 ),
514  scale( s2_tilde, omega1 ),
515  1.0 );
516  }
517 
518  void startSequentialIteration( integer const & iter,
519  DomainPartition & domain ) override
520  {
521  if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
522  {
523  if( iter == 0 )
524  {
525  recordAverageMeanTotalStressIncrement( domain, m_s1 );
526  }
527  else
528  {
529  m_s0 = m_s1;
530  m_s1 = m_s2;
531  m_s1_tilde = m_s2_tilde;
532  m_omega0 = m_omega1;
533  }
534  }
535  }
536 
537  void finishSequentialIteration( integer const & iter,
538  DomainPartition & domain ) override
539  {
540  if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
541  {
542  if( iter == 0 )
543  {
544  m_s2 = m_s2_tilde;
545  m_omega1 = 1.0;
546  }
547  else
548  {
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 );
552  }
553  }
554  }
555 
556  virtual void mapSolutionBetweenSolvers( DomainPartition & domain, integer const solverType ) override
557  {
559 
561  if( solverType == static_cast< integer >( SolverType::Flow ) )
562  {
563  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
564  MeshLevel & mesh,
565  arrayView1d< string const > const & regionNames )
566  {
567 
568  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
569  auto & subRegion )
570  {
571  // update the bulk density
572  // TODO: ideally, we would not recompute the bulk density, but a more general "rhs" containing the body force and the
573  // pressure/temperature terms
574  updateBulkDensity( subRegion );
575  } );
576  } );
577  }
578 
580  if( solverType == static_cast< integer >( SolverType::SolidMechanics )
581  && !m_performStressInitialization ) // do not update during poromechanics initialization
582  {
583  // compute the average of the mean total stress increment over quadrature points
584  averageMeanTotalStressIncrement( domain );
585 
586  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
587  MeshLevel & mesh,
588  arrayView1d< string const > const & regionNames )
589  {
590 
591  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
592  auto & subRegion )
593  {
594  // update the porosity after a change in displacement (after mechanics solve)
595  // or a change in pressure/temperature (after a flow solve)
596  flowSolver()->updatePorosityAndPermeability( subRegion );
597  // update bulk density to reflect porosity change into mechanics
598  updateBulkDensity( subRegion );
599  } );
600  } );
601  }
602 
603  // needed to perform nonlinear acceleration
604  if( solverType == static_cast< integer >( SolverType::SolidMechanics ) &&
605  this->getNonlinearSolverParameters().m_nonlinearAccelerationType== NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
606  {
607  recordAverageMeanTotalStressIncrement( domain, m_s2_tilde );
608  }
609  }
610 
616  {
617  this->template forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
618  MeshLevel & mesh,
619  arrayView1d< string const > const & regionNames )
620  {
621  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
622  auto & subRegion )
623  {
624  // get the solid model (to access stress increment)
625  string const solidName = subRegion.template getReference< string >( viewKeyStruct::porousMaterialNamesString() );
626  constitutive::CoupledSolidBase & solid = this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
627 
628  arrayView2d< real64 const > const meanTotalStressIncrement_k = solid.getMeanTotalStressIncrement_k();
629  arrayView1d< real64 > const averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
630 
631  finiteElement::FiniteElementBase & subRegionFE =
632  subRegion.template getReference< finiteElement::FiniteElementBase >( solidMechanicsSolver()->getDiscretizationName() );
633 
634  // determine the finite element type
635  finiteElement::FiniteElementDispatchHandler< BASE_FE_TYPES >::
636  dispatch3D( subRegionFE, [&] ( auto const finiteElement )
637  {
638  using FE_TYPE = decltype( finiteElement );
639 
640  // call the factory and launch the kernel
641  AverageOverQuadraturePoints1DKernelFactory::
642  createAndLaunch< CellElementSubRegion,
643  FE_TYPE,
644  parallelDevicePolicy<> >( mesh.getNodeManager(),
645  mesh.getEdgeManager(),
646  mesh.getFaceManager(),
647  subRegion,
648  finiteElement,
649  meanTotalStressIncrement_k,
650  averageMeanTotalStressIncrement_k );
651  } );
652  } );
653  } );
654  }
655 
656  virtual void updateBulkDensity( ElementSubRegionBase & subRegion ) = 0;
657 
658  virtual void validateNonlinearAcceleration() override
659  {
660  if( MpiWrapper::commSize( MPI_COMM_GEOS ) > 1 )
661  {
662  GEOS_ERROR( "Nonlinear acceleration is not implemented for MPI runs" );
663  }
664  }
665 
668 
671 
673  stabilization::StabilizationType m_stabilizationType;
674 
677 
680 
682  array1d< real64 > m_s0; // Accelerated averageMeanTotalStresIncrement @ outer iteration v ( two iterations ago )
683  array1d< real64 > m_s1; // Accelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( previous iteration )
684  array1d< real64 > m_s1_tilde; // Unaccelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( previous iteration )
685  array1d< real64 > m_s2; // Accelerated averageMeanTotalStresIncrement @ outer iteration v + 2 ( current iteration )
686  array1d< real64 > m_s2_tilde; // Unaccelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( current iteration )
687  real64 m_omega0; // Old Aitken relaxation factor
688  real64 m_omega1; // New Aitken relaxation factor
689 
690 };
691 
692 } /* namespace geos */
693 
694 #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< 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.
Definition: Group.hpp:1540
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:1343
void setRestartFlags(RestartFlags flags)
Set flags that control restart output of this group.
Definition: Group.hpp:1424
string const & getName() const
Get group name.
Definition: Group.hpp:1329
Group & setSizedFromParent(int val)
Set whether this wrapper is resized when its parent is resized.
Definition: Group.hpp:1409
DataContext const & getWrapperDataContext(KEY key) const
Definition: Group.hpp:1354
virtual void initializePostInitialConditionsPreSubGroups()
Called by InitializePostInitialConditions() prior to initializing sub-Groups.
Definition: Group.hpp:1552
@ 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.
Definition: DataTypes.hpp:180
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:401
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.
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:176
Exception class used to report errors in user input.
Definition: Logger.hpp:502
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.
Definition: Group.hpp:1442