|
[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
#include <target_boundary_functional.h>
Public Member Functions | |
| TargetBoundaryFunctional (std::shared_ptr< DGBase< dim, real >> dg_input, const dealii::LinearAlgebra::distributed::Vector< real > &target_solution, const bool uses_solution_values=true, const bool uses_solution_gradient=false) | |
| Constructor. | |
| template<typename real2 > | |
| real2 | evaluate_volume_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, real2 > &, const dealii::Point< dim, real2 > &, const std::array< real2, nstate > &, const std::array< real, nstate > &, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &) const |
| Zero out the default inverse target volume functional. | |
| real | evaluate_volume_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &physics, const dealii::Point< dim, real > &phys_coord, const std::array< real, nstate > &soln_at_q, const std::array< real, nstate > &target_soln_at_q, const std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_at_q, const std::array< dealii::Tensor< 1, dim, real >, nstate > &target_soln_grad_at_q) const override |
| non-template functions to override the template classes | |
| FadFadType | evaluate_volume_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, FadFadType > &physics, const dealii::Point< dim, FadFadType > &phys_coord, const std::array< FadFadType, nstate > &soln_at_q, const std::array< real, nstate > &target_soln_at_q, const std::array< dealii::Tensor< 1, dim, FadFadType >, nstate > &soln_grad_at_q, const std::array< dealii::Tensor< 1, dim, FadFadType >, nstate > &target_soln_grad_at_q) const override |
| non-template functions to override the template classes | |
Public Member Functions inherited from PHiLiP::TargetFunctional< dim, nstate, real > | |
| TargetFunctional (std::shared_ptr< DGBase< dim, real >> dg_input, const bool uses_solution_values=true, const bool uses_solution_gradient=true) | |
| Constructor. More... | |
| TargetFunctional (std::shared_ptr< DGBase< dim, real >> dg_input, const dealii::LinearAlgebra::distributed::Vector< real > &target_solution, const bool uses_solution_values=true, const bool uses_solution_gradient=true) | |
| Constructor. More... | |
| TargetFunctional (std::shared_ptr< DGBase< dim, real >> dg_input, std::shared_ptr< Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< Sacado::Fad::DFad< real >>> > _physics_fad_fad, const bool uses_solution_values=true, const bool uses_solution_gradient=true) | |
| Constructor. More... | |
| TargetFunctional (std::shared_ptr< DGBase< dim, real >> dg_input, const dealii::LinearAlgebra::distributed::Vector< real > &target_solution, std::shared_ptr< Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< Sacado::Fad::DFad< real >>> > _physics_fad_fad, const bool uses_solution_values=true, const bool uses_solution_gradient=true) | |
| Constructor. More... | |
| virtual real | evaluate_functional (const bool compute_dIdW=false, const bool compute_dIdX=false, const bool compute_d2I=false) override |
| Evaluates the functional derivative with respect to the solution variable. More... | |
| dealii::LinearAlgebra::distributed::Vector< real > | evaluate_dIdw_finiteDifferences (DGBase< dim, real > &dg, const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &physics, const double stepsize) |
| dealii::LinearAlgebra::distributed::Vector< real > | evaluate_dIdX_finiteDifferences (DGBase< dim, real > &dg, const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &physics, const double stepsize) |
Public Member Functions inherited from PHiLiP::Functional< dim, nstate, real > | |
| virtual | ~Functional ()=default |
| Destructor. | |
| Functional (std::shared_ptr< PHiLiP::DGBase< dim, real, dealii::parallel::distributed::Triangulation< dim > >> _dg, const bool _uses_solution_values=true, const bool _uses_solution_gradient=true) | |
| Functional (std::shared_ptr< PHiLiP::DGBase< dim, real, dealii::parallel::distributed::Triangulation< dim > >> _dg, std::shared_ptr< PHiLiP::Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< Sacado::Fad::DFad< real >> >> _physics_fad_fad, const bool _uses_solution_values=true, const bool _uses_solution_gradient=true) | |
| void | set_state (const dealii::LinearAlgebra::distributed::Vector< real > &solution_set) |
| void | set_geom (const dealii::LinearAlgebra::distributed::Vector< real > &volume_nodes_set) |
| dealii::LinearAlgebra::distributed::Vector< real > | evaluate_dIdw_finiteDifferences (DGBase< dim, real, dealii::parallel::distributed::Triangulation< dim > > &dg, const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &physics, const double stepsize) |
| dealii::LinearAlgebra::distributed::Vector< real > | evaluate_dIdX_finiteDifferences (DGBase< dim, real, dealii::parallel::distributed::Triangulation< dim > > &dg, const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &physics, const double stepsize) |
Private Types | |
| using | FadType = Sacado::Fad::DFad< real > |
| Sacado AD type for first derivatives. | |
| using | FadFadType = Sacado::Fad::DFad< FadType > |
| Sacado AD type that allows 2nd derivatives. | |
Additional Inherited Members | |
Public Types inherited from PHiLiP::TargetFunctional< dim, nstate, real > | |
| using | FadType = Sacado::Fad::DFad< real > |
| Sacado AD type for first derivatives. | |
| using | FadFadType = Sacado::Fad::DFad< FadType > |
| Sacado AD type that allows 2nd derivatives. | |
Public Attributes inherited from PHiLiP::Functional< dim, nstate, real > | |
| std::shared_ptr< DGBase< dim, real, dealii::parallel::distributed::Triangulation< dim > > > | dg |
| Smart pointer to DGBase. | |
| real | current_functional_value |
| Store the functional value from the last time evaluate_functional() was called. | |
| dealii::LinearAlgebra::distributed::Vector< real > | dIdw |
| Vector for storing the derivatives with respect to each solution DoF. | |
| dealii::LinearAlgebra::distributed::Vector< real > | dIdX |
| Vector for storing the derivatives with respect to each grid DoF. | |
| std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > | d2IdWdW |
| Sparse matrix for storing the functional partial second derivatives. | |
| std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > | d2IdWdX |
| Sparse matrix for storing the functional partial second derivatives. | |
| std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > | d2IdXdX |
| Sparse matrix for storing the functional partial second derivatives. | |
Protected Member Functions inherited from PHiLiP::TargetFunctional< dim, nstate, real > | |
| virtual real | evaluate_volume_cell_functional (const Physics::PhysicsBase< dim, nstate, real > &physics, const std::vector< real > &soln_coeff, const std::vector< real > &target_soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< real > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const dealii::Quadrature< dim > &volume_quadrature) const |
| Corresponding real function to evaluate a cell's volume functional. | |
| virtual Sacado::Fad::DFad< Sacado::Fad::DFad< real > > | evaluate_volume_cell_functional (const Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< Sacado::Fad::DFad< real >>> &physics_fad_fad, const std::vector< Sacado::Fad::DFad< Sacado::Fad::DFad< real >> > &soln_coeff, const std::vector< real > &target_soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< Sacado::Fad::DFad< Sacado::Fad::DFad< real >> > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const dealii::Quadrature< dim > &volume_quadrature) const |
| Corresponding FadFadType function to evaluate a cell's volume functional. | |
| virtual real | evaluate_boundary_cell_functional (const Physics::PhysicsBase< dim, nstate, real > &physics, const unsigned int boundary_id, const std::vector< real > &soln_coeff, const std::vector< real > &target_soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< real > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const dealii::Quadrature< dim-1 > &face_quadrature, const unsigned int face_number) const |
| Corresponding real function to evaluate a cell's face functional. | |
| virtual Sacado::Fad::DFad< Sacado::Fad::DFad< real > > | evaluate_boundary_cell_functional (const Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< Sacado::Fad::DFad< real >>> &physics_fad_fad, const unsigned int boundary_id, const std::vector< Sacado::Fad::DFad< Sacado::Fad::DFad< real >> > &soln_coeff, const std::vector< real > &target_soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< Sacado::Fad::DFad< Sacado::Fad::DFad< real >> > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const dealii::Quadrature< dim-1 > &face_quadrature, const unsigned int face_number) const |
| Corresponding FadFadType function to evaluate a cell's face functional. | |
| virtual real | evaluate_boundary_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &, const unsigned int, const dealii::Point< dim, real > &, const dealii::Tensor< 1, dim, real > &, const std::array< real, nstate > &soln_at_q, const std::array< real, nstate > &target_soln_at_q, const std::array< dealii::Tensor< 1, dim, real >, nstate > &, const std::array< dealii::Tensor< 1, dim, real >, nstate > &) const |
| Virtual function for computation of cell face functional term. More... | |
| virtual FadFadType | evaluate_boundary_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, FadFadType > &, const unsigned int, const dealii::Point< dim, FadFadType > &, const dealii::Tensor< 1, dim, FadFadType > &, const std::array< FadFadType, nstate > &soln_at_q, const std::array< real, nstate > &target_soln_at_q, const std::array< dealii::Tensor< 1, dim, FadFadType >, nstate > &, const std::array< dealii::Tensor< 1, dim, FadFadType >, nstate > &) const |
| Virtual function for Sacado computation of cell face functional term and derivatives. More... | |
Protected Member Functions inherited from PHiLiP::Functional< dim, nstate, real > | |
| void | allocate_derivatives (const bool compute_dIdW, const bool compute_dIdX, const bool compute_d2I) |
| Allocate and setup the derivative vectors/matrices. More... | |
| void | allocate_dIdX (dealii::LinearAlgebra::distributed::Vector< real > &dIdX) const |
| Allocate and setup the derivative dIdX vector. More... | |
| void | set_derivatives (const bool compute_dIdW, const bool compute_dIdX, const bool compute_d2I, const Sacado::Fad::DFad< Sacado::Fad::DFad< real >> volume_local_sum, std::vector< dealii::types::global_dof_index > cell_soln_dofs_indices, std::vector< dealii::types::global_dof_index > cell_metric_dofs_indices) |
| Set the derivative vectors/matrices. More... | |
| void | need_compute (bool &compute_value, bool &compute_dIdW, bool &compute_dIdX, bool &compute_d2I) |
| Checks which derivatives actually need to be recomputed. More... | |
| real2 | evaluate_volume_cell_functional (const Physics::PhysicsBase< dim, nstate, real2 > &physics, const std::vector< real2 > &soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< real2 > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const dealii::Quadrature< dim > &volume_quadrature) const |
| Templated function to evaluate a cell's volume functional. | |
| virtual real | evaluate_volume_cell_functional (const Physics::PhysicsBase< dim, nstate, real > &physics, const std::vector< real > &soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< real > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const dealii::Quadrature< dim > &volume_quadrature) const |
| Corresponding real function to evaluate a cell's volume functional. | |
| virtual Sacado::Fad::DFad< Sacado::Fad::DFad< real > > | evaluate_volume_cell_functional (const Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< Sacado::Fad::DFad< real >>> &physics_fad_fad, const std::vector< Sacado::Fad::DFad< Sacado::Fad::DFad< real >> > &soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< Sacado::Fad::DFad< Sacado::Fad::DFad< real >> > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const dealii::Quadrature< dim > &volume_quadrature) const |
| Corresponding FadFadType function to evaluate a cell's volume functional. | |
| real2 | evaluate_boundary_cell_functional (const Physics::PhysicsBase< dim, nstate, real2 > &physics, const unsigned int boundary_id, const std::vector< real2 > &soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< real2 > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const unsigned int face_number, const dealii::Quadrature< dim-1 > &face_quadrature) const |
| Templated function to evaluate a cell's boundary functional. More... | |
| virtual real | evaluate_boundary_cell_functional (const Physics::PhysicsBase< dim, nstate, real > &physics, const unsigned int boundary_id, const std::vector< real > &soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< real > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const unsigned int face_number, const dealii::Quadrature< dim-1 > &face_quadrature) const |
| Corresponding real function to evaluate a cell's boundary functional. | |
| virtual Sacado::Fad::DFad< Sacado::Fad::DFad< real > > | evaluate_boundary_cell_functional (const Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< Sacado::Fad::DFad< real >>> &physics_fad_fad, const unsigned int boundary_id, const std::vector< Sacado::Fad::DFad< Sacado::Fad::DFad< real >> > &soln_coeff, const dealii::FESystem< dim > &fe_solution, const std::vector< Sacado::Fad::DFad< Sacado::Fad::DFad< real >> > &coords_coeff, const dealii::FESystem< dim > &fe_metric, const unsigned int face_number, const dealii::Quadrature< dim-1 > &face_quadrature) const |
| Corresponding FadFadType function to evaluate a cell's boundary functional. | |
| virtual real | evaluate_volume_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &, const dealii::Point< dim, real > &, const std::array< real, nstate > &, const std::array< dealii::Tensor< 1, dim, real >, nstate > &) const |
| Virtual function for computation of cell volume functional term. More... | |
| virtual FadFadType | evaluate_volume_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, FadFadType > &, const dealii::Point< dim, FadFadType > &, const std::array< FadFadType, nstate > &, const std::array< dealii::Tensor< 1, dim, FadFadType >, nstate > &) const |
| Virtual function for Sacado computation of cell volume functional term and derivatives. More... | |
| virtual real | evaluate_boundary_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &, const unsigned int, const dealii::Point< dim, real > &, const dealii::Tensor< 1, dim, real > &, const std::array< real, nstate > &, const std::array< dealii::Tensor< 1, dim, real >, nstate > &) const |
| Virtual function for computation of cell boundary functional term. More... | |
| virtual FadFadType | evaluate_boundary_integrand (const PHiLiP::Physics::PhysicsBase< dim, nstate, FadFadType > &, const unsigned int, const dealii::Point< dim, FadFadType > &, const dealii::Tensor< 1, dim, FadFadType > &, const std::array< FadFadType, nstate > &, const std::array< dealii::Tensor< 1, dim, FadFadType >, nstate > &) const |
| Virtual function for Sacado computation of cell boundary functional term and derivatives. More... | |
Protected Attributes inherited from PHiLiP::TargetFunctional< dim, nstate, real > | |
| const dealii::LinearAlgebra::distributed::Vector< real > | target_solution |
| Solution used to evaluate target functional. | |
Protected Attributes inherited from PHiLiP::Functional< dim, nstate, real > | |
| std::shared_ptr< Physics::PhysicsBase< dim, nstate, FadFadType > > | physics_fad_fad |
| Physics that should correspond to the one in DGBase. | |
| const dealii::UpdateFlags | volume_update_flags |
| Update flags needed at volume points. | |
| const dealii::UpdateFlags | face_update_flags |
| Update flags needed at face points. | |
| const bool | uses_solution_values |
| Will evaluate solution values at quadrature points. | |
| const bool | uses_solution_gradient |
| Will evaluate solution gradient at quadrature points. | |
| dealii::ConditionalOStream | pcout |
| Parallel std::cout that only outputs on mpi_rank==0. | |
Target boundary values. Simply zero out the default volume contribution.
Definition at line 9 of file target_boundary_functional.h.