PolyFEM
Loading...
Searching...
No Matches
AdjointNLProblem.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
8
9#include <Eigen/Core>
10
11#include <memory>
12#include <vector>
13#include <fstream>
14
15namespace polyfem
16{
17 class State;
18}
19
20namespace polyfem::solver
21{
22 class AdjointForm;
23
25 {
26 public:
27 AdjointNLProblem(std::shared_ptr<AdjointForm> form,
28 const VariableToSimulationGroup &variables_to_simulation,
29 const std::vector<std::shared_ptr<State>> &all_states,
30 const std::vector<std::shared_ptr<DiffCache>> &all_diff_caches,
31 const json &args);
32
33 AdjointNLProblem(std::shared_ptr<AdjointForm> form,
34 const std::vector<std::shared_ptr<AdjointForm>> &stopping_conditions,
35 const VariableToSimulationGroup &variables_to_simulation,
36 const std::vector<std::shared_ptr<State>> &all_states,
37 const std::vector<std::shared_ptr<DiffCache>> &all_diff_caches,
38 const json &args);
39
40 double value(const Eigen::VectorXd &x) override;
41
42 void gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) override;
43 void hessian(const Eigen::VectorXd &x, StiffnessMatrix &hessian) override;
44 void save_to_file(const int iter_num, const Eigen::VectorXd &x0);
45 bool is_step_valid(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
46 bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
47 double max_step_size(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
48
49 void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
50 void line_search_end() override;
51 void post_step(const polysolve::nonlinear::PostStepData &data) override;
52 bool stop(const TVector &x) override;
53
54 void solution_changed(const Eigen::VectorXd &new_x) override;
55 bool after_line_search_custom_operation(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override;
56 void solve_pde();
57
58 private:
59 std::shared_ptr<AdjointForm> form_;
61 std::vector<std::shared_ptr<State>> all_states_;
62 std::vector<std::shared_ptr<DiffCache>> all_diff_caches_;
63 std::vector<bool> active_state_mask;
64 Eigen::VectorXd cur_grad;
65 Eigen::VectorXd curr_x;
66
68 std::ofstream solution_ostream;
69
72
74 std::vector<int> solve_in_order;
75
76 int save_iter = 0;
77
78 std::vector<std::shared_ptr<AdjointForm>> stopping_conditions_; // if all the stopping conditions are non-positive, stop the optimization
79 };
80} // namespace polyfem::solver
int x
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_
std::vector< std::shared_ptr< DiffCache > > all_diff_caches_
void hessian(const Eigen::VectorXd &x, StiffnessMatrix &hessian) override
void post_step(const polysolve::nonlinear::PostStepData &data) override
void solution_changed(const Eigen::VectorXd &new_x) override
bool is_step_collision_free(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
VariableToSimulationGroup variables_to_simulation_
double value(const Eigen::VectorXd &x) override
std::vector< std::shared_ptr< State > > all_states_
std::shared_ptr< AdjointForm > form_
void line_search_begin(const Eigen::VectorXd &x0, const Eigen::VectorXd &x1) override
A collection of VariableToSimulation.
nlohmann::json json
Definition Common.hpp:9
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24