GEOS
MsrsbLevelBuilderCoupled.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2018-2019 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2019 The Board of Trustees of the Leland Stanford Junior University
7  * Copyright (c) 2018-2019 Total, S.A
8  * Copyright (c) 2019- GEOS/GEOSX Contributors
9  * All right reserved
10  *
11  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
12  * ------------------------------------------------------------------------------------------------------------
13  */
14 
19 #ifndef GEOS_LINEARALGEBRA_MULTISCALE_MSRSBLEVELBUILDERCOUPLED_HPP
20 #define GEOS_LINEARALGEBRA_MULTISCALE_MSRSBLEVELBUILDERCOUPLED_HPP
21 
22 #include "MsrsbLevelBuilder.hpp"
24 
25 namespace geos
26 {
27 namespace multiscale
28 {
29 
34 template< typename LAI >
36 {
37 public:
38 
41 
43  using Vector = typename Base::Vector;
44 
46  using Matrix = typename Base::Matrix;
47 
53  explicit MsrsbLevelBuilderCoupled( string name, LinearSolverParameters params );
54 
61  virtual void initializeFineLevel( DomainPartition & domain,
63  MPI_Comm const & comm ) override;
64 
71  Matrix const & fineMatrix ) override;
72 
78  virtual bool updateProlongation( Matrix const & fineMatrix ) override;
79 
80  virtual std::unique_ptr< PreconditionerBase< LAI > > makeCoarseSolver() const override;
81 
82 private:
83 
84  using Base::m_params;
85  using Base::m_name;
87  using Base::m_restriction;
88  using Base::m_matrix;
89  using Base::m_dofManager;
90  using Base::m_preSmoother;
92  using Base::m_fineLevel;
93 
94  void initializeCommon( DomainPartition & domain, MPI_Comm const & comm );
95 
96  void createSmoothers();
97 
98  void buildProlongationStructure( DofManager const & fineDofManager );
99 
101  std::vector< string > m_fields;
102 
104  std::vector< Matrix > m_selectors;
105 
107  std::vector< std::unique_ptr< MsrsbLevelBuilder< LAI > > > m_builders;
108 
110  std::vector< CRSMatrix< real64, globalIndex > > m_prolongationBlocks;
111 
113  CRSMatrix< real64, globalIndex > m_localProlongation;
114 };
115 
116 } // namespace multiscale
117 
118 } // namespace geos
119 
120 #endif //GEOS_LINEARALGEBRA_MULTISCALE_MSRSBLEVELBUILDERCOUPLED_HPP
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Definition: DofManager.hpp:45
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
Degree-of-freedom manager that works with multiscale mesh levels.
Definition: DofManager.hpp:37
Base class for level builder implementations.
typename LAI::ParallelVector Vector
Alias for vector type.
Base class for MsRSB level builders.
multiscale::DofManager m_dofManager
DofManager for the matrix.
std::unique_ptr< Operator > m_restriction
Restriction (kept as abstract operator to allow for memory efficiency, e.g. when R = P^T)
LinearSolverParameters m_params
Linear solver top-level parameters.
MsrsbLevelBuilderBase const * m_fineLevel
Pointer to the fine level.
std::unique_ptr< PreconditionerBase< LAI > > m_postSmoother
Post-smoothing operator.
typename Base::Matrix Matrix
Alias for matrix type.
std::unique_ptr< PreconditionerBase< LAI > > m_preSmoother
Pre-smoothing operator.
typename Base::Vector Vector
Alias for vector type.
MsRSB level builder for coupled problems.
virtual void initializeFineLevel(DomainPartition &domain, geos::DofManager const &dofManager, MPI_Comm const &comm) override
Initialize the finest level (level 0).
typename Base::Matrix Matrix
Alias for matrix type.
virtual std::unique_ptr< PreconditionerBase< LAI > > makeCoarseSolver() const override
Instantiate coarsest level solver.
MsrsbLevelBuilderCoupled(string name, LinearSolverParameters params)
Constructor.
virtual bool updateProlongation(Matrix const &fineMatrix) override
Update current level's prolongation using a new previous level matrix.
virtual void initializeCoarseLevel(LevelBuilderBase< LAI > &fineLevel, Matrix const &fineMatrix) override
Initialize a coarse level (levels 1 and above).
LvArray::CRSMatrix< T, COL_INDEX, INDEX_TYPE, LvArray::ChaiBuffer > CRSMatrix
Alias for CRS Matrix class.
Definition: DataTypes.hpp:305
Set of parameters for a linear solver or preconditioner.