[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
GCL_test_invariant_curl.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 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 
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 
151 template<int dim>
152 dealii::Point<dim> CurvManifold<dim>::push_forward(const dealii::Point<dim> &chart_point) const
153 {
154  const double pi = atan(1)*4.0;
155 
156  dealii::Point<dim> x_ref;
157  dealii::Point<dim> x_phys;
158  for(int idim=0; idim<dim; idim++)
159  x_ref[idim] = chart_point[idim];
160  double beta = 1.0/10.0;
161  double alpha = 1.0/10.0;
162  if(dim==2){
163  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]);
164  x_phys[1] = x_ref[1] + beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
165  }
166  else{
167  x_phys[0] =x_ref[0] + alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
168  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]));
169  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]));
170  }
171  return dealii::Point<dim> (x_phys); // Trigonometric
172 }
173 
174 template<int dim>
175 dealii::DerivativeForm<1,dim,dim> CurvManifold<dim>::push_forward_gradient(const dealii::Point<dim> &chart_point) const
176 {
177  const double pi = atan(1)*4.0;
178  dealii::DerivativeForm<1, dim, dim> dphys_dref;
179  double beta = 1.0/10.0;
180  double alpha = 1.0/10.0;
181  dealii::Point<dim> x_ref;
182  for(int idim=0; idim<dim; idim++){
183  x_ref[idim] = chart_point[idim];
184  }
185 
186  if(dim==2){
187  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]);
188  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]);
189 
190  dphys_dref[1][0] = beta*2.0*pi*std::cos(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
191  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]);
192  }
193  else{
194  dphys_dref[0][0] = 1.0;
195  dphys_dref[0][1] = - alpha*pi*std::sin(pi*x_ref[1]);
196  dphys_dref[0][2] = - alpha*pi*std::sin(pi*x_ref[2]);
197 
198  dphys_dref[1][0] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[0]);
199  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]));
200  dphys_dref[1][2] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[2]);
201  dphys_dref[2][0] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[0]);
202  dphys_dref[2][1] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[1]);
203  dphys_dref[2][2] = 1.0;
204  }
205 
206  return dphys_dref;
207 }
208 
209 template<int dim>
210 std::unique_ptr<dealii::Manifold<dim,dim> > CurvManifold<dim>::clone() const
211 {
212  return std::make_unique<CurvManifold<dim>>();
213 }
214 
215 template <int dim>
216 static dealii::Point<dim> warp (const dealii::Point<dim> &p)
217 {
218  const double pi = atan(1)*4.0;
219  dealii::Point<dim> q = p;
220 
221  double beta =1.0/10.0;
222  double alpha =1.0/10.0;
223  if (dim == 2){
224  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]);
225  q[dim-1] =p[dim-1] + beta*std::sin(2.0 * pi * (p[dim-2])) * std::cos(pi /2.0 * p[dim-1]);
226  }
227  if(dim==3){
228  q[0] =p[0] + alpha*(std::cos(pi * p[2]) + std::cos(pi * p[1]));
229  q[1] =p[1] + alpha*exp(1.0-p[1])*(std::sin(pi * p[0]) + std::sin(pi* p[2]));
230  q[2] =p[2] + 1.0/20.0*( std::sin(2.0 * pi * p[0]) + std::sin(2.0 * pi * p[1]));
231  }
232 
233  return q;
234 }
235 
236 /****************************
237  * End of Curvilinear Grid
238  * ***************************/
239 
240 int main (int argc, char * argv[])
241 {
242 
243  dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
244  using real = double;
245  using namespace PHiLiP;
246  std::cout << std::setprecision(std::numeric_limits<long double>::digits10 + 1) << std::scientific;
247  const int dim = PHILIP_DIM;
248  const int nstate = 1;
249  dealii::ParameterHandler parameter_handler;
251  dealii::ConditionalOStream pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)==0);
252 
253  PHiLiP::Parameters::AllParameters all_parameters_new;
254  all_parameters_new.parse_parameters (parameter_handler);
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  for(unsigned int poly_degree = 2; poly_degree<5; poly_degree++){
283  unsigned int grid_degree = poly_degree;
284 
285  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);
286  dg->allocate_system ();
287 
288  dealii::QGaussLobatto<1> grid_quad(grid_degree +1);
289  const dealii::FE_DGQ<1> fe_grid(grid_degree);
290  const dealii::FESystem<1,1> fe_sys_grid(fe_grid, nstate);
291  dealii::QGauss<1> flux_quad(poly_degree +1);
292  dealii::QGauss<0> flux_quad_face(poly_degree +1);
293 
294  PHiLiP::OPERATOR::mapping_shape_functions<dim,2*dim,real> mapping_basis(nstate,poly_degree,grid_degree);
295  mapping_basis.build_1D_shape_functions_at_grid_nodes(fe_sys_grid, grid_quad);
296  mapping_basis.build_1D_shape_functions_at_flux_nodes(fe_sys_grid, flux_quad, flux_quad_face);
297 
298  const unsigned int n_quad_pts = pow(poly_degree+1,dim);
299 
300  const dealii::FESystem<dim> &fe_metric = (dg->high_order_grid->fe_system);
301  const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
302  auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
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  pcout<<" degree "<<grid_degree<<" metric dofs "<<n_metric_dofs<<std::endl;
307  std::vector<dealii::types::global_dof_index> current_metric_dofs_indices(n_metric_dofs);
308  metric_cell->get_dof_indices (current_metric_dofs_indices);
309  std::array<std::vector<real>,dim> mapping_support_points;
310  for(int idim=0; idim<dim; idim++){
311  mapping_support_points[idim].resize(n_metric_dofs/dim);
312  }
313  dealii::QGaussLobatto<dim> vol_GLL(grid_degree +1);
314  for (unsigned int igrid_node = 0; igrid_node< n_metric_dofs/dim; ++igrid_node) {
315  for (unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
316  const real val = (dg->high_order_grid->volume_nodes[current_metric_dofs_indices[idof]]);
317  const unsigned int istate = fe_metric.system_to_component_index(idof).first;
318  mapping_support_points[istate][igrid_node] += val * fe_metric.shape_value_component(idof,vol_GLL.point(igrid_node),istate);
319  }
320  }
321 
322  PHiLiP::OPERATOR::metric_operators<real,dim,2*dim> metric_oper(nstate,poly_degree,grid_degree);
323  metric_oper.build_volume_metric_operators(
324  n_quad_pts, n_metric_dofs/dim,
325  mapping_support_points,
326  mapping_basis,
327  true);
328 
329  std::array<std::vector<real>,dim> GCL;
330  for(int idim=0; idim<dim; idim++){
331  GCL[idim].resize(n_quad_pts);
332  }
333 
334  const dealii::FE_DGQArbitraryNodes<1> fe_poly(flux_quad);
335  const dealii::FESystem<1,1> fe_sys_poly(fe_poly, nstate);
337  flux_basis_quad.build_1D_gradient_state_operator(fe_sys_poly, flux_quad);
338  flux_basis_quad.build_1D_volume_state_operator(fe_sys_poly, flux_quad);
339  for(int idim=0; idim<dim; idim++){
340  flux_basis_quad.divergence_matrix_vector_mult(metric_oper.metric_cofactor_vol[idim], GCL[idim],
341  flux_basis_quad.oneD_vol_state_operator[0],
342  flux_basis_quad.oneD_vol_state_operator[0],
343  flux_basis_quad.oneD_vol_state_operator[0],
344  flux_basis_quad.oneD_grad_state_operator[0],
345  flux_basis_quad.oneD_grad_state_operator[0],
346  flux_basis_quad.oneD_grad_state_operator[0]);
347  }
348 
349  for(int idim=0; idim<dim; idim++){
350  // printf("\n GCL for derivative x_%d \n", idim);
351  for(unsigned int idof=0; idof<n_quad_pts; idof++){
352  // printf(" %.16g \n", GCL[idim][idof]);
353  if( std::abs(GCL[idim][idof]) > max_GCL){
354  max_GCL = std::abs(GCL[idim][idof]);
355  }
356  }
357  }
358 
359  }
360 
361  }//end poly degree loop
362  const double max_GCL_mpi= (dealii::Utilities::MPI::max(max_GCL, MPI_COMM_WORLD));
363 
364  if( max_GCL_mpi > 1e-10){
365  pcout<<" Metrics Do NOT Satisfy GCL Condition\n"<<std::endl;
366  return 1;
367  }
368  else{
369  pcout<<" Metrics Satisfy GCL Condition\n"<<std::endl;
370  return 0;
371  }
372 }
373 
374 //}//end PHiLiP namespace
375 
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
The basis functions separated by nstate with n shape functions.
Definition: operators.h:1316
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.