PolyFEM
Loading...
Searching...
No Matches
AdjointForm.hpp
Go to the documentation of this file.
1#pragma once
2
5
6namespace polyfem::solver
7{
8 class AdjointForm : public Form
9 {
10 public:
11 AdjointForm(const VariableToSimulationGroup &variable_to_simulations) : variable_to_simulations_(variable_to_simulations) {}
12 virtual ~AdjointForm() {}
13
14 virtual std::string name() const override { return "adjoint"; }
15 void enable_energy_print(const std::string &print_energy_keyword);
16
17 double value(const Eigen::VectorXd &x) const override;
18 virtual void solution_changed(const Eigen::VectorXd &new_x) override;
19
21
22 virtual Eigen::MatrixXd compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const;
23 virtual void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const;
24 virtual void first_derivative(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override;
25 virtual Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const;
26
27 // not used functions from base class
28 virtual void update_quantities(const double t, const Eigen::VectorXd &x) final override;
29 virtual void init_lagging(const Eigen::VectorXd &x) final override;
30 virtual void update_lagging(const Eigen::VectorXd &x, const int iter_num) final override;
31
32 protected:
33 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override;
34 virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const final override;
35
37
38 enum class PrintStage
39 {
43 };
44
47 };
48
49 class StaticForm : public AdjointForm
50 {
51 public:
53 virtual ~StaticForm() = default;
54
55 virtual std::string name() const override { return "static"; }
56
57 Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const final override;
58 void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override;
59 double value_unweighted(const Eigen::VectorXd &x) const final override;
60
61 virtual Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const = 0;
62 virtual Eigen::VectorXd compute_adjoint_rhs_step_prev(const int time_step, const Eigen::VectorXd &x, const State &state) const;
63 virtual void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const = 0;
64 virtual double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const = 0;
65 virtual void solution_changed(const Eigen::VectorXd &new_x) final override;
66 virtual void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x) {}
67 virtual bool depends_on_step_prev() const final { return depends_on_step_prev_; }
68
69 protected:
71 };
72
74 {
75 public:
76 MaxStressForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const json &args) : StaticForm(variable_to_simulations), state_(state)
77 {
78 auto tmp_ids = args["volume_selection"].get<std::vector<int>>();
79 interested_ids_ = std::set(tmp_ids.begin(), tmp_ids.end());
80 }
81
82 std::string name() const override { return "max_stress"; }
83
84 Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override;
85 double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override;
86 void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
87
88 private:
89 std::set<int> interested_ids_;
90 const State &state_;
91 };
92} // namespace polyfem::solver
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
const VariableToSimulationGroup & get_variable_to_simulations() const
virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override
Compute the first derivative of the value wrt x.
virtual void solution_changed(const Eigen::VectorXd &new_x) override
Update cached fields upon a change in the solution.
virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const final override
Compute the second derivative of the value wrt x.
virtual void init_lagging(const Eigen::VectorXd &x) final override
Initialize lagged fields TODO: more than one step.
virtual void first_derivative(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override
Compute the first derivative of the value wrt x multiplied with the weigth.
virtual Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const
double value(const Eigen::VectorXd &x) const override
Compute the value of the form multiplied with the weigth.
virtual void update_lagging(const Eigen::VectorXd &x, const int iter_num) final override
Update lagged fields.
virtual void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
virtual Eigen::MatrixXd compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const
virtual void update_quantities(const double t, const Eigen::VectorXd &x) final override
Update time-dependent fields.
void enable_energy_print(const std::string &print_energy_keyword)
AdjointForm(const VariableToSimulationGroup &variable_to_simulations)
virtual std::string name() const override
const VariableToSimulationGroup variable_to_simulations_
MaxStressForm(const VariableToSimulationGroup &variable_to_simulations, const State &state, const json &args)
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const override
std::string name() const override
double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override
void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override
virtual void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const =0
Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const final override
virtual Eigen::VectorXd compute_adjoint_rhs_step_prev(const int time_step, const Eigen::VectorXd &x, const State &state) const
virtual void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x)
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override
virtual std::string name() const override
virtual Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const =0
double value_unweighted(const Eigen::VectorXd &x) const final override
Compute the value of the form.
virtual ~StaticForm()=default
virtual double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const =0
virtual bool depends_on_step_prev() const final
virtual void solution_changed(const Eigen::VectorXd &new_x) final override
Update cached fields upon a change in the solution.
A collection of VariableToSimulation.
nlohmann::json json
Definition Common.hpp:9
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:22