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> 36 #include <deal.II/base/convergence_table.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 "parameters/all_parameters.h" 46 #include "parameters/parameters.h" 47 #include "dg/dg_base.hpp" 48 #include <deal.II/grid/manifold_lib.h> 49 #include <deal.II/fe/mapping_q.h> 50 #include "dg/dg_factory.hpp" 51 #include "operators/operators.h" 54 const double TOLERANCE = 1E-6;
59 class CurvManifold:
public dealii::ChartManifold<dim,dim,dim> {
60 virtual dealii::Point<dim> pull_back(
const dealii::Point<dim> &space_point)
const override;
61 virtual dealii::Point<dim> push_forward(
const dealii::Point<dim> &chart_point)
const override;
62 virtual dealii::DerivativeForm<1,dim,dim> push_forward_gradient(
const dealii::Point<dim> &chart_point)
const override;
64 virtual std::unique_ptr<dealii::Manifold<dim,dim> > clone()
const override;
71 const double pi = atan(1)*4.0;
72 dealii::Point<dim> x_ref;
73 dealii::Point<dim> x_phys;
74 for(
int idim=0; idim<dim; idim++){
75 x_ref[idim] = space_point[idim];
76 x_phys[idim] = space_point[idim];
78 dealii::Vector<double>
function(dim);
79 dealii::FullMatrix<double> derivative(dim);
80 double beta =1.0/10.0;
81 double alpha =1.0/10.0;
85 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]);
86 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]);
89 function[0] = x_ref[0] - x_phys[0] +alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
90 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]));
91 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]));
96 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]);
97 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]);
99 derivative[1][0] = beta*2.0*pi*std::cos(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
100 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]);
103 derivative[0][0] = 1.0;
104 derivative[0][1] = - alpha*pi*std::sin(pi*x_ref[1]);
105 derivative[0][2] = - alpha*pi*std::sin(pi*x_ref[2]);
107 derivative[1][0] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[0]);
108 derivative[1][1] = 1.0 -alpha*exp(1.0-x_ref[1])*(std::sin(pi*x_ref[0])+std::sin(pi*x_ref[2]));
109 derivative[1][2] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[2]);
110 derivative[2][0] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[0]);
111 derivative[2][1] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[1]);
112 derivative[2][2] = 1.0;
115 dealii::FullMatrix<double> Jacobian_inv(dim);
116 Jacobian_inv.invert(derivative);
117 dealii::Vector<double> Newton_Step(dim);
118 Jacobian_inv.vmult(Newton_Step,
function);
119 for(
int idim=0; idim<dim; idim++){
120 x_ref[idim] -= Newton_Step[idim];
123 for(
int idim=0; idim<dim; idim++){
124 if(std::abs(
function[idim]) < 1e-15)
130 std::vector<double> function_check(dim);
132 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]);
133 function_check[1] = x_ref[1] + beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
136 function_check[0] = x_ref[0] +alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
137 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]));
138 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]));
140 std::vector<double> error(dim);
141 for(
int idim=0; idim<dim; idim++)
142 error[idim] = std::abs(function_check[idim] - x_phys[idim]);
143 if (error[0] > 1e-13) {
144 std::cout <<
"Large error " << error[0] << std::endl;
145 for(
int idim=0;idim<dim; idim++)
146 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;
155 const double pi = atan(1)*4.0;
157 dealii::Point<dim> x_ref;
158 dealii::Point<dim> x_phys;
159 for(
int idim=0; idim<dim; idim++)
160 x_ref[idim] = chart_point[idim];
161 double beta = 1.0/10.0;
162 double alpha = 1.0/10.0;
164 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]);
165 x_phys[1] = x_ref[1] + beta*std::sin(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
168 x_phys[0] =x_ref[0] + alpha*(std::cos(pi * x_ref[2]) + std::cos(pi * x_ref[1]));
169 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]));
170 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]));
172 return dealii::Point<dim> (x_phys);
178 const double pi = atan(1)*4.0;
179 dealii::DerivativeForm<1, dim, dim> dphys_dref;
180 double beta = 1.0/10.0;
181 double alpha = 1.0/10.0;
182 dealii::Point<dim> x_ref;
183 for(
int idim=0; idim<dim; idim++){
184 x_ref[idim] = chart_point[idim];
188 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]);
189 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]);
191 dphys_dref[1][0] = beta*2.0*pi*std::cos(2.0*pi*(x_ref[0]))*std::cos(pi/2.0*x_ref[1]);
192 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]);
195 dphys_dref[0][0] = 1.0;
196 dphys_dref[0][1] = - alpha*pi*std::sin(pi*x_ref[1]);
197 dphys_dref[0][2] = - alpha*pi*std::sin(pi*x_ref[2]);
199 dphys_dref[1][0] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[0]);
200 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]));
201 dphys_dref[1][2] = alpha*pi*exp(1.0-x_ref[1])*std::cos(pi*x_ref[2]);
202 dphys_dref[2][0] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[0]);
203 dphys_dref[2][1] = 1.0/10.0*pi*std::cos(2.0*pi*x_ref[1]);
204 dphys_dref[2][2] = 1.0;
213 return std::make_unique<CurvManifold<dim>>();
217 static dealii::Point<dim> warp (
const dealii::Point<dim> &p)
219 const double pi = atan(1)*4.0;
220 dealii::Point<dim> q = p;
222 double beta =1.0/10.0;
223 double alpha =1.0/10.0;
225 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]);
226 q[dim-1] =p[dim-1] + beta*std::sin(2.0 * pi * (p[dim-2])) * std::cos(pi /2.0 * p[dim-1]);
229 q[0] =p[0] + alpha*(std::cos(pi * p[2]) + std::cos(pi * p[1]));
230 q[1] =p[1] + alpha*exp(1.0-p[1])*(std::sin(pi * p[0]) + std::sin(pi* p[2]));
231 q[2] =p[2] + 1.0/20.0*( std::sin(2.0 * pi * p[0]) + std::sin(2.0 * pi * p[1]));
241 int main (
int argc,
char * argv[])
243 dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
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);
257 const int dim_check = 1;
261 const bool colorize =
true;
262 dealii::ConvergenceTable convergence_table;
263 const unsigned int igrid_start = 2;
264 const unsigned int n_grids = 5;
265 std::array<double,n_grids> grid_size;
266 std::array<double,n_grids> soln_error;
267 std::array<double,n_grids> soln_error_inf;
269 for(
unsigned int poly_degree = 2; poly_degree<5; poly_degree++){
270 unsigned int grid_degree = poly_degree;
271 for(
unsigned int igrid=igrid_start; igrid<n_grids; ++igrid){
272 pcout<<
" Grid Index"<<igrid<<std::endl;
274 using Triangulation = dealii::parallel::distributed::Triangulation<dim>;
275 std::shared_ptr<Triangulation> grid = std::make_shared<Triangulation>(
277 typename dealii::Triangulation<dim>::MeshSmoothing(
278 dealii::Triangulation<dim>::smoothing_on_refinement |
279 dealii::Triangulation<dim>::smoothing_on_coarsening));
280 dealii::GridGenerator::hyper_cube (*grid, left, right, colorize);
281 grid->refine_global(igrid);
282 pcout<<
" made grid for Index"<<igrid<<std::endl;
286 dealii::GridTools::transform (&warp<dim>, *grid);
290 unsigned int manifold_id=0;
291 grid->reset_all_manifolds();
292 grid->set_all_manifold_ids(manifold_id);
293 grid->set_manifold ( manifold_id, curv_manifold );
297 OPERATOR::OperatorsBaseState<dim,real,nstate,2*dim> operators(&all_parameters_new, poly_degree, poly_degree);
300 dg->allocate_system ();
302 dealii::IndexSet locally_owned_dofs;
303 dealii::IndexSet ghost_dofs;
304 dealii::IndexSet locally_relevant_dofs;
305 locally_owned_dofs = dg->dof_handler.locally_owned_dofs();
306 dealii::DoFTools::extract_locally_relevant_dofs(dg->dof_handler, ghost_dofs);
307 locally_relevant_dofs = ghost_dofs;
308 ghost_dofs.subtract_set(locally_owned_dofs);
309 dealii::LinearAlgebra::distributed::Vector<double> solution_deriv;
310 solution_deriv.reinit(locally_owned_dofs, ghost_dofs, MPI_COMM_WORLD);
313 const unsigned int max_dofs_per_cell = dg->dof_handler.get_fe_collection().max_dofs_per_cell();
314 std::vector<dealii::types::global_dof_index> current_dofs_indices(max_dofs_per_cell);
315 const unsigned int n_dofs_cell = dg->operators->fe_collection_basis[poly_degree].dofs_per_cell;
316 const unsigned int n_quad_pts = dg->operators->volume_quadrature_collection[poly_degree].size();
318 const dealii::FESystem<dim> &fe_metric = (dg->high_order_grid->fe_system);
319 const unsigned int n_metric_dofs = fe_metric.dofs_per_cell;
320 auto metric_cell = dg->high_order_grid->dof_handler_grid.begin_active();
321 for (
auto current_cell = dg->dof_handler.begin_active(); current_cell!=dg->dof_handler.end(); ++current_cell, ++metric_cell) {
322 if (!current_cell->is_locally_owned())
continue;
324 std::vector<dealii::types::global_dof_index> current_metric_dofs_indices(n_metric_dofs);
325 metric_cell->get_dof_indices (current_metric_dofs_indices);
326 std::array<std::vector<real>,dim> mapping_support_points;
327 std::vector<std::vector<real>> phys_quad_pts(dim);
328 for(
int idim=0; idim<dim; idim++){
329 mapping_support_points[idim].resize(n_metric_dofs/dim);
330 phys_quad_pts[idim].resize(n_quad_pts);
332 dealii::QGaussLobatto<dim> vol_GLL(grid_degree +1);
333 for (
unsigned int igrid_node = 0; igrid_node< n_metric_dofs/dim; ++igrid_node) {
334 for (
unsigned int idof = 0; idof< n_metric_dofs; ++idof) {
335 const real val = (dg->high_order_grid->volume_nodes[current_metric_dofs_indices[idof]]);
336 const unsigned int istate = fe_metric.system_to_component_index(idof).first;
337 mapping_support_points[istate][igrid_node] += val * fe_metric.shape_value_component(idof,vol_GLL.point(igrid_node),istate);
340 std::vector<dealii::FullMatrix<real>> metric_cofactor(n_quad_pts);
341 std::vector<real> determinant_Jacobian(n_quad_pts);
342 for(
unsigned int iquad=0;iquad<n_quad_pts; iquad++){
343 metric_cofactor[iquad].reinit(dim, dim);
345 operators.build_local_vol_metric_cofactor_matrix_and_det_Jac(grid_degree, poly_degree, n_quad_pts, n_metric_dofs/dim, mapping_support_points, determinant_Jacobian, metric_cofactor);
348 std::array<std::array<dealii::FullMatrix<real>,dim>,nstate> physical_gradient;
349 for(
unsigned int istate=0; istate<nstate; istate++){
350 for(
int idim=0; idim<dim; idim++){
351 physical_gradient[istate][idim].reinit(n_quad_pts, n_quad_pts);
354 operators.get_Jacobian_scaled_physical_gradient(
true, operators.gradient_flux_basis[poly_degree], metric_cofactor, n_quad_pts, physical_gradient);
357 current_dofs_indices.resize(n_dofs_cell);
358 current_cell->get_dof_indices (current_dofs_indices);
359 for(
int idim=0; idim<dim; idim++){
360 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
361 phys_quad_pts[idim][iquad] = 0.0;
362 for(
unsigned int jquad=0; jquad<n_metric_dofs/dim; jquad++){
363 phys_quad_pts[idim][iquad] += operators.mapping_shape_functions_vol_flux_nodes[grid_degree][poly_degree][iquad][jquad]
364 * mapping_support_points[idim][jquad];
368 std::vector<real> soln(n_quad_pts);
369 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
371 for (
int idim=0; idim<dim; idim++){
372 exact *= exp(-(phys_quad_pts[idim][iquad])*(phys_quad_pts[idim][iquad]));
377 dealii::Vector<real> soln_derivative_x(n_quad_pts);
378 for(
int istate=0; istate<nstate; istate++){
379 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
380 soln_derivative_x[iquad]=0.0;
381 for(
unsigned int idof=0; idof<n_quad_pts; idof++){
382 soln_derivative_x[iquad] += physical_gradient[istate][dim_check][iquad][idof] * soln[idof];
384 soln_derivative_x[iquad] /= determinant_Jacobian[iquad];
387 for(
unsigned int idof=0; idof<n_dofs_cell; idof++){
388 solution_deriv[current_dofs_indices[idof]] = 0.0;
389 for(
unsigned int iquad=0; iquad<n_quad_pts; iquad++){
390 solution_deriv[current_dofs_indices[idof]] += operators.vol_projection_operator[poly_degree][idof][iquad]
391 * soln_derivative_x[iquad];
398 pcout<<
"OOA here"<<std::endl;
399 double l2error = 0.0;
400 double linf_error = 0.0;
401 int overintegrate = 4;
402 dealii::QGauss<dim> quad_extra(poly_degree+1+overintegrate);
403 dealii::FEValues<dim,dim> fe_values_extra(*(dg->high_order_grid->mapping_fe_field), dg->operators->fe_collection_basis[poly_degree], quad_extra,
404 dealii::update_values | dealii::update_JxW_values |
405 dealii::update_jacobians |
406 dealii::update_quadrature_points | dealii::update_inverse_jacobians);
407 const unsigned int n_quad_pts_extra = fe_values_extra.n_quadrature_points;
408 std::vector<dealii::types::global_dof_index> dofs_indices (fe_values_extra.dofs_per_cell);
409 dealii::Vector<real> soln_at_q(n_quad_pts_extra);
410 for (
auto current_cell = dg->dof_handler.begin_active(); current_cell!=dg->dof_handler.end(); ++current_cell) {
411 if (!current_cell->is_locally_owned())
continue;
412 fe_values_extra.reinit(current_cell);
413 dofs_indices.resize(fe_values_extra.dofs_per_cell);
414 current_cell->get_dof_indices (dofs_indices);
416 for (
unsigned int iquad=0; iquad<n_quad_pts_extra; ++iquad){
417 soln_at_q[iquad] = 0.0;
418 for (
unsigned int idof=0; idof<fe_values_extra.dofs_per_cell; ++idof){
419 soln_at_q[iquad] += solution_deriv[dofs_indices[idof]] * fe_values_extra.shape_value_component(idof, iquad, 0);
421 const dealii::Point<dim> qpoint = (fe_values_extra.quadrature_point(iquad));
423 for(
int idim=0; idim<dim; idim++){
424 uexact_x *= exp(-((qpoint[idim]) * (qpoint[idim])));
426 uexact_x *= - 2.0 * qpoint[dim_check];
427 l2error += pow(soln_at_q[iquad] - uexact_x, 2) * fe_values_extra.JxW(iquad);
428 double inf_temp = std::abs(soln_at_q[iquad]-uexact_x);
429 if(inf_temp > linf_error){
430 linf_error = inf_temp;
435 pcout<<
"got OOA here"<<std::endl;
437 const unsigned int n_global_active_cells = grid->n_global_active_cells();
438 const double l2error_mpi_sum = std::sqrt(dealii::Utilities::MPI::sum(l2error, MPI_COMM_WORLD));
439 const double linferror_mpi= (dealii::Utilities::MPI::max(linf_error, MPI_COMM_WORLD));
441 const unsigned int n_dofs = dg->dof_handler.n_dofs();
442 const double dx = 1.0/pow(n_dofs,(1.0/dim));
443 grid_size[igrid] = dx;
444 soln_error[igrid] = l2error_mpi_sum;
445 soln_error_inf[igrid] = linferror_mpi;
447 convergence_table.add_value(
"p", poly_degree);
448 convergence_table.add_value(
"cells", n_global_active_cells);
449 convergence_table.add_value(
"DoFs", n_dofs);
450 convergence_table.add_value(
"dx", dx);
451 convergence_table.add_value(
"soln_L2_error", l2error_mpi_sum);
452 convergence_table.add_value(
"soln_Linf_error", linferror_mpi);
454 pcout <<
" Grid size h: " << dx
455 <<
" L2-soln_error: " << l2error_mpi_sum
456 <<
" Linf-soln_error: " << linferror_mpi
459 if (igrid > igrid_start) {
460 const double slope_soln_err = log(soln_error[igrid]/soln_error[igrid-1])
461 / log(grid_size[igrid]/grid_size[igrid-1]);
462 const double slope_soln_err_inf = log(soln_error_inf[igrid]/soln_error_inf[igrid-1])
463 / log(grid_size[igrid]/grid_size[igrid-1]);
464 pcout <<
"From grid " << igrid-1
465 <<
" to grid " << igrid
466 <<
" dimension: " << dim
467 <<
" polynomial degree p: " << poly_degree
469 <<
" solution_error1 " << soln_error[igrid-1]
470 <<
" solution_error2 " << soln_error[igrid]
471 <<
" slope " << slope_soln_err
472 <<
" solution_error1_inf " << soln_error_inf[igrid-1]
473 <<
" solution_error2_inf " << soln_error_inf[igrid]
474 <<
" slope " << slope_soln_err_inf
479 const int igrid = n_grids-1;
480 const double slope_soln_err = log(soln_error[igrid]/soln_error[igrid-1])
481 / log(grid_size[igrid]/grid_size[igrid-1]);
482 if(std::abs(slope_soln_err-poly_degree)>0.05){
485 pcout <<
" ********************************************" 487 <<
" Convergence rates for p = " << poly_degree
489 <<
" ********************************************" 491 convergence_table.evaluate_convergence_rates(
"soln_L2_error",
"cells", dealii::ConvergenceTable::reduction_rate_log2, dim);
492 convergence_table.evaluate_convergence_rates(
"soln_Linf_error",
"cells", dealii::ConvergenceTable::reduction_rate_log2, dim);
493 convergence_table.set_scientific(
"dx",
true);
494 convergence_table.set_scientific(
"soln_L2_error",
true);
495 convergence_table.set_scientific(
"soln_Linf_error",
true);
496 if (pcout.is_active()) convergence_table.write_text(pcout.get_stream());
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.
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.
bool use_curvilinear_split_form
Flag to use curvilinear metric split form.