GEOSX
SinglePhasePoromechanics.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2018-2020 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2020 The Board of Trustees of the Leland Stanford Junior University
7  * Copyright (c) 2018-2020 TotalEnergies
8  * Copyright (c) 2019- GEOSX Contributors
9  * All rights reserved
10  *
11  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
12  * ------------------------------------------------------------------------------------------------------------
13  */
14 
19 #ifndef GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRSINGLEPHASEPOROMECHANICS_HPP_
20 #define GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRSINGLEPHASEPOROMECHANICS_HPP_
21 
23 
24 namespace geos
25 {
26 
27 namespace hypre
28 {
29 
30 namespace mgr
31 {
32 
48 {
49 public:
50 
55  : MGRStrategyBase( 4 )
56  {
57  m_labels[0].push_back( 3 );
58 
59  setupLabels();
60 
61  // Level 0
63  m_levelFRelaxIters[0] = 1;
68  }
69 
76  HyprePrecWrapper & precond,
77  HypreMGRData & mgrData )
78  {
79  setReduction( precond, mgrData );
80 
81  GEOS_LAI_CHECK_ERROR( HYPRE_MGRSetPMaxElmts( precond.ptr, 0 ));
82 
83  // Configure the BoomerAMG solver used as F-relaxation for the first level
84  setMechanicsFSolver( precond, mgrData );
85 
86  // Configure the BoomerAMG solver used as mgr coarse solver for the pressure reduced system
87  setPressureAMG( mgrData.coarseSolver );
88  }
89 };
90 
91 } // namespace mgr
92 
93 } // namespace hypre
94 
95 } // namespace geos
96 
97 #endif /*GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRSINGLEPHASEPOROMECHANICS_HPP_*/
Helper struct for strategies that provides some basic parameter arrays needed by MGR.
Definition: HypreMGR.hpp:73
MGRCoarseGridMethod m_levelCoarseGridMethod[numLevels]
Coarse grid method for each level.
Definition: HypreMGR.hpp:90
MGRInterpolationType m_levelInterpType[numLevels]
Interpolation type for each level.
Definition: HypreMGR.hpp:88
std::vector< HYPRE_Int > m_labels[numLevels]
Dof labels kept at each level.
Definition: HypreMGR.hpp:82
void setupLabels()
Call this after populating lv_cindexes.
Definition: HypreMGR.hpp:119
HYPRE_Int m_levelFRelaxIters[numLevels]
Number of F-relaxation iterations for each level.
Definition: HypreMGR.hpp:87
void setPressureAMG(HyprePrecWrapper &solver)
Set up BoomerAMG to perform the solve for the pressure system.
Definition: HypreMGR.hpp:197
void setReduction(HyprePrecWrapper &precond, HypreMGRData &mgrData)
Helper function that sets the reduction features common to all mgr strategies.
Definition: HypreMGR.hpp:133
MGRFRelaxationType m_levelFRelaxType[numLevels]
F-relaxation type for each level.
Definition: HypreMGR.hpp:86
MGRGlobalSmootherType m_levelGlobalSmootherType[numLevels]
Global smoother type for each level.
Definition: HypreMGR.hpp:91
MGRRestrictionType m_levelRestrictType[numLevels]
Restriction type for each level.
Definition: HypreMGR.hpp:89
void setMechanicsFSolver(HyprePrecWrapper &precond, HypreMGRData &mgrData)
Set up BoomerAMG to perform the mechanics F-solve for the first F-relaxation.
Definition: HypreMGR.hpp:256
void setup(LinearSolverParameters::MGR const &, HyprePrecWrapper &precond, HypreMGRData &mgrData)
Setup the MGR strategy.
SinglePhasePoromechanics(arrayView1d< int const > const &)
Constructor.
#define GEOS_LAI_CHECK_ERROR(call)
Definition: common.hpp:117
@ amgVCycle
Full AMG VCycle solver.
@ none
no global smoothing is performed (default)
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:220
Container for hypre preconditioner auxiliary data for MGR.
Definition: HypreMGR.hpp:37
HyprePrecWrapper coarseSolver
MGR coarse solver pointer and functions.
Definition: HypreMGR.hpp:39
Container for hypre preconditioner function pointers.
Definition: HypreUtils.hpp:53
HYPRE_Solver ptr
pointer to preconditioner
Definition: HypreUtils.hpp:63
Multigrid reduction parameters.