[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
operators.h
1 #ifndef __OPERATORS_H__
2 #define __OPERATORS_H__
3 
4 #include <deal.II/base/conditional_ostream.h>
5 #include <deal.II/base/parameter_handler.h>
6 
7 #include <deal.II/base/qprojector.h>
8 
9 #include <deal.II/grid/tria.h>
10 
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>
16 
17 #include <deal.II/dofs/dof_handler.h>
18 
19 #include <deal.II/hp/q_collection.h>
20 #include <deal.II/hp/mapping_collection.h>
21 #include <deal.II/hp/fe_values.h>
22 
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>
27 
28 #include <Epetra_RowMatrixTransposer.h>
29 #include <AztecOO.h>
30 
31 #include "ADTypes.hpp"
32 #include <Sacado.hpp>
33 #include <CoDiPack/include/codi.hpp>
34 
35 #include "parameters/all_parameters.h"
36 #include "parameters/parameters.h"
37 
38 namespace PHiLiP {
39 namespace OPERATOR {
40 
42 
52 template <int dim, int n_faces, typename real>
54 {
55 public:
57  virtual ~OperatorsBase() = default;
58 
61  const int nstate_input,//number of states input
62  const unsigned int max_degree_input,//max poly degree for operators
63  const unsigned int grid_degree_input);//max grid degree for operators
64 
66  const unsigned int max_degree;
68  const unsigned int max_grid_degree;
70  const int nstate;
71 
72 protected:
74  unsigned int max_grid_degree_check;
75 
76 public:
77 
79  dealii::FullMatrix<double> tensor_product(
80  const dealii::FullMatrix<double> &basis_x,
81  const dealii::FullMatrix<double> &basis_y,
82  const dealii::FullMatrix<double> &basis_z);
83 
85 
91  dealii::FullMatrix<double> tensor_product_state(
92  const int nstate,
93  const dealii::FullMatrix<double> &basis_x,
94  const dealii::FullMatrix<double> &basis_y,
95  const dealii::FullMatrix<double> &basis_z);
96 
98  double compute_factorial(double n);
99 
101  virtual void matrix_vector_mult(
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;
110  virtual void inner_product(
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;
119 protected:
120 
121  const MPI_Comm mpi_communicator;
122  dealii::ConditionalOStream pcout;
123 };//End of OperatorsBase
124 
126 /* All dim-sized operators are constructed by their one dimensional equivalent, and we use
127 * sum factorization to perform their operations.
128 * Note that we assume tensor product elements in this operators class.
129 */
130 template<int dim, int n_faces, typename real>
131 class SumFactorizedOperators : public OperatorsBase<dim,n_faces,real>
132 {
133 public:
136  const int nstate_input,
137  const unsigned int max_degree_input,
138  const unsigned int grid_degree_input);
139 
141 
147  void matrix_vector_mult(
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;
156 
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);
175 
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);
182 
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);
199 
201 
203  void inner_product(
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;
212 
213 
215 
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);//the only direction that isn't identity
225 
226 
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);
237 
239 
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,//the only direction that isn't identity
252  const std::vector<real> &weights,//vector storing diagonal entries for case not identity
253  const int direction);//direction for the derivative that corresponds to basis
254 
255 
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,//vector of non-zero row indices
261  std::vector<std::array<unsigned int,dim>> &columns);//vector of non-zero column indices
262 
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,//vector of non-zero row indices
268  const std::vector<std::array<unsigned int,dim>> &columns,//vector of non-zero column indices
269  const dealii::FullMatrix<double> &basis,//1D dense basis
270  const std::vector<double> &weights,//Diagonal weights
271  std::array<dealii::FullMatrix<double>,dim> &basis_sparse);//Sparse basis
272 
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,//vector of non-zero row indices
278  std::vector<unsigned int> &columns,//vector of non-zero column indices
279  const int dim_not_zero);//ref direction face is on
280 
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,//vector of non-zero row indices
286  const std::vector<unsigned int> &columns,//vector of non-zero column indices
287  const dealii::FullMatrix<double> &basis,//1D dense basis
288  const std::vector<double> &weights,//Diagonal weights
289  dealii::FullMatrix<double> &basis_sparse,//Sparse basis
290  const int dim_not_zero);//ref direction face is on
291 
293 
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);
302 
304  /* This is for the case where the operator of size dim is the dyadic product of
305  * the same 1D operator in each direction
306  */
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);
314 
316 
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,//only 2 faces in 1D
327  const dealii::FullMatrix<double> &basis_vol,
328  const bool adding = false,
329  const double factor = 1.0);
330 
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,//only 2 faces in 1D
338  const dealii::FullMatrix<double> &basis_vol,
339  const bool adding = false,
340  const double factor = 1.0);
341 
342 
344 
347  void Hadamard_product(
348  const dealii::FullMatrix<real> &input_mat1,
349  const dealii::FullMatrix<real> &input_mat2,
350  dealii::FullMatrix<real> &output_mat);
351 
352 //protected:
353 public:
355  dealii::FullMatrix<double> oneD_vol_operator;
356 
358 
360  std::array<dealii::FullMatrix<double>,2> oneD_surf_operator;
361 
363  dealii::FullMatrix<double> oneD_grad_operator;
364 
366  std::array<dealii::FullMatrix<double>,2> oneD_surf_grad_operator;
367 
368 };//End of SumFactorizedOperators Class
369 
370 /************************************************************************
371 *
372 * VOLUME OPERATORS
373 *
374 ************************************************************************/
375 
377 /* This class stores the basis functions evaluated at volume and facet
378 * cubature nodes, as well as it's gradient in REFERENCE space.
379 */
380 template<int dim, int n_faces, typename real>
381 class basis_functions : public SumFactorizedOperators<dim,n_faces,real>
382 {
383 public:
386  const int nstate_input,
387  const unsigned int max_degree_input,
388  const unsigned int grid_degree_input);
389 
391  unsigned int current_degree;
392 
394  void build_1D_volume_operator(
395  const dealii::FESystem<1,1> &finite_element,
396  const dealii::Quadrature<1> &quadrature);
397 
399  void build_1D_gradient_operator(
400  const dealii::FESystem<1,1> &finite_element,
401  const dealii::Quadrature<1> &quadrature);
402 
404  void build_1D_surface_operator(
405  const dealii::FESystem<1,1> &finite_element,
406  const dealii::Quadrature<0> &quadrature);
407 
409  void build_1D_surface_gradient_operator(
410  const dealii::FESystem<1,1> &finite_element,
411  const dealii::Quadrature<0> &quadrature);
412 };
413 
415 template<int dim, int n_faces, typename real>
416 class vol_integral_basis : public SumFactorizedOperators<dim,n_faces,real>
417 {
418 public:
421  const int nstate_input,
422  const unsigned int max_degree_input,
423  const unsigned int grid_degree_input);
424 
426  unsigned int current_degree;
427 
429  void build_1D_volume_operator(
430  const dealii::FESystem<1,1> &finite_element,
431  const dealii::Quadrature<1> &quadrature);
432 };
433 
435 template<int dim, int n_faces, typename real>
436 class local_mass : public SumFactorizedOperators<dim,n_faces,real>
437 {
438 public:
440  local_mass (
441  const int nstate_input,
442  const unsigned int max_degree_input,
443  const unsigned int grid_degree_input);
444 
446  unsigned int current_degree;
447 
449  void build_1D_volume_operator(
450  const dealii::FESystem<1,1> &finite_element,
451  const dealii::Quadrature<1> &quadrature); //override;
452 
454 
457  dealii::FullMatrix<double> build_dim_mass_matrix(
458  const int nstate,
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);
463 };
464 
466 
471 template<int dim, int n_faces, typename real>
472 class local_basis_stiffness : public SumFactorizedOperators<dim,n_faces,real>
473 {
474 public:
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);
481 
483  unsigned int current_degree;
484 
487 
489  void build_1D_volume_operator(
490  const dealii::FESystem<1,1> &finite_element,
491  const dealii::Quadrature<1> &quadrature);
492 
494  dealii::FullMatrix<double> oneD_skew_symm_vol_oper;
495 };
496 
498 template<int dim, int n_faces, typename real>
500 {
501 public:
504  const int nstate_input,
505  const unsigned int max_degree_input,
506  const unsigned int grid_degree_input);
507 
509  unsigned int current_degree;
510 
512  void build_1D_volume_operator(
513  const dealii::FESystem<1,1> &finite_element,
514  const dealii::Quadrature<1> &quadrature);
515 };
516 
518 template<int dim, int n_faces, typename real>
519 class derivative_p : public SumFactorizedOperators<dim,n_faces,real>
520 {
521 public:
523  derivative_p (
524  const int nstate_input,
525  const unsigned int max_degree_input,
526  const unsigned int grid_degree_input);
527 
529  unsigned int current_degree;
530 
532  void build_1D_volume_operator(
533  const dealii::FESystem<1,1> &finite_element,
534  const dealii::Quadrature<1> &quadrature);
535 };
536 
538 template<int dim, int n_faces, typename real>
540 {
541 public:
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);
549 
552 
555 
557  unsigned int current_degree;
558 
560  double FR_param;
561 
563  /* This parameter recovers Huynh, Hung T. "A flux reconstruction approach to high-order schemes including discontinuous Galerkin methods." 18th AIAA computational fluid dynamics conference. 2007.
564  */
565  void get_Huynh_g2_parameter (
566  const unsigned int curr_cell_degree,
567  double &c);
568 
570 
572  void get_spectral_difference_parameter (
573  const unsigned int curr_cell_degree,
574  double &c);
575 
577 
579  void get_c_negative_FR_parameter (
580  const unsigned int curr_cell_degree,
581  double &c);
582 
584 
587  void get_c_negative_divided_by_two_FR_parameter (
588  const unsigned int curr_cell_degree,
589  double &c);
590 
592 
596  void get_c_plus_parameter (
597  const unsigned int curr_cell_degree,
598  double &c);
599 
601 
608  void get_FR_correction_parameter (
609  const unsigned int curr_cell_degree,
610  double &c);
611 
613 
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,
620  const double c,
621  dealii::FullMatrix<double> &Flux_Reconstruction_operator);
622 
624  void build_1D_volume_operator(
625  const dealii::FESystem<1,1> &finite_element,
626  const dealii::Quadrature<1> &quadrature);
627 
629 
634  dealii::FullMatrix<double> build_dim_Flux_Reconstruction_operator(
635  const dealii::FullMatrix<double> &local_Mass_Matrix,
636  const int nstate,
637  const unsigned int n_dofs);
638 
640 
646  dealii::FullMatrix<double> build_dim_Flux_Reconstruction_operator_directly(
647  const int nstate,
648  const unsigned int n_dofs,
649  dealii::FullMatrix<double> &pth_deriv,
650  dealii::FullMatrix<double> &mass_matrix);
651 };
652 
654 
658 template<int dim, int n_faces, typename real>
660 {
661 public:
664  const int nstate_input,
665  const unsigned int max_degree_input,
666  const unsigned int grid_degree_input,
667  const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_aux_input);
668 
670  unsigned int current_degree;
671 
674 
676  double FR_param_aux;
677 
679 
686  void get_FR_aux_correction_parameter (
687  const unsigned int curr_cell_degree,
688  double &k);
689 
691  void build_1D_volume_operator(
692  const dealii::FESystem<1,1> &finite_element,
693  const dealii::Quadrature<1> &quadrature);
694 };
695 
697 template<int dim, int n_faces, typename real>
698 class vol_projection_operator : public SumFactorizedOperators<dim,n_faces,real>
699 {
700 public:
703  const int nstate_input,
704  const unsigned int max_degree_input,
705  const unsigned int grid_degree_input);
706 
708  unsigned int current_degree;
709 
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);
715 
717  void build_1D_volume_operator(
718  const dealii::FESystem<1,1> &finite_element,
719  const dealii::Quadrature<1> &quadrature);
720 };
721 
723 template<int dim, int n_faces, typename real>
724 class vol_projection_operator_FR : public vol_projection_operator<dim,n_faces,real>
725 {
726 public:
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);
735 
738 
741 
744 
746  unsigned int current_degree;
747 
749  void build_1D_volume_operator(
750  const dealii::FESystem<1,1> &finite_element,
751  const dealii::Quadrature<1> &quadrature);
752 
754  dealii::FullMatrix<double> oneD_transpose_vol_operator;
755 };
756 
758 template<int dim, int n_faces, typename real>
760 {
761 public:
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);
769 
771  unsigned int current_degree;
772 
775 
778 
780  void build_1D_volume_operator(
781  const dealii::FESystem<1,1> &finite_element,
782  const dealii::Quadrature<1> &quadrature);
783 
785  dealii::FullMatrix<double> oneD_transpose_vol_operator;
786 };
787 
789 template<int dim, int n_faces, typename real>
790 class FR_mass_inv : public SumFactorizedOperators<dim,n_faces,real>
791 {
792 public:
794  FR_mass_inv (
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);
800 
802  unsigned int current_degree;
803 
806 
809 
811  void build_1D_volume_operator(
812  const dealii::FESystem<1,1> &finite_element,
813  const dealii::Quadrature<1> &quadrature);
814 };
816 template<int dim, int n_faces, typename real>
817 class FR_mass_inv_aux : public SumFactorizedOperators<dim,n_faces,real>
818 {
819 public:
822  const int nstate_input,
823  const unsigned int max_degree_input,
824  const unsigned int grid_degree_input,
826 
828  unsigned int current_degree;
829 
832 
834  void build_1D_volume_operator(
835  const dealii::FESystem<1,1> &finite_element,
836  const dealii::Quadrature<1> &quadrature);
837 };
839 template<int dim, int n_faces, typename real>
840 class FR_mass : public SumFactorizedOperators<dim,n_faces,real>
841 {
842 public:
844  FR_mass (
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);
850 
852  unsigned int current_degree;
853 
856 
859 
861  void build_1D_volume_operator(
862  const dealii::FESystem<1,1> &finite_element,
863  const dealii::Quadrature<1> &quadrature);
864 };
865 
867 template<int dim, int n_faces, typename real>
868 class FR_mass_aux : public SumFactorizedOperators<dim,n_faces,real>
869 {
870 public:
872  FR_mass_aux (
873  const int nstate_input,
874  const unsigned int max_degree_input,
875  const unsigned int grid_degree_input,
877 
879  unsigned int current_degree;
880 
883 
885  void build_1D_volume_operator(
886  const dealii::FESystem<1,1> &finite_element,
887  const dealii::Quadrature<1> &quadrature);
888 };
889 
891 
898 template <int dim, int n_faces, typename real>
899 class vol_integral_gradient_basis : public SumFactorizedOperators<dim,n_faces,real>
900 {
901 public:
904  const int nstate_input,
905  const unsigned int max_degree_input,
906  const unsigned int grid_degree_input);
907 
909  unsigned int current_degree;
910 
912  void build_1D_gradient_operator(
913  const dealii::FESystem<1,1> &finite_element,
914  const dealii::Quadrature<1> &quadrature);
915 };
916 
917 /************************************************************************
918 *
919 * SURFACE OPERATORS
920 *
921 ************************************************************************/
922 
923 
925 
931 template<int dim, int n_faces, typename real>
932 class face_integral_basis : public SumFactorizedOperators<dim,n_faces,real>
933 {
934 public:
937  const int nstate_input,
938  const unsigned int max_degree_input,
939  const unsigned int grid_degree_input);
940 
942  unsigned int current_degree;
943 
945  void build_1D_surface_operator(
946  const dealii::FESystem<1,1> &finite_element,
947  const dealii::Quadrature<0> &face_quadrature);
948 };
949 
951 
956 template<int dim, int n_faces, typename real>
957 class lifting_operator : public SumFactorizedOperators<dim,n_faces,real>
958 {
959 public:
962  const int nstate_input,
963  const unsigned int max_degree_input,
964  const unsigned int grid_degree_input);
965 
967  unsigned int current_degree;
968 
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);
975 
977 
979  void build_1D_volume_operator(
980  const dealii::FESystem<1,1> &finite_element,
981  const dealii::Quadrature<1> &face_quadrature);
982 
984  void build_1D_surface_operator(
985  const dealii::FESystem<1,1> &finite_element,
986  const dealii::Quadrature<0> &face_quadrature);
987 };
988 
990 
997 template<int dim, int n_faces, typename real>
998 class lifting_operator_FR : public lifting_operator<dim,n_faces,real>
999 {
1000 public:
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);
1008 
1010  unsigned int current_degree;
1011 
1014 
1017 
1019 
1021  void build_1D_volume_operator(
1022  const dealii::FESystem<1,1> &finite_element,
1023  const dealii::Quadrature<1> &face_quadrature);
1024 
1026  void build_1D_surface_operator(
1027  const dealii::FESystem<1,1> &finite_element,
1028  const dealii::Quadrature<0> &face_quadrature);
1029 };
1030 
1031 
1032 /************************************************************************
1033 *
1034 * METRIC MAPPING OPERATORS
1035 *
1036 ************************************************************************/
1037 
1038 
1040 
1045 template<int dim, int n_faces, typename real>
1046 class mapping_shape_functions: public SumFactorizedOperators<dim,n_faces,real>
1047 {
1048 public:
1051  const int nstate_input,
1052  const unsigned int max_degree_input,
1053  const unsigned int grid_degree_input);
1054 
1056  unsigned int current_degree;
1057 
1059  unsigned int current_grid_degree;
1060 
1063 
1066 
1068 
1075  void build_1D_shape_functions_at_grid_nodes(
1076  const dealii::FESystem<1,1> &finite_element,
1077  const dealii::Quadrature<1> &quadrature);
1078 
1080 
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);
1087 
1089 
1093  void build_1D_shape_functions_at_volume_flux_nodes(
1094  const dealii::FESystem<1,1> &finite_element,
1095  const dealii::Quadrature<1> &quadrature);
1096 
1097 };
1098 
1099 /*****************************************************************************
1100 *
1101 * METRIC OPERTAORS TO BE CALLED ON-THE-FLY
1102 *
1103 *****************************************************************************/
1105 template <typename real, int dim, int n_faces>
1106 class metric_operators: public SumFactorizedOperators<dim,n_faces,real>
1107 {
1108 public:
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);
1117 
1119  const bool store_Jacobian;
1120 
1123 
1126 
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);
1132 
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);
1138 
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);
1144 
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);
1151 
1153  void build_determinant_volume_metric_Jacobian(
1154  const unsigned int n_quad_pts,//number volume quad pts
1155  const unsigned int n_metric_dofs,//dofs of metric basis. NOTE: this is the number of mapping support points
1156  const std::array<std::vector<real>,dim> &mapping_support_points,
1158 
1160 
1164  void build_volume_metric_operators(
1165  const unsigned int n_quad_pts,//number volume quad pts
1166  const unsigned int n_metric_dofs,//dofs of metric basis. NOTE: this is the number of mapping support points
1167  const std::array<std::vector<real>,dim> &mapping_support_points,
1169  const bool use_invariant_curl_form = false);
1170 
1172 
1176  void build_facet_metric_operators(
1177  const unsigned int iface,
1178  const unsigned int n_quad_pts,//number facet quad pts
1179  const unsigned int n_metric_dofs,//dofs of metric basis. NOTE: this is the number of mapping support points
1180  const std::array<std::vector<real>,dim> &mapping_support_points,
1182  const bool use_invariant_curl_form = false);
1183 
1185  dealii::Tensor<2,dim,std::vector<real>> metric_cofactor_vol;
1186 
1188  dealii::Tensor<2,dim,std::vector<real>> metric_cofactor_surf;
1189 
1191  std::vector<real> det_Jac_vol;
1192 
1194  std::vector<real> det_Jac_surf;
1195 
1197  dealii::Tensor<2,dim,std::vector<real>> metric_Jacobian_vol_cubature;
1198 
1200  dealii::Tensor<1,dim,std::vector<real>> flux_nodes_vol;
1201 
1203  std::array<dealii::Tensor<1,dim,std::vector<real>>,n_faces> flux_nodes_surf;
1204 
1205 protected:
1206 
1208 
1211  void build_metric_Jacobian(
1212  const unsigned int n_quad_pts,//the dim sized n_quad_pts, NOT the 1D
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);
1221 
1223 
1228  void build_determinant_metric_Jacobian(
1229  const unsigned int n_quad_pts,//number volume 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);
1238 
1240  void build_local_metric_cofactor_matrix(
1241  const unsigned int n_quad_pts,//number volume quad pts
1242  const unsigned int n_metric_dofs,//dofs of metric basis. NOTE: this is the number of mapping support points
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);
1258 
1260 
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);
1299 };
1300 
1301 /************************************************************
1302 *
1303 * SUMFACTORIZED STATE
1304 *
1305 ************************************************************/
1306 
1308 
1311 template <int dim, int nstate, int n_faces, typename real>
1312 class SumFactorizedOperatorsState : public SumFactorizedOperators<dim,n_faces,real>
1313 {
1314 public:
1317  const unsigned int max_degree_input,
1318  const unsigned int grid_degree_input);
1319 
1321  std::array<dealii::FullMatrix<double>,nstate> oneD_vol_state_operator;
1322 
1324  std::array<std::array<dealii::FullMatrix<double>,2>,nstate> oneD_surf_state_operator;
1325 
1327  std::array<dealii::FullMatrix<double>,nstate> oneD_grad_state_operator;
1328 
1330  std::array<std::array<dealii::FullMatrix<double>,2>,nstate> oneD_surf_grad_state_operator;
1331 
1332 };//end of OperatorsBaseState Class
1333 
1335 template <int dim, int nstate, int n_faces, typename real>
1336 class basis_functions_state : public SumFactorizedOperatorsState<dim,nstate,n_faces,real>
1337 {
1338 public:
1341  const unsigned int max_degree_input,
1342  const unsigned int grid_degree_input);
1343 
1345  unsigned int current_degree;
1346 
1348  void build_1D_volume_state_operator(
1349  const dealii::FESystem<1,1> &finite_element,
1350  const dealii::Quadrature<1> &quadrature);
1351 
1353  void build_1D_gradient_state_operator(
1354  const dealii::FESystem<1,1> &finite_element,
1355  const dealii::Quadrature<1> &quadrature);
1356 
1358  void build_1D_surface_state_operator(
1359  const dealii::FESystem<1,1> &finite_element,
1360  const dealii::Quadrature<0> &face_quadrature);
1361 };
1362 
1364 
1367 template <int dim, int nstate, int n_faces, typename real>
1368 class flux_basis_functions_state : public SumFactorizedOperatorsState<dim,nstate,n_faces,real>
1369 {
1370 public:
1373  const unsigned int max_degree_input,
1374  const unsigned int grid_degree_input);
1375 
1377  unsigned int current_degree;
1378 
1380  virtual void build_1D_volume_state_operator(
1381  const dealii::FESystem<1,1> &finite_element,
1382  const dealii::Quadrature<1> &quadrature);
1383 
1385  void build_1D_gradient_state_operator(
1386  const dealii::FESystem<1,1> &finite_element,
1387  const dealii::Quadrature<1> &quadrature);
1388 
1390  void build_1D_surface_state_operator(
1391  const dealii::FESystem<1,1> &finite_element,
1392  const dealii::Quadrature<0> &face_quadrature);
1393 };
1394 
1396 
1402 template <int dim, int nstate, int n_faces, typename real>
1403 class local_flux_basis_stiffness : public flux_basis_functions_state<dim,nstate,n_faces,real>
1404 {
1405 public:
1408  const unsigned int max_degree_input,
1409  const unsigned int grid_degree_input);
1410 
1412  unsigned int current_degree;
1413 
1415  void build_1D_volume_state_operator(
1416  const dealii::FESystem<1,1> &finite_element,//pass the finite element of the TEST FUNCTION
1417  const dealii::Quadrature<1> &quadrature);
1418 };
1419 
1420 }
1421 }
1422 
1423 #endif
1424 
The FLUX basis functions separated by nstate with n shape functions.
Definition: operators.h:1368
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:746
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:1345
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
Definition: operators.h:554
basis_functions< dim, n_faces, real > mapping_shape_functions_grid_nodes
Object of mapping shape functions evaluated at grid nodes.
Definition: operators.h:1062
dealii::Tensor< 2, dim, std::vector< real > > metric_cofactor_vol
The volume metric cofactor matrix.
Definition: operators.h:1185
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
Definition: operators.h:740
std::array< dealii::FullMatrix< double >, 2 > oneD_surf_operator
Stores the one dimensional surface operator.
Definition: operators.h:360
The metric independent inverse of the FR mass matrix .
Definition: operators.h:790
bool store_transpose
Flag is store transpose operator.
Definition: operators.h:737
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:426
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.
Definition: operators.h:899
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:828
const MPI_Comm mpi_communicator
MPI communicator.
Definition: operators.h:121
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:557
Sum Factorization derived class.
Definition: operators.h:131
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
Definition: operators.h:858
virtual ~OperatorsBase()=default
Destructor.
unsigned int current_grid_degree
Stores the degree of the current grid degree.
Definition: operators.h:1059
const unsigned int max_grid_degree
Max grid degree.
Definition: operators.h:68
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:879
const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_type
Flux reconstruction parameter type.
Definition: operators.h:777
dealii::FullMatrix< double > oneD_grad_operator
Stores the one dimensional gradient operator.
Definition: operators.h:363
dealii::ConditionalOStream pcout
Parallel std::cout that only outputs on mpi_rank==0.
Definition: operators.h:122
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.
Definition: operators.cpp:55
The metric independent FR mass matrix for auxiliary equation .
Definition: operators.h:868
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.
Definition: operators.cpp:106
bool store_transpose
Flag is store transpose operator.
Definition: operators.h:774
const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_type
Flux reconstruction parameter type.
Definition: operators.h:831
const int nstate
Number of states.
Definition: operators.h:70
double compute_factorial(double n)
Standard function to compute factorial of a number.
Definition: operators.cpp:173
Projection operator corresponding to basis functions onto -norm for auxiliary equation.
Definition: operators.h:759
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:909
Files for the baseline physics.
Definition: ADTypes.hpp:10
double FR_param_aux
Flux reconstruction paramater value.
Definition: operators.h:676
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:1056
const unsigned int max_degree
Max polynomial degree.
Definition: operators.h:66
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:942
std::array< dealii::FullMatrix< double >, 2 > oneD_surf_grad_operator
Stores the one dimensional surface gradient operator.
Definition: operators.h:366
dealii::Tensor< 2, dim, std::vector< real > > metric_cofactor_surf
The facet metric cofactor matrix, for ONE face.
Definition: operators.h:1188
dealii::FullMatrix< double > oneD_transpose_vol_operator
Stores the transpose of the operator for fast weight-adjusted solves.
Definition: operators.h:785
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:1010
In order to have all state operators be arrays of array, we template by dim, type, nstate, and number of faces.
Definition: operators.h:1312
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:967
Flux_Reconstruction
Type of correction in Flux Reconstruction.
dealii::FullMatrix< double > oneD_skew_symm_vol_oper
Skew-symmetric volume operator .
Definition: operators.h:494
const bool store_Jacobian
Flag if store metric Jacobian at flux nodes.
Definition: operators.h:1119
-th order modal derivative of basis fuctions, ie/
Definition: operators.h:519
That is Quadrature Weights multiplies with basis_at_vol_cubature.
Definition: operators.h:416
const bool store_surf_flux_nodes
Flag if store metric Jacobian at flux nodes.
Definition: operators.h:1125
This is the solution basis , the modal differential opertaor commonly seen in DG defined as ...
Definition: operators.h:499
ESFR correction matrix without jac dependence.
Definition: operators.h:539
std::vector< real > det_Jac_vol
The determinant of the metric Jacobian at volume cubature nodes.
Definition: operators.h:1191
Local mass matrix without jacobian dependence.
Definition: operators.h:436
The DG lifting operator is defined as the operator that lifts inner products of polynomials of some o...
Definition: operators.h:957
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
Definition: operators.h:743
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:391
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:771
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:1412
dealii::FullMatrix< double > oneD_transpose_vol_operator
Stores the transpose of the operator for fast weight-adjusted solves.
Definition: operators.h:754
dealii::Tensor< 2, dim, std::vector< real > > metric_Jacobian_vol_cubature
Stores the metric Jacobian at flux nodes.
Definition: operators.h:1197
The ESFR lifting operator.
Definition: operators.h:998
const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_type
Flux reconstruction parameter type.
Definition: operators.h:882
double FR_param
Flux reconstruction paramater value.
Definition: operators.h:560
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
Definition: operators.h:808
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:670
Base metric operators class that stores functions used in both the volume and on surface.
Definition: operators.h:1106
unsigned int max_grid_degree_check
Check to see if the metrics used are a higher order then the initialized grid.
Definition: operators.h:74
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
Definition: operators.h:855
ESFR correction matrix for AUX EQUATION without jac dependence.
Definition: operators.h:659
The mapping shape functions evaluated at the desired nodes (facet set included in volume grid nodes f...
Definition: operators.h:1046
const bool store_vol_flux_nodes
Flag if store metric Jacobian at flux nodes.
Definition: operators.h:1122
OperatorsBase(const int nstate_input, const unsigned int max_degree_input, const unsigned int grid_degree_input)
Constructor.
Definition: operators.cpp:42
The basis functions separated by nstate with n shape functions.
Definition: operators.h:1336
std::array< dealii::FullMatrix< double >, nstate > oneD_vol_state_operator
Stores the one dimensional volume operator.
Definition: operators.h:1321
The metric independent FR mass matrix .
Definition: operators.h:840
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
Definition: operators.h:551
std::array< std::array< dealii::FullMatrix< double >, 2 >, nstate > oneD_surf_grad_state_operator
Stores the one dimensional surface gradient operator.
Definition: operators.h:1330
dealii::FullMatrix< double > oneD_vol_operator
Stores the one dimensional volume operator.
Definition: operators.h:355
const double FR_user_specified_correction_parameter_value
User specified flux recontruction correction parameter value.
Definition: operators.h:1016
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.
Definition: operators.h:802
std::vector< real > det_Jac_surf
The determinant of the metric Jacobian at facet cubature nodes.
Definition: operators.h:1194
basis_functions< dim, n_faces, real > mapping_shape_functions_flux_nodes
Object of mapping shape functions evaluated at flux nodes.
Definition: operators.h:1065
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:708
const Parameters::AllParameters::Flux_Reconstruction_Aux FR_param_aux_type
Flux reconstruction parameter type.
Definition: operators.h:673
The metric independent inverse of the FR mass matrix for auxiliary equation .
Definition: operators.h:817
std::array< std::array< dealii::FullMatrix< double >, 2 >, nstate > oneD_surf_state_operator
Stores the one dimensional surface operator.
Definition: operators.h:1324
Projection operator corresponding to basis functions onto -norm.
Definition: operators.h:724
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:446
const bool store_skew_symmetric_form
Flag to store the skew symmetric form .
Definition: operators.h:486
std::array< dealii::Tensor< 1, dim, std::vector< real > >, n_faces > flux_nodes_surf
Stores the physical facet flux nodes.
Definition: operators.h:1203
"Stiffness" operator used in DG Strong form.
Definition: operators.h:1403
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
Definition: operators.h:1013
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:852
std::array< dealii::FullMatrix< double >, nstate > oneD_grad_state_operator
Stores the one dimensional gradient operator.
Definition: operators.h:1327
Operator base class.
Definition: operators.h:53
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:509
const Parameters::AllParameters::Flux_Reconstruction FR_param_type
Flux reconstruction parameter type.
Definition: operators.h:805
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:1377
dealii::Tensor< 1, dim, std::vector< real > > flux_nodes_vol
Stores the physical volume flux nodes.
Definition: operators.h:1200
The surface integral of test functions.
Definition: operators.h:932
unsigned int current_degree
Stores the degree of the current poly degree.
Definition: operators.h:529
Projection operator corresponding to basis functions onto M-norm (L2).
Definition: operators.h:698
Local stiffness matrix without jacobian dependence.
Definition: operators.h:472
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.
Definition: operators.h:483