1 #include <Epetra_RowMatrixTransposer.h>     3 #include <deal.II/base/conditional_ostream.h>     5 #include <deal.II/distributed/tria.h>     6 #include <deal.II/grid/grid_generator.h>     7 #include <deal.II/grid/grid_tools.h>     9 #include <deal.II/numerics/vector_tools.h>    11 #include "physics/physics_factory.h"    12 #include "parameters/all_parameters.h"    13 #include "dg/dg_factory.hpp"    14 #include "functional/functional.h"    16 #include <deal.II/lac/full_matrix.h>    18 #include <deal.II/lac/solver_bicgstab.h>    19 #include <deal.II/lac/solver_cg.h>    20 #include <deal.II/lac/solver_gmres.h>    21 #include <deal.II/lac/solver_minres.h>    23 #include <deal.II/lac/block_sparsity_pattern.h>    25 #include <deal.II/lac/precondition.h>    27 #include <deal.II/lac/trilinos_block_sparse_matrix.h>    28 #include <deal.II/lac/trilinos_precondition.h>    29 #include <deal.II/lac/trilinos_solver.h>    30 #include <deal.II/lac/trilinos_sparsity_pattern.h>    31 #include <deal.II/lac/trilinos_sparse_matrix.h>    33 #include <deal.II/lac/trilinos_parallel_block_vector.h>    34 #include <deal.II/lac/la_parallel_block_vector.h>    37 #include <deal.II/lac/packaged_operation.h>    38 #include <deal.II/lac/trilinos_linear_operator.h>    41 #include <deal.II/lac/solver_control.h>    44 #include <deal.II/lac/read_write_vector.h>    46 const double STEPSIZE = 1e-7;
    47 const double TOLERANCE = 1e-6;
    50     using Triangulation = dealii::Triangulation<PHILIP_DIM>;
    52     using Triangulation = dealii::parallel::distributed::Triangulation<PHILIP_DIM>;
    55 std::pair<unsigned int, double>
    57     dealii::TrilinosWrappers::BlockSparseMatrix &matrix_A,
    58     dealii::LinearAlgebra::distributed::BlockVector<double> &right_hand_side,
    59     dealii::LinearAlgebra::distributed::BlockVector<double> &solution,
    64     using block_vec = dealii::LinearAlgebra::distributed::BlockVector<double>;
    65     const auto kkt_operator = dealii::block_operator<block_vec> (matrix_A);
    66     dealii::SolverControl solver_control(100000, 1.e-15, 
true, 
true);
    67     const unsigned int  max_n_tmp_vectors = 2000;
    68     const bool  right_preconditioning = 
false;
    69     const bool  use_default_residual = 
true;
    70     const bool  force_re_orthogonalization = 
false;
    71     dealii::SolverGMRES<block_vec>::AdditionalData add_data( max_n_tmp_vectors, right_preconditioning, use_default_residual, force_re_orthogonalization);
    72     dealii::SolverGMRES<block_vec> solver_gmres(solver_control, add_data);
    73     auto kkt_inv = inverse_operator(kkt_operator, solver_gmres, dealii::PreconditionIdentity());
    74     solution = kkt_inv * right_hand_side;
    80 template <
int dim, 
int nstate, 
typename real>
    96     template <
