32 template<
typename T >
54 index_t numberOfKnotsToBeInserted = 1,
55 index_t multiplicityOfKnotsToBeInserted = 1,
75 options.
askInt(
"Levels", 3 ),
76 options.
askInt(
"NumberOfKnotsToBeInserted", 1 ),
77 options.
askInt(
"MultiplicityOfKnotsToBeInserted", 1 ),
123 options.
askInt(
"Levels", 3 ),
124 options.
askInt(
"DegreesOfFreedom", 0 ),
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 );
146 m_transferMatrices.clear();
155 { o =
give(m_mBases);
return *
this; }
159 {
return m_transferMatrices; }
163 { o =
give(m_transferMatrices);
return *
this; }
166 std::vector< gsMultiBasis<T> > m_mBases;
167 std::vector< gsSparseMatrix<T, RowMajor> > m_transferMatrices;
172 #ifndef GISMO_BUILD_LIB
173 #include GISMO_HPP_HEADER(gsGridHierarchy.hpp)
Knot vector for B-splines.
const std::vector< gsSparseMatrix< T, RowMajor > > & getTransferMatrices() const
Get the vector of transfer matrices (by reference)
Definition: gsGridHierarchy.h:158
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 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
void clear()
Reset the object (to save memory)
Definition: gsGridHierarchy.h:143
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
S give(S &x)
Definition: gsMemory.h:266
#define index_t
Definition: gsConfig.h:32
Grid Hierarchy.
Definition: gsGridHierarchy.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
Provides the gsDofMapper class for re-indexing DoFs.
Provides a list of labeled parameters/options that can be set and accessed easily.
Holds a set of patch-wise bases and their topology information.
Definition: gsMultiBasis.h:36
const std::vector< gsMultiBasis< T > > & getMultiBases() const
Get the vector of gsMultiBasis objects (by reference)
Definition: gsGridHierarchy.h:150
static gsOptionList defaultOptions()
Get the default options.
Definition: gsGridHierarchy.h:130
gsGridHierarchy & moveTransferMatricesTo(std::vector< gsSparseMatrix< T, RowMajor > > &o)
Get the vector of transfer matrices.
Definition: gsGridHierarchy.h:162
Class containing a set of boundary conditions.
Definition: gsBoundaryConditions.h:341
index_t askInt(const std::string &label, const index_t &value=0) const
Reads value for option label from options.
Definition: gsOptionList.cpp:117
gsGridHierarchy & moveMultiBasesTo(std::vector< gsMultiBasis< T > > &o)
Get the vector of gsMultiBasis objects.
Definition: gsGridHierarchy.h:154
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
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