[P]arallel [Hi]gh-order [Li]brary for [P]DEs
Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
|
DGBase is independent of the number of state variables. More...
#include <dg_base.hpp>
Public Types | |
using | Triangulation = MeshType |
using | MassiveCollectionTuple = std::tuple< dealii::hp::FECollection< dim >, dealii::hp::QCollection< dim >, dealii::hp::QCollection< dim-1 >, dealii::hp::FECollection< dim >, dealii::hp::FECollection< 1 >, dealii::hp::FECollection< 1 >, dealii::hp::FECollection< 1 >, dealii::hp::QCollection< 1 > > |
Makes for cleaner doxygen documentation. | |
Public Member Functions | |
virtual | ~DGBase ()=default |
Destructor. | |
DGBase (const int nstate_input, 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) | |
Principal constructor that will call delegated constructor. More... | |
void | reinit () |
Reinitializes the DG object after a change of triangulation. More... | |
DGBase (const int nstate_input, 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, const MassiveCollectionTuple collection_tuple) | |
Delegated constructor that initializes collections. More... | |
void | set_high_order_grid (std::shared_ptr< HighOrderGrid< dim, real, MeshType >> new_high_order_grid) |
Sets the associated high order grid with the provided one. | |
void | set_all_cells_fe_degree (const unsigned int degree) |
Refers to a collection Mappings, which represents the high-order grid. More... | |
unsigned int | get_max_fe_degree () |
Gets the maximum value of currently active FE degree. | |
unsigned int | get_min_fe_degree () |
Gets the minimum value of currently active FE degree. | |
dealii::Point< dim > | coordinates_of_highest_refined_cell (bool check_for_p_refined_cell=false) |
Returns the coordinates of the most refined cell. | |
virtual void | allocate_system (const bool compute_dRdW=true, const bool compute_dRdX=true, const bool compute_d2R=true) |
Allocates the system. More... | |
void | time_scale_solution_update (dealii::LinearAlgebra::distributed::Vector< double > &solution_update, const real CFL) const |
Scales a solution update with the appropriate maximum time step. More... | |
void | time_scaled_mass_matrices (const real scale) |
void | reinit_operators_for_cell_residual_loop (const unsigned int poly_degree_int, const unsigned int poly_degree_ext, const unsigned int grid_degree, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis_int, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis_ext, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis_int, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis_ext, OPERATOR::local_basis_stiffness< dim, 2 *dim, real > &flux_basis_stiffness, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_int, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_ext, OPERATOR::mapping_shape_functions< dim, 2 *dim, real > &mapping_basis) |
Builds needed operators for cell residual loop. | |
void | reinit_operators_for_mass_matrix (const bool Cartesian_element, const unsigned int poly_degree, const unsigned int grid_degree, OPERATOR::mapping_shape_functions< dim, 2 *dim, real > &mapping_basis, OPERATOR::basis_functions< dim, 2 *dim, real > &basis, OPERATOR::local_mass< dim, 2 *dim, real > &reference_mass_matrix, OPERATOR::local_Flux_Reconstruction_operator< dim, 2 *dim, real > &reference_FR, OPERATOR::local_Flux_Reconstruction_operator_aux< dim, 2 *dim, real > &reference_FR_aux, OPERATOR::derivative_p< dim, 2 *dim, real > &deriv_p) |
Builds needed operators to compute mass matrices/inverses efficiently. | |
void | evaluate_mass_matrices (bool do_inverse_mass_matrix=false) |
Allocates and evaluates the mass matrices for the entire grid. | |
void | evaluate_local_metric_dependent_mass_matrix_and_set_in_global_mass_matrix (const bool Cartesian_element, const bool do_inverse_mass_matrix, const unsigned int poly_degree, const unsigned int curr_grid_degree, const unsigned int n_quad_pts, const unsigned int n_dofs_cell, const std::vector< dealii::types::global_dof_index > dofs_indices, OPERATOR::metric_operators< real, dim, 2 *dim > &metric_oper, OPERATOR::basis_functions< dim, 2 *dim, real > &basis, OPERATOR::local_mass< dim, 2 *dim, real > &reference_mass_matrix, OPERATOR::local_Flux_Reconstruction_operator< dim, 2 *dim, real > &reference_FR, OPERATOR::local_Flux_Reconstruction_operator_aux< dim, 2 *dim, real > &reference_FR_aux, OPERATOR::derivative_p< dim, 2 *dim, real > &deriv_p) |
Evaluates the metric dependent local mass matrices and inverses, then sets them in the global matrices. | |
void | apply_inverse_global_mass_matrix (const dealii::LinearAlgebra::distributed::Vector< double > &input_vector, dealii::LinearAlgebra::distributed::Vector< double > &output_vector, const bool use_auxiliary_eq=false) |
Applies the inverse of the local metric dependent mass matrices when the global is not stored. More... | |
void | apply_global_mass_matrix (const dealii::LinearAlgebra::distributed::Vector< double > &input_vector, dealii::LinearAlgebra::distributed::Vector< double > &output_vector, const bool use_auxiliary_eq=false, const bool use_unmodified_mass_matrix=false) |
Applies the local metric dependent mass matrices when the global is not stored. More... | |
std::vector< real > | evaluate_time_steps (const bool exact_time_stepping) |
Evaluates the maximum stable time step. More... | |
void | add_mass_matrices (const real scale) |
Add mass matrices to the system scaled by a factor (likely time-step) More... | |
void | add_time_scaled_mass_matrices () |
Add time scaled mass matrices to the system. More... | |
double | get_residual_l2norm () const |
Returns the L2-norm of the right_hand_side vector. | |
double | get_residual_linfnorm () const |
Returns the Linf-norm of the right_hand_side vector. | |
unsigned int | n_dofs () const |
Number of degrees of freedom. | |
void | set_anisotropic_flags () |
Set anisotropic flags based on jump indicator. More... | |
template<typename real2 > | |
real2 | discontinuity_sensor (const dealii::Quadrature< dim > &volume_quadrature, const std::vector< real2 > &soln_coeff_high, const dealii::FiniteElement< dim, dim > &fe_high, const std::vector< real2 > &jac_det) |
void | set_dual (const dealii::LinearAlgebra::distributed::Vector< real > &dual_input) |
Sets the stored dual variables used to compute the dual dotted with the residual Hessians. | |
dealii::SparsityPattern | get_dRdX_sparsity_pattern () |
Evaluate SparsityPattern of dRdX. | |
dealii::SparsityPattern | get_dRdW_sparsity_pattern () |
Evaluate SparsityPattern of dRdW. | |
dealii::SparsityPattern | get_d2RdWdW_sparsity_pattern () |
Evaluate SparsityPattern of the residual Hessian dual.d2RdWdW. | |
dealii::SparsityPattern | get_d2RdXdX_sparsity_pattern () |
Evaluate SparsityPattern of the residual Hessian dual.d2RdXdX. | |
dealii::SparsityPattern | get_d2RdWdX_sparsity_pattern () |
Evaluate SparsityPattern of the residual Hessian dual.d2RdXdW. | |
dealii::SparsityPattern | get_dRdXs_sparsity_pattern () |
Evaluate SparsityPattern of dRdXs. | |
dealii::SparsityPattern | get_d2RdXsdXs_sparsity_pattern () |
Evaluate SparsityPattern of the residual Hessian dual.d2RdXsdXs. | |
dealii::SparsityPattern | get_d2RdWdXs_sparsity_pattern () |
Evaluate SparsityPattern of the residual Hessian dual.d2RdXsdW. | |
dealii::TrilinosWrappers::SparseMatrix | get_dRdX_finite_differences (dealii::SparsityPattern dRdX_sparsity_pattern) |
Evaluate dRdX using finite-differences. | |
void | initialize_manufactured_solution () |
Virtual function defined in DG. | |
void | output_results_vtk (const unsigned int cycle, const double current_time=0.0) |
Output solution. | |
void | output_face_results_vtk (const unsigned int cycle, const double current_time=0.0) |
Output Euler face solution. | |
void | assemble_residual (const bool compute_dRdW=false, const bool compute_dRdX=false, const bool compute_d2R=false, const double CFL_mass=0.0) |
Main loop of the DG class. More... | |
template<typename DoFCellAccessorType1 , typename DoFCellAccessorType2 > | |
void | assemble_cell_residual (const DoFCellAccessorType1 ¤t_cell, const DoFCellAccessorType2 ¤t_metric_cell, const bool compute_dRdW, const bool compute_dRdX, const bool compute_d2R, dealii::hp::FEValues< dim, dim > &fe_values_collection_volume, dealii::hp::FEFaceValues< dim, dim > &fe_values_collection_face_int, dealii::hp::FEFaceValues< dim, dim > &fe_values_collection_face_ext, dealii::hp::FESubfaceValues< dim, dim > &fe_values_collection_subface, dealii::hp::FEValues< dim, dim > &fe_values_collection_volume_lagrange, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis_int, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis_ext, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis_int, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis_ext, OPERATOR::local_basis_stiffness< dim, 2 *dim, real > &flux_basis_stiffness, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_int, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_ext, OPERATOR::mapping_shape_functions< dim, 2 *dim, real > &mapping_basis, const bool compute_auxiliary_right_hand_side, dealii::LinearAlgebra::distributed::Vector< double > &rhs, std::array< dealii::LinearAlgebra::distributed::Vector< double >, dim > &rhs_aux) |
Used in assemble_residual(). More... | |
void | set_current_time (const real current_time_input) |
Sets the current time within DG to be used for unsteady source terms. | |
void | allocate_auxiliary_equation () |
Allocates the auxiliary equations' variables and right hand side (primarily for Strong form diffusive) | |
virtual void | assemble_auxiliary_residual ()=0 |
Asembles the auxiliary equations' residuals and solves. | |
virtual void | allocate_dual_vector ()=0 |
Allocate the dual vector for optimization. More... | |
void | update_artificial_dissipation_discontinuity_sensor () |
Update discontinuity sensor. More... | |
virtual void | allocate_model_variables ()=0 |
Allocate the necessary variables declared in src/physics/model.h. | |
virtual void | update_model_variables ()=0 |
Update the necessary variables declared in src/physics/model.h. | |
virtual void | set_use_auxiliary_eq ()=0 |
Set use_auxiliary_eq flag. | |
Public Attributes | |
const Parameters::AllParameters *const | all_parameters |
Pointer to all parameters. | |
const int | nstate |
Number of state variables. More... | |
const unsigned int | initial_degree |
Initial polynomial degree assigned during constructor. | |
const unsigned int | max_degree |
Maximum degree used for p-refi1nement. More... | |
const unsigned int | max_grid_degree |
Maximum grid degree used for hp-refi1nement. More... | |
std::shared_ptr< Triangulation > | triangulation |
Mesh. | |
dealii::SparsityPattern | sparsity_pattern |
Sparsity pattern used on the system_matrix. More... | |
dealii::SparsityPattern | mass_sparsity_pattern |
Sparsity pattern used on the system_matrix. More... | |
dealii::TrilinosWrappers::SparseMatrix | time_scaled_global_mass_matrix |
Global mass matrix divided by the time scales. More... | |
dealii::TrilinosWrappers::SparseMatrix | global_mass_matrix |
Global mass matrix. More... | |
dealii::TrilinosWrappers::SparseMatrix | global_inverse_mass_matrix |
Global inverser mass matrix. More... | |
dealii::TrilinosWrappers::SparseMatrix | global_mass_matrix_auxiliary |
Global auxiliary mass matrix. More... | |
dealii::TrilinosWrappers::SparseMatrix | global_inverse_mass_matrix_auxiliary |
Global inverse of the auxiliary mass matrix. | |
dealii::TrilinosWrappers::SparseMatrix | system_matrix |
dealii::TrilinosWrappers::SparseMatrix | system_matrix_transpose |
std::unique_ptr< Epetra_RowMatrixTransposer > | epetra_rowmatrixtransposer_dRdW |
Epetra_RowMatrixTransposer used to transpose the system_matrix. | |
dealii::TrilinosWrappers::SparseMatrix | dRdXv |
dealii::TrilinosWrappers::SparseMatrix | d2RdWdW |
dealii::TrilinosWrappers::SparseMatrix | d2RdXdX |
dealii::TrilinosWrappers::SparseMatrix | d2RdWdX |
dealii::LinearAlgebra::distributed::Vector< double > | right_hand_side |
Residual of the current solution. More... | |
dealii::IndexSet | locally_owned_dofs |
Locally own degrees of freedom. | |
dealii::IndexSet | ghost_dofs |
Locally relevant ghost degrees of freedom. | |
dealii::IndexSet | locally_relevant_dofs |
Union of locally owned degrees of freedom and relevant ghost degrees of freedom. | |
dealii::IndexSet | locally_owned_dofs_grid |
Locally own degrees of freedom for the grid. | |
dealii::IndexSet | ghost_dofs_grid |
Locally relevant ghost degrees of freedom for the grid. | |
dealii::IndexSet | locally_relevant_dofs_grid |
dealii::LinearAlgebra::distributed::Vector< double > | solution |
Current modal coefficients of the solution. More... | |
std::array< dealii::LinearAlgebra::distributed::Vector< double >, dim > | auxiliary_right_hand_side |
The auxiliary equations' right hand sides. | |
std::array< dealii::LinearAlgebra::distributed::Vector< double >, dim > | auxiliary_solution |
The auxiliary equations' solution. | |
dealii::Vector< double > | cell_volume |
Time it takes for the maximum wavespeed to cross the cell domain. More... | |
dealii::Vector< double > | max_dt_cell |
Time it takes for the maximum wavespeed to cross the cell domain. More... | |
dealii::Vector< double > | reduced_mesh_weights |
dealii::Vector< double > | artificial_dissipation_coeffs |
Artificial dissipation in each cell. | |
dealii::Vector< double > | artificial_dissipation_se |
Artificial dissipation error ratio sensor in each cell. | |
dealii::LinearAlgebra::distributed::Vector< real > | dual |
Current optimization dual variables corresponding to the residual constraints also known as the adjoint. More... | |
bool | update_artificial_diss |
const dealii::hp::FECollection< dim > | fe_collection |
Finite Element Collection for p-finite-element to represent the solution. More... | |
dealii::hp::QCollection< dim > | volume_quadrature_collection |
Finite Element Collection to represent the high-order grid. More... | |
dealii::hp::QCollection< dim-1 > | face_quadrature_collection |
Quadrature used to evaluate face integrals. | |
const dealii::hp::FECollection< dim > | fe_collection_lagrange |
Lagrange basis used in strong form. More... | |
const dealii::hp::FECollection< 1 > | oneD_fe_collection |
1D Finite Element Collection for p-finite-element to represent the solution More... | |
const dealii::hp::FECollection< 1 > | oneD_fe_collection_1state |
1D Finite Element Collection for p-finite-element to represent the solution for a single state. More... | |
const dealii::hp::FECollection< 1 > | oneD_fe_collection_flux |
1D collocated flux basis used in strong form More... | |
dealii::hp::QCollection< 1 > | oneD_quadrature_collection |
1D quadrature to generate Lagrange polynomials for the sake of flux interpolation. | |
dealii::QGauss< 0 > | oneD_face_quadrature |
1D surface quadrature is always one single point for all poly degrees. | |
dealii::DoFHandler< dim > | dof_handler |
Finite Element Collection to represent the high-order grid. More... | |
std::shared_ptr< HighOrderGrid< dim, real, MeshType > > | high_order_grid |
High order grid that will provide the MappingFEField. | |
double | assemble_residual_time |
Computational time for assembling residual. | |
bool | freeze_artificial_dissipation |
Flag to freeze artificial dissipation. | |
double | max_artificial_dissipation_coeff |
Stores maximum artificial dissipation while assembling the residual. | |
bool | use_auxiliary_eq |
Flag for using the auxiliary equation. | |
Protected Member Functions | |
virtual void | assemble_volume_term_and_build_operators (typename dealii::DoFHandler< dim >::active_cell_iterator cell, const dealii::types::global_dof_index current_cell_index, const std::vector< dealii::types::global_dof_index > &cell_dofs_indices, const std::vector< dealii::types::global_dof_index > &metric_dof_indices, const unsigned int poly_degree, const unsigned int grid_degree, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis, OPERATOR::local_basis_stiffness< dim, 2 *dim, real > &flux_basis_stiffness, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_int, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_ext, OPERATOR::metric_operators< real, dim, 2 *dim > &metric_oper, OPERATOR::mapping_shape_functions< dim, 2 *dim, real > &mapping_basis, std::array< std::vector< real >, dim > &mapping_support_points, dealii::hp::FEValues< dim, dim > &fe_values_collection_volume, dealii::hp::FEValues< dim, dim > &fe_values_collection_volume_lagrange, const dealii::FESystem< dim, dim > ¤t_fe_ref, dealii::Vector< real > &local_rhs_int_cell, std::vector< dealii::Tensor< 1, dim, real >> &local_auxiliary_RHS, const bool compute_auxiliary_right_hand_side, const bool compute_dRdW, const bool compute_dRdX, const bool compute_d2R)=0 |
Builds the necessary operators/fe values and assembles volume residual. | |
virtual void | assemble_boundary_term_and_build_operators (typename dealii::DoFHandler< dim >::active_cell_iterator cell, const dealii::types::global_dof_index current_cell_index, const unsigned int iface, const unsigned int boundary_id, const real penalty, const std::vector< dealii::types::global_dof_index > &cell_dofs_indices, const std::vector< dealii::types::global_dof_index > &metric_dof_indices, const unsigned int poly_degree, const unsigned int grid_degree, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis, OPERATOR::local_basis_stiffness< dim, 2 *dim, real > &flux_basis_stiffness, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_int, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_ext, OPERATOR::metric_operators< real, dim, 2 *dim > &metric_oper, OPERATOR::mapping_shape_functions< dim, 2 *dim, real > &mapping_basis, std::array< std::vector< real >, dim > &mapping_support_points, dealii::hp::FEFaceValues< dim, dim > &fe_values_collection_face_int, const dealii::FESystem< dim, dim > ¤t_fe_ref, dealii::Vector< real > &local_rhs_int_cell, std::vector< dealii::Tensor< 1, dim, real >> &local_auxiliary_RHS, const bool compute_auxiliary_right_hand_side, const bool compute_dRdW, const bool compute_dRdX, const bool compute_d2R)=0 |
Builds the necessary operators/fe values and assembles boundary residual. | |
virtual void | assemble_face_term_and_build_operators (typename dealii::DoFHandler< dim >::active_cell_iterator cell, typename dealii::DoFHandler< dim >::active_cell_iterator neighbor_cell, const dealii::types::global_dof_index current_cell_index, const dealii::types::global_dof_index neighbor_cell_index, const unsigned int iface, const unsigned int neighbor_iface, const real penalty, const std::vector< dealii::types::global_dof_index > ¤t_dofs_indices, const std::vector< dealii::types::global_dof_index > &neighbor_dofs_indices, const std::vector< dealii::types::global_dof_index > ¤t_metric_dofs_indices, const std::vector< dealii::types::global_dof_index > &neighbor_metric_dofs_indices, const unsigned int poly_degree_int, const unsigned int poly_degree_ext, const unsigned int grid_degree_int, const unsigned int grid_degree_ext, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis_int, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis_ext, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis_int, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis_ext, OPERATOR::local_basis_stiffness< dim, 2 *dim, real > &flux_basis_stiffness, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_int, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_ext, OPERATOR::metric_operators< real, dim, 2 *dim > &metric_oper_int, OPERATOR::metric_operators< real, dim, 2 *dim > &metric_oper_ext, OPERATOR::mapping_shape_functions< dim, 2 *dim, real > &mapping_basis, std::array< std::vector< real >, dim > &mapping_support_points, dealii::hp::FEFaceValues< dim, dim > &fe_values_collection_face_int, dealii::hp::FEFaceValues< dim, dim > &fe_values_collection_face_ext, dealii::Vector< real > ¤t_cell_rhs, dealii::Vector< real > &neighbor_cell_rhs, std::vector< dealii::Tensor< 1, dim, real >> ¤t_cell_rhs_aux, dealii::LinearAlgebra::distributed::Vector< double > &rhs, std::array< dealii::LinearAlgebra::distributed::Vector< double >, dim > &rhs_aux, const bool compute_auxiliary_right_hand_side, const bool compute_dRdW, const bool compute_dRdX, const bool compute_d2R)=0 |
Builds the necessary operators/fe values and assembles face residual. | |
virtual void | assemble_subface_term_and_build_operators (typename dealii::DoFHandler< dim >::active_cell_iterator cell, typename dealii::DoFHandler< dim >::active_cell_iterator neighbor_cell, const dealii::types::global_dof_index current_cell_index, const dealii::types::global_dof_index neighbor_cell_index, const unsigned int iface, const unsigned int neighbor_iface, const unsigned int neighbor_i_subface, const real penalty, const std::vector< dealii::types::global_dof_index > ¤t_dofs_indices, const std::vector< dealii::types::global_dof_index > &neighbor_dofs_indices, const std::vector< dealii::types::global_dof_index > ¤t_metric_dofs_indices, const std::vector< dealii::types::global_dof_index > &neighbor_metric_dofs_indices, const unsigned int poly_degree_int, const unsigned int poly_degree_ext, const unsigned int grid_degree_int, const unsigned int grid_degree_ext, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis_int, OPERATOR::basis_functions< dim, 2 *dim, real > &soln_basis_ext, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis_int, OPERATOR::basis_functions< dim, 2 *dim, real > &flux_basis_ext, OPERATOR::local_basis_stiffness< dim, 2 *dim, real > &flux_basis_stiffness, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_int, OPERATOR::vol_projection_operator< dim, 2 *dim, real > &soln_basis_projection_oper_ext, OPERATOR::metric_operators< real, dim, 2 *dim > &metric_oper_int, OPERATOR::metric_operators< real, dim, 2 *dim > &metric_oper_ext, OPERATOR::mapping_shape_functions< dim, 2 *dim, real > &mapping_basis, std::array< std::vector< real >, dim > &mapping_support_points, dealii::hp::FEFaceValues< dim, dim > &fe_values_collection_face_int, dealii::hp::FESubfaceValues< dim, dim > &fe_values_collection_subface, dealii::Vector< real > ¤t_cell_rhs, dealii::Vector< real > &neighbor_cell_rhs, std::vector< dealii::Tensor< 1, dim, real >> ¤t_cell_rhs_aux, dealii::LinearAlgebra::distributed::Vector< double > &rhs, std::array< dealii::LinearAlgebra::distributed::Vector< double >, dim > &rhs_aux, const bool compute_auxiliary_right_hand_side, const bool compute_dRdW, const bool compute_dRdX, const bool compute_d2R)=0 |
Builds the necessary operators/fe values and assembles subface residual. | |
virtual void | assemble_volume_term_derivatives (typename dealii::DoFHandler< dim >::active_cell_iterator cell, const dealii::types::global_dof_index current_cell_index, const dealii::FEValues< dim, dim > &, const dealii::FESystem< dim, dim > &fe, const dealii::Quadrature< dim > &quadrature, const std::vector< dealii::types::global_dof_index > &metric_dof_indices, const std::vector< dealii::types::global_dof_index > &soln_dof_indices, dealii::Vector< real > &local_rhs_cell, const dealii::FEValues< dim, dim > &, const bool compute_dRdW, const bool compute_dRdX, const bool compute_d2R)=0 |
Evaluate the integral over the cell volume and the specified derivatives. More... | |
virtual void | assemble_boundary_term_derivatives (typename dealii::DoFHandler< dim >::active_cell_iterator cell, const dealii::types::global_dof_index current_cell_index, const unsigned int face_number, const unsigned int boundary_id, const dealii::FEFaceValuesBase< dim, dim > &fe_values_boundary, const real penalty, const dealii::FESystem< dim, dim > &fe, const dealii::Quadrature< dim-1 > &quadrature, const std::vector< dealii::types::global_dof_index > &metric_dof_indices, const std::vector< dealii::types::global_dof_index > &soln_dof_indices, dealii::Vector< real > &local_rhs_cell, const bool compute_dRdW, const bool compute_dRdX, const bool compute_d2R)=0 |
Evaluate the integral over the cell edges that are on domain boundaries and the specified derivatives. More... | |
virtual void | assemble_face_term_derivatives (typename dealii::DoFHandler< dim >::active_cell_iterator cell, const dealii::types::global_dof_index current_cell_index, const dealii::types::global_dof_index neighbor_cell_index, const std::pair< unsigned int, int > face_subface_int, const std::pair< unsigned int, int > face_subface_ext, const typename dealii::QProjector< dim >::DataSetDescriptor face_data_set_int, const typename dealii::QProjector< dim >::DataSetDescriptor face_data_set_ext, const dealii::FEFaceValuesBase< dim, dim > &, const dealii::FEFaceValuesBase< dim, dim > &, const real penalty, const dealii::FESystem< dim, dim > &fe_int, const dealii::FESystem< dim, dim > &fe_ext, const dealii::Quadrature< dim-1 > &face_quadrature, const std::vector< dealii::types::global_dof_index > &metric_dof_indices_int, const std::vector< dealii::types::global_dof_index > &metric_dof_indices_ext, const std::vector< dealii::types::global_dof_index > &soln_dof_indices_int, const std::vector< dealii::types::global_dof_index > &soln_dof_indices_ext, dealii::Vector< real > &local_rhs_int_cell, dealii::Vector< real > &local_rhs_ext_cell, const bool compute_dRdW, const bool compute_dRdX, const bool compute_d2R)=0 |
Evaluate the integral over the internal cell edges and its specified derivatives. More... | |
virtual void | assemble_volume_term_explicit (typename dealii::DoFHandler< dim >::active_cell_iterator cell, const dealii::types::global_dof_index current_cell_index, const dealii::FEValues< dim, dim > &fe_values_volume, const std::vector< dealii::types::global_dof_index > ¤t_dofs_indices, const std::vector< dealii::types::global_dof_index > &metric_dof_indices, const unsigned int poly_degree, const unsigned int grid_degree, dealii::Vector< real > ¤t_cell_rhs, const dealii::FEValues< dim, dim > &fe_values_lagrange)=0 |
Evaluate the integral over the cell volume. | |
Protected Attributes | |
real | current_time |
The current time set in set_current_time() | |
const dealii::FE_Q< dim > | fe_q_artificial_dissipation |
Continuous distribution of artificial dissipation. | |
dealii::DoFHandler< dim > | dof_handler_artificial_dissipation |
Degrees of freedom handler for C0 artificial dissipation. | |
dealii::LinearAlgebra::distributed::Vector< double > | artificial_dissipation_c0 |
Artificial dissipation coefficients. | |
const dealii::UpdateFlags | volume_update_flags |
Update flags needed at volume points. More... | |
const dealii::UpdateFlags | face_update_flags |
Update flags needed at face points. More... | |
const dealii::UpdateFlags | neighbor_face_update_flags = dealii::update_values | dealii::update_gradients | dealii::update_quadrature_points | dealii::update_JxW_values |
Update flags needed at neighbor' face points. More... | |
MPI_Comm | mpi_communicator |
MPI communicator. | |
dealii::ConditionalOStream | pcout |
Parallel std::cout that only outputs on mpi_rank==0. | |
Private Member Functions | |
virtual void | allocate_second_derivatives () |
Allocates the second derivatives. More... | |
virtual void | allocate_dRdX () |
Allocates the residual derivatives w.r.t the volume nodes. More... | |
void | allocate_artificial_dissipation () |
Allocates variables of artificial dissipation. More... | |
template<typename DoFCellAccessorType > | |
real | evaluate_penalty_scaling (const DoFCellAccessorType &cell, const int iface, const dealii::hp::FECollection< dim > fe_collection) const |
template<typename DoFCellAccessorType1 , typename DoFCellAccessorType2 > | |
bool | current_cell_should_do_the_work (const DoFCellAccessorType1 ¤t_cell, const DoFCellAccessorType2 &neighbor_cell) const |
In the case that two cells have the same coarseness, this function decides if the current cell should perform the work. More... | |
MassiveCollectionTuple | create_collection_tuple (const unsigned int max_degree, const int nstate, const Parameters::AllParameters *const parameters_input) const |
Used in the delegated constructor. More... | |
Private Attributes | |
dealii::LinearAlgebra::distributed::Vector< double > | solution_dRdW |
dealii::LinearAlgebra::distributed::Vector< double > | volume_nodes_dRdW |
double | CFL_mass_dRdW |
CFL used to add mass matrix in the optimization FlowConstraints class. | |
dealii::LinearAlgebra::distributed::Vector< double > | solution_dRdX |
dealii::LinearAlgebra::distributed::Vector< double > | volume_nodes_dRdX |
dealii::LinearAlgebra::distributed::Vector< double > | solution_d2R |
dealii::LinearAlgebra::distributed::Vector< double > | volume_nodes_d2R |
dealii::LinearAlgebra::distributed::Vector< double > | dual_d2R |
DGBase is independent of the number of state variables.
This base class allows the use of arrays to efficiently allocate the data structures through std::array in the derived class DGBaseState. This class is the one being returned by the DGFactory and is the main interface for a user to call its main functions such as "assemble_residual".
Discretizes the problem
\[ \frac{\partial \mathbf{u}}{\partial t} + \boldsymbol\nabla \cdot ( \mathbf{F}_{conv}(\mathbf{u}) + \mathbf{F}_{diss}(\mathbf{u},\boldsymbol\nabla\mathbf{u}) ) = \mathbf{q} \]
Definition at line 82 of file dg_base.hpp.
using PHiLiP::DGBase< dim, real, MeshType >::Triangulation = MeshType |
Triangulation to store the grid. In 1D, dealii::Triangulation<dim> is used. In 2D, 3D, dealii::parallel::distributed::Triangulation<dim> is used.
Definition at line 89 of file dg_base.hpp.
PHiLiP::DGBase< dim, real, MeshType >::DGBase | ( | const int | nstate_input, |
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 | ||
) |
Principal constructor that will call delegated constructor.
Will initialize mapping, fe_dg, all_parameters, volume_quadrature, and face_quadrature from DGBase. The it will new some FEValues that will be used to retrieve the finite element values at physical locations.
Passes create_collection_tuple() to the delegated constructor.
Definition at line 60 of file dg_base.cpp.
PHiLiP::DGBase< dim, real, MeshType >::DGBase | ( | const int | nstate_input, |
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, | ||
const MassiveCollectionTuple | collection_tuple | ||
) |
Delegated constructor that initializes collections.
Since a function is used to generate multiple different objects, a delegated constructor is used to unwrap the tuple and initialize the collections.
The tuple is built from create_collection_tuple().
Definition at line 71 of file dg_base.cpp.
void PHiLiP::DGBase< dim, real, MeshType >::add_mass_matrices | ( | const real | scale | ) |
Add mass matrices to the system scaled by a factor (likely time-step)
Although straightforward, this has not been tested yet. Will be required for accurate time-stepping or nonlinear problems
Definition at line 2802 of file dg_base.cpp.
void PHiLiP::DGBase< dim, real, MeshType >::add_time_scaled_mass_matrices | ( | ) |
Add time scaled mass matrices to the system.
For pseudotime-stepping where the scaling depends on wavespeed and cell-size.
Definition at line 2807 of file dg_base.cpp.
|
private |
Allocates variables of artificial dissipation.
It is called by allocate_system() when artificial dissipation is set to true in the parameters file.
Definition at line 1971 of file dg_base.cpp.
|
privatevirtual |
Allocates the residual derivatives w.r.t the volume nodes.
Is called when assembling the residual's second derivatives, and is currently empty due to being cleared by the allocate_system().
Definition at line 2016 of file dg_base.cpp.
|
pure virtual |
Allocate the dual vector for optimization.
Currently only used in weak form.
Implemented in PHiLiP::DGWeak< dim, nstate, real, MeshType >, and PHiLiP::DGStrong< dim, nstate, real, MeshType >.
|
privatevirtual |
Allocates the second derivatives.
Is called when assembling the residual's second derivatives, and is currently empty due to being cleared by the allocate_system().
Definition at line 1990 of file dg_base.cpp.
|
virtual |
Allocates the system.
Must be done after setting the mesh and before assembling the system.
Definition at line 1866 of file dg_base.cpp.
void PHiLiP::DGBase< dim, real, MeshType >::apply_global_mass_matrix | ( | const dealii::LinearAlgebra::distributed::Vector< double > & | input_vector, |
dealii::LinearAlgebra::distributed::Vector< double > & | output_vector, | ||
const bool | use_auxiliary_eq = false , |
||
const bool | use_unmodified_mass_matrix = false |
||
) |
Applies the local metric dependent mass matrices when the global is not stored.
We use matrix-free methods to apply the local mass matrix on-the-fly in each cell using sum-factorization techniques. use_unmodified_mass_matrix flag allows the unmodified mass matrix to be used for FR, i.e., use M rather than M+K. For example, if this function is used for an inner product <a,b>, setting use_unmodified_mass_matrix = true
will result in an M norm (L2) inner product, whereas setting use_unmodified_mass_matrix = false
will result in an M+K norm (broken Sobolev) inner product.
Definition at line 2625 of file dg_base.cpp.
void PHiLiP::DGBase< dim, real, MeshType >::apply_inverse_global_mass_matrix | ( | const dealii::LinearAlgebra::distributed::Vector< double > & | input_vector, |
dealii::LinearAlgebra::distributed::Vector< double > & | output_vector, | ||
const bool | use_auxiliary_eq = false |
||
) |
Applies the inverse of the local metric dependent mass matrices when the global is not stored.
We use matrix-free methods to apply the inverse of the local mass matrix on-the-fly in each cell using sum-factorization techniques.
Definition at line 2452 of file dg_base.cpp.
|
protectedpure virtual |
Evaluate the integral over the cell edges that are on domain boundaries and the specified derivatives.
Compute both the right-hand side and the corresponding block of dRdW, dRdX, and/or d2R.
Implemented in PHiLiP::DGWeak< dim, nstate, real, MeshType >, and PHiLiP::DGStrong< dim, nstate, real, MeshType >.
template void PHiLiP::DGBase< dim, real, MeshType >::assemble_cell_residual< dealii::TriaActiveIterator< dealii::DoFCellAccessor< PHILIP_DIM, PHILIP_DIM, false > >, dealii::TriaActiveIterator< dealii::DoFCellAccessor< PHILIP_DIM, PHILIP_DIM, false > > > | ( | const DoFCellAccessorType1 & | current_cell, |
const DoFCellAccessorType2 & | current_metric_cell, | ||
const bool | compute_dRdW, | ||
const bool | compute_dRdX, | ||
const bool | compute_d2R, | ||
dealii::hp::FEValues< dim, dim > & | fe_values_collection_volume, | ||
dealii::hp::FEFaceValues< dim, dim > & | fe_values_collection_face_int, | ||
dealii::hp::FEFaceValues< dim, dim > & | fe_values_collection_face_ext, | ||
dealii::hp::FESubfaceValues< dim, dim > & | fe_values_collection_subface, | ||
dealii::hp::FEValues< dim, dim > & | fe_values_collection_volume_lagrange, | ||
OPERATOR::basis_functions< dim, 2 *dim, real > & | soln_basis_int, | ||
OPERATOR::basis_functions< dim, 2 *dim, real > & | soln_basis_ext, | ||
OPERATOR::basis_functions< dim, 2 *dim, real > & | flux_basis_int, | ||
OPERATOR::basis_functions< dim, 2 *dim, real > & | flux_basis_ext, | ||
OPERATOR::local_basis_stiffness< dim, 2 *dim, real > & | flux_basis_stiffness, | ||
OPERATOR::vol_projection_operator< dim, 2 *dim, real > & | soln_basis_projection_oper_int, | ||
OPERATOR::vol_projection_operator< dim, 2 *dim, real > & | soln_basis_projection_oper_ext, | ||
OPERATOR::mapping_shape_functions< dim, 2 *dim, real > & | mapping_basis, | ||
const bool | compute_auxiliary_right_hand_side, | ||
dealii::LinearAlgebra::distributed::Vector< double > & | rhs, | ||
std::array< dealii::LinearAlgebra::distributed::Vector< double >, dim > & | rhs_aux | ||
) |
Used in assemble_residual().
IMPORTANT: This does not fully compute the cell residual since it might not perform the work on all the faces. All the active cells must be traversed to ensure that the right hand side is correct.
Definition at line 454 of file dg_base.cpp.
|
protectedpure virtual |
Evaluate the integral over the internal cell edges and its specified derivatives.
Compute both the right-hand side and the block of the Jacobian. This adds the contribution to both cell's residual and effectively computes 4 block contributions to dRdX blocks.
Implemented in PHiLiP::DGWeak< dim, nstate, real, MeshType >, and PHiLiP::DGStrong< dim, nstate, real, MeshType >.
void PHiLiP::DGBase< dim, real, MeshType >::assemble_residual | ( | const bool | compute_dRdW = false , |
const bool | compute_dRdX = false , |
||
const bool | compute_d2R = false , |
||
const double | CFL_mass = 0.0 |
||
) |
Main loop of the DG class.
Evaluates the right-hand-side \( \mathbf{R(\mathbf{u}}) \) of the system
\[ \frac{\partial \mathbf{u}}{\partial t} = \mathbf{R(\mathbf{u}}) = - \boldsymbol\nabla \cdot ( \mathbf{F}_{conv}(\mathbf{u}) + \mathbf{F}_{diss}(\mathbf{u},\boldsymbol\nabla\mathbf{u}) ) + \mathbf{q} \]
As well as sets the
\[ \mathbf{\text{system_matrix}} = \frac{\partial \mathbf{R}}{\partial \mathbf{u}} \]
It loops over all the cells, evaluates the volume contributions, then loops over the faces of the current cell. Four scenarios may happen
< FEValues of volume.
< FEValues of interior face.
< FEValues of exterior face.
< FEValues of subface.
Definition at line 1091 of file dg_base.cpp.
|
protectedpure virtual |
Evaluate the integral over the cell volume and the specified derivatives.
Compute both the right-hand side and the corresponding block of dRdW, dRdX, and/or d2R.
Implemented in PHiLiP::DGWeak< dim, nstate, real, MeshType >, and PHiLiP::DGStrong< dim, nstate, real, MeshType >.
|
private |
Used in the delegated constructor.
The main reason we use this weird function is because all of the above objects need to be looped with the various p-orders. This function allows us to do this in a single function instead of having like 6 different functions to initialize each of them.
Definition at line 142 of file dg_base.cpp.
|
private |
In the case that two cells have the same coarseness, this function decides if the current cell should perform the work.
In the case the neighbor is a ghost cell, we let the processor with the lower rank do the work on that face. We cannot use the cell->index() because the index is relative to the distributed triangulation. Therefore, the cell index of a ghost cell might be different to the physical cell index even if they refer to the same cell.
For a locally owned neighbor cell, cell with lower index does work or if both cells have same index, then cell at the lower level does the work See https://www.dealii.org/developer/doxygen/deal.II/classTriaAccessorBase.html#a695efcbe84fefef3e4c93ee7bdb446ad
Definition at line 418 of file dg_base.cpp.
template RadFadType PHiLiP::DGBase< dim, real, MeshType >::discontinuity_sensor< RadFadType > | ( | const dealii::Quadrature< dim > & | volume_quadrature, |
const std::vector< real2 > & | soln_coeff_high, | ||
const dealii::FiniteElement< dim, dim > & | fe_high, | ||
const std::vector< real2 > & | jac_det | ||
) |
Discontinuity sensor with 4 parameters, based on projecting to p-1.
Definition at line 2932 of file dg_base.cpp.
|
private |
Evaluate the average penalty term at the face. For a cell with solution of degree p, and Hausdorff measure h, which represents the element dimension orthogonal to the face, the penalty term is given by p*(p+1)/h .
Definition at line 398 of file dg_base.cpp.
std::vector<real> PHiLiP::DGBase< dim, real, MeshType >::evaluate_time_steps | ( | const bool | exact_time_stepping | ) |
Evaluates the maximum stable time step.
If exact_time_stepping = true, use the same time step for the entire solution NOT YET IMPLEMENTED
void PHiLiP::DGBase< dim, real, MeshType >::reinit | ( | ) |
Reinitializes the DG object after a change of triangulation.
Calls respective function for high-order-grid and initializes dof_handler again. Also resets all fe_degrees to intial_degree set during constructor.
Definition at line 112 of file dg_base.cpp.
void PHiLiP::DGBase< dim, real, MeshType >::set_all_cells_fe_degree | ( | const unsigned int | degree | ) |
Refers to a collection Mappings, which represents the high-order grid.
Since we are interested in performing mesh movement for optimization purposes, this is not a constant member variables.
Definition at line 293 of file dg_base.cpp.
void PHiLiP::DGBase< dim, real, MeshType >::set_anisotropic_flags | ( | ) |
Set anisotropic flags based on jump indicator.
Some cells must have already been tagged for refinement through some other indicator
void PHiLiP::DGBase< dim, real, MeshType >::time_scale_solution_update | ( | dealii::LinearAlgebra::distributed::Vector< double > & | solution_update, |
const real | CFL | ||
) | const |
Scales a solution update with the appropriate maximum time step.
Used for steady state solutions using the explicit ODE solver.
Definition at line 265 of file dg_base.cpp.
void PHiLiP::DGBase< dim, real, MeshType >::time_scaled_mass_matrices | ( | const real | scale | ) |
Evaluate the time_scaled_global_mass_matrix such that the maximum time step cell-wise is taken into account.
Definition at line 2812 of file dg_base.cpp.
void PHiLiP::DGBase< dim, real, MeshType >::update_artificial_dissipation_discontinuity_sensor | ( | ) |
dealii::Vector<double> PHiLiP::DGBase< dim, real, MeshType >::cell_volume |
Time it takes for the maximum wavespeed to cross the cell domain.
Uses evaluate_CFL() which would be defined in the subclasses. This is because DGBase isn't templated on nstate and therefore, can't use the Physics to compute maximum wavespeeds.
Definition at line 450 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::d2RdWdW |
System matrix corresponding to the second derivatives of the right_hand_side with respect to the solution
Definition at line 360 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::d2RdWdX |
System matrix corresponding to the mixed second derivatives of the right_hand_side with respect to the solution and the volume volume_nodes
Definition at line 368 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::d2RdXdX |
System matrix corresponding to the second derivatives of the right_hand_side with respect to the volume volume_nodes
Definition at line 364 of file dg_base.hpp.
dealii::DoFHandler<dim> PHiLiP::DGBase< dim, real, MeshType >::dof_handler |
Finite Element Collection to represent the high-order grid.
This is a collection of FESystems. Unfortunately, deal.II doesn't have a working hp Mapping FE field. Therefore, every grid/cell will use the maximal polynomial mapping regardless of the solution order.Degrees of freedom handler
Definition at line 652 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::dRdXv |
System matrix corresponding to the derivative of the right_hand_side with respect to the volume volume_nodes Xv
Definition at line 356 of file dg_base.hpp.
dealii::LinearAlgebra::distributed::Vector<real> PHiLiP::DGBase< dim, real, MeshType >::dual |
Current optimization dual variables corresponding to the residual constraints also known as the adjoint.
This is used to evaluate the dot-product between the dual and the 2nd derivatives of the residual since storing the 2nd order partials of the residual is a very large 3rd order tensor.
Definition at line 479 of file dg_base.hpp.
|
private |
Dual variables to compute d2R last Will be used to avoid recomputing d2R.
Definition at line 442 of file dg_base.hpp.
|
protected |
Update flags needed at face points.
Definition at line 873 of file dg_base.hpp.
const dealii::hp::FECollection<dim> PHiLiP::DGBase< dim, real, MeshType >::fe_collection |
Finite Element Collection for p-finite-element to represent the solution.
This is a collection of FESystems
Definition at line 597 of file dg_base.hpp.
const dealii::hp::FECollection<dim> PHiLiP::DGBase< dim, real, MeshType >::fe_collection_lagrange |
Lagrange basis used in strong form.
This is a collection of scalar Lagrange bases
Definition at line 614 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::global_inverse_mass_matrix |
Global inverser mass matrix.
Should be block diagonal where each block contains the inverse mass matrix of each cell.
Definition at line 332 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::global_mass_matrix |
Global mass matrix.
Should be block diagonal where each block contains the mass matrix of each cell.
Definition at line 329 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::global_mass_matrix_auxiliary |
Global auxiliary mass matrix.
Note that it has a mass matrix in each dimension since the auxiliary variable is a tensor of size dim. We use the same matrix in each dim.
Definition at line 336 of file dg_base.hpp.
dealii::IndexSet PHiLiP::DGBase< dim, real, MeshType >::locally_relevant_dofs_grid |
Union of locally owned degrees of freedom and relevant ghost degrees of freedom for the grid
Definition at line 404 of file dg_base.hpp.
dealii::SparsityPattern PHiLiP::DGBase< dim, real, MeshType >::mass_sparsity_pattern |
Sparsity pattern used on the system_matrix.
Not sure we need to store it.
Definition at line 321 of file dg_base.hpp.
const unsigned int PHiLiP::DGBase< dim, real, MeshType >::max_degree |
Maximum degree used for p-refi1nement.
This is known through the constructor parameters. DGBase cannot use nstate as a compile-time known.
Definition at line 104 of file dg_base.hpp.
dealii::Vector<double> PHiLiP::DGBase< dim, real, MeshType >::max_dt_cell |
Time it takes for the maximum wavespeed to cross the cell domain.
Uses evaluate_CFL() which would be defined in the subclasses. This is because DGBase isn't templated on nstate and therefore, can't use the Physics to compute maximum wavespeeds.
Definition at line 457 of file dg_base.hpp.
const unsigned int PHiLiP::DGBase< dim, real, MeshType >::max_grid_degree |
Maximum grid degree used for hp-refi1nement.
This is known through the constructor parameters. DGBase cannot use nstate as a compile-time known.
Definition at line 109 of file dg_base.hpp.
|
protected |
Update flags needed at neighbor' face points.
NOTE: With hp-adaptation, might need to query neighbor's quadrature points depending on the order of the cells.
Definition at line 877 of file dg_base.hpp.
const int PHiLiP::DGBase< dim, real, MeshType >::nstate |
Number of state variables.
This is known through the constructor parameters. DGBase cannot use nstate as a compile-time known.
Definition at line 96 of file dg_base.hpp.
const dealii::hp::FECollection<1> PHiLiP::DGBase< dim, real, MeshType >::oneD_fe_collection |
1D Finite Element Collection for p-finite-element to represent the solution
This is a collection of FESystems for 1D.
Definition at line 620 of file dg_base.hpp.
const dealii::hp::FECollection<1> PHiLiP::DGBase< dim, real, MeshType >::oneD_fe_collection_1state |
1D Finite Element Collection for p-finite-element to represent the solution for a single state.
This is a collection of FESystems for 1D. Since each state is represented by the same polynomial degree, for the RHS, we only need to store the 1D basis functions for a single state.
Definition at line 627 of file dg_base.hpp.
const dealii::hp::FECollection<1> PHiLiP::DGBase< dim, real, MeshType >::oneD_fe_collection_flux |
1D collocated flux basis used in strong form
This is a collection of collocated Lagrange bases for 1D.
Definition at line 630 of file dg_base.hpp.
dealii::LinearAlgebra::distributed::Vector<double> PHiLiP::DGBase< dim, real, MeshType >::right_hand_side |
Residual of the current solution.
Weak form.
The right-hand side sends all the term to the side of the source term.
Given
\[ \frac{\partial \mathbf{u}}{\partial t} + \boldsymbol\nabla \cdot ( \mathbf{F}_{conv}(\mathbf{u}) + \mathbf{F}_{diss}(\mathbf{u},\boldsymbol\nabla\mathbf{u}) ) = \mathbf{q} \]
The right-hand side is given by
\[ \mathbf{\text{rhs}} = - \boldsymbol\nabla \cdot ( \mathbf{F}_{conv}(\mathbf{u}) + \mathbf{F}_{diss}(\mathbf{u},\boldsymbol\nabla\mathbf{u}) ) + \mathbf{q} \]
It is important to note that the \(\mathbf{F}_{diss}\) is positive in the DG formulation. Therefore, the PhysicsBase class should have a negative when considering stable applications of diffusion.
Definition at line 396 of file dg_base.hpp.
dealii::LinearAlgebra::distributed::Vector<double> PHiLiP::DGBase< dim, real, MeshType >::solution |
Current modal coefficients of the solution.
Note that the current processor has read-access to all locally_relevant_dofs and has write-access to all locally_owned_dofs
Definition at line 409 of file dg_base.hpp.
|
private |
Modal coefficients of the solution used to compute d2R last Will be used to avoid recomputing d2R.
Definition at line 436 of file dg_base.hpp.
|
private |
Modal coefficients of the solution used to compute dRdW last Will be used to avoid recomputing dRdW.
Definition at line 419 of file dg_base.hpp.
|
private |
Modal coefficients of the solution used to compute dRdX last Will be used to avoid recomputing dRdX.
Definition at line 429 of file dg_base.hpp.
dealii::SparsityPattern PHiLiP::DGBase< dim, real, MeshType >::sparsity_pattern |
Sparsity pattern used on the system_matrix.
Not sure we need to store it.
Definition at line 317 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::system_matrix |
System matrix corresponding to the derivative of the right_hand_side with respect to the solution
Definition at line 343 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::system_matrix_transpose |
System matrix corresponding to the derivative of the right_hand_side with respect to the solution TRANSPOSED.
Definition at line 347 of file dg_base.hpp.
dealii::TrilinosWrappers::SparseMatrix PHiLiP::DGBase< dim, real, MeshType >::time_scaled_global_mass_matrix |
Global mass matrix divided by the time scales.
Should be block diagonal where each block contains the scaled mass matrix of each cell.
Definition at line 325 of file dg_base.hpp.
|
private |
Modal coefficients of the grid nodes used to compute d2R last Will be used to avoid recomputing d2R.
Definition at line 439 of file dg_base.hpp.
|
private |
Modal coefficients of the grid nodes used to compute dRdW last Will be used to avoid recomputing dRdW.
Definition at line 422 of file dg_base.hpp.
|
private |
Modal coefficients of the grid nodes used to compute dRdX last Will be used to avoid recomputing dRdX.
Definition at line 432 of file dg_base.hpp.
dealii::hp::QCollection<dim> PHiLiP::DGBase< dim, real, MeshType >::volume_quadrature_collection |
Finite Element Collection to represent the high-order grid.
This is a collection of FESystems. Unfortunately, deal.II doesn't have a working hp Mapping FE field. Therefore, every grid/cell will use the maximal polynomial mapping regardless of the solution order.Quadrature used to evaluate volume integrals.
Definition at line 608 of file dg_base.hpp.
|
protected |
Update flags needed at volume points.
Definition at line 870 of file dg_base.hpp.