GEOSX
CompositionalMultiphaseHybridFVM.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_HYPREMGRCOMPOSITIONALMULTIPHASEHYBRIDFVM_HPP_
20 #define GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRCOMPOSITIONALMULTIPHASEHYBRIDFVM_HPP_
21 
23 
24 namespace geos
25 {
26 
27 namespace hypre
28 {
29 
30 namespace mgr
31 {
32 
50 {
51 public:
56  explicit CompositionalMultiphaseHybridFVM( arrayView1d< int const > const & numComponentsPerField )
57  : MGRStrategyBase( LvArray::integerConversion< HYPRE_Int >( numComponentsPerField[0] + numComponentsPerField[1] ) )
58  {
59  // Level 0: eliminate the last density of the cell-centered block
60  m_labels[0].resize( m_numBlocks - 1 );
61  std::iota( m_labels[0].begin(), m_labels[0].begin() + m_numBlocks - 2, 0 );
62  m_labels[0][m_numBlocks - 2] = m_numBlocks - 1;
63  // Level 1: eliminate remaining densities of the cell-centered block
64  m_labels[1].resize( 2 );
65  m_labels[1][0] = 0;
66  m_labels[1][1] = m_numBlocks - 1;
67  // Level 2: eliminate reservoir cell-centered pressure
68  m_labels[2].resize( 1 );
69  m_labels[2][0] = m_numBlocks - 1;
70 
71  setupLabels();
72 
73  // Level 0
75  m_levelFRelaxIters[0] = 1;
80 
81  // Level 1
87 
88  // Level 2
90  m_levelFRelaxIters[2] = 1;
96  }
97 
104  HyprePrecWrapper & precond,
105  HypreMGRData & mgrData )
106  {
107  setReduction( precond, mgrData );
108 
109  // Configure the BoomerAMG solver used as mgr coarse solver for the pressure reduced system
110  setPressureAMG( mgrData.coarseSolver );
111  }
112 };
113 
114 } // namespace mgr
115 
116 } // namespace hypre
117 
118 } // namespace geos
119 
120 #endif /*GEOS_LINEARALGEBRA_INTERFACES_HYPREMGRCOMPOSITIONALMULTIPHASEHYBRIDFVM_HPP_*/
Labels description stored in point_marker_array 0 = pressure 1 = density ... = ......
void setup(LinearSolverParameters::MGR const &, HyprePrecWrapper &precond, HypreMGRData &mgrData)
Setup the MGR strategy.
CompositionalMultiphaseHybridFVM(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
HYPRE_Int m_numBlocks
Number of different matrix blocks treated separately.
Definition: HypreMGR.hpp:80
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
@ 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.