[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
surface_GCL_Superparametric_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 template <int dim>
57 class CurvManifold: public dealii::ChartManifold<dim,dim,dim> {
58  virtual dealii::Point<dim> pull_back(const dealii::Point<dim> &space_point) const override;
59  virtual dealii::Point<dim> push_forward(const dealii::Point<dim> &chart_point) const override;
60  virtual dealii::DerivativeForm<1,dim,dim> push_forward_gradient(const dealii::Point<dim> &chart_point) const override;
61 
62  virtual std::unique_ptr<dealii::Manifold<dim,dim> > clone() const override;
63 };
64 
65 template<int dim>
66 dealii::Point<dim> CurvManifold<dim>::pull_back(const dealii::Point<dim> &space_point) const
67 {
68  using namespace PHiLiP;
69  const double pi = atan(1)*4.0;
70  dealii::Point<dim> x_ref;
71  dealii::Point<dim> x_phys;
72  for(int idim=0; idim<dim; idim++){
73  x_ref[idim] = space_point[idim];
74  x_phys[idim] = space_point[idim];
75  }
76  dealii::Vector<double> function(dim);
77  dealii::FullMatrix<double> derivative(dim);
78  double beta =1.0/10.0;
79  double alpha =1.0/10.0;
80  int flag =0;
81  while(flag != dim){
82  if(dim==2){
83  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]);
84  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]);
85  }
86  else{
87  function[0] = x_ref[0] - x_phys[0] +alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
88  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]));
89  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]));
90  }
91  if(dim==2){
92  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]);
93  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]);
94 
95  derivative[1][0] = beta*2.0*pi*std::cos(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
96  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]);
97  }
98  else{
99  derivative[0][0] = 1.0;
100  derivative[0][1] = - alpha*pi*std::sin(pi*x_ref[1]);
101  derivative[0][2] = - alpha*pi*std::sin(pi*x_ref[2]);
102 
103  derivative[1][0] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[0]);
104  derivative[1][1] = 1.0 -alpha*exp(1.0-x_ref[1])*(std::sin(pi*x_ref[0])+std::sin(pi*x_ref[2]));
105  derivative[1][2] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[2]);
106  derivative[2][0] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[0]);
107  derivative[2][1] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[1]);
108  derivative[2][2] = 1.0;
109  }
110 
111  dealii::FullMatrix<double> Jacobian_inv(dim);
112  Jacobian_inv.invert(derivative);
113  dealii::Vector<double> Newton_Step(dim);
114  Jacobian_inv.vmult(Newton_Step, function);
115  for(int idim=0; idim<dim; idim++){
116  x_ref[idim] -= Newton_Step[idim];
117  }
118  flag=0;
119  for(int idim=0; idim<dim; idim++){
120  if(std::abs(function[idim]) < 1e-15)
121  flag++;
122  }
123  if(flag == dim)
124  break;
125  }
126  std::vector<double> function_check(dim);
127  if(dim==2){
128  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]);
129  function_check[1] = x_ref[1] + beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
130  }
131  else{
132  function_check[0] = x_ref[0] +alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
133  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]));
134  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]));
135  }
136  std::vector<double> error(dim);
137  for(int idim=0; idim<dim; idim++)
138  error[idim] = std::abs(function_check[idim] - x_phys[idim]);
139  if (error[0] > 1e-13) {
140  std::cout << "Large error " << error[0] << std::endl;
141  for(int idim=0;idim<dim; idim++)
142  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;
143  }
144 
145  return x_ref;
146 }
147 
148 template<int dim>
149 dealii::Point<dim> CurvManifold<dim>::push_forward(const dealii::Point<dim> &chart_point) const
150 {
151  const double pi = atan(1)*4.0;
152 
153  dealii::Point<dim> x_ref;
154  dealii::Point<dim> x_phys;
155  for(int idim=0; idim<dim; idim++)
156  x_ref[idim] = chart_point[idim];
157  double beta = 1.0/10.0;
158  double alpha = 1.0/10.0;
159  if(dim==2){
160  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]);
161  x_phys[1] = x_ref[1] + beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
162  }
163  else{
164  x_phys[0] =x_ref[0] + alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
165  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]));
166  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]));
167  }
168  return dealii::Point<dim> (x_phys); // Trigonometric
169 }
170 
171 template<int dim>
172 dealii::DerivativeForm<1,dim,dim> CurvManifold<dim>::push_forward_gradient(const dealii::Point<dim> &chart_point) const
173 {
174  const double pi = atan(1)*4.0;
175  dealii::DerivativeForm<1, dim, dim> dphys_dref;
176  double beta = 1.0/10.0;
177  double alpha = 1.0/10.0;
178  dealii::Point<dim> x_ref;
179  for(int idim=0; idim<dim; idim++){
180  x_ref[idim] = chart_point[idim];
181  }
182 
183  if(dim==2){
184  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]);
185  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]);
186 
187  dphys_dref[1][0] = beta*2.0*pi*std::cos(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
188  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]);
189  }
190  else{
191  dphys_dref[0][0] = 1.0;
192  dphys_dref[0][1] = - alpha*pi*std::sin(pi*x_ref[1]);
193  dphys_dref[0][2] = - alpha*pi*std::sin(pi*x_ref[2]);
194 
195  dphys_dref[1][0] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[0]);
196  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]));
197  dphys_dref[1][2] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[2]);
198  dphys_dref[2][0] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[0]);
199  dphys_dref[2][1] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[1]);
200  dphys_dref[2][2] = 1.0;
201  }
202 
203  return dphys_dref;
204 }
205 
206 template<int dim>
207 std::unique_ptr<dealii::Manifold<dim,dim> > CurvManifold<dim>::clone() const
208 {
209  return std::make_unique<CurvManifold<dim>>();
210 }
211 
212 template <int dim>
213 static dealii::Point<dim> warp (const dealii::Point<dim> &p)
214 {
215  const double pi = atan(1)*4.0;
216  dealii::Point<dim> q = p;
217 
218  double beta =1.0/10.0;
219  double alpha =1.0/10.0;
220  if (dim == 2){
221  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]);
222  q[dim-1] =p[dim-1] + beta*std::sin(2.0 * pi * (p[dim-2])) * std::cos(pi /2.0 * p[dim-1]);
223  }
224  if(dim==3){
225  q[0] =p[0] + alpha*(std::cos(pi * p[2]) + std::cos(pi * p[1]));
226  q[1] =p[1] + alpha*exp(1.0-p[1])*(std::sin(pi * p[0]) + std::sin(pi* p[2]));
227  q[2] =p[2] + 1.0/20.0*( std::sin(2.0 * pi * p[0]) + std::sin(2.0 * pi * p[1]));
228  }
229 
230  return q;
231 }
232 
233 /****************************
234  * End of Curvilinear Grid
235  * ***************************/
236 
237 int main (int argc, char * argv[])
238 {
239  dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
240  using real = double;
241  using namespace PHiLiP;
242  std::cout << std::setprecision(std::numeric_limits<long double>::digits10 + 1) << std::scientific;
243  const int dim = PHILIP_DIM;
244  const int nstate = 1;
245  dealii::ParameterHandler parameter_handler;
247  dealii::ConditionalOStream pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)==0);
248 
249  PHiLiP::Parameters::AllParameters all_parameters_new;
250  all_parameters_new.parse_parameters (parameter_handler);
251 
252  // all_parameters_new.flux_nodes_type = Parameters::AllParameters::FluxNodes::GLL;
253 
254  //unsigned int poly_degree = 3;
255  double left = 0.0;
256  double right = 1.0;
257  const bool colorize = true;
258  //Generate a standard grid
259  using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
260  std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
261  MPI_COMM_WORLD,
262  typename dealii::Triangulation<dim>::MeshSmoothing(
263  dealii::Triangulation<dim>::smoothing_on_refinement |
264  dealii::Triangulation<dim>::smoothing_on_coarsening));
265  dealii::GridGenerator::hyper_cube (*grid, left, right, colorize);
266  grid->refine_global(0);
267 
268  //Warp the grid
269  //IF WANT NON-WARPED GRID COMMENT UNTIL SAYS "NOT COMMENT"
270  dealii::GridTools::transform (&warp<dim>, *grid);
271 
272  // Assign a manifold to have curved geometry
273  const CurvManifold<dim> curv_manifold;
274  unsigned int manifold_id=0; // top face, see GridGenerator::hyper_rectangle, colorize=true
275  grid->reset_all_manifolds();
276  grid->set_all_manifold_ids(manifold_id);
277  grid->set_manifold ( manifold_id, curv_manifold );
278  //"END COMMENT" TO NOT WARP GRID
279 
280  double surf_int = 0.0;
281  for(unsigned int poly_degree = 2; poly_degree<6; poly_degree++){
282  unsigned int grid_degree = poly_degree + 1;
283  //setup DG
284  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);
285  dg->allocate_system ();
286 
287  const dealii::FESystem<dim> &fe_metric = (dg->high_order_grid->fe_system);
288  const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
289  auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
290 
291  dealii::QGaussLobatto<1> grid_quad(grid_degree +1);
292  const dealii::FE_DGQ<1> fe_grid(grid_degree);
293  const dealii::FESystem<1,1> fe_sys_grid(fe_grid, nstate);
294  dealii::QGauss<1> flux_quad(poly_degree +1);
295  dealii::QGauss<0> flux_quad_face(poly_degree +1);
296 
297  PHiLiP::OPERATOR::mapping_shape_functions<dim,2*dim,real> mapping_basis(nstate,poly_degree,grid_degree);
298  mapping_basis.build_1D_shape_functions_at_grid_nodes(fe_sys_grid, grid_quad);
299  mapping_basis.build_1D_shape_functions_at_flux_nodes(fe_sys_grid, flux_quad, flux_quad_face);
300 
301  PHiLiP::OPERATOR::metric_operators<real,dim,2*dim> metric_oper(nstate, poly_degree, grid_degree);
302 
303  for (auto current_cell = dg->dof_handler.begin_active(); current_cell!=dg->dof_handler.end(); ++current_cell, ++metric_cell) {
304  if (!current_cell->is_locally_owned()) continue;
305 
306  std::vector<dealii::types::global_dof_index> current_metric_dofs_indices(n_metric_dofs);
307  metric_cell->get_dof_indices (current_metric_dofs_indices);
308  std::array<std::vector<real>,dim> mapping_support_points;
309  for(int idim=0; idim<dim; idim++){
310  mapping_support_points[idim].resize(n_metric_dofs/dim);
311  }
312  dealii::QGaussLobatto<dim> vol_GLL(grid_degree +1);
313  for (unsigned int igrid_node = 0; igrid_node< n_metric_dofs/dim; ++igrid_node) {
314  for (unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
315  const real val = (dg->high_order_grid->volume_nodes[current_metric_dofs_indices[idof]]);
316  const unsigned int istate = fe_metric.system_to_component_index(idof).first;
317  mapping_support_points[istate][igrid_node] += val * fe_metric.shape_value_component(idof,vol_GLL.point(igrid_node),istate);
318  }
319  }
320 
321  const unsigned int n_quad_face_pts = dg->face_quadrature_collection[poly_degree].size();
322  const std::vector<real> &quad_weights = dg->face_quadrature_collection[poly_degree].get_weights ();
323  for (unsigned int iface=0; iface < dealii::GeometryInfo<dim>::faces_per_cell; ++iface) {
324 
325  metric_oper.build_facet_metric_operators(
326  iface,
327  n_quad_face_pts, n_metric_dofs/dim,
328  mapping_support_points,
329  mapping_basis,
330  false);
331 
332  const dealii::Tensor<1,dim,real> unit_normal_int = dealii::GeometryInfo<dim>::unit_normal_vector[iface];
333  std::vector<dealii::Tensor<1,dim,real>> normals_int(n_quad_face_pts);
334  for(unsigned int iquad=0; iquad<n_quad_face_pts; iquad++){
335  for(unsigned int idim=0; idim<dim; idim++){
336  normals_int[iquad][idim] = 0.0;
337  for(int idim2=0; idim2<dim; idim2++){
338  normals_int[iquad][idim] += unit_normal_int[idim2] * metric_oper.metric_cofactor_surf[idim][idim2][iquad];//\hat{n}^r * C_m^T
339  }
340  }
341  }
342  for(unsigned int iquad=0; iquad<n_quad_face_pts; iquad++){
343  for(int idim=0; idim<dim; idim++){
344  surf_int += 1.0 * quad_weights[iquad] * normals_int[iquad][idim] * 1.0;
345  }
346  }
347  }//end of face loop
348  }//end of cell loop
349  }//end poly degree loop
350  const double surf_int_mpi= (dealii::Utilities::MPI::max(surf_int, MPI_COMM_WORLD));
351 
352  if(std::abs(surf_int_mpi) > 1e-13){
353  pcout<<" Metrics Do NOT Satisfy GCL Condition\n"<<std::endl;
354  return 1;
355  }
356  else{
357  pcout<<" Metrics Satisfy GCL Condition\n"<<std::endl;
358  return 0;
359  }
360 }//end of main
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
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
The mapping shape functions evaluated at the desired nodes (facet set included in volume grid nodes f...
Definition: operators.h:1026
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.
virtual dealii::Point< dim > push_forward(const dealii::Point< dim > &chart_point) const override
See dealii::Manifold.
static void declare_parameters(dealii::ParameterHandler &prm)
Declare parameters that can be set as inputs and set up the default options.