Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
AdjointTools.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
4#include <set>
5
6namespace polyfem
7{
8 class State;
9 class IntegrableFunctional;
10
11 namespace solver
12 {
13 class PeriodicMeshToMesh;
14 }
15} // namespace polyfem
16
17namespace polyfem::solver
18{
31
33 {
34 Volume,
35 Surface,
37 };
38
39 namespace AdjointTools
40 {
42 const State &state,
43 const IntegrableFunctional &j,
44 const Eigen::MatrixXd &solution,
45 const std::set<int> &interested_ids, // either body id or surface id
46 const SpatialIntegralType spatial_integral_type,
47 const int cur_step = 0);
48 void dJ_du_step(
49 const State &state,
50 const IntegrableFunctional &j,
51 const Eigen::MatrixXd &solution,
52 const std::set<int> &interested_ids,
53 const SpatialIntegralType spatial_integral_type,
54 const int cur_step,
55 Eigen::VectorXd &term);
57 const State &state,
58 const Eigen::MatrixXd &solution,
59 const IntegrableFunctional &j,
60 const std::set<int> &interested_ids, // either body id or surface id
61 const SpatialIntegralType spatial_integral_type,
62 Eigen::VectorXd &term,
63 const int cur_time_step);
65 const State &state,
66 const Eigen::MatrixXd &sol,
67 const Eigen::MatrixXd &adjoint,
68 Eigen::VectorXd &one_form);
69
70 // The periodic simulation result may not be differentiable wrt.
71 // every vertex of the periodic mesh, since a pair of periodic
72 // vertices cannot move freely, in which case
73 // dJ_periodic_shape_adjoint_term() instead of
74 // dJ_shape_homogenization_adjoint_term() should be used: When it
75 // computes shape derivatives, it considers the pair of periodic
76 // vertices as only one degree of freedom.
78 const State &state,
79 const Eigen::MatrixXd &sol,
80 const Eigen::MatrixXd &adjoint,
81 Eigen::VectorXd &one_form);
83 const State &state,
84 const PeriodicMeshToMesh &periodic_mesh_map,
85 const Eigen::VectorXd &periodic_mesh_representation,
86 const Eigen::MatrixXd &sol,
87 const Eigen::MatrixXd &adjoint,
88 Eigen::VectorXd &one_form);
89
91 const State &state,
92 const Eigen::MatrixXd &adjoint_nu,
93 const Eigen::MatrixXd &adjoint_p,
94 Eigen::VectorXd &one_form);
96 const State &state,
97 const Eigen::MatrixXd &sol,
98 const Eigen::MatrixXd &adjoint,
99 Eigen::VectorXd &one_form);
101 const State &state,
102 const Eigen::MatrixXd &adjoint_nu,
103 const Eigen::MatrixXd &adjoint_p,
104 Eigen::VectorXd &one_form);
106 const State &state,
107 const Eigen::MatrixXd &adjoint_nu,
108 const Eigen::MatrixXd &adjoint_p,
109 Eigen::VectorXd &one_form);
111 const State &state,
112 const Eigen::MatrixXd &adjoint_nu,
113 const Eigen::MatrixXd &adjoint_p,
114 Eigen::VectorXd &one_form);
116 const State &state,
117 const Eigen::MatrixXd &adjoint_nu,
118 const Eigen::MatrixXd &adjoint_p,
119 Eigen::VectorXd &one_form);
121 const State &state,
122 const Eigen::MatrixXd &adjoint,
123 Eigen::VectorXd &one_form);
125 const State &state,
126 const Eigen::MatrixXd &adjoint_nu,
127 const Eigen::MatrixXd &adjoint_p,
128 Eigen::VectorXd &one_form);
130 const State &state,
131 const std::vector<int> &boundary_ids,
132 const Eigen::MatrixXd &sol,
133 const Eigen::MatrixXd &adjoint,
134 Eigen::VectorXd &one_form);
136 const State &state,
137 const std::vector<int> &boundary_ids,
138 const Eigen::MatrixXd &adjoint_nu,
139 const Eigen::MatrixXd &adjoint_p,
140 Eigen::VectorXd &one_form);
141
142 Eigen::VectorXd map_primitive_to_node_order(
143 const State &state,
144 const Eigen::VectorXd &primitives);
145 Eigen::VectorXd map_node_to_primitive_order(
146 const State &state,
147 const Eigen::VectorXd &nodes);
148
149 Eigen::MatrixXd edge_normal_gradient(
150 const Eigen::MatrixXd &V);
151 Eigen::MatrixXd face_normal_gradient(
152 const Eigen::MatrixXd &V);
153
154 Eigen::MatrixXd edge_velocity_divergence(
155 const Eigen::MatrixXd &V);
156 Eigen::MatrixXd face_velocity_divergence(
157 const Eigen::MatrixXd &V);
158
159 double triangle_jacobian(const Eigen::VectorXd &v1, const Eigen::VectorXd &v2, const Eigen::VectorXd &v3);
160 double tet_determinant(const Eigen::VectorXd &v1, const Eigen::VectorXd &v2, const Eigen::VectorXd &v3, const Eigen::VectorXd &v4);
161 void scaled_jacobian(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, Eigen::VectorXd &quality);
162 bool is_flipped(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F);
163 }; // namespace AdjointTools
164} // namespace polyfem::solver
int V
main class that contains the polyfem solver and all its state
Definition State.hpp:79
double tet_determinant(const Eigen::VectorXd &v1, const Eigen::VectorXd &v2, const Eigen::VectorXd &v3, const Eigen::VectorXd &v4)
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_periodic_shape_adjoint_term(const State &state, const PeriodicMeshToMesh &periodic_mesh_map, const Eigen::VectorXd &periodic_mesh_representation, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
Eigen::MatrixXd face_velocity_divergence(const Eigen::MatrixXd &V)
void scaled_jacobian(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F, Eigen::VectorXd &quality)
void dJ_dirichlet_static_adjoint_term(const State &state, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
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_friction_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 std::vector< int > &boundary_ids, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, 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)
void dJ_material_transient_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_shape_homogenization_adjoint_term(const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
void dJ_shape_transient_adjoint_term(const State &state, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
double triangle_jacobian(const Eigen::VectorXd &v1, const Eigen::VectorXd &v2, const Eigen::VectorXd &v3)
void dJ_material_static_adjoint_term(const State &state, 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_shape_static_adjoint_term(const State &state, const Eigen::MatrixXd &sol, const Eigen::MatrixXd &adjoint, Eigen::VectorXd &one_form)
void dJ_damping_transient_adjoint_term(const State &state, const Eigen::MatrixXd &adjoint_nu, const Eigen::MatrixXd &adjoint_p, Eigen::VectorXd &one_form)
bool is_flipped(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)
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)