PolyFEM
Loading...
Searching...
No Matches
NavierStokesSolver.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
7
8#include <polysolve/linear/Solver.hpp>
9
11
12#include <memory>
13
14namespace polyfem
15{
16 namespace solver
17 {
19 {
20 public:
22
23 void minimize(const int n_bases,
24 const int n_pressure_bases,
25 const std::vector<basis::ElementBases> &bases,
26 const std::vector<basis::ElementBases> &pressure_bases,
27 const std::vector<basis::ElementBases> &gbases,
28 const assembler::Assembler &velocity_stokes_assembler,
29 assembler::NavierStokesVelocity &velocity_assembler,
30 const assembler::MixedAssembler &mixed_assembler,
31 const assembler::Assembler &pressure_assembler,
32 const assembler::AssemblyValsCache &ass_vals_cache,
33 const assembler::AssemblyValsCache &pressure_ass_vals_cache,
34 const std::vector<int> &boundary_nodes,
35 const bool use_avg_pressure,
36 const int problem_dim,
37 const bool is_volume,
38 const Eigen::MatrixXd &rhs,
39 Eigen::VectorXd &x);
40
41 void get_info(json &params)
42 {
43 params = solver_info;
44 }
45
46 int error_code() const { return 0; }
47
48 private:
49 int minimize_aux(
50 const bool is_picard,
51 const std::vector<int> &skipping,
52 const int n_bases,
53 const int n_pressure_bases,
54 const std::vector<basis::ElementBases> &bases,
55 const std::vector<basis::ElementBases> &gbases,
56 assembler::NavierStokesVelocity &velocity_assembler,
57 const assembler::AssemblyValsCache &ass_vals_cache,
58 const std::vector<int> &boundary_nodes,
59 const bool use_avg_pressure,
60 const int problem_dim,
61 const bool is_volume,
62 const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness,
63 const Eigen::VectorXd &rhs, const double grad_norm,
64 std::unique_ptr<polysolve::linear::Solver> &solver, double &nlres_norm,
65 Eigen::VectorXd &x);
66
68
69 double gradNorm;
71
73
74 json internal_solver = json::array();
75
80
81 bool has_nans(const polyfem::StiffnessMatrix &hessian);
82 };
83 } // namespace solver
84} // namespace polyfem
int x
Caches basis evaluation and geometric mapping at every element.
void minimize(const int n_bases, const int n_pressure_bases, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &pressure_bases, const std::vector< basis::ElementBases > &gbases, const assembler::Assembler &velocity_stokes_assembler, assembler::NavierStokesVelocity &velocity_assembler, const assembler::MixedAssembler &mixed_assembler, const assembler::Assembler &pressure_assembler, const assembler::AssemblyValsCache &ass_vals_cache, const assembler::AssemblyValsCache &pressure_ass_vals_cache, const std::vector< int > &boundary_nodes, const bool use_avg_pressure, const int problem_dim, const bool is_volume, const Eigen::MatrixXd &rhs, Eigen::VectorXd &x)
int minimize_aux(const bool is_picard, const std::vector< int > &skipping, const int n_bases, const int n_pressure_bases, const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, assembler::NavierStokesVelocity &velocity_assembler, const assembler::AssemblyValsCache &ass_vals_cache, const std::vector< int > &boundary_nodes, const bool use_avg_pressure, const int problem_dim, const bool is_volume, const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness, const Eigen::VectorXd &rhs, const double grad_norm, std::unique_ptr< polysolve::linear::Solver > &solver, double &nlres_norm, Eigen::VectorXd &x)
bool has_nans(const polyfem::StiffnessMatrix &hessian)
nlohmann::json json
Definition Common.hpp:9
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22