GEOS
ThermalCompositionalMultiphaseFVM.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_HYPREMGRTHERMALCOMPOSITIONALMULTIPHASEFVM_HPP_
21 #define GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRTHERMALCOMPOSITIONALMULTIPHASEFVM_HPP_
22 
24 
25 namespace geos
26 {
27 
28 namespace hypre
29 {
30 
31 namespace mgr
32 {
33 
51 {
52 public:
57  explicit ThermalCompositionalMultiphaseFVM( arrayView1d< int const > const & numComponentsPerField )
58  : MGRStrategyBase( LvArray::integerConversion< HYPRE_Int >( numComponentsPerField[0] ) )
59  {
60  // Level 0: eliminate last density which corresponds to the volume constraint equation
61  m_labels[0].resize( m_numBlocks - 2 );
62  std::iota( m_labels[0].begin(), m_labels[0].end(), 0 );
63  m_labels[0].push_back( m_numBlocks-1 ); // keep temperature
64  // Level 1: eliminate the other densities
65  m_labels[1].push_back( 0 ); // keep pressure
66  m_labels[1].push_back( m_numBlocks-1 ); // keep temperature
67 
68  setupLabels();
69 
71  m_levelFRelaxIters[0] = 1;
72  m_levelInterpType[0] = MGRInterpolationType::jacobi; // Diagonal scaling (Jacobi)
77 
81  m_levelCoarseGridMethod[1] = MGRCoarseGridMethod::cprLikeBlockDiag; // Non-Galerkin Quasi-IMPES CPR
84  }
85 
92  HyprePrecWrapper & precond,
93  HypreMGRData & mgrData )
94  {
95  setReduction( precond, mgrData );
96 
97  // Configure the BoomerAMG solver used as mgr coarse solver for the pressure/temperature reduced system
99  }
100 };
101 
102 } // namespace mgr
103 
104 } // namespace hypre
105 
106 } // namespace geos
107 
108 #endif /*GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRCOMPOSITIONALMULTIPHASEFVM_HPP_*/
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
HYPRE_Int m_numBlocks
Number of different matrix blocks treated separately.
Definition: HypreMGR.hpp:81
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 setPressureTemperatureAMG(HyprePrecWrapper &solver)
Set up BoomerAMG to perform the solve for the pressure/temperature system.
Definition: HypreMGR.hpp:232
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 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
ThermalCompositionalMultiphaseFVM(arrayView1d< int const > const &numComponentsPerField)
Constructor.
void setup(LinearSolverParameters::MGR const &, HyprePrecWrapper &precond, HypreMGRData &mgrData)
Setup the MGR strategy.
@ galerkin
Galerkin coarse grid computation using RAP.
@ none
no F-relaxation if performed
@ 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.