PolyFEM
Loading...
Searching...
No Matches
Optimizations.cpp
Go to the documentation of this file.
2
3#include <polyfem/State.hpp>
4#include <polyfem/Common.hpp>
5
20
23
25
28
29#include <polysolve/nonlinear/BoxConstraintSolver.hpp>
30
31#include <jse/jse.h>
32
33#include <Eigen/Core>
34
35#include <memory>
36#include <algorithm>
37#include <fstream>
38#include <vector>
39#include <string>
40#include <set>
41#include <stdexcept>
42
44{
45 NLOHMANN_JSON_SERIALIZE_ENUM(
46 spdlog::level::level_enum,
47 {{spdlog::level::level_enum::trace, "trace"},
48 {spdlog::level::level_enum::debug, "debug"},
49 {spdlog::level::level_enum::info, "info"},
50 {spdlog::level::level_enum::warn, "warning"},
51 {spdlog::level::level_enum::err, "error"},
52 {spdlog::level::level_enum::critical, "critical"},
53 {spdlog::level::level_enum::off, "off"},
54 {spdlog::level::level_enum::trace, 0},
55 {spdlog::level::level_enum::debug, 1},
56 {spdlog::level::level_enum::info, 2},
57 {spdlog::level::level_enum::warn, 3},
58 {spdlog::level::level_enum::err, 3},
59 {spdlog::level::level_enum::critical, 4},
60 {spdlog::level::level_enum::off, 5}})
61}
62
63namespace polyfem::solver
64{
65
66 std::shared_ptr<polysolve::nonlinear::Solver> AdjointOptUtils::make_nl_solver(const json &solver_params, const json &linear_solver_params, const double characteristic_length)
67 {
68 auto names = polysolve::nonlinear::Solver::available_solvers();
69 if (std::find(names.begin(), names.end(), solver_params["solver"]) != names.end())
70 return polysolve::nonlinear::Solver::create(solver_params, linear_solver_params, characteristic_length, adjoint_logger());
71
72 names = polysolve::nonlinear::BoxConstraintSolver::available_solvers();
73 if (std::find(names.begin(), names.end(), solver_params["solver"]) != names.end())
74 return polysolve::nonlinear::BoxConstraintSolver::create(solver_params, linear_solver_params, characteristic_length, adjoint_logger());
75
76 log_and_throw_adjoint_error("Invalid nonlinear solver name!");
77 }
78
79 Eigen::VectorXd AdjointOptUtils::inverse_evaluation(const json &args, const int ndof, const std::vector<int> &variable_sizes, VariableToSimulationGroup &var2sim)
80 {
81 Eigen::VectorXd x;
82 x.setZero(ndof);
83 int accumulative = 0;
84 int var = 0;
85 for (const auto &arg : args)
86 {
87 const auto &arg_initial = arg["initial"];
88 Eigen::VectorXd tmp(variable_sizes[var]);
89 if (arg_initial.is_array() && arg_initial.size() > 0)
90 {
91 tmp = arg_initial;
92 x.segment(accumulative, tmp.size()) = tmp;
93 }
94 else if (arg_initial.is_number())
95 {
96 tmp.setConstant(arg_initial.get<double>());
97 x.segment(accumulative, tmp.size()) = tmp;
98 }
99 else // arg["initial"] is empty array
100 x += var2sim.data[var]->inverse_eval();
101
102 accumulative += tmp.size();
103 var++;
104 }
105
106 return x;
107 }
108
110 {
111 state.assemble_rhs();
112 state.assemble_mass_mat();
113 Eigen::MatrixXd sol, pressure;
114 state.solve_problem(sol, pressure);
115 }
116
117 void apply_objective_json_spec(json &args, const json &rules)
118 {
119 if (args.is_array())
120 {
121 for (auto &arg : args)
122 apply_objective_json_spec(arg, rules);
123 }
124 else if (args.is_object())
125 {
126 jse::JSE jse;
127 const bool valid_input = jse.verify_json(args, rules);
128
129 if (!valid_input)
130 {
131 logger().error("invalid objective json:\n{}", jse.log2str());
132 throw std::runtime_error("Invald objective json file");
133 }
134
135 args = jse.inject_defaults(args, rules);
136
137 for (auto &it : args.items())
138 {
139 if (it.key().find("objective") != std::string::npos)
140 apply_objective_json_spec(it.value(), rules);
141 }
142 }
143 }
144
145 json AdjointOptUtils::apply_opt_json_spec(const json &input_args, bool strict_validation)
146 {
147 json args_in = input_args;
148
149 // CHECK validity json
150 json rules;
151 jse::JSE jse;
152 {
153 jse.strict = strict_validation;
154 std::ifstream file(POLYFEM_OPT_INPUT_SPEC);
155
156 if (file.is_open())
157 file >> rules;
158 else
159 {
160 logger().error("unable to open {} rules", POLYFEM_OPT_INPUT_SPEC);
161 throw std::runtime_error("Invald spec file");
162 }
163
164 jse.include_directories.push_back(POLYFEM_JSON_SPEC_DIR);
165 jse.include_directories.push_back(POLYSOLVE_JSON_SPEC_DIR);
166 rules = jse.inject_include(rules);
167
168 // polysolve::linear::Solver::apply_default_solver(rules, "/solver/linear");
169 }
170
171 // polysolve::linear::Solver::select_valid_solver(args_in["solver"]["linear"], logger());
172
173 const bool valid_input = jse.verify_json(args_in, rules);
174
175 if (!valid_input)
176 {
177 logger().error("invalid input json:\n{}", jse.log2str());
178 throw std::runtime_error("Invald input json file");
179 }
180
181 json args = jse.inject_defaults(args_in, rules);
182
183 json obj_rules;
184 {
185 const std::string polyfem_objective_spec = POLYFEM_OBJECTIVE_INPUT_SPEC;
186 std::ifstream file(polyfem_objective_spec);
187
188 if (file.is_open())
189 file >> obj_rules;
190 else
191 {
192 logger().error("unable to open {} rules", polyfem_objective_spec);
193 throw std::runtime_error("Invald spec file");
194 }
195 }
196 apply_objective_json_spec(args["functionals"], obj_rules);
197
198 if (args.contains("stopping_conditions"))
199 apply_objective_json_spec(args["stopping_conditions"], obj_rules);
200
201 return args;
202 }
203
204 int AdjointOptUtils::compute_variable_size(const json &args, const std::vector<std::shared_ptr<State>> &states)
205 {
206 if (args["number"].is_number())
207 {
208 return args["number"].get<int>();
209 }
210 else if (args["number"].is_null() && args["initial"].size() > 0)
211 {
212 return args["initial"].size();
213 }
214 else if (args["number"].is_object())
215 {
216 auto selection = args["number"];
217 if (selection.contains("surface_selection"))
218 {
219 auto surface_selection = selection["surface_selection"].get<std::vector<int>>();
220 auto state_id = selection["state"];
221 std::set<int> node_ids = {};
222 for (const auto &surface : surface_selection)
223 {
224 std::vector<int> ids;
225 compute_surface_node_ids(*states[state_id], surface, ids);
226 for (const auto &i : ids)
227 node_ids.insert(i);
228 }
229 return node_ids.size() * states[state_id]->mesh->dimension();
230 }
231 else if (selection.contains("volume_selection"))
232 {
233 auto volume_selection = selection["volume_selection"].get<std::vector<int>>();
234 auto state_id = selection["state"];
235 std::set<int> node_ids = {};
236 for (const auto &volume : volume_selection)
237 {
238 std::vector<int> ids;
239 compute_volume_node_ids(*states[state_id], volume, ids);
240 for (const auto &i : ids)
241 node_ids.insert(i);
242 }
243
244 if (selection["exclude_boundary_nodes"])
245 {
246 std::vector<int> ids;
247 compute_total_surface_node_ids(*states[state_id], ids);
248 for (const auto &i : ids)
249 node_ids.erase(i);
250 }
251
252 return node_ids.size() * states[state_id]->mesh->dimension();
253 }
254 }
255
256 log_and_throw_adjoint_error("Incorrect specification for parameters.");
257 return -1;
258 }
259} // namespace polyfem::solver
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:113
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:1461
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:1585
void solve_problem(Eigen::MatrixXd &sol, Eigen::MatrixXd &pressure, UserPostStepCallback user_post_step={}, const InitialConditionOverride *ic_override=nullptr)
solves the problems
Definition State.cpp:1657
A collection of VariableToSimulation.
std::vector< std::shared_ptr< VariableToSimulation > > data
void apply_objective_json_spec(json &args, const json &rules)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:44
void compute_volume_node_ids(const State &state, const int volume_selection, std::vector< int > &node_ids)
void compute_surface_node_ids(const State &state, const int surface_selection, std::vector< int > &node_ids)
spdlog::logger & adjoint_logger()
Retrieves the current logger for adjoint.
Definition Logger.cpp:30
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:79
void compute_total_surface_node_ids(const State &state, std::vector< int > &node_ids)
static std::shared_ptr< polysolve::nonlinear::Solver > make_nl_solver(const json &solver_params, const json &linear_solver_params, const double characteristic_length)
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 int compute_variable_size(const json &args, const std::vector< std::shared_ptr< State > > &states)