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