Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Optimizations.cpp
Go to the documentation of this file.
1#include "Optimizations.hpp"
2
4#include <jse/jse.h>
5
7
8#include <polyfem/State.hpp>
18
21
23
27
28#include <polysolve/nonlinear/BoxConstraintSolver.hpp>
29
30namespace spdlog::level
31{
33 spdlog::level::level_enum,
34 {{spdlog::level::level_enum::trace, "trace"},
35 {spdlog::level::level_enum::debug, "debug"},
36 {spdlog::level::level_enum::info, "info"},
37 {spdlog::level::level_enum::warn, "warning"},
38 {spdlog::level::level_enum::err, "error"},
39 {spdlog::level::level_enum::critical, "critical"},
40 {spdlog::level::level_enum::off, "off"},
41 {spdlog::level::level_enum::trace, 0},
42 {spdlog::level::level_enum::debug, 1},
43 {spdlog::level::level_enum::info, 2},
44 {spdlog::level::level_enum::warn, 3},
45 {spdlog::level::level_enum::err, 3},
46 {spdlog::level::level_enum::critical, 4},
47 {spdlog::level::level_enum::off, 5}})
48}
49
50namespace polyfem::solver
51{
52 namespace
53 {
54 bool load_json(const std::string &json_file, json &out)
55 {
56 std::ifstream file(json_file);
57
58 if (!file.is_open())
59 return false;
60
61 file >> out;
62
63 out["root_path"] = json_file;
64
65 return true;
66 }
67 } // namespace
68
69 std::shared_ptr<polysolve::nonlinear::Solver> AdjointOptUtils::make_nl_solver(const json &solver_params, const json &linear_solver_params, const double characteristic_length)
70 {
71 auto names = polysolve::nonlinear::Solver::available_solvers();
72 if (std::find(names.begin(), names.end(), solver_params["solver"]) != names.end())
73 return polysolve::nonlinear::Solver::create(solver_params, linear_solver_params, characteristic_length, adjoint_logger());
74
75 names = polysolve::nonlinear::BoxConstraintSolver::available_solvers();
76 if (std::find(names.begin(), names.end(), solver_params["solver"]) != names.end())
77 return polysolve::nonlinear::BoxConstraintSolver::create(solver_params, linear_solver_params, characteristic_length, adjoint_logger());
78
79 log_and_throw_adjoint_error("Invalid nonlinear solver name!");
80 }
81
82 std::shared_ptr<AdjointForm> AdjointOptUtils::create_simple_form(const std::string &obj_type, const std::string &param_type, const std::shared_ptr<State> &state, const json &args)
83 {
84 std::shared_ptr<AdjointForm> obj;
87 if (obj_type == "compliance")
88 obj = std::make_shared<ComplianceForm>(var2sim, *state, args);
89 else if (obj_type == "acceleration")
90 obj = std::make_shared<AccelerationForm>(var2sim, *state, args);
91 else if (obj_type == "kinetic")
92 obj = std::make_shared<AccelerationForm>(var2sim, *state, args);
93 else if (obj_type == "position")
94 obj = std::make_shared<PositionForm>(var2sim, *state, args);
95 else if (obj_type == "stress")
96 obj = std::make_shared<StressForm>(var2sim, *state, args);
97 else if (obj_type == "stress_norm")
98 obj = std::make_shared<StressNormForm>(var2sim, *state, args);
99 else if (obj_type == "dirichlet_energy")
100 obj = std::make_shared<DirichletEnergyForm>(var2sim, *state, args);
101 else if (obj_type == "elastic_energy")
102 obj = std::make_shared<ElasticEnergyForm>(var2sim, *state, args);
103 else if (obj_type == "max_stress")
104 obj = std::make_shared<MaxStressForm>(var2sim, *state, args);
105 else if (obj_type == "volume")
106 obj = std::make_shared<VolumeForm>(var2sim, *state, args);
107 else if (obj_type == "min_jacobian")
108 obj = std::make_shared<MinJacobianForm>(var2sim, *state);
109 else if (obj_type == "AMIPS")
110 obj = std::make_shared<AMIPSForm>(var2sim, *state);
111 else if (obj_type == "layer_thickness")
112 obj = std::make_shared<LayerThicknessForm>(var2sim, *state, args["boundary_ids"].get<std::vector<int>>(), args["dhat"]);
113 else if (obj_type == "layer_thickness_log")
114 obj = std::make_shared<LayerThicknessForm>(var2sim, *state, args["boundary_ids"].get<std::vector<int>>(), args["dhat"], true, args["dmin"]);
115 else if (obj_type == "function-target")
116 {
117 std::shared_ptr<TargetForm> tmp = std::make_shared<TargetForm>(var2sim, *state, args);
118 tmp->set_reference(args["target_function"], args["target_function_gradient"]);
119 obj = tmp;
120 }
121 else if (obj_type == "boundary_smoothing")
122 {
123 if (args["surface_selection"].is_array())
124 obj = std::make_shared<BoundarySmoothingForm>(var2sim, *state, args["scale_invariant"], args["power"], args["surface_selection"].get<std::vector<int>>());
125 else
126 obj = std::make_shared<BoundarySmoothingForm>(var2sim, *state, args["scale_invariant"], args["power"], std::vector<int>{args["surface_selection"].get<int>()});
127 }
128 else
129 log_and_throw_adjoint_error("Invalid simple form type {}!", obj_type);
130
131 return obj;
132 }
133
134 std::shared_ptr<AdjointForm> AdjointOptUtils::create_form(const json &args, const VariableToSimulationGroup &var2sim, const std::vector<std::shared_ptr<State>> &states)
135 {
136 std::shared_ptr<AdjointForm> obj;
137 if (args.is_array())
138 {
139 std::vector<std::shared_ptr<AdjointForm>> forms;
140 for (const auto &arg : args)
141 forms.push_back(create_form(arg, var2sim, states));
142
143 obj = std::make_shared<SumCompositeForm>(var2sim, forms);
144 }
145 else
146 {
147 const std::string type = args["type"];
148 if (type == "transient_integral")
149 {
150 std::shared_ptr<StaticForm> static_obj = std::dynamic_pointer_cast<StaticForm>(create_form(args["static_objective"], var2sim, states));
151 if (!static_obj)
152 log_and_throw_adjoint_error("Transient integral objective must have a static objective!");
153 const auto &state = states[args["state"]];
154 obj = std::make_shared<TransientForm>(var2sim, state->args["time"]["time_steps"], state->args["time"]["dt"], args["integral_type"], args["steps"].get<std::vector<int>>(), static_obj);
155 }
156 else if (type == "proxy_transient_integral")
157 {
158 std::shared_ptr<StaticForm> static_obj = std::dynamic_pointer_cast<StaticForm>(create_form(args["static_objective"], var2sim, states));
159 if (!static_obj)
160 log_and_throw_adjoint_error("Transient integral objective must have a static objective!");
161 if (args["steps"].size() == 0)
162 log_and_throw_adjoint_error("ProxyTransientForm requires non-empty \"steps\"!");
163 const auto &state = states[args["state"]];
164 obj = std::make_shared<ProxyTransientForm>(var2sim, state->args["time"]["time_steps"], state->args["time"]["dt"], args["integral_type"], args["steps"].get<std::vector<int>>(), static_obj);
165 }
166 else if (type == "power")
167 {
168 std::shared_ptr<AdjointForm> obj_aux = create_form(args["objective"], var2sim, states);
169 obj = std::make_shared<PowerForm>(obj_aux, args["power"]);
170 }
171 else if (type == "divide")
172 {
173 std::shared_ptr<AdjointForm> obj1 = create_form(args["objective"][0], var2sim, states);
174 std::shared_ptr<AdjointForm> obj2 = create_form(args["objective"][1], var2sim, states);
175 std::vector<std::shared_ptr<AdjointForm>> objs({obj1, obj2});
176 obj = std::make_shared<DivideForm>(objs);
177 }
178 else if (type == "plus-const")
179 {
180 obj = std::make_shared<PlusConstCompositeForm>(create_form(args["objective"], var2sim, states), args["value"]);
181 }
182 else if (type == "compliance")
183 {
184 obj = std::make_shared<ComplianceForm>(var2sim, *(states[args["state"]]), args);
185 }
186 else if (type == "acceleration")
187 {
188 obj = std::make_shared<AccelerationForm>(var2sim, *(states[args["state"]]), args);
189 }
190 else if (type == "kinetic")
191 {
192 obj = std::make_shared<AccelerationForm>(var2sim, *(states[args["state"]]), args);
193 }
194 else if (type == "target")
195 {
196 std::shared_ptr<TargetForm> tmp = std::make_shared<TargetForm>(var2sim, *(states[args["state"]]), args);
197 auto reference_cached = args["reference_cached_body_ids"].get<std::vector<int>>();
198 tmp->set_reference(states[args["target_state"]], std::set(reference_cached.begin(), reference_cached.end()));
199 obj = tmp;
200 }
201 else if (type == "displacement-target")
202 {
203 std::shared_ptr<TargetForm> tmp = std::make_shared<TargetForm>(var2sim, *(states[args["state"]]), args);
204 Eigen::VectorXd target_displacement;
205 target_displacement.setZero(states[args["state"]]->mesh->dimension());
206 if (target_displacement.size() != args["target_displacement"].size())
207 log_and_throw_error("Target displacement shape must match the dimension of the simulation");
208 for (int i = 0; i < target_displacement.size(); ++i)
209 target_displacement(i) = args["target_displacement"][i].get<double>();
210 if (args["active_dimension"].size() > 0)
211 {
212 if (target_displacement.size() != args["active_dimension"].size())
213 log_and_throw_error("Active dimension shape must match the dimension of the simulation");
214 std::vector<bool> active_dimension_mask(args["active_dimension"].size());
215 for (int i = 0; i < args["active_dimension"].size(); ++i)
216 active_dimension_mask[i] = args["active_dimension"][i].get<bool>();
217 tmp->set_active_dimension(active_dimension_mask);
218 }
219 tmp->set_reference(target_displacement);
220 obj = tmp;
221 }
222 else if (type == "center-target")
223 {
224 obj = std::make_shared<BarycenterTargetForm>(var2sim, args, states[args["state"]], states[args["target_state"]]);
225 }
226 else if (type == "sdf-target")
227 {
228 std::shared_ptr<SDFTargetForm> tmp = std::make_shared<SDFTargetForm>(var2sim, *(states[args["state"]]), args);
229 double delta = args["delta"].get<double>();
230 if (!states[args["state"]]->mesh->is_volume())
231 {
232 int dim = 2;
233 Eigen::MatrixXd control_points(args["control_points"].size(), dim);
234 for (int i = 0; i < control_points.rows(); ++i)
235 for (int j = 0; j < control_points.cols(); ++j)
236 control_points(i, j) = args["control_points"][i][j].get<double>();
237 Eigen::VectorXd knots(args["knots"].size());
238 for (int i = 0; i < knots.size(); ++i)
239 knots(i) = args["knots"][i].get<double>();
240 tmp->set_bspline_target(control_points, knots, delta);
241 }
242 else
243 {
244 int dim = 3;
245 Eigen::MatrixXd control_points_grid(args["control_points_grid"].size(), dim);
246 for (int i = 0; i < control_points_grid.rows(); ++i)
247 for (int j = 0; j < control_points_grid.cols(); ++j)
248 control_points_grid(i, j) = args["control_points_grid"][i][j].get<double>();
249 Eigen::VectorXd knots_u(args["knots_u"].size());
250 for (int i = 0; i < knots_u.size(); ++i)
251 knots_u(i) = args["knots_u"][i].get<double>();
252 Eigen::VectorXd knots_v(args["knots_v"].size());
253 for (int i = 0; i < knots_v.size(); ++i)
254 knots_v(i) = args["knots_v"][i].get<double>();
255 tmp->set_bspline_target(control_points_grid, knots_u, knots_v, delta);
256 }
257
258 obj = tmp;
259 }
260 else if (type == "mesh-target")
261 {
262 std::shared_ptr<MeshTargetForm> tmp = std::make_shared<MeshTargetForm>(var2sim, *(states[args["state"]]), args);
263 double delta = args["delta"].get<double>();
264 Eigen::MatrixXd V;
265 Eigen::MatrixXi E, F;
266 bool read = polyfem::io::OBJReader::read(args["mesh_path"], V, E, F);
267 if (!read)
268 log_and_throw_error(fmt::format("Could not read mesh! {}", args["mesh"]));
269 tmp->set_surface_mesh_target(V, F, delta);
270 obj = tmp;
271 }
272 else if (type == "function-target")
273 {
274 std::shared_ptr<TargetForm> tmp = std::make_shared<TargetForm>(var2sim, *(states[args["state"]]), args);
275 tmp->set_reference(args["target_function"], args["target_function_gradient"]);
276 obj = tmp;
277 }
278 else if (type == "node-target")
279 {
280 obj = std::make_shared<NodeTargetForm>(*(states[args["state"]]), var2sim, args);
281 }
282 else if (type == "min-dist-target")
283 {
284 obj = std::make_shared<MinTargetDistForm>(var2sim, args["steps"], args["target"], args, states[args["state"]]);
285 }
286 else if (type == "position")
287 {
288 obj = std::make_shared<PositionForm>(var2sim, *(states[args["state"]]), args);
289 }
290 else if (type == "stress")
291 {
292 obj = std::make_shared<StressForm>(var2sim, *(states[args["state"]]), args);
293 }
294 else if (type == "stress_norm")
295 {
296 obj = std::make_shared<StressNormForm>(var2sim, *(states[args["state"]]), args);
297 }
298 else if (type == "dirichlet_energy")
299 {
300 obj = std::make_shared<DirichletEnergyForm>(var2sim, *(states[args["state"]]), args);
301 }
302 else if (type == "elastic_energy")
303 {
304 obj = std::make_shared<ElasticEnergyForm>(var2sim, *(states[args["state"]]), args);
305 }
306 else if (type == "quadratic_contact_force_norm")
307 {
308 obj = std::make_shared<ProxyContactForceForm>(var2sim, *(states[args["state"]]), args["dhat"], true, args);
309 }
310 else if (type == "log_contact_force_norm")
311 {
312 obj = std::make_shared<ProxyContactForceForm>(var2sim, *(states[args["state"]]), args["dhat"], false, args);
313 }
314 else if (type == "max_stress")
315 {
316 obj = std::make_shared<MaxStressForm>(var2sim, *(states[args["state"]]), args);
317 }
318 else if (type == "volume")
319 {
320 obj = std::make_shared<VolumeForm>(var2sim, *(states[args["state"]]), args);
321 }
322 else if (type == "soft_constraint")
323 {
324 std::vector<std::shared_ptr<AdjointForm>> forms({create_form(args["objective"], var2sim, states)});
325 Eigen::VectorXd bounds = args["soft_bound"];
326 obj = std::make_shared<InequalityConstraintForm>(forms, bounds, args["power"]);
327 }
328 else if (type == "min_jacobian")
329 {
330 obj = std::make_shared<MinJacobianForm>(var2sim, *(states[args["state"]]));
331 }
332 else if (type == "AMIPS")
333 {
334 obj = std::make_shared<AMIPSForm>(var2sim, *(states[args["state"]]));
335 }
336 else if (type == "boundary_smoothing")
337 {
338 if (args["surface_selection"].is_array())
339 obj = std::make_shared<BoundarySmoothingForm>(var2sim, *(states[args["state"]]), args["scale_invariant"], args["power"], args["surface_selection"].get<std::vector<int>>());
340 else
341 obj = std::make_shared<BoundarySmoothingForm>(var2sim, *(states[args["state"]]), args["scale_invariant"], args["power"], std::vector<int>{args["surface_selection"].get<int>()});
342 }
343 else if (type == "collision_barrier")
344 {
345 obj = std::make_shared<CollisionBarrierForm>(var2sim, *(states[args["state"]]), args["dhat"]);
346 }
347 else if (type == "layer_thickness")
348 {
349 obj = std::make_shared<LayerThicknessForm>(var2sim, *(states[args["state"]]), args["boundary_ids"].get<std::vector<int>>(), args["dhat"]);
350 }
351 else if (type == "layer_thickness_log")
352 {
353 obj = std::make_shared<LayerThicknessForm>(var2sim, *(states[args["state"]]), args["boundary_ids"].get<std::vector<int>>(), args["dhat"], true, args["dmin"]);
354 }
355 else if (type == "deformed_collision_barrier")
356 {
357 obj = std::make_shared<DeformedCollisionBarrierForm>(var2sim, *(states[args["state"]]), args["dhat"]);
358 }
359 else if (type == "parametrized_product")
360 {
361 std::vector<std::shared_ptr<Parametrization>> map_list;
362 for (const auto &arg : args["parametrization"])
363 map_list.push_back(create_parametrization(arg, states, {}));
364 obj = std::make_shared<ParametrizedProductForm>(CompositeParametrization(std::move(map_list)));
365 }
366 else
367 log_and_throw_adjoint_error("Objective not implemented!");
368
369 obj->set_weight(args["weight"]);
370 if (args["print_energy"].get<std::string>() != "")
371 obj->enable_energy_print(args["print_energy"]);
372 }
373
374 return obj;
375 }
376
377 std::shared_ptr<Parametrization> AdjointOptUtils::create_parametrization(const json &args, const std::vector<std::shared_ptr<State>> &states, const std::vector<int> &variable_sizes)
378 {
379 std::shared_ptr<Parametrization> map;
380 const std::string type = args["type"];
381 if (type == "per-body-to-per-elem")
382 {
383 map = std::make_shared<PerBody2PerElem>(*(states[args["state"]]->mesh));
384 }
385 else if (type == "per-body-to-per-node")
386 {
387 map = std::make_shared<PerBody2PerNode>(*(states[args["state"]]->mesh), states[args["state"]]->bases, states[args["state"]]->n_bases);
388 }
389 else if (type == "E-nu-to-lambda-mu")
390 {
391 map = std::make_shared<ENu2LambdaMu>(args["is_volume"]);
392 }
393 else if (type == "slice")
394 {
395 if (args["from"] != -1 || args["to"] != -1)
396 map = std::make_shared<SliceMap>(args["from"], args["to"], args["last"]);
397 else if (args["parameter_index"] != -1)
398 {
399 int idx = args["parameter_index"].get<int>();
400 int from, to, last;
401 int cumulative = 0;
402 for (int i = 0; i < variable_sizes.size(); ++i)
403 {
404 if (i == idx)
405 {
406 from = cumulative;
407 to = from + variable_sizes[i];
408 }
409 cumulative += variable_sizes[i];
410 }
411 last = cumulative;
412 map = std::make_shared<SliceMap>(from, to, last);
413 }
414 else
415 log_and_throw_adjoint_error("Incorrect spec for SliceMap!");
416 }
417 else if (type == "exp")
418 {
419 map = std::make_shared<ExponentialMap>(args["from"], args["to"]);
420 }
421 else if (type == "scale")
422 {
423 map = std::make_shared<Scaling>(args["value"]);
424 }
425 else if (type == "power")
426 {
427 map = std::make_shared<PowerMap>(args["power"]);
428 }
429 else if (type == "append-values")
430 {
431 Eigen::VectorXd vals = args["values"];
432 map = std::make_shared<InsertConstantMap>(vals, args["start"]);
433 }
434 else if (type == "append-const")
435 {
436 map = std::make_shared<InsertConstantMap>(args["size"], args["value"], args["start"]);
437 }
438 else if (type == "linear-filter")
439 {
440 map = std::make_shared<LinearFilter>(*(states[args["state"]]->mesh), args["radius"]);
441 }
442 else if (type == "bounded-biharmonic-weights")
443 {
444 map = std::make_shared<BoundedBiharmonicWeights2Dto3D>(args["num_control_vertices"], args["num_vertices"], *states[args["state"]], args["allow_rotations"]);
445 }
446 else if (type == "scalar-velocity-parametrization")
447 {
448 map = std::make_shared<ScalarVelocityParametrization>(args["start_val"], args["dt"]);
449 }
450 else
451 log_and_throw_adjoint_error("Unkown parametrization!");
452
453 return map;
454 }
455
456 std::unique_ptr<VariableToSimulation> AdjointOptUtils::create_variable_to_simulation(const json &args, const std::vector<std::shared_ptr<State>> &states, const std::vector<int> &variable_sizes)
457 {
458 const std::string type = args["type"];
459
460 std::vector<std::shared_ptr<State>> cur_states;
461 if (args["state"].is_array())
462 for (int i : args["state"])
463 cur_states.push_back(states[i]);
464 else
465 cur_states.push_back(states[args["state"]]);
466
467 std::vector<std::shared_ptr<Parametrization>> map_list;
468 for (const auto &arg : args["composition"])
469 map_list.push_back(create_parametrization(arg, states, variable_sizes));
470
471 std::unique_ptr<VariableToSimulation> var2sim = VariableToSimulation::create(type, cur_states, CompositeParametrization(std::move(map_list)));
472 var2sim->set_output_indexing(args);
473
474 return var2sim;
475 }
476
477 Eigen::VectorXd AdjointOptUtils::inverse_evaluation(const json &args, const int ndof, const std::vector<int> &variable_sizes, VariableToSimulationGroup &var2sim)
478 {
479 Eigen::VectorXd x;
480 x.setZero(ndof);
481 int accumulative = 0;
482 int var = 0;
483 for (const auto &arg : args)
484 {
485 const auto &arg_initial = arg["initial"];
486 Eigen::VectorXd tmp(variable_sizes[var]);
487 if (arg_initial.is_array() && arg_initial.size() > 0)
488 {
489 tmp = arg_initial;
490 x.segment(accumulative, tmp.size()) = tmp;
491 }
492 else if (arg_initial.is_number())
493 {
494 tmp.setConstant(arg_initial.get<double>());
495 x.segment(accumulative, tmp.size()) = tmp;
496 }
497 else // arg["initial"] is empty array
498 x += var2sim[var]->inverse_eval();
499
500 accumulative += tmp.size();
501 var++;
502 }
503
504 return x;
505 }
506
507 std::shared_ptr<State> AdjointOptUtils::create_state(const json &args, CacheLevel level, const size_t max_threads)
508 {
509 std::shared_ptr<State> state = std::make_shared<State>();
510 state->set_max_threads(max_threads);
511
512 json in_args = args;
513 in_args["solver"]["max_threads"] = max_threads;
514 if (!args.contains("output") || !args["output"].contains("log") || !args["output"]["log"].contains("level"))
515 {
516 const json tmp = R"({
517 "output": {
518 "log": {
519 "level": "error"
520 }
521 }
522 })"_json;
523
524 in_args.merge_patch(tmp);
525 }
526
527 state->optimization_enabled = level;
528 state->init(in_args, true);
529 state->load_mesh();
530 Eigen::MatrixXd sol, pressure;
531 state->build_basis();
532 state->assemble_rhs();
533 state->assemble_mass_mat();
534
535 return state;
536 }
537
538 std::vector<std::shared_ptr<State>> AdjointOptUtils::create_states(const json &state_args, const CacheLevel &level, const size_t max_threads)
539 {
540 std::vector<std::shared_ptr<State>> states(state_args.size());
541 int i = 0;
542 for (const json &args : state_args)
543 {
544 json cur_args;
545 if (!load_json(args["path"], cur_args))
546 log_and_throw_adjoint_error("Can't find json for State {}", i);
547
548 states[i++] = AdjointOptUtils::create_state(cur_args, level, max_threads);
549 }
550 return states;
551 }
552
554 {
555 state.assemble_rhs();
556 state.assemble_mass_mat();
557 Eigen::MatrixXd sol, pressure;
558 state.solve_problem(sol, pressure);
559 }
560
561 void apply_objective_json_spec(json &args, const json &rules)
562 {
563 if (args.is_array())
564 {
565 for (auto &arg : args)
566 apply_objective_json_spec(arg, rules);
567 }
568 else if (args.is_object())
569 {
570 jse::JSE jse;
571 const bool valid_input = jse.verify_json(args, rules);
572
573 if (!valid_input)
574 {
575 logger().error("invalid objective json:\n{}", jse.log2str());
576 throw std::runtime_error("Invald objective json file");
577 }
578
579 args = jse.inject_defaults(args, rules);
580
581 for (auto &it : args.items())
582 {
583 if (it.key().find("objective") != std::string::npos)
584 apply_objective_json_spec(it.value(), rules);
585 }
586 }
587 }
588
589 json AdjointOptUtils::apply_opt_json_spec(const json &input_args, bool strict_validation)
590 {
591 json args_in = input_args;
592
593 // CHECK validity json
594 json rules;
595 jse::JSE jse;
596 {
597 jse.strict = strict_validation;
598 std::ifstream file(POLYFEM_OPT_INPUT_SPEC);
599
600 if (file.is_open())
601 file >> rules;
602 else
603 {
604 logger().error("unable to open {} rules", POLYFEM_OPT_INPUT_SPEC);
605 throw std::runtime_error("Invald spec file");
606 }
607
608 jse.include_directories.push_back(POLYFEM_JSON_SPEC_DIR);
609 jse.include_directories.push_back(POLYSOLVE_JSON_SPEC_DIR);
610 rules = jse.inject_include(rules);
611
612 // polysolve::linear::Solver::apply_default_solver(rules, "/solver/linear");
613 }
614
615 // polysolve::linear::Solver::select_valid_solver(args_in["solver"]["linear"], logger());
616
617 const bool valid_input = jse.verify_json(args_in, rules);
618
619 if (!valid_input)
620 {
621 logger().error("invalid input json:\n{}", jse.log2str());
622 throw std::runtime_error("Invald input json file");
623 }
624
625 json args = jse.inject_defaults(args_in, rules);
626
627 json obj_rules;
628 {
629 const std::string polyfem_objective_spec = POLYFEM_OBJECTIVE_INPUT_SPEC;
630 std::ifstream file(polyfem_objective_spec);
631
632 if (file.is_open())
633 file >> obj_rules;
634 else
635 {
636 logger().error("unable to open {} rules", polyfem_objective_spec);
637 throw std::runtime_error("Invald spec file");
638 }
639 }
640 apply_objective_json_spec(args["functionals"], obj_rules);
641
642 if (args.contains("stopping_conditions"))
643 apply_objective_json_spec(args["stopping_conditions"], obj_rules);
644
645 return args;
646 }
647
648 int AdjointOptUtils::compute_variable_size(const json &args, const std::vector<std::shared_ptr<State>> &states)
649 {
650 if (args["number"].is_number())
651 {
652 return args["number"].get<int>();
653 }
654 else if (args["number"].is_null() && args["initial"].size() > 0)
655 {
656 return args["initial"].size();
657 }
658 else if (args["number"].is_object())
659 {
660 auto selection = args["number"];
661 if (selection.contains("surface_selection"))
662 {
663 auto surface_selection = selection["surface_selection"].get<std::vector<int>>();
664 auto state_id = selection["state"];
665 std::set<int> node_ids = {};
666 for (const auto &surface : surface_selection)
667 {
668 std::vector<int> ids;
669 states[state_id]->compute_surface_node_ids(surface, ids);
670 for (const auto &i : ids)
671 node_ids.insert(i);
672 }
673 return node_ids.size() * states[state_id]->mesh->dimension();
674 }
675 else if (selection.contains("volume_selection"))
676 {
677 auto volume_selection = selection["volume_selection"].get<std::vector<int>>();
678 auto state_id = selection["state"];
679 std::set<int> node_ids = {};
680 for (const auto &volume : volume_selection)
681 {
682 std::vector<int> ids;
683 states[state_id]->compute_volume_node_ids(volume, ids);
684 for (const auto &i : ids)
685 node_ids.insert(i);
686 }
687
688 if (selection["exclude_boundary_nodes"])
689 {
690 std::vector<int> ids;
691 states[state_id]->compute_total_surface_node_ids(ids);
692 for (const auto &i : ids)
693 node_ids.erase(i);
694 }
695
696 return node_ids.size() * states[state_id]->mesh->dimension();
697 }
698 }
699
700 log_and_throw_adjoint_error("Incorrect specification for parameters.");
701 return -1;
702 }
703} // namespace polyfem::solver
int V
ElementAssemblyValues vals
Definition Assembler.cpp:22
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
void assemble_mass_mat()
assemble mass, step 4 of solve build mass matrix based on defined basis modifies mass (and maybe more...
Definition State.cpp:1487
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure)
solves the problems
Definition State.cpp:1683
void assemble_rhs()
compute rhs, step 3 of solve build rhs vector based on defined basis and given rhs of the problem mod...
Definition State.cpp:1611
static bool read(const std::string obj_file_name, std::vector< std::vector< double > > &V, std::vector< std::vector< double > > &TC, std::vector< std::vector< double > > &N, std::vector< std::vector< int > > &F, std::vector< std::vector< int > > &FTC, std::vector< std::vector< int > > &FN, std::vector< std::vector< int > > &L)
Read a mesh from an ascii obj file.
Definition OBJReader.cpp:32
A collection of VariableToSimulation.
static std::unique_ptr< VariableToSimulation > create(const std::string &type, const std::vector< std::shared_ptr< State > > &states, CompositeParametrization &&parametrization)
bool load_json(const std::string &json_file, json &out)
Definition main.cpp:26
NLOHMANN_JSON_SERIALIZE_ENUM(CollisionProxyTessellation, {{CollisionProxyTessellation::REGULAR, "regular"}, {CollisionProxyTessellation::IRREGULAR, "irregular"}})
void apply_objective_json_spec(json &args, const json &rules)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
spdlog::logger & adjoint_logger()
Retrieves the current logger for adjoint.
Definition Logger.cpp:28
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
void log_and_throw_error(const std::string &msg)
Definition Logger.cpp:71
static std::shared_ptr< AdjointForm > create_form(const json &args, const VariableToSimulationGroup &var2sim, const std::vector< std::shared_ptr< State > > &states)
static std::shared_ptr< AdjointForm > create_simple_form(const std::string &obj_type, const std::string &param_type, const std::shared_ptr< State > &state, const json &args)
static std::shared_ptr< polysolve::nonlinear::Solver > make_nl_solver(const json &solver_params, const json &linear_solver_params, const double characteristic_length)
static std::unique_ptr< VariableToSimulation > create_variable_to_simulation(const json &args, const std::vector< std::shared_ptr< State > > &states, const std::vector< int > &variable_sizes)
static json apply_opt_json_spec(const json &input_args, bool strict_validation)
static Eigen::VectorXd inverse_evaluation(const json &args, const int ndof, const std::vector< int > &variable_sizes, VariableToSimulationGroup &var2sim)
static void solve_pde(State &state)
static std::shared_ptr< Parametrization > create_parametrization(const json &args, const std::vector< std::shared_ptr< State > > &states, const std::vector< int > &variable_sizes)
static int compute_variable_size(const json &args, const std::vector< std::shared_ptr< State > > &states)
static std::vector< std::shared_ptr< State > > create_states(const json &state_args, const CacheLevel &level, const size_t max_threads)
static std::shared_ptr< State > create_state(const json &args, CacheLevel level, const size_t max_threads)