1 #ifndef __OPERATORS_H__ 2 #define __OPERATORS_H__ 4 #include <deal.II/base/conditional_ostream.h> 5 #include <deal.II/base/parameter_handler.h> 7 #include <deal.II/base/qprojector.h> 9 #include <deal.II/grid/tria.h> 11 #include <deal.II/fe/fe_dgq.h> 12 #include <deal.II/fe/fe_dgp.h> 13 #include <deal.II/fe/fe_system.h> 14 #include <deal.II/fe/mapping_fe_field.h> 15 #include <deal.II/fe/fe_q.h> 17 #include <deal.II/dofs/dof_handler.h> 19 #include <deal.II/hp/q_collection.h> 20 #include <deal.II/hp/mapping_collection.h> 21 #include <deal.II/hp/fe_values.h> 23 #include <deal.II/lac/vector.h> 24 #include <deal.II/lac/sparsity_pattern.h> 25 #include <deal.II/lac/trilinos_sparse_matrix.h> 26 #include <deal.II/lac/trilinos_vector.h> 28 #include <Epetra_RowMatrixTransposer.h> 31 #include "ADTypes.hpp" 33 #include <CoDiPack/include/codi.hpp> 35 #include "parameters/all_parameters.h" 36 #include "parameters/parameters.h" 52 template <
int dim,
int n_faces,
typename real>
61 const int nstate_input,
62 const unsigned int max_degree_input,
63 const unsigned int grid_degree_input);
80 const dealii::FullMatrix<double> &basis_x,
81 const dealii::FullMatrix<double> &basis_y,
82 const dealii::FullMatrix<double> &basis_z);
93 const dealii::FullMatrix<double> &basis_x,
94 const dealii::FullMatrix<double> &basis_y,
95 const dealii::FullMatrix<double> &basis_z);
102 const std::vector<real> &input_vect,
103 std::vector<real> &output_vect,
104 const dealii::FullMatrix<double> &basis_x,
105 const dealii::FullMatrix<double> &basis_y,
106 const dealii::FullMatrix<double> &basis_z,
107 const bool adding =
false,
108 const double factor = 1.0) = 0;
111 const std::vector<real> &input_vect,
112 const std::vector<real> &weight_vect,
113 std::vector<real> &output_vect,
114 const dealii::FullMatrix<double> &basis_x,
115 const dealii::FullMatrix<double> &basis_y,
116 const dealii::FullMatrix<double> &basis_z,
117 const bool adding =
false,
118 const double factor = 1.0) = 0;
130 template<
int dim,
int n_faces,
typename real>
136 const int nstate_input,
137 const unsigned int max_degree_input,
138 const unsigned int grid_degree_input);
148 const std::vector<real> &input_vect,
149 std::vector<real> &output_vect,
150 const dealii::FullMatrix<double> &basis_x,
151 const dealii::FullMatrix<double> &basis_y,
152 const dealii::FullMatrix<double> &basis_z,
153 const bool adding =
false,
154 const double factor = 1.0)
override;
166 void divergence_matrix_vector_mult(
167 const dealii::Tensor<1,dim,std::vector<real>> &input_vect,
168 std::vector<real> &output_vect,
169 const dealii::FullMatrix<double> &basis_x,
170 const dealii::FullMatrix<double> &basis_y,
171 const dealii::FullMatrix<double> &basis_z,
172 const dealii::FullMatrix<double> &gradient_basis_x,
173 const dealii::FullMatrix<double> &gradient_basis_y,
174 const dealii::FullMatrix<double> &gradient_basis_z);
177 void divergence_matrix_vector_mult_1D(
178 const dealii::Tensor<1,dim,std::vector<real>> &input_vect,
179 std::vector<real> &output_vect,
180 const dealii::FullMatrix<double> &basis,
181 const dealii::FullMatrix<double> &gradient_basis);
184 void gradient_matrix_vector_mult(
185 const std::vector<real> &input_vect,
186 dealii::Tensor<1,dim,std::vector<real>> &output_vect,
187 const dealii::FullMatrix<double> &basis_x,
188 const dealii::FullMatrix<double> &basis_y,
189 const dealii::FullMatrix<double> &basis_z,
190 const dealii::FullMatrix<double> &gradient_basis_x,
191 const dealii::FullMatrix<double> &gradient_basis_y,
192 const dealii::FullMatrix<double> &gradient_basis_z);
194 void gradient_matrix_vector_mult_1D(
195 const std::vector<real> &input_vect,
196 dealii::Tensor<1,dim,std::vector<real>> &output_vect,
197 const dealii::FullMatrix<double> &basis,
198 const dealii::FullMatrix<double> &gradient_basis);
204 const std::vector<real> &input_vect,
205 const std::vector<real> &weight_vect,
206 std::vector<real> &output_vect,
207 const dealii::FullMatrix<double> &basis_x,
208 const dealii::FullMatrix<double> &basis_y,
209 const dealii::FullMatrix<double> &basis_z,
210 const bool adding =
false,
211 const double factor = 1.0)
override;
219 void divergence_two_pt_flux_Hadamard_product(
220 const dealii::Tensor<1,dim,dealii::FullMatrix<real>> &input_mat,
221 std::vector<real> &output_vect,
222 const std::vector<real> &weights,
223 const dealii::FullMatrix<double> &basis,
224 const double scaling = 2.0);
228 void surface_two_pt_flux_Hadamard_product(
229 const dealii::FullMatrix<real> &input_mat,
230 std::vector<real> &output_vect_vol,
231 std::vector<real> &output_vect_surf,
232 const std::vector<real> &weights,
233 const std::array<dealii::FullMatrix<double>,2> &surf_basis,
234 const unsigned int iface,
235 const unsigned int dim_not_zero,
236 const double scaling = 2.0);
248 void two_pt_flux_Hadamard_product(
249 const dealii::FullMatrix<real> &input_mat,
250 dealii::FullMatrix<real> &output_mat,
251 const dealii::FullMatrix<double> &basis,
252 const std::vector<real> &weights,
253 const int direction);
257 void sum_factorized_Hadamard_sparsity_pattern(
258 const unsigned int rows_size,
259 const unsigned int columns_size,
260 std::vector<std::array<unsigned int,dim>> &rows,
261 std::vector<std::array<unsigned int,dim>> &columns);
264 void sum_factorized_Hadamard_basis_assembly(
265 const unsigned int rows_size_1D,
266 const unsigned int columns_size_1D,
267 const std::vector<std::array<unsigned int,dim>> &rows,
268 const std::vector<std::array<unsigned int,dim>> &columns,
269 const dealii::FullMatrix<double> &basis,
270 const std::vector<double> &weights,
271 std::array<dealii::FullMatrix<double>,dim> &basis_sparse);
274 void sum_factorized_Hadamard_surface_sparsity_pattern(
275 const unsigned int rows_size,
276 const unsigned int columns_size,
277 std::vector<unsigned int> &rows,
278 std::vector<unsigned int> &columns,
279 const int dim_not_zero);
282 void sum_factorized_Hadamard_surface_basis_assembly(
283 const unsigned int rows_size,
284 const unsigned int columns_size_1D,
285 const std::vector<unsigned int> &rows,
286 const std::vector<unsigned int> &columns,
287 const dealii::FullMatrix<double> &basis,
288 const std::vector<double> &weights,
289 dealii::FullMatrix<double> &basis_sparse,
290 const int dim_not_zero);
296 void matrix_vector_mult_1D(
297 const std::vector<real> &input_vect,
298 std::vector<real> &output_vect,
299 const dealii::FullMatrix<double> &basis_x,
300 const bool adding =
false,
301 const double factor = 1.0);
307 void inner_product_1D(
308 const std::vector<real> &input_vect,
309 const std::vector<real> &weight_vect,
310 std::vector<real> &output_vect,
311 const dealii::FullMatrix<double> &basis_x,
312 const bool adding =
false,
313 const double factor = 1.0);
322 void matrix_vector_mult_surface_1D(
323 const unsigned int face_number,
324 const std::vector<real> &input_vect,
325 std::vector<real> &output_vect,
326 const std::array<dealii::FullMatrix<double>,2> &basis_surf,
327 const dealii::FullMatrix<double> &basis_vol,
328 const bool adding =
false,
329 const double factor = 1.0);
332 void inner_product_surface_1D(
333 const unsigned int face_number,
334 const std::vector<real> &input_vect,
335 const std::vector<real> &weight_vect,
336 std::vector<real> &output_vect,
337 const std::array<dealii::FullMatrix<double>,2> &basis_surf,
338 const dealii::FullMatrix<double> &basis_vol,
339 const bool adding =
false,
340 const double factor = 1.0);
347 void Hadamard_product(
348 const dealii::FullMatrix<real> &input_mat1,
349 const dealii::FullMatrix<real> &input_mat2,
350 dealii::FullMatrix<real> &output_mat);
380 template<
int dim,
int n_faces,
typename real>
386 const int nstate_input,
387 const unsigned int max_degree_input,
388 const unsigned int grid_degree_input);
394 void build_1D_volume_operator(
395 const dealii::FESystem<1,1> &finite_element,
396 const dealii::Quadrature<1> &quadrature);
399 void build_1D_gradient_operator(
400 const dealii::FESystem<1,1> &finite_element,
401 const dealii::Quadrature<1> &quadrature);
404 void build_1D_surface_operator(
405 const dealii::FESystem<1,1> &finite_element,
406 const dealii::Quadrature<0> &quadrature);
409 void build_1D_surface_gradient_operator(
410 const dealii::FESystem<1,1> &finite_element,
411 const dealii::Quadrature<0> &quadrature);
415 template<
int dim,
int n_faces,
typename real>
421 const int nstate_input,
422 const unsigned int max_degree_input,
423 const unsigned int grid_degree_input);
429 void build_1D_volume_operator(
430 const dealii::FESystem<1,1> &finite_element,
431 const dealii::Quadrature<1> &quadrature);
435 template<
int dim,
int n_faces,
typename real>
441 const int nstate_input,
442 const unsigned int max_degree_input,
443 const unsigned int grid_degree_input);
449 void build_1D_volume_operator(
450 const dealii::FESystem<1,1> &finite_element,
451 const dealii::Quadrature<1> &quadrature);
457 dealii::FullMatrix<double> build_dim_mass_matrix(
459 const unsigned int n_dofs,
const unsigned int n_quad_pts,
461 const std::vector<double> &det_Jac,
462 const std::vector<double> &quad_weights);
471 template<
int dim,
int n_faces,
typename real>
477 const int nstate_input,
478 const unsigned int max_degree_input,
479 const unsigned int grid_degree_input,
480 const bool store_skew_symmetric_form_input =
false);
489 void build_1D_volume_operator(
490 const dealii::FESystem<1,1> &finite_element,
491 const dealii::Quadrature<1> &quadrature);
498 template<
int dim,
int n_faces,
typename real>
504 const int nstate_input,
505 const unsigned int max_degree_input,
506 const unsigned int grid_degree_input);
512 void build_1D_volume_operator(
513 const dealii::FESystem<1,1> &finite_element,
514 const dealii::Quadrature<1> &quadrature);
518 template<
int dim,
int n_faces,
typename real>
524 const int nstate_input,
525 const unsigned int max_degree_input,
526 const unsigned int grid_degree_input);
532 void build_1D_volume_operator(
533 const dealii::FESystem<1,1> &finite_element,
534 const dealii::Quadrature<1> &quadrature);
538 template<
int dim,
int n_faces,
typename real>
544 const int nstate_input,
545 const unsigned int max_degree_input,
546 const unsigned int grid_degree_input,
548 const double FR_user_specified_correction_parameter_value_input=0.0);
565 void get_Huynh_g2_parameter (
566 const unsigned int curr_cell_degree,
572 void get_spectral_difference_parameter (
573 const unsigned int curr_cell_degree,
579 void get_c_negative_FR_parameter (
580 const unsigned int curr_cell_degree,
587 void get_c_negative_divided_by_two_FR_parameter (
588 const unsigned int curr_cell_degree,
596 void get_c_plus_parameter (
597 const unsigned int curr_cell_degree,
608 void get_FR_correction_parameter (
609 const unsigned int curr_cell_degree,
616 void build_local_Flux_Reconstruction_operator(
617 const dealii::FullMatrix<double> &local_Mass_Matrix,
618 const dealii::FullMatrix<double> &pth_derivative,
619 const unsigned int n_dofs,
621 dealii::FullMatrix<double> &Flux_Reconstruction_operator);
624 void build_1D_volume_operator(
625 const dealii::FESystem<1,1> &finite_element,
626 const dealii::Quadrature<1> &quadrature);
634 dealii::FullMatrix<double> build_dim_Flux_Reconstruction_operator(
635 const dealii::FullMatrix<double> &local_Mass_Matrix,
637 const unsigned int n_dofs);
646 dealii::FullMatrix<double> build_dim_Flux_Reconstruction_operator_directly(
648 const unsigned int n_dofs,
649 dealii::FullMatrix<double> &pth_deriv,
650 dealii::FullMatrix<double> &mass_matrix);
658 template<
int dim,
int n_faces,
typename real>
664 const int nstate_input,
665 const unsigned int max_degree_input,
666 const unsigned int grid_degree_input,
686 void get_FR_aux_correction_parameter (
687 const unsigned int curr_cell_degree,
691 void build_1D_volume_operator(
692 const dealii::FESystem<1,1> &finite_element,
693 const dealii::Quadrature<1> &quadrature);
697 template<
int dim,
int n_faces,
typename real>
703 const int nstate_input,
704 const unsigned int max_degree_input,
705 const unsigned int grid_degree_input);
711 void compute_local_vol_projection_operator(
712 const dealii::FullMatrix<double> &norm_matrix_inverse,
713 const dealii::FullMatrix<double> &integral_vol_basis,
714 dealii::FullMatrix<double> &volume_projection);
717 void build_1D_volume_operator(
718 const dealii::FESystem<1,1> &finite_element,
719 const dealii::Quadrature<1> &quadrature);
723 template<
int dim,
int n_faces,
typename real>
729 const int nstate_input,
730 const unsigned int max_degree_input,
731 const unsigned int grid_degree_input,
733 const double FR_user_specified_correction_parameter_value_input=0.0,
734 const bool store_transpose_input =
false);
749 void build_1D_volume_operator(
750 const dealii::FESystem<1,1> &finite_element,
751 const dealii::Quadrature<1> &quadrature);
758 template<
int dim,
int n_faces,
typename real>
764 const int nstate_input,
765 const unsigned int max_degree_input,
766 const unsigned int grid_degree_input,
768 const bool store_transpose_input =
false);
780 void build_1D_volume_operator(
781 const dealii::FESystem<1,1> &finite_element,
782 const dealii::Quadrature<1> &quadrature);
789 template<
int dim,
int n_faces,
typename real>
795 const int nstate_input,
796 const unsigned int max_degree_input,
797 const unsigned int grid_degree_input,
799 const double FR_user_specified_correction_parameter_value_input=0.0);
811 void build_1D_volume_operator(
812 const dealii::FESystem<1,1> &finite_element,
813 const dealii::Quadrature<1> &quadrature);
816 template<
int dim,
int n_faces,
typename real>
822 const int nstate_input,
823 const unsigned int max_degree_input,
824 const unsigned int grid_degree_input,
834 void build_1D_volume_operator(
835 const dealii::FESystem<1,1> &finite_element,
836 const dealii::Quadrature<1> &quadrature);
839 template<
int dim,
int n_faces,
typename real>
845 const int nstate_input,
846 const unsigned int max_degree_input,
847 const unsigned int grid_degree_input,
849 const double FR_user_specified_correction_parameter_value_input=0.0);
861 void build_1D_volume_operator(
862 const dealii::FESystem<1,1> &finite_element,
863 const dealii::Quadrature<1> &quadrature);
867 template<
int dim,
int n_faces,
typename real>
873 const int nstate_input,
874 const unsigned int max_degree_input,
875 const unsigned int grid_degree_input,
885 void build_1D_volume_operator(
886 const dealii::FESystem<1,1> &finite_element,
887 const dealii::Quadrature<1> &quadrature);
898 template <
int dim,
int n_faces,
typename real>
904 const int nstate_input,
905 const unsigned int max_degree_input,
906 const unsigned int grid_degree_input);
912 void build_1D_gradient_operator(
913 const dealii::FESystem<1,1> &finite_element,
914 const dealii::Quadrature<1> &quadrature);
931 template<
int dim,
int n_faces,
typename real>
937 const int nstate_input,
938 const unsigned int max_degree_input,
939 const unsigned int grid_degree_input);
945 void build_1D_surface_operator(
946 const dealii::FESystem<1,1> &finite_element,
947 const dealii::Quadrature<0> &face_quadrature);
956 template<
int dim,
int n_faces,
typename real>
962 const int nstate_input,
963 const unsigned int max_degree_input,
964 const unsigned int grid_degree_input);
970 void build_local_surface_lifting_operator (
971 const unsigned int n_dofs,
972 const dealii::FullMatrix<double> &norm_matrix,
973 const dealii::FullMatrix<double> &face_integral,
974 dealii::FullMatrix<double> &lifting);
979 void build_1D_volume_operator(
980 const dealii::FESystem<1,1> &finite_element,
981 const dealii::Quadrature<1> &face_quadrature);
984 void build_1D_surface_operator(
985 const dealii::FESystem<1,1> &finite_element,
986 const dealii::Quadrature<0> &face_quadrature);
997 template<
int dim,
int n_faces,
typename real>
1003 const int nstate_input,
1004 const unsigned int max_degree_input,
1005 const unsigned int grid_degree_input,
1007 const double FR_user_specified_correction_parameter_value_input=0.0);
1021 void build_1D_volume_operator(
1022 const dealii::FESystem<1,1> &finite_element,
1023 const dealii::Quadrature<1> &face_quadrature);
1026 void build_1D_surface_operator(
1027 const dealii::FESystem<1,1> &finite_element,
1028 const dealii::Quadrature<0> &face_quadrature);
1045 template<
int dim,
int n_faces,
typename real>
1051 const int nstate_input,
1052 const unsigned int max_degree_input,
1053 const unsigned int grid_degree_input);
1075 void build_1D_shape_functions_at_grid_nodes(
1076 const dealii::FESystem<1,1> &finite_element,
1077 const dealii::Quadrature<1> &quadrature);
1083 void build_1D_shape_functions_at_flux_nodes(
1084 const dealii::FESystem<1,1> &finite_element,
1085 const dealii::Quadrature<1> &quadrature,
1086 const dealii::Quadrature<0> &face_quadrature);
1093 void build_1D_shape_functions_at_volume_flux_nodes(
1094 const dealii::FESystem<1,1> &finite_element,
1095 const dealii::Quadrature<1> &quadrature);
1105 template <
typename real,
int dim,
int n_faces>
1111 const int nstate_input,
1112 const unsigned int max_degree_input,
1113 const unsigned int grid_degree_input,
1114 const bool store_vol_flux_nodes_input =
false,
1115 const bool store_surf_flux_nodes_input =
false,
1116 const bool store_Jacobian_input =
false);
1128 void transform_physical_to_reference(
1129 const dealii::Tensor<1,dim,real> &phys,
1130 const dealii::Tensor<2,dim,real> &metric_cofactor,
1131 dealii::Tensor<1,dim,real> &ref);
1134 void transform_reference_to_physical(
1135 const dealii::Tensor<1,dim,real> &ref,
1136 const dealii::Tensor<2,dim,real> &metric_cofactor,
1137 dealii::Tensor<1,dim,real> &phys);
1140 void transform_physical_to_reference_vector(
1141 const dealii::Tensor<1,dim,std::vector<real>> &phys,
1142 const dealii::Tensor<2,dim,std::vector<real>> &metric_cofactor,
1143 dealii::Tensor<1,dim,std::vector<real>> &ref);
1146 void transform_reference_unit_normal_to_physical_unit_normal(
1147 const unsigned int n_quad_pts,
1148 const dealii::Tensor<1,dim,real> &ref,
1149 const dealii::Tensor<2,dim,std::vector<real>> &metric_cofactor,
1150 std::vector<dealii::Tensor<1,dim,real>> &phys);
1153 void build_determinant_volume_metric_Jacobian(
1154 const unsigned int n_quad_pts,
1155 const unsigned int n_metric_dofs,
1156 const std::array<std::vector<real>,dim> &mapping_support_points,
1164 void build_volume_metric_operators(
1165 const unsigned int n_quad_pts,
1166 const unsigned int n_metric_dofs,
1167 const std::array<std::vector<real>,dim> &mapping_support_points,
1169 const bool use_invariant_curl_form =
false);
1176 void build_facet_metric_operators(
1177 const unsigned int iface,
1178 const unsigned int n_quad_pts,
1179 const unsigned int n_metric_dofs,
1180 const std::array<std::vector<real>,dim> &mapping_support_points,
1182 const bool use_invariant_curl_form =
false);
1211 void build_metric_Jacobian(
1212 const unsigned int n_quad_pts,
1213 const std::array<std::vector<real>,dim> &mapping_support_points,
1214 const dealii::FullMatrix<double> &basis_x_flux_nodes,
1215 const dealii::FullMatrix<double> &basis_y_flux_nodes,
1216 const dealii::FullMatrix<double> &basis_z_flux_nodes,
1217 const dealii::FullMatrix<double> &grad_basis_x_flux_nodes,
1218 const dealii::FullMatrix<double> &grad_basis_y_flux_nodes,
1219 const dealii::FullMatrix<double> &grad_basis_z_flux_nodes,
1220 std::vector<dealii::Tensor<2,dim,real>> &local_Jac);
1228 void build_determinant_metric_Jacobian(
1229 const unsigned int n_quad_pts,
1230 const std::array<std::vector<real>,dim> &mapping_support_points,
1231 const dealii::FullMatrix<double> &basis_x_flux_nodes,
1232 const dealii::FullMatrix<double> &basis_y_flux_nodes,
1233 const dealii::FullMatrix<double> &basis_z_flux_nodes,
1234 const dealii::FullMatrix<double> &grad_basis_x_flux_nodes,
1235 const dealii::FullMatrix<double> &grad_basis_y_flux_nodes,
1236 const dealii::FullMatrix<double> &grad_basis_z_flux_nodes,
1237 std::vector<real> &det_metric_Jac);
1240 void build_local_metric_cofactor_matrix(
1241 const unsigned int n_quad_pts,
1242 const unsigned int n_metric_dofs,
1243 const std::array<std::vector<real>,dim> &mapping_support_points,
1244 const dealii::FullMatrix<double> &basis_x_grid_nodes,
1245 const dealii::FullMatrix<double> &basis_y_grid_nodes,
1246 const dealii::FullMatrix<double> &basis_z_grid_nodes,
1247 const dealii::FullMatrix<double> &basis_x_flux_nodes,
1248 const dealii::FullMatrix<double> &basis_y_flux_nodes,
1249 const dealii::FullMatrix<double> &basis_z_flux_nodes,
1250 const dealii::FullMatrix<double> &grad_basis_x_grid_nodes,
1251 const dealii::FullMatrix<double> &grad_basis_y_grid_nodes,
1252 const dealii::FullMatrix<double> &grad_basis_z_grid_nodes,
1253 const dealii::FullMatrix<double> &grad_basis_x_flux_nodes,
1254 const dealii::FullMatrix<double> &grad_basis_y_flux_nodes,
1255 const dealii::FullMatrix<double> &grad_basis_z_flux_nodes,
1256 dealii::Tensor<2,dim,std::vector<real>> &metric_cofactor,
1257 const bool use_invariant_curl_form =
false);
1281 void compute_local_3D_cofactor(
1282 const unsigned int n_metric_dofs,
1283 const unsigned int n_quad_pts,
1284 const std::array<std::vector<real>,dim> &mapping_support_points,
1285 const dealii::FullMatrix<double> &basis_x_grid_nodes,
1286 const dealii::FullMatrix<double> &basis_y_grid_nodes,
1287 const dealii::FullMatrix<double> &basis_z_grid_nodes,
1288 const dealii::FullMatrix<double> &basis_x_flux_nodes,
1289 const dealii::FullMatrix<double> &basis_y_flux_nodes,
1290 const dealii::FullMatrix<double> &basis_z_flux_nodes,
1291 const dealii::FullMatrix<double> &grad_basis_x_grid_nodes,
1292 const dealii::FullMatrix<double> &grad_basis_y_grid_nodes,
1293 const dealii::FullMatrix<double> &grad_basis_z_grid_nodes,
1294 const dealii::FullMatrix<double> &grad_basis_x_flux_nodes,
1295 const dealii::FullMatrix<double> &grad_basis_y_flux_nodes,
1296 const dealii::FullMatrix<double> &grad_basis_z_flux_nodes,
1297 dealii::Tensor<2,dim,std::vector<real>> &metric_cofactor,
1298 const bool use_invariant_curl_form =
false);
1311 template <
int dim,
int nstate,
int n_faces,
typename real>
1317 const unsigned int max_degree_input,
1318 const unsigned int grid_degree_input);
1335 template <
int dim,
int nstate,
int n_faces,
typename real>
1341 const unsigned int max_degree_input,
1342 const unsigned int grid_degree_input);
1348 void build_1D_volume_state_operator(
1349 const dealii::FESystem<1,1> &finite_element,
1350 const dealii::Quadrature<1> &quadrature);
1353 void build_1D_gradient_state_operator(
1354 const dealii::FESystem<1,1> &finite_element,
1355 const dealii::Quadrature<1> &quadrature);
1358 void build_1D_surface_state_operator(
1359 const dealii::FESystem<1,1> &finite_element,
1360 const dealii::Quadrature<0> &face_quadrature);
1367 template <
int dim,
int nstate,
int n_faces,
typename real>
1373 const unsigned int max_degree_input,
1374 const unsigned int grid_degree_input);
1380 virtual void build_1D_volume_state_operator(
1381 const dealii::FESystem<1,1> &finite_element,
1382 const dealii::Quadrature<1> &quadrature);
1385 void build_1D_gradient_state_operator(
1386 const dealii::FESystem<1,1> &finite_element,
1387 const dealii::Quadrature<1> &quadrature);
1390 void build_1D_surface_state_operator(
1391 const dealii::FESystem<1,1> &finite_element,
1392 const dealii::Quadrature<0> &face_quadrature);
1402 template <
int dim,
int nstate,
int n_faces,
typename real>
1408 const unsigned int max_degree_input,
1409 const unsigned int grid_degree_input);
1415 void build_1D_volume_state_operator(
1416 const dealii::FESystem<1,1> &finite_element,
1417 const dealii::Quadrature<1> &quadrature);
The FLUX basis functions separated by nstate with n shape functions.
unsigned int current_degree
Stores the degree of the current poly degree.
unsigned int current_degree
Stores the degree of the current poly degree.
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
basis_functions< dim, n_faces, real > mapping_shape_functions_grid_nodes
Object of mapping shape functions evaluated at grid nodes.
dealii::Tensor< 2, dim, std::vector< real > > metric_cofactor_vol
The volume metric cofactor matrix.
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
std::array< dealii::FullMatrix< double >, 2 > oneD_surf_operator
Stores the one dimensional surface operator.
The metric independent inverse of the FR mass matrix .
bool store_transpose
Flag is store transpose operator.
unsigned int current_degree
Stores the degree of the current poly degree.
virtual void inner_product(const std::vector< real > &input_vect, const std::vector< real > &weight_vect, std::vector< real > &output_vect, const dealii::FullMatrix< double > &basis_x, const dealii::FullMatrix< double > &basis_y, const dealii::FullMatrix< double > &basis_z, const bool adding=false, const double factor=1.0)=0
virtual function to be defined.
The integration of gradient of solution basis.
unsigned int current_degree
Stores the degree of the current poly degree.
const MPI_Comm mpi_communicator
MPI communicator.
unsigned int current_degree
Stores the degree of the current poly degree.
Sum Factorization derived class.
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
virtual ~OperatorsBase()=default
Destructor.
unsigned int current_grid_degree
Stores the degree of the current grid degree.
const unsigned int max_grid_degree
Max grid degree.
unsigned int current_degree
Stores the degree of the current poly degree.
const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_type
Flux reconstruction parameter type.
dealii::FullMatrix< double > oneD_grad_operator
Stores the one dimensional gradient operator.
dealii::ConditionalOStream pcout
Parallel std::cout that only outputs on mpi_rank==0.
dealii::FullMatrix< double > tensor_product(const dealii::FullMatrix< double > &basis_x, const dealii::FullMatrix< double > &basis_y, const dealii::FullMatrix< double > &basis_z)
Returns the tensor product of matrices passed.
The metric independent FR mass matrix for auxiliary equation .
dealii::FullMatrix< double > tensor_product_state(const int nstate, const dealii::FullMatrix< double > &basis_x, const dealii::FullMatrix< double > &basis_y, const dealii::FullMatrix< double > &basis_z)
Returns the tensor product of matrices passed, but makes it sparse diagonal by state.
bool store_transpose
Flag is store transpose operator.
const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_type
Flux reconstruction parameter type.
const int nstate
Number of states.
double compute_factorial(double n)
Standard function to compute factorial of a number.
Projection operator corresponding to basis functions onto -norm for auxiliary equation.
unsigned int current_degree
Stores the degree of the current poly degree.
Files for the baseline physics.
double FR_param_aux
Flux reconstruction paramater value.
unsigned int current_degree
Stores the degree of the current poly degree.
const unsigned int max_degree
Max polynomial degree.
unsigned int current_degree
Stores the degree of the current poly degree.
std::array< dealii::FullMatrix< double >, 2 > oneD_surf_grad_operator
Stores the one dimensional surface gradient operator.
dealii::Tensor< 2, dim, std::vector< real > > metric_cofactor_surf
The facet metric cofactor matrix, for ONE face.
dealii::FullMatrix< double > oneD_transpose_vol_operator
Stores the transpose of the operator for fast weight-adjusted solves.
unsigned int current_degree
Stores the degree of the current poly degree.
unsigned int current_degree
Stores the degree of the current poly degree.
Flux_Reconstruction
Type of correction in Flux Reconstruction.
dealii::FullMatrix< double > oneD_skew_symm_vol_oper
Skew-symmetric volume operator .
const bool store_Jacobian
Flag if store metric Jacobian at flux nodes.
-th order modal derivative of basis fuctions, ie/
That is Quadrature Weights multiplies with basis_at_vol_cubature.
const bool store_surf_flux_nodes
Flag if store metric Jacobian at flux nodes.
This is the solution basis , the modal differential opertaor commonly seen in DG defined as ...
ESFR correction matrix without jac dependence.
std::vector< real > det_Jac_vol
The determinant of the metric Jacobian at volume cubature nodes.
Local mass matrix without jacobian dependence.
The DG lifting operator is defined as the operator that lifts inner products of polynomials of some o...
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
unsigned int current_degree
Stores the degree of the current poly degree.
unsigned int current_degree
Stores the degree of the current poly degree.
unsigned int current_degree
Stores the degree of the current poly degree.
dealii::FullMatrix< double > oneD_transpose_vol_operator
Stores the transpose of the operator for fast weight-adjusted solves.
dealii::Tensor< 2, dim, std::vector< real > > metric_Jacobian_vol_cubature
Stores the metric Jacobian at flux nodes.
The ESFR lifting operator.
const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_type
Flux reconstruction parameter type.
double FR_param
Flux reconstruction paramater value.
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
unsigned int current_degree
Stores the degree of the current poly degree.
Base metric operators class that stores functions used in both the volume and on surface.
unsigned int max_grid_degree_check
Check to see if the metrics used are a higher order then the initialized grid.
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
ESFR correction matrix for AUX EQUATION without jac dependence.
The mapping shape functions evaluated at the desired nodes (facet set included in volume grid nodes f...
const bool store_vol_flux_nodes
Flag if store metric Jacobian at flux nodes.
OperatorsBase(const int nstate_input, const unsigned int max_degree_input, const unsigned int grid_degree_input)
Constructor.
The basis functions separated by nstate with n shape functions.
The metric independent FR mass matrix .
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
Flux_Reconstruction_Aux
Type of correction in Flux Reconstruction for the auxiliary variables.
unsigned int current_degree
Stores the degree of the current poly degree.
std::vector< real > det_Jac_surf
The determinant of the metric Jacobian at facet cubature nodes.
basis_functions< dim, n_faces, real > mapping_shape_functions_flux_nodes
Object of mapping shape functions evaluated at flux nodes.
unsigned int current_degree
Stores the degree of the current poly degree.
const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_aux_type
Flux reconstruction parameter type.
The metric independent inverse of the FR mass matrix for auxiliary equation .
Projection operator corresponding to basis functions onto -norm.
unsigned int current_degree
Stores the degree of the current poly degree.
const bool store_skew_symmetric_form
Flag to store the skew symmetric form .
std::array< dealii::Tensor< 1, dim, std::vector< real > >, n_faces > flux_nodes_surf
Stores the physical facet flux nodes.
"Stiffness" operator used in DG Strong form.
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
unsigned int current_degree
Stores the degree of the current poly degree.
unsigned int current_degree
Stores the degree of the current poly degree.
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
unsigned int current_degree
Stores the degree of the current poly degree.
dealii::Tensor< 1, dim, std::vector< real > > flux_nodes_vol
Stores the physical volume flux nodes.
The surface integral of test functions.
unsigned int current_degree
Stores the degree of the current poly degree.
Projection operator corresponding to basis functions onto M-norm (L2).
Local stiffness matrix without jacobian dependence.
virtual void matrix_vector_mult(const std::vector< real > &input_vect, std::vector< real > &output_vect, const dealii::FullMatrix< double > &basis_x, const dealii::FullMatrix< double > &basis_y, const dealii::FullMatrix< double > &basis_z, const bool adding=false, const double factor=1.0)=0
virtual function to be defined.
unsigned int current_degree
Stores the degree of the current poly degree.