G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsFlowTerms.h
Go to the documentation of this file.
1
12#pragma once
13#include <gismo.h>
14
15namespace gismo
16{
17
20template <class T>
22{
23
24public: // *** Smart pointers ***
25
26 typedef memory::shared_ptr<gsFlowTerm> Ptr;
27 typedef memory::unique_ptr<gsFlowTerm> uPtr;
28
29
30protected: // *** Class members ***
31
32 unsigned m_geoFlags, m_testFunFlags, m_shapeFunFlags; // evaluation flags
33 gsVector<T> m_coeff;
34
35
36public: // *** Constructor/destructor ***
37
39 {
40 m_geoFlags = 0;
41 m_testFunFlags = 0;
42 m_shapeFunFlags = 0;
43 }
44
45 virtual ~gsFlowTerm()
46 {}
47
48
49public: // *** Member functions ***
50
57 virtual void assemble(const gsMapData<T>& mapData, const gsVector<T>& quWeights, const std::vector< gsMatrix<T> >& testFunData, const std::vector< gsMatrix<T> >& shapeFunData, gsMatrix<T>& localMat)
59
66 virtual void assemble(const gsMapData<T>& mapData, const gsVector<T>& quWeights, const std::vector< gsMatrix<T> >& testFunData, const std::vector< gsMatrix<T> >& shapeFunData, std::vector< gsMatrix<T> >& localMat)
68
69 void updateEvalFlags(unsigned& geoFlags, unsigned& testFunFlags, unsigned& shapeFunFlags)
70 {
71 geoFlags |= m_geoFlags;
72 testFunFlags |= m_testFunFlags;
73 shapeFunFlags |= m_shapeFunFlags;
74 }
75
76
77protected: // *** Member functions ***
78
79 void setConstCoeff(T value)
80 {
81 m_coeff.resize(1);
82 m_coeff << value;
83 }
84
94 virtual void evalCoeff(const gsMapData<T>& mapData)
95 { setConstCoeff(1.0); }
96
97 virtual gsVector<T> getCoeffGeoMapProduct(const gsMapData<T>& mapData);
98
99};
100
101// ===================================================================================================================
102
105template <class T>
107{
108
109public: // *** Smart pointers ***
110
111 typedef memory::shared_ptr<gsFlowTermNonlin> Ptr;
112 typedef memory::unique_ptr<gsFlowTermNonlin> uPtr;
113
114public: // *** Type definitions ***
115
116 typedef gsFlowTerm<T> Base;
117
118
119protected: // *** Class members ***
120
121 gsField<T> m_currentSolU;
122 bool m_isCurrentSolSet;
123 gsMatrix<T> m_solUVals;
124
125
126public: // *** Constructor/destructor ***
127
129 { }
130
131
132public: // *** Member functions ***
133
134 void setCurrentSolution(std::vector<gsField<T> >& solutions)
135 {
136 m_currentSolU = solutions.front();
137 m_isCurrentSolSet = true;
138 }
139
140 void setCurrentSolution(gsField<T>& solution)
141 {
142 m_currentSolU = solution;
143 m_isCurrentSolSet = true;
144 }
145
146
147protected: // *** Member functions ***
148
149 virtual void computeCoeffSolU(const gsMapData<T>& mapData)
150 {
151 GISMO_ASSERT(m_isCurrentSolSet, "No velocity solution set in gsFlowTermNonlin.");
152
153 m_solUVals.resize(mapData.dim.first, mapData.points.cols());
154 m_solUVals = m_currentSolU.value(mapData.points, mapData.patchId);
155 }
156
157};
158
159// ===================================================================================================================
160// ===================================================================================================================
161
164template <class T>
166{
167
168public: // *** Constructor/destructor ***
169
171 {
172 this->m_geoFlags = NEED_MEASURE;
173 this->m_testFunFlags = NEED_VALUE;
174 this->m_shapeFunFlags = NEED_VALUE;
175 }
176
177
178public: // *** Member functions ***
179
180 virtual void assemble(const gsMapData<T>& mapData, const gsVector<T>& quWeights, const std::vector< gsMatrix<T> >& testFunData, const std::vector< gsMatrix<T> >& shapeFunData, gsMatrix<T>& localMat);
181
182};
183
184// ===================================================================================================================
185
188template <class T>
190{
191
192protected: // *** Class members ***
193
194 real_t m_timeStep;
195
196
197public: // *** Constructor/destructor ***
198
199 gsFlowTerm_TimeDiscr(real_t timeStep) :
200 m_timeStep(timeStep)
201 { }
202
203
204protected: // *** Member functions ***
205
206 virtual void evalCoeff(const gsMapData<T>& mapData0)
207 { this->setConstCoeff(1./m_timeStep); }
208
209};
210
211// ===================================================================================================================
212
215template <class T>
217{
218
219public: // *** Constructor/destructor ***
220
222 {
223 this->m_geoFlags = NEED_MEASURE | NEED_GRAD_TRANSFORM;
224 this->m_testFunFlags = NEED_DERIV;
225 this->m_shapeFunFlags = NEED_DERIV;
226 }
227
228
229public: // *** Member functions ***
230
231 virtual void assemble(const gsMapData<T>& mapData, const gsVector<T>& quWeights, const std::vector< gsMatrix<T> >& testFunData, const std::vector< gsMatrix<T> >& shapeFunData, gsMatrix<T>& localMat);
232
233};
234
235// ===================================================================================================================
236
239template <class T>
241{
242
243protected: // *** Class members ***
244
245 real_t m_viscosity;
246
247public: // *** Constructor/destructor ***
248
249 gsFlowTerm_Diffusion(real_t viscosity) :
250 m_viscosity(viscosity)
251 { }
252
253
254protected: // *** Member functions ***
255
256 virtual void evalCoeff(const gsMapData<T>& mapData)
257 { this->setConstCoeff(m_viscosity); }
258
259};
260
261// ===================================================================================================================
262// ===================================================================================================================
263
266template <class T>
267class gsFlowTerm_rhs : public gsFlowTerm<T>
268{
269
270protected: // *** Class members ***
271
272 const gsFunction<T>* m_pRhsFun;
273 gsMatrix<T> m_rhsVals;
274
275public: // *** Constructor/destructor ***
276
277 gsFlowTerm_rhs(const gsFunction<T>* pRhsFun)
278 {
279 this->m_geoFlags = NEED_VALUE | NEED_MEASURE;
280 this->m_testFunFlags = NEED_VALUE;
281 m_pRhsFun = pRhsFun;
282 }
283
284
285public: // *** Member functions ***
286
287 virtual void assemble(const gsMapData<T>& mapData, const gsVector<T>& quWeights, const std::vector< gsMatrix<T> >& testFunData, const std::vector< gsMatrix<T> >& shapeFunData, gsMatrix<T>& localMat);
288
289};
290
291
292} // namespace gismo
293
294#ifndef GISMO_BUILD_LIB
295#include GISMO_HPP_HEADER(gsFlowTerms.hpp)
296#endif
A scalar of vector field defined on a m_parametric geometry.
Definition gsField.h:55
gsMatrix< T > value(const gsMatrix< T > &u, int i=0) const
Evaluation of the field at points u.
Definition gsField.h:122
A class computing nonlinear terms of the weak formulation appearing in incompressible flow problems.
Definition gsFlowTerms.h:107
A class for integrals of the form: viscosity * test function gradient * shape function gradient.
Definition gsFlowTerms.h:241
virtual void evalCoeff(const gsMapData< T > &mapData)
Evaluate the term coefficient.
Definition gsFlowTerms.h:256
A class for integrals of the form: test function gradient * shape function gradient.
Definition gsFlowTerms.h:217
virtual void assemble(const gsMapData< T > &mapData, const gsVector< T > &quWeights, const std::vector< gsMatrix< T > > &testFunData, const std::vector< gsMatrix< T > > &shapeFunData, gsMatrix< T > &localMat)
Assemble the current local matrix.
Definition gsFlowTerms.hpp:56
A class for integrals of the form: (1 / time step) * test function value * shape function value.
Definition gsFlowTerms.h:190
virtual void evalCoeff(const gsMapData< T > &mapData0)
Evaluate the term coefficient.
Definition gsFlowTerms.h:206
A class for integrals of the form: test function value * shape function value.
Definition gsFlowTerms.h:166
virtual void assemble(const gsMapData< T > &mapData, const gsVector< T > &quWeights, const std::vector< gsMatrix< T > > &testFunData, const std::vector< gsMatrix< T > > &shapeFunData, gsMatrix< T > &localMat)
Assemble the current local matrix.
Definition gsFlowTerms.hpp:37
A class for integrals of the form: test function value * rhs function value.
Definition gsFlowTerms.h:268
virtual void assemble(const gsMapData< T > &mapData, const gsVector< T > &quWeights, const std::vector< gsMatrix< T > > &testFunData, const std::vector< gsMatrix< T > > &shapeFunData, gsMatrix< T > &localMat)
Assemble the current local matrix.
Definition gsFlowTerms.hpp:80
A class computing individual terms of the weak formulation appearing in incompressible flow problems.
Definition gsFlowTerms.h:22
virtual void assemble(const gsMapData< T > &mapData, const gsVector< T > &quWeights, const std::vector< gsMatrix< T > > &testFunData, const std::vector< gsMatrix< T > > &shapeFunData, gsMatrix< T > &localMat)
Assemble the current local matrix.
Definition gsFlowTerms.h:57
virtual void assemble(const gsMapData< T > &mapData, const gsVector< T > &quWeights, const std::vector< gsMatrix< T > > &testFunData, const std::vector< gsMatrix< T > > &shapeFunData, std::vector< gsMatrix< T > > &localMat)
Assemble the current local matrices.
Definition gsFlowTerms.h:66
virtual void evalCoeff(const gsMapData< T > &mapData)
Evaluate the term coefficient.
Definition gsFlowTerms.h:94
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
the gsMapData is a cache of pre-computed function (map) values.
Definition gsFuncData.h:349
gsMatrix< T > points
input (parametric) points
Definition gsFuncData.h:372
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
Main header to be included by clients using the G+Smo library.
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
The G+Smo namespace, containing all definitions for the library.
@ NEED_VALUE
Value of the object.
Definition gsForwardDeclarations.h:72
@ NEED_DERIV
Gradient of the object.
Definition gsForwardDeclarations.h:73
@ NEED_GRAD_TRANSFORM
Gradient transformation matrix.
Definition gsForwardDeclarations.h:77
@ NEED_MEASURE
The density of the measure pull back.
Definition gsForwardDeclarations.h:76