PolyFEM
Loading...
Searching...
No Matches
StateSolveLinear.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2
5
8
12#include <polysolve/linear/FEMSolver.hpp>
13
15
16#include <unsupported/Eigen/SparseExtra>
18
19namespace polyfem
20{
21 using namespace mesh;
22 using namespace time_integrator;
23 using namespace utils;
24 using namespace solver;
25 using namespace io;
26
28 {
29 igl::Timer timer;
30 timer.start();
31 logger().info("Assembling stiffness mat...");
32 assert(assembler->is_linear());
33
34 if (mixed_assembler != nullptr)
35 {
36
37 StiffnessMatrix velocity_stiffness, mixed_stiffness, pressure_stiffness;
38 assembler->assemble(mesh->is_volume(), n_bases, bases, geom_bases(), ass_vals_cache, 0, velocity_stiffness);
40 pressure_assembler->assemble(mesh->is_volume(), n_pressure_bases, pressure_bases, geom_bases(), pressure_ass_vals_cache, 0, pressure_stiffness);
41
42 const int problem_dim = problem->is_scalar() ? 1 : mesh->dimension();
43
45 velocity_stiffness, mixed_stiffness, pressure_stiffness,
46 stiffness);
47 }
48 else
49 {
50 assembler->assemble(mesh->is_volume(), n_bases, bases, geom_bases(), ass_vals_cache, 0, stiffness);
51 }
52
53 timer.stop();
54 timings.assembling_stiffness_mat_time = timer.getElapsedTime();
55 logger().info(" took {}s", timings.assembling_stiffness_mat_time);
56
57 stats.nn_zero = stiffness.nonZeros();
58 stats.num_dofs = stiffness.rows();
59 stats.mat_size = (long long)stiffness.rows() * (long long)stiffness.cols();
60 logger().info("sparsity: {}/{}", stats.nn_zero, stats.mat_size);
61
62 const std::string full_mat_path = args["output"]["data"]["full_mat"];
63 if (!full_mat_path.empty())
64 {
65 Eigen::saveMarket(stiffness, full_mat_path);
66 }
67 }
68
70 const std::unique_ptr<polysolve::linear::Solver> &solver,
72 Eigen::VectorXd &b,
73 const bool compute_spectrum,
74 Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
75 {
76 assert(assembler->is_linear() && !is_contact_enabled());
77 assert(solve_data.rhs_assembler != nullptr);
78
79 const int problem_dim = problem->is_scalar() ? 1 : mesh->dimension();
80 const int full_size = A.rows();
81 int precond_num = problem_dim * n_bases;
82
83 std::vector<int> boundary_nodes_tmp;
84 if (has_periodic_bc())
85 {
86 boundary_nodes_tmp = periodic_bc->full_to_periodic(boundary_nodes);
87 precond_num = periodic_bc->full_to_periodic(A);
88 Eigen::MatrixXd tmp = periodic_bc->full_to_periodic(b, true);
89 b = tmp;
90 }
91 else
92 boundary_nodes_tmp = boundary_nodes;
93
94 Eigen::VectorXd x;
96 {
97 auto A_tmp = A;
98 prefactorize(*solver, A, boundary_nodes_tmp, precond_num, args["output"]["data"]["stiffness_mat"]);
99 dirichlet_solve_prefactorized(*solver, A_tmp, b, boundary_nodes_tmp, x);
100 }
101 else
102 {
103 stats.spectrum = dirichlet_solve(
104 *solver, A, b, boundary_nodes_tmp, x, precond_num, args["output"]["data"]["stiffness_mat"], compute_spectrum,
105 assembler->is_fluid(), use_avg_pressure);
106 }
107 if (has_periodic_bc())
108 {
109 sol = periodic_bc->periodic_to_full(full_size, x);
110 if (args["/boundary_conditions/periodic_boundary/force_zero_mean"_json_pointer].get<bool>())
111 {
112 Eigen::VectorXd integral = Evaluator::integrate_function(bases, geom_bases(), sol, mesh->dimension(), problem_dim);
113 double area = Evaluator::integrate_function(bases, geom_bases(), Eigen::VectorXd::Ones(n_bases), mesh->dimension(), 1)(0);
114 for (int d = 0; d < problem_dim; d++)
115 sol(Eigen::seqN(d, n_bases, problem_dim), 0).array() -= integral(d) / area;
116 }
117 }
118 else
119 sol = x; // Explicit copy because sol is a MatrixXd (with one column)
120
121 solver->get_info(stats.solver_info);
122
123 const auto error = (A * x - b).norm();
124
125 if (error > 1e-4)
126 logger().error("Solver error: {}", error);
127 else
128 logger().debug("Solver error: {}", error);
129
130 if (mixed_assembler != nullptr)
131 sol_to_pressure(sol, pressure);
132 }
133
134 void State::solve_linear(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
135 {
136 assert(!problem->is_time_dependent());
137 assert(assembler->is_linear() && !is_contact_enabled());
138
139 // --------------------------------------------------------------------
141 lin_solver_cached.reset();
142
144 polysolve::linear::Solver::create(args["solver"]["linear"], logger());
145 logger().info("{}...", lin_solver_cached->name());
146
147 // --------------------------------------------------------------------
148
151 (assembler->name() != "Bilaplacian") ? local_neumann_boundary : std::vector<LocalBoundary>(), rhs);
152
155
156 Eigen::VectorXd b = rhs;
157
158 // --------------------------------------------------------------------
159
160 solve_linear(lin_solver_cached, A, b, args["output"]["advanced"]["spectrum"], sol, pressure);
161 }
162
163 void State::init_linear_solve(Eigen::MatrixXd &sol, const double t)
164 {
165 assert(sol.cols() == 1);
166 assert(assembler->is_linear() && !is_contact_enabled()); // linear
167
168 if (mixed_assembler != nullptr)
169 return;
170
171 const int ndof = n_bases * mesh->dimension();
172
173 solve_data.elastic_form = std::make_shared<ElasticForm>(
176 t, problem->is_time_dependent() ? args["time"]["dt"].get<double>() : 0.0,
177 mesh->is_volume());
178
179 solve_data.body_form = std::make_shared<BodyForm>(
183 mass_matrix_assembler->density(),
184 /*apply_DBC=*/true, /*is_formulation_mixed=*/false, problem->is_time_dependent());
185 solve_data.body_form->update_quantities(t, sol);
186
187 solve_data.inertia_form = nullptr;
188 solve_data.damping_form = nullptr;
189 if (problem->is_time_dependent())
190 {
192 solve_data.inertia_form = std::make_shared<InertiaForm>(mass, *solve_data.time_integrator);
193 }
194
195 solve_data.contact_form = nullptr;
196 solve_data.friction_form = nullptr;
197
199 // Initialize time integrator
200 if (problem->is_time_dependent() && assembler->is_tensor())
201 {
202 POLYFEM_SCOPED_TIMER("Initialize time integrator");
203
204 Eigen::MatrixXd solution, velocity, acceleration;
205 initial_solution(solution); // Reload this because we need all previous solutions
206 solution.col(0) = sol; // Make sure the current solution is the same as `sol`
207 assert(solution.rows() == sol.size());
208 initial_velocity(velocity);
209 assert(velocity.rows() == sol.size());
210 initial_acceleration(acceleration);
211 assert(acceleration.rows() == sol.size());
212
213 const double dt = args["time"]["dt"];
214 solve_data.time_integrator->init(solution, velocity, acceleration, dt);
215 }
217 }
218
219 void State::solve_transient_linear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
220 {
221 assert(sol.cols() == 1);
222 assert(problem->is_time_dependent());
223 assert(assembler->is_linear() && !is_contact_enabled());
224 assert(solve_data.rhs_assembler != nullptr);
225
226 const bool is_scalar_or_mixed = problem->is_scalar() || mixed_assembler != nullptr;
227
228 // --------------------------------------------------------------------
229
230 auto solver =
231 polysolve::linear::Solver::create(args["solver"]["linear"], logger());
232 logger().info("{}...", solver->name());
233
234 // --------------------------------------------------------------------
235
236 std::shared_ptr<ImplicitTimeIntegrator> time_integrator;
237 if (is_scalar_or_mixed)
238 {
239 time_integrator = std::make_shared<BDF>();
240 time_integrator->set_parameters(args["time"]);
241 time_integrator->init(sol, Eigen::VectorXd::Zero(sol.size()), Eigen::VectorXd::Zero(sol.size()), dt);
242 }
243 else
244 {
245 Eigen::MatrixXd solution, velocity, acceleration;
246 initial_solution(solution); // Reload this because we need all previous solutions
247 solution.col(0) = sol; // Make sure the current solution is the same as `sol`
248 assert(solution.rows() == sol.size());
249 initial_velocity(velocity);
250 assert(velocity.rows() == sol.size());
251 initial_acceleration(acceleration);
252 assert(acceleration.rows() == sol.size());
253
254 time_integrator = ImplicitTimeIntegrator::construct_time_integrator(args["time"]["integrator"]);
255 time_integrator->init(solution, velocity, acceleration, dt);
256 }
257
258 // --------------------------------------------------------------------
259
260 const int n_b_samples = n_boundary_samples();
261
263 {
264 log_and_throw_adjoint_error("Transient linear problems are not differentiable yet!");
265 cache_transient_adjoint_quantities(0, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
266 }
267
268 Eigen::MatrixXd current_rhs = rhs;
269
270 StiffnessMatrix stiffness;
271 build_stiffness_mat(stiffness);
272
273 // --------------------------------------------------------------------
274 // TODO rebuild stiffnes if material are time dept
275 for (int t = 1; t <= time_steps; ++t)
276 {
277 const double time = t0 + t * dt;
278
280 Eigen::VectorXd b;
281 bool compute_spectrum = args["output"]["advanced"]["spectrum"];
282
283 if (is_scalar_or_mixed)
284 {
285 solve_data.rhs_assembler->compute_energy_grad(
287 current_rhs);
288
290 local_boundary, boundary_nodes, n_b_samples, local_neumann_boundary, current_rhs, sol, time);
291
292 if (mixed_assembler != nullptr)
293 {
294 // divergence free
295 int fluid_offset = use_avg_pressure ? (assembler->is_fluid() ? 1 : 0) : 0;
296 current_rhs
297 .block(
298 current_rhs.rows() - n_pressure_bases - use_avg_pressure, 0,
299 n_pressure_bases + use_avg_pressure, current_rhs.cols())
300 .setZero();
301 }
302
303 std::shared_ptr<BDF> bdf = std::dynamic_pointer_cast<BDF>(time_integrator);
304 A = mass / bdf->beta_dt() + stiffness;
305 b = (mass * bdf->weighted_sum_x_prevs()) / bdf->beta_dt();
306 for (int i : boundary_nodes)
307 b[i] = 0;
308 b += current_rhs;
309
310 compute_spectrum &= t == time_steps;
311 }
312 else
313 {
314 solve_data.rhs_assembler->assemble(mass_matrix_assembler->density(), current_rhs, time);
315
316 current_rhs *= -1;
317
319 std::vector<LocalBoundary>(), std::vector<int>(), n_b_samples, local_neumann_boundary, current_rhs, sol, time);
320
321 current_rhs *= time_integrator->acceleration_scaling();
322 current_rhs += mass * time_integrator->x_tilde();
323
325 local_boundary, boundary_nodes, n_b_samples, std::vector<LocalBoundary>(), current_rhs, sol, time);
326
327 A = stiffness * time_integrator->acceleration_scaling() + mass;
328 b = current_rhs;
329
330 compute_spectrum &= t == 1;
331 }
332
333 solve_linear(solver, A, b, compute_spectrum, sol, pressure);
334
336 {
337 log_and_throw_adjoint_error("Transient linear problems are not differentiable yet!");
338 cache_transient_adjoint_quantities(t, sol, Eigen::MatrixXd::Zero(mesh->dimension(), mesh->dimension()));
339 }
340
341 time_integrator->update_quantities(sol);
342
343 save_timestep(time, t, t0, dt, sol, pressure);
344 logger().info("{}/{} t={}", t, time_steps, time);
345 }
346
347 time_integrator->save_state(resolve_output_path(args["output"]["data"]["state"]));
348 }
349} // namespace polyfem
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
io::OutRuntimeData timings
runtime statistics
Definition State.hpp:603
std::shared_ptr< utils::PeriodicBoundary > periodic_bc
periodic BC and periodic mesh utils
Definition State.hpp:389
int n_bases
number of bases
Definition State.hpp:178
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
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
std::string resolve_output_path(const std::string &path) const
Resolve output path relative to output_dir if the path is not absolute.
solver::CacheLevel optimization_enabled
Definition State.hpp:667
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
void initial_solution(Eigen::MatrixXd &solution) const
Load or compute the initial solution.
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
bool has_periodic_bc() const
Definition State.hpp:390
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
void initial_velocity(Eigen::MatrixXd &velocity) const
Load or compute the initial velocity.
io::OutStatsData stats
Other statistics.
Definition State.hpp:605
void initial_acceleration(Eigen::MatrixXd &acceleration) const
Load or compute the initial acceleration.
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
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 init_linear_solve(Eigen::MatrixXd &sol, const double t=1.0)
initialize the linear solve
int ndof() const
Definition State.hpp:673
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::unique_ptr< polysolve::linear::Solver > lin_solver_cached
Definition State.hpp:671
void solve_linear(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves a linear problem
bool is_contact_enabled() const
does the simulation has contact
Definition State.hpp:574
std::vector< mesh::LocalBoundary > local_neumann_boundary
mapping from elements to nodes for neumann boundary conditions
Definition State.hpp:435
void solve_transient_linear(const int time_steps, const double t0, const double dt, Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves transient linear problem
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
Definition State.hpp:159
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:207
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:373
json solver_info
information of the solver, eg num iteration, time, errors, etc the informations varies depending on t...
Definition OutData.hpp:398
Eigen::Vector4d spectrum
spectrum of the stiffness matrix, enable only if POLYSOLVE_WITH_SPECTRA is ON (off by default)
Definition OutData.hpp:394
long long nn_zero
non zeros and sytem matrix size num dof is the total dof in the system
Definition OutData.hpp:412
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:42
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20