9 #include <deal.II/distributed/solution_transfer.h> 11 #include "testing/tests.h" 14 #include <deal.II/base/parameter_handler.h> 15 #include <deal.II/base/tensor.h> 16 #include <deal.II/numerics/vector_tools.h> 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> 24 #include <deal.II/dofs/dof_handler.h> 25 #include <deal.II/dofs/dof_tools.h> 26 #include <deal.II/dofs/dof_renumbering.h> 28 #include <deal.II/dofs/dof_accessor.h> 30 #include <deal.II/lac/vector.h> 31 #include <deal.II/lac/dynamic_sparsity_pattern.h> 32 #include <deal.II/lac/sparse_matrix.h> 34 #include <deal.II/meshworker/dof_info.h> 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> 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" 52 const double TOLERANCE = 1E-6;
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;
63 virtual std::unique_ptr<dealii::Manifold<dim,dim> > clone()
const override;
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];
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;
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]);
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]));
95 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]);
96 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]);
98 derivative[1][0] = beta*2.0*pi*std::cos(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
99 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]);
102 derivative[0][0] = 1.0;
103 derivative[0][1] = - alpha*pi*std::sin(pi*x_ref[1]);
104 derivative[0][2] = - alpha*pi*std::sin(pi*x_ref[2]);
106 derivative[1][0] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[0]);
107 derivative[1][1] = 1.0 -alpha*exp(1.0-x_ref[1])*(std::sin(pi*x_ref[0])+std::sin(pi*x_ref[2]));
108 derivative[1][2] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[2]);
109 derivative[2][0] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[0]);
110 derivative[2][1] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[1]);
111 derivative[2][2] = 1.0;
114 dealii::FullMatrix<double> Jacobian_inv(dim);
115 Jacobian_inv.invert(derivative);
116 dealii::Vector<double> Newton_Step(dim);
117 Jacobian_inv.vmult(Newton_Step,
function);
118 for(
int idim=0; idim<dim; idim++){
119 x_ref[idim] -= Newton_Step[idim];
122 for(
int idim=0; idim<dim; idim++){
123 if(std::abs(
function[idim]) < 1e-15)
129 std::vector<double> function_check(dim);
131 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]);
132 function_check[1] = x_ref[1] + beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
135 function_check[0] = x_ref[0] +alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
136 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]));
137 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]));
139 std::vector<double> error(dim);
140 for(
int idim=0; idim<dim; idim++)
141 error[idim] = std::abs(function_check[idim] - x_phys[idim]);
142 if (error[0] > 1e-13) {
143 std::cout <<
"Large error " << error[0] << std::endl;
144 for(
int idim=0;idim<dim; idim++)
145 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;
154 const double pi = atan(1)*4.0;
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;
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]);
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]));
171 return dealii::Point<dim> (x_phys);
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];
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]);
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]);
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]);
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;
212 return std::make_unique<CurvManifold<dim>>();
216 static dealii::Point<dim> warp (
const dealii::Point<dim> &p)
218 const double pi = atan(1)*4.0;
219 dealii::Point<dim> q = p;
221 double beta =1.0/10.0;
222 double alpha =1.0/10.0;
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]);
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]));
240 int main (
int argc,
char * argv[])
242 dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
245 std::cout << std::setprecision(std::numeric_limits<long double>::digits10 + 1) << std::scientific;
246 const int dim = PHILIP_DIM;
247 const int nstate = 1;
248 dealii::ParameterHandler parameter_handler;
250 dealii::ConditionalOStream pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD)==0);
257 const bool colorize =
true;
259 using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
260 std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
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(3);
270 dealii::GridTools::transform (&warp<dim>, *grid);
273 unsigned int manifold_id=0;
274 grid->reset_all_manifolds();
275 grid->set_all_manifold_ids(manifold_id);
276 grid->set_manifold ( manifold_id, curv_manifold );
278 double max_GCL = 0.0;
279 for(
unsigned int poly_degree = 2; poly_degree<5; poly_degree++){
280 unsigned int grid_degree = poly_degree;
283 dg->allocate_system ();
285 dealii::QGaussLobatto<1> grid_quad(grid_degree +1);
286 const dealii::FE_DGQ<1> fe_grid(grid_degree);
287 const dealii::FESystem<1,1> fe_sys_grid(fe_grid, nstate);
288 dealii::QGauss<1> flux_quad(poly_degree +1);
289 dealii::QGauss<0> flux_quad_face(poly_degree +1);
292 mapping_basis.build_1D_shape_functions_at_grid_nodes(fe_sys_grid, grid_quad);
293 mapping_basis.build_1D_shape_functions_at_flux_nodes(fe_sys_grid, flux_quad, flux_quad_face);
295 const unsigned int n_quad_pts = pow(poly_degree+1,dim);
297 const dealii::FESystem<dim> &fe_metric = (dg->high_order_grid->fe_system);
298 const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
299 auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
300 for (
auto current_cell = dg->dof_handler.begin_active(); current_cell!=dg->dof_handler.end(); ++current_cell, ++metric_cell) {
301 if (!current_cell->is_locally_owned())
continue;
303 std::vector<dealii::types::global_dof_index> current_metric_dofs_indices(n_metric_dofs);
304 metric_cell->get_dof_indices (current_metric_dofs_indices);
305 std::array<std::vector<real>,dim> mapping_support_points;
306 for(
int idim=0; idim<dim; idim++){
307 mapping_support_points[idim].resize(n_metric_dofs/dim);
309 dealii::QGaussLobatto<dim> vol_GLL(grid_degree +1);
310 for (
unsigned int igrid_node = 0; igrid_node< n_metric_dofs/dim; ++igrid_node) {
311 for (
unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
312 const real val = (dg->high_order_grid->volume_nodes[current_metric_dofs_indices[idof]]);
313 const unsigned int istate = fe_metric.system_to_component_index(idof).first;
314 mapping_support_points[istate][igrid_node] += val * fe_metric.shape_value_component(idof,vol_GLL.point(igrid_node),istate);
319 metric_oper.build_volume_metric_operators(
320 n_quad_pts, n_metric_dofs/dim,
321 mapping_support_points,
325 std::array<std::vector<real>,dim> GCL;
326 for(
int idim=0; idim<dim; idim++){
327 GCL[idim].resize(n_quad_pts);
330 const dealii::FE_DGQArbitraryNodes<1> fe_poly(flux_quad);
331 const dealii::FESystem<1,1> fe_sys_poly(fe_poly, nstate);
333 flux_basis_quad.build_1D_gradient_state_operator(fe_sys_poly, flux_quad);
334 flux_basis_quad.build_1D_volume_state_operator(fe_sys_poly, flux_quad);
335 for(
int idim=0; idim<dim; idim++){
336 flux_basis_quad.divergence_matrix_vector_mult(metric_oper.metric_cofactor_vol[idim], GCL[idim],
337 flux_basis_quad.oneD_vol_state_operator[0],
338 flux_basis_quad.oneD_vol_state_operator[0],
339 flux_basis_quad.oneD_vol_state_operator[0],
340 flux_basis_quad.oneD_grad_state_operator[0],
341 flux_basis_quad.oneD_grad_state_operator[0],
342 flux_basis_quad.oneD_grad_state_operator[0]);
345 for(
int idim=0; idim<dim; idim++){
347 for(
unsigned int idof=0; idof<n_quad_pts; idof++){
349 if( std::abs(GCL[idim][idof]) > max_GCL){
350 max_GCL = std::abs(GCL[idim][idof]);
357 const double max_GCL_mpi= (dealii::Utilities::MPI::max(max_GCL, MPI_COMM_WORLD));
358 if(max_GCL_mpi > 1e-10){
359 pcout<<
" Metrics Do NOT Satisfy GCL Condition\n"<<std::endl;
363 pcout<<
" Metrics Satisfy GCL Condition\n"<<std::endl;
virtual dealii::Point< dim > pull_back(const dealii::Point< dim > &space_point) const override
See dealii::Manifold.
Files for the baseline physics.
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.
Base metric operators class that stores functions used in both the volume and on surface.
The mapping shape functions evaluated at the desired nodes (facet set included in volume grid nodes f...
The basis functions separated by nstate with n shape functions.
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.