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 
29 #include "constitutive/solid/CoupledSolidBase.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 {
42  None,
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 ),
85  {
87  setApplyDefaultValue( 0 ).
89  setDescription( "Flag indicating whether the problem is thermal or not. Set isThermal=\"1\" to enable the thermal coupling" );
90 
93  setDescription( "StabilizationType. Options are:\n" +
94  toString( stabilization::StabilizationType::None ) + "- Add no stabilization to mass equation \n" +
95  toString( stabilization::StabilizationType::Global ) + "- Add jump stabilization to all faces \n" +
96  toString( stabilization::StabilizationType::Local ) + "- Add jump stabilization on interior of macro elements" );
97 
99  setRTTypeName( rtTypes::CustomTypes::groupNameRefArray ).
100  setInputFlag( dataRepository::InputFlags::OPTIONAL ).
101  setDescription( "Regions where stabilization is applied." );
102 
104  setApplyDefaultValue( 1.0 ).
105  setInputFlag( dataRepository::InputFlags::OPTIONAL ).
106  setDescription( "Constant multiplier of stabilization strength" );
107 
108  LinearSolverParameters & linearSolverParameters = this->m_linearSolverParameters.get();
109  linearSolverParameters.dofsPerNode = 3;
110  linearSolverParameters.multiscale.label = "poro";
111  }
112 
113  virtual void postInputInitialization() override
114  {
116 
118  GEOS_FMT( "The attribute `{}` of solid mechanics solver `{}` must be `{}`",
119  SolidMechanicsLagrangianFEM::viewKeyStruct::timeIntegrationOptionString(),
120  this->solidMechanicsSolver()->getName(),
123 
124  setMGRStrategy();
125  }
126 
127  virtual void setConstitutiveNamesCallSuper( ElementSubRegionBase & subRegion ) const override final
128  {
130 
131  this->template setConstitutiveName< constitutive::CoupledSolidBase >( subRegion,
132  viewKeyStruct::porousMaterialNamesString(), "coupled solid" );
133 
134  // This is needed by the way the surface generator currently does things.
135  this->template setConstitutiveName< constitutive::PorosityBase >( subRegion,
136  constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString(), "porosity" );
137 
138  if( dynamic_cast< SurfaceElementSubRegion * >( &subRegion ) )
139  {
140  this->template setConstitutiveName< constitutive::HydraulicApertureBase >( subRegion,
141  viewKeyStruct::hydraulicApertureRelationNameString(), "hydraulic aperture" );
142  }
143  }
144 
145  virtual void initializePreSubGroups() override
146  {
148 
149  GEOS_THROW_IF( this->m_isThermal && !this->flowSolver()->isThermal(),
150  GEOS_FMT( "The attribute `{}` of the flow solver must be thermal since the poromechanics solver is thermal",
151  this->flowSolver()->getName() ),
152  InputError, this->flowSolver()->getDataContext() );
153 
154  GEOS_THROW_IF( m_stabilizationType == stabilization::StabilizationType::Local,
156  ": Local stabilization has been temporarily disabled",
158 
159  DomainPartition & domain = this->template getGroupByPath< DomainPartition >( "/Problem/domain" );
160 
161  this->template forDiscretizationOnMeshTargets<>( domain.getMeshBodies(), [&] ( string const &,
162  MeshLevel & mesh,
163  string_array const & regionNames )
164  {
165  ElementRegionManager & elementRegionManager = mesh.getElemManager();
166  elementRegionManager.forElementSubRegions< CellElementSubRegion >( regionNames,
167  [&]( localIndex const,
168  ElementSubRegionBase & subRegion )
169  {
170  if( subRegion.hasField< fields::poromechanics::bulkDensity >() )
171  {
172  // get the solid model to know the number of quadrature points and resize the bulk density
173  constitutive::CoupledSolidBase const & solid =
174  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion,
175  subRegion.getReference< string >( viewKeyStruct::porousMaterialNamesString() ) );
176  subRegion.getField< fields::poromechanics::bulkDensity >().resizeDimension< 1 >( solid.getDensity().size( 1 ) );
177  }
178  } );
179  } );
180  }
181 
182  virtual void registerDataOnMesh( dataRepository::Group & meshBodies ) override
183  {
184  Base::registerDataOnMesh( meshBodies );
185 
186  if( this->getNonlinearSolverParameters().m_couplingType == NonlinearSolverParameters::CouplingType::Sequential )
187  {
188  // to let the solid mechanics solver that there is a pressure and temperature RHS in the mechanics solve
189  solidMechanicsSolver()->enableFixedStressPoromechanicsUpdate();
190  // to let the flow solver that saving pressure_k and temperature_k is necessary (for the fixed-stress porosity terms)
191  flowSolver()->enableFixedStressPoromechanicsUpdate();
192  }
193 
194  if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
195  {
196  flowSolver()->enableJumpStabilization();
197  }
198 
199  this->forDiscretizationOnMeshTargets( meshBodies, [&] ( string const &,
200  MeshLevel & mesh,
201  string_array const & regionNames )
202  {
203  ElementRegionManager & elemManager = mesh.getElemManager();
204 
205  elemManager.forElementSubRegions< CellElementSubRegion >( regionNames,
206  [&]( localIndex const,
207  ElementSubRegionBase & subRegion )
208  {
209  subRegion.registerWrapper< string >( viewKeyStruct::porousMaterialNamesString() ).
210  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
211  setRestartFlags( dataRepository::RestartFlags::NO_WRITE ).
212  setSizedFromParent( 0 );
213 
214  // This is needed by the way the surface generator currently does things.
215  subRegion.registerWrapper< string >( constitutive::CoupledSolidBase::viewKeyStruct::porosityModelNameString() ).
216  setPlotLevel( dataRepository::PlotLevel::NOPLOT ).
217  setRestartFlags( dataRepository::RestartFlags::NO_WRITE ).
218  setSizedFromParent( 0 );
219 
220  // register the bulk density for use in the solid mechanics solver
221  // ideally we would resize it here as well, but the solid model name is not available yet (see below)
222  subRegion.registerField< fields::poromechanics::bulkDensity >( this->getName() );
223 
224  if( m_stabilizationType == stabilization::StabilizationType::Global || m_stabilizationType == stabilization::StabilizationType::Local )
225  {
226  subRegion.registerField< fields::flow::macroElementIndex >( this->getName());
227  subRegion.registerField< fields::flow::elementStabConstant >( this->getName());
228  }
229  } );
230  } );
231  }
232 
233  virtual void implicitStepSetup( real64 const & time_n,
234  real64 const & dt,
235  DomainPartition & domain ) override
236  {
237  flowSolver()->setKeepVariablesConstantDuringInitStep( m_performStressInitialization );
238 
239  if( this->m_stabilizationType == stabilization::StabilizationType::Global || this->m_stabilizationType == stabilization::StabilizationType::Local )
240  {
241  this->updateStabilizationParameters( domain );
242  }
243 
244  Base::implicitStepSetup( time_n, dt, domain );
245  }
246 
247  virtual void setupDofs( DomainPartition const & domain,
248  DofManager & dofManager ) const override
249  {
250  // note that the order of operations matters a lot here (for instance for the MGR labels)
251  // we must set up dofs for solid mechanics first, and then for flow
252  // that's the reason why this function is here and not in CoupledSolvers.hpp
253  solidMechanicsSolver()->setupDofs( domain, dofManager );
254  flowSolver()->setupDofs( domain, dofManager );
255  this->setupCoupling( domain, dofManager );
256  }
257 
258  virtual void setupCoupling( DomainPartition const & GEOS_UNUSED_PARAM( domain ),
259  DofManager & dofManager ) const override
260  {
261  dofManager.addCoupling( fields::solidMechanics::totalDisplacement::key(),
262  getFlowDofKey(),
264  }
265 
266  virtual void assembleSystem( real64 const time,
267  real64 const dt,
268  DomainPartition & domain,
269  DofManager const & dofManager,
270  CRSMatrixView< real64, globalIndex const > const & localMatrix,
271  arrayView1d< real64 > const & localRhs ) override
272  {
274 
275  // Steps 1 and 2: compute element-based terms (mechanics and local flow terms)
276  assembleElementBasedTerms( time,
277  dt,
278  domain,
279  dofManager,
280  localMatrix,
281  localRhs );
282 
283  // Step 3: compute the fluxes (face-based contributions)
284 
285  if( m_stabilizationType == stabilization::StabilizationType::Global ||
286  m_stabilizationType == stabilization::StabilizationType::Local )
287  {
288  this->flowSolver()->assembleStabilizedFluxTerms( dt,
289  domain,
290  dofManager,
291  localMatrix,
292  localRhs );
293  }
294  else
295  {
296  this->flowSolver()->assembleFluxTerms( dt,
297  domain,
298  dofManager,
299  localMatrix,
300  localRhs );
301  }
302  }
303 
304  virtual void assembleElementBasedTerms( real64 const time_n,
305  real64 const dt,
306  DomainPartition & domain,
307  DofManager const & dofManager,
308  CRSMatrixView< real64, globalIndex const > const & localMatrix,
309  arrayView1d< real64 > const & localRhs ) = 0;
310 
311  virtual bool checkSequentialConvergence( integer const cycleNumber,
312  integer const iter,
313  real64 const & time_n,
314  real64 const & dt,
315  DomainPartition & domain ) override
316  {
317  // always force outer loop for initialization
318  auto & subcycling = this->getNonlinearSolverParameters().m_subcyclingOption;
319  auto const subcycling_orig = subcycling;
320  if( m_performStressInitialization )
321  subcycling = 1;
322  bool isConverged = Base::checkSequentialConvergence( cycleNumber, iter, time_n, dt, domain );
323 
324  // restore original
325  subcycling = subcycling_orig;
326 
327  return isConverged;
328  }
329 
334  MECHANICS_SOLVER * solidMechanicsSolver() const
335  {
336  return std::get< toUnderlying( SolverType::SolidMechanics ) >( m_solvers );
337  }
338 
343  FLOW_SOLVER * flowSolver() const
344  {
345  return std::get< toUnderlying( SolverType::Flow ) >( m_solvers );
346  }
347 
348  /*
349  * @brief Utility function to set the stress initialization flag
350  * @param[in] performStressInitialization true if the solver has to initialize stress, false otherwise
351  */
352  void setStressInitialization( bool const performStressInitialization )
353  {
354  m_performStressInitialization = performStressInitialization;
355  solidMechanicsSolver()->setStressInitialization( performStressInitialization );
356  }
357 
359  {
361  constexpr static char const * porousMaterialNamesString() { return "porousMaterialNames"; }
362 
364  constexpr static char const * isThermalString() { return "isThermal"; }
365 
367  constexpr static char const * stabilizationTypeString() {return "stabilizationType"; }
368 
370  constexpr static const char * stabilizationRegionNamesString() {return "stabilizationRegionNames"; }
371 
373  constexpr static const char * stabilizationMultiplierString() {return "stabilizationMultiplier"; }
374 
376  static constexpr char const * hydraulicApertureRelationNameString() {return "hydraulicApertureRelationName"; }
377 
378  };
379 
380  void updateStabilizationParameters( DomainPartition & domain ) const
381  {
382  // Step 1: we loop over the regions where stabilization is active and collect their name
383  set< string > regionFilter;
384  for( string const & regionName : m_stabilizationRegionNames )
385  {
386  regionFilter.insert( regionName );
387  }
388 
389  // Step 2: loop over target regions of solver, and tag the elements belonging to the stabilization regions
390  this->template forDiscretizationOnMeshTargets<>( domain.getMeshBodies(), [&] ( string const &,
391  MeshLevel & mesh,
392  string_array const & targetRegionNames )
393  {
394  //keep only target regions in filter
395  string_array filteredTargetRegionNames;
396  filteredTargetRegionNames.reserve( targetRegionNames.size() );
397 
398  for( string const & targetRegionName : targetRegionNames )
399  {
400  if( regionFilter.count( targetRegionName ) )
401  {
402  filteredTargetRegionNames.emplace_back( targetRegionName );
403  }
404  }
405 
406  // Loop over elements and update stabilization constant
407  mesh.getElemManager().forElementSubRegions( filteredTargetRegionNames, [&]( localIndex const,
408  ElementSubRegionBase & subRegion )
409  {
410  arrayView1d< integer > const macroElementIndex = subRegion.getField< fields::flow::macroElementIndex >();
411  arrayView1d< real64 > const elementStabConstant = subRegion.getField< fields::flow::elementStabConstant >();
412 
413  constitutive::CoupledSolidBase const & porousSolid =
414  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion,
415  subRegion.getReference< string >( viewKeyStruct::porousMaterialNamesString() ) );
416 
417  arrayView1d< real64 const > const bulkModulus = porousSolid.getBulkModulus();
418  arrayView1d< real64 const > const shearModulus = porousSolid.getShearModulus();
419  arrayView1d< real64 const > const biotCoefficient = porousSolid.getBiotCoefficient();
420 
421  real64 const stabilizationMultiplier = m_stabilizationMultiplier;
422 
423  forAll< parallelDevicePolicy<> >( subRegion.size(), [bulkModulus,
424  shearModulus,
425  biotCoefficient,
426  stabilizationMultiplier,
427  macroElementIndex,
428  elementStabConstant] GEOS_HOST_DEVICE ( localIndex const ei )
429 
430  {
431  real64 const bM = bulkModulus[ei];
432  real64 const sM = shearModulus[ei];
433  real64 const bC = biotCoefficient[ei];
434 
435  macroElementIndex[ei] = 1;
436  elementStabConstant[ei] = stabilizationMultiplier * 9.0 * (bC * bC) / (32.0 * (10.0 * sM / 3.0 + bM));
437  } );
438  } );
439  } );
440 
441  }
442 
443  void updateBulkDensity( DomainPartition & domain )
444  {
445  this->template forDiscretizationOnMeshTargets<>( domain.getMeshBodies(), [&]( string const &,
446  MeshLevel & mesh,
447  string_array const & regionNames )
448  {
449  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
450  auto & subRegion )
451  {
452  // update the bulk density
453  // TODO: ideally, we would not recompute the bulk density, but a more general "rhs" containing the body force and the
454  // pressure/temperature terms
455  updateBulkDensity( subRegion );
456  } );
457  } );
458  }
459 
460 protected:
461 
463  {
464  Base::initializePostInitialConditionsPreSubGroups();
465 
466  string_array const & poromechanicsTargetRegionNames =
467  this->template getReference< string_array >( PhysicsSolverBase::viewKeyStruct::targetRegionsString() );
468  string_array const & solidMechanicsTargetRegionNames =
469  this->solidMechanicsSolver()->template getReference< string_array >( PhysicsSolverBase::viewKeyStruct::targetRegionsString() );
470  string_array const & flowTargetRegionNames =
471  this->flowSolver()->template getReference< string_array >( PhysicsSolverBase::viewKeyStruct::targetRegionsString() );
472  for( size_t i = 0; i < poromechanicsTargetRegionNames.size(); ++i )
473  {
474  GEOS_THROW_IF( std::find( solidMechanicsTargetRegionNames.begin(), solidMechanicsTargetRegionNames.end(),
475  poromechanicsTargetRegionNames[i] )
476  == solidMechanicsTargetRegionNames.end(),
477  GEOS_FMT( "Region {} must be a target region of {}",
478  poromechanicsTargetRegionNames[i], this->solidMechanicsSolver()->getName() ),
479  InputError, this->getDataContext(), this->solidMechanicsSolver()->getDataContext());
480  GEOS_THROW_IF( std::find( flowTargetRegionNames.begin(), flowTargetRegionNames.end(), poromechanicsTargetRegionNames[i] )
481  == flowTargetRegionNames.end(),
482  GEOS_FMT( "Region `{}` must be a target region of `{}`",
483  poromechanicsTargetRegionNames[i], this->flowSolver()->getCatalogName() ),
484  InputError, this->getDataContext(), this->flowSolver()->getDataContext() );
485  }
486  }
487 
488  template< typename TYPE_LIST,
489  typename KERNEL_WRAPPER,
490  typename ... PARAMS >
491  real64 assemblyLaunch( MeshLevel & mesh,
492  DofManager const & dofManager,
493  string_array const & regionNames,
494  string const & materialNamesString,
495  CRSMatrixView< real64, globalIndex const > const & localMatrix,
496  arrayView1d< real64 > const & localRhs,
497  real64 const dt,
498  PARAMS && ... params )
499  {
501 
502  NodeManager const & nodeManager = mesh.getNodeManager();
503 
504  string const dofKey = dofManager.getKey( fields::solidMechanics::totalDisplacement::key() );
505  arrayView1d< globalIndex const > const & dofNumber = nodeManager.getReference< globalIndex_array >( dofKey );
506 
507  real64 const gravityVectorData[3] = LVARRAY_TENSOROPS_INIT_LOCAL_3( this->gravityVector() );
508 
509  KERNEL_WRAPPER kernelWrapper( dofNumber,
510  dofManager.rankOffset(),
511  localMatrix,
512  localRhs,
513  dt,
514  gravityVectorData,
515  std::forward< PARAMS >( params )... );
516 
517  return finiteElement::
518  regionBasedKernelApplication< parallelDevicePolicy< >,
519  TYPE_LIST >( mesh,
520  regionNames,
521  this->solidMechanicsSolver()->getDiscretizationName(),
522  materialNamesString,
523  kernelWrapper );
524  }
525 
526  /* Implementation of Nonlinear Acceleration (Aitken) of averageMeanTotalStressIncrement */
527 
528  void recordAverageMeanTotalStressIncrement( DomainPartition & domain,
529  array1d< real64 > & averageMeanTotalStressIncrement )
530  {
531  averageMeanTotalStressIncrement.resize( 0 );
532  this->forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
533  MeshLevel & mesh,
534  string_array const & regionNames )
535  {
536  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
537  auto & subRegion )
538  {
539  // get the solid model (to access stress increment)
540  string const & solidName =
541  subRegion.template getReference< string >( viewKeyStruct::porousMaterialNamesString() );
542  constitutive::CoupledSolidBase & solid =
543  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
544 
545  arrayView1d< const real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
546  for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
547  {
548  averageMeanTotalStressIncrement.emplace_back( averageMeanTotalStressIncrement_k[k] );
549  }
550  } );
551  } );
552  }
553 
554  void applyAcceleratedAverageMeanTotalStressIncrement( DomainPartition & domain,
555  array1d< real64 > & averageMeanTotalStressIncrement )
556  {
557  integer i = 0;
558  this->forDiscretizationOnMeshTargets( domain.getMeshBodies(), [&]( string const &,
559  MeshLevel & mesh,
560  string_array const & regionNames )
561  {
562  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
563  auto & subRegion )
564  {
565  // get the solid model (to access stress increment)
566  string const & solidName = subRegion.template getReference< string >( viewKeyStruct::porousMaterialNamesString() );
567  constitutive::CoupledSolidBase & solid =
568  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
569  auto & porosityModel = dynamic_cast< constitutive::BiotPorosity const & >( solid.getBasePorosityModel());
570  arrayView1d< real64 > const & averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
571  for( localIndex k = 0; k < localIndex( averageMeanTotalStressIncrement_k.size()); k++ )
572  {
573  porosityModel.updateAverageMeanTotalStressIncrement( k, averageMeanTotalStressIncrement[i] );
574  i++;
575  }
576  } );
577  } );
578  }
579 
580  real64 computeAitkenRelaxationFactor( array1d< real64 > const & s0,
581  array1d< real64 > const & s1,
582  array1d< real64 > const & s1_tilde,
583  array1d< real64 > const & s2_tilde,
584  real64 const omega0 )
585  {
586  array1d< real64 > r1 = axpy( s1_tilde, s0, -1.0 );
587  array1d< real64 > r2 = axpy( s2_tilde, s1, -1.0 );
588 
589  // diff = r2 - r1
590  array1d< real64 > diff = axpy( r2, r1, -1.0 );
591 
592  real64 const denom = dot( diff, diff );
593  real64 const numer = dot( r1, diff );
594 
595  real64 omega1 = 1.0;
596  if( !isZero( denom ))
597  {
598  omega1 = -1.0 * omega0 * numer / denom;
599  }
600  return omega1;
601  }
602 
603  array1d< real64 > computeUpdate( array1d< real64 > const & s1,
604  array1d< real64 > const & s2_tilde,
605  real64 const omega1 )
606  {
607  return axpy( scale( s1, 1.0 - omega1 ),
608  scale( s2_tilde, omega1 ),
609  1.0 );
610  }
611 
612  void startSequentialIteration( integer const & iter,
613  DomainPartition & domain ) override
614  {
615  if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
616  {
617  if( iter == 0 )
618  {
619  recordAverageMeanTotalStressIncrement( domain, m_s1 );
620  }
621  else
622  {
623  m_s0 = m_s1;
624  m_s1 = m_s2;
625  m_s1_tilde = m_s2_tilde;
626  m_omega0 = m_omega1;
627  }
628  }
629  }
630 
631  void finishSequentialIteration( integer const & iter,
632  DomainPartition & domain ) override
633  {
634  if( this->getNonlinearSolverParameters().m_nonlinearAccelerationType == NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
635  {
636  if( iter == 0 )
637  {
638  m_s2 = m_s2_tilde;
639  m_omega1 = 1.0;
640  }
641  else
642  {
643  m_omega1 = computeAitkenRelaxationFactor( m_s0, m_s1, m_s1_tilde, m_s2_tilde, m_omega0 );
644  m_s2 = computeUpdate( m_s1, m_s2_tilde, m_omega1 );
645  applyAcceleratedAverageMeanTotalStressIncrement( domain, m_s2 );
646  }
647  }
648  }
649 
650  virtual void mapSolutionBetweenSolvers( DomainPartition & domain, integer const solverType ) override
651  {
653 
655  if( solverType == static_cast< integer >( SolverType::Flow ) )
656  {
657  updateBulkDensity( domain );
658  }
659 
661  if( solverType == static_cast< integer >( SolverType::SolidMechanics )
662  && !m_performStressInitialization ) // do not update during poromechanics initialization
663  {
664  // compute the average of the mean total stress increment over quadrature points
665  averageMeanTotalStressIncrement( domain );
666 
667  this->template forDiscretizationOnMeshTargets<>( domain.getMeshBodies(), [&]( string const &,
668  MeshLevel & mesh,
669  string_array const & regionNames )
670  {
671 
672 
673  mesh.getElemManager().forElementSubRegions< CellElementSubRegion, SurfaceElementSubRegion >( regionNames, [&]( localIndex const,
674  auto & subRegion )
675  {
676  // update the porosity after a change in displacement (after mechanics solve)
677  // or a change in pressure/temperature (after a flow solve)
678  flowSolver()->updatePorosityAndPermeability( subRegion );
679  } );
680  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
681  auto & subRegion )
682 
683  {
684  // update bulk density to reflect porosity change into mechanics
685  updateBulkDensity( subRegion );
686  } );
687  } );
688  }
689 
690  // needed to perform nonlinear acceleration
691  if( solverType == static_cast< integer >( SolverType::SolidMechanics ) &&
692  this->getNonlinearSolverParameters().m_nonlinearAccelerationType== NonlinearSolverParameters::NonlinearAccelerationType::Aitken )
693  {
694  recordAverageMeanTotalStressIncrement( domain, m_s2_tilde );
695  }
696  }
697 
703  {
704  this->template forDiscretizationOnMeshTargets<>( domain.getMeshBodies(), [&]( string const &,
705  MeshLevel & mesh,
706  string_array const & regionNames )
707  {
708  mesh.getElemManager().forElementSubRegions< CellElementSubRegion >( regionNames, [&]( localIndex const,
709  auto & subRegion )
710  {
711  // get the solid model (to access stress increment)
712  string const & solidName = subRegion.template getReference< string >( viewKeyStruct::porousMaterialNamesString() );
713  constitutive::CoupledSolidBase & solid =
714  this->template getConstitutiveModel< constitutive::CoupledSolidBase >( subRegion, solidName );
715 
716  arrayView2d< real64 const > const meanTotalStressIncrement_k = solid.getMeanTotalStressIncrement_k();
717  arrayView1d< real64 > const averageMeanTotalStressIncrement_k = solid.getAverageMeanTotalStressIncrement_k();
718 
719  finiteElement::FiniteElementBase & subRegionFE =
720  subRegion.template getReference< finiteElement::FiniteElementBase >( solidMechanicsSolver()->getDiscretizationName() );
721 
722  // determine the finite element type
723  finiteElement::FiniteElementDispatchHandler< BASE_FE_TYPES >::
724  dispatch3D( subRegionFE, [&] ( auto const finiteElement )
725  {
726  using FE_TYPE = decltype( finiteElement );
727 
728  // call the factory and launch the kernel
729  AverageOverQuadraturePoints1DKernelFactory::
730  createAndLaunch< CellElementSubRegion,
731  FE_TYPE,
732  parallelDevicePolicy<> >( mesh.getNodeManager(),
733  mesh.getEdgeManager(),
734  mesh.getFaceManager(),
735  subRegion,
736  finiteElement,
737  meanTotalStressIncrement_k,
738  averageMeanTotalStressIncrement_k );
739  } );
740  } );
741  } );
742  }
743 
744  virtual void setMGRStrategy() = 0;
745 
746  virtual void updateBulkDensity( ElementSubRegionBase & subRegion ) = 0;
747 
748  virtual string getFlowDofKey() const = 0;
749 
750  virtual void validateNonlinearAcceleration() override
751  {
752  if( MpiWrapper::commSize( MPI_COMM_GEOS ) > 1 )
753  {
754  GEOS_ERROR( "Nonlinear acceleration is not implemented for MPI runs" );
755  }
756  }
757 
760 
763 
765  stabilization::StabilizationType m_stabilizationType;
766 
769 
772 
774  array1d< real64 > m_s0; // Accelerated averageMeanTotalStresIncrement @ outer iteration v ( two iterations ago )
775  array1d< real64 > m_s1; // Accelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( previous iteration )
776  array1d< real64 > m_s1_tilde; // Unaccelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( previous iteration )
777  array1d< real64 > m_s2; // Accelerated averageMeanTotalStresIncrement @ outer iteration v + 2 ( current iteration )
778  array1d< real64 > m_s2_tilde; // Unaccelerated averageMeanTotalStresIncrement @ outer iteration v + 1 ( current iteration )
779  real64 m_omega0; // Old Aitken relaxation factor
780  real64 m_omega1; // New Aitken relaxation factor
781 
782 };
783 
784 } /* namespace geos */
785 
786 #endif //GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSSOLVER_HPP_
#define GEOS_UNUSED_PARAM(X)
Mark an unused argument and silence compiler warnings.
Definition: GeosxMacros.hpp:72
#define GEOS_ERROR(...)
Raise a hard error and terminate the program.
Definition: Logger.hpp:204
#define GEOS_THROW_IF(COND, MSG,...)
Conditionally raise a hard error and terminate the program.
Definition: Logger.hpp:263
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
virtual void initializePreSubGroups() override
Called by Initialize() prior to initializing sub-Groups.
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
globalIndex rankOffset(string const &fieldName) const
void addCoupling(string const &rowFieldName, string const &colFieldName, Connector connectivity, stdVector< FieldSupport > const &regions={}, bool symmetric=true)
Add coupling between two fields.
@ Elem
connectivity is element (like in finite elements)
string const & getKey(string const &fieldName) const
Return the key used to record the field in the DofManager.
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
NodeManager const & getNodeManager() const
Get the node manager.
Definition: MeshLevel.hpp:155
ElementRegionManager const & getElemManager() const
Get the element region manager.
Definition: MeshLevel.hpp:207
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
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 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.
virtual void assembleSystem(real64 const time, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
function to assemble the linear system matrix and rhs
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.
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
string const & getName() const
Get group name.
Definition: Group.hpp:1329
GEOS_DECLTYPE_AUTO_RETURN getReference(LOOKUP_TYPE const &lookup) const
Look up a wrapper and get reference to wrapped object.
Definition: Group.hpp:1273
DataContext const & getWrapperDataContext(KEY key) const
Definition: Group.hpp:1354
@ NOPLOT
Do not ever write to plot file.
@ OPTIONAL
Optional in input.
@ NO_WRITE
Do not write into restart.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:179
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
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:309
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.
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:1446