|
[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
#include <adjoint.h>
Public Types | |
| enum | AdjointStateEnum { coarse, fine } |
| For storing the current state in the adjoint. More... | |
Public Member Functions | |
| Adjoint (std::shared_ptr< DGBase< dim, real, MeshType > > _dg, std::shared_ptr< Functional< dim, nstate, real, MeshType > > _functional, std::shared_ptr< Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< real >> > _physics) | |
| Constructor. More... | |
| void | reinit () |
| Reinitialize Adjoint with the same pointers. More... | |
| void | convert_to_state (AdjointStateEnum state) |
| Converts the adjoint to specified state. More... | |
| void | coarse_to_fine () |
| Projects the problem to a p-enriched space. More... | |
| void | fine_to_coarse () |
| Return the problem to the original solution and polynomial distribution. More... | |
| dealii::LinearAlgebra::distributed::Vector< real > | fine_grid_adjoint () |
| Computes the fine grid adjoint. More... | |
| dealii::LinearAlgebra::distributed::Vector< real > | coarse_grid_adjoint () |
| Computes the coarse grid adjoint. More... | |
| dealii::Vector< real > | dual_weighted_residual () |
| compute the Dual Weighted Residual (DWR) More... | |
| void | output_results_vtk (const unsigned int cycle) |
| Outputs the current solution and adjoint values. More... | |
Public Attributes | |
| std::shared_ptr< DGBase< dim, real, MeshType > > | dg |
| DG class pointer. | |
| std::shared_ptr< Functional< dim, nstate, real, MeshType > > | functional |
| Functional class pointer. | |
| std::shared_ptr< Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< real > > > | physics |
| Problem physics (for calling the functional class) | |
| std::shared_ptr< MeshType > | triangulation |
| Grid. | |
| dealii::LinearAlgebra::distributed::Vector< real > | solution_coarse |
| original solution | |
| dealii::LinearAlgebra::distributed::Vector< real > | dIdw_fine |
| functional derivative (on the fine grid) | |
| dealii::LinearAlgebra::distributed::Vector< real > | dIdw_coarse |
| functional derivative (on the coarse grid) | |
| dealii::LinearAlgebra::distributed::Vector< real > | adjoint_fine |
| fine grid adjoint ( \(\psi_h\)) | |
| dealii::LinearAlgebra::distributed::Vector< real > | adjoint_coarse |
| coarse grid adjoint ( \(\psi_H\)) | |
| dealii::Vector< real > | dual_weighted_residual_fine |
| dual weighted residual More... | |
| dealii::Vector< real > | coarse_fe_index |
| Original FE_index distribution. | |
| AdjointStateEnum | adjoint_state |
| Current adjoint state. | |
Protected Attributes | |
| MPI_Comm | mpi_communicator |
| MPI communicator. | |
| dealii::ConditionalOStream | pcout |
| Parallel std::cout that only outputs on mpi_rank==0. | |
Adjoint class.
This class computes the discrete adjoint of the system based on a functional of interest and a computed DG solution. Uses the Sacado functions Functional::evaluate_dIdw() and DGBase::assemble_residual() to generate and solve the discrete adjoint system
\[ \left( \frac{\partial \mathbf{R}}{\partial \mathbf{u}} \right)^T \psi + \left(\frac{\partial \mathcal{J}}{\partial \mathbf{u}}\right)^T = \mathbf{0} \]
Includes functions for solving both the coarse and fine \(p\)-enriched adjoint problems. Subscripts \(H\) and \(h\) are used to denote coarse and fine grid variables respectively.
| enum PHiLiP::Adjoint::AdjointStateEnum |
| PHiLiP::Adjoint< dim, nstate, real, MeshType >::Adjoint | ( | std::shared_ptr< DGBase< dim, real, MeshType > > | _dg, |
| std::shared_ptr< Functional< dim, nstate, real, MeshType > > | _functional, | ||
| std::shared_ptr< Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< real >> > | _physics | ||
| ) |
Constructor.
Initializes the Adjoint as being in the AdjointEnum::coarse state. Also stores the current solution and distribution of polynomial orders for the mesh for converting back to coarse state after refinement.
Definition at line 30 of file adjoint.cpp.
| dealii::LinearAlgebra::distributed::Vector< real > PHiLiP::Adjoint< dim, nstate, real, MeshType >::coarse_grid_adjoint | ( | ) |
Computes the coarse grid adjoint.
Reverts the state to the coarse grid (if needed) and solves for Adjoint::adjoint_coarse from
\[ \left(\left. \frac{\partial \mathbf{R}_H}{\partial \mathbf{u}} \right|_{\mathbf{u}_H}\right)^T \psi_H + \left(\left. \frac{\partial \mathcal{J}_H}{\partial \mathbf{u}} \right|_{\mathbf{u}_H}\right)^T=\mathbf{0} \]
Definition at line 187 of file adjoint.cpp.
| void PHiLiP::Adjoint< dim, nstate, real, MeshType >::coarse_to_fine | ( | ) |
Projects the problem to a p-enriched space.
Raises the FE_index on each cell and transfers the coarse solution to a fine solution (stored in DGBase::solution)
Definition at line 93 of file adjoint.cpp.
| void PHiLiP::Adjoint< dim, nstate, real, MeshType >::convert_to_state | ( | AdjointStateEnum | state | ) |
Converts the adjoint to specified state.
Calls the functions coarse_to_fine() or fine_to_coarse() if the Adjoint::adjoint_state is different than the input state
Definition at line 78 of file adjoint.cpp.
| dealii::Vector< real > PHiLiP::Adjoint< dim, nstate, real, MeshType >::dual_weighted_residual | ( | ) |
compute the Dual Weighted Residual (DWR)
Computes Adjoint::dual_weighted_resiudal_fine ( \(\eta\)) on the fine grid. This value should be zero on the coarse grid due to Galerkin Orthogonality. It is calculated from
\[ \eta = \mathbf{R}_h(\mathbf{u}_h^H)^T \psi_h \]
Uses Adjoint::adjoint_fine and should only be called after fine_grid_adjoint().
Definition at line 216 of file adjoint.cpp.
| dealii::LinearAlgebra::distributed::Vector< real > PHiLiP::Adjoint< dim, nstate, real, MeshType >::fine_grid_adjoint | ( | ) |
Computes the fine grid adjoint.
Converts the state to a refined grid (if needed) and solves for Adjoint::adjoint_fine from
\[ \left(\left. \frac{\partial \mathbf{R}_h}{\partial \mathbf{u}} \right|_{\mathbf{u}_h^H}\right)^T \psi_h + \left(\left. \frac{\partial \mathcal{J}_h}{\partial \mathbf{u}} \right|_{\mathbf{u}_h^H}\right)^T=\mathbf{0} \]
where \(\mathbf{u}_h^H\) is the projected solution on the fine grid.
Definition at line 157 of file adjoint.cpp.
| void PHiLiP::Adjoint< dim, nstate, real, MeshType >::fine_to_coarse | ( | ) |
Return the problem to the original solution and polynomial distribution.
Copies the values that were stored in solution_coarse and Adjoint::coarse_fe_index at intilization
Definition at line 136 of file adjoint.cpp.
| void PHiLiP::Adjoint< dim, nstate, real, MeshType >::output_results_vtk | ( | const unsigned int | cycle | ) |
Outputs the current solution and adjoint values.
Similar to DGBase::output_results_vtk() but will also include the adjoint and dIdw related to the current adjoint state. Will also output Adjoint::dual_weighted_residual_fine if currenly on the fine grid.
Definition at line 249 of file adjoint.cpp.
| void PHiLiP::Adjoint< dim, nstate, real, MeshType >::reinit | ( | ) |
Reinitialize Adjoint with the same pointers.
Sets adjoint_state to AdjointEnum::coarse and stores the current solution and polynomial order distribution
Definition at line 53 of file adjoint.cpp.
| dealii::Vector<real> PHiLiP::Adjoint< dim, nstate, real, MeshType >::dual_weighted_residual_fine |