GEOSX
SinglePhaseReservoirFVM.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_HYPREMGRSINGLEPHASERESERVOIRFVM_HPP_
20 #define GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRSINGLEPHASERESERVOIRFVM_HPP_
21 
23 
24 namespace geos
25 {
26 
27 namespace hypre
28 {
29 
30 namespace mgr
31 {
32 
48 {
49 public:
54  : MGRStrategyBase( LvArray::integerConversion< HYPRE_Int >( 3 ) )
55  {
56  // Level 0: eliminate the well variables, and just keep the cell-centered pressures
57  m_labels[0].push_back( 0 );
58 
59  setupLabels();
60 
61  // Level 0
63  m_levelFRelaxIters[0] = 1;
68  }
69 
76  void setup( LinearSolverParameters::MGR const & mgrParams,
77  HyprePrecWrapper & precond,
78  HypreMGRData & mgrData )
79  {
80  // if the wells are shut, using Gaussian elimination as F-relaxation for the well block is an overkill
81  // in that case, we just use Jacobi
82  if( mgrParams.areWellsShut )
83  {
85  }
86 
87  setReduction( precond, mgrData );
88 
89  // Configure the BoomerAMG solver used as mgr coarse solver for the pressure reduced system
90  setPressureAMG( mgrData.coarseSolver );
91  }
92 };
93 
94 } // namespace mgr
95 
96 } // namespace hypre
97 
98 } // namespace geos
99 
100 #endif /*GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRSINGLEPHASERESERVOIRFVM_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
SinglePhaseReservoirFVM(arrayView1d< int const > const &)
Constructor.
void setup(LinearSolverParameters::MGR const &mgrParams, HyprePrecWrapper &precond, HypreMGRData &mgrData)
Setup the MGR strategy.
@ galerkin
Galerkin coarse grid computation using RAP.
@ gsElimWInverse
Direct Inversion with Gaussian Elimination (OK for larger systems)
@ 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
Multigrid reduction parameters.