PolyFEM
Loading...
Searching...
No Matches
Stokes.hpp
Go to the documentation of this file.
1#pragma once
2
5
7
8namespace polyfem::assembler
9{
10 // stokes local assembler for velocity
12 {
13 public:
15
17
18 VectorNd compute_rhs(const AutodiffHessianPt &pt) const override;
19 // res is R^{dimĀ²}
20 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>
21 assemble(const LinearAssemblerData &data) const override;
22
23 void add_multimaterial(const int index, const json &params, const Units &units) override;
24
25 const GenericMatParam &viscosity() const { return viscosity_; }
26
27 virtual std::string name() const override { return "Stokes"; }
28 std::map<std::string, ParamFunc> parameters() const override;
29
30 bool is_fluid() const override { return true; }
31 bool is_tensor() const override { return true; }
32
33 private:
35 };
36
37 // stokes mixed assembler (velocity phi and pressure psi)
39 {
40 public:
41 std::string name() const override { return "StokesMixed"; }
42
43 // res is R^{dim}
44 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 3, 1>
45 assemble(const MixedAssemblerData &data) const override;
46
47 inline int rows() const override { return size(); }
48 inline int cols() const override { return 1; }
49 };
50
51 // pressure only for stokes is zero
53 {
54 public:
56
57 std::string name() const override { return "StokesPressure"; }
58 std::map<std::string, ParamFunc> parameters() const override { return std::map<std::string, ParamFunc>(); }
59
60 // res is R^{dimĀ²}
61 Eigen::Matrix<double, Eigen::Dynamic, 1, 0, 9, 1>
62 assemble(const LinearAssemblerData &data) const override
63 {
64 return Eigen::Matrix<double, 1, 1>::Zero(1, 1);
65 }
66
67 bool is_fluid() const override { return true; }
68
69 void set_size(const int) override { size_ = 1; }
70 };
71
73 {
74 public:
75 OperatorSplitting() = default;
76 ~OperatorSplitting() = default;
77
78 std::string name() const override { return "OperatorSplitting"; }
79 };
80} // namespace polyfem::assembler
assemble matrix based on the local assembler local assembler is eg Laplace, LinearElasticity etc
void assemble(const bool is_volume, const int n_basis, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, StiffnessMatrix &stiffness, const bool is_mass=false) const override
assembles the stiffness matrix for the given basis the bilinear form (local assembler) is encoded by ...
std::string name() const override
Definition Stokes.hpp:78
int rows() const override
Definition Stokes.hpp:47
std::string name() const override
Definition Stokes.hpp:41
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > assemble(const MixedAssemblerData &data) const override
Definition Stokes.cpp:68
int cols() const override
Definition Stokes.hpp:48
std::string name() const override
Definition Stokes.hpp:57
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
local assembly function that defines the bilinear form (LHS) computes and returns a single local stif...
Definition Stokes.hpp:62
void set_size(const int) override
Definition Stokes.hpp:69
bool is_fluid() const override
Definition Stokes.hpp:67
std::map< std::string, ParamFunc > parameters() const override
Definition Stokes.hpp:58
bool is_fluid() const override
Definition Stokes.hpp:30
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
local assembly function that defines the bilinear form (LHS) computes and returns a single local stif...
Definition Stokes.cpp:17
virtual std::string name() const override
Definition Stokes.hpp:27
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
Definition Stokes.cpp:38
bool is_tensor() const override
Definition Stokes.hpp:31
std::map< std::string, ParamFunc > parameters() const override
Definition Stokes.cpp:56
void add_multimaterial(const int index, const json &params, const Units &units) override
Definition Stokes.cpp:9
const GenericMatParam & viscosity() const
Definition Stokes.hpp:25
Used for test only.
Eigen::Matrix< AutodiffScalarHessian, Eigen::Dynamic, 1, 0, 3, 1 > AutodiffHessianPt
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:11
nlohmann::json json
Definition Common.hpp:9