[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
rol_objective.cpp
1 #include "rol_to_dealii_vector.hpp"
2 #include "rol_objective.hpp"
3 
4 #include <deal.II/optimization/rol/vector_adaptor.h>
5 
6 #include "global_counter.hpp"
7 
8 namespace PHiLiP {
9 
10 template <int dim, int nstate>
12  Functional<dim,nstate,double> &_functional,
13  std::shared_ptr<BaseParameterization<dim>> _design_parameterization,
14  std::shared_ptr<dealii::TrilinosWrappers::SparseMatrix> precomputed_dXvdXp)
15  : functional(_functional)
16  , design_parameterization(_design_parameterization)
17 {
18  Assert(functional.dg->high_order_grid == design_parameterization->high_order_grid,
19  dealii::ExcMessage("Functional and DesignParameterization do not point to the same high order grid."));
20  design_parameterization->initialize_design_variables(design_var);
21  const unsigned int n_design_variables = design_parameterization->get_number_of_design_variables();
22 
23  if (precomputed_dXvdXp) {
24  if (precomputed_dXvdXp->m() == functional.dg->high_order_grid->volume_nodes.size() && precomputed_dXvdXp->n() == n_design_variables) {
25  dXvdXp.copy_from(*precomputed_dXvdXp);
26  }
27  } else {
28  design_parameterization->compute_dXv_dXp(dXvdXp);
29  }
30 }
31 
32 
33 template <int dim, int nstate>
35  const ROL::Vector<double> &des_var_sim,
36  const ROL::Vector<double> &des_var_ctl,
37  bool /*flag*/, int /*iter*/)
38 {
40 
42  design_parameterization->update_mesh_from_design_variables(dXvdXp, design_var);
43 }
44 
45 
46 template <int dim, int nstate>
48  const ROL::Vector<double> &des_var_sim,
49  const ROL::Vector<double> &des_var_ctl,
50  double &tol )
51 {
52  // Tolerance tends to not be used except in the case of a reduced objective function.
53  // In that scenario, tol is the constraint norm.
54  // If the flow has not converged (>1e-5 or is nan), simply return a high functional.
55  // This is likely happening in the linesearch while optimizing in the reduced-space.
56  if (tol > 1e-5 || std::isnan(tol)) return 1e200;
57  update(des_var_sim, des_var_ctl);
58 
59  const bool compute_dIdW = false;
60  const bool compute_dIdX = false;
61  const bool compute_d2I = false;
62  return functional.evaluate_functional( compute_dIdW, compute_dIdX, compute_d2I );
63 }
64 
65 template <int dim, int nstate>
67  ROL::Vector<double> &gradient_sim,
68  const ROL::Vector<double> &des_var_sim,
69  const ROL::Vector<double> &des_var_ctl,
70  double &/*tol*/ )
71 {
72  update(des_var_sim, des_var_ctl);
73 
74  const bool compute_dIdW = true;
75  const bool compute_dIdX = false;
76  const bool compute_d2I = false;
77  functional.evaluate_functional( compute_dIdW, compute_dIdX, compute_d2I );
78  auto &dIdW = ROL_vector_to_dealii_vector_reference(gradient_sim);
79  dIdW = functional.dIdw;
80 }
81 
82 template <int dim, int nstate>
84  ROL::Vector<double> &gradient_ctl,
85  const ROL::Vector<double> &des_var_sim,
86  const ROL::Vector<double> &des_var_ctl,
87  double &/*tol*/ )
88 {
89  update(des_var_sim, des_var_ctl);
90 
91  const bool compute_dIdW = false, compute_dIdX = true, compute_d2I = false;
92  functional.evaluate_functional( compute_dIdW, compute_dIdX, compute_d2I );
93 
94  const auto &dIdXv = functional.dIdX;
95 
96  auto &dealii_output = ROL_vector_to_dealii_vector_reference(gradient_ctl);
97  dXvdXp.Tvmult(dealii_output, dIdXv);
98 
99  //n_vmult += 1;
100 
101  // auto dIdXvs = dIdXv;
102  // {
103  // dealii::LinearAlgebra::distributed::Vector<double> dummy_vector(functional.dg->high_order_grid->surface_nodes);
104  // MeshMover::LinearElasticity<dim, double, dealii::LinearAlgebra::distributed::Vector<double>, dealii::DoFHandler<dim>>
105  // meshmover(*(functional.dg->high_order_grid->triangulation),
106  // functional.dg->high_order_grid->initial_mapping_fe_field,
107  // functional.dg->high_order_grid->dof_handler_grid,
108  // functional.dg->high_order_grid->surface_to_volume_indices,
109  // dummy_vector);
110  // meshmover.apply_dXvdXvs_transpose(dIdXv, dIdXvs);
111  // }
112 
113  // auto &dIdXp = ROL_vector_to_dealii_vector_reference(gradient_ctl);
114  // {
115  // dealii::TrilinosWrappers::SparseMatrix dXvsdXp;
116  // ffd.get_dXvsdXp (*(functional.dg->high_order_grid), ffd_design_variables_indices_dim, dXvsdXp);
117  // dXvsdXp.Tvmult(dIdXp, dIdXvs);
118  // }
119 
120 }
121 
122 template <int dim, int nstate>
124  ROL::Vector<double> &output_vector,
125  const ROL::Vector<double> &input_vector,
126  const ROL::Vector<double> &des_var_sim,
127  const ROL::Vector<double> &des_var_ctl,
128  double &/*tol*/ )
129 {
130  update(des_var_sim, des_var_ctl);
131 
132  const bool compute_dIdW = false;
133  const bool compute_dIdX = false;
134  const bool compute_d2I = true;
135  functional.evaluate_functional( compute_dIdW, compute_dIdX, compute_d2I );
136 
137  const auto &dealii_input = ROL_vector_to_dealii_vector_reference(input_vector);
138  auto &hv = ROL_vector_to_dealii_vector_reference(output_vector);
139 
140  functional.d2IdWdW->vmult(hv, dealii_input);
141 
142  //n_vmult += 1;
143 }
144 
145 template <int dim, int nstate>
147  ROL::Vector<double> &output_vector,
148  const ROL::Vector<double> &input_vector,
149  const ROL::Vector<double> &des_var_sim,
150  const ROL::Vector<double> &des_var_ctl,
151  double &/*tol*/ )
152 {
153  update(des_var_sim, des_var_ctl);
154 
155  const auto &dealii_input = ROL_vector_to_dealii_vector_reference(input_vector);
156 
157  // auto dXvsdXp_input = functional.dg->high_order_grid->volume_nodes;
158  // {
159  // dealii::TrilinosWrappers::SparseMatrix dXvsdXp;
160  // ffd.get_dXvsdXp (functional.dg->high_order_grid, ffd_design_variables_indices_dim, dXvsdXp);
161  // dXvsdXp.vmult(dXvsdXp_input, dealii_input);
162  // }
163 
164  // auto dXvdXp_input = dXvsdXp_input;
165  // {
166  // dealii::LinearAlgebra::distributed::Vector<double> dummy_vector(functional.dg->high_order_grid->surface_nodes);
167  // MeshMover::LinearElasticity<dim, double, dealii::LinearAlgebra::distributed::Vector<double>, dealii::DoFHandler<dim>>
168  // meshmover(*(functional.dg->high_order_grid->triangulation),
169  // functional.dg->high_order_grid->initial_mapping_fe_field,
170  // functional.dg->high_order_grid->dof_handler_grid,
171  // functional.dg->high_order_grid->surface_to_volume_indices,
172  // dummy_vector);
173  // meshmover.apply_dXvdXvs(dXvsdXp_input, dXvdXp_input);
174  // }
175 
176  // auto &d2IdWdXp_input = ROL_vector_to_dealii_vector_reference(output_vector);
177  // {
178  // const bool compute_dIdW = false, compute_dIdX = false, compute_d2I = true;
179  // functional.evaluate_functional( compute_dIdW, compute_dIdX, compute_d2I );
180  // functional.d2IdWdX.vmult(d2IdWdXp_input, dXvdXp_input);
181  // }
182 
183  auto dXvdXp_input = functional.dg->high_order_grid->volume_nodes;
184  dXvdXp.vmult(dXvdXp_input, dealii_input);
185 
186  auto &dealii_output = ROL_vector_to_dealii_vector_reference(output_vector);
187  {
188  const bool compute_dIdW = false, compute_dIdX = false, compute_d2I = true;
189  functional.evaluate_functional( compute_dIdW, compute_dIdX, compute_d2I );
190  functional.d2IdWdX->vmult(dealii_output, dXvdXp_input);
191  }
192 
193  //n_vmult += 2;
194 
195 }
196 
197 template <int dim, int nstate>
199  ROL::Vector<double> &output_vector,
200  const ROL::Vector<double> &input_vector,
201  const ROL::Vector<double> &des_var_sim,
202  const ROL::Vector<double> &des_var_ctl,
203  double &/*tol*/ )
204 {
205  update(des_var_sim, des_var_ctl);
206 
207  const bool compute_dIdW = false;
208  const bool compute_dIdX = false;
209  const bool compute_d2I = true;
210  functional.evaluate_functional( compute_dIdW, compute_dIdX, compute_d2I );
211 
212  const auto &dealii_input = ROL_vector_to_dealii_vector_reference(input_vector);
213 
214  auto d2IdXdW_input = functional.dg->high_order_grid->volume_nodes;
215  functional.d2IdWdX->Tvmult(d2IdXdW_input, dealii_input);
216 
217  // auto d2IdXvsdW_input = functional.dg->high_order_grid->volume_nodes;
218  // {
219  // dealii::LinearAlgebra::distributed::Vector<double> dummy_vector(functional.dg->high_order_grid->surface_nodes);
220  // MeshMover::LinearElasticity<dim, double, dealii::LinearAlgebra::distributed::Vector<double>, dealii::DoFHandler<dim>>
221  // meshmover(*(functional.dg->high_order_grid->triangulation),
222  // functional.dg->high_order_grid->initial_mapping_fe_field,
223  // functional.dg->high_order_grid->dof_handler_grid,
224  // functional.dg->high_order_grid->surface_to_volume_indices,
225  // dummy_vector);
226  // meshmover.apply_dXvdXvs_transpose(d2IdXdW_input, d2IdXvsdW_input);
227  // }
228 
229  // auto &d2IdXpdW_input = ROL_vector_to_dealii_vector_reference(output_vector);
230  // {
231  // dealii::TrilinosWrappers::SparseMatrix dXvsdXp;
232  // ffd.get_dXvsdXp (*(functional.dg->high_order_grid), ffd_design_variables_indices_dim, dXvsdXp);
233  // dXvsdXp.Tvmult(d2IdXpdW_input, d2IdXvsdW_input);
234  // }
235 
236  auto &dealii_output = ROL_vector_to_dealii_vector_reference(output_vector);
237  dXvdXp.Tvmult(dealii_output, d2IdXdW_input);
238 
239  //n_vmult += 2;
240 }
241 
242 template <int dim, int nstate>
244  ROL::Vector<double> &output_vector,
245  const ROL::Vector<double> &input_vector,
246  const ROL::Vector<double> &des_var_sim,
247  const ROL::Vector<double> &des_var_ctl,
248  double &/*tol*/ )
249 {
250  update(des_var_sim, des_var_ctl);
251 
252 
253  const auto &dealii_input = ROL_vector_to_dealii_vector_reference(input_vector);
254 
255  // dealii::TrilinosWrappers::SparseMatrix dXvsdXp;
256  // ffd.get_dXvsdXp (*(functional.dg->high_order_grid), ffd_design_variables_indices_dim, dXvsdXp);
257 
258  // auto dXvsdXp_input = functional.dg->high_order_grid->volume_nodes;
259  // {
260  // dXvsdXp.vmult(dXvsdXp_input, dealii_input);
261  // }
262 
263  // auto dXvdXp_input = dXvsdXp_input;
264  // {
265  // dealii::LinearAlgebra::distributed::Vector<double> dummy_vector(functional.dg->high_order_grid->surface_nodes);
266  // MeshMover::LinearElasticity<dim, double, dealii::LinearAlgebra::distributed::Vector<double>, dealii::DoFHandler<dim>>
267  // meshmover(*(functional.dg->high_order_grid->triangulation),
268  // functional.dg->high_order_grid->initial_mapping_fe_field,
269  // functional.dg->high_order_grid->dof_handler_grid,
270  // functional.dg->high_order_grid->surface_to_volume_indices,
271  // dummy_vector);
272  // meshmover.apply_dXvdXvs(dXvsdXp_input, dXvdXp_input);
273  // }
274 
275  auto dXvdXp_input = functional.dg->high_order_grid->volume_nodes;
276  dXvdXp.vmult(dXvdXp_input, dealii_input);
277 
278  auto d2IdXdXp_input = functional.dg->high_order_grid->volume_nodes;
279  {
280  const bool compute_dIdW = false, compute_dIdX = false, compute_d2I = true;
281  functional.evaluate_functional( compute_dIdW, compute_dIdX, compute_d2I );
282  functional.d2IdXdX->vmult(d2IdXdXp_input, dXvdXp_input);
283  }
284 
285  //auto d2IdXvsdXp_input = functional.dg->high_order_grid->volume_nodes;
286  //{
287  // dealii::LinearAlgebra::distributed::Vector<double> dummy_vector(functional.dg->high_order_grid->surface_nodes);
288  // MeshMover::LinearElasticity<dim, double, dealii::LinearAlgebra::distributed::Vector<double>, dealii::DoFHandler<dim>>
289  // meshmover(*(functional.dg->high_order_grid->triangulation),
290  // functional.dg->high_order_grid->initial_mapping_fe_field,
291  // functional.dg->high_order_grid->dof_handler_grid,
292  // functional.dg->high_order_grid->surface_to_volume_indices,
293  // dummy_vector);
294  // meshmover.apply_dXvdXvs_transpose(d2IdXdXp_input, d2IdXvsdXp_input);
295  //}
296 
297  //auto &d2IdXpdXp_input = ROL_vector_to_dealii_vector_reference(output_vector);
298  //{
299  // dXvsdXp.Tvmult(d2IdXpdXp_input, d2IdXvsdXp_input);
300  //}
301 
302  auto &dealii_output = ROL_vector_to_dealii_vector_reference(output_vector);
303  dXvdXp.Tvmult(dealii_output, d2IdXdXp_input);
304 
305  //n_vmult += 3;
306 }
307 
308 template class ROLObjectiveSimOpt <PHILIP_DIM,1>;
309 template class ROLObjectiveSimOpt <PHILIP_DIM,2>;
310 template class ROLObjectiveSimOpt <PHILIP_DIM,3>;
311 template class ROLObjectiveSimOpt <PHILIP_DIM,4>;
312 template class ROLObjectiveSimOpt <PHILIP_DIM,5>;
313 
314 } // PHiLiP namespace
std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > d2IdXdX
Sparse matrix for storing the functional partial second derivatives.
Definition: functional.h:128
dealii::LinearAlgebra::distributed::Vector< real > dIdX
Vector for storing the derivatives with respect to each grid DoF.
Definition: functional.h:121
void hessVec_22(ROL::Vector< double > &output_vector, const ROL::Vector< double > &input_vector, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override
Applies the functional Hessian w. r. t. the control variables onto a vector.
Interface between the ROL::Objective_SimOpt PHiLiP::Functional.
void set_state(const dealii::LinearAlgebra::distributed::Vector< real > &solution_set)
Definition: functional.cpp:376
dealii::LinearAlgebra::distributed::Vector< real > dIdw
Vector for storing the derivatives with respect to each solution DoF.
Definition: functional.h:119
std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > d2IdWdW
Sparse matrix for storing the functional partial second derivatives.
Definition: functional.h:124
void hessVec_11(ROL::Vector< double > &output_vector, const ROL::Vector< double > &input_vector, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override
Applies the functional Hessian w. r. t. the simulation variables onto a vector.
std::shared_ptr< BaseParameterization< dim > > design_parameterization
Design parameterization to link design variables with volume nodes.
void gradient_1(ROL::Vector< double > &g, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override
Returns the gradient w. r. t. the simulation variables of the Functional object.
Files for the baseline physics.
Definition: ADTypes.hpp:10
Abstract class for design parameterization. Objective function and the constraints take this class&#39;s ...
void gradient_2(ROL::Vector< double > &g, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override
Returns the gradient w. r. t. the control variables of the Functional object.
void update(const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, bool flag=true, int iter=-1) override
Update the simulation and control variables.
dealii::TrilinosWrappers::SparseMatrix dXvdXp
Stored mesh sensitivity evaluated at initialization.
void hessVec_12(ROL::Vector< double > &output_vector, const ROL::Vector< double > &input_vector, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override
Applies the functional Hessian w. r. t. the simulation and control variables onto a vector...
void hessVec_21(ROL::Vector< double > &output_vector, const ROL::Vector< double > &input_vector, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override
Applies the functional Hessian w. r. t. the control and simulation variables onto a vector...
Functional< dim, nstate, double > & functional
Functional to be evaluated.
ROLObjectiveSimOpt(Functional< dim, nstate, double > &_functional, std::shared_ptr< BaseParameterization< dim >> _design_parameterization, std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > precomputed_dXvdXp=nullptr)
Constructor.
dealii::LinearAlgebra::distributed::Vector< double > design_var
Design variables.
std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > d2IdWdX
Sparse matrix for storing the functional partial second derivatives.
Definition: functional.h:126
std::shared_ptr< DGBase< dim, real, MeshType > > dg
Smart pointer to DGBase.
Definition: functional.h:50
virtual real evaluate_functional(const bool compute_dIdW=false, const bool compute_dIdX=false, const bool compute_d2I=false)
Evaluates the functional derivative with respect to the solution variable.
Definition: functional.cpp:795
const dealii::LinearAlgebra::distributed::Vector< double > & ROL_vector_to_dealii_vector_reference(const ROL::Vector< double > &x)
Access the read-write deali.II Vector stored within the ROL::Vector.
double value(const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override
Returns the value of the Functional object.