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> 40 #include <deal.II/lac/read_write_vector.h> 42 const double STEPSIZE = 1e-7;
43 const double TOLERANCE = 1e-6;
46 using Triangulation = dealii::Triangulation<PHILIP_DIM>;
48 using Triangulation = dealii::parallel::distributed::Triangulation<PHILIP_DIM>;
51 std::pair<unsigned int, double>
53 dealii::TrilinosWrappers::BlockSparseMatrix &matrix_A,
54 dealii::LinearAlgebra::distributed::BlockVector<double> &right_hand_side,
55 dealii::LinearAlgebra::distributed::BlockVector<double> &solution,
60 using trilinos_vector_type = dealii::TrilinosWrappers::MPI::Vector;
61 using vector_type = dealii::LinearAlgebra::distributed::Vector<double>;
64 vector_type &_rhs1 = right_hand_side.block(0);
65 vector_type &_rhs2 = right_hand_side.block(1);
67 dealii::IndexSet rhs1_locally_owned = _rhs1.locally_owned_elements();
68 dealii::IndexSet rhs1_ghost = _rhs1.get_partitioner()->ghost_indices();
70 dealii::IndexSet rhs2_locally_owned = _rhs2.locally_owned_elements();
71 dealii::IndexSet rhs2_ghost = _rhs2.get_partitioner()->ghost_indices();
74 trilinos_vector_type rhs1(rhs1_locally_owned);
75 trilinos_vector_type rhs2(rhs2_locally_owned);
76 rhs1_locally_owned.print(std::cout);
77 std::cout << rhs1.size() << std::endl;
80 trilinos_vector_type soln1(_rhs1.locally_owned_elements());
81 trilinos_vector_type soln2(_rhs2.locally_owned_elements());
84 dealii::LinearAlgebra::ReadWriteVector<double> rw_vector(rhs1_locally_owned);
85 rw_vector.import(_rhs1, dealii::VectorOperation::insert);
86 rhs1.import(rw_vector, dealii::VectorOperation::insert);
89 dealii::LinearAlgebra::ReadWriteVector<double> rw_vector(rhs2_locally_owned);
90 rw_vector.import(_rhs2, dealii::VectorOperation::insert);
91 rhs2.import(rw_vector, dealii::VectorOperation::insert);
97 using matrix_type = dealii::TrilinosWrappers::SparseMatrix;
98 using payload_type = dealii::TrilinosWrappers::internal::LinearOperatorImplementation::TrilinosPayload;
101 const auto &L11 = matrix_A.block(0,0);
102 const auto &L12 = matrix_A.block(0,1);
103 const auto &L21 = matrix_A.block(1,0);
104 const auto &L22 = matrix_A.block(1,1);
106 const auto op_L11 = dealii::linear_operator<trilinos_vector_type,trilinos_vector_type,payload_type>(L11);
107 const auto op_L12 = dealii::linear_operator<trilinos_vector_type,trilinos_vector_type,payload_type>(L12);
108 const auto op_L21 = dealii::linear_operator<trilinos_vector_type,trilinos_vector_type,payload_type>(L21);
109 const auto op_L22 = dealii::linear_operator<trilinos_vector_type,trilinos_vector_type,payload_type>(L22);
111 dealii::ReductionControl reduction_control_L11(2000, 1.0e-15, 1.0e-12);
112 dealii::SolverGMRES<trilinos_vector_type> solver_L11(reduction_control_L11);
114 dealii::TrilinosWrappers::PreconditionILU preconditioner_L11;
116 preconditioner_L11.initialize(L11);
119 const auto op_L11_inv = dealii::inverse_operator(op_L11, solver_L11, preconditioner_L11);
120 const auto op_Schur = op_L22 - op_L21 * op_L11_inv * op_L12;
122 const auto op_preconditioner_L11 = dealii::linear_operator<trilinos_vector_type,trilinos_vector_type,payload_type>(L11,preconditioner_L11);
123 const auto op_approxSchur = op_L22 - op_L21 * op_preconditioner_L11 * op_L12;
125 const trilinos_vector_type schur_rhs = rhs2 - op_L21 * op_L11_inv * rhs1;
126 std::cout << reduction_control_L11.last_step() <<
" GMRES iterations to solve L11inv*rhs1." << std::endl;
134 dealii::TrilinosWrappers::SparseMatrix approxSchur;
135 const auto L11_rows = L11.locally_owned_range_indices();
136 trilinos_vector_type L11_diag_inv(L11_rows);
137 for (
auto row = L11_rows.begin(); row != L11_rows.end(); ++row) {
138 L11_diag_inv[*row] = 1.0/L11.diag_element(*row);
140 L21.mmult(approxSchur, L12, L11_diag_inv);
142 approxSchur.add(1.0,L22);
143 dealii::TrilinosWrappers::PreconditionILU preconditioner_Schur;
144 const unsigned int ilu_fill = 20, overlap = 1;
145 const double ilu_atol = 1e-5, ilu_rtol = 1e-2;
146 preconditioner_Schur.initialize(approxSchur, dealii::TrilinosWrappers::PreconditionILU::AdditionalData(ilu_fill,ilu_atol,ilu_rtol,overlap));
149 dealii::SolverControl solver_control_Schur(2000, 1.e-15,
true);
151 dealii::SolverGMRES<trilinos_vector_type> solver_Schur(solver_control_Schur,
152 dealii::SolverGMRES<trilinos_vector_type>::AdditionalData (2000,
false,
true,
false) );
153 const auto op_Schur_inv = dealii::inverse_operator(op_Schur, solver_Schur, preconditioner_Schur);
155 soln2 = op_Schur_inv * schur_rhs;
156 std::cout << solver_control_Schur.last_step() <<
" GMRES iterations to obtain convergence." << std::endl;
157 soln1 = op_L11_inv * (rhs1 - op_L12 * soln2);
160 dealii::LinearAlgebra::ReadWriteVector<double> rw_vector;
161 rw_vector.reinit(soln1);
162 solution.block(0).import(rw_vector, dealii::VectorOperation::insert);
165 dealii::LinearAlgebra::ReadWriteVector<double> rw_vector;
166 rw_vector.reinit(soln2);
167 solution.block(1).import(rw_vector, dealii::VectorOperation::insert);
205 const trilinos_vector_type r1 = op_L11 * soln1 + op_L12 * soln2 - rhs1;
206 const trilinos_vector_type r2 = op_L21 * soln1 + op_L22 * soln2 - rhs2;
208 std::cout<<
"r1 norm: "<<r1.l2_norm()<<std::endl;
210 std::cout<<
"r2 norm: "<<r2.l2_norm()<<std::endl;
234 template <
int dim,
int nstate,
typename real>
249 template <
typename real2>
252 const dealii::Point<dim,real2> &phys_coord,
253 const std::array<real2,nstate> &soln_at_q,
254 const std::array<dealii::Tensor<1,dim,real2>,nstate> &)
const 258 for (
int istate=0; istate<nstate; ++istate) {
260 l2error += std::pow(soln_at_q[istate] - uexact, 2);
268 template<
typename real2>
272 const dealii::Point<dim,real2> &phys_coord,
273 const dealii::Tensor<1,dim,real2> &,
274 const std::array<real2,nstate> &soln_at_q,
275 const std::array<dealii::Tensor<1,dim,real2>,nstate> &)
const 278 for (
int istate=0; istate<nstate; ++istate) {
280 l1error += std::abs(soln_at_q[istate] - uexact);
289 const unsigned int boundary_id,
290 const dealii::Point<dim,real> &phys_coord,
291 const dealii::Tensor<1,dim,real> &normal,
292 const std::array<real,nstate> &soln_at_q,
293 const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_at_q)
const override 295 return evaluate_boundary_integrand<real>(
307 const unsigned int boundary_id,
308 const dealii::Point<dim,FadFadType> &phys_coord,
309 const dealii::Tensor<1,dim,FadFadType> &normal,
310 const std::array<FadFadType,nstate> &soln_at_q,
311 const std::array<dealii::Tensor<1,dim,FadFadType>,nstate> &soln_grad_at_q)
const override 313 return evaluate_boundary_integrand<FadFadType>(
326 const dealii::Point<dim,real> &phys_coord,
327 const std::array<real,nstate> &soln_at_q,
328 const std::array<dealii::Tensor<1,dim,real>,nstate> &soln_grad_at_q)
const override 330 return evaluate_volume_integrand<>(physics, phys_coord, soln_at_q, soln_grad_at_q);
335 const dealii::Point<dim,FadFadType> &phys_coord,
336 const std::array<FadFadType,nstate> &soln_at_q,
337 const std::array<dealii::Tensor<1,dim,FadFadType>,nstate> &soln_grad_at_q)
const override 339 return evaluate_volume_integrand<>(physics, phys_coord, soln_at_q, soln_grad_at_q);
345 template <
int dim,
int nstate>
348 dealii::LinearAlgebra::distributed::Vector<double> solution_no_ghost;
354 int main(
int argc,
char *argv[])
357 const int dim = PHILIP_DIM;
358 const int nstate = 1;
359 int fail_bool =
false;
362 dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
363 const int this_mpi_process = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
364 dealii::ConditionalOStream
pcout(std::cout, this_mpi_process==0);
367 dealii::ParameterHandler parameter_handler;
373 const unsigned poly_degree = 1;
376 std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
380 typename dealii::Triangulation<dim>::MeshSmoothing(
381 dealii::Triangulation<dim>::smoothing_on_refinement |
382 dealii::Triangulation<dim>::smoothing_on_coarsening));
384 const unsigned int n_refinements = 2;
387 const bool colorize =
true;
389 dealii::GridGenerator::hyper_cube(*grid, left, right, colorize);
390 grid->refine_global(n_refinements);
391 const double random_factor = 0.2;
392 const bool keep_boundary =
false;
393 if (random_factor > 0.0) dealii::GridTools::distort_random (random_factor, *grid, keep_boundary);
395 pcout <<
"Grid generated and refined" << std::endl;
399 pcout <<
"dg created" << std::endl;
401 dg->allocate_system();
402 pcout <<
"dg allocated" << std::endl;
404 const int n_refine = 2;
405 for (
int i=0; i<n_refine;i++) {
406 dg->high_order_grid->prepare_for_coarsening_and_refinement();
407 grid->prepare_coarsening_and_refinement();
408 unsigned int icell = 0;
409 for (
auto cell = grid->begin_active(); cell!=grid->end(); ++cell) {
411 if (!cell->is_locally_owned())
continue;
412 if (icell < grid->n_global_active_cells()/2) {
413 cell->set_refine_flag();
416 grid->execute_coarsening_and_refinement();
417 bool mesh_out = (i==n_refine-1);
418 dg->high_order_grid->execute_coarsening_and_refinement(mesh_out);
420 dg->allocate_system ();
424 pcout <<
"Physics created" << std::endl;
427 initialize_perturbed_solution(*dg, *physics_double);
428 pcout <<
"solution initialized" << std::endl;
431 pcout << std::endl <<
"Starting Hessian AD... " << std::endl;
433 const bool compute_dIdW =
false, compute_dIdX =
false, compute_d2I =
true;
434 double functional_value = functional.
evaluate_functional(compute_dIdW, compute_dIdX, compute_d2I);
435 (void) functional_value;
438 bool compute_dRdW, compute_dRdX, compute_d2R;
440 pcout <<
"Evaluating RHS only to use as dual variables..." << std::endl;
441 compute_dRdW =
false; compute_dRdX =
false, compute_d2R =
false;
442 dg->assemble_residual(compute_dRdW, compute_dRdX, compute_d2R);
443 dealii::LinearAlgebra::distributed::Vector<double> dummy_dual(dg->right_hand_side);
444 dg->set_dual(dummy_dual);
446 pcout <<
"Evaluating RHS with d2R..." << std::endl;
447 compute_dRdW =
true; compute_dRdX =
false, compute_d2R =
false;
448 dg->assemble_residual(compute_dRdW, compute_dRdX, compute_d2R);
449 compute_dRdW =
false; compute_dRdX =
true, compute_d2R =
false;
450 dg->assemble_residual(compute_dRdW, compute_dRdX, compute_d2R);
451 compute_dRdW =
false; compute_dRdX =
false, compute_d2R =
true;
452 dg->assemble_residual(compute_dRdW, compute_dRdX, compute_d2R);
453 dealii::LinearAlgebra::distributed::Vector<double> rhs_d2R(dg->right_hand_side);
456 dealii::TrilinosWrappers::SparseMatrix dRdW_transpose;
458 Epetra_CrsMatrix *transpose_CrsMatrix;
459 Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&dg->system_matrix.trilinos_matrix()));
460 epmt.CreateTranspose(
false, transpose_CrsMatrix);
461 dRdW_transpose.reinit(*transpose_CrsMatrix);
464 dealii::TrilinosWrappers::SparseMatrix dRdX_transpose;
466 Epetra_CrsMatrix *transpose_CrsMatrix;
467 Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&dg->dRdXv.trilinos_matrix()));
468 epmt.CreateTranspose(
false, transpose_CrsMatrix);
469 dRdX_transpose.reinit(*transpose_CrsMatrix);
472 dealii::TrilinosWrappers::SparseMatrix d2RdXdW;
474 Epetra_CrsMatrix *transpose_CrsMatrix;
475 Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&dg->d2RdWdX.trilinos_matrix()));
476 epmt.CreateTranspose(
false, transpose_CrsMatrix);
477 d2RdXdW.reinit(*transpose_CrsMatrix);
479 dealii::TrilinosWrappers::SparseMatrix d2IdXdW;
481 Epetra_CrsMatrix *transpose_CrsMatrix;
482 Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&functional.
d2IdWdX->trilinos_matrix()));
483 epmt.CreateTranspose(
false, transpose_CrsMatrix);
484 d2IdXdW.reinit(*transpose_CrsMatrix);
488 functional.
d2IdWdW->add(1.0,dg->d2RdWdW);
489 functional.
d2IdWdX->add(1.0,dg->d2RdWdX);
490 d2IdXdW.add(1.0,d2RdXdW);
491 functional.
d2IdXdX->add(1.0,dg->d2RdXdX);
514 dealii::TrilinosWrappers::BlockSparseMatrix kkt_hessian;
515 kkt_hessian.reinit(3,3);
516 kkt_hessian.block(0, 0).copy_from( *functional.
d2IdWdW);
517 kkt_hessian.block(0, 1).copy_from( *functional.
d2IdWdX);
518 kkt_hessian.block(0, 2).copy_from( dRdW_transpose);
520 kkt_hessian.block(1, 0).copy_from( d2IdXdW);
521 kkt_hessian.block(1, 1).copy_from( *functional.
d2IdXdX);
522 kkt_hessian.block(1, 2).copy_from( dRdX_transpose);
524 kkt_hessian.block(2, 0).copy_from( dg->system_matrix);
525 kkt_hessian.block(2, 1).copy_from( dg->dRdXv);
526 dealii::TrilinosWrappers::SparsityPattern zero_sparsity_pattern(dg->locally_owned_dofs, MPI_COMM_WORLD, 0);
532 zero_sparsity_pattern.compress();
533 kkt_hessian.block(2, 2).reinit(zero_sparsity_pattern);
535 kkt_hessian.collect_sizes();
537 pcout <<
"kkt_hessian.frobenius_norm() " << kkt_hessian.frobenius_norm() << std::endl;
547 dealii::LinearAlgebra::distributed::BlockVector<double> block_vector(3);
548 block_vector.block(0) = dg->solution;
549 block_vector.block(1) = dg->high_order_grid->volume_nodes;
550 block_vector.block(2) = dummy_dual;
551 dealii::LinearAlgebra::distributed::BlockVector<double> Hv(3);
552 dealii::LinearAlgebra::distributed::BlockVector<double> Htv(3);
553 Hv.reinit(block_vector);
554 Htv.reinit(block_vector);
556 kkt_hessian.vmult(Hv,block_vector);
557 kkt_hessian.Tvmult(Htv,block_vector);
561 const double vector_norm = Hv.l2_norm();
562 const double vector_abs_diff = Htv.l2_norm();
563 const double vector_rel_diff = vector_abs_diff / vector_norm;
565 const double tol = 1e-11;
567 <<
" vector_abs_diff: " << vector_abs_diff
568 <<
" vector_rel_diff: " << vector_rel_diff
570 <<
" vector_abs_diff: " << vector_abs_diff
571 <<
" vector_rel_diff: " << vector_rel_diff
573 if (vector_abs_diff > tol && vector_rel_diff > tol) fail_bool =
true;
575 const int n_mpi_processes = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
576 if (n_mpi_processes == 1) {
577 dealii::FullMatrix<double> fullA(kkt_hessian.m());
578 fullA.copy_from(kkt_hessian);
579 const int n_digits = 8;
580 if (pcout.is_active()) fullA.print_formatted(pcout.get_stream(), n_digits,
true, n_digits+7,
"0", 1., 0.);
583 dealii::deallog.depth_console(3);
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.