[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
DualWeightedResidualError class. More...
#include <mesh_error_estimate.h>
Public Types | |
enum | SolutionRefinementStateEnum { coarse, fine } |
For storing the current refinement state of the solution. More... | |
Public Member Functions | |
DualWeightedResidualError (std::shared_ptr< DGBase< dim, real, MeshType > > dg_input) | |
Constructor. More... | |
void | reinit () |
Reinitializes member variables of DualWeightedResidualError. More... | |
void | convert_dgsolution_to_coarse_or_fine (SolutionRefinementStateEnum refinement_state) |
Converts DG solution to the 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 () |
Computes the Dual Weighted Residual (DWR) More... | |
dealii::Vector< real > | compute_cellwise_errors () override |
Computes dual weighted residual error in each cell, by integrating over all quadrature points. Overwrites the virtual function in MeshErrorEstimateBase. | |
real | total_dual_weighted_residual_error () |
Computes the sum of dual weighted residual error over all the cells in the domain. | |
void | output_results_vtk (const unsigned int cycle) |
Outputs the current solution and adjoint values. More... | |
dealii::LinearAlgebra::distributed::Vector< real > | compute_adjoint (dealii::LinearAlgebra::distributed::Vector< real > &derivative_functional_wrt_solution, dealii::LinearAlgebra::distributed::Vector< real > &adjoint_variable) |
Solves the adjoint equation. | |
![]() | |
MeshErrorEstimateBase (std::shared_ptr< DGBase< dim, real, MeshType > > dg_input) | |
Constructor. | |
virtual | ~MeshErrorEstimateBase ()=default |
Virtual Destructor. | |
Public Attributes | |
std::shared_ptr< Functional< dim, nstate, real, MeshType > > | functional |
Functional class pointer. | |
dealii::LinearAlgebra::distributed::Vector< real > | solution_coarse |
original solution | |
dealii::LinearAlgebra::distributed::Vector< real > | derivative_functional_wrt_solution_fine |
functional derivative (on the fine grid) | |
dealii::LinearAlgebra::distributed::Vector< real > | derivative_functional_wrt_solution_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 (. | |
Original FE_index distribution dealii::Vector< real > | coarse_fe_index |
Current refinement state of the solution SolutionRefinementStateEnum | solution_refinement_state |
![]() | |
std::shared_ptr< DGBase< dim, real, MeshType > > | dg |
Pointer to DGBase. | |
Protected Attributes | |
MPI_Comm | mpi_communicator |
MPI communicator. | |
dealii::ConditionalOStream | pcout |
Parallel std::cout that only outputs on mpi_rank==0. | |
DualWeightedResidualError 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. Reference: Venditti AND Darmafol, "Adjoint Error Estimation and Grid Adaptation for Functional Outputs: Application to Quasi-One-Dimensional Flow". Journal of Computational Physics 164, 1 (2000), 204–227.
Definition at line 89 of file mesh_error_estimate.h.
enum PHiLiP::DualWeightedResidualError::SolutionRefinementStateEnum |
For storing the current refinement state of the solution.
Enumerator | |
---|---|
coarse | Initial state. |
fine | Refined state. |
Definition at line 94 of file mesh_error_estimate.h.
|
explicit |
Constructor.
Initializes the solution as being in the SolutionRefinementStateEnum::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 70 of file mesh_error_estimate.cpp.
dealii::LinearAlgebra::distributed::Vector< real > PHiLiP::DualWeightedResidualError< dim, nstate, real, MeshType >::coarse_grid_adjoint | ( | ) |
Computes the coarse grid adjoint.
Reverts the state to the coarse grid (if needed) and solves for DualWeightedResidualError::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} \]
Eq(7) from Venditti and Darmafol (2000), cited above, using coarse-grid solution.
Definition at line 270 of file mesh_error_estimate.cpp.
void PHiLiP::DualWeightedResidualError< 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 174 of file mesh_error_estimate.cpp.
void PHiLiP::DualWeightedResidualError< dim, nstate, real, MeshType >::convert_dgsolution_to_coarse_or_fine | ( | SolutionRefinementStateEnum | refinement_state | ) |
Converts DG solution to the specified state.
Calls the functions coarse_to_fine() or fine_to_coarse() if the DualWeightedResidualError::solution_refinement_state is different than the input state
Definition at line 149 of file mesh_error_estimate.cpp.
dealii::Vector< real > PHiLiP::DualWeightedResidualError< dim, nstate, real, MeshType >::dual_weighted_residual | ( | ) |
Computes the Dual Weighted Residual (DWR)
Computes DualWeightedResidualError::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 DualWeightedResidualError::adjoint_fine and should only be called after fine_grid_adjoint(). Eq(6) from Venditti and Darmafol (2000), cited above.
Definition at line 309 of file mesh_error_estimate.cpp.
dealii::LinearAlgebra::distributed::Vector< real > PHiLiP::DualWeightedResidualError< dim, nstate, real, MeshType >::fine_grid_adjoint | ( | ) |
Computes the fine grid adjoint.
Converts the state to a refined grid (if needed) and solves for DualWeightedResidualError::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. Eq(7) from Venditti and Darmafol (2000), cited above.
Definition at line 260 of file mesh_error_estimate.cpp.
void PHiLiP::DualWeightedResidualError< 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 DualWeightedResidualError::coarse_fe_index at intilization
Definition at line 231 of file mesh_error_estimate.cpp.
void PHiLiP::DualWeightedResidualError< 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 derivative_functional_wrt_solution related to the current adjoint state. Will also output DualWeightedResidualError::dual_weighted_residual_fine if currenly on the fine grid.
Definition at line 344 of file mesh_error_estimate.cpp.
void PHiLiP::DualWeightedResidualError< dim, nstate, real, MeshType >::reinit | ( | ) |
Reinitializes member variables of DualWeightedResidualError.
Sets solution_refinement_state to SolutionRefinementStateEnum::coarse and stores the current solution and polynomial order distribution
Definition at line 121 of file mesh_error_estimate.cpp.