[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
Use DGBase as a Simulation Constraint ROL::Constraint_SimOpt. More...
#include <flow_constraints.hpp>
Public Member Functions | |
FlowConstraints (std::shared_ptr< DGBase< dim, double >> &_dg, std::shared_ptr< BaseParameterization< dim >> _design_parameterization, std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > precomputed_dXvdXp=nullptr) | |
Constructor. | |
~FlowConstraints () | |
Destructor. | |
void | update_1 (const ROL::Vector< double > &des_var_sim, bool flag=true, int iter=-1) |
Update the simulation variables. | |
void | update_2 (const ROL::Vector< double > &des_var_ctl, bool flag=true, int iter=-1) |
Update the control variables. More... | |
void | solve (ROL::Vector< double > &constraint_values, ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override |
Solve the Simulation Constraint and returns the constraints values. More... | |
void | value (ROL::Vector< double > &constraint_values, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &) override |
Returns the current constraint residual given a set of control and simulation variables. | |
void | applyJacobian_1 (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 Jacobian of the Constraints w. r. t. the simulation variables onto a vector. | |
void | applyAdjointJacobian_1 (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 Jacobian of the Constraints w. r. t. the simulation variables onto a vector. | |
void | applyInverseJacobian_1 (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 inverse Jacobian of the Constraints w. r. t. the simulation variables onto a vector. | |
int | construct_JacobianPreconditioner_1 (const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl) |
Constructs the Jacobian preconditioner. | |
void | destroy_JacobianPreconditioner_1 () |
Frees Jacobian preconditioner from memory;. | |
void | applyInverseJacobianPreconditioner_1 (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 &) |
Applies the inverse Jacobian preconditioner. More... | |
int | construct_AdjointJacobianPreconditioner_1 (const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl) |
Constructs the Adjoint Jacobian preconditioner. | |
void | destroy_AdjointJacobianPreconditioner_1 () |
Frees adjoint Jacobian preconditioner from memory;. | |
void | applyInverseAdjointJacobianPreconditioner_1 (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 &) |
Applies the inverse Adjoint Jacobian preconditioner. More... | |
void | applyInverseAdjointJacobian_1 (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 adjoint Jacobian of the Constraints w. r. t. the simulation variables onto a vector. | |
void | applyJacobian_2 (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 Jacobian of the Constraints w. r. t. the control variables onto a vector. | |
void | applyAdjointJacobian_2 (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 Jacobian of the Constraints w. r. t. the control variables onto a vector. | |
void | applyAdjointHessian_11 (ROL::Vector< double > &output_vector, const ROL::Vector< double > &dual, const ROL::Vector< double > &input_vector, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &tol) override |
Applies the adjoint of the Hessian of the constraints w. r. t. the simulation variables onto a vector. More... | |
void | applyAdjointHessian_12 (ROL::Vector< double > &output_vector, const ROL::Vector< double > &dual, const ROL::Vector< double > &input_vector, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &tol) override |
Applies the adjoint of the Hessian of the constraints w. r. t. the simulation variables onto a vector. More... | |
void | applyAdjointHessian_21 (ROL::Vector< double > &output_vector, const ROL::Vector< double > &dual, const ROL::Vector< double > &input_vector, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &tol) override |
Applies the adjoint of the Hessian of the constraints w. r. t. the simulation variables onto a vector. More... | |
void | applyAdjointHessian_22 (ROL::Vector< double > &output_vector, const ROL::Vector< double > &dual, const ROL::Vector< double > &input_vector, const ROL::Vector< double > &des_var_sim, const ROL::Vector< double > &des_var_ctl, double &tol) override |
Applies the adjoint of the Hessian of the constraints w. r. t. the simulation variables onto a vector. More... | |
Public Attributes | |
dealii::TrilinosWrappers::SparseMatrix | dXvdXp |
Stores the mesh sensitivities. | |
double | flow_CFL_ |
Regularization of the constraint by adding flow_CFL_ times the mass matrix. | |
Protected Attributes | |
int | i_out = 1000 |
ID used when outputting the flow solution. | |
int | iupdate = 9000 |
ID used when outputting the flow solution. | |
Private Attributes | |
const int | mpi_rank |
MPI rank used for printing. | |
const bool | i_print |
Whether the current processor should print or not. | |
std::shared_ptr< DGBase< dim, double > > | dg |
Smart pointer to DGBase. | |
std::shared_ptr< BaseParameterization< dim > > | design_parameterization |
Parameterization which links design variables to the volume nodes. | |
Parameters::LinearSolverParam | linear_solver_param |
Linear solver parameters. More... | |
dealii::LinearAlgebra::distributed::Vector< double > | design_var |
Design variables values. | |
Ifpack_Preconditioner * | jacobian_prec |
Jacobian preconditioner. More... | |
Ifpack_Preconditioner * | adjoint_jacobian_prec |
Adjoint Jacobian preconditioner. More... | |
Use DGBase as a Simulation Constraint ROL::Constraint_SimOpt.
The simulation variables will be the DoFs stored in the DGBase::solution.
The control variables will be some of the the FFD points/directions. The given ffd_design_variables_indices_dim
point to the points/directions used as design variables.
Definition at line 26 of file flow_constraints.hpp.
|
override |
Applies the adjoint of the Hessian of the constraints w. r. t. the simulation variables onto a vector.
More specifically, apply
\[ \mathbf{v}_{output} = \left( \sum_i \psi_i \frac{\partial^2 R_i}{\partial u \partial u} \right)^T \mathbf{v}_{input} \]
onto the input_vector
to obtain the output_vector
Definition at line 526 of file flow_constraints.cpp.
|
override |
Applies the adjoint of the Hessian of the constraints w. r. t. the simulation variables onto a vector.
More specifically, apply
\[ \mathbf{v}_{output} = \left( \sum_i \psi_i \frac{\partial^2 R_i}{\partial u \partial x} \right)^T \mathbf{v}_{input} \]
onto the input_vector
to obtain the output_vector
Definition at line 553 of file flow_constraints.cpp.
|
override |
Applies the adjoint of the Hessian of the constraints w. r. t. the simulation variables onto a vector.
More specifically, apply
\[ \mathbf{v}_{output} = \left( \sum_i \psi_i \frac{\partial^2 R_i}{\partial x \partial u} \right)^T \mathbf{v}_{input} \]
onto the input_vector
to obtain the output_vector
Definition at line 608 of file flow_constraints.cpp.
|
override |
Applies the adjoint of the Hessian of the constraints w. r. t. the simulation variables onto a vector.
More specifically, apply
\[ \mathbf{v}_{output} = \left( \sum_i \psi_i \frac{\partial^2 R_i}{\partial x \partial x} \right)^T \mathbf{v}_{input} \]
onto the input_vector
to obtain the output_vector
[out] | output_vector | Resulting vector \( \mathbf{v}_{output} \in \mathbb{R}^{N_{ctl}} \) |
[in] | dual | Lagrange multiplier associated with constraint vector \( \boldsymbol{\psi} \in \mathbb{R}^{N_{residual}} \) |
[in] | input_vector | Input vector \( \mathbf{v}_{input} \in \mathbb{R}^{N_{ctl}} \) |
[in] | des_var_sim | Simulation variables \( \mathbf{v}_{input} \in \mathbb{R}^{N_{flow}} \) |
[in] | des_var_ctl | Simulation variables \( \mathbf{v}_{input} \in \mathbb{R}^{N_{ctl}} \) |
[in] | tol | Tolerance, not used. From virtual ROL function. |
Definition at line 667 of file flow_constraints.cpp.
void PHiLiP::FlowConstraints< dim >::applyInverseAdjointJacobianPreconditioner_1 | ( | 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 & | |||
) |
Applies the inverse Adjoint Jacobian preconditioner.
construct_AdjointJacobianPreconditioner_1 needs to be called beforehand.
Definition at line 348 of file flow_constraints.cpp.
void PHiLiP::FlowConstraints< dim >::applyInverseJacobianPreconditioner_1 | ( | 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 & | |||
) |
Applies the inverse Jacobian preconditioner.
construct_JacobianPreconditioner_1 needs to be called beforehand.
Definition at line 316 of file flow_constraints.cpp.
|
override |
Solve the Simulation Constraint and returns the constraints values.
In this case, we use the ODESolver to solve the steady state solution of the flow given the geometry.
Definition at line 96 of file flow_constraints.cpp.
void PHiLiP::FlowConstraints< dim >::update_2 | ( | const ROL::Vector< double > & | des_var_ctl, |
bool | flag = true , |
||
int | iter = -1 |
||
) |
Update the control variables.
Update FFD from design variables and then deforms the mesh.
Definition at line 81 of file flow_constraints.cpp.
|
private |
Adjoint Jacobian preconditioner.
Currently uses ILUT
Definition at line 51 of file flow_constraints.hpp.
|
private |
|
private |
Linear solver parameters.
Currently set such that the linear systems are fully solved
Definition at line 41 of file flow_constraints.hpp.