PolyFEM
Loading...
Searching...
No Matches
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#ifdef POLYFEM_WITH_ITR
25#include <wmtk/utils/Logger.hpp>
26#endif
27
30
31#include <sstream>
32
33namespace spdlog::level
34{
36 spdlog::level::level_enum,
37 {{spdlog::level::level_enum::trace, "trace"},
38 {spdlog::level::level_enum::debug, "debug"},
39 {spdlog::level::level_enum::info, "info"},
40 {spdlog::level::level_enum::warn, "warning"},
41 {spdlog::level::level_enum::err, "error"},
42 {spdlog::level::level_enum::critical, "critical"},
43 {spdlog::level::level_enum::off, "off"},
44 {spdlog::level::level_enum::trace, 0},
45 {spdlog::level::level_enum::debug, 1},
46 {spdlog::level::level_enum::info, 2},
47 {spdlog::level::level_enum::warn, 3},
48 {spdlog::level::level_enum::err, 3},
49 {spdlog::level::level_enum::critical, 4},
50 {spdlog::level::level_enum::off, 5}})
51}
52
53namespace polyfem
54{
55 using namespace problem;
56 using namespace utils;
57
59 {
60 using namespace polysolve;
61
63
65 }
66
68 const std::string &log_file,
69 const spdlog::level::level_enum log_level,
70 const spdlog::level::level_enum file_log_level,
71 const bool is_quiet)
72 {
73 std::vector<spdlog::sink_ptr> sinks;
74
75 if (!is_quiet)
76 {
77 console_sink_ = std::make_shared<spdlog::sinks::stdout_color_sink_mt>();
78 sinks.emplace_back(console_sink_);
79 }
80
81 if (!log_file.empty())
82 {
83 file_sink_ = std::make_shared<spdlog::sinks::basic_file_sink_mt>(log_file, /*truncate=*/true);
84 // Set the file sink separately from the console so it can save all messages
85 file_sink_->set_level(file_log_level);
86 sinks.push_back(file_sink_);
87 }
88
89 init_logger(sinks, log_level);
90 spdlog::flush_every(std::chrono::seconds(3));
91 }
92
93 void State::init_logger(std::ostream &os, const spdlog::level::level_enum log_level)
94 {
95 std::vector<spdlog::sink_ptr> sinks;
96 sinks.emplace_back(std::make_shared<spdlog::sinks::ostream_sink_mt>(os, false));
97 init_logger(sinks, log_level);
98 }
99
101 const std::vector<spdlog::sink_ptr> &sinks,
102 const spdlog::level::level_enum log_level)
103 {
104 set_logger(std::make_shared<spdlog::logger>("polyfem", sinks.begin(), sinks.end()));
106
107 ipc::set_logger(std::make_shared<spdlog::logger>("ipctk", sinks.begin(), sinks.end()));
108
109#ifdef POLYFEM_WITH_ITR
110 wmtk::set_logger(std::make_shared<spdlog::logger>("wmtk", sinks.begin(), sinks.end()));
111#endif
112
113 // Set the logger at the lowest level, so all messages are passed to the sinks
114 logger().set_level(spdlog::level::trace);
115 ipc::logger().set_level(spdlog::level::trace);
116#ifdef POLYFEM_WITH_ITR
117 wmtk::logger().set_level(spdlog::level::trace);
118#endif
119
120 set_log_level(log_level);
121 }
122
123 void State::set_log_level(const spdlog::level::level_enum log_level)
124 {
125 spdlog::set_level(log_level);
126 if (console_sink_)
127 {
128 // Set only the level of the console
129 console_sink_->set_level(log_level); // Shared by all loggers
130 }
131 else
132 {
133 // Set the level of all sinks
134 logger().set_level(log_level);
135 ipc::logger().set_level(log_level);
136 }
137 }
138
139 void State::init(const json &p_args_in, const bool strict_validation)
140 {
141 json args_in = p_args_in; // mutable copy
142
143 has_constraints_ = p_args_in.contains("constraints");
144
145 apply_common_params(args_in);
146
147 // CHECK validity json
148 json rules;
149 jse::JSE jse;
150 {
151 jse.strict = strict_validation;
152 const std::string polyfem_input_spec = POLYFEM_INPUT_SPEC;
153 std::ifstream file(polyfem_input_spec);
154
155 if (file.is_open())
156 file >> rules;
157 else
158 {
159 logger().error("unable to open {} rules", polyfem_input_spec);
160 throw std::runtime_error("Invalid spec file");
161 }
162
163 jse.include_directories.push_back(POLYFEM_JSON_SPEC_DIR);
164 jse.include_directories.push_back(POLYSOLVE_JSON_SPEC_DIR);
165 rules = jse.inject_include(rules);
166
167 polysolve::linear::Solver::apply_default_solver(rules, "/solver/linear");
168 polysolve::linear::Solver::apply_default_solver(rules, "/solver/adjoint_linear");
169 }
170
171 polysolve::linear::Solver::select_valid_solver(args_in["solver"]["linear"], logger());
172 if (args_in["solver"]["adjoint_linear"].is_null())
173 args_in["solver"]["adjoint_linear"] = args_in["solver"]["linear"];
174 else
175 polysolve::linear::Solver::select_valid_solver(args_in["solver"]["adjoint_linear"], logger());
176
177 // Use the /solver/nonlinear settings as the default for /solver/augmented_lagrangian/nonlinear
178 if (args_in.contains("/solver/nonlinear"_json_pointer))
179 {
180 if (args_in.contains("/solver/augmented_lagrangian/nonlinear"_json_pointer))
181 {
182 assert(args_in["solver"]["augmented_lagrangian"]["nonlinear"].is_object());
183 // Merge the augmented lagrangian settings into the nonlinear settings,
184 // and then replace the augmented lagrangian settings with the merged settings.
185 json nonlinear = args_in["solver"]["nonlinear"]; // copy
186 nonlinear.merge_patch(args_in["solver"]["augmented_lagrangian"]["nonlinear"]);
187 args_in["solver"]["augmented_lagrangian"]["nonlinear"] = nonlinear;
188 }
189 else
190 {
191 // Copy the nonlinear settings to the augmented_lagrangian settings
192 args_in["solver"]["augmented_lagrangian"]["nonlinear"] = args_in["solver"]["nonlinear"];
193 }
194 }
195 const bool valid_input = jse.verify_json(args_in, rules);
196
197 if (!valid_input)
198 {
199 logger().error("invalid input json:\n{}", jse.log2str());
200 throw std::runtime_error("Invalid input json file");
201 }
202 // end of check
203
204 this->args = jse.inject_defaults(args_in, rules);
205 units.init(this->args["units"]);
206
207 if (!args_in.contains("/space/advanced/bc_method"_json_pointer) && this->args["space"]["basis_type"] != "Lagrange")
208 {
209 logger().warn("Setting bc method to lsq for non-Lagrange basis");
210 this->args["space"]["advanced"]["bc_method"] = "lsq";
211 }
212
213 // Save output directory and resolve output paths dynamically
214 const std::string output_dir = resolve_input_path(this->args["output"]["directory"]);
215 if (!output_dir.empty())
216 {
217 std::filesystem::create_directories(output_dir);
218 }
219 this->output_dir = output_dir;
220
221 std::string out_path_log = this->args["output"]["log"]["path"];
222 if (!out_path_log.empty())
223 {
224 out_path_log = resolve_output_path(out_path_log);
225 }
226
227 for (auto &path : this->args["constraints"]["hard"])
228 {
229 path = resolve_input_path(path);
230 }
231 for (auto &path : this->args["constraints"]["soft"])
232 {
233 path["data"] = resolve_input_path(path["data"]);
234 }
235
237 out_path_log,
238 this->args["output"]["log"]["level"],
239 this->args["output"]["log"]["file_level"],
240 this->args["output"]["log"]["quiet"]);
241
242 logger().info("Saving output to {}", output_dir);
243
244 const unsigned int thread_in = this->args["solver"]["max_threads"];
245 set_max_threads(thread_in);
246
247 has_dhat = args_in["contact"].contains("dhat");
248
249 init_time();
250
251 if (is_contact_enabled())
252 {
253 if (args["solver"]["contact"]["friction_iterations"] == 0)
254 {
255 logger().info("specified friction_iterations is 0; disabling friction");
256 args["contact"]["friction_coefficient"] = 0.0;
257 }
258 else if (args["solver"]["contact"]["friction_iterations"] < 0)
259 {
260 args["solver"]["contact"]["friction_iterations"] = std::numeric_limits<int>::max();
261 }
262 if (args["contact"]["friction_coefficient"] == 0.0)
263 {
264 args["solver"]["contact"]["friction_iterations"] = 0;
265 }
266 }
267 else
268 {
269 args["solver"]["contact"]["friction_iterations"] = 0;
270 args["contact"]["friction_coefficient"] = 0;
271 args["contact"]["periodic"] = false;
272 }
273
274 const std::string formulation = this->formulation();
276 assert(assembler->name() == formulation);
277 mass_matrix_assembler = std::make_shared<assembler::Mass>();
278 pure_mass_matrix_assembler = std::make_shared<assembler::HRZMass>();
280
281 if (!other_name.empty())
282 {
285 }
286
287 if (args["solver"]["advanced"]["check_inversion"] == "Conservative")
288 {
289 if (auto elastic_assembler = std::dynamic_pointer_cast<assembler::ElasticityAssembler>(assembler))
290 elastic_assembler->set_use_robust_jacobian();
291 }
292
293 if (!args.contains("preset_problem"))
294 {
295 if (!assembler->is_tensor())
296 problem = std::make_shared<assembler::GenericScalarProblem>("GenericScalar");
297 else
298 problem = std::make_shared<assembler::GenericTensorProblem>("GenericTensor");
299
300 problem->clear();
301 if (!args["time"].is_null())
302 {
303 const auto tmp = R"({"is_time_dependent": true})"_json;
304 problem->set_parameters(tmp, root_path());
305 }
306 // important for the BC
307
308 auto bc = args["boundary_conditions"];
309 bc["root_path"] = root_path();
310 problem->set_parameters(bc, root_path());
311 problem->set_parameters(args["initial_conditions"], root_path());
312
313 problem->set_parameters(args["output"], root_path());
314 }
315 else
316 {
317 if (args["preset_problem"]["type"] == "Kernel")
318 {
319 problem = std::make_shared<KernelProblem>("Kernel", *assembler);
320 problem->clear();
321 KernelProblem &kprob = *dynamic_cast<KernelProblem *>(problem.get());
322 }
323 else
324 {
325 problem = ProblemFactory::factory().get_problem(args["preset_problem"]["type"]);
326 problem->clear();
327 }
328 // important for the BC
329 problem->set_parameters(args["preset_problem"], root_path());
330 }
331
332 problem->set_units(*assembler, units);
333 }
334
335 void State::set_max_threads(const int max_threads)
336 {
337 NThread::get().set_num_threads(max_threads);
338 }
339
341 {
342 if (!is_param_valid(args, "time"))
343 return;
344
345 const double t0 = Units::convert(args["time"]["t0"], units.time());
346 double tend, dt;
347 int time_steps;
348
349 // from "tend", "dt", "time_steps" only two can be used at a time
350 const int num_valid = is_param_valid(args["time"], "tend")
351 + is_param_valid(args["time"], "dt")
352 + is_param_valid(args["time"], "time_steps");
353 if (num_valid < 2)
354 {
355 log_and_throw_error("Exactly two of (tend, dt, time_steps) must be specified");
356 }
357 else if (num_valid == 2)
358 {
359 if (is_param_valid(args["time"], "tend"))
360 {
361 tend = Units::convert(args["time"]["tend"], units.time());
362 assert(tend > t0);
363 if (is_param_valid(args["time"], "dt"))
364 {
365 dt = Units::convert(args["time"]["dt"], units.time());
366 assert(dt > 0);
367 time_steps = int(ceil((tend - t0) / dt));
368 assert(time_steps > 0);
369 }
370 else if (is_param_valid(args["time"], "time_steps"))
371 {
372 time_steps = args["time"]["time_steps"];
373 assert(time_steps > 0);
374 dt = (tend - t0) / time_steps;
375 assert(dt > 0);
376 }
377 else
378 {
379 throw std::runtime_error("This code should be unreachable!");
380 }
381 }
382 else if (is_param_valid(args["time"], "dt"))
383 {
384 // tend is already confirmed to be invalid, so time_steps must be valid
385 assert(is_param_valid(args["time"], "time_steps"));
386
387 dt = Units::convert(args["time"]["dt"], units.time());
388 assert(dt > 0);
389
390 time_steps = args["time"]["time_steps"];
391 assert(time_steps > 0);
392
393 tend = t0 + time_steps * dt;
394 }
395 else
396 {
397 // tend and dt are already confirmed to be invalid
398 throw std::runtime_error("This code should be unreachable!");
399 }
400 }
401 else if (num_valid == 3)
402 {
403 tend = Units::convert(args["time"]["tend"], units.time());
404 dt = Units::convert(args["time"]["dt"], units.time());
405 time_steps = args["time"]["time_steps"];
406
407 // Check that all parameters agree
408 if (abs(t0 + dt * time_steps - tend) > 1e-12)
409 {
410 log_and_throw_error("Exactly two of (tend, dt, time_steps) must be specified");
411 }
412 }
413
414 // Store these for use later
415 args["time"]["tend"] = tend;
416 args["time"]["dt"] = dt;
417 args["time"]["time_steps"] = time_steps;
418
420
421 logger().info("t0={}, dt={}, tend={}", t0, dt, tend);
422 }
423
424 void State::set_materials(std::vector<std::shared_ptr<assembler::Assembler>> &assemblers) const
425 {
426 const int size = (assembler->is_tensor() || assembler->is_fluid()) ? mesh->dimension() : 1;
427 for (auto &a : assemblers)
428 a->set_size(size);
429
430 if (!utils::is_param_valid(args, "materials"))
431 return;
432
433 std::vector<int> body_ids(mesh->n_elements());
434 for (int i = 0; i < mesh->n_elements(); ++i)
435 body_ids[i] = mesh->get_body_id(i);
436
437 for (auto &a : assemblers)
438 a->set_materials(body_ids, args["materials"], units, root_path());
439 }
440
442 {
443 const int size = (this->assembler->is_tensor() || this->assembler->is_fluid()) ? this->mesh->dimension() : 1;
444 assembler.set_size(size);
445
446 if (!utils::is_param_valid(args, "materials"))
447 return;
448
449 std::vector<int> body_ids(mesh->n_elements());
450 for (int i = 0; i < mesh->n_elements(); ++i)
451 body_ids[i] = mesh->get_body_id(i);
452
453 assembler.set_materials(body_ids, args["materials"], units, root_path());
454 }
455
456} // namespace polyfem
std::string formulation() const
return the formulation (checks if the problem is scalar or not and deals with multiphysics)
Definition State.cpp:332
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:190
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:738
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:196
std::shared_ptr< assembler::Mass > mass_matrix_assembler
Definition State.hpp:192
bool has_dhat
stores if input json contains dhat
Definition State.hpp:731
std::unique_ptr< mesh::Mesh > mesh
current mesh, it can be a Mesh2D or Mesh3D
Definition State.hpp:594
std::shared_ptr< assembler::HRZMass > pure_mass_matrix_assembler
Definition State.hpp:193
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:204
bool has_constraints_
Definition State.hpp:538
spdlog::sink_ptr file_sink_
Definition State.hpp:178
spdlog::sink_ptr console_sink_
logger sink to stdout
Definition State.hpp:177
State()
Constructor.
Definition StateInit.cpp:58
json args
main input arguments containing all defaults
Definition State.hpp:136
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:67
bool is_contact_enabled() const
does the simulation have contact
Definition State.hpp:708
std::shared_ptr< assembler::MixedAssembler > mixed_assembler
Definition State.hpp:195
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)
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:44
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:62
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:73