Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
StateInit.cpp
Go to the documentation of this file.
1#include <polyfem/State.hpp>
2
6
9
14
16
17#include <jse/jse.h>
18
19#include <spdlog/sinks/stdout_color_sinks.h>
20#include <spdlog/sinks/basic_file_sink.h>
21#include <spdlog/sinks/ostream_sink.h>
22
23#include <ipc/utils/logger.hpp>
24#include <wmtk/utils/Logger.hpp>
25
28
29#include <sstream>
30
31namespace spdlog::level
32{
34 spdlog::level::level_enum,
35 {{spdlog::level::level_enum::trace, "trace"},
36 {spdlog::level::level_enum::debug, "debug"},
37 {spdlog::level::level_enum::info, "info"},
38 {spdlog::level::level_enum::warn, "warning"},
39 {spdlog::level::level_enum::err, "error"},
40 {spdlog::level::level_enum::critical, "critical"},
41 {spdlog::level::level_enum::off, "off"},
42 {spdlog::level::level_enum::trace, 0},
43 {spdlog::level::level_enum::debug, 1},
44 {spdlog::level::level_enum::info, 2},
45 {spdlog::level::level_enum::warn, 3},
46 {spdlog::level::level_enum::err, 3},
47 {spdlog::level::level_enum::critical, 4},
48 {spdlog::level::level_enum::off, 5}})
49}
50
51namespace polyfem
52{
53 using namespace problem;
54 using namespace utils;
55
57 {
58 using namespace polysolve;
59
61
63 }
64
66 const std::string &log_file,
67 const spdlog::level::level_enum log_level,
68 const spdlog::level::level_enum file_log_level,
69 const bool is_quiet)
70 {
71 std::vector<spdlog::sink_ptr> sinks;
72
73 if (!is_quiet)
74 {
75 console_sink_ = std::make_shared<spdlog::sinks::stdout_color_sink_mt>();
76 sinks.emplace_back(console_sink_);
77 }
78
79 if (!log_file.empty())
80 {
81 file_sink_ = std::make_shared<spdlog::sinks::basic_file_sink_mt>(log_file, /*truncate=*/true);
82 // Set the file sink separately from the console so it can save all messages
83 file_sink_->set_level(file_log_level);
84 sinks.push_back(file_sink_);
85 }
86
87 init_logger(sinks, log_level);
88 spdlog::flush_every(std::chrono::seconds(3));
89 }
90
91 void State::init_logger(std::ostream &os, const spdlog::level::level_enum log_level)
92 {
93 std::vector<spdlog::sink_ptr> sinks;
94 sinks.emplace_back(std::make_shared<spdlog::sinks::ostream_sink_mt>(os, false));
95 init_logger(sinks, log_level);
96 }
97
99 const std::vector<spdlog::sink_ptr> &sinks,
100 const spdlog::level::level_enum log_level)
101 {
102 set_logger(std::make_shared<spdlog::logger>("polyfem", sinks.begin(), sinks.end()));
104
105 ipc::set_logger(std::make_shared<spdlog::logger>("ipctk", sinks.begin(), sinks.end()));
106
107 wmtk::set_logger(std::make_shared<spdlog::logger>("wmtk", sinks.begin(), sinks.end()));
108
109 // Set the logger at the lowest level, so all messages are passed to the sinks
110 logger().set_level(spdlog::level::trace);
111 ipc::logger().set_level(spdlog::level::trace);
112 wmtk::logger().set_level(spdlog::level::trace);
113
114 set_log_level(log_level);
115 }
116
117 void State::set_log_level(const spdlog::level::level_enum log_level)
118 {
119 spdlog::set_level(log_level);
120 if (console_sink_)
121 {
122 // Set only the level of the console
123 console_sink_->set_level(log_level); // Shared by all loggers
124 }
125 else
126 {
127 // Set the level of all sinks
128 logger().set_level(log_level);
129 ipc::logger().set_level(log_level);
130 }
131 }
132
133 void State::init(const json &p_args_in, const bool strict_validation)
134 {
135 json args_in = p_args_in; // mutable copy
136
137 apply_common_params(args_in);
138
139 // CHECK validity json
140 json rules;
141 jse::JSE jse;
142 {
143 jse.strict = strict_validation;
144 const std::string polyfem_input_spec = POLYFEM_INPUT_SPEC;
145 std::ifstream file(polyfem_input_spec);
146
147 if (file.is_open())
148 file >> rules;
149 else
150 {
151 logger().error("unable to open {} rules", polyfem_input_spec);
152 throw std::runtime_error("Invald spec file");
153 }
154
155 jse.include_directories.push_back(POLYFEM_JSON_SPEC_DIR);
156 jse.include_directories.push_back(POLYSOLVE_JSON_SPEC_DIR);
157 rules = jse.inject_include(rules);
158
159 polysolve::linear::Solver::apply_default_solver(rules, "/solver/linear");
160 polysolve::linear::Solver::apply_default_solver(rules, "/solver/adjoint_linear");
161 }
162
163 polysolve::linear::Solver::select_valid_solver(args_in["solver"]["linear"], logger());
164 if (args_in["solver"]["adjoint_linear"].is_null())
165 args_in["solver"]["adjoint_linear"] = args_in["solver"]["linear"];
166 else
167 polysolve::linear::Solver::select_valid_solver(args_in["solver"]["adjoint_linear"], logger());
168
169 // Use the /solver/nonlinear settings as the default for /solver/augmented_lagrangian/nonlinear
170 if (args_in.contains("/solver/nonlinear"_json_pointer))
171 {
172 if (args_in.contains("/solver/augmented_lagrangian/nonlinear"_json_pointer))
173 {
174 assert(args_in["solver"]["augmented_lagrangian"]["nonlinear"].is_object());
175 // Merge the augmented lagrangian settings into the nonlinear settings,
176 // and then replace the augmented lagrangian settings with the merged settings.
177 json nonlinear = args_in["solver"]["nonlinear"]; // copy
178 nonlinear.merge_patch(args_in["solver"]["augmented_lagrangian"]["nonlinear"]);
179 args_in["solver"]["augmented_lagrangian"]["nonlinear"] = nonlinear;
180 }
181 else
182 {
183 // Copy the nonlinear settings to the augmented_lagrangian settings
184 args_in["solver"]["augmented_lagrangian"]["nonlinear"] = args_in["solver"]["nonlinear"];
185 }
186 }
187 const bool valid_input = jse.verify_json(args_in, rules);
188
189 if (!valid_input)
190 {
191 logger().error("invalid input json:\n{}", jse.log2str());
192 throw std::runtime_error("Invald input json file");
193 }
194 // end of check
195
196 this->args = jse.inject_defaults(args_in, rules);
197 units.init(this->args["units"]);
198
199 if (!args_in.contains("/space/advanced/bc_method"_json_pointer) && this->args["space"]["basis_type"] != "Lagrange")
200 {
201 logger().warn("Setting bc method to lsq for non-Lagrange basis");
202 this->args["space"]["advanced"]["bc_method"] = "lsq";
203 }
204
205 // Save output directory and resolve output paths dynamically
206 const std::string output_dir = resolve_input_path(this->args["output"]["directory"]);
207 if (!output_dir.empty())
208 {
209 std::filesystem::create_directories(output_dir);
210 }
211 this->output_dir = output_dir;
212
213 std::string out_path_log = this->args["output"]["log"]["path"];
214 if (!out_path_log.empty())
215 {
216 out_path_log = resolve_output_path(out_path_log);
217 }
218
220 out_path_log,
221 this->args["output"]["log"]["level"],
222 this->args["output"]["log"]["file_level"],
223 this->args["output"]["log"]["quiet"]);
224
225 logger().info("Saving output to {}", output_dir);
226
227 const unsigned int thread_in = this->args["solver"]["max_threads"];
228 set_max_threads(thread_in);
229
230 has_dhat = args_in["contact"].contains("dhat");
231
232 init_time();
233
234 if (is_contact_enabled())
235 {
236 if (args["solver"]["contact"]["friction_iterations"] == 0)
237 {
238 logger().info("specified friction_iterations is 0; disabling friction");
239 args["contact"]["friction_coefficient"] = 0.0;
240 }
241 else if (args["solver"]["contact"]["friction_iterations"] < 0)
242 {
243 args["solver"]["contact"]["friction_iterations"] = std::numeric_limits<int>::max();
244 }
245 if (args["contact"]["friction_coefficient"] == 0.0)
246 {
247 args["solver"]["contact"]["friction_iterations"] = 0;
248 }
249 }
250 else
251 {
252 args["solver"]["contact"]["friction_iterations"] = 0;
253 args["contact"]["friction_coefficient"] = 0;
254 args["contact"]["periodic"] = false;
255 }
256
257 const std::string formulation = this->formulation();
259 assert(assembler->name() == formulation);
260 mass_matrix_assembler = std::make_shared<assembler::Mass>();
262
263 if (!other_name.empty())
264 {
267 }
268
269 if (args["solver"]["advanced"]["check_inversion"] == "Conservative")
270 {
271 if (auto elastic_assembler = std::dynamic_pointer_cast<assembler::ElasticityAssembler>(assembler))
272 elastic_assembler->set_use_robust_jacobian();
273 }
274
275 if (!args.contains("preset_problem"))
276 {
277 if (!assembler->is_tensor())
278 problem = std::make_shared<assembler::GenericScalarProblem>("GenericScalar");
279 else
280 problem = std::make_shared<assembler::GenericTensorProblem>("GenericTensor");
281
282 problem->clear();
283 if (!args["time"].is_null())
284 {
285 const auto tmp = R"({"is_time_dependent": true})"_json;
286 problem->set_parameters(tmp);
287 }
288 // important for the BC
289
290 auto bc = args["boundary_conditions"];
291 bc["root_path"] = root_path();
292 problem->set_parameters(bc);
293 problem->set_parameters(args["initial_conditions"]);
294
295 problem->set_parameters(args["output"]);
296 }
297 else
298 {
299 if (args["preset_problem"]["type"] == "Kernel")
300 {
301 problem = std::make_shared<KernelProblem>("Kernel", *assembler);
302 problem->clear();
303 KernelProblem &kprob = *dynamic_cast<KernelProblem *>(problem.get());
304 }
305 else
306 {
307 problem = ProblemFactory::factory().get_problem(args["preset_problem"]["type"]);
308 problem->clear();
309 }
310 // important for the BC
311 problem->set_parameters(args["preset_problem"]);
312 }
313
314 problem->set_units(*assembler, units);
315
317 {
318 if (is_contact_enabled())
319 {
320 if (!args["contact"]["use_convergent_formulation"])
321 {
322 args["contact"]["use_convergent_formulation"] = true;
323 logger().info("Use convergent formulation for differentiable contact...");
324 }
325 if (args["/solver/contact/barrier_stiffness"_json_pointer].is_string())
326 {
327 logger().error("Only constant barrier stiffness is supported in differentiable contact!");
328 }
329 }
330
331 if (args.contains("boundary_conditions") && args["boundary_conditions"].contains("rhs"))
332 {
333 json rhs = args["boundary_conditions"]["rhs"];
334 if ((rhs.is_array() && rhs.size() > 0 && rhs[0].is_string()) || rhs.is_string())
335 logger().error("Only constant rhs over space is supported in differentiable code!");
336 }
337 }
338 }
339
340 void State::set_max_threads(const int max_threads)
341 {
342 NThread::get().set_num_threads(max_threads);
343 }
344
346 {
347 if (!is_param_valid(args, "time"))
348 return;
349
350 const double t0 = Units::convert(args["time"]["t0"], units.time());
351 double tend, dt;
352 int time_steps;
353
354 // from "tend", "dt", "time_steps" only two can be used at a time
355 const int num_valid = is_param_valid(args["time"], "tend")
356 + is_param_valid(args["time"], "dt")
357 + is_param_valid(args["time"], "time_steps");
358 if (num_valid < 2)
359 {
360 log_and_throw_error("Exactly two of (tend, dt, time_steps) must be specified");
361 }
362 else if (num_valid == 2)
363 {
364 if (is_param_valid(args["time"], "tend"))
365 {
366 tend = Units::convert(args["time"]["tend"], units.time());
367 assert(tend > t0);
368 if (is_param_valid(args["time"], "dt"))
369 {
370 dt = Units::convert(args["time"]["dt"], units.time());
371 assert(dt > 0);
372 time_steps = int(ceil((tend - t0) / dt));
373 assert(time_steps > 0);
374 }
375 else if (is_param_valid(args["time"], "time_steps"))
376 {
377 time_steps = args["time"]["time_steps"];
378 assert(time_steps > 0);
379 dt = (tend - t0) / time_steps;
380 assert(dt > 0);
381 }
382 else
383 {
384 throw std::runtime_error("This code should be unreachable!");
385 }
386 }
387 else if (is_param_valid(args["time"], "dt"))
388 {
389 // tend is already confirmed to be invalid, so time_steps must be valid
390 assert(is_param_valid(args["time"], "time_steps"));
391
392 dt = Units::convert(args["time"]["dt"], units.time());
393 assert(dt > 0);
394
395 time_steps = args["time"]["time_steps"];
396 assert(time_steps > 0);
397
398 tend = t0 + time_steps * dt;
399 }
400 else
401 {
402 // tend and dt are already confirmed to be invalid
403 throw std::runtime_error("This code should be unreachable!");
404 }
405 }
406 else if (num_valid == 3)
407 {
408 tend = Units::convert(args["time"]["tend"], units.time());
409 dt = Units::convert(args["time"]["dt"], units.time());
410 time_steps = args["time"]["time_steps"];
411
412 // Check that all parameters agree
413 if (abs(t0 + dt * time_steps - tend) > 1e-12)
414 {
415 log_and_throw_error("Exactly two of (tend, dt, time_steps) must be specified");
416 }
417 }
418
419 // Store these for use later
420 args["time"]["tend"] = tend;
421 args["time"]["dt"] = dt;
422 args["time"]["time_steps"] = time_steps;
423
425
426 logger().info("t0={}, dt={}, tend={}", t0, dt, tend);
427 }
428
429 void State::set_materials(std::vector<std::shared_ptr<assembler::Assembler>> &assemblers) const
430 {
431 const int size = (assembler->is_tensor() || assembler->is_fluid()) ? mesh->dimension() : 1;
432 for (auto &a : assemblers)
433 a->set_size(size);
434
435 if (!utils::is_param_valid(args, "materials"))
436 return;
437
438 if (!args["materials"].is_array() && args["materials"]["type"] == "AMIPSAutodiff")
439 {
440 json transform_params = {};
441 transform_params["canonical_transformation"] = json::array();
442 if (!mesh->is_volume())
443 {
444 Eigen::MatrixXd regular_tri(3, 3);
445 regular_tri << 0, 0, 1,
446 1, 0, 1,
447 1. / 2., std::sqrt(3) / 2., 1;
448 regular_tri.transposeInPlace();
449 Eigen::MatrixXd regular_tri_inv = regular_tri.inverse();
450
451 const auto &mesh2d = *dynamic_cast<mesh::Mesh2D *>(mesh.get());
452 for (int e = 0; e < mesh->n_elements(); e++)
453 {
454 Eigen::MatrixXd transform;
455 mesh2d.compute_face_jacobian(e, regular_tri_inv, transform);
456 transform_params["canonical_transformation"].push_back(json({
457 {
458 transform(0, 0),
459 transform(0, 1),
460 },
461 {
462 transform(1, 0),
463 transform(1, 1),
464 },
465 }));
466 }
467 }
468 else
469 {
470 Eigen::MatrixXd regular_tet(4, 4);
471 regular_tet << 0, 0, 0, 1,
472 1, 0, 0, 1,
473 1. / 2., std::sqrt(3) / 2., 0, 1,
474 1. / 2., 1. / 2. / std::sqrt(3), std::sqrt(3) / 2., 1;
475 regular_tet.transposeInPlace();
476 Eigen::MatrixXd regular_tet_inv = regular_tet.inverse();
477
478 const auto &mesh3d = *dynamic_cast<mesh::Mesh3D *>(mesh.get());
479 for (int e = 0; e < mesh->n_elements(); e++)
480 {
481 Eigen::MatrixXd transform;
482 mesh3d.compute_cell_jacobian(e, regular_tet_inv, transform);
483 transform_params["canonical_transformation"].push_back(json({
484 {
485 transform(0, 0),
486 transform(0, 1),
487 transform(0, 2),
488 },
489 {
490 transform(1, 0),
491 transform(1, 1),
492 transform(1, 2),
493 },
494 {
495 transform(2, 0),
496 transform(2, 1),
497 transform(2, 2),
498 },
499 }));
500 }
501 }
502 transform_params["solve_displacement"] = true;
503 assembler->set_materials({}, transform_params, units);
504
505 return;
506 }
507
508 std::vector<int> body_ids(mesh->n_elements());
509 for (int i = 0; i < mesh->n_elements(); ++i)
510 body_ids[i] = mesh->get_body_id(i);
511
512 for (auto &a : assemblers)
513 a->set_materials(body_ids, args["materials"], units);
514 }
515
517 {
518 const int size = (this->assembler->is_tensor() || this->assembler->is_fluid()) ? this->mesh->dimension() : 1;
519 assembler.set_size(size);
520
521 if (!utils::is_param_valid(args, "materials"))
522 return;
523
524 std::vector<int> body_ids(mesh->n_elements());
525 for (int i = 0; i < mesh->n_elements(); ++i)
526 body_ids[i] = mesh->get_body_id(i);
527
528 assembler.set_materials(body_ids, args["materials"], units);
529 }
530
531} // namespace polyfem
std::string formulation() const
return the formulation (checks if the problem is scalar or not and deals with multiphysics)
Definition State.cpp:328
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 ".")
std::shared_ptr< assembler::Assembler > assembler
assemblers
Definition State.hpp:155
std::string resolve_output_path(const std::string &path) const
Resolve output path relative to output_dir if the path is not absolute.
std::string output_dir
Directory for output files.
Definition State.hpp:573
solver::CacheLevel optimization_enabled
Definition State.hpp:647
void set_materials(std::vector< std::shared_ptr< assembler::Assembler > > &assemblers) const
set the material and the problem dimension
void set_max_threads(const int max_threads=std::numeric_limits< int >::max())
std::shared_ptr< assembler::Assembler > pressure_assembler
Definition State.hpp:160
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:157
bool has_dhat
stores if input json contains dhat
Definition State.hpp:566
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:466
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
spdlog::sink_ptr file_sink_
Definition State.hpp:143
spdlog::sink_ptr console_sink_
logger sink to stdout
Definition State.hpp:142
State()
Constructor.
Definition StateInit.cpp:56
json args
main input arguments containing all defaults
Definition State.hpp:101
void set_log_level(const spdlog::level::level_enum log_level)
change log level
void init_time()
initialize time settings if args contains "time"
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.
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
bool is_contact_enabled() const
does the simulation has contact
Definition State.hpp:551
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
Definition State.hpp:159
Eigen::MatrixXd rhs
System right-hand side.
Definition State.hpp:207
void init(const json &json)
Definition Units.cpp:9
double characteristic_length() const
Definition Units.hpp:22
static double convert(const json &val, const std::string &unit_type)
Definition Units.cpp:31
const std::string & time() const
Definition Units.hpp:21
virtual bool is_tensor() const
virtual void set_size(const int size)
Definition Assembler.hpp:64
virtual bool is_fluid() const
static std::shared_ptr< MixedAssembler > make_mixed_assembler(const std::string &formulation)
static std::string other_assembler_name(const std::string &formulation)
static std::shared_ptr< Assembler > make_assembler(const std::string &formulation)
void compute_face_jacobian(const int el_id, const Eigen::MatrixXd &reference_map, Eigen::MatrixXd &jacobian) const
Definition Mesh2D.cpp:101
void compute_cell_jacobian(const int el_id, const Eigen::MatrixXd &reference_map, Eigen::MatrixXd &jacobian) const
Definition Mesh3D.cpp:424
static const ProblemFactory & factory()
std::shared_ptr< assembler::Problem > get_problem(const std::string &problem) const
void set_logger(spdlog::logger &logger)
static GeogramUtils & instance()
static NThread & get()
Definition par_for.hpp:19
void set_num_threads(const int max_threads)
Definition par_for.hpp:27
NLOHMANN_JSON_SERIALIZE_ENUM(CollisionProxyTessellation, {{CollisionProxyTessellation::REGULAR, "regular"}, {CollisionProxyTessellation::IRREGULAR, "irregular"}})
bool is_param_valid(const json &params, const std::string &key)
Determine if a key exists and is non-null in a json object.
void apply_common_params(json &args)
Definition JSONUtils.cpp:14
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
nlohmann::json json
Definition Common.hpp:9
void set_logger(std::shared_ptr< spdlog::logger > p_logger)
Setup a logger object to be used by Polyfem.
Definition Logger.cpp:60
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71