typename real2>
    99                 const dealii::Point<dim,real2> &phys_coord,
   100                 const std::array<real2,nstate> &soln_at_q,
   101                 const std::array<dealii::Tensor<1,dim,real2>,nstate> &)
 const   104         for (
int istate=0; istate<nstate; ++istate) {
   106             l2error += std::pow(soln_at_q[istate] - uexact, 2);
   113     template<
typename real2>
   117         const dealii::Point<dim,real2> &phys_coord,
   118         const dealii::Tensor<1,dim,real2> &,
   119         const std::array<real2,nstate> &soln_at_q,
   120         const std::array<dealii::Tensor<1,dim,real2>,nstate> &)
 const   123         for (
int istate=0; istate<nstate; ++istate) {
   125             l1error += std::abs(soln_at_q[istate] - uexact);
   134         const unsigned int boundary_id,
   135         const dealii::Point<dim,real> &phys_coord,
   136         const dealii::Tensor<1,dim,real> &normal,
   137         const std::array<real,nstate> &soln_at_q,
   138         const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_at_q)
 const override   140         return evaluate_boundary_integrand<real>(
   152         const unsigned int boundary_id,
   153         const dealii::Point<dim,FadFadType> &phys_coord,
   154         const dealii::Tensor<1,dim,FadFadType> &normal,
   155         const std::array<FadFadType,nstate> &soln_at_q,
   156         const std::array<dealii::Tensor<1,dim,FadFadType>,nstate> &soln_grad_at_q)
 const override   158         return evaluate_boundary_integrand<FadFadType>(
   172             const dealii::Point<dim,real> &phys_coord,
   173             const std::array<real,nstate> &soln_at_q,
   174             const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_at_q)
 const override   176         return evaluate_volume_integrand<>(physics, phys_coord, soln_at_q, soln_grad_at_q);
   181             const dealii::Point<dim,FadFadType> &phys_coord,
   182             const std::array<FadFadType,nstate> &soln_at_q,
   183             const std::array<dealii::Tensor<1,dim,FadFadType>,nstate> &soln_grad_at_q)
 const override   185         return evaluate_volume_integrand<>(physics, phys_coord, soln_at_q, soln_grad_at_q);
   191 template <
int dim, 
int nstate>
   194     dealii::LinearAlgebra::distributed::Vector<double> solution_no_ghost;
   200 int main(
int argc, 
char *argv[])
   203     const int dim = PHILIP_DIM;
   204     const int nstate = 1;
   205     int fail_bool = 
false;
   208     dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
   209     const int this_mpi_process = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
   210     dealii::ConditionalOStream 
pcout(std::cout, this_mpi_process==0);
   213     dealii::ParameterHandler parameter_handler;
   219     const unsigned poly_degree = 1;
   222     std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
   226         typename dealii::Triangulation<dim>::MeshSmoothing(
   227             dealii::Triangulation<dim>::smoothing_on_refinement |
   228             dealii::Triangulation<dim>::smoothing_on_coarsening));
   230     const unsigned int n_refinements = 2;
   233     const bool colorize = 
true;
   235     dealii::GridGenerator::hyper_cube(*grid, left, right, colorize);
   236     grid->refine_global(n_refinements);
   237     const double random_factor = 0.2;
   238     const bool keep_boundary = 
false;
   239     if (random_factor > 0.0) dealii::GridTools::distort_random (random_factor, *grid, keep_boundary);
   241     pcout << 
"Grid generated and refined" << std::endl;
   245     pcout << 
"dg created" << std::endl;
   247     dg->allocate_system();
   248     pcout << 
"dg allocated" << std::endl;
   250     const int n_refine = 2;
   251     for (
int i=0; i<n_refine;i++) {
   252         dg->high_order_grid->prepare_for_coarsening_and_refinement();
   253         grid->prepare_coarsening_and_refinement();
   254         unsigned int icell = 0;
   255         for (
auto cell = grid->begin_active(); cell!=grid->end(); ++cell) {
   257             if (!cell->is_locally_owned()) 
continue;
   258             if (icell < grid->n_global_active_cells()/2) {
   259                 cell->set_refine_flag();
   262         grid->execute_coarsening_and_refinement();
   263         bool mesh_out = (i==n_refine-1);
   264         dg->high_order_grid->execute_coarsening_and_refinement(mesh_out);
   266     dg->allocate_system ();
   270     pcout << 
"Physics created" << std::endl;
   273     initialize_perturbed_solution(*dg, *physics_double);
   274     pcout << 
"solution initialized" << std::endl;
   277     pcout << std::endl << 
"Starting Hessian AD... " << std::endl;
   279     const bool compute_dIdW = 
false, compute_dIdX = 
false, compute_d2I = 
true;
   280     double functional_value = functional.
evaluate_functional(compute_dIdW, compute_dIdX, compute_d2I);
   281     (void) functional_value;
   284     bool compute_dRdW, compute_dRdX, compute_d2R;
   286     pcout << 
"Evaluating RHS only to use as dual variables..." << std::endl;
   287     compute_dRdW = 
false; compute_dRdX = 
false, compute_d2R = 
false;
   288     dg->assemble_residual(compute_dRdW, compute_dRdX, compute_d2R);
   289     dealii::LinearAlgebra::distributed::Vector<double> dummy_dual(dg->right_hand_side);
   290     dg->set_dual(dummy_dual);
   292     pcout << 
"Evaluating RHS with d2R..." << std::endl;
   293     compute_dRdW = 
true; compute_dRdX = 
false, compute_d2R = 
false;
   294     dg->assemble_residual(compute_dRdW, compute_dRdX, compute_d2R);
   295     compute_dRdW = 
false; compute_dRdX = 
true, compute_d2R = 
false;
   296     dg->assemble_residual(compute_dRdW, compute_dRdX, compute_d2R);
   297     compute_dRdW = 
false; compute_dRdX = 
false, compute_d2R = 
true;
   298     dg->assemble_residual(compute_dRdW, compute_dRdX, compute_d2R);
   299     dealii::LinearAlgebra::distributed::Vector<double> rhs_d2R(dg->right_hand_side);
   302     dealii::TrilinosWrappers::SparseMatrix dRdW_transpose;
   304         Epetra_CrsMatrix *transpose_CrsMatrix;
   305         Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&dg->system_matrix.trilinos_matrix()));
   306         epmt.CreateTranspose(
false, transpose_CrsMatrix);
   307         dRdW_transpose.reinit(*transpose_CrsMatrix);
   310     dealii::TrilinosWrappers::SparseMatrix dRdX_transpose;
   312         Epetra_CrsMatrix *transpose_CrsMatrix;
   313         Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&dg->dRdXv.trilinos_matrix()));
   314         epmt.CreateTranspose(
false, transpose_CrsMatrix);
   315         dRdX_transpose.reinit(*transpose_CrsMatrix);
   318     dealii::TrilinosWrappers::SparseMatrix d2RdXdW;
   320         Epetra_CrsMatrix *transpose_CrsMatrix;
   321         Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&dg->d2RdWdX.trilinos_matrix()));
   322         epmt.CreateTranspose(
false, transpose_CrsMatrix);
   323         d2RdXdW.reinit(*transpose_CrsMatrix);
   325     dealii::TrilinosWrappers::SparseMatrix d2IdXdW;
   327         Epetra_CrsMatrix *transpose_CrsMatrix;
   328         Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&functional.
d2IdWdX->trilinos_matrix()));
   329         epmt.CreateTranspose(
false, transpose_CrsMatrix);
   330         d2IdXdW.reinit(*transpose_CrsMatrix);
   334     functional.
d2IdWdW->add(1.0,dg->d2RdWdW);
   335     functional.
d2IdWdX->add(1.0,dg->d2RdWdX);
   336     d2IdXdW.add(1.0,d2RdXdW);
   337     functional.
