1 #ifndef __DIFFUSION_EXACT_ADJOINT_H__     2 #define __DIFFUSION_EXACT_ADJOINT_H__     4 #include <deal.II/base/tensor.h>     5 #include <deal.II/base/types.h>     9 #include "dg/dg_base.hpp"    10 #include "functional/functional.h"    11 #include "parameters/all_parameters.h"    12 #include "physics/convection_diffusion.h"    13 #include "physics/manufactured_solution.h"    14 #include "physics/physics.h"    53 template <
int dim, 
typename real>
    58     using dealii::Function<dim,real>::value;
    59     using dealii::Function<dim,real>::gradient;
    60     using dealii::Function<dim,real>::hessian;
    63     real 
value (
const dealii::Point<dim,real> &pos, 
const unsigned int istate = 0) 
const override;
    66     dealii::Tensor<1,dim,real> 
gradient (
const dealii::Point<dim,real> &pos, 
const unsigned int istate = 0) 
const override;
    69     dealii::SymmetricTensor<2,dim,real> 
hessian (
const dealii::Point<dim,real> &, 
const unsigned int )
 const     71         return dealii::SymmetricTensor<2,dim,real>();
    76 template <
int dim, 
typename real>
    80     using dealii::Function<dim,real>::value;
    81     using dealii::Function<dim,real>::gradient;
    82     using dealii::Function<dim,real>::hessian;
    85     real 
value (
const dealii::Point<dim,real> &pos, 
const unsigned int istate = 0) 
const override;
    88     dealii::Tensor<1,dim,real> 
gradient (
const dealii::Point<dim,real> &pos, 
const unsigned int istate = 0) 
const override;
    91     dealii::SymmetricTensor<2,dim,real> 
hessian (
const dealii::Point<dim,real> &, 
const unsigned int )
 const override    93         return dealii::SymmetricTensor<2,dim,real>();
    98 template <
int dim, 
int nstate, 
typename real>
   115         const bool                                              convection, 
   116         const bool                                              diffusion,
   118             Physics::ConvectionDiffusion<dim,
nstate,real>::ConvectionDiffusion(
   122                 default_diffusion_tensor(),
   123                 Parameters::ManufacturedSolutionParam::get_default_advection_vector(),
   124                 default_diffusion_coefficient(),
   125                 manufactured_solution_function)
   137         dealii::Tensor<2,3,double> diffusion_tensor;
   138         diffusion_tensor[0][0] = 1;
   139         if constexpr(dim>=2) {
   140             diffusion_tensor[0][1] = 0.0;
   141             diffusion_tensor[1][0] = 0.0;
   142             diffusion_tensor[1][1] = 1.0;
   144         if constexpr(dim>=3) {
   145             diffusion_tensor[0][2] = 0.0;
   146             diffusion_tensor[1][2] = 0.0;
   147             diffusion_tensor[2][0] = 0.0;
   148             diffusion_tensor[2][1] = 0.0;
   149             diffusion_tensor[2][2] = 1.0;
   151         return diffusion_tensor;
   155     virtual real objective_function(
   156         const dealii::Point<dim,real> &pos) 
const = 0;
   160 template <
int dim, 
int nstate, 
typename real>
   177         const bool convection, 
   178         const bool diffusion): 
   187     std::array<real,nstate> source_term (
   188         const dealii::Point<dim,real> &pos,
   189         const std::array<real,nstate> &solution,
   190         const real current_time,
   191         const dealii::types::global_dof_index cell_index) 
const override;
   194     real objective_function(
   195         const dealii::Point<dim,real> &pos) 
const override;
   199 template <
int dim, 
int nstate, 
typename real>
   216         const bool convection, 
   217         const bool diffusion): 
   226     std::array<real,nstate> source_term (
   227         const dealii::Point<dim,real> &pos,
   228         const std::array<real,nstate> &solution,
   229         const real current_time,
   230         const dealii::types::global_dof_index cell_index) 
const override;
   233     real objective_function(
   234         const dealii::Point<dim,real> &pos) 
const override;
   238 template <
int dim, 
int nstate, 
typename real>
   248         const bool uses_solution_values = 
