G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsINSVisitors.h
Go to the documentation of this file.
1
12#pragma once
13
16
17namespace gismo
18{
19
20// ===================================================================================================================
21// VELOCITY-VELOCITY VISITORS
22
23template <class T, int MatOrder>
24class gsINSVisitorUU : public gsFlowVisitor<T, MatOrder>
25{
26
27public:
28 typedef gsFlowVisitor<T, MatOrder> Base;
29
30
31protected: // *** Base class members ***
32
33 using Base::m_params;
34 using Base::m_patchID;
35 using Base::m_testUnkID;
36 using Base::m_shapeUnkID;
37 using Base::m_dofMappers;
38 using Base::m_testFunActives;
39 using Base::m_shapeFunActives;
40 using Base::m_localMat;
41
42public: // *** Constructor/destructor ***
43
44 gsINSVisitorUU() {}
45
46 gsINSVisitorUU(const gsFlowSolverParams<T>& params) : Base(params)
47 { }
48
49
50protected: // *** Member functions ***
51
52 virtual void defineTestShapeUnknowns()
53 {
54 m_testUnkID = 0; // velocity
55 m_shapeUnkID = 0; // velocity
56 }
57
58public: // *** Member functions ***
59
60 virtual void localToGlobal(const std::vector<gsMatrix<T> >& eliminatedDofs, gsSparseMatrix<T, MatOrder>& globalMat, gsMatrix<T>& globalRhs);
61
62};
63
64// ===================================================================================================================
65
66template <class T, int MatOrder>
67class gsINSVisitorUUlin : public gsINSVisitorUU<T, MatOrder>
68{
69
70public:
71 typedef gsINSVisitorUU<T, MatOrder> Base;
72
73
74protected: // *** Base class members ***
75
76 using Base::m_params;
77 using Base::m_terms;
78
79
80public: // *** Constructor/destructor ***
81
82 gsINSVisitorUUlin() {}
83
84 gsINSVisitorUUlin(const gsFlowSolverParams<T>& params) :
85 Base(params)
86 { }
87
88
89protected: // *** Member functions ***
90
91 virtual void defineTerms()
92 {
93 m_terms.push_back( new gsFlowTerm_Diffusion<T>(m_params.getPde().viscosity()) );
94
95 // if(m_params.options().getSwitch("unsteady"))
96 // m_terms.push_back( new gsFlowTermTimeDiscr<T>(m_params.options().getReal("timeStep")) );
97
98 // ... other terms, e.g. from stabilizations
99 }
100
101};
102
103// ===================================================================================================================
104
105template <class T, int MatOrder>
106class gsINSVisitorUUnonlin : public gsINSVisitorUU<T, MatOrder>
107{
108
109public:
110 typedef gsINSVisitorUU<T, MatOrder> Base;
111
112
113protected: // *** Base class members ***
114
115 using Base::m_params;
116 using Base::m_terms;
117
118
119public: // *** Constructor/destructor ***
120
121 gsINSVisitorUUnonlin() {}
122
123 gsINSVisitorUUnonlin(const gsFlowSolverParams<T>& params) : Base(params)
124 { }
125
126
127protected: // *** Member functions ***
128
129 virtual void defineTerms()
130 {
131 m_terms.push_back( new typename gsFlowVisitor<T, MatOrder>::ConvectionTerm() );
132
133 // ... other terms, e.g. from stabilizations
134 }
135
136};
137
138// ===================================================================================================================
139
140template <class T, int MatOrder>
141class gsINSVisitorUUmass : public gsINSVisitorUU<T, MatOrder>
142{
143
144public:
145 typedef gsINSVisitorUU<T, MatOrder> Base;
146
147
148protected: // *** Base class members ***
149
150 using Base::m_params;
151 using Base::m_terms;
152
153
154public: // *** Constructor/destructor ***
155
156 gsINSVisitorUUmass() {}
157
158 gsINSVisitorUUmass(const gsFlowSolverParams<T>& params) :
159 Base(params)
160 { }
161
162
163protected: // *** Member functions ***
164
165 virtual void defineTerms()
166 {
167 m_terms.push_back( new gsFlowTerm_ValVal<T>() );
168 }
169
170};
171
172// ===================================================================================================================
173
174template <class T, int MatOrder>
175class gsINSVisitorUUtimeDiscr : public gsINSVisitorUU<T, MatOrder>
176{
177
178public:
179 typedef gsINSVisitorUU<T, MatOrder> Base;
180
181
182protected: // *** Base class members ***
183
184 using Base::m_params;
185 using Base::m_terms;
186
187
188public: // *** Constructor/destructor ***
189
190 gsINSVisitorUUtimeDiscr() {}
191
192 gsINSVisitorUUtimeDiscr(const gsFlowSolverParams<T>& params) :
193 Base(params)
194 { }
195
196
197protected: // *** Member functions ***
198
199 virtual void defineTerms()
200 {
201 m_terms.push_back( new gsFlowTerm_TimeDiscr<T>(m_params.options().getReal("timeStep")) );
202 }
203
204};
205
206
207// ===================================================================================================================
208// ===================================================================================================================
209
210// VELOCITY-PRESSURE VISITORS
211template <class T, int MatOrder>
212class gsINSVisitorPU : public gsFlowVisitorVectorValued<T, MatOrder> // order: shape, test
213{
214
215public:
216 typedef gsFlowVisitorVectorValued<T, MatOrder> Base;
217
218
219protected: // *** Base class members ***
220
221 using Base::m_locMatVec;
222 using Base::m_params;
223 using Base::m_patchID;
224 using Base::m_testUnkID;
225 using Base::m_shapeUnkID;
226 using Base::m_dofMappers;
227 using Base::m_testFunActives;
228 using Base::m_shapeFunActives;
229 using Base::m_terms;
230
231public: // *** Constructor/destructor ***
232
233 gsINSVisitorPU() {}
234
235 gsINSVisitorPU(const gsFlowSolverParams<T>& params) : Base(params)
236 { }
237
238
239protected: // *** Member functions ***
240
241 virtual void defineTerms()
242 {
243 m_terms.push_back( new gsINSTerm_PvalUdiv<T>() );
244 }
245
246 virtual void defineTestShapeUnknowns()
247 {
248 m_testUnkID = 0; // velocity
249 m_shapeUnkID = 1; // pressure
250 }
251
252public: // *** Member functions ***
253
254 virtual void localToGlobal(const std::vector<gsMatrix<T> >& eliminatedDofs, gsSparseMatrix<T, MatOrder>& globalMat, gsMatrix<T>& globalRhs);
255
256};
257
258// ===================================================================================================================
259
260template <class T, int MatOrder>
261class gsINSVisitorPU_withUPrhs : public gsINSVisitorPU<T, MatOrder> // order: shape, test
262{
263
264public:
265 typedef gsINSVisitorPU<T, MatOrder> Base;
266
267
268protected: // *** Base class members ***
269
270 using Base::m_locMatVec;
271 using Base::m_params;
272 using Base::m_patchID;
273 using Base::m_testUnkID;
274 using Base::m_shapeUnkID;
275 using Base::m_dofMappers;
276 using Base::m_testFunActives;
277 using Base::m_shapeFunActives;
278
279public: // *** Constructor/destructor ***
280
281 gsINSVisitorPU_withUPrhs() {}
282
283 gsINSVisitorPU_withUPrhs(const gsFlowSolverParams<T>& params) : Base(params)
284 { }
285
286
287public: // *** Member functions ***
288
289 virtual void localToGlobal(const std::vector<gsMatrix<T> >& eliminatedDofs, gsSparseMatrix<T, MatOrder>& globalMat, gsMatrix<T>& globalRhs);
290
291};
292
293// ===================================================================================================================
294// ===================================================================================================================
295
296// PRESSURE-VELOCITY VISITORS
297template <class T, int MatOrder>
298class gsINSVisitorUP : public gsFlowVisitorVectorValued<T, MatOrder> // order: shape, test
299{
300
301public:
302 typedef gsFlowVisitorVectorValued<T, MatOrder> Base;
303
304
305protected: // *** Base class members ***
306
307 using Base::m_locMatVec;
308 using Base::m_params;
309 using Base::m_patchID;
310 using Base::m_testUnkID;
311 using Base::m_shapeUnkID;
312 using Base::m_dofMappers;
313 using Base::m_testFunActives;
314 using Base::m_shapeFunActives;
315 using Base::m_terms;
316
317public: // *** Constructor/destructor ***
318
319 gsINSVisitorUP() {}
320
321 gsINSVisitorUP(const gsFlowSolverParams<T>& params) : Base(params)
322 { }
323
324
325protected: // *** Member functions ***
326
327 virtual void defineTerms()
328 {
329 m_terms.push_back( new gsINSTerm_UdivPval<T>() );
330 }
331
332 virtual void defineTestShapeUnknowns()
333 {
334 m_testUnkID = 1; // pressure
335 m_shapeUnkID = 0; // velocity
336 }
337
338public: // *** Member functions ***
339
340 virtual void localToGlobal(const std::vector<gsMatrix<T> >& eliminatedDofs, gsSparseMatrix<T, MatOrder>& globalMat, gsMatrix<T>& globalRhs);
341
342};
343
344// ===================================================================================================================
345// ===================================================================================================================
346
347// PRESSURE-PRESSURE VISITORS
348template <class T, int MatOrder>
349class gsINSVisitorPP : public gsFlowVisitor<T, MatOrder>
350{
351
352public:
353 typedef gsFlowVisitor<T, MatOrder> Base;
354
355
356protected: // *** Base class members ***
357
358 using Base::m_params;
359 using Base::m_patchID;
360 using Base::m_testUnkID;
361 using Base::m_shapeUnkID;
362 using Base::m_dofMappers;
363 using Base::m_testFunActives;
364 using Base::m_shapeFunActives;
365 using Base::m_localMat;
366
367
368public: // *** Constructor/destructor ***
369
370 gsINSVisitorPP() {}
371
372 gsINSVisitorPP(const gsFlowSolverParams<T>& params) : Base(params)
373 { }
374
375
376protected: // *** Member functions ***
377
378 virtual void defineTestShapeUnknowns()
379 {
380 m_testUnkID = 1; // pressure
381 m_shapeUnkID = 1; // pressure
382 }
383
384public: // *** Member functions ***
385
386 virtual void localToGlobal(const std::vector<gsMatrix<T> >& eliminatedDofs, gsSparseMatrix<T, MatOrder>& globalMat, gsMatrix<T>& globalRhs);
387
388};
389
390// ===================================================================================================================
391
392// template <class T, int MatOrder>
393// class gsINSVisitorPPlin : public gsINSVisitorPP<T, MatOrder>
394// {
395
396// public:
397// typedef gsINSVisitorPP<T, MatOrder> Base;
398
399
400// protected: // *** Base class members ***
401
402// using Base::m_params;
403// using Base::m_terms;
404
405
406// public: // *** Constructor/destructor ***
407
408// gsINSVisitorPPlin() {}
409
410// gsINSVisitorPPlin(const gsFlowSolverParams<T>& params) : Base(params)
411// { }
412
413
414// protected: // *** Member functions ***
415
416// virtual void defineTerms()
417// {
418// // no default pressure-pressure terms in the Navier-Stokes eqns
419// // optionally stabilization for inf-sup unstable discretizations
420// }
421
422// };
423
424// // ===================================================================================================================
425
426// template <class T, int MatOrder>
427// class gsINSVisitorPPnonlin : public gsINSVisitorPP<T, MatOrder>
428// {
429
430// public:
431// typedef gsINSVisitorPP<T, MatOrder> Base;
432
433
434// protected: // *** Base class members ***
435
436// using Base::m_params;
437// using Base::m_terms;
438
439
440// public: // *** Constructor/destructor ***
441
442// gsINSVisitorPPnonlin() {}
443
444// gsINSVisitorPPnonlin(const gsFlowSolverParams<T>& params) : Base(params)
445// { }
446
447
448// protected: // *** Member functions ***
449
450// virtual void defineTerms()
451// {
452// // no default pressure-pressure terms in the Navier-Stokes eqns
453// // optionally stabilization for inf-sup unstable discretizations
454// }
455
456// };
457
458// ===================================================================================================================
459
460template <class T, int MatOrder>
461class gsINSVisitorPPmass : public gsINSVisitorPP<T, MatOrder>
462{
463
464public:
465 typedef gsINSVisitorPP<T, MatOrder> Base;
466
467
468protected: // *** Base class members ***
469
470 using Base::m_params;
471 using Base::m_terms;
472
473
474public: // *** Constructor/destructor ***
475
476 gsINSVisitorPPmass() {}
477
478 gsINSVisitorPPmass(const gsFlowSolverParams<T>& params) : Base(params)
479 { }
480
481
482protected: // *** Member functions ***
483
484 virtual void defineTerms()
485 {
486 m_terms.push_back( new gsFlowTerm_ValVal<T>() );
487 }
488
489};
490
491// ===================================================================================================================
492
493template <class T, int MatOrder>
494class gsINSVisitorPPlaplace : public gsINSVisitorPP<T, MatOrder>
495{
496
497public:
498 typedef gsINSVisitorPP<T, MatOrder> Base;
499
500
501protected: // *** Base class members ***
502
503 using Base::m_params;
504 using Base::m_terms;
505
506
507public: // *** Constructor/destructor ***
508
509 gsINSVisitorPPlaplace() {}
510
511 gsINSVisitorPPlaplace(const gsFlowSolverParams<T>& params) : Base(params)
512 { }
513
514
515protected: // *** Member functions ***
516
517 virtual void defineTerms()
518 {
519 m_terms.push_back( new gsFlowTerm_GradGrad<T>() );
520 }
521
522};
523
524// ===================================================================================================================
525
526template <class T, int MatOrder>
527class gsINSVisitorPPconvection : public gsINSVisitorPP<T, MatOrder>
528{
529
530public:
531 typedef gsINSVisitorPP<T, MatOrder> Base;
532
533protected: // *** Base class members ***
534
535 using Base::m_params;
536 using Base::m_terms;
537
538
539public: // *** Constructor/destructor ***
540
541 gsINSVisitorPPconvection() {}
542
543 gsINSVisitorPPconvection(const gsFlowSolverParams<T>& params) : Base(params)
544 { }
545
546
547protected: // *** Member functions ***
548
549 virtual void defineTerms()
550 {
551 m_terms.push_back( new gsINSTerm_UsolGradVal<T>() );
552 }
553
554};
555
556// ===================================================================================================================
557
558// template <class T, int MatOrder>
559// class gsINSVisitorPP_PCDrobinBC : public gsINSVisitorPP<T, MatOrder>
560// {
561
562// public:
563// typedef gsINSVisitorPP<T, MatOrder> Base;
564
565// protected: // *** Base class members ***
566
567// using Base::m_params;
568// using Base::m_terms;
569
570
571// public: // *** Constructor/destructor ***
572
573// gsINSVisitorPP_PCDrobinBC() {}
574
575// gsINSVisitorPP_PCDrobinBC(const gsFlowSolverParams<T>& params) : Base(params)
576// { }
577
578
579// protected: // *** Member functions ***
580
581// virtual void defineTerms()
582// {
583// // TODO
584// }
585
586// };
587
588// ===================================================================================================================
589// ===================================================================================================================
590
591// RHS VISITORS
592
593template <class T, int MatOrder>
594class gsINSVisitorRhsU : public gsFlowVisitor<T, MatOrder>
595{
596
597public:
598 typedef gsFlowVisitor<T, MatOrder> Base;
599
600
601protected: // *** Class members ***
602
603 const gsFunction<T>* m_pRhsFun;
604
605
606protected: // *** Base class members ***
607
608 using Base::m_params;
609 using Base::m_terms;
610 using Base::m_patchID;
611 using Base::m_testUnkID;
612 using Base::m_shapeUnkID;
613 using Base::m_dofMappers;
614 using Base::m_testFunActives;
615 using Base::m_localMat;
616 using Base::m_mapData;
617 using Base::m_quWeights;
618 using Base::m_testFunData;
619 using Base::m_shapeFunData;
620
621
622public: // *** Constructor/destructor ***
623
624 gsINSVisitorRhsU() {}
625
626 gsINSVisitorRhsU(const gsFlowSolverParams<T>& params):
627 Base(params), m_pRhsFun(params.getPde().force())
628 {
629 GISMO_ASSERT(m_pRhsFun->targetDim() == m_params.getPde().dim(), "Wrong RHS function passed into gsINSRhsU.");
630 }
631
632
633protected: // *** Member functions ***
634
635 virtual void defineTestShapeUnknowns()
636 {
637 m_testUnkID = 0; // velocity
638 m_shapeUnkID = 0; // velocity (not needed here)
639 }
640
641 virtual void defineTerms()
642 {
643 m_terms.push_back( new gsFlowTerm_rhs<T>(m_pRhsFun) );
644 }
645
646public: // *** Member functions ***
647
648 virtual void assemble();
649
650 virtual void localToGlobal(gsMatrix<T>& globalRhs);
651
652};
653
654// ===================================================================================================================
655
656template <class T, int MatOrder>
657class gsINSVisitorRhsP : public gsFlowVisitor<T, MatOrder>
658{
659
660public:
661 typedef gsFlowVisitor<T, MatOrder> Base;
662
663
664protected: // *** Class members ***
665
666 const gsFunction<T>* m_pRhsFun;
667
668
669protected: // *** Base class members ***
670
671 using Base::m_params;
672 using Base::m_terms;
673 using Base::m_patchID;
674 using Base::m_testUnkID;
675 using Base::m_shapeUnkID;
676 using Base::m_dofMappers;
677 using Base::m_testFunActives;
678 using Base::m_localMat;
679 using Base::m_mapData;
680 using Base::m_quWeights;
681 using Base::m_testFunData;
682 using Base::m_shapeFunData;
683
684
685public: // *** Constructor/destructor ***
686
687 gsINSVisitorRhsP() {}
688
689 gsINSVisitorRhsP(const gsFlowSolverParams<T>& params):
690 Base(params), m_pRhsFun(params.getPde().source())
691 {
692 GISMO_ASSERT(m_pRhsFun == NULL || m_pRhsFun->targetDim() == 1, "Wrong RHS function passed into gsINSRhsP.");
693 }
694
695
696protected: // *** Member functions ***
697
698 virtual void defineTestShapeUnknowns()
699 {
700 m_testUnkID = 1; // pressure
701 m_shapeUnkID = 1; // pressure (not needed here)
702 }
703
704 virtual void defineTerms()
705 {
706 m_terms.push_back( new gsFlowTerm_rhs<T>(m_pRhsFun) );
707 }
708
709public: // *** Member functions ***
710
711 virtual void assemble();
712
713 virtual void localToGlobal(gsMatrix<T>& globalRhs);
714
715};
716
717} // namespace gismo
718
719#ifndef GISMO_BUILD_LIB
720#include GISMO_HPP_HEADER(gsINSVisitors.hpp)
721#endif
A class that holds all parameters needed by the incompressible flow solver.
Definition gsFlowSolverParams.h:34
A class for integrals of the form: viscosity * test function gradient * shape function gradient.
Definition gsFlowTerms.h:241
A class for integrals of the form: test function gradient * shape function gradient.
Definition gsFlowTerms.h:217
A class for integrals of the form: (1 / time step) * test function value * shape function value.
Definition gsFlowTerms.h:190
A class for integrals of the form: test function value * shape function value.
Definition gsFlowTerms.h:166
A class for integrals of the form: test function value * rhs function value.
Definition gsFlowTerms.h:268
Base class for incompressible flow visitors.
Definition gsFlowVisitors.h:29
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
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
virtual short_t targetDim() const
Dimension of the target space.
Definition gsFunctionSet.h:595
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 shape function divergence * pressure test function value.
Definition gsINSTerms.h:53
A class for integrals of the form: velocity solution * shape function gradient * test function value.
Definition gsINSTerms.h:77
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
The G+Smo namespace, containing all definitions for the library.