PolyFEM
Loading...
Searching...
No Matches
AdjointNLProblem.cpp
Go to the documentation of this file.
2
8#include <polyfem/State.hpp>
9#include <igl/boundary_facets.h>
10#include <igl/writeOBJ.h>
11
15
16#include <list>
17#include <stack>
18
19namespace polyfem::solver
20{
21 namespace
22 {
23 using namespace std;
24 // Class to represent a graph
25 class Graph
26 {
27 int V; // No. of vertices'
28
29 // adjacency lists
30 vector<list<int>> adj;
31
32 // A function used by topologicalSort
33 void topologicalSortUtil(int v, vector<bool> &visited, stack<int> &Stack);
34
35 public:
36 Graph(int V); // Constructor
37
38 // function to add an edge to graph
39 void addEdge(int v, int w);
40
41 // prints a Topological Sort of the complete graph
42 vector<int> topologicalSort();
43 };
44
45 Graph::Graph(int V)
46 {
47 this->V = V;
48 adj.resize(V);
49 }
50
51 void Graph::addEdge(int v, int w)
52 {
53 adj[v].push_back(w); // Add w to v’s list.
54 }
55
56 // A recursive function used by topologicalSort
57 void Graph::topologicalSortUtil(int v, vector<bool> &visited,
58 stack<int> &Stack)
59 {
60 // Mark the current node as visited.
61 visited[v] = true;
62
63 // Recur for all the vertices adjacent to this vertex
64 list<int>::iterator i;
65 for (i = adj[v].begin(); i != adj[v].end(); ++i)
66 if (!visited[*i])
67 topologicalSortUtil(*i, visited, Stack);
68
69 // Push current vertex to stack which stores result
70 Stack.push(v);
71 }
72
73 // The function to do Topological Sort. It uses recursive
74 // topologicalSortUtil()
75 vector<int> Graph::topologicalSort()
76 {
77 stack<int> Stack;
78
79 // Mark all the vertices as not visited
80 vector<bool> visited(V, false);
81
82 // Call the recursive helper function to store Topological
83 // Sort starting from all vertices one by one
84 for (int i = 0; i < V; i++)
85 if (visited[i] == false)
86 topologicalSortUtil(i, visited, Stack);
87
88 // Print contents of stack
89 vector<int> sorted;
90 while (Stack.empty() == false)
91 {
92 sorted.push_back(Stack.top());
93 Stack.pop();
94 }
95
96 return sorted;
97 }
98 } // namespace
99
100 AdjointNLProblem::AdjointNLProblem(std::shared_ptr<AdjointForm> form, const VariableToSimulationGroup &variables_to_simulation, const std::vector<std::shared_ptr<State>> &all_states, const json &args)
101 : FullNLProblem({form}),
102 form_(form),
103 variables_to_simulation_(variables_to_simulation),
104 all_states_(all_states),
105 save_freq(args["output"]["save_frequency"]),
106 solve_in_parallel(args["solver"]["advanced"]["solve_in_parallel"])
107 {
108 cur_grad.setZero(0);
109
110 if (args["output"]["solution"] != "")
111 {
112 solution_ostream.open(args["output"]["solution"].get<std::string>(), std::ofstream::out);
113 if (!solution_ostream.is_open())
114 adjoint_logger().error("Cannot open solution file for writing!");
115 }
116
117 solve_in_order.clear();
118 {
119 Graph G(all_states.size());
120 for (int k = 0; k < all_states.size(); k++)
121 {
122 auto &arg = args["states"][k];
123 if (arg["initial_guess"].get<int>() >= 0)
124 G.addEdge(arg["initial_guess"].get<int>(), k);
125 }
126
127 solve_in_order = G.topologicalSort();
128 }
129
130 active_state_mask.assign(all_states_.size(), false);
131 for (int i = 0; i < all_states_.size(); i++)
132 {
133 for (const auto &v2sim : variables_to_simulation_)
134 {
135 for (const auto &state : v2sim->get_states())
136 {
137 if (all_states_[i].get() == state.get())
138 {
139 active_state_mask[i] = true;
140 break;
141 }
142 }
143 }
144 }
145 }
146
147 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)
148 {
149 stopping_conditions_ = stopping_conditions;
150 }
151
152 void AdjointNLProblem::hessian(const Eigen::VectorXd &x, StiffnessMatrix &hessian)
153 {
154 log_and_throw_adjoint_error("Hessian not supported!");
155 }
156
157 double AdjointNLProblem::value(const Eigen::VectorXd &x)
158 {
159 return form_->value(x);
160 }
161
162 void AdjointNLProblem::gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv)
163 {
164 if (cur_grad.size() == x.size())
165 gradv = cur_grad;
166 else
167 {
168 gradv.setZero(x.size());
169
170 {
171 POLYFEM_SCOPED_TIMER("adjoint solve");
172 for (int i = 0; i < all_states_.size(); i++)
173 all_states_[i]->solve_adjoint_cached(form_->compute_reduced_adjoint_rhs(x, *all_states_[i])); // caches inside state
174 }
175
176 {
177 POLYFEM_SCOPED_TIMER("gradient assembly");
178 form_->first_derivative(x, gradv);
179 }
180
181 cur_grad = gradv;
182 }
183 }
184
185 bool AdjointNLProblem::is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
186 {
187 return form_->is_step_valid(x0, x1);
188 }
189
190 bool AdjointNLProblem::is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
191 {
192 return form_->is_step_collision_free(x0, x1);
193 }
194
195 double AdjointNLProblem::max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
196 {
197 return form_->max_step_size(x0, x1);
198 }
199
200 void AdjointNLProblem::line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1)
201 {
202 form_->line_search_begin(x0, x1);
203 }
204
206 {
207 form_->line_search_end();
208 }
209
210 void AdjointNLProblem::post_step(const polysolve::nonlinear::PostStepData &data)
211 {
212 save_to_file(save_iter++, data.x);
213
214 form_->post_step(data);
215 }
216
217 void AdjointNLProblem::save_to_file(const int iter_num, const Eigen::VectorXd &x0)
218 {
219 int id = 0;
220
221 if (solution_ostream.is_open())
222 {
223 adjoint_logger().debug("Save solution at iteration {} to file...", iter_num);
224 solution_ostream << iter_num << ": " << std::setprecision(16) << x0.transpose() << std::endl;
225 solution_ostream.flush();
226 }
227
228 if (iter_num % save_freq != 0)
229 return;
230 adjoint_logger().info("Saving iteration {}", iter_num);
231 for (const auto &state : all_states_)
232 {
233 bool save_vtu = true;
234 bool save_rest_mesh = true;
235
236 std::string vis_mesh_path = state->resolve_output_path(fmt::format("opt_state_{:d}_iter_{:d}.vtu", id, iter_num));
237 std::string rest_mesh_path = state->resolve_output_path(fmt::format("opt_state_{:d}_iter_{:d}.obj", id, iter_num));
238 id++;
239
240 if (!save_vtu)
241 continue;
242 adjoint_logger().debug("Save final vtu to file {} ...", vis_mesh_path);
243
244 double tend = state->args.value("tend", 1.0);
245 double dt = 1;
246 if (!state->args["time"].is_null())
247 dt = state->args["time"]["dt"];
248
249 Eigen::MatrixXd sol = state->diff_cached.u(-1);
250
251 state->out_geom.save_vtu(
252 vis_mesh_path,
253 *state,
254 sol,
255 Eigen::MatrixXd::Zero(state->n_pressure_bases, 1),
256 tend, dt,
257 io::OutGeometryData::ExportOptions(state->args, state->mesh->is_linear(), state->problem->is_scalar(), state->solve_export_to_file),
258 state->is_contact_enabled(),
259 state->solution_frames);
260
261 if (!save_rest_mesh)
262 continue;
263 adjoint_logger().debug("Save rest mesh to file {} ...", rest_mesh_path);
264
265 // If shape opt, save rest meshes as well
266 Eigen::MatrixXd V;
267 Eigen::MatrixXi F;
268 state->get_vertices(V);
269 state->get_elements(F);
270 if (state->mesh->dimension() == 3)
271 F = igl::boundary_facets<Eigen::MatrixXi, Eigen::MatrixXi>(F);
272
273 io::OBJWriter::write(rest_mesh_path, V, F);
274 }
275 }
276
277 void AdjointNLProblem::solution_changed(const Eigen::VectorXd &newX)
278 {
279 bool need_rebuild_basis = false;
280
281 // update to new parameter and check if the new parameter is valid to solve
282 for (const auto &v : variables_to_simulation_)
283 {
284 v->update(newX);
285 if (v->get_parameter_type() == ParameterType::Shape)
286 need_rebuild_basis = true;
287 }
288
289 if (need_rebuild_basis)
290 {
291 for (const auto &state : all_states_)
292 state->build_basis();
293 }
294
295 // solve PDE
296 solve_pde();
297
298 form_->solution_changed(newX);
299 }
300
302 {
304 {
305 adjoint_logger().info("Run simulations in parallel...");
306
307 utils::maybe_parallel_for(all_states_.size(), [&](int start, int end, int thread_id) {
308 for (int i = start; i < end; i++)
309 {
310 auto state = all_states_[i];
311 if (active_state_mask[i] || state->diff_cached.size() == 0)
312 {
313 state->assemble_rhs();
314 state->assemble_mass_mat();
315 Eigen::MatrixXd sol, pressure; // solution is also cached in state
316 state->solve_problem(sol, pressure);
317 }
318 }
319 });
320 }
321 else
322 {
323 adjoint_logger().info("Run simulations in serial...");
324
325 Eigen::MatrixXd sol, pressure; // solution is also cached in state
326 for (int i : solve_in_order)
327 {
328 auto state = all_states_[i];
329 if (active_state_mask[i] || state->diff_cached.size() == 0)
330 {
331 state->assemble_rhs();
332 state->assemble_mass_mat();
333
334 state->solve_problem(sol, pressure);
335 }
336 }
337 }
338
339 cur_grad.resize(0);
340 }
341
342 bool AdjointNLProblem::stop(const TVector &x)
343 {
344 if (stopping_conditions_.size() == 0)
345 return false;
346
347 for (auto &obj : stopping_conditions_)
348 {
349 obj->solution_changed(x);
350 if (obj->value(x) > 0)
351 return false;
352 }
353 return true;
354 }
355
356} // namespace polyfem::solver
vector< list< int > > adj
int V
int x
#define POLYFEM_SCOPED_TIMER(...)
Definition Timer.hpp:10
static bool write(const std::string &path, const Eigen::MatrixXd &v, const Eigen::MatrixXi &e, const Eigen::MatrixXi &f)
Definition OBJWriter.cpp:18
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.
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:28
nlohmann::json json
Definition Common.hpp:9
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20