GEOS
CompositionalMultiphaseReservoirFVM.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 Total, S.A
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 
20 #ifndef GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRCOMPOSITIONALMULTIPHASERESERVOIRFVM_HPP_
21 #define GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRCOMPOSITIONALMULTIPHASERESERVOIRFVM_HPP_
22 
24 
25 namespace geos
26 {
27 
28 namespace hypre
29 {
30 
31 namespace mgr
32 {
33 
55 {
56 public:
57 
62  explicit CompositionalMultiphaseReservoirFVM( arrayView1d< int const > const & numComponentsPerField )
63  : MGRStrategyBase( LvArray::integerConversion< HYPRE_Int >( numComponentsPerField[0] + numComponentsPerField[1] ) )
64  {
65  HYPRE_Int const numResLabels = LvArray::integerConversion< HYPRE_Int >( numComponentsPerField[0] );
66 
67  // Level 0: eliminate the well block
68  m_labels[0].resize( numResLabels );
69  std::iota( m_labels[0].begin(), m_labels[0].end(), 0 );
70  // Level 1: eliminate the last density of the reservoir block
71  m_labels[1].resize( numResLabels - 1 );
72  std::iota( m_labels[1].begin(), m_labels[1].end(), 0 );
73  // Level 2: eliminate the rest of the densities
74  m_labels[2].push_back( 0 );
75 
76  setupLabels();
77 
78  // level 0
80  m_levelFRelaxIters[0] = 1;
85 
86  // level 1
88  m_levelFRelaxIters[1] = 1;
93 
94  // level 2
101  }
102 
109  void setup( LinearSolverParameters::MGR const & mgrParams,
110  HyprePrecWrapper & precond,
111  HypreMGRData & mgrData )
112  {
113  // if the wells are shut, using Gaussian elimination as F-relaxation for the well block is an overkill
114  // in that case, we just use Jacobi
115  if( mgrParams.areWellsShut )
116  {
118  }
119 
120  setReduction( precond, mgrData );
121 
122  // Configure the BoomerAMG solver used as mgr coarse solver for the pressure reduced system
123  setPressureAMG( mgrData.coarseSolver );
124  }
125 };
126 
127 } // namespace mgr
128 
129 } // namespace hypre
130 
131 } // namespace geos
132 
133 #endif /*GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRCOMPOSITIONALMULTIPHASERESERVOIRFVM_HPP_*/
CompositionalMultiphaseReservoirFVM(arrayView1d< int const > const &numComponentsPerField)
Constructor.
void setup(LinearSolverParameters::MGR const &mgrParams, HyprePrecWrapper &precond, HypreMGRData &mgrData)
Setup the MGR strategy.
Helper struct for strategies that provides some basic parameter arrays needed by MGR.
Definition: HypreMGR.hpp:74
MGRCoarseGridMethod m_levelCoarseGridMethod[numLevels]
Coarse grid method for each level.
Definition: HypreMGR.hpp:91
MGRInterpolationType m_levelInterpType[numLevels]
Interpolation type for each level.
Definition: HypreMGR.hpp:89
std::vector< HYPRE_Int > m_labels[numLevels]
Dof labels kept at each level.
Definition: HypreMGR.hpp:83
void setupLabels()
Call this after populating lv_cindexes.
Definition: HypreMGR.hpp:121
HYPRE_Int m_levelFRelaxIters[numLevels]
Number of F-relaxation iterations for each level.
Definition: HypreMGR.hpp:88
void setPressureAMG(HyprePrecWrapper &solver)
Set up BoomerAMG to perform the solve for the pressure system.
Definition: HypreMGR.hpp:204
void setReduction(HyprePrecWrapper &precond, HypreMGRData &mgrData)
Helper function that sets the reduction features common to all mgr strategies.
Definition: HypreMGR.hpp:135
MGRFRelaxationType m_levelFRelaxType[numLevels]
F-relaxation type for each level.
Definition: HypreMGR.hpp:87
HYPRE_Int m_levelGlobalSmootherIters[numLevels]
Number of global smoother iterations for each level.
Definition: HypreMGR.hpp:93
MGRGlobalSmootherType m_levelGlobalSmootherType[numLevels]
Global smoother type for each level.
Definition: HypreMGR.hpp:92
MGRRestrictionType m_levelRestrictType[numLevels]
Restriction type for each level.
Definition: HypreMGR.hpp:90
@ 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)
@ ilu0
incomplete LU factorization
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
Container for hypre preconditioner auxiliary data for MGR.
Definition: HypreMGR.hpp:38
HyprePrecWrapper coarseSolver
MGR coarse solver pointer and functions.
Definition: HypreMGR.hpp:40
Container for hypre preconditioner function pointers.
Definition: HypreUtils.hpp:54
Multigrid reduction parameters.