Reference documentation for deal.II version 8.5.1
output.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2015 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii__mesh_worker_output_h
18 #define dealii__mesh_worker_output_h
19 
20 #include <deal.II/meshworker/dof_info.h>
21 #include <deal.II/base/smartpointer.h>
22 #include <deal.II/base/utilities.h>
23 #include <deal.II/lac/block_vector.h>
24 #include <deal.II/base/mg_level_object.h>
25 
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 namespace MeshWorker
30 {
31  namespace Assembler
32  {
33 
58  {
59  public:
63  GnuplotPatch();
64 
74  void initialize (const unsigned int n_points,
75  const unsigned int n_vectors);
76 
81  void initialize_stream (std::ostream &stream);
82 
90  template <int dim>
91  void initialize_info(DoFInfo<dim> &info, bool face);
92 
96  template<int dim>
97  void assemble(const DoFInfo<dim> &info);
98 
102  template<int dim>
103  void assemble(const DoFInfo<dim> &info1,
104  const DoFInfo<dim> &info2);
105 
106  private:
111  template<typename T>
112  void write(const T &t) const;
113 
119  void write_endl () const;
120 
124  unsigned int n_vectors;
128  unsigned int n_points;
129 
133  std::ostream *os;
134  };
135 
136 //----------------------------------------------------------------------//
137 
138  template <typename T>
139  inline void
140  GnuplotPatch::write(const T &d) const
141  {
142  if (os == 0)
143  deallog << d;
144  else
145  (*os) << d;
146  }
147 
148 
149  inline void
151  {
152  if (os == 0)
153  deallog << std::endl;
154  else
155  (*os) << std::endl;
156  }
157 
158 
159  inline
161  :
162  os(0)
163  {}
164 
165 
166  inline void
167  GnuplotPatch::initialize (const unsigned int np,
168  const unsigned int nv)
169  {
170  n_vectors = nv;
171  n_points = np;
172  }
173 
174 
175  inline void
176  GnuplotPatch::initialize_stream (std::ostream &stream)
177  {
178  os = &stream;
179  }
180 
181 
182  template <int dim>
183  inline void
185  {
186  if (face)
187  info.initialize_quadrature(Utilities::fixed_power<dim-1>(n_points), n_vectors+dim);
188  else
189  info.initialize_quadrature(Utilities::fixed_power<dim>(n_points), n_vectors+dim);
190  }
191 
192 
193  template <int dim>
194  inline void
196  {
197  const unsigned int np = info.n_quadrature_points();
198  const unsigned int nv = info.n_quadrature_values();
199  const unsigned int patch_dim = (info.face_number == numbers::invalid_unsigned_int)
200  ? dim : (dim-1);
201  const unsigned int row_length = n_points;
202  // If patches are 1D, end the
203  // patch after a row, else end
204  // it after a square
205  const unsigned int row_length2 = (patch_dim==1) ? row_length : (row_length*row_length);
206 
207 // AssertDimension(np, Utilities::fixed_power<dim>(n_points));
208  AssertDimension(nv, n_vectors+dim);
209 
210 
211  for (unsigned int k=0; k<np; ++k)
212  {
213  if (k % row_length == 0)
214  write_endl();
215  if (k % row_length2 == 0)
216  write_endl();
217 
218  for (unsigned int i=0; i<nv; ++i)
219  {
220  write(info.quadrature_value(k,i));
221  write('\t');
222  }
223  write_endl();
224  }
225  }
226 
227 
228  template <int dim>
229  inline void
231  {
232  assemble(info1);
233  assemble(info2);
234  }
235  }
236 }
237 
238 DEAL_II_NAMESPACE_CLOSE
239 
240 #endif
static const unsigned int invalid_unsigned_int
Definition: types.h:170
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
unsigned int n_quadrature_values() const
void assemble(const DoFInfo< dim > &info)
Definition: output.h:195
void initialize_quadrature(unsigned int np, unsigned int nv)
number & quadrature_value(unsigned int k, unsigned int i)
void write(const T &t) const
Definition: output.h:140
void initialize_info(DoFInfo< dim > &info, bool face)
Definition: output.h:184
void initialize(const unsigned int n_points, const unsigned int n_vectors)
Definition: output.h:167
void initialize_stream(std::ostream &stream)
Definition: output.h:176
unsigned int face_number
Definition: dof_info.h:83
unsigned int n_quadrature_points() const