PolyFEM
Loading...
Searching...
No Matches
StateSolveLinear.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2
6
9
13
16
18
19#include <Eigen/Core>
20#include <unsupported/Eigen/SparseExtra>
21#include <spdlog/fmt/fmt.h>
22#include <polysolve/linear/FEMSolver.hpp>
23
24#include <cassert>
25#include <memory>
26#include <string>
27#include <vector>
28
29namespace polyfem
30{
31 using namespace mesh;
32 using namespace time_integrator;
33 using namespace utils;
34 using namespace solver;
35 using namespace io;
36
38 {
39 igl::Timer timer;
40 timer.start();
41 logger().info("Assembling stiffness mat...");
42 assert(assembler->is_linear());
43
44 if (mixed_assembler != nullptr)
45 {
46
47 StiffnessMatrix velocity_stiffness, mixed_stiffness, pressure_stiffness;
48 assembler->assemble(mesh->is_volume(), n_bases, bases, geom_bases(), ass_vals_cache, 0, velocity_stiffness);
50 pressure_assembler->assemble(mesh->is_volume(), n_pressure_bases, pressure_bases, geom_bases(), pressure_ass_vals_cache, 0, pressure_stiffness);
51
52 const int problem_dim = problem->is_scalar() ? 1 : mesh->dimension();
53
55 velocity_stiffness, mixed_stiffness, pressure_stiffness,
56 stiffness);
57 }
58 else
59 {
60 assembler->assemble(mesh->is_volume(), n_bases, bases, geom_bases(), ass_vals_cache, 0, stiffness);
61 }
62
63 timer.stop();
64 timings.assembling_stiffness_mat_time = timer.getElapsedTime();
65 logger().info(" took {}s", timings.assembling_stiffness_mat_time);
66
67 stats.nn_zero = stiffness.nonZeros();
68 stats.num_dofs = stiffness.rows();
69 stats.mat_size = (long long)stiffness.rows() * (long long)stiffness.cols();
70 logger().info("sparsity: {}/{}", stats.nn_zero, stats.mat_size);
71
72 const std::string full_mat_path = args["output"]["data"]["full_mat"];
73 if (!full_mat_path.empty())
74 {
75 Eigen::saveMarket(stiffness, full_mat_path);
76 }
77 }
78
80 int step,
81 const std::unique_ptr<polysolve::linear::Solver> &solver,
83 Eigen::VectorXd &b,
84 const bool compute_spectrum,
85 Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step)
86 {
87 assert(assembler->is_linear() && !is_contact_enabled());
88 assert(solve_data.rhs_assembler != nullptr);
89
90 const int problem_dim = problem->is_scalar() ? 1 : mesh->dimension();
91 const int full_size = A.rows();
92 int precond_num = problem_dim * n_bases;
93
94 std::vector<int> boundary_nodes_tmp;
95 if (has_periodic_bc())
96 {
97 boundary_nodes_tmp = periodic_bc->full_to_periodic(boundary_nodes);
98 precond_num = periodic_bc->full_to_periodic(A);
99 Eigen::MatrixXd tmp = periodic_bc->full_to_periodic(b, true);
100 b = tmp;
101 }
102 else
103 boundary_nodes_tmp = boundary_nodes;
104
105 Eigen::VectorXd x;
106 stats.spectrum = dirichlet_solve(
107 *solver,
108 A,
109 b,
110 boundary_nodes_tmp,
111 x,
112 precond_num,
113 args["output"]["data"]["stiffness_mat"],
114 compute_spectrum,
115 assembler->is_fluid(),
117
118 if (has_periodic_bc())
119 {
120 sol = periodic_bc->periodic_to_full(full_size, x);
121 if (args["/boundary_conditions/periodic_boundary/force_zero_mean"_json_pointer].get<bool>())
122 {
123 Eigen::VectorXd integral = Evaluator::integrate_function(bases, geom_bases(), sol, mesh->dimension(), problem_dim);
124 double area = Evaluator::integrate_function(bases, geom_bases(), Eigen::VectorXd::Ones(n_bases), mesh->dimension(), 1)(0);
125 for (int d = 0; d < problem_dim; d++)
126 sol(Eigen::seqN(d, n_bases, problem_dim), 0).array() -= integral(d) / area;
127 }
128 }
129 else
130 sol = x; // Explicit copy because sol is a MatrixXd (with one column)
131
132 solver->get_info(stats.solver_info);
133
134 const auto error = (A * x - b).norm();
135
136 if (error > 1e-4)
137 logger().error("Solver error: {}", error);
138 else
139 logger().debug("Solver error: {}", error);
140
141 if (mixed_assembler != nullptr)
142 sol_to_pressure(sol, pressure);
143
144 if (user_post_step)
145 {
146 user_post_step(step, *this, sol, nullptr, nullptr);
147 }
148 }
149
150 void State::solve_linear(int step, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step)
151 {
152 assert(!problem->is_time_dependent());
153 assert(assembler->is_linear() && !is_contact_enabled());
154
155 // --------------------------------------------------------------------
156
158 polysolve::linear::Solver::create(args["solver"]["linear"], logger());
159 logger().info("{}...", static_linear_solver_cache->name());
160
161 // --------------------------------------------------------------------
162
165 (assembler->name() != "Bilaplacian") ? local_neumann_boundary : std::vector<LocalBoundary>(), rhs);
166
169
170 Eigen::VectorXd b = rhs;
171
172 // --------------------------------------------------------------------
173
174 solve_linear(step, static_linear_solver_cache, A, b, args["output"]["advanced"]["spectrum"], sol, pressure, user_post_step);
175 }
176
177 void State::init_linear_solve(Eigen::MatrixXd &sol, const double t, const InitialConditionOverride *ic_override)
178 {
179 assert(sol.cols() == 1);
180 assert(assembler->is_linear() && !is_contact_enabled()); // linear
181
182 if (mixed_assembler != nullptr)
183 return;
184
185 const int ndof = n_bases * mesh->dimension();
186
187 solve_data.elastic_form = std::make_shared<ElasticForm>(
190 t, problem->is_time_dependent() ? args["time"]["dt"].get<double>() : 0.0,
191 mesh->is_volume());
192
193 solve_data.body_form = std::make_shared<BodyForm>(
197 mass_matrix_assembler->density(),
198 /*is_formulation_mixed=*/false, problem->is_time_dependent());
199 solve_data.body_form->update_quantities(t, sol);
200
201 solve_data.inertia_form = nullptr;
202 solve_data.damping_form = nullptr;
203 if (problem->is_time_dependent())
204 {
206 solve_data.inertia_form = std::make_shared<InertiaForm>(mass, *solve_data.time_integrator);
207 }
208
209 solve_data.contact_form = nullptr;
210 solve_data.friction_form = nullptr;
211
213 // Initialize time integrator
214 if (problem->is_time_dependent() && assembler->is_tensor())
215 {
216 POLYFEM_SCOPED_TIMER("Initialize time integrator");
217
218 Eigen::MatrixXd solution, velocity, acceleration;
219 initial_solution(solution, ic_override); // Reload this because we need all previous solutions
220 solution.col(0) = sol; // Make sure the current solution is the same as `sol`
221 assert(solution.rows() == sol.size());
222 initial_velocity(velocity, ic_override);
223 assert(velocity.rows() == sol.size());
224 initial_acceleration(acceleration, ic_override);
225 assert(acceleration.rows() == sol.size());
226
227 if (solution.cols() != velocity.cols() || solution.cols() != acceleration.cols())
228 {
230 "Incompatible initial-condition history for transient solve: "
231 "solution has {} columns, velocity has {}, acceleration has {}.",
232 solution.cols(), velocity.cols(), acceleration.cols());
233 }
234
235 const double dt = args["time"]["dt"];
236 solve_data.time_integrator->init(solution, velocity, acceleration, dt);
237 }
239 }
240
241 void State::solve_transient_linear(const int time_steps,
242 const double t0,
243 const double dt,
244 Eigen::MatrixXd &sol,
245 Eigen::MatrixXd &pressure,
246 UserPostStepCallback user_post_step,
247 const InitialConditionOverride *ic_override)
248 {
249 assert(sol.cols() == 1);
250 assert(problem->is_time_dependent());
251 assert(assembler->is_linear() && !is_contact_enabled());
252 assert(solve_data.rhs_assembler != nullptr);
253
254 const bool is_scalar_or_mixed = problem->is_scalar() || mixed_assembler != nullptr;
255
256 // --------------------------------------------------------------------
257
258 auto solver =
259 polysolve::linear::Solver::create(args["solver"]["linear"], logger());
260 logger().info("{}...", solver->name());
261
262 // --------------------------------------------------------------------
263
264 std::shared_ptr<ImplicitTimeIntegrator> time_integrator;
265 if (is_scalar_or_mixed)
266 {
267 time_integrator = std::make_shared<BDF>();
268 time_integrator->set_parameters(args["time"]);
269 time_integrator->init(sol, Eigen::VectorXd::Zero(sol.size()), Eigen::VectorXd::Zero(sol.size()), dt);
270 }
271 else
272 {
273 Eigen::MatrixXd solution, velocity, acceleration;
274 initial_solution(solution, ic_override); // Reload this because we need all previous solutions
275 solution.col(0) = sol; // Make sure the current solution is the same as `sol`
276 assert(solution.rows() == sol.size());
277 initial_velocity(velocity, ic_override);
278 assert(velocity.rows() == sol.size());
279 initial_acceleration(acceleration, ic_override);
280 assert(acceleration.rows() == sol.size());
281
282 if (solution.cols() != velocity.cols() || solution.cols() != acceleration.cols())
283 {
285 "Incompatible initial-condition history for transient solve: "
286 "solution has {} columns, velocity has {}, acceleration has {}.",
287 solution.cols(), velocity.cols(), acceleration.cols());
288 }
289
290 time_integrator = ImplicitTimeIntegrator::construct_time_integrator(args["time"]["integrator"]);
291 time_integrator->init(solution, velocity, acceleration, dt);
292 }
293
294 // --------------------------------------------------------------------
295
296 const QuadratureOrders &n_b_samples = n_boundary_samples();
297
298 // Step 0.
299 if (user_post_step)
300 {
301 user_post_step(0, *this, sol, nullptr, nullptr);
302 }
303
304 Eigen::MatrixXd current_rhs = rhs;
305
306 StiffnessMatrix stiffness;
307 build_stiffness_mat(stiffness);
308
309 // --------------------------------------------------------------------
310 // TODO rebuild stiffnes if material are time dept
311 for (int t = 1; t <= time_steps; ++t)
312 {
313 const double time = t0 + t * dt;
314
316 Eigen::VectorXd b;
317 bool compute_spectrum = args["output"]["advanced"]["spectrum"];
318
319 if (is_scalar_or_mixed)
320 {
321 solve_data.rhs_assembler->compute_energy_grad(
323 current_rhs);
324
326 local_boundary, boundary_nodes, n_b_samples, local_neumann_boundary, current_rhs, sol, time);
327
328 if (mixed_assembler != nullptr)
329 {
330 // divergence free
331 int fluid_offset = use_avg_pressure ? (assembler->is_fluid() ? 1 : 0) : 0;
332 current_rhs
333 .block(
334 current_rhs.rows() - n_pressure_bases - use_avg_pressure, 0,
335 n_pressure_bases + use_avg_pressure, current_rhs.cols())
336 .setZero();
337 }
338
339 std::shared_ptr<BDF> bdf = std::dynamic_pointer_cast<BDF>(time_integrator);
340 A = mass / bdf->beta_dt() + stiffness;
341 b = (mass * bdf->weighted_sum_x_prevs()) / bdf->beta_dt();
342 for (int i : boundary_nodes)
343 b[i] = 0;
344 b += current_rhs;
345
346 compute_spectrum &= t == time_steps;
347 }
348 else
349 {
350 solve_data.rhs_assembler->assemble(mass_matrix_assembler->density(), current_rhs, time);
351
352 current_rhs *= -1;
353
355 std::vector<LocalBoundary>(), std::vector<int>(), n_b_samples, local_neumann_boundary, current_rhs, sol, time);
356
357 current_rhs *= time_integrator->acceleration_scaling();
358 current_rhs += mass * time_integrator->x_tilde();
359
361 local_boundary, boundary_nodes, n_b_samples, std::vector<LocalBoundary>(), current_rhs, sol, time);
362
363 A = stiffness * time_integrator->acceleration_scaling() + mass;
364 b = current_rhs;
365
366 compute_spectrum &= t == 1;
367 }
368
369 solve_linear(t, solver, A, b, compute_spectrum, sol, pressure, user_post_step);
370
371 time_integrator->update_quantities(sol);
372
373 save_timestep(time, t, t0, dt, sol, pressure);
374
375 const std::string &state_path = resolve_output_path(fmt::format(args["output"]["data"]["state"], t));
376 if (!state_path.empty())
377 time_integrator->save_state(state_path);
378
379 logger().info("{}/{} t={}", t, time_steps, time);
380 }
381 }
382} // namespace polyfem
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
Runtime override for initial-condition histories.
Definition State.hpp:90
io::OutRuntimeData timings
runtime statistics
Definition State.hpp:733
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
assembler::AssemblyValsCache ass_vals_cache
used to store assembly values for small problems
Definition State.hpp:230
const std::vector< basis::ElementBases > & geom_bases() const
Get a constant reference to the geometry mapping bases.
Definition State.hpp:257
void sol_to_pressure(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
splits the solution in solution and pressure for mixed problems
Definition State.cpp:375
void solve_linear(int step, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={})
solves a linear problem
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:189
void initial_velocity(Eigen::MatrixXd &velocity, const InitialConditionOverride *ic_override=nullptr) const
Load or compute the initial velocity.
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.
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
std::shared_ptr< assembler::Assembler > pressure_assembler
Definition State.hpp:194
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:191
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:587
StiffnessMatrix mass
Mass matrix, it is computed only for time dependent problems.
Definition State.hpp:236
bool has_periodic_bc() const
Definition State.hpp:500
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:202
json args
main input arguments containing all defaults
Definition State.hpp:135
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::unique_ptr< polysolve::linear::Solver > static_linear_solver_cache
Linear solver instance from the most recent static linear solve.
Definition State.hpp:382
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.
assembler::AssemblyValsCache pressure_ass_vals_cache
used to store assembly values for pressure for small problems
Definition State.hpp:233
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
QuadratureOrders n_boundary_samples() const
quadrature used for projecting boundary conditions
Definition State.hpp:297
int ndof() const
Definition State.hpp:309
std::vector< int > boundary_nodes
list of boundary nodes
Definition State.hpp:548
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:375
bool is_contact_enabled() const
does the simulation have contact
Definition State.hpp:699
std::vector< mesh::LocalBoundary > local_neumann_boundary
mapping from elements to nodes for neumann boundary conditions
Definition State.hpp:556
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::MixedAssembler > mixed_assembler
Definition State.hpp:193
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:241
static void merge_mixed_matrices(const int n_bases, const int n_pressure_bases, const int problem_dim, const bool add_average, const StiffnessMatrix &velocity_stiffness, const StiffnessMatrix &mixed_stiffness, const StiffnessMatrix &pressure_stiffness, StiffnessMatrix &stiffness)
utility to merge 3 blocks of mixed matrices, A=velocity_stiffness, B=mixed_stiffness,...
static Eigen::VectorXd integrate_function(const std::vector< basis::ElementBases > &bases, const std::vector< basis::ElementBases > &gbases, const Eigen::MatrixXd &fun, const int dim, const int actual_dim)
double assembling_stiffness_mat_time
time to assembly
Definition OutData.hpp:357
json solver_info
information of the solver, eg num iteration, time, errors, etc the informations varies depending on t...
Definition OutData.hpp:382
Eigen::Vector4d spectrum
spectrum of the stiffness matrix, enable only if POLYSOLVE_WITH_SPECTRA is ON (off by default)
Definition OutData.hpp:378
long long nn_zero
non zeros and sytem matrix size num dof is the total dof in the system
Definition OutData.hpp:396
std::shared_ptr< solver::FrictionForm > friction_form
std::shared_ptr< solver::InertiaForm > inertia_form
void update_dt()
updates the dt inside the different forms
std::shared_ptr< solver::BodyForm > body_form
std::shared_ptr< solver::ContactForm > contact_form
std::shared_ptr< solver::ElasticForm > damping_form
std::shared_ptr< solver::ElasticForm > elastic_form
std::shared_ptr< time_integrator::ImplicitTimeIntegrator > time_integrator
std::shared_ptr< assembler::RhsAssembler > rhs_assembler
static std::shared_ptr< ImplicitTimeIntegrator > construct_time_integrator(const json &params)
Factory method for constructing implicit time integrators from the name of the integrator.
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
std::array< int, 2 > QuadratureOrders
Definition Types.hpp:19
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
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24