Amandus: Simulations based on multilevel Schwarz methods
visualize_solution.h
Go to the documentation of this file.
1 /**********************************************************************
2  * Copyright (C) 2011 - 2014 by the authors
3  * Distributed under the MIT License
4  *
5  * See the files AUTHORS and LICENSE in the project root directory
6  **********************************************************************/
7 #include <deal.II/base/function.h>
8 #include <deal.II/base/parameter_handler.h>
9 #include <deal.II/base/quadrature.h>
10 #include <deal.II/dofs/dof_handler.h>
11 #include <deal.II/dofs/dof_tools.h>
12 #include <deal.II/lac/constraint_matrix.h>
13 #include <deal.II/lac/vector.h>
14 #include <deal.II/numerics/data_out.h>
15 #include <deal.II/numerics/vector_tools.h>
16 
17 #include <fstream>
18 #include <sstream>
19 #include <string>
20 #include <vector>
21 
22 namespace debug
23 {
24 template <int dim>
25 void
26 output_solution(const dealii::Function<dim>& function, const dealii::DoFHandler<dim>& dof_handler,
27  const dealii::Quadrature<dim>& quadrature, dealii::ParameterHandler& param)
28 {
29  dealii::Vector<double> projection(dof_handler.n_dofs());
30 
31  dealii::ConstraintMatrix constraints;
32  dealii::DoFTools::make_hanging_node_constraints(dof_handler, constraints);
33  constraints.close();
34 
35  dealii::VectorTools::project(dof_handler, constraints, quadrature, function, projection);
36 
37  std::vector<std::string> names;
38  std::string basename = "component";
39  for (unsigned int component = 0; component < function.n_components; ++component)
40  {
41  std::ostringstream name;
42  name << basename << component;
43  names.push_back(name.str());
44  }
45 
46  dealii::DataOut<dim> data_out;
47  data_out.parse_parameters(param);
48 
49  data_out.attach_dof_handler(dof_handler);
50  data_out.add_data_vector(projection, names);
51  data_out.build_patches(dof_handler.get_fe().tensor_degree());
52 
53  std::ostringstream filename;
54  filename << "projected_solution" << data_out.default_suffix();
55  std::ofstream output(filename.str().c_str());
56  data_out.write(output);
57 }
58 
59 template <int dim>
60 void
61 output_fe_data(const dealii::DoFHandler<dim>& dofh, const dealii::Vector<double> dof_vector,
62  const std::string& name, dealii::ParameterHandler& param)
63 {
64  std::vector<std::string> names;
65  for (unsigned int c = 0; c < dealii::DoFTools::n_components(dofh); ++c)
66  {
67  std::ostringstream c_name;
68  c_name << name << "component" << c;
69  names.push_back(c_name.str());
70  }
71 
72  dealii::DataOut<dim> data_out;
73  data_out.parse_parameters(param);
74 
75  data_out.attach_dof_handler(dofh);
76  data_out.add_data_vector(dof_vector, names);
77  data_out.build_patches(dofh.get_fe().tensor_degree());
78 
79  std::ostringstream filename;
80  filename << name << data_out.default_suffix();
81  std::ofstream output(filename.str().c_str());
82  data_out.write(output);
83 }
84 
85 template <int dim>
86 void
87 output_errors(dealii::BlockVector<double> errors, const dealii::Triangulation<dim>& tria,
88  unsigned int s, std::string prefix)
89 {
90  dealii::DataOut<dim> data_out;
91  data_out.attach_triangulation(tria);
92 
93  for (std::size_t i = 0; i < errors.n_blocks(); ++i)
94  {
95  std::ostringstream name;
96  name << prefix << "_error" << i;
97  data_out.add_data_vector(errors.block(i), name.str(), dealii::DataOut<dim>::type_cell_data);
98  }
99  data_out.build_patches();
100 
101  std::ostringstream filename;
102  filename << prefix << "_errors-" << s << ".vtk";
103  std::ofstream out_file(filename.str().c_str());
104  data_out.write_vtk(out_file);
105 }
106 }
void output_fe_data(const dealii::DoFHandler< dim > &dofh, const dealii::Vector< double > dof_vector, const std::string &name, dealii::ParameterHandler &param)
Definition: visualize_solution.h:61
void output_errors(dealii::BlockVector< double > errors, const dealii::Triangulation< dim > &tria, unsigned int s, std::string prefix)
Definition: visualize_solution.h:87
void output_solution(const dealii::Function< dim > &function, const dealii::DoFHandler< dim > &dof_handler, const dealii::Quadrature< dim > &quadrature, dealii::ParameterHandler &param)
Definition: visualize_solution.h:26
Definition: visualize_solution.h:22