PolyFEM
Loading...
Searching...
No Matches
AdjointTools.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
5
6#include <Eigen/Core>
7
8#include <vector>
9#include <set>
10
11namespace polyfem
12{
13 class State;
14 class IntegrableFunctional;
15
16 namespace solver
17 {
18 class PeriodicMeshToMesh;
19 }
20} // namespace polyfem
21
22namespace polyfem::solver
23{
36
38 {
39 Volume,
40 Surface,
42 };
43
44 namespace AdjointTools
45 {
47 const State &state,
48 const IntegrableFunctional &j,
49 const Eigen::MatrixXd &solution,
50 const std::set<int> &interested_ids, // either body id or surface id
51 const SpatialIntegralType spatial_integral_type,
52 const int cur_step = 0);
53 void dJ_du_step(
54 const State &state,
55 const IntegrableFunctional &j,
56 const Eigen::MatrixXd &solution,
57 const std::set<int> &interested_ids,
58 const SpatialIntegralType spatial_integral_type,
59 const int cur_step,
60 Eigen::VectorXd &term);
62 const State &state,
63 const Eigen::MatrixXd &solution,
64 const IntegrableFunctional &j,
65 const std::set<int> &interested_ids, // either body id or surface id
66 const SpatialIntegralType spatial_integral_type,
67 Eigen::VectorXd &term,
68 const int cur_time_step);
70 const State &state,
71 const DiffCache &diff_cache,
72 const Eigen::MatrixXd &sol,
73 const Eigen::MatrixXd &adjoint,
74 Eigen::VectorXd &one_form);
75
76 // The periodic simulation result may not be differentiable wrt.
77 // every vertex of the periodic mesh, since a pair of periodic
78 // vertices cannot move freely, in which case
79 // dJ_periodic_shape_adjoint_term() instead of
80 // dJ_shape_homogenization_adjoint_term() should be used: When it
81 // computes shape derivatives, it considers the pair of periodic
82 // vertices as only one degree of freedom.
84 const State &state,
85 const DiffCache &diff_cache,
86 const Eigen::MatrixXd &sol,
87 const Eigen::MatrixXd &adjoint,
88 Eigen::VectorXd &one_form);
90 const State &state,
91 const DiffCache &diff_cache,
92 const PeriodicMeshToMesh &periodic_mesh_map,
93 const Eigen::VectorXd &periodic_mesh_representation,
94 const Eigen::MatrixXd &sol,
95 const Eigen::MatrixXd &adjoint,
96 Eigen::VectorXd &one_form);
97
99 const State &state,
100 const DiffCache &diff_cache,
101 const Eigen::MatrixXd &adjoint_nu,
102 const Eigen::MatrixXd &adjoint_p,
103 Eigen::VectorXd &one_form);
105 const State &state,
106 const Eigen::MatrixXd &sol,
107 const Eigen::MatrixXd &adjoint,
108 Eigen::VectorXd &one_form);
110 const State &state,
111 const DiffCache &diff_cache,
112 const Eigen::MatrixXd &adjoint_nu,
113 const Eigen::MatrixXd &adjoint_p,
114 Eigen::VectorXd &one_form);
116 const State &state,
117 const DiffCache &diff_cache,
118 const Eigen::MatrixXd &adjoint_nu,
119 const Eigen::MatrixXd &adjoint_p,
120 Eigen::VectorXd &one_form);
122 const State &state,
123 const DiffCache &diff_cache,
124 const Eigen::MatrixXd &adjoint_nu,
125 const Eigen::MatrixXd &adjoint_p,
126 Eigen::VectorXd &one_form);
128 const State &state,
129 const Eigen::MatrixXd &adjoint_nu,
130 const Eigen::MatrixXd &adjoint_p,
131 Eigen::VectorXd &one_form);
133 const State &state,
134 const DiffCache &diff_cache,
135 const Eigen::MatrixXd &adjoint,
136 Eigen::VectorXd &one_form);
138 const State &state,
139 const Eigen::MatrixXd &adjoint_nu,
140 const Eigen::MatrixXd &adjoint_p,
141 Eigen::VectorXd &one_form);
143 const State &state,
144 const std::vector<int> &boundary_ids,
145 const Eigen::MatrixXd &sol,
146 const Eigen::MatrixXd &adjoint,
147 Eigen::VectorXd &one_form);
149 const State &state,
150 const DiffCache &diff_cache,
151 const std::vector<int> &boundary_ids,
152 const Eigen::MatrixXd &adjoint_nu,
153 const Eigen::MatrixXd &adjoint_p,
154 Eigen::VectorXd &one_form);
155
156 Eigen::VectorXd map_primitive_to_node_order(
157 const State &state,
158 const Eigen::VectorXd &primitives);
159 Eigen::VectorXd map_node_to_primitive_order(
160 const State &state,
161 const Eigen::VectorXd &nodes);
162
163 Eigen::MatrixXd edge_normal_gradient(
164 const Eigen::MatrixXd &V);
165 Eigen::MatrixXd face_normal_gradient(
166 const Eigen::MatrixXd &V);
167
168 Eigen::MatrixXd edge_velocity_divergence(
169 const Eigen::MatrixXd &V);
170 Eigen::MatrixXd face_velocity_divergence(
171 const Eigen::MatrixXd &V);
172
173 void scaled_jacobian(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, Eigen::VectorXd &quality);
174 }; // namespace AdjointTools
175} // namespace polyfem::solver
int V
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
void dJ_pressure_static_adjoint_term(const State &state, const std::vector< int > &boundary_ids, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
void dJ_friction_transient_adjoint_term(const State &state, const DiffCache &diff_cache, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
Eigen::MatrixXd face_velocity_divergence(const Eigen::MatrixXd &V)
void dJ_shape_homogenization_adjoint_term(const State &state, const DiffCache &diff_cache, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
void scaled_jacobian(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, Eigen::VectorXd &quality)
double integrate_objective(const State &state, const IntegrableFunctional &j, const Eigen::MatrixXd &solution, const std::set< int > &interested_ids, const SpatialIntegralType spatial_integral_type, const int cur_step=0)
Eigen::MatrixXd face_normal_gradient(const Eigen::MatrixXd &V)
void dJ_shape_static_adjoint_term(const State &state, const DiffCache &diff_cache, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
void dJ_initial_condition_adjoint_term(const State &state, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
Eigen::MatrixXd edge_velocity_divergence(const Eigen::MatrixXd &V)
Eigen::MatrixXd edge_normal_gradient(const Eigen::MatrixXd &V)
void compute_shape_derivative_functional_term(const State &state, const Eigen::MatrixXd &solution, const IntegrableFunctional &j, const std::set< int > &interested_ids, const SpatialIntegralType spatial_integral_type, Eigen::VectorXd &term, const int cur_time_step)
Eigen::VectorXd map_primitive_to_node_order(const State &state, const Eigen::VectorXd &primitives)
void dJ_material_transient_adjoint_term(const State &state, const DiffCache &diff_cache, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
void dJ_dirichlet_static_adjoint_term(const State &state, const DiffCache &diff_cache, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
void dJ_shape_transient_adjoint_term(const State &state, const DiffCache &diff_cache, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
void dJ_material_static_adjoint_term(const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
void dJ_periodic_shape_adjoint_term(const State &state, const DiffCache &diff_cache, const PeriodicMeshToMesh &periodic_mesh_map, const Eigen::VectorXd &periodic_mesh_representation, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
void dJ_dirichlet_transient_adjoint_term(const State &state, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
void dJ_pressure_transient_adjoint_term(const State &state, const DiffCache &diff_cache, const std::vector< int > &boundary_ids, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
void dJ_damping_transient_adjoint_term(const State &state, const DiffCache &diff_cache, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
Eigen::VectorXd map_node_to_primitive_order(const State &state, const Eigen::VectorXd &nodes)
void dJ_du_step(const State &state, const IntegrableFunctional &j, const Eigen::MatrixXd &solution, const std::set< int > &interested_ids, const SpatialIntegralType spatial_integral_type, const int cur_step, Eigen::VectorXd &term)