23 class gsBiharmonicExprAssembler
35 gsBiharmonicExprAssembler(
const gsMultiPatch<T> & mp,
36 const gsMultiBasis<T> & mb,
37 const gsFunctionSet<T> & force,
38 const gsBoundaryConditions<T> & bcs
42 gsBiharmonicExprAssembler() { }
45 gsBiharmonicExprAssembler(
const gsBiharmonicExprAssembler& other )
51 gsBiharmonicExprAssembler( gsBiharmonicExprAssembler&& other )
53 operator=(
give(other));
57 gsBiharmonicExprAssembler& operator= (
const gsBiharmonicExprAssembler& other );
60 gsBiharmonicExprAssembler& operator= ( gsBiharmonicExprAssembler&& other );
64 void _defaultOptions();
70 void _setup(
const expr::gsFeSpace<T> & u);
75 void setSpaceBasis(
const gsFunctionSet<T> & spaceBasis);
89 const gsSparseMatrix<T> & matrix()
const {
return m_assembler.matrix(); }
90 const gsMatrix<T> & rhs()
const {
return m_assembler.rhs(); }
92 T penalty()
const {
return m_penalty; };
93 index_t numDofs()
const {
return m_assembler.numDofs(); };
96 T l2error(gsMatrix<T> & solVector,
const gsFunctionSet<T> & exact);
99 T h1error(gsMatrix<T> & solVector,
const gsFunctionSet<T> & exact);
102 T h2error(gsMatrix<T> & solVector,
const gsFunctionSet<T> & exact);
105 std::tuple<T,T,T> errors(gsMatrix<T> & solVector,
const gsFunctionSet<T> & exact);
108 T interfaceError(gsMatrix<T> & solVector,
const gsFunctionSet<T> & exact);
111 gsOptionList & options() {
return m_options;}
114 void setOptions(gsOptionList & options);
117 void constructSolution(gsMatrix<T> & solVector);
120 typename gsFunctionSet<T>::Ptr getSolution()
const;
125 void _setMapperForBiharmonic(
const gsBoundaryConditions<T> & bc,
126 const gsFunctionSet<T> & bb2,
127 gsDofMapper & mapper);
129 void _getDirichletNeumannValuesL2Projection(
const gsMultiPatch<T> & mp,
130 const gsMultiBasis<T> & dbasis,
131 const gsBoundaryConditions<T> & bc,
132 const gsFunctionSet<T> & bb2,
133 const expr::gsFeSpace<T> & u);
135 void _computeStabilityParameter(
const gsMultiPatch<T> & mp,
136 const gsMultiBasis<T> & dbasis,
137 gsMatrix<T> & mu_interfaces);
141 typedef typename gsExprAssembler<T>::geometryMap geometryMap;
142 typedef typename gsExprAssembler<T>::space space;
143 typedef typename gsExprAssembler<T>::solution solution;
144 typedef typename gsExprAssembler<T>::element element;
148 mutable bool m_second;
151 gsExprAssembler<T> m_assembler;
152 gsExprEvaluator<T> m_evaluator;
154 gsMultiPatch<T> m_patches;
155 mutable gsMultiBasis<T> m_basis;
156 const gsFunctionSet<T> * m_spaceBasis;
157 const gsFunctionSet<T> * m_force;
158 gsBoundaryConditions<T> m_bcs;
160 mutable gsMappedSpline<2,T> m_mspline;
161 mutable gsMultiPatch<T> m_sol;
163 mutable gsOptionList m_options;
170 #ifdef GISMO_WITH_PYBIND11
175 void pybind11_init_gsBiharmonicExprAssembler(pybind11::module &m);
177 #endif // GISMO_WITH_PYBIND11
186 #ifndef GISMO_BUILD_LIB
187 #include GISMO_HPP_HEADER(gsBiharmonicExprAssembler.hpp)
S give(S &x)
Definition: gsMemory.h:266
#define index_t
Definition: gsConfig.h:32
Generic expressions evaluator.
Generic expressions matrix assembly.