From 1c4be894f28598f0a4aefc6b278a3a7d766ec2c5 Mon Sep 17 00:00:00 2001 From: hartmann Date: Wed, 12 Sep 2001 12:14:09 +0000 Subject: [PATCH] same changes. output of data in tables. output of grids in gnuplot files. git-svn-id: https://svn.dealii.org/trunk@4966 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/examples/step-10/step-10.cc | 112 +++++++++++++++++++++------- 1 file changed, 85 insertions(+), 27 deletions(-) diff --git a/deal.II/examples/step-10/step-10.cc b/deal.II/examples/step-10/step-10.cc index 4b36f1de1d..12f35c6f52 100644 --- a/deal.II/examples/step-10/step-10.cc +++ b/deal.II/examples/step-10/step-10.cc @@ -2,20 +2,24 @@ /* Author: Wolfgang Bangerth, University of Heidelberg, 1999 */ #include -#include -#include +#include #include #include #include +#include +#include +#include +#include #include #include -#include #include -#include #include +#include +static const long double pi=3.141592653589793238462643; + template void compute_pi_by_area () @@ -23,9 +27,10 @@ void compute_pi_by_area () std::cout << "Computation of Pi by the area:" << std::endl << "==============================" << std::endl; + const QGauss4 quadrature; for (unsigned int order=1; order<5; ++order) { - cout << "Order = " << order << endl; + cout << "Order = " << order << std::endl; Triangulation triangulation; GridGenerator::hyper_ball (triangulation); @@ -37,29 +42,39 @@ void compute_pi_by_area () DoFHandler dof_handler (triangulation); - for (unsigned int refinement=0; refinement<5; ++refinement) + FEValues fe_values (mapping, fe, quadrature, update_JxW_values); + ConvergenceTable table; + + for (unsigned int refinement=0; refinement<6; + ++refinement, triangulation.refine_global (1)) { - triangulation.refine_global (1); - dof_handler.distribute_dofs (fe); - - QGauss4 quadrature; - FEValues fe_values (mapping, fe, quadrature, update_JxW_values); + table.add_value("cells", triangulation.n_active_cells()); + + dof_handler.distribute_dofs (fe); typename DoFHandler::active_cell_iterator cell = dof_handler.begin_active(), endc = dof_handler.end(); - double area = 0; + long double area = 0; for (; cell!=endc; ++cell) { fe_values.reinit (cell); for (unsigned int i=0; i (area)); + table.add_value("error", fabs(area-pi)); }; - std::cout << std::endl; + + table.omit_column_from_convergence_rate_evaluation("cells"); + table.omit_column_from_convergence_rate_evaluation("eval.pi"); + table.evaluate_all_convergence_rates( + ConvergenceTable::reduction_rate_log2); + + table.set_precision("eval.pi", 16); + table.set_scientific("error", true); + + table.write_text(cout); }; }; @@ -71,6 +86,7 @@ void compute_pi_by_perimeter () std::cout << "Computation of Pi by the perimeter:" << std::endl << "===================================" << std::endl; + const QGauss4 quadrature; for (unsigned int order=1; order<5; ++order) { cout << "Order = " << order << endl; @@ -84,19 +100,21 @@ void compute_pi_by_perimeter () const FE_Q fe (1); DoFHandler dof_handler (triangulation); - - for (unsigned int refinement=0; refinement<5; ++refinement) + + FEFaceValues fe_face_values (mapping, fe, quadrature, update_JxW_values); + ConvergenceTable table; + + for (unsigned int refinement=0; refinement<6; + ++refinement, triangulation.refine_global (1)) { - triangulation.refine_global (1); + table.add_value("cells", triangulation.n_active_cells()); + dof_handler.distribute_dofs (fe); - QGauss4 quadrature; - FEFaceValues fe_face_values (mapping, fe, quadrature, update_JxW_values); - typename DoFHandler::active_cell_iterator cell = dof_handler.begin_active(), endc = dof_handler.end(); - double perimeter = 0; + long double perimeter = 0; for (; cell!=endc; ++cell) for (unsigned int face_no=0; face_no::faces_per_cell; ++face_no) if (cell->face(face_no)->at_boundary()) @@ -105,21 +123,61 @@ void compute_pi_by_perimeter () for (unsigned int i=0; i (perimeter/2.)); + table.add_value("error", fabs(perimeter/2.-pi)); }; - std::cout << std::endl; + + table.omit_column_from_convergence_rate_evaluation("cells"); + table.omit_column_from_convergence_rate_evaluation("eval.pi"); + table.evaluate_all_convergence_rates( + ConvergenceTable::reduction_rate_log2); + + table.set_precision("eval.pi", 16); + table.set_scientific("error", true); + + table.write_text(cout); }; }; +template +void gnuplot_output() +{ + cout << "Output of grids as gnuplot file:" << std::endl; + Triangulation triangulation; + GridGenerator::hyper_ball (triangulation); + + static const HyperBallBoundary boundary; + triangulation.set_boundary (0, boundary); + + GridOut grid_out; + // on boundary faces plot 30 + // additional points per face. + GridOutFlags::Gnuplot gnuplot_flags(false, 30); + grid_out.set_flags(gnuplot_flags); + + + for (unsigned int order=1; order<4; ++order) + { + cout << "Order = " << order << std::endl; + + const MappingQ mapping (order); + string filename="ball_mapping_q"; + filename += ('0'+order); + filename += ".dat"; + ofstream gnuplot_file(filename.c_str()); + + grid_out.write_gnuplot(triangulation, gnuplot_file, &mapping); + } +} int main () { cout.precision (16); + gnuplot_output<2>(); + compute_pi_by_area<2> (); compute_pi_by_perimeter<2> (); -- 2.39.5