Loading [MathJax]/extensions/tex2jax.js
PolyFEM
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
VariableToSimulation.hpp
Go to the documentation of this file.
1#pragma once
2
6
7#include <iostream>
8
9namespace polyfem::solver
10{
13 {
14 public:
15 VariableToSimulation(const std::vector<std::shared_ptr<State>> &states, const CompositeParametrization &parametrization) : states_(states), parametrization_(parametrization) { assert(states.size() > 0); }
16 VariableToSimulation(const std::shared_ptr<State> &state, const CompositeParametrization &parametrization) : states_({state}), parametrization_(parametrization) {}
18
19 static std::unique_ptr<VariableToSimulation> create(const std::string &type, const std::vector<std::shared_ptr<State>> &states, CompositeParametrization &&parametrization);
20
21 inline virtual void update(const Eigen::VectorXd &x)
22 {
24 }
25
26 virtual std::string name() const = 0;
27
28 inline int n_states() const { return states_.size(); }
29 inline const std::vector<std::shared_ptr<State>> &get_states() const { return states_; }
30
32 virtual ParameterType get_parameter_type() const = 0;
33 virtual Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const = 0;
34 virtual Eigen::VectorXd inverse_eval();
35
36 virtual void set_output_indexing(const json &args);
37 Eigen::VectorXi get_output_indexing(const Eigen::VectorXd &x) const;
38
39 virtual Eigen::VectorXd apply_parametrization_jacobian(const Eigen::VectorXd &term, const Eigen::VectorXd &x) const;
40
41 protected:
42 virtual void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices);
43 const std::vector<std::shared_ptr<State>> states_;
45
46 Eigen::VectorXi output_indexing_; // if a derived class overrides apply_parametrization_jacobian(term, x), this is not necessarily used.
47 };
48
51 {
52 public:
53 using ValueType = std::shared_ptr<VariableToSimulation>;
54
56 virtual ~VariableToSimulationGroup() = default;
57
58 void init(const json &args, const std::vector<std::shared_ptr<State>> &states, const std::vector<int> &variable_sizes);
59
62 inline void update(const Eigen::VectorXd &x)
63 {
64 for (auto &v2s : L)
65 v2s->update(x);
66 }
67
71 void compute_state_variable(const ParameterType type, const State *state_ptr, const Eigen::VectorXd &x, Eigen::VectorXd &state_variable) const;
72
76 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const;
77
84 virtual Eigen::VectorXd apply_parametrization_jacobian(const ParameterType type, const State *state_ptr, const Eigen::VectorXd &x, const std::function<Eigen::VectorXd()> &grad) const;
85
86 typedef std::vector<ValueType>::const_iterator const_iterator;
87
88 inline ValueType &operator[](size_t i) { return L[i]; }
89 inline const ValueType &operator[](size_t i) const { return L[i]; }
90 inline const_iterator begin() const { return L.begin(); }
91 inline const_iterator end() const { return L.end(); }
92 inline void push_back(const ValueType &v2s) { L.push_back(v2s); }
93 inline void clear() { L.clear(); }
94
95 private:
96 std::vector<ValueType> L;
97 };
98
99 // state variable dof = dim * n_vertices
101 {
102 public:
105
106 std::string name() const override { return "shape"; }
107
109
110 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
111 Eigen::VectorXd inverse_eval() override;
112
113 void set_output_indexing(const json &args) override;
114
115 protected:
116 virtual void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override;
117 };
118
119 // To optimize per element elastic parameters
120 // state variable dof = 2 * n_elements
122 {
123 public:
126
127 std::string name() const override { return "elastic"; }
128
130
131 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
132 Eigen::VectorXd inverse_eval() override;
133
134 protected:
135 void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override;
136 };
137
138 // To optimize the friction constant
139 // state variable dof = 1
141 {
142 public:
145
146 std::string name() const override { return "friction"; }
147
149
150 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
151 Eigen::VectorXd inverse_eval() override;
152
153 protected:
154 void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override;
155 };
156
157 // To optimize the damping constant psi and phi
158 // state variable dof = 2
160 {
161 public:
164
165 std::string name() const override { return "damping"; }
166
168
169 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
170 Eigen::VectorXd inverse_eval() override;
171
172 protected:
173 void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override;
174 };
175
176 // To optimize the per node initial condition
177 // state variable dof = dim * n_bases
179 {
180 public:
183
184 std::string name() const override { return "initial"; }
185
187
188 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
189 Eigen::VectorXd inverse_eval() override;
190
191 protected:
192 void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override;
193 };
194
195 // To optimize the per node Dirichlet values, only work in transient simulations
196 // state variable dof = dim * n_time_steps * n_dirichlet_nodes
198 {
199 public:
202
203 std::string name() const override { return "dirichlet"; }
204
205 void set_dirichlet_boundaries(const std::vector<int> &dirichlet_boundaries)
206 {
207 dirichlet_boundaries_ = dirichlet_boundaries;
208 }
209
211
212 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
213 Eigen::VectorXd inverse_eval() override;
214
215 void set_output_indexing(const json &args) override;
216
217 protected:
218 void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override;
219
220 private:
221 std::string variable_to_string(const Eigen::VectorXd &variable);
222
223 std::vector<int> dirichlet_boundaries_;
224 };
225
227 {
228 public:
231
232 std::string name() const override { return "dirichlet-nodes"; }
233
234 void set_dirichlet_nodes(const Eigen::VectorXi &dirichlet_nodes)
235 {
236 dirichlet_nodes_ = dirichlet_nodes;
237 }
238
240
241 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
242 Eigen::VectorXd inverse_eval() override;
243
244 void set_output_indexing(const json &args) override;
245
246 protected:
247 void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override;
248
249 private:
250 std::string variable_to_string(const Eigen::VectorXd &variable);
251
252 Eigen::VectorXi dirichlet_nodes_;
253 };
254
255 // To optimize the per node pressure boundaries
256 // Each pressure boundary will have the same value
257 // state variable dof = dim * n_time_steps
259 {
260 public:
263
264 std::string name() const override { return "pressure"; }
265
266 void set_pressure_boundaries(const std::vector<int> &pressure_boundaries)
267 {
268 pressure_boundaries_ = pressure_boundaries;
269 }
270
272
273 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
274 Eigen::VectorXd inverse_eval() override;
275
276 void set_output_indexing(const json &args) override;
277
278 protected:
279 void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override;
280
281 private:
282 std::string variable_to_string(const Eigen::VectorXd &variable);
283
284 std::vector<int> pressure_boundaries_;
285 };
286
287 // state variable dof = dim * n_vertices - periodic dof
289 {
290 public:
293
294 std::string name() const override { return "periodic-shape"; }
295
297
298 Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override;
299 Eigen::VectorXd inverse_eval() override;
300
301 void update(const Eigen::VectorXd &x) override;
302
303 Eigen::VectorXd apply_parametrization_jacobian(const Eigen::VectorXd &term, const Eigen::VectorXd &x) const override;
304
305 protected:
306 std::unique_ptr<PeriodicMeshToMesh> periodic_mesh_map;
308 };
309} // namespace polyfem::solver
int x
main class that contains the polyfem solver and all its state
Definition State.hpp:79
Eigen::VectorXd eval(const Eigen::VectorXd &x) const override
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override
void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override
void set_dirichlet_nodes(const Eigen::VectorXi &dirichlet_nodes)
std::string variable_to_string(const Eigen::VectorXd &variable)
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
std::string variable_to_string(const Eigen::VectorXd &variable)
void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override
void set_dirichlet_boundaries(const std::vector< int > &dirichlet_boundaries)
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
Eigen::VectorXd apply_parametrization_jacobian(const Eigen::VectorXd &term, const Eigen::VectorXd &x) const override
std::unique_ptr< PeriodicMeshToMesh > periodic_mesh_map
void set_pressure_boundaries(const std::vector< int > &pressure_boundaries)
std::string variable_to_string(const Eigen::VectorXd &variable)
void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
virtual void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices) override
void set_output_indexing(const json &args) override
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const override
ParameterType get_parameter_type() const override
A collection of VariableToSimulation.
std::shared_ptr< VariableToSimulation > ValueType
void compute_state_variable(const ParameterType type, const State *state_ptr, const Eigen::VectorXd &x, Eigen::VectorXd &state_variable) const
Evaluate the variable to simulations and overwrite the state_variable based on x.
Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const
Computes the sum of adjoint terms for all VariableToSimulation.
virtual Eigen::VectorXd apply_parametrization_jacobian(const ParameterType type, const State *state_ptr, const Eigen::VectorXd &x, const std::function< Eigen::VectorXd()> &grad) const
Maps the partial gradient wrt.
const ValueType & operator[](size_t i) const
void init(const json &args, const std::vector< std::shared_ptr< State > > &states, const std::vector< int > &variable_sizes)
void update(const Eigen::VectorXd &x)
Update parameters in simulators.
std::vector< ValueType >::const_iterator const_iterator
Maps the optimization variable to the state variable.
virtual Eigen::VectorXd apply_parametrization_jacobian(const Eigen::VectorXd &term, const Eigen::VectorXd &x) const
VariableToSimulation(const std::vector< std::shared_ptr< State > > &states, const CompositeParametrization &parametrization)
virtual ParameterType get_parameter_type() const =0
VariableToSimulation(const std::shared_ptr< State > &state, const CompositeParametrization &parametrization)
virtual Eigen::VectorXd compute_adjoint_term(const Eigen::VectorXd &x) const =0
const std::vector< std::shared_ptr< State > > states_
virtual void update_state(const Eigen::VectorXd &state_variable, const Eigen::VectorXi &indices)
virtual void set_output_indexing(const json &args)
virtual void update(const Eigen::VectorXd &x)
CompositeParametrization & get_parametrization()
virtual std::string name() const =0
Eigen::VectorXi get_output_indexing(const Eigen::VectorXd &x) const
static std::unique_ptr< VariableToSimulation > create(const std::string &type, const std::vector< std::shared_ptr< State > > &states, CompositeParametrization &&parametrization)
const std::vector< std::shared_ptr< State > > & get_states() const
nlohmann::json json
Definition Common.hpp:9