PolyFEM
Loading...
Searching...
No Matches
StateSolveNavierStokes.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2
7
14
15namespace polyfem
16{
17 using namespace solver;
18 using namespace time_integrator;
19
20 void State::solve_navier_stokes(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
21 {
22 assert(!problem->is_time_dependent());
23 assert(assembler->name() == "NavierStokes");
24
25 assert(solve_data.rhs_assembler != nullptr);
28
29 std::shared_ptr<assembler::Assembler> velocity_stokes_assembler = std::make_shared<assembler::StokesVelocity>();
30 set_materials(*velocity_stokes_assembler);
31
32 Eigen::VectorXd x;
33 solver::NavierStokesSolver ns_solver(args["solver"]);
36 geom_bases(),
37 *velocity_stokes_assembler,
38 *dynamic_cast<assembler::NavierStokesVelocity *>(assembler.get()),
45 mesh->dimension(),
46 mesh->is_volume(),
47 rhs, x);
48
49 sol = x;
50 sol_to_pressure(sol, pressure);
51 }
52
53 void State::solve_transient_navier_stokes_split(const int time_steps, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
54 {
55 assert(assembler->name() == "OperatorSplitting" && problem->is_time_dependent());
56
57 Eigen::MatrixXd local_pts;
58 auto &gbases = geom_bases();
59 if (mesh->dimension() == 2)
60 {
61 if (gbases[0].bases.size() == 3)
62 autogen::p_nodes_2d(args["space"]["discr_order"], local_pts);
63 else
64 autogen::q_nodes_2d(args["space"]["discr_order"], local_pts);
65 }
66 else
67 {
68 if (gbases[0].bases.size() == 4)
69 autogen::p_nodes_3d(args["space"]["discr_order"], local_pts);
70 else
71 autogen::q_nodes_3d(args["space"]["discr_order"], local_pts);
72 }
73
74 std::vector<int> bnd_nodes;
75 bnd_nodes.reserve(boundary_nodes.size() / mesh->dimension());
76 for (auto it = boundary_nodes.begin(); it != boundary_nodes.end(); it++)
77 {
78 if (!(*it % mesh->dimension()))
79 continue;
80 bnd_nodes.push_back(*it / mesh->dimension());
81 }
82
83 const int dim = mesh->dimension();
84 const int n_el = int(bases.size()); // number of elements
85 const int shape = gbases[0].bases.size(); // number of geometry vertices in an element
86
87 auto fluid_assembler = std::dynamic_pointer_cast<assembler::OperatorSplitting>(assembler);
88 if (!fluid_assembler)
89 log_and_throw_error("Invalid assembler {}!", assembler->name());
90 // TODO
91 double viscosity_ = fluid_assembler->viscosity()(0, 0, 0, 0, 0);
92 assert(viscosity_ >= 0);
93
94 logger().info("Matrices assembly...");
95 StiffnessMatrix stiffness_viscosity, mixed_stiffness, velocity_mass, stiffness;
96
97 build_stiffness_mat(stiffness);
98 // coefficient matrix of viscosity
99 assembler::Laplacian lapl_assembler;
100 lapl_assembler.set_size(1);
101 lapl_assembler.assemble(mesh->is_volume(), n_bases, bases, gbases, ass_vals_cache, 0, stiffness_viscosity);
102 mass_matrix_assembler->set_size(1);
103 mass_matrix_assembler->assemble(mesh->is_volume(), n_bases, bases, gbases, mass_ass_vals_cache, 0, mass, true);
104
105 // coefficient matrix of pressure projection
106 lapl_assembler.assemble(mesh->is_volume(), n_pressure_bases, pressure_bases, gbases, pressure_ass_vals_cache, 0, stiffness);
107
108 // matrix used to calculate divergence of velocity
109 mixed_assembler->assemble(mesh->is_volume(), n_pressure_bases, n_bases, pressure_bases, bases, gbases,
110 pressure_ass_vals_cache, ass_vals_cache, 0, mixed_stiffness);
111 mass_matrix_assembler->set_size(mesh->dimension());
112 mass_matrix_assembler->assemble(mesh->is_volume(), n_bases, bases, gbases, mass_ass_vals_cache, 0, velocity_mass, true);
113 mixed_stiffness = mixed_stiffness.transpose();
114 logger().info("Matrices assembly ends!");
115
118 stiffness_viscosity, stiffness, velocity_mass, dt, viscosity_, args["solver"]["linear"]);
119
120 /* initialize solution */
121 pressure = Eigen::MatrixXd::Zero(n_pressure_bases, 1);
122
123 const int n_b_samples = n_boundary_samples();
124 for (int t = 1; t <= time_steps; t++)
125 {
126 double time = t * dt;
127 logger().info("{}/{} steps, t={}s", t, time_steps, time);
128
129 /* advection */
130 logger().info("Advection...");
131 if (args["space"]["advanced"]["use_particle_advection"])
132 ss.advection_FLIP(*mesh, gbases, bases, sol, dt, local_pts);
133 else
134 ss.advection(*mesh, gbases, bases, sol, dt, local_pts);
135 logger().info("Advection finished!");
136
137 /* apply boundary condition */
139 local_boundary, boundary_nodes, n_b_samples, local_neumann_boundary, sol, Eigen::MatrixXd(), time);
140
141 /* viscosity */
142 logger().info("Solving diffusion...");
143 if (viscosity_ > 0)
144 ss.solve_diffusion_1st(mass, bnd_nodes, sol);
145 logger().info("Diffusion solved!");
146
147 /* external force */
148 ss.external_force(*mesh, *assembler, gbases, bases, dt, sol, local_pts, problem, time);
149
150 /* incompressibility */
151 logger().info("Pressure projection...");
152 ss.solve_pressure(mixed_stiffness, pressure_boundary_nodes, sol, pressure);
153
154 ss.projection(n_bases, gbases, bases, pressure_bases, local_pts, pressure, sol);
155 // ss.projection(velocity_mass, mixed_stiffness, boundary_nodes, sol, pressure);
156 logger().info("Pressure projection finished!");
157
158 pressure = pressure / dt;
159
160 /* apply boundary condition */
162 local_boundary, boundary_nodes, n_b_samples, local_neumann_boundary, sol, Eigen::MatrixXd(), time);
163
164 /* export to vtu */
165 save_timestep(time, t, 0, dt, sol, pressure);
166 }
167 }
168
169 void State::solve_transient_navier_stokes(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
170 {
171 assert(assembler->name() == "NavierStokes" && problem->is_time_dependent());
172
173 const auto &gbases = geom_bases();
174 Eigen::MatrixXd current_rhs = rhs;
175
176 StiffnessMatrix velocity_mass;
177 // TODO mass might not be constant
178 mass_matrix_assembler->assemble(mesh->is_volume(), n_bases, bases, gbases, mass_ass_vals_cache, 0, velocity_mass, true);
179
180 StiffnessMatrix velocity_stiffness, mixed_stiffness, pressure_stiffness;
181
182 Eigen::VectorXd prev_sol;
183
184 BDF time_integrator;
185 if (args["time"]["integrator"].is_object() && args["time"]["integrator"]["type"] == "BDF")
186 time_integrator.set_parameters(args["time"]["integrator"]);
187 time_integrator.init(sol, Eigen::VectorXd::Zero(sol.size()), Eigen::VectorXd::Zero(sol.size()), dt);
188
189 std::shared_ptr<assembler::Assembler> velocity_stokes_assembler = std::make_shared<assembler::StokesVelocity>();
190 set_materials(*velocity_stokes_assembler);
191
192 mixed_assembler->assemble(mesh->is_volume(), n_pressure_bases, n_bases, pressure_bases, bases, gbases,
193 pressure_ass_vals_cache, ass_vals_cache, 0, mixed_stiffness);
195 pressure_stiffness);
196
197 solver::TransientNavierStokesSolver ns_solver(args["solver"]);
198 const int n_larger = n_pressure_bases + (use_avg_pressure ? 1 : 0);
199
200 const int n_b_samples = n_boundary_samples();
201
202 for (int t = 1; t <= time_steps; ++t)
203 {
204 double time = t0 + t * dt;
205 double current_dt = dt;
206
207 velocity_stokes_assembler->assemble(mesh->is_volume(), n_bases, bases, gbases, ass_vals_cache, time, velocity_stiffness);
208
209 logger().info("{}/{} steps, dt={}s t={}s", t, time_steps, current_dt, time);
210
211 prev_sol = time_integrator.weighted_sum_x_prevs();
212 solve_data.rhs_assembler->compute_energy_grad(
213 local_boundary, boundary_nodes, mass_matrix_assembler->density(), n_b_samples, local_neumann_boundary, rhs, time, current_rhs);
215 local_boundary, boundary_nodes, n_b_samples, local_neumann_boundary, current_rhs, Eigen::MatrixXd(), time);
216
217 const int prev_size = current_rhs.size();
218 if (prev_size != rhs.size())
219 {
220 current_rhs.conservativeResize(prev_size + n_larger, current_rhs.cols());
221 current_rhs.block(prev_size, 0, n_larger, current_rhs.cols()).setZero();
222 }
223
224 Eigen::VectorXd tmp_sol;
225 ns_solver.minimize(
227 time,
228 bases, geom_bases(),
229 *dynamic_cast<assembler::NavierStokesVelocity *>(assembler.get()),
233 mesh->dimension(),
234 mesh->is_volume(),
235 sqrt(time_integrator.acceleration_scaling()), prev_sol, velocity_stiffness, mixed_stiffness,
236 pressure_stiffness, velocity_mass, current_rhs, tmp_sol);
237 sol = tmp_sol;
238 time_integrator.update_quantities(sol.topRows(n_bases * mesh->dimension()));
239 sol_to_pressure(sol, pressure);
240
241 save_timestep(time, t, t0, dt, sol, pressure);
242 }
243 }
244} // namespace polyfem
int x
int n_bases
number of bases
Definition State.hpp:178
int n_pressure_bases
number of pressure bases
Definition State.hpp:180
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:196
int n_boundary_samples() const
quadrature used for projecting boundary conditions
Definition State.hpp:263
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:223
void sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
splits the solution in solution and pressure for mixed problems
Definition State.cpp:445
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:155
bool use_avg_pressure
use average pressure for stokes problem to fix the additional dofs, true by default if false,...
Definition State.hpp:211
void set_materials(std::vector< std::shared_ptr< assembler::Assembler > > &assemblers) const
set the material and the problem dimension
std::vector< int > pressure_boundary_nodes
list of neumann boundary nodes
Definition State.hpp:429
void save_timestep(const double time, const int t, const double t0, const double dt, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves a timestep
std::vector< basis::ElementBases > pressure_bases
FE pressure bases for mixed elements, the size is #elements.
Definition State.hpp:173
std::shared_ptr< assembler::Assembler > pressure_assembler
Definition State.hpp:160
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:157
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
StiffnessMatrix mass
Mass matrix, it is computed only for time dependent problems.
Definition State.hpp:202
void solve_navier_stokes(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves a navier stokes
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
json args
main input arguments containing all defaults
Definition State.hpp:101
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
void solve_transient_navier_stokes_split(const int time_steps, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves transient navier stokes with operator splitting
assembler::AssemblyValsCache pressure_ass_vals_cache
used to store assembly values for pressure for small problems
Definition State.hpp:199
void build_stiffness_mat(StiffnessMatrix &stiffness)
utility that builds the stiffness matrix and collects stats, used only for linear problems
std::vector< mesh::LocalBoundary > local_boundary
mapping from elements to nodes for dirichlet boundary conditions
Definition State.hpp:433
void solve_transient_navier_stokes(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves transient navier stokes with FEM
assembler::AssemblyValsCache mass_ass_vals_cache
Definition State.hpp:197
std::vector< int > boundary_nodes
list of boundary nodes
Definition State.hpp:427
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:324
std::vector< mesh::LocalBoundary > local_neumann_boundary
mapping from elements to nodes for neumann boundary conditions
Definition State.hpp:435
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
Definition State.hpp:159
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:207
virtual void set_size(const int size)
Definition Assembler.hpp:64
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 9, 1 > assemble(const LinearAssemblerData &data) const override
computes local stiffness matrix (1x1) for bases i,j where i,j is passed in through data ie integral o...
Definition Laplacian.cpp:13
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)
void solve_pressure(const StiffnessMatrix &mixed_stiffness, const std::vector< int > &pressure_boundary_nodes, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
void projection(const StiffnessMatrix &velocity_mass, const StiffnessMatrix &mixed_stiffness, const std::vector< int > &boundary_nodes_, Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
void external_force(const mesh::Mesh &mesh, const assembler::Assembler &assembler, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, const double dt, Eigen::MatrixXd &sol, const Eigen::MatrixXd &local_pts, const std::shared_ptr< assembler::Problem > problem, const double time)
void advection(const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, Eigen::MatrixXd &sol, const double dt, const Eigen::MatrixXd &local_pts, const int order=1, const int RK=1)
void advection_FLIP(const mesh::Mesh &mesh, const std::vector< basis::ElementBases > &gbases, const std::vector< basis::ElementBases > &bases, Eigen::MatrixXd &sol, const double dt, const Eigen::MatrixXd &local_pts, const int order=1)
void solve_diffusion_1st(const StiffnessMatrix &mass, const std::vector< int > &bnd_nodes, Eigen::MatrixXd &sol)
std::shared_ptr< assembler::RhsAssembler > rhs_assembler
void minimize(const int n_bases, const int n_pressure_bases, const double t, 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 double beta_dt, const Eigen::VectorXd &prev_sol, const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness, const StiffnessMatrix &velocity_mass1, const Eigen::MatrixXd &rhs, Eigen::VectorXd &x)
Backward Differential Formulas.
Definition BDF.hpp:14
void set_parameters(const json &params) override
Set the number of steps parameters from a json object.
Definition BDF.cpp:14
Eigen::VectorXd weighted_sum_x_prevs() const
Compute the weighted sum of the previous solutions.
Definition BDF.cpp:49
void update_quantities(const Eigen::VectorXd &x) override
Update the time integration quantities (i.e., , , and ).
Definition BDF.cpp:75
double acceleration_scaling() const override
Compute the acceleration scaling used to scale forces when integrating a second order ODE.
Definition BDF.cpp:110
virtual void init(const Eigen::MatrixXd &x_prevs, const Eigen::MatrixXd &v_prevs, const Eigen::MatrixXd &a_prevs, double dt)
Initialize the time integrator with the previous values for , , and .
void q_nodes_2d(const int q, Eigen::MatrixXd &val)
void p_nodes_2d(const int p, Eigen::MatrixXd &val)
void p_nodes_3d(const int p, Eigen::MatrixXd &val)
void q_nodes_3d(const int q, Eigen::MatrixXd &val)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20