[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
target_boundary_functional.h
1 #include "target_functional.h"
2 
3 namespace PHiLiP {
4 
8 template <int dim, int nstate, typename real>
9 class TargetBoundaryFunctional : public TargetFunctional<dim, nstate, real>
10 {
11  using FadType = Sacado::Fad::DFad<real>;
12  using FadFadType = Sacado::Fad::DFad<FadType>;
13 
15 
20 
21 public:
24  std::shared_ptr<DGBase<dim,real>> dg_input,
25  const dealii::LinearAlgebra::distributed::Vector<real> &target_solution,
26  const bool uses_solution_values = true,
27  const bool uses_solution_gradient = false)
28  : TargetFunctional<dim,nstate,real>(dg_input, target_solution, uses_solution_values, uses_solution_gradient)
29  {}
30 
32  template <typename real2>
35  const dealii::Point<dim,real2> &/*phys_coord*/,
36  const std::array<real2,nstate> &,//soln_at_q,
37  const std::array<real,nstate> &,//target_soln_at_q,
38  const std::array<dealii::Tensor<1,dim,real2>,nstate> &/*soln_grad_at_q*/,
39  const std::array<dealii::Tensor<1,dim,real2>,nstate> &/*target_soln_grad_at_q*/) const
40  {
41  real2 l2error = 0;
42 
43  return l2error;
44  }
45 
49  const dealii::Point<dim,real> &phys_coord,
50  const std::array<real,nstate> &soln_at_q,
51  const std::array<real,nstate> &target_soln_at_q,
52  const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_at_q,
53  const std::array<dealii::Tensor<1,dim,real>,nstate> &target_soln_grad_at_q) const override
54  {
55  return evaluate_volume_integrand<>(physics, phys_coord, soln_at_q, target_soln_at_q, soln_grad_at_q, target_soln_grad_at_q);
56  }
60  const dealii::Point<dim,FadFadType> &phys_coord,
61  const std::array<FadFadType,nstate> &soln_at_q,
62  const std::array<real,nstate> &target_soln_at_q,
63  const std::array<dealii::Tensor<1,dim,FadFadType>,nstate> &soln_grad_at_q,
64  const std::array<dealii::Tensor<1,dim,FadFadType>,nstate> &target_soln_grad_at_q) const override
65  {
66  return evaluate_volume_integrand<>(physics, phys_coord, soln_at_q, target_soln_at_q, soln_grad_at_q, target_soln_grad_at_q);
67  }
68 };
69 
70 } // PHiLiP namespace
Sacado::Fad::DFad< real > FadType
Sacado AD type for first derivatives.
const bool uses_solution_gradient
Will evaluate solution gradient at quadrature points.
Definition: functional.h:315
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
TargetFunctional base class.
const bool uses_solution_values
Will evaluate solution values at quadrature points.
Definition: functional.h:314
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived.
Definition: physics.h:34
Files for the baseline physics.
Definition: ADTypes.hpp:10
Sacado::Fad::DFad< FadType > FadFadType
Sacado AD type that allows 2nd derivatives.
const dealii::LinearAlgebra::distributed::Vector< real > target_solution
Solution used to evaluate target functional.
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
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.
Functional base class.
Definition: functional.h:43
DGBase is independent of the number of state variables.
Definition: dg_base.hpp:82
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.