35#include <polysolve/linear/Solver.hpp>
38#include <Eigen/Sparse>
42#include <unordered_map>
44#include <spdlog/sinks/basic_file_sink.h>
46#include <ipc/collision_mesh.hpp>
47#include <ipc/utils/logger.hpp>
59 class ViscousDampingPrev;
90 void set_max_threads(
const int max_threads = std::numeric_limits<int>::max());
95 void init(
const json &
args,
const bool strict_validation);
113 const std::string &log_file,
114 const spdlog::level::level_enum log_level,
115 const spdlog::level::level_enum file_log_level,
116 const bool is_quiet);
121 void init_logger(std::ostream &os,
const spdlog::level::level_enum log_level);
125 void set_log_level(
const spdlog::level::level_enum log_level);
130 std::string
get_log(
const Eigen::MatrixXd &sol)
132 std::stringstream ss;
139 void init_logger(
const std::vector<spdlog::sink_ptr> &sinks,
const spdlog::level::level_enum log_level);
155 std::shared_ptr<assembler::Assembler>
assembler =
nullptr;
171 std::vector<basis::ElementBases>
bases;
185 std::map<int, Eigen::MatrixXd>
polys;
187 std::map<int, std::pair<Eigen::MatrixXd, Eigen::MatrixXi>>
polys_3d;
245 const std::vector<basis::ElementBases> &
bases,
255 const std::vector<basis::ElementBases> &bases_)
const;
266 const int n_b_samples_j =
args[
"space"][
"advanced"][
"n_boundary_samples"];
267 const int gdiscr_order =
mesh->orders().size() <= 0 ? 1 :
mesh->orders().maxCoeff();
268 const int discr_order = std::max(
disc_orders.maxCoeff(), gdiscr_order);
286 void set_materials(std::vector<std::shared_ptr<assembler::Assembler>> &assemblers)
const;
299 void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure);
303 void solve(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
307 logger().error(
"Load the mesh first!");
328 void init_solve(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure);
348 void solve_transient_linear(
const int time_steps,
const double t0,
const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure);
374 void solve_linear(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure);
386 std::shared_ptr<polysolve::nonlinear::Solver>
make_nl_solver(
bool for_al)
const;
392 return args[
"boundary_conditions"][
"periodic_boundary"][
"enabled"].get<
bool>();
403 const std::unique_ptr<polysolve::linear::Solver> &solver,
406 const bool compute_spectrum,
407 Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure);
423 std::unordered_map<int, std::array<bool, 3>>
466 std::unique_ptr<mesh::Mesh>
mesh;
475 void load_mesh(
bool non_conforming =
false,
476 const std::vector<std::string> &names = std::vector<std::string>(),
477 const std::vector<Eigen::MatrixXi> &cells = std::vector<Eigen::MatrixXi>(),
478 const std::vector<Eigen::MatrixXd> &vertices = std::vector<Eigen::MatrixXd>());
485 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);
491 void load_mesh(
const Eigen::MatrixXd &
V,
const Eigen::MatrixXi &
F,
bool non_conforming =
false)
508 bool remesh(
const double time,
const double dt, Eigen::MatrixXd &sol);
526 const std::vector<basis::ElementBases> &
bases,
527 const std::vector<basis::ElementBases> &
geom_bases,
553 return args[
"contact"][
"enabled"];
561 return (
args[
"boundary_conditions"][
"pressure_boundary"].size() > 0)
562 || (
args[
"boundary_conditions"][
"pressure_cavity"].size() > 0);
593 void export_data(
const Eigen::MatrixXd &sol,
const Eigen::MatrixXd &pressure);
602 void save_timestep(
const double time,
const int t,
const double t0,
const double dt,
const Eigen::MatrixXd &sol,
const Eigen::MatrixXd &pressure);
609 void save_subsolve(
const int i,
const int t,
const Eigen::MatrixXd &sol,
const Eigen::MatrixXd &pressure);
614 void save_json(
const Eigen::MatrixXd &sol, std::ostream &out);
618 void save_json(
const Eigen::MatrixXd &sol);
636 std::string
resolve_input_path(
const std::string &path,
const bool only_if_exists =
false)
const;
655 const int actual_dim =
problem->is_scalar() ? 1 :
mesh->dimension();
672 if (
problem->is_time_dependent())
710 void solve_homogenization(
const int time_steps,
const double t0,
const double dt, Eigen::MatrixXd &sol);
713 return args[
"boundary_conditions"][
"periodic_boundary"][
"linear_displacement_offset"].size() > 0;
main class that contains the polyfem solver and all its state
void get_vertices(Eigen::MatrixXd &vertices) const
Eigen::MatrixXd initial_vel_update
io::OutRuntimeData timings
runtime statistics
assembler::MacroStrainValue macro_strain_constraint
std::shared_ptr< utils::PeriodicBoundary > periodic_bc
periodic BC and periodic mesh utils
int n_bases
number of bases
void cache_transient_adjoint_quantities(const int current_step, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &disp_grad)
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
int n_boundary_samples() const
quadrature used for projecting boundary conditions
void set_mesh_vertex(int v_id, const Eigen::VectorXd &vertex)
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 build_polygonal_basis()
builds bases for polygons, called inside build_basis
void solve(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves the problem, call other methods
mesh::Obstacle obstacle
Obstacles used in collisions.
std::shared_ptr< assembler::Assembler > assembler
assemblers
void compute_errors(const Eigen::MatrixXd &sol)
computes all errors
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.
solver::CacheLevel optimization_enabled
void solve_homogenization_step(Eigen::MatrixXd &sol, const int t=0, bool adaptive_initial_weight=false)
In Elasticity PDE, solve for "min W(disp_grad + \grad u)" instead of "min W(\grad u)".
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
Eigen::MatrixXd get_adjoint_mat(int type) 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_transient_tensor_nonlinear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol)
solves transient tensor nonlinear problem
void initial_solution(Eigen::MatrixXd &solution) const
Load or compute the initial solution.
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
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
Eigen::MatrixXd solve_static_adjoint(const Eigen::MatrixXd &adjoint_rhs) const
void compute_force_jacobian(const Eigen::MatrixXd &sol, const Eigen::MatrixXd &disp_grad, StiffnessMatrix &hessian)
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_homogenization(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol)
std::vector< int > dirichlet_nodes
per node dirichlet
bool has_periodic_bc() const
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves the problems
void build_basis()
builds the bases step 2 of solve modifies bases, pressure_bases, geom_bases_, boundary_nodes,...
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
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.
std::vector< io::SolutionFrame > solution_frames
saves the frames in a vector instead of VTU
solver::DiffCache diff_cached
void initial_velocity(Eigen::MatrixXd &velocity) const
Load or compute the initial velocity.
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
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"
void initial_acceleration(Eigen::MatrixXd &acceleration) const
Load or compute the initial acceleration.
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
void solve_adjoint_cached(const Eigen::MatrixXd &rhs)
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
Eigen::MatrixXd solve_adjoint(const Eigen::MatrixXd &rhs) const
void get_elements(Eigen::MatrixXi &elements) const
int n_geom_bases
number of geometric bases
assembler::AssemblyValsCache pressure_ass_vals_cache
used to store assembly values for pressure for small problems
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 iso_parametric() const
check if using iso parametric bases
std::shared_ptr< assembler::RhsAssembler > build_rhs_assembler() const
build a RhsAssembler for the problem
void compute_volume_node_ids(const int volume_selection, std::vector< int > &node_ids) const
std::vector< mesh::LocalBoundary > total_local_boundary
mapping from elements to nodes for all mesh
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
void init_solve(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
initialize solver
void init_linear_solve(Eigen::MatrixXd &sol, const double t=1.0)
initialize the linear solve
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
void solve_tensor_nonlinear(Eigen::MatrixXd &sol, const int t=0, const bool init_lagging=true)
solves nonlinear problems
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
std::unique_ptr< polysolve::linear::Solver > lin_solver_cached
std::shared_ptr< polyfem::mesh::MeshNodes > geom_mesh_nodes
void solve_linear(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves a linear problem
void init_nonlinear_tensor_solve(Eigen::MatrixXd &sol, const double t=1.0, const bool init_time_integrator=true)
initialize the nonlinear solver
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 has contact
Eigen::MatrixXd solve_transient_adjoint(const Eigen::MatrixXd &adjoint_rhs) const
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
void compute_force_jacobian_prev(const int force_step, const int sol_step, StiffnessMatrix &hessian_prev) const
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
bool solve_export_to_file
flag to decide if exporting the time dependent solution to files or save it in the solution_frames ar...
void solve_transient_linear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves transient linear problem
void compute_total_surface_node_ids(std::vector< int > &node_ids) const
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.
void compute_surface_node_ids(const int surface_selection, std::vector< int > &node_ids) const
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
Eigen::MatrixXd initial_sol_update
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)
StiffnessMatrix basis_nodes_to_gbasis_nodes
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
const Eigen::MatrixXd & adjoint_mat() const
class to store time stepping data
spdlog::logger & logger()
Retrieves the current logger.
void log_and_throw_adjoint_error(const std::string &msg)
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix