deal.II version GIT relicensing-2291-g668cd86249 2024-12-24 11:30:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
output.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2010 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
16#ifndef dealii_mesh_worker_output_h
17#define dealii_mesh_worker_output_h
18
19#include <deal.II/base/config.h>
20
24
26
28
29
31
32namespace MeshWorker
33{
34 namespace Assembler
35 {
57 {
58 public:
63
73 void
74 initialize(const unsigned int n_points, const unsigned int n_vectors);
75
80 void
81 initialize_stream(std::ostream &stream);
82
90 template <int dim>
91 void
92 initialize_info(DoFInfo<dim> &info, bool face);
93
97 template <int dim>
98 void
99 assemble(const DoFInfo<dim> &info);
100
104 template <int dim>
105 void
106 assemble(const DoFInfo<dim> &info1, const DoFInfo<dim> &info2);
107
108 private:
113 template <typename T>
114 void
115 write(const T &t) const;
116
122 void
123 write_endl() const;
124
128 unsigned int n_vectors;
132 unsigned int n_points;
133
137 std::ostream *os;
138 };
139
140 //----------------------------------------------------------------------//
141
142 template <typename T>
143 inline void
144 GnuplotPatch::write(const T &d) const
145 {
146 if (os == nullptr)
147 deallog << d;
148 else
149 (*os) << d;
150 }
151
152
153 inline void
155 {
156 if (os == nullptr)
157 deallog << std::endl;
158 else
159 (*os) << std::endl;
160 }
161
162
164 : n_vectors(numbers::invalid_unsigned_int)
165 , n_points(numbers::invalid_unsigned_int)
166 , os(nullptr)
167 {}
168
169
170 inline void
171 GnuplotPatch::initialize(const unsigned int np, const unsigned int nv)
172 {
173 n_vectors = nv;
174 n_points = np;
175 }
176
177
178 inline void
180 {
181 os = &stream;
182 }
183
184
185 template <int dim>
186 inline void
188 {
189 if (face)
190 info.initialize_quadrature(Utilities::fixed_power<dim - 1>(n_points),
191 n_vectors + dim);
192 else
193 info.initialize_quadrature(Utilities::fixed_power<dim>(n_points),
194 n_vectors + dim);
195 }
196
197
198 template <int dim>
199 inline void
201 {
202 const unsigned int np = info.n_quadrature_points();
203 const unsigned int nv = info.n_quadrature_values();
204 const unsigned int patch_dim =
205 (info.face_number == numbers::invalid_unsigned_int) ? dim : (dim - 1);
206 const unsigned int row_length = n_points;
207 // If patches are 1d, end the
208 // patch after a row, else end
209 // it after a square
210 const unsigned int row_length2 =
211 (patch_dim == 1) ? row_length : (row_length * row_length);
212
213 // AssertDimension(np, Utilities::fixed_power<dim>(n_points));
214 AssertDimension(nv, n_vectors + dim);
215
216
217 for (unsigned int k = 0; k < np; ++k)
218 {
219 if (k % row_length == 0)
220 write_endl();
221 if (k % row_length2 == 0)
222 write_endl();
223
224 for (unsigned int i = 0; i < nv; ++i)
225 {
226 write(info.quadrature_value(k, i));
227 write('\t');
228 }
229 write_endl();
230 }
231 }
232
233
234 template <int dim>
235 inline void
237 {
238 assemble(info1);
239 assemble(info2);
240 }
241 } // namespace Assembler
242} // namespace MeshWorker
243
245
246#endif
void write(const T &t) const
Definition output.h:144
void initialize_stream(std::ostream &stream)
Definition output.h:179
void assemble(const DoFInfo< dim > &info)
Definition output.h:200
void initialize_info(DoFInfo< dim > &info, bool face)
Definition output.h:187
void initialize(const unsigned int n_points, const unsigned int n_vectors)
Definition output.h:171
unsigned int face_number
Definition dof_info.h:89
unsigned int n_quadrature_values() const
number & quadrature_value(const unsigned int k, const unsigned int i)
void initialize_quadrature(const unsigned int np, const unsigned int nv)
unsigned int n_quadrature_points() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#define AssertDimension(dim1, dim2)
LogStream deallog
Definition logstream.cc:36
static const unsigned int invalid_unsigned_int
Definition types.h:220