34#include <polysolve/linear/Solver.hpp>
37#include <Eigen/Sparse>
39#include <spdlog/sinks/basic_file_sink.h>
41#include <ipc/collision_mesh.hpp>
42#include <ipc/utils/logger.hpp>
46#include <unordered_map>
67 class ViscousDampingPrev;
86 using UserPostStepCallback = std::function<void(
int step,
State &state,
const Eigen::MatrixXd &sol,
const Eigen::MatrixXd *disp_grad,
const Eigen::MatrixXd *pressure)>;
124 void set_max_threads(
const int max_threads = std::numeric_limits<int>::max());
129 void init(
const json &
args,
const bool strict_validation);
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);
155 void init_logger(std::ostream &os,
const spdlog::level::level_enum log_level);
159 void set_log_level(
const spdlog::level::level_enum log_level);
164 std::string
get_log(
const Eigen::MatrixXd &sol)
166 std::stringstream ss;
173 void init_logger(
const std::vector<spdlog::sink_ptr> &sinks,
const spdlog::level::level_enum log_level);
189 std::shared_ptr<assembler::Assembler>
assembler =
nullptr;
205 std::vector<basis::ElementBases>
bases;
219 std::map<int, Eigen::MatrixXd>
polys;
221 std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>>
polys_3d;
279 const std::vector<basis::ElementBases> &
bases,
289 const std::vector<basis::ElementBases> &bases_)
const;
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);
306 return {{n_b_samples, n_b_samples}};
311 const int actual_dim =
problem->is_scalar() ? 1 :
mesh->dimension();
329 void set_materials(std::vector<std::shared_ptr<assembler::Assembler>> &assemblers)
const;
345 Eigen::MatrixXd &pressure,
347 const InitialConditionOverride *ic_override =
nullptr);
355 Eigen::MatrixXd &pressure,
361 logger().error(
"Load the mesh first!");
389 Eigen::MatrixXd &pressure,
400 Eigen::MatrixXd &sol,
401 Eigen::MatrixXd &pressure,
414 Eigen::MatrixXd &sol,
415 Eigen::MatrixXd &pressure,
429 Eigen::MatrixXd &sol,
430 Eigen::MatrixXd &pressure,
432 const InitialConditionOverride *ic_override =
nullptr);
444 Eigen::MatrixXd &sol,
446 const InitialConditionOverride *ic_override =
nullptr);
453 const double t = 1.0,
454 const bool init_time_integrator =
true,
455 const InitialConditionOverride *ic_override =
nullptr);
461 const double t = 1.0,
462 const InitialConditionOverride *ic_override =
nullptr);
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;
496 std::shared_ptr<polysolve::nonlinear::Solver>
make_nl_solver(
bool for_al)
const;
502 return args[
"boundary_conditions"][
"periodic_boundary"][
"enabled"].get<
bool>();
516 const std::unique_ptr<polysolve::linear::Solver> &solver,
519 const bool compute_spectrum,
544 std::unordered_map<int, std::array<bool, 3>>
587 std::unique_ptr<mesh::Mesh>
mesh;
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>());
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);
612 void load_mesh(
const Eigen::MatrixXd &
V,
const Eigen::MatrixXi &
F,
bool non_conforming =
false)
629 bool remesh(
const double time,
const double dt, Eigen::MatrixXd &sol);
634 vertices.setZero(
mesh->n_vertices(),
mesh->dimension());
636 for (
int v = 0; v <
mesh->n_vertices(); v++)
637 vertices.row(v) =
mesh->point(v);
643 assert(
mesh->is_simplicial());
648 int dim =
mesh->dimension();
649 elements.setZero(gbases.size(), dim + 1);
650 for (
int e = 0; e < gbases.size(); e++)
653 for (
const auto &gbs : gbases[e].bases)
654 elements(e, i++) = node_to_primitive_map[gbs.global()[0].index];
674 const std::vector<basis::ElementBases> &
bases,
675 const std::vector<basis::ElementBases> &
geom_bases,
701 return args[
"contact"][
"enabled"];
709 return args[
"contact"][
"adhesion"][
"adhesion_enabled"];
717 return (
args[
"boundary_conditions"][
"pressure_boundary"].size() > 0)
718 || (
args[
"boundary_conditions"][
"pressure_cavity"].size() > 0);
745 void export_data(
const Eigen::MatrixXd &sol,
const Eigen::MatrixXd &pressure);
754 void save_timestep(
const double time,
const int t,
const double t0,
const double dt,
const Eigen::MatrixXd &sol,
const Eigen::MatrixXd &pressure);
761 void save_subsolve(
const int i,
const int t,
const Eigen::MatrixXd &sol,
const Eigen::MatrixXd &pressure);
766 void save_json(
const Eigen::MatrixXd &sol, std::ostream &out);
770 void save_json(
const Eigen::MatrixXd &sol);
788 std::string
resolve_input_path(
const std::string &path,
const bool only_if_exists =
false)
const;
813 return args[
"boundary_conditions"][
"periodic_boundary"][
"linear_displacement_offset"].size() > 0;
Runtime override for initial-condition histories.
Eigen::MatrixXd solution
ndof by H (history size) matrix where each col is solution from a time step.
bool is_empty() const
Returns true when no quantity is overridden.
Eigen::MatrixXd acceleration
ndof by H (history size) matrix where each col is acceleration from a time step.
Eigen::MatrixXd velocity
ndof by H (history size) matrix where each col is velocity from a time step.
main class that contains the polyfem solver and all its state
void get_vertices(Eigen::MatrixXd &vertices) const
Gather geometry vertices into a dense matrix.
io::OutRuntimeData timings
runtime statistics
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
std::shared_ptr< utils::PeriodicBoundary > periodic_bc
periodic BC and periodic mesh utils
int n_bases
number of bases
int n_pressure_bases
number of pressure bases
std::string formulation() const
return the formulation (checks if the problem is scalar or not and deals with multiphysics)
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
bool is_adhesion_enabled() const
does the simulation have adhesion
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
double starting_max_edge_length
std::vector< mesh::LocalBoundary > local_pressure_boundary
mapping from elements to nodes for pressure boundary conditions
void sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
splits the solution in solution and pressure for mixed problems
std::shared_ptr< assembler::ViscousDamping > damping_assembler
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
void assemble_mass_mat()
assemble mass, step 4 of solve build mass matrix based on defined basis modifies mass (and maybe more...
Eigen::VectorXi in_node_to_node
Inpute nodes (including high-order) to polyfem nodes, only for isoparametric.
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
mesh::Obstacle obstacle
Obstacles used in collisions.
std::shared_ptr< assembler::Assembler > assembler
assemblers
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.
void build_periodic_collision_mesh()
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
bool use_avg_pressure
use average pressure for stokes problem to fix the additional dofs, true by default if false,...
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
ipc::CollisionMesh collision_mesh
IPC collision mesh.
std::string output_dir
Directory for output files.
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
std::vector< basis::ElementBases > geom_bases_
Geometric mapping bases, if the elements are isoparametric, this list is empty.
Eigen::VectorXi in_primitive_to_primitive
maps in vertices/edges/faces/cells to polyfem vertices/edges/faces/cells
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,...
std::map< int, basis::InterfaceData > poly_edge_to_data
nodes on the boundary of polygonal elements, used for harmonic bases
std::vector< int > pressure_boundary_nodes
list of neumann boundary nodes
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.
void save_json(const Eigen::MatrixXd &sol, std::ostream &out)
saves the output statistic to a stream
bool is_homogenization() const
void load_mesh(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, bool non_conforming=false)
loads the mesh from V and F,
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
std::shared_ptr< assembler::Assembler > pressure_assembler
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
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
void build_node_mapping()
build the mapping from input nodes to polyfem nodes
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
bool has_dhat
stores if input json contains dhat
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
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
void reset_mesh()
Resets the mesh.
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.
StiffnessMatrix mass
Mass matrix, it is computed only for time dependent problems.
void solve(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves the problem, call other methods
std::vector< int > dirichlet_nodes
per node dirichlet
bool has_periodic_bc() const
void build_basis()
builds the bases step 2 of solve modifies bases, pressure_bases, geom_bases_, boundary_nodes,...
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
std::vector< RowVectorNd > dirichlet_nodes_position
spdlog::sink_ptr file_sink_
spdlog::sink_ptr console_sink_
logger sink to stdout
std::map< int, std::pair< Eigen::MatrixXd, Eigen::MatrixXi > > polys_3d
polyhedra, used since poly have no geom mapping
std::vector< int > node_to_primitive() const
json args
main input arguments containing all defaults
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.
std::vector< RowVectorNd > neumann_nodes_position
void assemble_rhs()
compute rhs, step 3 of solve build rhs vector based on defined basis and given rhs of the problem mod...
double min_boundary_edge_length
std::unique_ptr< polysolve::linear::Solver > static_linear_solver_cache
Linear solver instance from the most recent static linear solve.
Eigen::VectorXi periodic_collision_mesh_to_basis
index mapping from periodic 2x2 collision mesh to FE periodic mesh
void init_time()
initialize time settings if args contains "time"
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
void initial_solution(Eigen::MatrixXd &solution, const InitialConditionOverride *ic_override=nullptr) const
Load or compute the initial solution.
Eigen::VectorXi disc_ordersq
void get_elements(Eigen::MatrixXi &elements) const
Gather geometry elements into a dense matrix.
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
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves the problems
assembler::AssemblyValsCache pressure_ass_vals_cache
used to store assembly values for pressure for small problems
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
double avg_mass
average system mass, used for contact with IPC
bool optimization_enabled
std::function< void(int, int, double, double)> time_callback
bool iso_parametric() const
check if using iso parametric bases
std::shared_ptr< assembler::RhsAssembler > build_rhs_assembler() const
build a RhsAssembler for the problem
std::vector< mesh::LocalBoundary > total_local_boundary
mapping from elements to nodes for all mesh
QuadratureOrders n_boundary_samples() const
quadrature used for projecting boundary conditions
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
bool is_obstacle_vertex(const size_t vi) const
checks if vertex is obstacle
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
solver::SolveData solve_data
timedependent stuff cached
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
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
std::shared_ptr< assembler::ViscousDampingPrev > damping_prev_assembler
bool is_contact_enabled() const
does the simulation have contact
void build_collision_mesh()
extracts the boundary mesh for collision, called in build_basis
Eigen::VectorXi disc_orders
vector of discretization orders, used when not all elements have the same degree, one per element
std::vector< mesh::LocalBoundary > local_neumann_boundary
mapping from elements to nodes for neumann boundary conditions
std::shared_ptr< assembler::PressureAssembler > build_pressure_assembler() const
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
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.
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
std::map< int, Eigen::MatrixXd > polys
polygons, used since poly have no geom mapping
Eigen::MatrixXd rhs
System right-hand side.
std::unordered_map< int, std::vector< mesh::LocalBoundary > > local_pressure_cavity
mapping from elements to nodes for pressure boundary conditions
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.
void compute_mesh_stats(const polyfem::mesh::Mesh &mesh)
compute stats (counts els type, mesh lenght, etc), step 1 of solve
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
static std::unique_ptr< Mesh > create(const std::string &path, const bool non_conforming=false)
factory to build the proper mesh
class to store time stepping data
spdlog::logger & logger()
Retrieves the current logger.
std::array< int, 2 > QuadratureOrders
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.
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix