PolyFEM
Loading...
Searching...
No Matches
AdjointNLProblem.cpp
Go to the documentation of this file.
2
10#include <polyfem/State.hpp>
12
13#include <list>
14#include <stack>
15
17{
18 namespace
19 {
20
21 Eigen::VectorXd get_updated_mesh_nodes(const VariableToSimulationGroup &variables_to_simulation, const std::shared_ptr<State> &curr_state, const Eigen::VectorXd &x)
22 {
23 Eigen::MatrixXd V;
24 curr_state->get_vertices(V);
25 Eigen::VectorXd X = utils::flatten(V);
26
27 for (auto &p : variables_to_simulation)
28 {
29 for (const auto &state : p->get_states())
30 if (state.get() != curr_state.get())
31 continue;
32 if (p->get_parameter_type() != ParameterType::Shape)
33 continue;
34 auto state_variable = p->get_parametrization().eval(x);
35 auto output_indexing = p->get_output_indexing(x);
36 for (int i = 0; i < output_indexing.size(); ++i)
37 X(output_indexing(i)) = state_variable(i);
38 }
39
40 return X;
41 }
42
43 using namespace std;
44 // Class to represent a graph
45 class Graph
46 {
47 int V; // No. of vertices'
48
49 // adjacency lists
50 vector<list<int>> adj;
51
52 // A function used by topologicalSort
53 void topologicalSortUtil(int v, vector<bool> &visited, stack<int> &Stack);
54
55 public:
56 Graph(int V); // Constructor
57
58 // function to add an edge to graph
59 void addEdge(int v, int w);
60
61 // prints a Topological Sort of the complete graph
62 vector<int> topologicalSort();
63 };
64
65 Graph::Graph(int V)
66 {
67 this->V = V;
68 adj.resize(V);
69 }
70
71 void Graph::addEdge(int v, int w)
72 {
73 adj[v].push_back(w); // Add w to v’s list.
74 }
75
76 // A recursive function used by topologicalSort
77 void Graph::topologicalSortUtil(int v, vector<bool> &visited,
78 stack<int> &Stack)
79 {
80 // Mark the current node as visited.
81 visited[v] = true;
82
83 // Recur for all the vertices adjacent to this vertex
84 list<int>::iterator i;
85 for (i = adj[v].begin(); i != adj[v].end(); ++i)
86 if (!visited[*i])
87 topologicalSortUtil(*i, visited, Stack);
88
89 // Push current vertex to stack which stores result
90 Stack.push(v);
91 }
92
93 // The function to do Topological Sort. It uses recursive
94 // topologicalSortUtil()
95 vector<int> Graph::topologicalSort()
96 {
97 stack<int> Stack;
98
99 // Mark all the vertices as not visited
100 vector<bool> visited(V, false);
101
102 // Call the recursive helper function to store Topological
103 // Sort starting from all vertices one by one
104 for (int i = 0; i < V; i++)
105 if (visited[i] == false)
106 topologicalSortUtil(i, visited, Stack);
107
108 // Print contents of stack
109 vector<int> sorted;
110 while (Stack.empty() == false)
111 {
112 sorted.push_back(Stack.top());
113 Stack.pop();
114 }
115
116 return sorted;
117 }
118 } // namespace
119
120 AdjointNLProblem::AdjointNLProblem(std::shared_ptr<AdjointForm> form, const VariableToSimulationGroup &variables_to_simulation, const std::vector<std::shared_ptr<State>> &all_states, const json &args)
121 : FullNLProblem({form}),
122 form_(form),
123 variables_to_simulation_(variables_to_simulation),
124 all_states_(all_states),
125 save_freq(args["output"]["save_frequency"]),
126 enable_slim(args["solver"]["advanced"]["enable_slim"]),
127 smooth_line_search(args["solver"]["advanced"]["smooth_line_search"]),
128 solve_in_parallel(args["solver"]["advanced"]["solve_in_parallel"])
129 {
130 cur_grad.setZero(0);
131
132 if (enable_slim && args["solver"]["nonlinear"]["advanced"]["apply_gradient_fd"] != "None")
133 adjoint_logger().warn("SLIM may affect the finite difference result!");
134
135 if (enable_slim && smooth_line_search)
136 adjoint_logger().warn("Both in-line-search SLIM and after-line-search SLIM are ON!");
137
138 if (args["output"]["solution"] != "")
139 {
140 solution_ostream.open(args["output"]["solution"].get<std::string>(), std::ofstream::out);
141 if (!solution_ostream.is_open())
142 adjoint_logger().error("Cannot open solution file for writing!");
143 }
144
145 solve_in_order.clear();
146 {
147 Graph G(all_states.size());
148 for (int k = 0; k < all_states.size(); k++)
149 {
150 auto &arg = args["states"][k];
151 if (arg["initial_guess"].get<int>() >= 0)
152 G.addEdge(arg["initial_guess"].get<int>(), k);
153 }
154
155 solve_in_order = G.topologicalSort();
156 }
157
158 active_state_mask.assign(all_states_.size(), false);
159 for (int i = 0; i < all_states_.size(); i++)
160 {
161 for (const auto &v2sim : variables_to_simulation_)
162 {
163 for (const auto &state : v2sim->get_states())
164 {
165 if (all_states_[i].get() == state.get())
166 {
167 active_state_mask[i] = true;
168 break;
169 }
170 }
171 }
172 }
173 }
174
175 AdjointNLProblem::AdjointNLProblem(std::shared_ptr<AdjointForm> form, const std::vector<std::shared_ptr<AdjointForm>> &stopping_conditions, const VariableToSimulationGroup &variables_to_simulation, const std::vector<std::shared_ptr<State>> &all_states, const json &args) : AdjointNLProblem(form, variables_to_simulation, all_states, args)
176 {
177 stopping_conditions_ = stopping_conditions;
178 }
179
180 void AdjointNLProblem::hessian(const Eigen::VectorXd &x, StiffnessMatrix &hessian)
181 {
182 log_and_throw_adjoint_error("Hessian not supported!");
183 }
184
185 double AdjointNLProblem::value(const Eigen::VectorXd &x)
186 {
187 return form_->value(x);
188 }
189
190 void AdjointNLProblem::gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv)
191 {
192 if (cur_grad.size() == x.size())
193 gradv = cur_grad;
194 else
195 {
196 gradv.setZero(x.size());
197
198 {
199 POLYFEM_SCOPED_TIMER("adjoint solve");
200 for (int i = 0; i < all_states_.size(); i++)
201 all_states_[i]->solve_adjoint_cached(form_->compute_reduced_adjoint_rhs(x, *all_states_[i])); // caches inside state
202 }
203
204 {
205 POLYFEM_SCOPED_TIMER("gradient assembly");
206 form_->first_derivative(x, gradv);
207 if (x.size() < 10)
208 {
209 adjoint_logger().trace("x {}", x.transpose());
210 adjoint_logger().trace("gradient {}", gradv.transpose());
211 }
212 }
213
214 cur_grad = gradv;
215 }
216 }
217
218 bool AdjointNLProblem::is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
219 {
220 bool need_rebuild_basis = false;
221
222 // update to new parameter and check if the new parameter is valid to solve
223 for (const auto &v : variables_to_simulation_)
224 if (v->get_parameter_type() == ParameterType::Shape || v->get_parameter_type() == ParameterType::PeriodicShape)
225 need_rebuild_basis = true;
226
227 if (need_rebuild_basis && smooth_line_search)
228 {
229 Eigen::MatrixXd X, V0, V1;
230 Eigen::MatrixXi F;
231
232 int state_count = 0;
233 for (auto state_ : all_states_)
234 {
235
237 get_updated_mesh_nodes(variables_to_simulation_, state_, x1),
238 state_->mesh->dimension());
239 state_->get_vertices(V0);
240 state_->get_elements(F);
241
242 Eigen::MatrixXd V_smooth;
243 bool slim_success = polyfem::mesh::apply_slim(V0, F, V1, V_smooth);
244 if (!slim_success)
245 {
246 adjoint_logger().info("SLIM failed, step not valid!");
247 return false;
248 }
249
250 V1 = V_smooth;
251
252 bool flipped = utils::is_flipped(V1, F);
253 if (flipped)
254 {
255 adjoint_logger().info("Found flipped element in LS, step not valid!");
256 return false;
257 }
258
259 state_count++;
260 }
261 }
262
263 return form_->is_step_valid(x0, x1);
264 }
265
266 bool AdjointNLProblem::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
267 {
268 return form_->is_step_collision_free(x0, x1);
269 }
270
271 double AdjointNLProblem::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
272 {
273 return form_->max_step_size(x0, x1);
274 }
275
276 void AdjointNLProblem::line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
277 {
278 form_->line_search_begin(x0, x1);
279 }
280
282 {
283 form_->line_search_end();
284 }
285
286 void AdjointNLProblem::post_step(const polysolve::nonlinear::PostStepData &data)
287 {
288 save_to_file(save_iter++, data.x);
289
290 form_->post_step(data);
291 }
292
293 void AdjointNLProblem::save_to_file(const int iter_num, const Eigen::VectorXd &x0)
294 {
295 int id = 0;
296
297 if (solution_ostream.is_open())
298 {
299 adjoint_logger().debug("Save solution at iteration {} to file...", iter_num);
300 solution_ostream << iter_num << ": " << std::setprecision(16) << x0.transpose() << std::endl;
301 solution_ostream.flush();
302 }
303
304 if (iter_num % save_freq != 0)
305 return;
306 adjoint_logger().info("Saving iteration {}", iter_num);
307 for (const auto &state : all_states_)
308 {
309 bool save_vtu = true;
310 bool save_rest_mesh = true;
311
312 std::string vis_mesh_path = state->resolve_output_path(fmt::format("opt_state_{:d}_iter_{:d}.vtu", id, iter_num));
313 std::string mesh_ext = state->mesh->is_volume() ? ".msh" : ".obj";
314 std::string rest_mesh_path = state->resolve_output_path(fmt::format("opt_state_{:d}_iter_{:d}" + mesh_ext, id, iter_num));
315 id++;
316
317 if (!save_vtu)
318 continue;
319 adjoint_logger().debug("Save final vtu to file {} ...", vis_mesh_path);
320
321 double tend = state->args.value("tend", 1.0);
322 double dt = 1;
323 if (!state->args["time"].is_null())
324 dt = state->args["time"]["dt"];
325
326 Eigen::MatrixXd sol = state->diff_cached.u(-1);
327
328 state->out_geom.save_vtu(
329 vis_mesh_path,
330 *state,
331 sol,
332 Eigen::MatrixXd::Zero(state->n_pressure_bases, 1),
333 tend, dt,
335 state->mesh->is_linear(),
336 state->mesh->has_prism(),
337 state->problem->is_scalar()),
338 state->is_contact_enabled());
339
340 if (!save_rest_mesh)
341 continue;
342 adjoint_logger().debug("Save rest mesh to file {} ...", rest_mesh_path);
343
344 // If shape opt, save rest meshes as well
345 Eigen::MatrixXd V;
346 Eigen::MatrixXi F;
347 state->get_vertices(V);
348 state->get_elements(F);
349 if (state->mesh->is_volume())
350 io::MshWriter::write(rest_mesh_path, V, F, state->mesh->get_body_ids(), true, false);
351 else
352 io::OBJWriter::write(rest_mesh_path, V, F);
353 }
354 }
355
356 void AdjointNLProblem::solution_changed(const Eigen::VectorXd &newX)
357 {
358 bool need_rebuild_basis = false;
359
360 // update to new parameter and check if the new parameter is valid to solve
361 for (const auto &v : variables_to_simulation_)
362 {
363 v->update(newX);
364 if (v->get_parameter_type() == ParameterType::Shape || v->get_parameter_type() == ParameterType::PeriodicShape)
365 need_rebuild_basis = true;
366 }
367
368 if (need_rebuild_basis)
369 {
370 for (const auto &state : all_states_)
371 state->build_basis();
372 }
373
374 // solve PDE
375 solve_pde();
376
377 form_->solution_changed(newX);
378
379 curr_x = newX;
380 }
381
382 bool AdjointNLProblem::after_line_search_custom_operation(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
383 {
384 if (!enable_slim)
385 return false;
386
387 for (const auto &v : variables_to_simulation_)
388 v->update(x0);
389
390 std::vector<Eigen::MatrixXd> V_old_list;
391 for (auto state : all_states_)
392 {
393 Eigen::MatrixXd V;
394 state->get_vertices(V);
395 V_old_list.push_back(V);
396 }
397
398 for (const auto &v : variables_to_simulation_)
399 v->update(x1);
400
401 // Apply slim to all states on a frequency
402 int state_num = 0;
403 for (auto state : all_states_)
404 {
405 Eigen::MatrixXd V_new, V_smooth;
406 Eigen::MatrixXi F;
407 state->get_vertices(V_new);
408 state->get_elements(F);
409
410 bool slim_success = polyfem::mesh::apply_slim(V_old_list[state_num++], F, V_new, V_smooth, 50);
411
412 if (!slim_success)
413 log_and_throw_adjoint_error("SLIM cannot succeed");
414
415 for (int i = 0; i < V_smooth.rows(); ++i)
416 state->set_mesh_vertex(i, V_smooth.row(i));
417 }
418 adjoint_logger().debug("SLIM succeeded!");
419
420 return true;
421 }
422
424 {
426 {
427 adjoint_logger().info("Run simulations in parallel...");
428
429 utils::maybe_parallel_for(all_states_.size(), [&](int start, int end, int thread_id) {
430 for (int i = start; i < end; i++)
431 {
432 auto state = all_states_[i];
433 if (active_state_mask[i] || state->diff_cached.size() == 0)
434 {
435 state->assemble_rhs();
436 state->assemble_mass_mat();
437 Eigen::MatrixXd sol, pressure; // solution is also cached in state
438 state->solve_problem(sol, pressure);
439 }
440 }
441 });
442 }
443 else
444 {
445 adjoint_logger().info("Run simulations in serial...");
446
447 Eigen::MatrixXd sol, pressure; // solution is also cached in state
448 for (int i : solve_in_order)
449 {
450 auto state = all_states_[i];
451 if (active_state_mask[i] || state->diff_cached.size() == 0)
452 {
453 state->assemble_rhs();
454 state->assemble_mass_mat();
455
456 state->solve_problem(sol, pressure);
457 }
458 }
459 }
460
461 cur_grad.resize(0);
462 }
463
464 bool AdjointNLProblem::stop(const TVector &x)
465 {
466 if (stopping_conditions_.size() == 0)
467 return false;
468
469 for (auto &obj : stopping_conditions_)
470 {
471 obj->solution_changed(x);
472 if (obj->value(x) > 0)
473 return false;
474 }
475 return true;
476 }
477
478} // namespace polyfem::solver
vector< list< int > > adj
int V
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
static void write(const std::string &path, const mesh::Mesh &mesh, const bool binary)
saves the mesh
Definition MshWriter.cpp:7
static bool write(const std::string &path, const Eigen::MatrixXd &v, const Eigen::MatrixXi &e, const Eigen::MatrixXi &f)
Definition OBJWriter.cpp:18
bool after_line_search_custom_operation(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
void gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) override
void save_to_file(const int iter_num, const Eigen::VectorXd &x0)
bool stop(const TVector &x) override
bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
std::vector< std::shared_ptr< AdjointForm > > stopping_conditions_
void hessian(const Eigen::VectorXd &x, StiffnessMatrix &hessian) override
void post_step(const polysolve::nonlinear::PostStepData &data) override
AdjointNLProblem(std::shared_ptr< AdjointForm > form, const VariableToSimulationGroup &variables_to_simulation, const std::vector< std::shared_ptr< State > > &all_states, const json &args)
void solution_changed(const Eigen::VectorXd &new_x) override
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
double value(const Eigen::VectorXd &x) override
std::vector< std::shared_ptr< State > > all_states_
std::shared_ptr< AdjointForm > form_
const VariableToSimulationGroup variables_to_simulation_
void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
A collection of VariableToSimulation.
bool apply_slim(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, const Eigen::MatrixXd &V_new, Eigen::MatrixXd &V_smooth, const int max_iters)
bool is_flipped(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
Determine if any simplex is inverted or collapses.
Eigen::MatrixXd unflatten(const Eigen::VectorXd &x, int dim)
Unflatten rowwises, so every dim elements in x become a row.
Eigen::VectorXd flatten(const Eigen::MatrixXd &X)
Flatten rowwises.
void maybe_parallel_for(int size, const std::function< void(int, int, int)> &partial_for)
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
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24