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/contact/HydraulicApertureBase.hpp"
30 #include "mesh/DomainPartition.hpp"
32 #include "codingUtilities/Utilities.hpp"
33 
34 namespace geos
35 {
36 
37 namespace stabilization
38 {
39 enum class StabilizationType : integer
40 {
41  None,
42  Global,
43  Local,
44 };
45 
46 ENUM_STRINGS( StabilizationType,
47  "None",
48  "Global",
49  "Local" );
50 }
51 
52 
53 template< typename FLOW_SOLVER, typename MECHANICS_SOLVER = SolidMechanicsLagrangianFEM >
54 class PoromechanicsSolver : public CoupledSolver< FLOW_SOLVER, MECHANICS_SOLVER >
55 {
56 public:
57 
59  using Base::m_solvers;
60  using Base::m_dofManager;
61  using Base::m_localMatrix;
62  using Base::m_rhs;
63  using Base::m_solution;
64 
65  enum class SolverType : integer
66  {
67  Flow = 0,
68  SolidMechanics = 1
69  };
70 
72  static string coupledSolverAttributePrefix() { return "poromechanics"; }
73 
79  PoromechanicsSolver( const string & name,
80  dataRepository::Group * const parent )
81  : Base( name, parent ),
82  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 
92  setDescription( "StabilizationType. Options are:\n" +
93  toString( stabilization::StabilizationType::None ) + "- Add no stabilization to mass equation \n" +
94  toString( stabilization::StabilizationType::Global ) + "- Add jump stabilization to all faces \n" +
95  toString( stabilization::StabilizationType::Local ) + "- Add jump stabilization on interior of macro elements" );
96 
98  setRTTypeName( rtTypes::CustomTypes::groupNameRefArray ).
100  setDescription( "Regions where stabilization is applied." );
101 
103  setApplyDefaultValue( 1.0 ).
104  setInputFlag( dataRepository::InputFlags::OPTIONAL ).
105  setDescription( "Constant multiplier of stabilization strength" );
106 
107  LinearSolverParameters & linearSolverParameters = this->m_linearSolverParameters.get();
108  linearSolverParameters.dofsPerNode = 3;
109  linearSolverParameters.multiscale.label = "poro";
110  }
111 
112  virtual void postInputInitialization() override
113  {
115 
116  GEOS_THROW_IF( this->m_isThermal && !this->flowSolver()->isThermal(),
117  GEOS_FMT( "{} {}: The attribute `{}` of the flow solver must be thermal since the poromechanics solver is thermal",
118  this->getCatalogName(), this->getName(), this->flowSolver()->getName() ),
119  InputError );
120 
122  GEOS_FMT( "{} {}: The attribute `{}` of solid mechanics solver `{}` must be `{}`",
123  this->getCatalogName(), this->getName(),
124  SolidMechanicsLagrangianFEM::viewKeyStruct::timeIntegrationOptionString(),
125  this->solidMechanicsSolver()->getName(),
127  InputError );
128  }
129 
130  virtual void setConstitutiveNamesCallSuper( ElementSubRegionBase & subRegion ) const override final
131  {
133 
134  this->template setConstitutiveName< constitutive::CoupledSolidBase >( subRegion,
135  viewKeyStruct::porousMaterialNamesString(), "coupled solid" );
136 
137  // This is needed by the way the surface generator currently does things.
138  this->template setConstitutiveName< constitutive::PorosityBase >( subRegion,
139  constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString(), "porosity" );
140 
141  if( dynamic_cast< SurfaceElementSubRegion * >( &subRegion ) )
142  {
143  this->template setConstitutiveName< constitutive::HydraulicApertureBase >( subRegion,
144  viewKeyStruct::hydraulicApertureRelationNameString(), "hydraulic aperture" );
145  }
146  }
147 
148  virtual void initializePreSubGroups() override
149  {
151 
152  GEOS_THROW_IF( m_stabilizationType == stabilization::StabilizationType::Local,
154  ": Local stabilization has been temporarily disabled",
155  InputError );
156 
157  DomainPartition & domain = this->template getGroupByPath< DomainPartition >( "/Problem/domain" );
158 
159  this->template forDiscretizationOnMeshTargets<>( domain.getMeshBodies(), [&] ( string const &,
160  MeshLevel & mesh,
161  string_array const & regionNames )
162  {
163  ElementRegionManager & elementRegionManager = mesh.getElemManager();
164  elementRegionManager.forElementSubRegions< CellElementSubRegion >( regionNames,
165  [&]( localIndex const,
166  ElementSubRegionBase & subRegion )
167  {
168  if( subRegion.hasField< fields::poromechanics::bulkDensity >() )
169  {
170  // get the solid model to know the number of quadrature points and resize the bulk density
171  constitutive::CoupledSolidBase const & solid =
172  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion,
173  subRegion.getReference< string >( viewKeyStruct::porousMaterialNamesString() ) );
174  subRegion.getField< fields::poromechanics::bulkDensity >().resizeDimension< 1 >( solid.getDensity().size( 1 ) );
175  }
176  } );
177  } );
178  }
179 
180  virtual void registerDataOnMesh( dataRepository::Group & meshBodies ) override
181  {
182  Base::registerDataOnMesh( meshBodies );
183 
184  if( this->getNonlinearSolverParameters().m_couplingType == NonlinearSolverParameters::CouplingType::Sequential )
185  {
186  // to let the solid mechanics solver that there is a pressure and temperature RHS in the mechanics solve
187  solidMechanicsSolver()->enableFixedStressPoromechanicsUpdate();
188  // to let the flow solver that saving pressure_k and temperature_k is necessary (for the fixed-stress porosity terms)
189  flowSolver()->enableFixedStressPoromechanicsUpdate();
190  }
191 
192  if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
193  {
194  flowSolver()->enableJumpStabilization();
195  }
196 
197  this->forDiscretizationOnMeshTargets( meshBodies, [&] ( string const &,
198  MeshLevel & mesh,
199  string_array const & regionNames )
200  {
201  ElementRegionManager & elemManager = mesh.getElemManager();
202 
203  elemManager.forElementSubRegions< CellElementSubRegion >( regionNames,
204  [&]( localIndex const,
205  ElementSubRegionBase & subRegion )
206  {
207  subRegion.registerWrapper< string >( viewKeyStruct::porousMaterialNamesString() ).
208  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
209  setRestartFlags( dataRepository::RestartFlags::NO_WRITE ).
210  setSizedFromParent( 0 );
211 
212  // This is needed by the way the surface generator currently does things.
213  subRegion.registerWrapper< string >( constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString() ).
214  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
215  setRestartFlags( dataRepository::RestartFlags::NO_WRITE ).
216  setSizedFromParent( 0 );
217 
218  // register the bulk density for use in the solid mechanics solver
219  // ideally we would resize it here as well, but the solid model name is not available yet (see below)
220  subRegion.registerField< fields::poromechanics::bulkDensity >( this->getName() );
221 
222  if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
223  {
224  subRegion.registerField< fields::flow::macroElementIndex >( this->getName());
225  subRegion.registerField< fields::flow::elementStabConstant >( this->getName());
226  }
227  } );
228  } );
229  }
230 
231  virtual void implicitStepSetup( real64 const & time_n,
232  real64 const & dt,
233  DomainPartition & domain ) override
234  {
235  flowSolver()->setKeepVariablesConstantDuringInitStep( m_performStressInitialization );
236 
237  if( this->m_stabilizationType == stabilization::StabilizationType::Global || this->m_stabilizationType == stabilization::StabilizationType::Local )
238  {
239  this->updateStabilizationParameters( domain );
240  }
241 
242  Base::implicitStepSetup( time_n, dt, domain );
243  }
244 
245  virtual void setupDofs( DomainPartition const & domain,
246  DofManager & dofManager ) const override
247  {
248  // note that the order of operations matters a lot here (for instance for the MGR labels)
249  // we must set up dofs for solid mechanics first, and then for flow
250  // that's the reason why this function is here and not in CoupledSolvers.hpp
251  solidMechanicsSolver()->setupDofs( domain, dofManager );
252  flowSolver()->setupDofs( domain, dofManager );
253  this->setupCoupling( domain, dofManager );
254  }
255 
256  virtual bool checkSequentialConvergence( integer const cycleNumber,
257  integer 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  bool isConverged = Base::checkSequentialConvergence( cycleNumber, iter, time_n, dt, domain );
268 
269  // restore original
270  subcycling = subcycling_orig;
271 
272  return isConverged;
273  }
274 
279  MECHANICS_SOLVER * solidMechanicsSolver() const
280  {
281  return std::get< toUnderlying( SolverType::SolidMechanics ) >( m_solvers );
282  }
283 
288  FLOW_SOLVER * flowSolver() const
289  {
290  return std::get< toUnderlying( SolverType::Flow ) >( m_solvers );
291  }
292 
293  /*
294  * @brief Utility function to set the stress initialization flag
295  * @param[in] performStressInitialization true if the solver has to initialize stress, false otherwise
296  */
297  void setStressInitialization( bool const performStressInitialization )
298  {
299  m_performStressInitialization = performStressInitialization;
300  solidMechanicsSolver()->setStressInitialization( performStressInitialization );
301  }
302 
304  {
306  constexpr static char const * porousMaterialNamesString() { return "porousMaterialNames"; }
307 
309  constexpr static char const * isThermalString() { return "isThermal"; }
310 
312  constexpr static char const * stabilizationTypeString() {return "stabilizationType"; }
313 
315  constexpr static const char * stabilizationRegionNamesString() {return "stabilizationRegionNames"; }
316 
318  constexpr static const char * stabilizationMultiplierString() {return "stabilizationMultiplier"; }
319 
321  static constexpr char const * hydraulicApertureRelationNameString() {return "hydraulicApertureRelationName"; }
322 
323  };
324 
325  void updateStabilizationParameters( DomainPartition & domain ) const
326  {
327  // Step 1: we loop over the regions where stabilization is active and collect their name
328  set< string > regionFilter;
329  for( string const & regionName : m_stabilizationRegionNames )
330  {
331  regionFilter.insert( regionName );
332  }
333 
334  // Step 2: loop over target regions of solver, and tag the elements belonging to the stabilization regions
335  this->template forDiscretizationOnMeshTargets<>( domain.getMeshBodies(), [&] ( string const &,
336  MeshLevel & mesh,
337  string_array const & targetRegionNames )
338  {
339  //keep only target regions in filter
340  string_array filteredTargetRegionNames;
341  filteredTargetRegionNames.reserve( targetRegionNames.size() );
342 
343  for( string const & targetRegionName : targetRegionNames )
344  {
345  if( regionFilter.count( targetRegionName ) )
346  {
347  filteredTargetRegionNames.emplace_back( targetRegionName );
348  }
349  }
350 
351  // Loop over elements and update stabilization constant
352  mesh.getElemManager().forElementSubRegions( filteredTargetRegionNames, [&]( localIndex const,
353  ElementSubRegionBase & subRegion )
354  {
355  arrayView1d< integer > const macroElementIndex = subRegion.getField< fields::flow::macroElementIndex >();
356  arrayView1d< real64 > const elementStabConstant = subRegion.getField< fields::flow::elementStabConstant >();
357 
358  constitutive::CoupledSolidBase const & porousSolid =
359  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion,
360  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 TYPE_LIST,
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( this->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  TYPE_LIST >( mesh,
439  regionNames,
440  this->solidMechanicsSolver()->getDiscretizationName(),
441  materialNamesString,
442  kernelWrapper );
443  }
444 
445  /* Implementation of Nonlinear Acceleration (Aitken) of averageMeanTotalStressIncrement */
446 
447  void recordAverageMeanTotalStressIncrement( DomainPartition & domain,
448  array1d< real64 > & averageMeanTotalStressIncrement )
449  {
450  averageMeanTotalStressIncrement.resize( 0 );
451  this->forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
452  MeshLevel & mesh,
453  string_array const & regionNames )
454  {
455  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
456  auto & subRegion )
457  {
458  // get the solid model (to access stress increment)
459  string const & solidName =
460  subRegion.template getReference< string >( viewKeyStruct::porousMaterialNamesString() );
461  constitutive::CoupledSolidBase & solid =
462  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
463 
464  arrayView1d< const real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
465  for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
466  {
467  averageMeanTotalStressIncrement.emplace_back( averageMeanTotalStressIncrement_k[k] );
468  }
469  } );
470  } );
471  }
472 
473  void applyAcceleratedAverageMeanTotalStressIncrement( DomainPartition & domain,
474  array1d< real64 > & averageMeanTotalStressIncrement )
475  {
476  integer i = 0;
477  this->forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
478  MeshLevel & mesh,
479  string_array const & regionNames )
480  {
481  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
482  auto & subRegion )
483  {
484  // get the solid model (to access stress increment)
485  string const & solidName = subRegion.template getReference< string >( viewKeyStruct::porousMaterialNamesString() );
486  constitutive::CoupledSolidBase & solid =
487  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
488  auto & porosityModel = dynamic_cast< constitutive::BiotPorosity const & >( solid.getBasePorosityModel());
489  arrayView1d< real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
490  for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
491  {
492  porosityModel.updateAverageMeanTotalStressIncrement( k, averageMeanTotalStressIncrement[i] );
493  i++;
494  }
495  } );
496  } );
497  }
498 
499  real64 computeAitkenRelaxationFactor( array1d< real64 > const & s0,
500  array1d< real64 > const & s1,
501  array1d< real64 > const & s1_tilde,
502  array1d< real64 > const & s2_tilde,
503  real64 const omega0 )
504  {
505  array1d< real64 > r1 = axpy( s1_tilde, s0, -1.0 );
506  array1d< real64 > r2 = axpy( s2_tilde, s1, -1.0 );
507 
508  // diff = r2 - r1
509  array1d< real64 > diff = axpy( r2, r1, -1.0 );
510 
511  real64 const denom = dot( diff, diff );
512  real64 const numer = dot( r1, diff );
513 
514  real64 omega1 = 1.0;
515  if( !isZero( denom ))
516  {
517  omega1 = -1.0 * omega0 * numer / denom;
518  }
519  return omega1;
520  }
521 
522  array1d< real64 > computeUpdate( array1d< real64 > const & s1,
523  array1d< real64 > const & s2_tilde,
524  real64 const omega1 )
525  {
526  return axpy( scale( s1, 1.0 - omega1 ),
527  scale( s2_tilde, omega1 ),
528  1.0 );
529  }
530 
531  void startSequentialIteration( integer const & iter,
532  DomainPartition & domain ) override
533  {
534  if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
535  {
536  if( iter == 0 )
537  {
538  recordAverageMeanTotalStressIncrement( domain, m_s1 );
539  }
540  else
541  {
542  m_s0 = m_s1;
543  m_s1 = m_s2;
544  m_s1_tilde = m_s2_tilde;
545  m_omega0 = m_omega1;
546  }
547  }
548  }
549 
550  void finishSequentialIteration( integer const & iter,
551  DomainPartition & domain ) override
552  {
553  if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
554  {
555  if( iter == 0 )
556  {
557  m_s2 = m_s2_tilde;
558  m_omega1 = 1.0;
559  }
560  else
561  {
562  m_omega1 = computeAitkenRelaxationFactor( m_s0, m_s1, m_s1_tilde, m_s2_tilde, m_omega0 );
563  m_s2 = computeUpdate( m_s1, m_s2_tilde, m_omega1 );
564  applyAcceleratedAverageMeanTotalStressIncrement( domain, m_s2 );
565  }
566  }
567  }
568 
569  virtual void mapSolutionBetweenSolvers( DomainPartition & domain, integer const solverType ) override
570  {
572 
574  if( solverType == static_cast< integer >( SolverType::Flow ) )
575  {
576  updateBulkDensity( domain );
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  string_array 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  string_array 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 =
627  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
628 
629  arrayView2d< real64 const > const meanTotalStressIncrement_k = solid.getMeanTotalStressIncrement_k();
630  arrayView1d< real64 > const averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
631 
632  finiteElement::FiniteElementBase & subRegionFE =
633  subRegion.template getReference< finiteElement::FiniteElementBase >( solidMechanicsSolver()->getDiscretizationName() );
634 
635  // determine the finite element type
636  finiteElement::FiniteElementDispatchHandler< BASE_FE_TYPES >::
637  dispatch3D( subRegionFE, [&] ( auto const finiteElement )
638  {
639  using FE_TYPE = decltype( finiteElement );
640 
641  // call the factory and launch the kernel
642  AverageOverQuadraturePoints1DKernelFactory::
643  createAndLaunch< CellElementSubRegion,
644  FE_TYPE,
645  parallelDevicePolicy<> >( mesh.getNodeManager(),
646  mesh.getEdgeManager(),
647  mesh.getFaceManager(),
648  subRegion,
649  finiteElement,
650  meanTotalStressIncrement_k,
651  averageMeanTotalStressIncrement_k );
652  } );
653  } );
654  } );
655  }
656 
657  virtual void updateBulkDensity( ElementSubRegionBase & subRegion ) = 0;
658 
659  virtual void validateNonlinearAcceleration() override
660  {
661  if( MpiWrapper::commSize( MPI_COMM_GEOS ) > 1 )
662  {
663  GEOS_ERROR( "Nonlinear acceleration is not implemented for MPI runs" );
664  }
665  }
666 
669 
672 
674  stabilization::StabilizationType m_stabilizationType;
675 
678 
681 
683  array1d< real64 > m_s0; // Accelerated averageMeanTotalStresIncrement @ outer iteration v ( two iterations ago )
684  array1d< real64 > m_s1; // Accelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( previous iteration )
685  array1d< real64 > m_s1_tilde; // Unaccelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( previous iteration )
686  array1d< real64 > m_s2; // Accelerated averageMeanTotalStresIncrement @ outer iteration v + 2 ( current iteration )
687  array1d< real64 > m_s2_tilde; // Unaccelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( current iteration )
688  real64 m_omega0; // Old Aitken relaxation factor
689  real64 m_omega1; // New Aitken relaxation factor
690 
691 };
692 
693 } /* namespace geos */
694 
695 #endif //GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSSOLVER_HPP_
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
Definition: Logger.hpp:157
#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.
virtual void postInputInitialization() override
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:45
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 string getCatalogName() const =0
virtual void setConstitutiveNamesCallSuper(ElementSubRegionBase &subRegion) const
This function sets constitutive name fields on an ElementSubRegionBase, and calls the base function i...
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.
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.
virtual void postInputInitialization() override
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 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.
string const & getName() const
Get group name.
Definition: Group.hpp:1331
DataContext const & getWrapperDataContext(KEY key) const
Definition: Group.hpp:1356
@ NOPLOT
Do not ever write to plot file.
@ OPTIONAL
Optional in input.
@ NO_WRITE
Do not write into restart.
stdVector< string > string_array
A 1-dimensional array of geos::string types.
Definition: DataTypes.hpp:361
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:370
std::set< T > set
A set of local indices.
Definition: DataTypes.hpp:262
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:84
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:175
ENUM_STRINGS(LinearSolverParameters::SolverType, "direct", "cg", "gmres", "fgmres", "bicgstab", "richardson", "preconditioner")
Declare strings associated with enumeration values.
Provides enum <-> string conversion facilities.
Exception class used to report errors in user input.
Definition: Logger.hpp:464
string label
User-displayed label of the scheme.
Set of parameters for a linear solver or preconditioner.
integer dofsPerNode
Dofs per node (or support location) for non-scalar problems.
struct geos::LinearSolverParameters::Multiscale multiscale
Multiscale preconditioner parameters.
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