G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsFlowVisitors.h
Go to the documentation of this file.
1
12#pragma once
13
14#include <unordered_map>
15
20
21namespace gismo
22{
23
27template <class T, int MatOrder>
29{
30
31protected: // *** Type definitions ***
32
36
37
38protected: // *** Class members ***
39
40 // constant and same for all visitors
41 index_t m_patchID;
42 gsFlowSolverParams<T> m_params; // pde, bases, assemblerOptions, options, precOptions
43 // pde members: dim, viscosity, rhs (f,g), gsBoundaryConditions, patches, unknownDim
44
45 // constant for individual visitors
46 index_t m_testUnkID, m_shapeUnkID; // used, e.g., to reference the corresponding mapper
47 unsigned m_geoFlags, m_testFunFlags, m_shapeFunFlags;
48 const gsBasis<T>* m_testBasisPtr;
49 const gsBasis<T>* m_shapeBasisPtr;
50 std::vector< gsDofMapper > m_dofMappers;
51 std::vector< gsFlowTerm<T>* > m_terms;
52
53 // updated repeatedly
54 gsMatrix<T> m_localMat;
55 gsMatrix<index_t> m_testFunActives, m_shapeFunActives;
56 gsMapData<T> m_mapData; // members: points, dim, patchID
57
58 // will be changed:
59 gsQuadRule<T> m_quRule;
60 gsMatrix<T> m_quNodes;
61 gsVector<T> m_quWeights;
62 std::vector< gsMatrix<T> > m_testFunData;
63 std::vector< gsMatrix<T> > m_shapeFunData;
64
65
66public: // *** Constructor/destructor ***
67
68 gsFlowVisitor() {}
69
70 gsFlowVisitor(const gsFlowSolverParams<T>& params) : m_params(params)
71 { }
72
74 {
75 deleteTerms();
76 }
77
78
79protected: // *** Member functions ***
80
81 void deleteTerms()
82 {
83 for(size_t i = 0; i < m_terms.size(); i++)
84 delete m_terms[i];
85
86 m_terms.clear();
87 }
88
90 virtual void defineTerms()
92
96
98 void gatherEvalFlags();
99
101 virtual void defineTestShapeBases()
102 {
103 m_testBasisPtr = &(m_params.getBases()[m_testUnkID].piece(m_patchID));
104 m_shapeBasisPtr = &(m_params.getBases()[m_shapeUnkID].piece(m_patchID));
105 }
106
108 virtual void setupQuadrature();
109 //{ GISMO_NO_IMPLEMENTATION }
110
116 void evalSingleFunData(const unsigned& basisFlags, const gsBasis<T>* basisPtr, const index_t funID, std::vector< gsMatrix<T> >& basisData);
117
123 void evalBasisData(const unsigned& basisFlags, const gsBasis<T>* basisPtr, gsMatrix<index_t>& activesUnique, std::vector< gsMatrix<T> >& basisData);
124
125
126public: // *** Member functions ***
127
129 void initialize();
130
133 void updateDofMappers(const std::vector<gsDofMapper>& mappers) { m_dofMappers = mappers; }
134
137 void initOnPatch(index_t patchID);
138
141 void setCurrentSolution(std::vector<gsField<T> >& solutions);
142
145 void setCurrentSolution(gsField<T>& solution);
146
149 void evaluate(index_t testFunID);
150
153 void evaluate(const gsDomainIterator<T>* domIt);
154
156 virtual void assemble();
157
162 virtual void localToGlobal(const std::vector<gsMatrix<T> >& eliminatedDofs, gsSparseMatrix<T, MatOrder>& globalMat, gsMatrix<T>& globalRhs)
164
167 virtual void localToGlobal(gsMatrix<T>& globalRhs)
169
170};
171
172// ===================================================================================================================
173
174template <class T, int MatOrder>
175class gsFlowVisitorVectorValued : public gsFlowVisitor<T, MatOrder>
176{
177
178public:
179 typedef gsFlowVisitor<T, MatOrder> Base;
180
181
182protected: // *** Class members ***
183
184 std::vector< gsMatrix<T> > m_locMatVec;
185
186
187protected: // *** Base class members ***
188
189 using Base::m_params;
190 using Base::m_mapData;
191 using Base::m_testFunActives;
192 using Base::m_shapeFunActives;
193 using Base::m_quWeights;
194 using Base::m_testFunData;
195 using Base::m_shapeFunData;
196 using Base::m_terms;
197
198public: // *** Constructor/destructor ***
199
200 gsFlowVisitorVectorValued() {}
201
202 gsFlowVisitorVectorValued(const gsFlowSolverParams<T>& params) : Base(params)
203 { }
204
205
206public: // *** Member functions ***
207
208 virtual void assemble();
209
210};
211
212// ===================================================================================================================
213
214// TODO
215
216// /// @brief Base class for incompressible flow boundary visitors.
217// /// @tparam T real number type
218// /// @tparam MatOrder sparse matrix storage order (ColMajor/RowMajor)
219// template <class T, int MatOrder>
220// class gsFlowVisitorBnd : public gsFlowVisitor<T, MatOrder>
221// {
222
223// public:
224// typedef gsFlowVisitor<T, MatOrder> Base;
225
226// protected: // *** Class members ***
227
228// boxSide m_side;
229
230// protected: // *** Base class members ***
231
232// //using Base::;
233
234// protected: // *** Class members ***
235
236
237
238// public: // *** Constructor/destructor ***
239
240// gsFlowVisitorBnd() {}
241
242// gsFlowVisitorBnd(const gsFlowSolverParams<T>& params) : Base(params)
243// { }
244
245
246// protected: // *** Member functions ***
247
248// /// @brief Set pointers to test and shape basis.
249// virtual void defineTestShapeBases()
250// {
251// m_testBasisPtr = &(m_params.getBases()[m_testUnkID].piece(m_patchID).boundaryBasis(m_side));
252// m_shapeBasisPtr = &(m_params.getBases()[m_shapeUnkID].piece(m_patchID).boundaryBasis(m_side));
253// }
254
255// /// @brief Setup the quadrature rule.
256// void setupQuadrature();
257
258
259// public: // *** Member functions ***
260
261// /// @brief Initialize the visitor on the given patch.
262// /// @param[in] patchID the patch number
263// void initOnPatchSide(index_t patchID, boxSide side)
264// {
265// m_side = side;
266// Base::initOnPatch(patchID);
267// }
268
269// /// @brief Evaluate basis data on the support of a given test function (used for row-by-row assembly).
270// /// @param[in] testFunID the local test function index on the current patch
271// void evaluate(index_t testFunID);
272
273// /// @brief Evaluate basis data on the current element (used for element-by-element assembly).
274// /// @param[in] domIt domain iterator pointing to the current element
275// void evaluate(const gsDomainIterator<T>* domIt);
276// };
277
278} // namespace gismo
279
280#ifndef GISMO_BUILD_LIB
281#include GISMO_HPP_HEADER(gsFlowVisitors.hpp)
282#endif
A basis represents a family of scalar basis functions defined over a common parameter domain.
Definition gsBasis.h:79
Class which enables iteration over all elements of a parameter domain.
Definition gsDomainIterator.h:68
A scalar of vector field defined on a m_parametric geometry.
Definition gsField.h:55
A class that holds all parameters needed by the incompressible flow solver.
Definition gsFlowSolverParams.h:34
std::vector< gsMultiBasis< T > > & getBases()
Returns a reference to the discretization bases.
Definition gsFlowSolverParams.h:155
A class for integrals of the form: test function value * shape function value.
Definition gsFlowTerms.h:166
Base class for incompressible flow visitors.
Definition gsFlowVisitors.h:29
void setCurrentSolution(std::vector< gsField< T > > &solutions)
Set all needed current solution fields.
Definition gsFlowVisitors.hpp:164
virtual void defineTestShapeBases()
Set pointers to test and shape basis.
Definition gsFlowVisitors.h:101
virtual void defineTestShapeUnknowns()
Set test and shape basis unknown IDs according to visitor type.
Definition gsFlowVisitors.h:94
void initialize()
Initialize the visitor.
Definition gsFlowVisitors.hpp:137
void evalBasisData(const unsigned &basisFlags, const gsBasis< T > *basisPtr, gsMatrix< index_t > &activesUnique, std::vector< gsMatrix< T > > &basisData)
Evaluate required data for the given basis.
Definition gsFlowVisitors.hpp:46
virtual void setupQuadrature()
Setup the quadrature rule.
Definition gsFlowVisitors.hpp:217
virtual void defineTerms()
Define terms to be evaluated according to visitor type and given parameters.
Definition gsFlowVisitors.h:90
void gatherEvalFlags()
Gather evaluation flags from all terms.
Definition gsFlowVisitors.hpp:19
void initOnPatch(index_t patchID)
Initialize the visitor on the given patch.
Definition gsFlowVisitors.hpp:150
virtual void localToGlobal(gsMatrix< T > &globalRhs)
Map local rhs vector to the global rhs vector.
Definition gsFlowVisitors.h:167
virtual void localToGlobal(const std::vector< gsMatrix< T > > &eliminatedDofs, gsSparseMatrix< T, MatOrder > &globalMat, gsMatrix< T > &globalRhs)
Map local matrix to the global matrix.
Definition gsFlowVisitors.h:162
void updateDofMappers(const std::vector< gsDofMapper > &mappers)
Update DoF mappers.
Definition gsFlowVisitors.h:133
virtual void assemble()
Assemble the local matrix.
Definition gsFlowVisitors.hpp:190
void evalSingleFunData(const unsigned &basisFlags, const gsBasis< T > *basisPtr, const index_t funID, std::vector< gsMatrix< T > > &basisData)
Evaluate required data for the given basis function.
Definition gsFlowVisitors.hpp:31
A class for integrals of the form: pressure shape function value * velocity test function divergence.
Definition gsINSTerms.h:23
A class for integrals of the form: velocity solution * shape function gradient * test function value.
Definition gsINSTerms.h:77
the gsMapData is a cache of pre-computed function (map) values.
Definition gsFuncData.h:349
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class representing a reference quadrature rule.
Definition gsQuadRule.h:29
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
A class that holds all parameters needed by the incompressible flow solver.
The G+Smo namespace, containing all definitions for the library.