3 #include <Epetra_RowMatrixTransposer.h>     4 #include <deal.II/distributed/shared_tria.h>     5 #include <deal.II/distributed/solution_transfer.h>     6 #include <deal.II/distributed/tria.h>     7 #include <deal.II/dofs/dof_tools.h>     8 #include <deal.II/fe/fe_q.h>     9 #include <deal.II/fe/fe_values.h>    10 #include <deal.II/grid/tria.h>    11 #include <deal.II/lac/la_parallel_vector.h>    12 #include <deal.II/numerics/data_out.h>    13 #include <deal.II/numerics/vector_tools.h>    19 #include "dg/dg_base.hpp"    20 #include "functional.h"    21 #include "linear_solver/linear_solver.h"    22 #include "parameters/all_parameters.h"    23 #include "physics/physics.h"    24 #include "post_processor/physics_post_processor.h"    29 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
    35     functional(_functional),
    37     triangulation(dg->triangulation),
    38     solution_coarse(dg->solution),
    40     mpi_communicator(MPI_COMM_WORLD),
    41     pcout(
std::cout, dealii::Utilities::MPI::this_mpi_process(mpi_communicator)==0)
    47     for(
auto cell = 
dg->dof_handler.begin_active(); cell != 
dg->dof_handler.end(); ++cell)
    48         if(cell->is_locally_owned())
    52 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
    64     for(
auto cell = 
dg->dof_handler.begin_active(); cell != 
dg->dof_handler.end(); ++cell)
    65         if(cell->is_locally_owned())
    69     dIdw_fine      = dealii::LinearAlgebra::distributed::Vector<real>();
    70     dIdw_coarse    = dealii::LinearAlgebra::distributed::Vector<real>();
    71     adjoint_fine   = dealii::LinearAlgebra::distributed::Vector<real>();
    72     adjoint_coarse = dealii::LinearAlgebra::distributed::Vector<real>();
    77 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
    85     if(
adjoint_state == AdjointStateEnum::coarse && state == AdjointStateEnum::fine) 
    88     if(
adjoint_state == AdjointStateEnum::fine && state == AdjointStateEnum::coarse)
    92 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
    95     dealii::IndexSet locally_owned_dofs, locally_relevant_dofs;
    96     locally_owned_dofs =  
dg->dof_handler.locally_owned_dofs();
    97     dealii::DoFTools::extract_locally_relevant_dofs(
dg->dof_handler, locally_relevant_dofs);
   103     using VectorType       = 
typename dealii::LinearAlgebra::distributed::Vector<double>;
   104     using DoFHandlerType   = 
typename dealii::DoFHandler<dim>;
   107     SolutionTransfer solution_transfer(
dg->dof_handler);
   108     solution_transfer.prepare_for_coarsening_and_refinement(
solution_coarse);
   110     dg->high_order_grid->prepare_for_coarsening_and_refinement();
   111     dg->triangulation->prepare_coarsening_and_refinement();
   113     for (
auto cell = 
dg->dof_handler.begin_active(); cell != 
dg->dof_handler.end(); ++cell)
   114         if (cell->is_locally_owned()) 
   115             cell->set_future_fe_index(cell->active_fe_index()+1);
   117     dg->triangulation->execute_coarsening_and_refinement();
   118     dg->high_order_grid->execute_coarsening_and_refinement();
   120     dg->allocate_system();
   121     dg->solution.zero_out_ghosts();
   123     if constexpr (std::is_same_v<
typename dealii::SolutionTransfer<dim,VectorType,DoFHandlerType>, 
   124                                  decltype(solution_transfer)>){
   127         solution_transfer.interpolate(
dg->solution);
   130     dg->solution.update_ghost_values();
   135 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
   138     dg->high_order_grid->prepare_for_coarsening_and_refinement();
   139     dg->triangulation->prepare_coarsening_and_refinement();
   141     for (
auto cell = 
dg->dof_handler.begin_active(); cell != 
dg->dof_handler.end(); ++cell)
   142         if (cell->is_locally_owned()) 
   145     dg->triangulation->execute_coarsening_and_refinement();
   146     dg->high_order_grid->execute_coarsening_and_refinement();
   148     dg->allocate_system();
   149     dg->solution.zero_out_ghosts();
   156 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
   163     const bool compute_dIdW = 
true, compute_dIdX = 
false;
   164     const real functional_value = 
functional->evaluate_functional(compute_dIdW,compute_dIdX);
   165     (void) functional_value;
   170     dg->assemble_residual(
true);
   171     dg->system_matrix *= -1.0;
   173     dealii::TrilinosWrappers::SparseMatrix system_matrix_transpose;
   174     Epetra_CrsMatrix *system_matrix_transpose_tril;
   176     Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&
dg->system_matrix.trilinos_matrix()));
   177     epmt.CreateTranspose(
false, system_matrix_transpose_tril);
   178     system_matrix_transpose.reinit(*system_matrix_transpose_tril,
true);
   179     delete system_matrix_transpose_tril;
   186 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
   193     const bool compute_dIdW = 
true, compute_dIdX = 
false;
   194     const real functional_value = 
functional->evaluate_functional(compute_dIdW,compute_dIdX);
   195     (void) functional_value;
   200     dg->assemble_residual(
true);
   201     dg->system_matrix *= -1.0;
   203     dealii::TrilinosWrappers::SparseMatrix system_matrix_transpose;
   204     Epetra_CrsMatrix *system_matrix_transpose_tril;
   206     Epetra_RowMatrixTransposer epmt(const_cast<Epetra_CrsMatrix *>(&
dg->system_matrix.trilinos_matrix()));
   207     epmt.CreateTranspose(
false, system_matrix_transpose_tril);
   208     system_matrix_transpose.reinit(*system_matrix_transpose_tril);
   215 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
   223     const unsigned int max_dofs_per_cell = 
dg->dof_handler.get_fe_collection().max_dofs_per_cell();
   224     std::vector<dealii::types::global_dof_index> current_dofs_indices(max_dofs_per_cell);
   227     for(
auto cell = 
dg->dof_handler.begin_active(); cell != 
dg->dof_handler.end(); ++cell){
   228         if(!cell->is_locally_owned()) 
continue;
   230         const unsigned int fe_index_curr_cell = cell->active_fe_index();
   231         const dealii::FESystem<dim,dim> ¤t_fe_ref = 
dg->fe_collection[fe_index_curr_cell];
   232         const unsigned int n_dofs_curr_cell = current_fe_ref.n_dofs_per_cell();
   234         current_dofs_indices.resize(n_dofs_curr_cell);
   235         cell->get_dof_indices(current_dofs_indices);
   238         for(
unsigned int idof = 0; idof < n_dofs_curr_cell; ++idof){
   239             dwr_cell += 
dg->right_hand_side[current_dofs_indices[idof]]*
adjoint_fine[current_dofs_indices[idof]];
   248 template <
int dim, 
int nstate, 
typename real, 
typename MeshType>
   251     dealii::DataOut<dim, dealii::DoFHandler<dim>> data_out;
   252     data_out.attach_dof_handler(
dg->dof_handler);
   255     data_out.add_data_vector(
dg->solution, *post_processor);
   257     dealii::Vector<float> subdomain(
dg->triangulation->n_active_cells());
   258     for (
unsigned int i = 0; i < subdomain.size(); ++i) {
   259         subdomain(i) = 
dg->triangulation->locally_owned_subdomain();
   261     data_out.add_data_vector(subdomain, 
"subdomain", dealii::DataOut_DoFData<dealii::DoFHandler<dim>,dim>::DataVectorType::type_cell_data);
   264     std::vector<unsigned int> active_fe_indices;
   265     dg->dof_handler.get_active_fe_indices(active_fe_indices);
   266     dealii::Vector<double> active_fe_indices_dealiivector(active_fe_indices.begin(), active_fe_indices.end());
   267     dealii::Vector<double> cell_poly_degree = active_fe_indices_dealiivector;
   269     data_out.add_data_vector(active_fe_indices_dealiivector, 
"PolynomialDegree", dealii::DataOut_DoFData<dealii::DoFHandler<dim>,dim>::DataVectorType::type_cell_data);
   271     std::vector<std::string> residual_names;
   272     for(
int s=0;s<nstate;++s) {
   273         std::string varname = 
"residual" + dealii::Utilities::int_to_string(s,1);
   274         residual_names.push_back(varname);
   277     data_out.add_data_vector(
dg->right_hand_side, residual_names, dealii::DataOut_DoFData<dealii::DoFHandler<dim>,dim>::DataVectorType::type_dof_data);
   280     std::vector<std::string> dIdw_names;
   281     for(
int s=0;s<nstate;++s) {
   282         std::string varname = 
"dIdw" + dealii::Utilities::int_to_string(s,1);
   283         dIdw_names.push_back(varname);
   286     std::vector<std::string> adjoint_names;
   287     for(
int s=0;s<nstate;++s) {
   288         std::string varname = 
"psi" + dealii::Utilities::int_to_string(s,1);
   289         adjoint_names.push_back(varname);
   294         data_out.add_data_vector(
dIdw_fine, dIdw_names, dealii::DataOut_DoFData<dealii::DoFHandler<dim>,dim>::DataVectorType::type_dof_data);
   295         data_out.add_data_vector(
adjoint_fine, adjoint_names, dealii::DataOut_DoFData<dealii::DoFHandler<dim>,dim>::DataVectorType::type_dof_data);
   297         data_out.add_data_vector(
dual_weighted_residual_fine, 
"DWR", dealii::DataOut_DoFData<dealii::DoFHandler<dim>,dim>::DataVectorType::type_cell_data);
   299         data_out.add_data_vector(
dIdw_coarse, dIdw_names, dealii::DataOut_DoFData<dealii::DoFHandler<dim>,dim>::DataVectorType::type_dof_data);
   300         data_out.add_data_vector(
adjoint_coarse, adjoint_names, dealii::DataOut_DoFData<dealii::DoFHandler<dim>,dim>::DataVectorType::type_dof_data);
   303     const int iproc = dealii::Utilities::MPI::this_mpi_process(
mpi_communicator);
   305     data_out.build_patches();
   308     std::string filename = 
"adjoint-" ;
   310     filename += dealii::Utilities::int_to_string(cycle, 4) + 
".";
   311     filename += dealii::Utilities::int_to_string(iproc, 4);
   313     std::ofstream output(filename);
   314     data_out.write_vtu(output);
   317         std::vector<std::string> filenames;
   318         for (
unsigned int iproc = 0; iproc < dealii::Utilities::MPI::n_mpi_processes(
mpi_communicator); ++iproc) {
   319             std::string fn = 
"adjoint-";
   324             fn += dealii::Utilities::int_to_string(dim, 1) + 
"D-";
   325             fn += dealii::Utilities::int_to_string(cycle, 4) + 
".";
   326             fn += dealii::Utilities::int_to_string(iproc, 4);
   328             filenames.push_back(fn);
   330         std::string master_fn = 
"adjoint-";
   332             master_fn += 
"fine-";
   334             master_fn += 
"coarse-";
   335         master_fn += dealii::Utilities::int_to_string(dim, 1) +
"D-";
   336         master_fn += dealii::Utilities::int_to_string(cycle, 4) + 
".pvtu";
   337         std::ofstream master_output(master_fn);
   338         data_out.write_pvtu_record(master_output, filenames);
 dealii::LinearAlgebra::distributed::Vector< real > dIdw_coarse
functional derivative (on the coarse grid) 
std::shared_ptr< Functional< dim, nstate, real, MeshType > > functional
Functional class pointer. 
void convert_to_state(AdjointStateEnum state)
Converts the adjoint to specified state. 
AdjointStateEnum adjoint_state
Current adjoint state. 
dealii::LinearAlgebra::distributed::Vector< real > fine_grid_adjoint()
Computes the fine grid adjoint. 
void fine_to_coarse()
Return the problem to the original solution and polynomial distribution. 
dealii::LinearAlgebra::distributed::Vector< real > adjoint_coarse
coarse grid adjoint ( ) 
Base class from which Advection, Diffusion, ConvectionDiffusion, and Euler is derived. 
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)
Adjoint(std::shared_ptr< DGBase< dim, real, MeshType > > _dg, std::shared_ptr< Functional< dim, nstate, real, MeshType > > _functional, std::shared_ptr< Physics::PhysicsBase< dim, nstate, Sacado::Fad::DFad< real >> > _physics)
Constructor. 
Files for the baseline physics. 
AdjointStateEnum
For storing the current state in the adjoint. 
dealii::Vector< real > dual_weighted_residual()
compute the Dual Weighted Residual (DWR) 
dealii::LinearAlgebra::distributed::Vector< real > coarse_grid_adjoint()
Computes the coarse grid adjoint. 
dealii::LinearAlgebra::distributed::Vector< real > dIdw_fine
functional derivative (on the fine grid) 
void reinit()
Reinitialize Adjoint with the same pointers. 
dealii::Vector< real > coarse_fe_index
Original FE_index distribution. 
MPI_Comm mpi_communicator
MPI communicator. 
static std::unique_ptr< dealii::DataPostprocessor< dim > > create_Postprocessor(const Parameters::AllParameters *const parameters_input)
Create the post-processor with the correct template parameters. 
dealii::LinearAlgebra::distributed::Vector< real > solution_coarse
original solution 
std::shared_ptr< DGBase< dim, real, MeshType > > dg
DG class pointer. 
void output_results_vtk(const unsigned int cycle)
Outputs the current solution and adjoint values. 
dealii::Vector< real > dual_weighted_residual_fine
dual weighted residual 
void coarse_to_fine()
Projects the problem to a p-enriched space. 
DGBase is independent of the number of state variables. 
dealii::LinearAlgebra::distributed::Vector< real > adjoint_fine
fine grid adjoint ( )