-
Notifications
You must be signed in to change notification settings - Fork 8
/
Copy pathgsNsAssembler.h
113 lines (83 loc) · 4.04 KB
/
gsNsAssembler.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
/** @file gsNsAssembler.h
@brief Provides matrix and rhs assebmly for stationary and transient incompressible
Stokes and Navier-Stokes equations.
This file is part of the G+Smo library.
This Source Code Form is subject to the terms of the Mozilla Public
License, v. 2.0. If a copy of the MPL was not distributed with this
file, You can obtain one at http://mozilla.org/MPL/2.0/.
Author(s):
A.Shamanskiy (2016 - ...., TU Kaiserslautern)
*/
#pragma once
#include <gsElasticity/gsBaseAssembler.h>
#include <gsElasticity/gsBaseUtils.h>
namespace gismo
{
/** @brief TODO: write.
*/
template <class T>
class gsNsAssembler : public gsBaseAssembler<T>
{
public:
typedef gsBaseAssembler<T> Base;
/// @brief Constructor
gsNsAssembler(const gsMultiPatch<T> & patches,
const gsMultiBasis<T> & basisVel,
const gsMultiBasis<T> & basisPres,
const gsBoundaryConditions<T> & bconditions,
const gsFunction<T> & body_force);
/// @brief Returns the list of default options for assembly
static gsOptionList defaultOptions();
/// @brief Refresh routine to set dof-mappers
virtual void refresh();
//--------------------- SYSTEM ASSEMBLY ----------------------------------//
using Base::assemble;
/// @brief Assembly of the linear system for the Stokes problem
/// @{
virtual void assemble(bool saveEliminationMatrix);
virtual void assemble() { assemble(false); };
/// @}
/// Assembles the tangential linear system for Newton's method given the current solution
/// in the form of free and fixed/Dirichelt degrees of freedom.
/// Returns the status of the assembly (for safe exit)
virtual bool assemble(const gsMatrix<T> & solutionVector,
const std::vector<gsMatrix<T> > & fixedDoFs);
/// Assembles the tangential linear system for Newton's method given the current solution
/// in the form of free and fixed/Dirichelt degrees of freedom.
virtual void assemble(const gsMultiPatch<T> & velocity, const gsMultiPatch<T> & pressure);
//--------------------- SOLUTION CONSTRUCTION ----------------------------------//
using Base::constructSolution;
/// @brief Construct velocity from computed solution vector and fixed degrees of freedom
virtual void constructSolution(const gsMatrix<T> & solVector,
const std::vector<gsMatrix<T> > & fixedDoFs,
gsMultiPatch<T> & velocity) const;
/// @brief Construct velocity and pressure from computed solution vector and fixed degrees of freedom
virtual void constructSolution(const gsMatrix<T> & solVector,
const std::vector<gsMatrix<T> > & fixedDoFs,
gsMultiPatch<T> & velocity, gsMultiPatch<T> & pressure) const;
/// @ brief Construct pressure from computed solution vector
virtual void constructPressure(const gsMatrix<T> & solVector,
const std::vector<gsMatrix<T> > & fixedDoFs,
gsMultiPatch<T> & pressure) const;
//--------------------- SPECIALS ----------------------------------//
/// compute forces acting on a given part of the boundary (drag and lift)
/// if split=true, returns pressure and velocity conrtibution separately
virtual gsMatrix<T> computeForce(const gsMultiPatch<T> & velocity, const gsMultiPatch<T> & pressure,
const std::vector<std::pair<index_t,boxSide> > & bdrySides, bool split = false) const;
protected:
/// a custom reserve function to allocate memory for the sparse matrix
virtual void reserve();
protected:
/// Dimension of the problem
/// parametric dim = physical dim = velocity dim
short_t m_dim;
using Base::m_pde_ptr;
using Base::m_bases;
using Base::m_ddof;
using Base::m_options;
using Base::m_system;
};
} // namespace gismo ends
#ifndef GISMO_BUILD_LIB
#include GISMO_HPP_HEADER(gsNsAssembler.hpp)
#endif