PolyFEM
Loading...
Searching...
No Matches
AdjointForm.cpp
Go to the documentation of this file.
1#include "AdjointForm.hpp"
5#include <polyfem/State.hpp>
7
8namespace polyfem::solver
9{
10 double AdjointForm::value(const Eigen::VectorXd &x) const
11 {
12 double val = Form::value(x);
14 {
15 logger().debug("[{}] {}", print_energy_keyword_, val);
17 }
18 return val;
19 }
20
21 void AdjointForm::enable_energy_print(const std::string &print_energy_keyword)
22 {
23 print_energy_keyword_ = print_energy_keyword;
25 }
26
32
33 void AdjointForm::second_derivative_unweighted(const Eigen::VectorXd &x, StiffnessMatrix &hessian) const
34 {
35 log_and_throw_adjoint_error("[{}] Second derivatives not implemented", name());
36 }
37
38 Eigen::MatrixXd AdjointForm::compute_reduced_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const
39 {
40 Eigen::MatrixXd rhs = compute_adjoint_rhs(x, state);
41 if (!state.problem->is_time_dependent() && !state.lin_solver_cached) // nonlinear static solve only
42 {
43 Eigen::MatrixXd reduced;
44 for (int i = 0; i < rhs.cols(); i++)
45 {
46 Eigen::VectorXd reduced_vec = state.solve_data.nl_problem->full_to_reduced_grad(rhs.col(i));
47 if (i == 0)
48 reduced.setZero(reduced_vec.rows(), rhs.cols());
49 reduced.col(i) = reduced_vec;
50 }
51 return reduced;
52 }
53 else
54 return rhs;
55 }
56
57 void AdjointForm::first_derivative(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
58 {
59 Eigen::VectorXd partial_grad;
60 compute_partial_gradient(x, partial_grad);
61 gradv = variable_to_simulations_.compute_adjoint_term(x) + partial_grad;
62 }
63
64 void AdjointForm::first_derivative_unweighted(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
65 {
66 log_and_throw_adjoint_error("first_derivative_unweighted cannot be defined for adjoint forms!");
67 }
68
69 void AdjointForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
70 {
71 gradv = Eigen::VectorXd::Zero(x.size());
72 }
73
74 Eigen::MatrixXd AdjointForm::compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const
75 {
76 return Eigen::MatrixXd::Zero(state.ndof(), state.diff_cached.size());
77 }
78
79 void AdjointForm::update_quantities(const double t, const Eigen::VectorXd &x)
80 {
81
82 }
83 void AdjointForm::init_lagging(const Eigen::VectorXd &x)
84 {
85
86 }
87 void AdjointForm::update_lagging(const Eigen::VectorXd &x, const int iter_num)
88 {
89
90 }
91 void AdjointForm::set_apply_DBC(const Eigen::VectorXd &x, bool apply_DBC)
92 {
93
94 }
95
96 double StaticForm::value_unweighted(const Eigen::VectorXd &x) const
97 {
98 return value_unweighted_step(0, x);
99 }
100
101 void StaticForm::compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
102 {
104 }
105
106 Eigen::VectorXd StaticForm::compute_adjoint_rhs_step_prev(const int time_step, const Eigen::VectorXd &x, const State &state) const
107 {
108 return Eigen::MatrixXd::Zero(state.ndof(), 1);
109 }
110
111 Eigen::MatrixXd StaticForm::compute_adjoint_rhs(const Eigen::VectorXd &x, const State &state) const
112 {
113 assert(!depends_on_step_prev());
114 Eigen::MatrixXd term = Eigen::MatrixXd::Zero(state.ndof(), state.diff_cached.size());
115 term.col(0) = compute_adjoint_rhs_step(0, x, state);
116
117 return term;
118 }
119
120 double MaxStressForm::value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const
121 {
122 const double t = state_.problem->is_time_dependent() ? time_step * state_.args["time"]["dt"].get<double>() + state_.args["time"]["t0"].get<double>() : 0;
123 Eigen::VectorXd max_stress;
124 max_stress.setZero(state_.bases.size());
125 utils::maybe_parallel_for(state_.bases.size(), [&](int start, int end, int thread_id) {
126 Eigen::MatrixXd local_vals;
127 assembler::ElementAssemblyValues vals;
128 for (int e = start; e < end; e++)
129 {
130 if (interested_ids_.size() != 0 && interested_ids_.find(state_.mesh->get_body_id(e)) == interested_ids_.end())
131 continue;
132
133 state_.ass_vals_cache.compute(e, state_.mesh->is_volume(), state_.bases[e], state_.geom_bases()[e], vals);
134 // std::vector<assembler::Assembler::NamedMatrix> result;
135 // state_.assembler->compute_tensor_value(e, state_.bases[e], state_.geom_bases()[e], vals.quadrature.points, state_.diff_cached.u(time_step), result);
136 std::dynamic_pointer_cast<assembler::ElasticityAssembler>(state_.assembler)->compute_stress_tensor(assembler::OutputData(t, e, state_.bases[e], state_.geom_bases()[e], vals.quadrature.points, state_.diff_cached.u(time_step)), ElasticityTensorType::PK1, local_vals);
137
138 Eigen::VectorXd stress_norms = local_vals.rowwise().norm();
139 max_stress(e) = std::max(max_stress(e), stress_norms.maxCoeff());
140 }
141 });
142
143 return max_stress.maxCoeff();
144 }
145 Eigen::VectorXd MaxStressForm::compute_adjoint_rhs_step(const int time_step, const Eigen::VectorXd &x, const State &state) const
146 {
147 log_and_throw_adjoint_error("[{}] Not differentiable!", name());
148 return Eigen::VectorXd();
149 }
150 void MaxStressForm::compute_partial_gradient_step(const int time_step, const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const
151 {
152 log_and_throw_adjoint_error("[{}] Not differentiable!", name());
153 }
154} // namespace polyfem::solver
double val
Definition Assembler.cpp:86
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
std::shared_ptr< assembler::Problem > problem
current problem, it contains rhs and bc
Definition State.hpp:168
json args
main input arguments containing all defaults
Definition State.hpp:101
solver::DiffCache diff_cached
Definition State.hpp:669
std::vector< basis::ElementBases > bases
FE bases, the size is #elements.
Definition State.hpp:171
int ndof() const
Definition State.hpp:673
solver::SolveData solve_data
timedependent stuff cached
Definition State.hpp:324
std::unique_ptr< polysolve::linear::Solver > lin_solver_cached
Definition State.hpp:671
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 set_apply_DBC(const Eigen::VectorXd &x, bool apply_DBC) final override
Set if the Dirichlet boundary conditions should be enforced.
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)
virtual std::string name() const override
const VariableToSimulationGroup variable_to_simulations_
virtual double value(const Eigen::VectorXd &x) const
Compute the value of the form multiplied with the weigth.
Definition Form.hpp:24
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
std::shared_ptr< solver::NLProblem > nl_problem
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
void compute_partial_gradient(const Eigen::VectorXd &x, Eigen::VectorXd &gradv) const final 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 double value_unweighted_step(const int time_step, const Eigen::VectorXd &x) const =0
virtual bool depends_on_step_prev() const final
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const
Computes the sum of adjoint terms for all VariableToSimulation.
void maybe_parallel_for(int size, const std::function< void(int, int, int)> &partial_for)
spdlog::logger & logger()
Retrieves the current logger.
Definition Logger.cpp:42
void log_and_throw_adjoint_error(const std::string &msg)
Definition Logger.cpp:77
Eigen::SparseMatrix< double, Eigen::ColMajor > StiffnessMatrix
Definition Types.hpp:20