GEOSX
CompositionalMultiphaseReservoirHybridFVM.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_HYPREMGRCOMPOSITIONALMULTIPHASERESERVOIRHYBRIDFVM_HPP_
20 #define GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRCOMPOSITIONALMULTIPHASERESERVOIRHYBRIDFVM_HPP_
21 
23 
24 namespace geos
25 {
26 
27 namespace hypre
28 {
29 
30 namespace mgr
31 {
32 
55 {
56 public:
61  explicit CompositionalMultiphaseReservoirHybridFVM( arrayView1d< int const > const & numComponentsPerField )
62  : MGRStrategyBase( LvArray::integerConversion< HYPRE_Int >( numComponentsPerField[0] + numComponentsPerField[1] + numComponentsPerField[2] ) )
63  {
64  HYPRE_Int const numResCellCenteredLabels = LvArray::integerConversion< HYPRE_Int >( numComponentsPerField[0] );
65  HYPRE_Int const numResFaceCenteredLabels = LvArray::integerConversion< HYPRE_Int >( numComponentsPerField[1] );
66  HYPRE_Int const numResLabels = numResCellCenteredLabels + numResFaceCenteredLabels;
67 
68  // Level 0: eliminate the well block
69  m_labels[0].resize( numResLabels );
70  std::iota( m_labels[0].begin(), m_labels[0].end(), 0 );
71  // Level 1: eliminate the last density of the reservoir block
72  m_labels[1].resize( numResLabels - 1 );
73  std::iota( m_labels[1].begin(), m_labels[1].begin() + numResCellCenteredLabels - 1, 0 );
74  m_labels[1][numResCellCenteredLabels-1] = numResCellCenteredLabels;
75  // Level 2: eliminate remaining densities of the reservoir block
76  m_labels[2].resize( 2 );
77  m_labels[2][0] = 0;
78  m_labels[2][1] = numResCellCenteredLabels;
79  // Level 3: eliminate reservoir cell centered pressure
80  m_labels[3].resize( 1 );
81  m_labels[3][0] = numResCellCenteredLabels;
82 
83  setupLabels();
84 
85  // Level 0
87  m_levelFRelaxIters[0] = 1;
92 
93  // Level 1
95  m_levelFRelaxIters[1] = 1;
100 
101  // Level 2
107 
108  // Level 3
110  m_levelFRelaxIters[3] = 1;
116  }
117 
124  void setup( LinearSolverParameters::MGR const & mgrParams,
125  HyprePrecWrapper & precond,
126  HypreMGRData & mgrData )
127  {
128  // if the wells are shut, using Gaussian elimination as F-relaxation for the well block is an overkill
129  // in that case, we just use Jacobi
130  if( mgrParams.areWellsShut )
131  {
133  }
134 
135  setReduction( precond, mgrData );
136 
137  // Configure the BoomerAMG solver used as mgr coarse solver for the pressure reduced system
138  setPressureAMG( mgrData.coarseSolver );
139  }
140 };
141 
142 } // namespace mgr
143 
144 } // namespace hypre
145 
146 } // namespace geos
147 
148 #endif /*GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRCOMPOSITIONALMULTIPHASERESERVOIRHYBRIDFVM_HPP_*/
void setup(LinearSolverParameters::MGR const &mgrParams, HyprePrecWrapper &precond, HypreMGRData &mgrData)
Setup the MGR strategy.
CompositionalMultiphaseReservoirHybridFVM(arrayView1d< int const > const &numComponentsPerField)
Constructor.
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
HYPRE_Int m_levelGlobalSmootherIters[numLevels]
Number of global smoother iterations for each level.
Definition: HypreMGR.hpp:92
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
@ galerkinRAI
Galerkin coarse grid computation with arbitrary classical restriction and injective prolongation.
@ galerkin
Galerkin coarse grid computation using RAP.
@ blockColLumped
Block column-lumped approximation.
@ none
no F-relaxation if performed
@ 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.