PolyFEM
Loading...
Searching...
No Matches
NavierStokes.hpp
Go to the documentation of this file.
1#pragma once
2
5
7
8// Navier-Stokes local assembler
9namespace polyfem::assembler
10{
12 {
13 public:
17
19
20 std::string name() const override { return "NavierStokes"; }
21 std::map<std::string, ParamFunc> parameters() const override;
22
23 // navier stokes is not energy based
24 double compute_energy(const NonLinearAssemblerData &data) const override
25 {
26 // not used, this formulation is gradient based!
27 assert(false);
28 return 0;
29 }
30
31 // res is R^{dimĀ²}
32 // pde
33 Eigen::VectorXd
34 assemble_gradient(const NonLinearAssemblerData &data) const override;
35
36 Eigen::MatrixXd
37 // gradient of pde, this returns full gradient or partil depending on the template
38 assemble_hessian(const NonLinearAssemblerData &data) const override;
39
40 // rhs for fabbricated solution
41 VectorNd compute_rhs(const AutodiffHessianPt &pt) const override;
42
43 // set viscosity
44 void add_multimaterial(const int index, const json &params, const Units &units) override;
45
46 bool is_fluid() const override { return true; }
47 bool is_tensor() const override { return true; }
48
49 void set_picard(const bool val) { full_gradient_ = !val; }
50
51 private:
53
54 // not full graidnet used for Picard iteration
55 bool full_gradient_ = true;
56
57 Eigen::MatrixXd compute_N(const NonLinearAssemblerData &data) const;
58 Eigen::MatrixXd compute_W(const NonLinearAssemblerData &data) const;
59 };
60
61} // namespace polyfem::assembler
double val
Definition Assembler.cpp:86
double assemble_energy(const bool is_volume, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev) const override
void assemble_gradient(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, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, Eigen::MatrixXd &rhs) const override
void assemble_hessian(const bool is_volume, const int n_basis, const bool project_to_psd, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const AssemblyValsCache &cache, const double t, const double dt, const Eigen::MatrixXd &displacement, const Eigen::MatrixXd &displacement_prev, utils::MatrixCache &mat_cache, StiffnessMatrix &grad) const override
Eigen::MatrixXd assemble_hessian(const NonLinearAssemblerData &data) const override
Eigen::VectorXd assemble_gradient(const NonLinearAssemblerData &data) const override
void add_multimaterial(const int index, const json &params, const Units &units) override
double compute_energy(const NonLinearAssemblerData &data) const override
VectorNd compute_rhs(const AutodiffHessianPt &pt) const override
std::map< std::string, ParamFunc > parameters() const override
Eigen::MatrixXd compute_N(const NonLinearAssemblerData &data) const
std::string name() const override
Eigen::MatrixXd compute_W(const NonLinearAssemblerData &data) const
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