[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
GCL_Collocated_test.cpp
1 #include <iomanip>
2 #include <cmath>
3 #include <limits>
4 #include <type_traits>
5 #include <math.h>
6 #include <iostream>
7 #include <stdlib.h>
8 
9 #include <deal.II/distributed/solution_transfer.h>
10 
11 #include "testing/tests.h"
12 
13 #include<fstream>
14 #include <deal.II/base/parameter_handler.h>
15 #include <deal.II/base/tensor.h>
16 #include <deal.II/numerics/vector_tools.h>
17 
18 #include <deal.II/grid/grid_generator.h>
19 #include <deal.II/grid/grid_refinement.h>
20 #include <deal.II/grid/grid_tools.h>
21 #include <deal.II/grid/grid_out.h>
22 #include <deal.II/grid/grid_in.h>
23 
24 #include <deal.II/dofs/dof_handler.h>
25 #include <deal.II/dofs/dof_tools.h>
26 #include <deal.II/dofs/dof_renumbering.h>
27 
28 #include <deal.II/dofs/dof_accessor.h>
29 
30 #include <deal.II/lac/vector.h>
31 #include <deal.II/lac/dynamic_sparsity_pattern.h>
32 #include <deal.II/lac/sparse_matrix.h>
33 
34 #include <deal.II/meshworker/dof_info.h>
35 
36 // Finally, we take our exact solution from the library as well as volume_quadrature
37 // and additional tools.
38 #include <deal.II/fe/mapping_q.h>
39 #include <deal.II/grid/manifold_lib.h>
40 #include <deal.II/numerics/data_out.h>
41 #include <deal.II/numerics/data_out_dof_data.h>
42 #include <deal.II/numerics/vector_tools.h>
43 #include <deal.II/numerics/vector_tools.templates.h>
44 
45 #include "dg/dg_base.hpp"
46 #include "dg/dg_factory.hpp"
47 #include "operators/operators.h"
48 #include "parameters/all_parameters.h"
49 #include "parameters/parameters.h"
50 //#include <GCL_test.h>
51 
52 const double TOLERANCE = 1E-6;
53 using namespace std;
54 //namespace PHiLiP {
55 
56 
57 template <int dim>
58 class CurvManifold: public dealii::ChartManifold<dim,dim,dim> {
59  virtual dealii::Point<dim> pull_back(const dealii::Point<dim> &space_point) const override;
60  virtual dealii::Point<dim> push_forward(const dealii::Point<dim> &chart_point) const override;
61  virtual dealii::DerivativeForm<1,dim,dim> push_forward_gradient(const dealii::Point<dim> &chart_point) const override;
62 
63  virtual std::unique_ptr<dealii::Manifold<dim,dim> > clone() const override;
64 };
65 
66 template<int dim>
67 dealii::Point<dim> CurvManifold<dim>::pull_back(const dealii::Point<dim> &space_point) const
68 {
69  using namespace PHiLiP;
70  const double pi = atan(1)*4.0;
71  dealii::Point<dim> x_ref;
72  dealii::Point<dim> x_phys;
73  for(int idim=0; idim<dim; idim++){
74  x_ref[idim] = space_point[idim];
75  x_phys[idim] = space_point[idim];
76  }
77  dealii::Vector<double> function(dim);
78  dealii::FullMatrix<double> derivative(dim);
79  double beta =1.0/10.0;
80  double alpha =1.0/10.0;
81  int flag =0;
82  while(flag != dim){
83  if(dim==2){
84  function[0] = x_ref[0] - x_phys[0] +beta*std::cos(pi/2.0*x_ref[0])*std::cos(3.0*pi/2.0*x_ref[1]);
85  function[1] = x_ref[1] - x_phys[1] +beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
86  }
87  else{
88  function[0] = x_ref[0] - x_phys[0] +alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
89  function[1] = x_ref[1] - x_phys[1] +alpha*exp(1.0-x_ref[1])*(std::sin(pi * x_ref[0]) + std::sin(pi* x_ref[2]));
90  function[2] = x_ref[2] - x_phys[2] +1.0/20.0*( std::sin(2.0 * pi * x_ref[0]) + std::sin(2.0 * pi * x_ref[1]));
91  }
92 
93  if(dim==2){
94  derivative[0][0] = 1.0 - beta* pi/2.0 * std::sin(pi/2.0*x_ref[0])*std::cos(3.0*pi/2.0*x_ref[1]);
95  derivative[0][1] = - beta*3.0 *pi/2.0 * std::cos(pi/2.0*x_ref[0])*std::sin(3.0*pi/2.0*x_ref[1]);
96 
97  derivative[1][0] = beta*2.0*pi*std::cos(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
98  derivative[1][1] = 1.0 -beta*pi/2.0*std::sin(2.0*pi*(x_ref[0]))*std::sin(pi/2.0*x_ref[1]);
99  }
100  else{
101  derivative[0][0] = 1.0;
102  derivative[0][1] = - alpha*pi*std::sin(pi*x_ref[1]);
103  derivative[0][2] = - alpha*pi*std::sin(pi*x_ref[2]);
104 
105  derivative[1][0] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[0]);
106  derivative[1][1] = 1.0 -alpha*exp(1.0-x_ref[1])*(std::sin(pi*x_ref[0])+std::sin(pi*x_ref[2]));
107  derivative[1][2] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[2]);
108  derivative[2][0] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[0]);
109  derivative[2][1] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[1]);
110  derivative[2][2] = 1.0;
111  }
112 
113  dealii::FullMatrix<double> Jacobian_inv(dim);
114  Jacobian_inv.invert(derivative);
115  dealii::Vector<double> Newton_Step(dim);
116  Jacobian_inv.vmult(Newton_Step, function);
117  for(int idim=0; idim<dim; idim++){
118  x_ref[idim] -= Newton_Step[idim];
119  }
120  flag=0;
121  for(int idim=0; idim<dim; idim++){
122  if(std::abs(function[idim]) < 1e-15)
123  flag++;
124  }
125  if(flag == dim)
126  break;
127  }
128  std::vector<double> function_check(dim);
129  if(dim==2){
130  function_check[0] = x_ref[0] + beta*std::cos(pi/2.0*x_ref[0])*std::cos(3.0*pi/2.0*x_ref[1]);
131  function_check[1] = x_ref[1] + beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
132  }
133  else{
134  function_check[0] = x_ref[0] +alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
135  function_check[1] = x_ref[1] +alpha*exp(1.0-x_ref[1])*(std::sin(pi * x_ref[0]) + std::sin(pi* x_ref[2]));
136  function_check[2] = x_ref[2] +1.0/20.0*( std::sin(2.0 * pi * x_ref[0]) + std::sin(2.0 * pi * x_ref[1]));
137  }
138  std::vector<double> error(dim);
139  for(int idim=0; idim<dim; idim++)
140  error[idim] = std::abs(function_check[idim] - x_phys[idim]);
141  if (error[0] > 1e-13) {
142  std::cout << "Large error " << error[0] << std::endl;
143  for(int idim=0;idim<dim; idim++)
144  std::cout << "dim " << idim << " xref " << x_ref[idim] << " x_phys " << x_phys[idim] << " function Check " << function_check[idim] << " Error " << error[idim] << " Flag " << flag << std::endl;
145  }
146 
147  return x_ref;
148 }
149 
150 template<int dim>
151 dealii::Point<dim> CurvManifold<dim>::push_forward(const dealii::Point<dim> &chart_point) const
152 {
153  const double pi = atan(1)*4.0;
154 
155  dealii::Point<dim> x_ref;
156  dealii::Point<dim> x_phys;
157  for(int idim=0; idim<dim; idim++)
158  x_ref[idim] = chart_point[idim];
159  double beta = 1.0/10.0;
160  double alpha = 1.0/10.0;
161  if(dim==2){
162  x_phys[0] = x_ref[0] + beta*std::cos(pi/2.0*x_ref[0])*std::cos(3.0*pi/2.0*x_ref[1]);
163  x_phys[1] = x_ref[1] + beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
164  }
165  else{
166  x_phys[0] =x_ref[0] + alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
167  x_phys[1] =x_ref[1] + alpha*exp(1.0-x_ref[1])*(std::sin(pi * x_ref[0]) + std::sin(pi* x_ref[2]));
168  x_phys[2] =x_ref[2] + 1.0/20.0*( std::sin(2.0 * pi * x_ref[0]) + std::sin(2.0 * pi * x_ref[1]));
169  }
170  return dealii::Point<dim> (x_phys); // Trigonometric
171 }
172 
173 template<int dim>
174 dealii::DerivativeForm<1,dim,dim> CurvManifold<dim>::push_forward_gradient(const dealii::Point<dim> &chart_point) const
175 {
176  const double pi = atan(1)*4.0;
177  dealii::DerivativeForm<1, dim, dim> dphys_dref;
178  double beta = 1.0/10.0;
179  double alpha = 1.0/10.0;
180  dealii::Point<dim> x_ref;
181  for(int idim=0; idim<dim; idim++){
182  x_ref[idim] = chart_point[idim];
183  }
184 
185  if(dim==2){
186  dphys_dref[0][0] = 1.0 - beta*pi/2.0 * std::sin(pi/2.0*x_ref[0])*std::cos(3.0*pi/2.0*x_ref[1]);
187  dphys_dref[0][1] = - beta*3.0*pi/2.0 * std::cos(pi/2.0*x_ref[0])*std::sin(3.0*pi/2.0*x_ref[1]);
188 
189  dphys_dref[1][0] = beta*2.0*pi*std::cos(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
190  dphys_dref[1][1] = 1.0 -beta*pi/2.0*std::sin(2.0*pi*(x_ref[0]))*std::sin(pi/2.0*x_ref[1]);
191  }
192  else{
193  dphys_dref[0][0] = 1.0;
194  dphys_dref[0][1] = - alpha*pi*std::sin(pi*x_ref[1]);
195  dphys_dref[0][2] = - alpha*pi*std::sin(pi*x_ref[2]);
196 
197  dphys_dref[1][0] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[0]);
198  dphys_dref[1][1] = 1.0 -alpha*exp(1.0-x_ref[1])*(std::sin(pi*x_ref[0])+std::sin(pi*x_ref[2]));
199  dphys_dref[1][2] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[2]);
200  dphys_dref[2][0] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[0]);
201  dphys_dref[2][1] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[1]);
202  dphys_dref[2][2] = 1.0;
203  }
204 
205  return dphys_dref;
206 }
207 
208 template<int dim>
209 std::unique_ptr<dealii::Manifold<dim,dim> > CurvManifold<dim>::clone() const
210 {
211  return std::make_unique<CurvManifold<dim>>();
212 }
213 
214 template <int dim>
215 static dealii::Point<dim> warp (const dealii::Point<dim> &p)
216 {
217  const double pi = atan(1)*4.0;
218  dealii::Point<dim> q = p;
219 
220  double beta =1.0/10.0;
221  double alpha =1.0/10.0;
222  if (dim == 2){
223  q[dim-2] =p[dim-2] + beta*std::cos(pi/2.0 * p[dim-2]) * std::cos(3.0 * pi/2.0 * p[dim-1]);
224  q[dim-1] =p[dim-1] + beta*std::sin(2.0 * pi * (p[dim-2])) * std::cos(pi /2.0 * p[dim-1]);
225  }
226  if(dim==3){
227  q[0] =p[0] + alpha*(std::cos(pi * p[2]) + std::cos(pi * p[1]));
228  q[1] =p[1] + alpha*exp(1.0-p[1])*(std::sin(pi * p[0]) + std::sin(pi* p[2]));
229  q[2] =p[2] + 1.0/20.0*( std::sin(2.0 * pi * p[0]) + std::sin(2.0 * pi * p[1]));
230  }
231 
232  return q;
233 }
234 
235 /****************************
236  * End of Curvilinear Grid
237  * ***************************/
238 
239 int main (int argc, char * argv[])
240 {
241  dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
242  using real = double;
243  using namespace PHiLiP;
244  std::cout << std::setprecision(std::numeric_limits<long double>::digits10 + 1) << std::scientific;
245  const int dim = PHILIP_DIM;
246  const int nstate = 1;
247  dealii::ParameterHandler parameter_handler;
249  dealii::ConditionalOStream pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)==0);
250 
251  PHiLiP::Parameters::AllParameters all_parameters_new;
252  all_parameters_new.parse_parameters (parameter_handler);
253 
254  all_parameters_new.flux_nodes_type = Parameters::AllParameters::FluxNodes::GLL;
255 
256  double left = 0.0;
257  double right = 1.0;
258  const bool colorize = true;
259  //Generate a standard grid
260 
261  using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
262  std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
263  MPI_COMM_WORLD,
264  typename dealii::Triangulation<dim>::MeshSmoothing(
265  dealii::Triangulation<dim>::smoothing_on_refinement |
266  dealii::Triangulation<dim>::smoothing_on_coarsening));
267  dealii::GridGenerator::hyper_cube (*grid, left, right, colorize);
268  grid->refine_global(0);
269 
270  //Warp the grid
271  //IF WANT NON-WARPED GRID COMMENT UNTIL SAYS "NOT COMMENT"
272  dealii::GridTools::transform (&warp<dim>, *grid);
273 
274  // Assign a manifold to have curved geometry
275  const CurvManifold<dim> curv_manifold;
276  unsigned int manifold_id=0; // top face, see GridGenerator::hyper_rectangle, colorize=true
277  grid->reset_all_manifolds();
278  grid->set_all_manifold_ids(manifold_id);
279  grid->set_manifold ( manifold_id, curv_manifold );
280  //"END COMMENT" TO NOT WARP GRID
281  double max_GCL = 0.0;
282 
283  // Loop through poly degree
284  for(unsigned int poly_degree = 2; poly_degree<5; poly_degree++){
285  unsigned int grid_degree = poly_degree;
286 
287  std::shared_ptr < PHiLiP::DGBase<dim, double> > dg = PHiLiP::DGFactory<dim,double>::create_discontinuous_galerkin(&all_parameters_new, poly_degree, poly_degree, grid_degree, grid);
288  dg->allocate_system ();
289 
290  dealii::QGaussLobatto<1> grid_quad(grid_degree +1);
291  const dealii::FE_DGQ<1> fe_grid(grid_degree);
292  const dealii::FESystem<1,1> fe_sys_grid(fe_grid, nstate);
293  dealii::QGauss<1> flux_quad(poly_degree +1);
294  dealii::QGauss<0> flux_quad_face(poly_degree +1);
295 
296  PHiLiP::OPERATOR::mapping_shape_functions<dim,2*dim,real> mapping_basis(nstate,poly_degree,grid_degree);
297  mapping_basis.build_1D_shape_functions_at_grid_nodes(fe_sys_grid, grid_quad);
298  mapping_basis.build_1D_shape_functions_at_flux_nodes(fe_sys_grid, flux_quad, flux_quad_face);
299 
300  const unsigned int n_quad_pts = pow(poly_degree+1,dim);
301 
302  const dealii::FESystem<dim> &fe_metric = (dg->high_order_grid->fe_system);
303  const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
304  auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
305  for (auto current_cell = dg->dof_handler.begin_active(); current_cell!=dg->dof_handler.end(); ++current_cell, ++metric_cell) {
306  if (!current_cell->is_locally_owned()) continue;
307 
308  pcout<<" degree "<<grid_degree<<" metric dofs "<<n_metric_dofs<<std::endl;
309  std::vector<dealii::types::global_dof_index> current_metric_dofs_indices(n_metric_dofs);
310  metric_cell->get_dof_indices (current_metric_dofs_indices);
311  std::array<std::vector<real>,dim> mapping_support_points;
312  for(int idim=0; idim<dim; idim++){
313  mapping_support_points[idim].resize(n_metric_dofs/dim);
314  }
315  dealii::QGaussLobatto<dim> vol_GLL(grid_degree +1);
316  for (unsigned int igrid_node = 0; igrid_node< n_metric_dofs/dim; ++igrid_node) {
317  for (unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
318  const real val = (dg->high_order_grid->volume_nodes[current_metric_dofs_indices[idof]]);
319  const unsigned int istate = fe_metric.system_to_component_index(idof).first;
320  mapping_support_points[istate][igrid_node] += val * fe_metric.shape_value_component(idof,vol_GLL.point(igrid_node),istate);
321  }
322  }
323 
324  PHiLiP::OPERATOR::metric_operators<real,dim,2*dim> metric_oper(nstate,poly_degree,grid_degree);
325  metric_oper.build_volume_metric_operators(
326  n_quad_pts, n_metric_dofs/dim,
327  mapping_support_points,
328  mapping_basis,
329  false);
330 
331  std::array<std::vector<real>,dim> GCL;
332  for(int idim=0; idim<dim; idim++){
333  GCL[idim].resize(n_quad_pts);
334  }
335 
336  const dealii::FE_DGQArbitraryNodes<1> fe_poly(flux_quad);
337  const dealii::FESystem<1,1> fe_sys_poly(fe_poly, nstate);
339  flux_basis_quad.build_1D_gradient_state_operator(fe_sys_poly, flux_quad);
340  flux_basis_quad.build_1D_volume_state_operator(fe_sys_poly, flux_quad);
341  for(int idim=0; idim<dim; idim++){
342  flux_basis_quad.divergence_matrix_vector_mult(metric_oper.metric_cofactor_vol[idim], GCL[idim],
343  flux_basis_quad.oneD_vol_state_operator[0],
344  flux_basis_quad.oneD_vol_state_operator[0],
345  flux_basis_quad.oneD_vol_state_operator[0],
346  flux_basis_quad.oneD_grad_state_operator[0],
347  flux_basis_quad.oneD_grad_state_operator[0],
348  flux_basis_quad.oneD_grad_state_operator[0]);
349  }
350 
351  for(int idim=0; idim<dim; idim++){
352  // printf("\n GCL for derivative x_%d \n", idim);
353  for(unsigned int idof=0; idof<n_quad_pts; idof++){
354  // printf(" %.16g \n", GCL[idim][idof]);
355  if(std::abs(GCL[idim][idof]) > max_GCL){
356  max_GCL = std::abs(GCL[idim][idof]);
357  }
358  }
359  }
360  } //end cell loop
361  } //end poly degree loop
362 
363  const double max_GCL_mpi= (dealii::Utilities::MPI::max(max_GCL, MPI_COMM_WORLD));
364 
365  if(max_GCL_mpi > 1e-10){
366  pcout<<" Metrics Do NOT Satisfy GCL Condition\n"<<std::endl;
367  return 1;
368  }
369  else{
370  pcout<<" Metrics Satisfy GCL Condition\n"<<std::endl;
371  return 0;
372  }
373 } // end of main
dealii::Tensor< 2, dim, std::vector< real > > metric_cofactor_vol
The volume metric cofactor matrix.
Definition: operators.h:1165
void build_1D_volume_state_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
Definition: operators.cpp:2760
virtual dealii::Point< dim > pull_back(const dealii::Point< dim > &space_point) const override
See dealii::Manifold.
Files for the baseline physics.
Definition: ADTypes.hpp:10
void build_1D_shape_functions_at_grid_nodes(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Constructs the volume operator and gradient operator.
Definition: operators.cpp:2155
Main parameter class that contains the various other sub-parameter classes.
static std::shared_ptr< DGBase< dim, real, MeshType > > create_discontinuous_galerkin(const Parameters::AllParameters *const parameters_input, const unsigned int degree, const unsigned int max_degree_input, const unsigned int grid_degree_input, const std::shared_ptr< Triangulation > triangulation_input)
Creates a derived object DG, but returns it as DGBase.
Definition: dg_factory.cpp:10
Base metric operators class that stores functions used in both the volume and on surface.
Definition: operators.h:1086
void build_1D_gradient_state_operator(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature)
Assembles the one dimensional operator.
Definition: operators.cpp:2784
The mapping shape functions evaluated at the desired nodes (facet set included in volume grid nodes f...
Definition: operators.h:1026
void divergence_matrix_vector_mult(const dealii::Tensor< 1, dim, 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 dealii::FullMatrix< double > &gradient_basis_x, const dealii::FullMatrix< double > &gradient_basis_y, const dealii::FullMatrix< double > &gradient_basis_z)
Computes the divergence using the sum factorization matrix-vector multiplication. ...
Definition: operators.cpp:381
The basis functions separated by nstate with n shape functions.
Definition: operators.h:1316
std::array< dealii::FullMatrix< double >, nstate > oneD_vol_state_operator
Stores the one dimensional volume operator.
Definition: operators.h:1301
virtual std::unique_ptr< dealii::Manifold< dim, dim > > clone() const override
See dealii::Manifold.
virtual dealii::DerivativeForm< 1, dim, dim > push_forward_gradient(const dealii::Point< dim > &chart_point) const override
See dealii::Manifold.
void parse_parameters(dealii::ParameterHandler &prm)
Retrieve parameters from dealii::ParameterHandler.
FluxNodes flux_nodes_type
Store selected FluxNodes from the input file.
void build_volume_metric_operators(const unsigned int n_quad_pts, const unsigned int n_metric_dofs, const std::array< std::vector< real >, dim > &mapping_support_points, mapping_shape_functions< dim, n_faces, real > &mapping_basis, const bool use_invariant_curl_form=false)
Builds the volume metric operators.
Definition: operators.cpp:2294
virtual dealii::Point< dim > push_forward(const dealii::Point< dim > &chart_point) const override
See dealii::Manifold.
std::array< dealii::FullMatrix< double >, nstate > oneD_grad_state_operator
Stores the one dimensional gradient operator.
Definition: operators.h:1307
static void declare_parameters(dealii::ParameterHandler &prm)
Declare parameters that can be set as inputs and set up the default options.
void build_1D_shape_functions_at_flux_nodes(const dealii::FESystem< 1, 1 > &finite_element, const dealii::Quadrature< 1 > &quadrature, const dealii::Quadrature< 0 > &face_quadrature)
Constructs the volume, gradient, surface, and surface gradient operator.
Definition: operators.cpp:2164