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