PolyFEM
Loading...
Searching...
No Matches
AdjointForm.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
8
9#include <memory>
10#include <utility>
11#include <set>
12#include <string>
13#include <vector>
14
15namespace polyfem::solver
16{
17 class AdjointForm : public Form
18 {
19 public:
20 AdjointForm(const VariableToSimulationGroup &variable_to_simulations) : variable_to_simulations_(variable_to_simulations) {}
21 virtual ~AdjointForm() {}
22
23 virtual std::string name() const override { return "adjoint"; }
24 void enable_energy_print(const std::string &print_energy_keyword);
25
26 double value(const Eigen::VectorXd &x) const override;
27 virtual void solution_changed(const Eigen::VectorXd &new_x) override;
28
30
31 virtual Eigen::MatrixXd compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const;
32 virtual void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const;
33 virtual void first_derivative(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override;
34 virtual Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const;
35
36 // not used functions from base class
37 virtual void update_quantities(const double t, const Eigen::VectorXd &x) final override;
38 virtual void init_lagging(const Eigen::VectorXd &x) final override;
39 virtual void update_lagging(const Eigen::VectorXd &x, const int iter_num) final override;
40
41 protected:
42 virtual void first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override;
43 virtual void second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const final override;
44
46
47 enum class PrintStage
48 {
52 };
53
56 };
57
58 class StaticForm : public AdjointForm
59 {
60 public:
62 virtual ~StaticForm() = default;
63
64 virtual std::string name() const override { return "static"; }
65
66 Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const final override;
67 void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final override;
68 double value_unweighted(const Eigen::VectorXd &x) const final override;
69
70 virtual Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const = 0;
71 virtual Eigen::VectorXd compute_adjoint_rhs_step_prev(const int time_step, const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const;
72 virtual void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const = 0;
73 virtual double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const = 0;
74 virtual void solution_changed(const Eigen::VectorXd &new_x) final override;
75 virtual void solution_changed_step(const int time_step, const Eigen::VectorXd &new_x) {}
76 virtual bool depends_on_step_prev() const final { return depends_on_step_prev_; }
77
78 protected:
80 };
81
83 {
84 public:
85 MaxStressForm(const VariableToSimulationGroup &variable_to_simulations, std::shared_ptr<const State> state, std::shared_ptr<const DiffCache> diff_cache, const json &args)
86 : StaticForm(variable_to_simulations), state_(std::move(state)), diff_cache_(std::move(diff_cache))
87 {
88 auto tmp_ids = args["volume_selection"].get<std::vector<int>>();
89 interested_ids_ = std::set(tmp_ids.begin(), tmp_ids.end());
90 }
91
92 std::string name() const override { return "max_stress"; }
93
94 Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const override;
95 double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const override;
96 void compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const override;
97
98 private:
99 std::set<int> interested_ids_;
100 std::shared_ptr<const State> state_;
101 std::shared_ptr<const DiffCache> diff_cache_;
102 };
103} // namespace polyfem::solver
int x
Storage for additional data required by differntial code.
Definition DiffCache.hpp:21
main class that contains the polyfem solver and all its state
Definition State.hpp:113
virtual Eigen::MatrixXd compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const
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 Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const
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.
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 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, std::shared_ptr< const State > state, std::shared_ptr< const DiffCache > diff_cache, const json &args)
Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const override
std::shared_ptr< const DiffCache > diff_cache_
std::shared_ptr< const State > state_
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
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
Eigen::MatrixXd compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const final override
virtual std::string name() const override
double value_unweighted(const Eigen::VectorXd &x) const final override
Compute the value of the form.
virtual Eigen::VectorXd compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const =0
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.
virtual Eigen::VectorXd compute_adjoint_rhs_step_prev(const int time_step, const Eigen::VectorXd &x, const State &state, const DiffCache &diff_cache) const
A collection of VariableToSimulation.
nlohmann::json json
Definition Common.hpp:9
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:24