G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsGridHierarchy.h
Go to the documentation of this file.
1
14#pragma once
15
16#include <gsCore/gsDofMapper.h>
18#include <gsIO/gsOptionList.h>
19
20namespace gismo
21{
22
32template< typename T >
34{
35
36public:
37
50 gsMultiBasis<T> mBasis,
51 const gsBoundaryConditions<T>& boundaryConditions,
52 const gsOptionList& assemblerOptions,
53 index_t levels,
54 index_t numberOfKnotsToBeInserted = 1,
55 index_t multiplicityOfKnotsToBeInserted = 1,
56 index_t unk = 0
57 );
58
66 gsMultiBasis<T> mBasis,
67 const gsBoundaryConditions<T>& boundaryConditions,
68 const gsOptionList& options
69 )
70 {
72 give(mBasis),
73 boundaryConditions,
74 options,
75 options.askInt( "Levels", 3 ),
76 options.askInt( "NumberOfKnotsToBeInserted", 1 ),
77 options.askInt( "MultiplicityOfKnotsToBeInserted", 1 ),
78 0
79 );
80 }
81
82
98 gsMultiBasis<T> mBasis,
99 const gsBoundaryConditions<T>& boundaryConditions,
100 const gsOptionList& assemblerOptions,
101 index_t levels,
102 index_t degreesOfFreedom = 0,
103 index_t unk = 0
104 );
105
114 gsMultiBasis<T> mBasis,
115 const gsBoundaryConditions<T>& boundaryConditions,
116 const gsOptionList& options
117 )
118 {
120 give(mBasis),
121 boundaryConditions,
122 options,
123 options.askInt( "Levels", 3 ),
124 options.askInt( "DegreesOfFreedom", 0 ),
125 0
126 );
127 }
128
131 {
132 gsOptionList opt;
133 opt.addInt( "DirichletStrategy", "Method for enforcement of Dirichlet BCs [11..14]", 11 );
134 opt.addInt( "InterfaceStrategy", "Method of treatment of patch interfaces [0..3]", 1 );
135 opt.addInt( "Levels", "Number of levels to be constructed in the grid hierarchy", 3 );
136 opt.addInt( "DegreesOfFreedom", "Number of dofs in the coarsest grid in the grid hierarchy (only buildByCoarsening)", 0 );
137 opt.addInt( "NumberOfKnotsToBeInserted", "The number of knots to be inserted (only buildByRefinement)", 1 );
138 opt.addInt( "MultiplicityOfKnotsToBeInserted", "The multiplicity of the knots to be inserted (only buildByRefinement)", 1 );
139 return opt;
140 }
141
143 void clear()
144 {
145 m_mBases.clear();
146 m_transferMatrices.clear();
147 }
148
150 const std::vector< gsMultiBasis<T> >& getMultiBases() const
151 { return m_mBases; }
152
155 { o = give(m_mBases); return *this; }
156
158 const std::vector< gsSparseMatrix<T, RowMajor> >& getTransferMatrices() const
159 { return m_transferMatrices; }
160
163 { o = give(m_transferMatrices); return *this; }
164
165private:
166 std::vector< gsMultiBasis<T> > m_mBases;
167 std::vector< gsSparseMatrix<T, RowMajor> > m_transferMatrices;
168};
169
170} // namespace gismo
171
172#ifndef GISMO_BUILD_LIB
173#include GISMO_HPP_HEADER(gsGridHierarchy.hpp)
174#endif
Class containing a set of boundary conditions.
Definition gsBoundaryConditions.h:342
Grid Hierarchy.
Definition gsGridHierarchy.h:34
static gsGridHierarchy buildByRefinement(gsMultiBasis< T > mBasis, const gsBoundaryConditions< T > &boundaryConditions, const gsOptionList &assemblerOptions, index_t levels, index_t numberOfKnotsToBeInserted=1, index_t multiplicityOfKnotsToBeInserted=1, index_t unk=0)
This function sets up a multigrid hierarchy by uniform refinement.
Definition gsGridHierarchy.hpp:27
gsGridHierarchy & moveMultiBasesTo(std::vector< gsMultiBasis< T > > &o)
Get the vector of gsMultiBasis objects.
Definition gsGridHierarchy.h:154
const std::vector< gsMultiBasis< T > > & getMultiBases() const
Get the vector of gsMultiBasis objects (by reference)
Definition gsGridHierarchy.h:150
static gsGridHierarchy buildByCoarsening(gsMultiBasis< T > mBasis, const gsBoundaryConditions< T > &boundaryConditions, const gsOptionList &options)
This function sets up a grid hierarchy by coarsening.
Definition gsGridHierarchy.h:113
static gsGridHierarchy buildByRefinement(gsMultiBasis< T > mBasis, const gsBoundaryConditions< T > &boundaryConditions, const gsOptionList &options)
This function sets up a multigrid hierarchy by uniform refinement.
Definition gsGridHierarchy.h:65
const std::vector< gsSparseMatrix< T, RowMajor > > & getTransferMatrices() const
Get the vector of transfer matrices (by reference)
Definition gsGridHierarchy.h:158
gsGridHierarchy & moveTransferMatricesTo(std::vector< gsSparseMatrix< T, RowMajor > > &o)
Get the vector of transfer matrices.
Definition gsGridHierarchy.h:162
void clear()
Reset the object (to save memory)
Definition gsGridHierarchy.h:143
static gsGridHierarchy buildByCoarsening(gsMultiBasis< T > mBasis, const gsBoundaryConditions< T > &boundaryConditions, const gsOptionList &assemblerOptions, index_t levels, index_t degreesOfFreedom=0, index_t unk=0)
This function sets up a grid hierarchy by coarsening.
Definition gsGridHierarchy.hpp:57
static gsOptionList defaultOptions()
Get the default options.
Definition gsGridHierarchy.h:130
Holds a set of patch-wise bases and their topology information.
Definition gsMultiBasis.h:37
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
void addInt(const std::string &label, const std::string &desc, const index_t &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:201
index_t askInt(const std::string &label, const index_t &value=0) const
Reads value for option label from options.
Definition gsOptionList.cpp:117
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
#define index_t
Definition gsConfig.h:32
Provides the gsDofMapper class for re-indexing DoFs.
Knot vector for B-splines.
Provides a list of labeled parameters/options that can be set and accessed easily.
The G+Smo namespace, containing all definitions for the library.
S give(S &x)
Definition gsMemory.h:266