true,
   249         const bool uses_solution_gradient = 
false)
   250     : 
PHiLiP::
Functional<dim,
nstate,real>(dg_input,_physics_fad_fad,uses_solution_values,uses_solution_gradient)
   254     template <
typename real2>
   255     real2 evaluate_volume_integrand(
   257         const dealii::Point<dim,real2> &phys_coord,
   258         const std::array<real2,nstate> &soln_at_q,
   259         const std::array<dealii::Tensor<1,dim,real2>,
nstate> &soln_grad_at_q) 
const;
   264         const dealii::Point<dim,real> &phys_coord,
   265         const std::array<real,nstate> &soln_at_q,
   266         const std::array<dealii::Tensor<1,dim,real>,
nstate> &soln_grad_at_q)
 const override   268         return evaluate_volume_integrand<>(physics, phys_coord, soln_at_q, soln_grad_at_q);
   274         const dealii::Point<dim,FadFadType> &phys_coord,
   275         const std::array<FadFadType,nstate> &soln_at_q,
   276         const std::array<dealii::Tensor<1,dim,FadFadType>,
nstate> &soln_grad_at_q)
 const override   278         return evaluate_volume_integrand<>(physics, phys_coord, soln_at_q, soln_grad_at_q);
   283 template <
int dim, 
int nstate>
   294     int run_test() 
const;
   298 double eval_avg_slope(std::vector<double> error, std::vector<double> grid_size, 
unsigned int n_grids);
   303 #endif //__DIFFUSION_EXACT_ADJOINT_H__ Functional that performs the inner product over the entire domain. 
Sacado::Fad::DFad< real > FadType
Sacado AD type for first derivatives. 
Sacado::Fad::DFad< FadType > FadFadType
Sacado AD type that allows 2nd derivatives. 
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived. 
Manufactured solution used for grid studies to check convergence orders. 
Convection-diffusion with linear advective and diffusive term. Derived from PhysicsBase. 
Files for the baseline physics. 
physics for the v variable 
diffusion_u(const Parameters::AllParameters *const parameters_input, const bool convection, const bool diffusion)
constructor 
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< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_at_q) const override
Non-template functions to override the template classes. 
static dealii::Tensor< 2, 3, double > default_diffusion_tensor()
Default diffusion tensor is set for isotropic diffusion, . 
Main parameter class that contains the various other sub-parameter classes. 
diffusion_objective(const Parameters::AllParameters *const parameters_input, const bool convection, const bool diffusion, std::shared_ptr< ManufacturedSolutionFunction< dim, real >> manufactured_solution_function)
constructor 
physics for the u variable 
manufactured solution for v 
dealii::SymmetricTensor< 2, dim, real > hessian(const dealii::Point< dim, real > &, const unsigned int) const override
Hessian of manufactured solution is unused but needed to make class concrete. 
dealii::SymmetricTensor< 2, dim, real > hessian(const dealii::Point< dim, real > &, const unsigned int) const
Hessian of manufactured solution is unused but needed to make class concrete. 
diffusion_v(const Parameters::AllParameters *const parameters_input, const bool convection, const bool diffusion)
constructor 
manufactured solution for u 
const unsigned int nstate
DiffusionFunctional(std::shared_ptr< PHiLiP::DGBase< dim, real >> dg_input, std::shared_ptr< PHiLiP::Physics::PhysicsBase< dim, nstate, FadFadType >> _physics_fad_fad, const bool uses_solution_values=true, const bool uses_solution_gradient=false)
Constructor. 
dealii::Tensor< 1, dim, real > gradient(const dealii::Point< dim, real > &pos, const unsigned int istate=0) const override
Gradient of the manufactured solution. 
DGBase is independent of the number of state variables. 
real value(const dealii::Point< dim, real > &pos, const unsigned int istate=0) const override
overriding the function for the value and gradient 
Base class of all the tests. 
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< dealii::Tensor< 1, dim, FadFadType >, nstate > &soln_grad_at_q) const override
Non-template functions to override the template classes. 
static double default_diffusion_coefficient()
negative one is used for the diffusion coefficient so that the problem becomes  
parent class to add the objective function directly to physics as a virtual class ...