[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.hpp
1 #ifndef __ROLOBJECTIVESIMOPT_H__
2 #define __ROLOBJECTIVESIMOPT_H__
3 
4 #include "ROL_Objective_SimOpt.hpp"
5 
6 #include "mesh/free_form_deformation.h"
7 
8 #include "functional/functional.h"
9 
10 #include "design_parameterization/base_parameterization.hpp"
11 
12 namespace PHiLiP {
13 
14 using ROL_Vector = ROL::Vector<double>;
15 
17 
22 template <int dim, int nstate>
23 class ROLObjectiveSimOpt : public ROL::Objective_SimOpt<double> {
24 private:
27 
29  std::shared_ptr<BaseParameterization<dim>> design_parameterization;
30 
32  dealii::LinearAlgebra::distributed::Vector<double> design_var;
33 
34 public:
35 
37  dealii::TrilinosWrappers::SparseMatrix dXvdXp;
38 
41  Functional<dim,nstate,double> &_functional,
42  std::shared_ptr<BaseParameterization<dim>> _design_parameterization,
43  std::shared_ptr<dealii::TrilinosWrappers::SparseMatrix> precomputed_dXvdXp = nullptr);
44 
45  using ROL::Objective_SimOpt<double>::value;
46  using ROL::Objective_SimOpt<double>::update;
47 
49  void update(
50  const ROL::Vector<double> &des_var_sim,
51  const ROL::Vector<double> &des_var_ctl,
52  bool flag = true,
53  int iter = -1) override;
54 
56  double value(
57  const ROL::Vector<double> &des_var_sim,
58  const ROL::Vector<double> &des_var_ctl,
59  double &/*tol*/ ) override;
60 
62  void gradient_1(
63  ROL::Vector<double> &g,
64  const ROL::Vector<double> &des_var_sim,
65  const ROL::Vector<double> &des_var_ctl,
66  double &/*tol*/ ) override;
67 
69  void gradient_2(
70  ROL::Vector<double> &g,
71  const ROL::Vector<double> &des_var_sim,
72  const ROL::Vector<double> &des_var_ctl,
73  double &/*tol*/ ) override;
74 
76 
81  void hessVec_11(
82  ROL::Vector<double> &output_vector,
83  const ROL::Vector<double> &input_vector,
84  const ROL::Vector<double> &des_var_sim,
85  const ROL::Vector<double> &des_var_ctl,
86  double &/*tol*/ ) override;
87 
89 
94  void hessVec_12(
95  ROL::Vector<double> &output_vector,
96  const ROL::Vector<double> &input_vector,
97  const ROL::Vector<double> &des_var_sim,
98  const ROL::Vector<double> &des_var_ctl,
99  double &/*tol*/ ) override;
100 
102 
107  void hessVec_21(
108  ROL::Vector<double> &output_vector,
109  const ROL::Vector<double> &input_vector,
110  const ROL::Vector<double> &des_var_sim,
111  const ROL::Vector<double> &des_var_ctl,
112  double &/*tol*/ ) override;
113 
115 
120  void hessVec_22(
121  ROL::Vector<double> &output_vector,
122  const ROL::Vector<double> &input_vector,
123  const ROL::Vector<double> &des_var_sim,
124  const ROL::Vector<double> &des_var_ctl,
125  double &/*tol*/ ) override;
126 
127 }; // ROLObjectiveSimOpt
128 
129 } // PHiLiP namespace
130 
131 #endif
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 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.
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.