d2IdXdX->add(1.0,dg->d2RdXdX);
   360     dealii::TrilinosWrappers::BlockSparseMatrix kkt_hessian;
   361     kkt_hessian.reinit(3,3);
   362     kkt_hessian.block(0, 0).copy_from( *functional.
d2IdWdW);
   363     kkt_hessian.block(0, 1).copy_from( *functional.
d2IdWdX);
   364     kkt_hessian.block(0, 2).copy_from( dRdW_transpose);
   366     kkt_hessian.block(1, 0).copy_from( d2IdXdW);
   367     kkt_hessian.block(1, 1).copy_from( *functional.
d2IdXdX);
   368     kkt_hessian.block(1, 2).copy_from( dRdX_transpose);
   370     kkt_hessian.block(2, 0).copy_from( dg->system_matrix);
   371     kkt_hessian.block(2, 1).copy_from( dg->dRdXv);
   372     dealii::TrilinosWrappers::SparsityPattern zero_sparsity_pattern(dg->locally_owned_dofs, MPI_COMM_WORLD, 0);
   378     zero_sparsity_pattern.compress();
   379     kkt_hessian.block(2, 2).reinit(zero_sparsity_pattern);
   381     kkt_hessian.collect_sizes();
   383     pcout << 
"kkt_hessian.frobenius_norm()  " << kkt_hessian.frobenius_norm() << std::endl;
   393     dealii::LinearAlgebra::distributed::BlockVector<double> block_vector(3);
   394     block_vector.block(0) = dg->solution;
   395     block_vector.block(1) = dg->high_order_grid->volume_nodes;
   396     block_vector.block(2) = dummy_dual;
   397     dealii::LinearAlgebra::distributed::BlockVector<double> Hv(3);
   398     dealii::LinearAlgebra::distributed::BlockVector<double> Htv(3);
   399     Hv.reinit(block_vector);
   400     Htv.reinit(block_vector);
   402     kkt_hessian.vmult(Hv,block_vector);
   403     kkt_hessian.Tvmult(Htv,block_vector);
   407     const double vector_norm = Hv.l2_norm();
   408     const double vector_abs_diff = Htv.l2_norm();
   409     const double vector_rel_diff = vector_abs_diff / vector_norm;
   411     const double tol = 1e-11;
   413                     << 
" vector_abs_diff: " << vector_abs_diff
   414                     << 
" vector_rel_diff: " << vector_rel_diff
   416                     << 
" vector_abs_diff: " << vector_abs_diff
   417                     << 
" vector_rel_diff: " << vector_rel_diff
   419     if (vector_abs_diff > tol && vector_rel_diff > tol) fail_bool = 
