PolyFEM
Loading...
Searching...
No Matches
State.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
4
5#include <polyfem/Units.hpp>
6
9
18
19#include <polyfem/mesh/Mesh.hpp>
23
25
33
34#include <polysolve/linear/Solver.hpp>
35
36#include <Eigen/Dense>
37#include <Eigen/Sparse>
38
39#include <spdlog/sinks/basic_file_sink.h>
40
41#include <ipc/collision_mesh.hpp>
42#include <ipc/utils/logger.hpp>
43
44#include <memory>
45#include <string>
46#include <unordered_map>
47#include <functional>
48#include <cassert>
49#include <map>
50#include <utility>
51#include <vector>
52#include <sstream>
53#include <utility>
54#include <algorithm>
55#include <cstddef>
56
57// Forward declaration
59{
60 class Solver;
61}
62
63namespace polyfem::assembler
64{
65 class Mass;
66 class ViscousDamping;
67 class ViscousDampingPrev;
68} // namespace polyfem::assembler
69
70namespace polyfem
71{
72 namespace mesh
73 {
74 class Mesh2D;
75 class Mesh3D;
76 } // namespace mesh
77
78 class State;
79
86 using UserPostStepCallback = std::function<void(int step, State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd *disp_grad, const Eigen::MatrixXd *pressure)>;
87
90 {
91 public:
94 Eigen::MatrixXd solution;
95
98 Eigen::MatrixXd velocity;
99
102 Eigen::MatrixXd acceleration;
103
105 bool is_empty() const
106 {
107 return solution.size() == 0 && velocity.size() == 0 && acceleration.size() == 0;
108 }
109 };
110
112 class State
113 {
114 public:
115 //---------------------------------------------------
116 //-----------------initialization--------------------
117 //---------------------------------------------------
118
119 ~State() = default;
121 State();
122
124 void set_max_threads(const int max_threads = std::numeric_limits<int>::max());
125
129 void init(const json &args, const bool strict_validation);
130
132 void init_time();
133
136
137 //---------------------------------------------------
138 //-----------------logger----------------------------
139 //---------------------------------------------------
140
146 void init_logger(
147 const std::string &log_file,
148 const spdlog::level::level_enum log_level,
149 const spdlog::level::level_enum file_log_level,
150 const bool is_quiet);
151
155 void init_logger(std::ostream &os, const spdlog::level::level_enum log_level);
156
159 void set_log_level(const spdlog::level::level_enum log_level);
160
164 std::string get_log(const Eigen::MatrixXd &sol)
165 {
166 std::stringstream ss;
167 save_json(sol, ss);
168 return ss.str();
169 }
170
171 private:
173 void init_logger(const std::vector<spdlog::sink_ptr> &sinks, const spdlog::level::level_enum log_level);
174
176 spdlog::sink_ptr console_sink_ = nullptr;
177 spdlog::sink_ptr file_sink_ = nullptr;
178
179 public:
180 //---------------------------------------------------
181 //-----------------assembly--------------------------
182 //---------------------------------------------------
183
185
187
189 std::shared_ptr<assembler::Assembler> assembler = nullptr;
190
191 std::shared_ptr<assembler::Mass> mass_matrix_assembler = nullptr;
192
193 std::shared_ptr<assembler::MixedAssembler> mixed_assembler = nullptr;
194 std::shared_ptr<assembler::Assembler> pressure_assembler = nullptr;
195
196 std::shared_ptr<assembler::PressureAssembler> elasticity_pressure_assembler = nullptr;
197
198 std::shared_ptr<assembler::ViscousDamping> damping_assembler = nullptr;
199 std::shared_ptr<assembler::ViscousDampingPrev> damping_prev_assembler = nullptr;
200
202 std::shared_ptr<assembler::Problem> problem;
203
205 std::vector<basis::ElementBases> bases;
207 std::vector<basis::ElementBases> pressure_bases;
209 std::vector<basis::ElementBases> geom_bases_;
210
217
219 std::map<int, Eigen::MatrixXd> polys;
221 std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>> polys_3d;
222
224 Eigen::VectorXi disc_orders, disc_ordersq;
225
227 std::shared_ptr<polyfem::mesh::MeshNodes> mesh_nodes, geom_mesh_nodes, pressure_mesh_nodes;
228
234
238 double avg_mass;
239
241 Eigen::MatrixXd rhs;
242
246
249 std::string formulation() const;
250
253 bool iso_parametric() const;
254
257 const std::vector<basis::ElementBases> &geom_bases() const
258 {
259 return iso_parametric() ? bases : geom_bases_;
260 }
261
266 void build_basis();
270 void assemble_rhs();
274 void assemble_mass_mat();
275
277 std::shared_ptr<assembler::RhsAssembler> build_rhs_assembler(
278 const int n_bases,
279 const std::vector<basis::ElementBases> &bases,
282 std::shared_ptr<assembler::RhsAssembler> build_rhs_assembler() const
283 {
285 }
286
287 std::shared_ptr<assembler::PressureAssembler> build_pressure_assembler(
288 const int n_bases_,
289 const std::vector<basis::ElementBases> &bases_) const;
290 std::shared_ptr<assembler::PressureAssembler> build_pressure_assembler() const
291 {
293 }
294
298 {
300 const int n_b_samples_j = args["space"]["advanced"]["n_boundary_samples"];
301 const int gdiscr_order = mesh->orders().size() <= 0 ? 1 : mesh->orders().maxCoeff();
302 const int discr_order = std::max(disc_orders.maxCoeff(), gdiscr_order);
303
304 const int n_b_samples = std::max(n_b_samples_j, AssemblerUtils::quadrature_order("Mass", discr_order, AssemblerUtils::BasisType::POLY, mesh->dimension()));
305 // todo prism
306 return {{n_b_samples, n_b_samples}};
307 }
308
309 int ndof() const
310 {
311 const int actual_dim = problem->is_scalar() ? 1 : mesh->dimension();
312 if (mixed_assembler == nullptr)
313 return actual_dim * n_bases;
314 else
315 return actual_dim * n_bases + n_pressure_bases;
316 }
317
318 private:
322 void sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure);
325
326 public:
329 void set_materials(std::vector<std::shared_ptr<assembler::Assembler>> &assemblers) const;
333
334 //---------------------------------------------------
335 //-----------------solver----------------------------
336 //---------------------------------------------------
337
338 public:
344 void solve_problem(Eigen::MatrixXd &sol,
345 Eigen::MatrixXd &pressure,
346 UserPostStepCallback user_post_step = {},
347 const InitialConditionOverride *ic_override = nullptr);
348
354 void solve(Eigen::MatrixXd &sol,
355 Eigen::MatrixXd &pressure,
356 UserPostStepCallback user_post_step = {},
357 const InitialConditionOverride *ic_override = nullptr)
358 {
359 if (!mesh)
360 {
361 logger().error("Load the mesh first!");
362 return;
363 }
365
366 build_basis();
367
368 assemble_rhs();
370
371 solve_problem(sol, pressure, user_post_step, ic_override);
372 }
373
376
382 std::unique_ptr<polysolve::linear::Solver> static_linear_solver_cache;
383
388 void init_solve(Eigen::MatrixXd &sol,
389 Eigen::MatrixXd &pressure,
390 const InitialConditionOverride *ic_override = nullptr);
391
398 void solve_transient_navier_stokes_split(const int time_steps,
399 const double dt,
400 Eigen::MatrixXd &sol,
401 Eigen::MatrixXd &pressure,
402 UserPostStepCallback user_post_step = {});
403
411 void solve_transient_navier_stokes(const int time_steps,
412 const double t0,
413 const double dt,
414 Eigen::MatrixXd &sol,
415 Eigen::MatrixXd &pressure,
416 UserPostStepCallback user_post_step = {});
417
426 void solve_transient_linear(const int time_steps,
427 const double t0,
428 const double dt,
429 Eigen::MatrixXd &sol,
430 Eigen::MatrixXd &pressure,
431 UserPostStepCallback user_post_step = {},
432 const InitialConditionOverride *ic_override = nullptr);
433
441 void solve_transient_tensor_nonlinear(const int time_steps,
442 const double t0,
443 const double dt,
444 Eigen::MatrixXd &sol,
445 UserPostStepCallback user_post_step = {},
446 const InitialConditionOverride *ic_override = nullptr);
447
452 void init_nonlinear_tensor_solve(Eigen::MatrixXd &sol,
453 const double t = 1.0,
454 const bool init_time_integrator = true,
455 const InitialConditionOverride *ic_override = nullptr);
456
460 void init_linear_solve(Eigen::MatrixXd &sol,
461 const double t = 1.0,
462 const InitialConditionOverride *ic_override = nullptr);
463
467 void initial_solution(Eigen::MatrixXd &solution, const InitialConditionOverride *ic_override = nullptr) const;
471 void initial_velocity(Eigen::MatrixXd &velocity, const InitialConditionOverride *ic_override = nullptr) const;
475 void initial_acceleration(Eigen::MatrixXd &acceleration, const InitialConditionOverride *ic_override = nullptr) const;
481 void solve_linear(int step, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step = {});
487 void solve_navier_stokes(int step, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step = {});
492 void solve_tensor_nonlinear(int step, Eigen::MatrixXd &sol, const bool init_lagging = true, UserPostStepCallback user_post_step = {});
493
496 std::shared_ptr<polysolve::nonlinear::Solver> make_nl_solver(bool for_al) const;
497
499 std::shared_ptr<utils::PeriodicBoundary> periodic_bc;
500 bool has_periodic_bc() const
501 {
502 return args["boundary_conditions"]["periodic_boundary"]["enabled"].get<bool>();
503 }
504
514 void solve_linear(
515 int step,
516 const std::unique_ptr<polysolve::linear::Solver> &solver,
518 Eigen::VectorXd &b,
519 const bool compute_spectrum,
520 Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step = {});
521
523 bool is_problem_linear() const { return assembler->is_linear() && !is_contact_enabled() && !is_pressure_enabled() && !has_constraints(); }
524
525 bool has_constraints() const
526 {
527 return has_constraints_;
528 }
529
530 private:
532
533 public:
536 void build_stiffness_mat(StiffnessMatrix &stiffness);
537
538 //---------------------------------------------------
539 //-----------------nodes flags-----------------------
540 //---------------------------------------------------
541
542 public:
544 std::unordered_map<int, std::array<bool, 3>>
545 boundary_conditions_ids(const std::string &bc_type) const;
546
548 std::vector<int> boundary_nodes;
550 std::vector<int> pressure_boundary_nodes;
552 std::vector<mesh::LocalBoundary> total_local_boundary;
554 std::vector<mesh::LocalBoundary> local_boundary;
556 std::vector<mesh::LocalBoundary> local_neumann_boundary;
558 std::vector<mesh::LocalBoundary> local_pressure_boundary;
560 std::unordered_map<int, std::vector<mesh::LocalBoundary>> local_pressure_cavity;
562 std::map<int, basis::InterfaceData> poly_edge_to_data;
564 std::vector<int> dirichlet_nodes;
565 std::vector<RowVectorNd> dirichlet_nodes_position;
567 std::vector<int> neumann_nodes;
568 std::vector<RowVectorNd> neumann_nodes_position;
569
571 Eigen::VectorXi in_node_to_node;
574
575 std::vector<int> primitive_to_node() const;
576 std::vector<int> node_to_primitive() const;
577
578 private:
580 void build_node_mapping();
581
582 //---------------------------------------------------
583 //-----------------Geometry--------------------------
584 //---------------------------------------------------
585 public:
587 std::unique_ptr<mesh::Mesh> mesh;
590
596 void load_mesh(bool non_conforming = false,
597 const std::vector<std::string> &names = std::vector<std::string>(),
598 const std::vector<Eigen::MatrixXi> &cells = std::vector<Eigen::MatrixXi>(),
599 const std::vector<Eigen::MatrixXd> &vertices = std::vector<Eigen::MatrixXd>());
600
606 void load_mesh(GEO::Mesh &meshin, const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &boundary_marker, bool non_conforming = false, bool skip_boundary_sideset = false);
607
612 void load_mesh(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, bool non_conforming = false)
613 {
614 mesh = mesh::Mesh::create(V, F, non_conforming);
615 load_mesh(non_conforming);
616 }
617
619 void reset_mesh();
620
622 void build_mesh_matrices(Eigen::MatrixXd &V, Eigen::MatrixXi &F);
623
629 bool remesh(const double time, const double dt, Eigen::MatrixXd &sol);
630
632 void get_vertices(Eigen::MatrixXd &vertices) const
633 {
634 vertices.setZero(mesh->n_vertices(), mesh->dimension());
635
636 for (int v = 0; v < mesh->n_vertices(); v++)
637 vertices.row(v) = mesh->point(v);
638 }
639
641 void get_elements(Eigen::MatrixXi &elements) const
642 {
643 assert(mesh->is_simplicial());
644
645 auto node_to_primitive_map = node_to_primitive();
646
647 const auto &gbases = geom_bases();
648 int dim = mesh->dimension();
649 elements.setZero(gbases.size(), dim + 1);
650 for (int e = 0; e < gbases.size(); e++)
651 {
652 int i = 0;
653 for (const auto &gbs : gbases[e].bases)
654 elements(e, i++) = node_to_primitive_map[gbs.global()[0].index];
655 }
656 }
657
658 //---------------------------------------------------
659 //-----------------IPC-------------------------------
660 //---------------------------------------------------
661
663 ipc::CollisionMesh collision_mesh;
664
666 ipc::CollisionMesh periodic_collision_mesh;
669
671 static void build_collision_mesh(
672 const mesh::Mesh &mesh,
673 const int n_bases,
674 const std::vector<basis::ElementBases> &bases,
675 const std::vector<basis::ElementBases> &geom_bases,
676 const std::vector<mesh::LocalBoundary> &total_local_boundary,
678 const json &args,
679 const std::function<std::string(const std::string &)> &resolve_input_path,
680 const Eigen::VectorXi &in_node_to_node,
681 ipc::CollisionMesh &collision_mesh);
682
686
690 bool is_obstacle_vertex(const size_t vi) const
691 {
692 // The obstalce vertices are at the bottom of the collision mesh vertices
693 return vi >= collision_mesh.full_num_vertices() - obstacle.n_vertices();
694 }
695
700 {
701 return args["contact"]["enabled"];
702 }
703
708 {
709 return args["contact"]["adhesion"]["adhesion_enabled"];
710 }
711
716 {
717 return (args["boundary_conditions"]["pressure_boundary"].size() > 0)
718 || (args["boundary_conditions"]["pressure_cavity"].size() > 0);
719 }
720
722 bool has_dhat = false;
723
724 //---------------------------------------------------
725 //-----------------OUTPUT----------------------------
726 //---------------------------------------------------
727 public:
729 std::string output_dir;
739
740 std::function<void(int, int, double, double)> time_callback = nullptr;
741
745 void export_data(const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure);
746
754 void save_timestep(const double time, const int t, const double t0, const double dt, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure);
755
761 void save_subsolve(const int i, const int t, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure);
762
766 void save_json(const Eigen::MatrixXd &sol, std::ostream &out);
767
770 void save_json(const Eigen::MatrixXd &sol);
771
773 void compute_errors(const Eigen::MatrixXd &sol);
774
777 void save_restart_json(const double t0, const double dt, const int t) const;
778
779 //-----------PATH management
782 std::string root_path() const;
783
788 std::string resolve_input_path(const std::string &path, const bool only_if_exists = false) const;
789
793 std::string resolve_output_path(const std::string &path) const;
794
795 //---------------------------------------------------
796 //-----------------differentiable--------------------
797 //---------------------------------------------------
798 public:
800
801 //---------------------------------------------------
802 //-----------------homogenization--------------------
803 //---------------------------------------------------
804 public:
806
808 void solve_homogenization_step(int step, Eigen::MatrixXd &sol, bool adaptive_initial_weight = false, UserPostStepCallback user_post_step = {}); // sol is the extended solution, i.e. [periodic fluctuation, macro strain]
809 void init_homogenization_solve(const double t);
810 void solve_homogenization(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, UserPostStepCallback user_post_step = {});
811 bool is_homogenization() const
812 {
813 return args["boundary_conditions"]["periodic_boundary"]["linear_displacement_offset"].size() > 0;
814 }
815 };
816
817} // namespace polyfem
int V
Runtime override for initial-condition histories.
Definition State.hpp:90
Eigen::MatrixXd solution
ndof by H (history size) matrix where each col is solution from a time step.
Definition State.hpp:94
bool is_empty() const
Returns true when no quantity is overridden.
Definition State.hpp:105
Eigen::MatrixXd acceleration
ndof by H (history size) matrix where each col is acceleration from a time step.
Definition State.hpp:102
Eigen::MatrixXd velocity
ndof by H (history size) matrix where each col is velocity from a time step.
Definition State.hpp:98
main class that contains the polyfem solver and all its state
Definition State.hpp:113
void get_vertices(Eigen::MatrixXd &vertices) const
Gather geometry vertices into a dense matrix.
Definition State.hpp:632
io::OutRuntimeData timings
runtime statistics
Definition State.hpp:733
void solve_navier_stokes(int step, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={})
solves a navier stokes
assembler::MacroStrainValue macro_strain_constraint
Definition State.hpp:805
std::shared_ptr< utils::PeriodicBoundary > periodic_bc
periodic BC and periodic mesh utils
Definition State.hpp:499
int n_bases
number of bases
Definition State.hpp:212
int n_pressure_bases
number of pressure bases
Definition State.hpp:214
std::string formulation() const
return the formulation (checks if the problem is scalar or not and deals with multiphysics)
Definition State.cpp:332
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:230
bool is_adhesion_enabled() const
does the simulation have adhesion
Definition State.hpp:707
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:257
double starting_max_edge_length
Definition State.hpp:737
std::vector< mesh::LocalBoundary > local_pressure_boundary
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:558
void sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
splits the solution in solution and pressure for mixed problems
Definition State.cpp:375
std::shared_ptr< assembler::ViscousDamping > damping_assembler
Definition State.hpp:198
void init(const json &args, const bool strict_validation)
initialize the polyfem solver with a json settings
std::string root_path() const
Get the root path for the state (e.g., args["root_path"] or ".")
bool is_pressure_enabled() const
does the simulation has pressure
Definition State.hpp:715
void assemble_mass_mat()
assemble mass, step 4 of solve build mass matrix based on defined basis modifies mass (and maybe more...
Definition State.cpp:1461
Eigen::VectorXi in_node_to_node
Inpute nodes (including high-order) to polyfem nodes, only for isoparametric.
Definition State.hpp:571
void solve_linear(int step, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={})
solves a linear problem
void build_polygonal_basis()
builds bases for polygons, called inside build_basis
Definition State.cpp:1029
mesh::Obstacle obstacle
Obstacles used in collisions.
Definition State.hpp:589
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:189
void compute_errors(const Eigen::MatrixXd &sol)
computes all errors
void initial_velocity(Eigen::MatrixXd &velocity, const InitialConditionOverride *ic_override=nullptr) const
Load or compute the initial velocity.
ipc::CollisionMesh periodic_collision_mesh
IPC collision mesh under periodic BC.
Definition State.hpp:666
void build_periodic_collision_mesh()
Definition State.cpp:1192
void save_subsolve(const int i, const int t, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves a subsolve when save_solve_sequence_debug is true
io::OutGeometryData out_geom
visualization stuff
Definition State.hpp:731
bool use_avg_pressure
use average pressure for stokes problem to fix the additional dofs, true by default if false,...
Definition State.hpp:245
std::string resolve_output_path(const std::string &path) const
Resolve output path relative to output_dir if the path is not absolute.
std::vector< int > neumann_nodes
per node neumann
Definition State.hpp:567
ipc::CollisionMesh collision_mesh
IPC collision mesh.
Definition State.hpp:663
std::string output_dir
Directory for output files.
Definition State.hpp:729
void load_mesh(bool non_conforming=false, const std::vector< std::string > &names=std::vector< std::string >(), const std::vector< Eigen::MatrixXi > &cells=std::vector< Eigen::MatrixXi >(), const std::vector< Eigen::MatrixXd > &vertices=std::vector< Eigen::MatrixXd >())
loads the mesh from the json arguments
Definition StateLoad.cpp:91
std::vector< basis::ElementBases > geom_bases_
Geometric mapping bases, if the elements are isoparametric, this list is empty.
Definition State.hpp:209
Eigen::VectorXi in_primitive_to_primitive
maps in vertices/edges/faces/cells to polyfem vertices/edges/faces/cells
Definition State.hpp:573
void set_materials(std::vector< std::shared_ptr< assembler::Assembler > > &assemblers) const
set the material and the problem dimension
std::string get_log(const Eigen::MatrixXd &sol)
gets the output log as json this is not what gets printed but more informative information,...
Definition State.hpp:164
std::map< int, basis::InterfaceData > poly_edge_to_data
nodes on the boundary of polygonal elements, used for harmonic bases
Definition State.hpp:562
std::vector< int > pressure_boundary_nodes
list of neumann boundary nodes
Definition State.hpp:550
void set_max_threads(const int max_threads=std::numeric_limits< int >::max())
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:207
void save_json(const Eigen::MatrixXd &sol, std::ostream &out)
saves the output statistic to a stream
bool is_homogenization() const
Definition State.hpp:811
void load_mesh(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, bool non_conforming=false)
loads the mesh from V and F,
Definition State.hpp:612
void solve_homogenization_step(int step, Eigen::MatrixXd &sol, bool adaptive_initial_weight=false, UserPostStepCallback user_post_step={})
In Elasticity PDE, solve for "min W(disp_grad + \grad u)" instead of "min W(\grad u)".
bool has_constraints() const
Definition State.hpp:525
std::shared_ptr< assembler::Assembler > pressure_assembler
Definition State.hpp:194
std::shared_ptr< polysolve::nonlinear::Solver > make_nl_solver(bool for_al) const
factory to create the nl solver depending on input
double starting_min_edge_length
Definition State.hpp:736
void solve_transient_tensor_nonlinear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves transient tensor nonlinear problem
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:191
void build_node_mapping()
build the mapping from input nodes to polyfem nodes
Definition State.cpp:249
bool remesh(const double time, const double dt, Eigen::MatrixXd &sol)
Remesh the FE space and update solution(s).
std::vector< int > primitive_to_node() const
Definition State.cpp:232
bool has_dhat
stores if input json contains dhat
Definition State.hpp:722
void solve_tensor_nonlinear(int step, Eigen::MatrixXd &sol, const bool init_lagging=true, UserPostStepCallback user_post_step={})
solves nonlinear problems
std::shared_ptr< polyfem::mesh::MeshNodes > pressure_mesh_nodes
Definition State.hpp:227
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:587
void reset_mesh()
Resets the mesh.
Definition StateLoad.cpp:22
std::unordered_map< int, std::array< bool, 3 > > boundary_conditions_ids(const std::string &bc_type) const
Construct a vector of boundary conditions ids with their dimension flags.
std::shared_ptr< polyfem::mesh::MeshNodes > mesh_nodes
Mapping from input nodes to FE nodes.
Definition State.hpp:227
StiffnessMatrix mass
Mass matrix, it is computed only for time dependent problems.
Definition State.hpp:236
void solve(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves the problem, call other methods
Definition State.hpp:354
std::vector< int > dirichlet_nodes
per node dirichlet
Definition State.hpp:564
bool has_periodic_bc() const
Definition State.hpp:500
void build_basis()
builds the bases step 2 of solve modifies bases, pressure_bases, geom_bases_, boundary_nodes,...
Definition State.cpp:507
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:202
bool has_constraints_
Definition State.hpp:531
std::vector< RowVectorNd > dirichlet_nodes_position
Definition State.hpp:565
spdlog::sink_ptr file_sink_
Definition State.hpp:177
spdlog::sink_ptr console_sink_
logger sink to stdout
Definition State.hpp:176
State()
Constructor.
Definition StateInit.cpp:56
std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > polys_3d
polyhedra, used since poly have no geom mapping
Definition State.hpp:221
std::vector< int > node_to_primitive() const
Definition State.cpp:239
~State()=default
json args
main input arguments containing all defaults
Definition State.hpp:135
void set_log_level(const spdlog::level::level_enum log_level)
change log level
void save_restart_json(const double t0, const double dt, const int t) const
Save a JSON sim file for restarting the simulation at time t.
void initial_acceleration(Eigen::MatrixXd &acceleration, const InitialConditionOverride *ic_override=nullptr) const
Load or compute the initial acceleration.
io::OutStatsData stats
Other statistics.
Definition State.hpp:735
std::vector< RowVectorNd > neumann_nodes_position
Definition State.hpp:568
void assemble_rhs()
compute rhs, step 3 of solve build rhs vector based on defined basis and given rhs of the problem mod...
Definition State.cpp:1585
double min_boundary_edge_length
Definition State.hpp:738
std::unique_ptr< polysolve::linear::Solver > static_linear_solver_cache
Linear solver instance from the most recent static linear solve.
Definition State.hpp:382
Eigen::VectorXi periodic_collision_mesh_to_basis
index mapping from periodic 2x2 collision mesh to FE periodic mesh
Definition State.hpp:668
void init_time()
initialize time settings if args contains "time"
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:205
void initial_solution(Eigen::MatrixXd &solution, const InitialConditionOverride *ic_override=nullptr) const
Load or compute the initial solution.
Eigen::VectorXi disc_ordersq
Definition State.hpp:224
void get_elements(Eigen::MatrixXi &elements) const
Gather geometry elements into a dense matrix.
Definition State.hpp:641
void solve_transient_navier_stokes_split(const int time_steps, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={})
solves transient navier stokes with operator splitting
int n_geom_bases
number of geometric bases
Definition State.hpp:216
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves the problems
Definition State.cpp:1657
assembler::AssemblyValsCache pressure_ass_vals_cache
used to store assembly values for pressure for small problems
Definition State.hpp:233
void init_solve(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, const InitialConditionOverride *ic_override=nullptr)
initialize solver
void solve_transient_linear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves transient linear problem
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:554
double avg_mass
average system mass, used for contact with IPC
Definition State.hpp:238
bool optimization_enabled
Definition State.hpp:799
std::function< void(int, int, double, double)> time_callback
Definition State.hpp:740
bool iso_parametric() const
check if using iso parametric bases
Definition State.cpp:466
std::shared_ptr< assembler::RhsAssembler > build_rhs_assembler() const
build a RhsAssembler for the problem
Definition State.hpp:282
std::vector< mesh::LocalBoundary > total_local_boundary
mapping from elements to nodes for all mesh
Definition State.hpp:552
QuadratureOrders n_boundary_samples() const
quadrature used for projecting boundary conditions
Definition State.hpp:297
int ndof() const
Definition State.hpp:309
void export_data(const Eigen::MatrixXd &sol, const Eigen::MatrixXd &pressure)
saves all data on the disk according to the input params
assembler::AssemblyValsCache mass_ass_vals_cache
Definition State.hpp:231
bool is_obstacle_vertex(const size_t vi) const
checks if vertex is obstacle
Definition State.hpp:690
std::string resolve_input_path(const std::string &path, const bool only_if_exists=false) const
Resolve input path relative to root_path() if the path is not absolute.
std::vector< int > boundary_nodes
list of boundary nodes
Definition State.hpp:548
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:375
void init_nonlinear_tensor_solve(Eigen::MatrixXd &sol, const double t=1.0, const bool init_time_integrator=true, const InitialConditionOverride *ic_override=nullptr)
initialize the nonlinear solver
void solve_transient_navier_stokes(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={})
solves transient navier stokes with FEM
std::shared_ptr< polyfem::mesh::MeshNodes > geom_mesh_nodes
Definition State.hpp:227
void init_logger(const std::string &log_file, const spdlog::level::level_enum log_level, const spdlog::level::level_enum file_log_level, const bool is_quiet)
initializing the logger
Definition StateInit.cpp:65
std::shared_ptr< assembler::ViscousDampingPrev > damping_prev_assembler
Definition State.hpp:199
bool is_contact_enabled() const
does the simulation have contact
Definition State.hpp:699
void build_collision_mesh()
extracts the boundary mesh for collision, called in build_basis
Definition State.cpp:1326
Eigen::VectorXi disc_orders
vector of discretization orders, used when not all elements have the same degree, one per element
Definition State.hpp:224
std::vector< mesh::LocalBoundary > local_neumann_boundary
mapping from elements to nodes for neumann boundary conditions
Definition State.hpp:556
std::shared_ptr< assembler::PressureAssembler > build_pressure_assembler() const
Definition State.hpp:290
void solve_homogenization(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, UserPostStepCallback user_post_step={})
void init_linear_solve(Eigen::MatrixXd &sol, const double t=1.0, const InitialConditionOverride *ic_override=nullptr)
initialize the linear solve
std::shared_ptr< assembler::PressureAssembler > elasticity_pressure_assembler
Definition State.hpp:196
void build_mesh_matrices(Eigen::MatrixXd &V, Eigen::MatrixXi &F)
Build the mesh matrices (vertices and elements) from the mesh using the bases node ordering.
bool is_problem_linear() const
Returns whether the system is linear. Collisions and pressure add nonlinearity to the problem.
Definition State.hpp:523
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
Definition State.hpp:193
std::map< int, Eigen::MatrixXd > polys
polygons, used since poly have no geom mapping
Definition State.hpp:219
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:241
std::unordered_map< int, std::vector< mesh::LocalBoundary > > local_pressure_cavity
mapping from elements to nodes for pressure boundary conditions
Definition State.hpp:560
void init_homogenization_solve(const double t)
static int quadrature_order(const std::string &assembler, const int basis_degree, const BasisType &b_type, const int dim)
utility for retrieving the needed quadrature order to precisely integrate the given form on the given...
Caches basis evaluation and geometric mapping at every element.
Utilies related to export of geometry.
Definition OutData.hpp:33
stores all runtime data
Definition OutData.hpp:348
all stats from polyfem
Definition OutData.hpp:375
void compute_mesh_stats(const polyfem::mesh::Mesh &mesh)
compute stats (counts els type, mesh lenght, etc), step 1 of solve
Definition OutData.cpp:3025
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:40
static std::unique_ptr< Mesh > create(const std::string &path, const bool non_conforming=false)
factory to build the proper mesh
Definition Mesh.cpp:173
class to store time stepping data
Definition SolveData.hpp:59
Used for test only.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
std::array< int, 2 > QuadratureOrders
Definition Types.hpp:19
nlohmann::json json
Definition Common.hpp:9
std::function< void(int step, State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd *disp_grad, const Eigen::MatrixXd *pressure)> UserPostStepCallback
User callback at the end of every solver step.
Definition State.hpp:86
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:13
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24