true;
   421     const int n_mpi_processes = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
   422     if (n_mpi_processes == 1) {
   423         dealii::FullMatrix<double> fullA(kkt_hessian.m());
   424         fullA.copy_from(kkt_hessian);
   425         const int n_digits = 8;
   426         if (pcout.is_active()) fullA.print_formatted(pcout.get_stream(), n_digits, 
true, n_digits+7, 
"0", 1., 0.);
   429     dealii::deallog.depth_console(3);
   436     dealii::LinearAlgebra::distributed::BlockVector<double> diff(3);
   437     diff.reinit(block_vector);
   438     kkt_hessian.vmult(diff,Htv);
   440     double diff_norm = diff.l2_norm();
 std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > d2IdXdX
Sparse matrix for storing the functional partial second derivatives. 
const bool uses_solution_gradient
Will evaluate solution gradient at quadrature points. 
Parameters related to the linear solver. 
const bool uses_solution_values
Will evaluate solution values at quadrature points. 
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived. 
real2 evaluate_volume_integrand(const PHiLiP::Physics::PhysicsBase< dim, nstate, real2 > &physics, const dealii::Point< dim, real2 > &phys_coord, const std::array< real2, nstate > &soln_at_q, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &) const
Templated volume integrand. 
std::pair< unsigned int, double > solve_linear(const dealii::TrilinosWrappers::SparseMatrix &system_matrix, dealii::LinearAlgebra::distributed::Vector< double > &right_hand_side, dealii::LinearAlgebra::distributed::Vector< double > &solution, const Parameters::LinearSolverParam ¶m)
virtual FadFadType evaluate_boundary_integrand(const PHiLiP::Physics::PhysicsBase< dim, nstate, FadFadType > &physics, const unsigned int boundary_id, const dealii::Point< dim, FadFadType > &phys_coord, const dealii::Tensor< 1, dim, FadFadType > &normal, const std::array< FadFadType, nstate > &soln_at_q, const std::array< dealii::Tensor< 1, dim, FadFadType >, nstate > &soln_grad_at_q) const override
Virtual function for Sacado computation of cell boundary functional term and derivatives. 
std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > d2IdWdW
Sparse matrix for storing the functional partial second derivatives. 
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 
Files for the baseline physics. 
Main parameter class that contains the various other sub-parameter classes. 
dealii::LinearAlgebra::distributed::Vector< double > solution
Current modal coefficients of the solution. 
dealii::IndexSet locally_owned_dofs
Locally own degrees of freedom. 
L2_Norm_Functional(std::shared_ptr< PHiLiP::DGBase< dim, real >> dg_input, const bool uses_solution_values=true, const bool uses_solution_gradient=false)
Constructor. 
Sacado::Fad::DFad< real > FadType
Sacado AD type for first derivatives. 
static std::shared_ptr< DGBase< dim, real, MeshType > > create_discontinuous_galerkin(const Parameters::AllParameters *const parameters_input, const unsigned int degree, const unsigned int max_degree_input, const unsigned int grid_degree_input, const std::shared_ptr< Triangulation > triangulation_input)
Creates a derived object DG, but returns it as DGBase. 
LinearSolverParam linear_solver_param
Contains parameters for linear solver. 
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 std::shared_ptr< PhysicsBase< dim, nstate, real > > create_Physics(const Parameters::AllParameters *const parameters_input, std::shared_ptr< ModelBase< dim, nstate, real > > model_input=nullptr)
Factory to return the correct physics given input file. 
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)
void parse_parameters(dealii::ParameterHandler &prm)
Retrieve parameters from dealii::ParameterHandler. 
dealii::DoFHandler< dim > dof_handler
Finite Element Collection to represent the high-order grid. 
Sacado::Fad::DFad< FadType > FadFadType
Sacado AD type that allows 2nd derivatives. 
std::shared_ptr< dealii::TrilinosWrappers::SparseMatrix > d2IdWdX
Sparse matrix for storing the functional partial second derivatives. 
static void declare_parameters(dealii::ParameterHandler &prm)
Declare parameters that can be set as inputs and set up the default options. 
virtual real evaluate_boundary_integrand(const PHiLiP::Physics::PhysicsBase< dim, nstate, real > &physics, const unsigned int boundary_id, const dealii::Point< dim, real > &phys_coord, const dealii::Tensor< 1, dim, real > &normal, const std::array< real, nstate > &soln_at_q, const std::array< dealii::Tensor< 1, dim, real >, nstate > &soln_grad_at_q) const override
Virtual function for computation of cell boundary functional term. 
std::shared_ptr< DGBase< dim, real, dealii::parallel::distributed::Triangulation< dim > > > dg
Smart pointer to DGBase. 
DGBase is independent of the number of state variables. 
std::shared_ptr< ManufacturedSolutionFunction< dim, real > > manufactured_solution_function
Manufactured solution function. 
virtual real evaluate_functional(const bool compute_dIdW=false, const bool compute_dIdX=false, const bool compute_d2I=false)
Evaluates the functional derivative with respect to the solution variable. 
real2 evaluate_boundary_integrand(const PHiLiP::Physics::PhysicsBase< dim, nstate, real2 > &physics, const unsigned int, const dealii::Point< dim, real2 > &phys_coord, const dealii::Tensor< 1, dim, real2 > &, const std::array< real2, nstate > &soln_at_q, const std::array< dealii::Tensor< 1, dim, real2 >, nstate > &) const
Virtual function for computation of cell boundary functional term. 
dealii::ConditionalOStream pcout
Parallel std::cout that only outputs on mpi_rank==0.