/usr/include/deal.II/numerics/data_out.h is in libdeal.ii-dev 8.5.1-3.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 | // ---------------------------------------------------------------------
//
// Copyright (C) 1999 - 2016 by the deal.II authors
//
// This file is part of the deal.II library.
//
// The deal.II library is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the deal.II distribution.
//
// ---------------------------------------------------------------------
#ifndef dealii__data_out_h
#define dealii__data_out_h
#include <deal.II/base/config.h>
#include <deal.II/numerics/data_out_dof_data.h>
#include <deal.II/base/std_cxx11/shared_ptr.h>
DEAL_II_NAMESPACE_OPEN
namespace internal
{
namespace DataOut
{
/**
* A derived class for use in the DataOut class. This is a class for the
* AdditionalData kind of data structure discussed in the documentation of
* the WorkStream context.
*/
template <int dim, int spacedim>
struct ParallelData : public ParallelDataBase<dim,spacedim>
{
ParallelData (const unsigned int n_datasets,
const unsigned int n_subdivisions,
const std::vector<unsigned int> &n_postprocessor_outputs,
const Mapping<dim,spacedim> &mapping,
const std::vector<std_cxx11::shared_ptr<dealii::hp::FECollection<dim,spacedim> > > &finite_elements,
const UpdateFlags update_flags,
const std::vector<std::vector<unsigned int> > &cell_to_patch_index_map);
std::vector<Point<spacedim> > patch_evaluation_points;
const std::vector<std::vector<unsigned int> > *cell_to_patch_index_map;
};
}
}
/**
* This class is the main class to provide output of data described by finite
* element fields defined on a collection of cells.
*
* This class is an actual implementation of the functionality proposed by the
* DataOut_DoFData class. It offers a function build_patches() that generates
* the data to be written in some graphics format. Most of the interface and
* an example of its use is described in the documentation of this base class.
*
* The only thing this class offers is the function build_patches() which
* loops over all cells of the triangulation stored by the
* attach_dof_handler() function of the base class (with the exception of
* cells of parallel::distributed::Triangulation objects that are not owned by
* the current processor) and converts the data on these to actual patches
* which are the objects that are later output by the functions of the base
* classes. You can give a parameter to the function which determines how many
* subdivisions in each coordinate direction are to be performed, i.e. of how
* many subcells each patch shall consist. Default is one, but you may want to
* choose a higher number for higher order elements, for example two for
* quadratic elements, three for cubic elements three, and so on. The purpose
* of this parameter is because most graphics programs do not allow to specify
* higher order polynomial functions in the file formats: only data at
* vertices can be plotted and is then shown as a bilinear interpolation
* within the interior of cells. This may be insufficient if you have higher
* order finite elements, and the only way to achieve better output is to
* subdivide each cell of the mesh into several cells for graphical output. Of
* course, what you get to see is still a bilinear interpolation on each cell
* of the output (where these cells are not subdivisions of the cells of the
* triangulation in use) due to the same limitations in output formats, but at
* least a bilinear interpolation of a higher order polynomial on a finer
* mesh.
*
* Note that after having called build_patches() once, you can call one or
* more of the write() functions of DataOutInterface. You can therefore output
* the same data in more than one format without having to rebuild the
* patches.
*
*
* <h3>User interface information</h3>
*
* The base classes of this class, DataOutBase, DataOutInterface and
* DataOut_DoFData offer several interfaces of their own. Refer to the
* DataOutBase class's documentation for a discussion of the different output
* formats presently supported, DataOutInterface for ways of selecting which
* format to use upon output at run-time and without the need to adapt your
* program when new formats become available, as well as for flags to
* determine aspects of output. The DataOut_DoFData() class's documentation
* has an example of using nodal data to generate output.
*
*
* <h3>Extensions</h3>
*
* By default, this class produces patches for all active cells. Sometimes,
* this is not what you want, maybe because they are simply too many (and too
* small to be seen individually) or because you only want to see a certain
* region of the domain (for example in parallel programs such as the step-18
* example program), or for some other reason.
*
* For this, internally build_patches() does not generate the sequence of
* cells to be converted into patches itself, but relies on the two functions
* first_cell() and next_cell(). By default, they return the first active
* cell, and the next active cell, respectively. Since they are @p virtual
* functions, you can write your own class derived from DataOut in which you
* overload these two functions to select other cells for output. This may,
* for example, include only cells that are in parts of a domain (e.g., if you
* don't care about the solution elsewhere, think for example a buffer region
* in which you attenuate outgoing waves in the Perfectly Matched Layer
* method) or if you don't want output to be generated at all levels of an
* adaptively refined mesh because this creates too much data (in this case,
* the set of cells returned by your implementations of first_cell() and
* next_cell() will include non-active cells, and DataOut::build_patches()
* will simply take interpolated values of the solution instead of the exact
* values on these cells children for output). Once you derive your own class,
* you would just create an object of this type instead of an object of type
* DataOut, and everything else will remain the same.
*
* The two functions are not constant, so you may store information within
* your derived class about the last accessed cell. This is useful if the
* information of the last cell which was accessed is not sufficient to
* determine the next one.
*
* There is one caveat, however: if you have cell data (in contrast to nodal,
* or dof, data) such as error indicators, then you must make sure that
* first_cell() and next_cell() only walk over active cells, since cell data
* cannot be interpolated to a coarser cell. If you do have cell data and use
* this pair of functions and they return a non-active cell, then an exception
* will be thrown.
*
* @pre This class only makes sense if the first template argument,
* <code>dim</code> equals the dimension of the DoFHandler type given as the
* second template argument, i.e., if <code>dim ==
* DoFHandlerType::dimension</code>. This redundancy is a historical relic
* from the time where the library had only a single DoFHandler class and this
* class consequently only a single template argument.
*
* @ingroup output
* @author Wolfgang Bangerth, 1999
*/
template <int dim, typename DoFHandlerType=DoFHandler<dim> >
class DataOut : public DataOut_DoFData<DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension>
{
public:
/**
* Typedef to the iterator type of the dof handler class under
* consideration.
*/
typedef typename DataOut_DoFData<DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension>::cell_iterator
cell_iterator;
typedef typename DataOut_DoFData<DoFHandlerType, DoFHandlerType::dimension, DoFHandlerType::space_dimension>::active_cell_iterator
active_cell_iterator;
/**
* Enumeration describing the part of the domain in which cells
* should be written with curved boundaries. In reality, no file
* format we are aware of really supports curved boundaries, but
* this can be emulated by plotting edges as a sequence of straight
* lines (and faces in 3d as a collection of bilinear patches) if
* DataOut::build_patches() is called with a number of subdivisions
* greater than 1.
*
* The elements of this enumeration then describe for which cells
* DataOut::build_patches() will query the manifold or boundary
* description for curved geometries.
*/
enum CurvedCellRegion
{
/**
* The geometry or boundary description
* will never be queried for curved geometries. This means that even
* if you have more than one subdivision per cell (see
* DataOut::build_patches() for what exactly this means) and even
* if the geometry really is curved, each cell will still be
* subdivided as if it was just a bi- or trilinear cell.
*/
no_curved_cells,
/**
* The geometry or boundary description
* will be queried for curved geometries for cells located
* at the boundary, i.e., for cells that have at least one
* face at the boundary. This is sufficient if you have not
* attached a manifold description to the interiors of cells
* but only to faces at the boundary.
*/
curved_boundary,
/**
* The geometry description will be
* queried for all cells and all faces, whether they are
* at the boundary or not. This option is appropriate if you
* have attached a manifold object to cells (not only to faces).
*/
curved_inner_cells
};
/**
* This is the central function of this class since it builds the list of
* patches to be written by the low-level functions of the base class. A
* patch is, in essence, some intermediate representation of the data on
* each cell of a triangulation and DoFHandler object that can then be used
* to write files in some format that is readable by visualization programs.
*
* You can find an overview of the use of this function in the general
* documentation of this class. An example is also provided in the
* documentation of this class's base class DataOut_DoFData.
*
* @param n_subdivisions A parameter that determines how many "patches" this
* function will build out of every cell. If you do not specify this value
* in calling, or provide the default value zero, then this is interpreted
* as DataOutInterface::default_subdivisions which most of the time will be
* equal to one (unless you have set it to something else). The purpose of
* this parameter is to subdivide each cell of the mesh into $2\times 2,
* 3\times 3, \ldots$ "patches" in 2d, and $2\times 2\times 2, 3\times
* 3\times 3, \ldots$ (if passed the value 2, 3, etc) where each patch
* represents the data from a regular subdivision of the cell into equal
* parts. Most of the times, this is not necessary and outputting one patch
* per cell is exactly what you want to plot the solution. That said, the
* data we write into files for visualization can only represent (bi-,
* tri)linear data on each cell, and most visualization programs can in fact
* only visualize this kind of data. That's good enough if you work with
* (bi-, tri)linear finite elements, in which case what you get to see is
* exactly what has been computed. On the other hand, if you work with (bi-,
* tri)quadratic elements, then what is written into the output file is just
* a (bi-, tri)linear interpolation onto the current mesh, i.e., only the
* values at the vertices. If this is not good enough, you can, for example,
* specify @p n_subdivisions equal to 2 to plot the solution on a once-
* refined mesh, or if set to 3, on a mesh where each cell is represented by
* 3-by-3 patches. On each of these smaller patches, given the limitations
* of output formats, the data is still linearly interpolated, but a linear
* interpolation of quadratic data on a finer mesh is still a better
* representation of the actual quadratic surface than on the original mesh.
* In other words, using this parameter can not help you plot the solution
* exactly, but it can get you closer if you use finite elements of higher
* polynomial degree.
*/
virtual void build_patches (const unsigned int n_subdivisions = 0);
/**
* Same as above, except that the additional first parameter defines a
* mapping that is to be used in the generation of output. If
* <tt>n_subdivisions>1</tt>, the points interior of subdivided patches
* which originate from cells at the boundary of the domain can be computed
* using the mapping, i.e., a higher order mapping leads to a representation
* of a curved boundary by using more subdivisions. Some mappings like
* MappingQEulerian result in curved cells in the interior of the domain.
* The same is true if you have attached a manifold description to
* the cells of a triangulation (see
* @ref manifold "Manifolds"
* for more information). However, there is no easy way to query the mapping
* or manifold whether it really does lead to curved cells.
* Thus the last argument @p curved_region takes one of three values
* resulting in no curved cells at all, curved cells at the boundary
* (default) or curved cells in the whole domain. For more information
* about these three options, see the CurvedCellRegion enum's
* description.
*
* Even for non-curved cells, the mapping argument can be used for
* Eulerian mappings (see class MappingQ1Eulerian) where a mapping is used
* not only to determine the position of points interior to a cell, but also
* of the vertices. It offers an opportunity to watch the solution on a
* deformed triangulation on which the computation was actually carried out,
* even if the mesh is internally stored in its undeformed configuration and
* the deformation is only tracked by an additional vector that holds the
* deformation of each vertex.
*
* @todo The @p mapping argument should be replaced by a
* hp::MappingCollection in case of a hp::DoFHandler.
*/
virtual void build_patches (const Mapping<DoFHandlerType::dimension, DoFHandlerType::space_dimension> &mapping,
const unsigned int n_subdivisions = 0,
const CurvedCellRegion curved_region = curved_boundary);
/**
* Return the first cell which we want output for. The default
* implementation returns the first active cell, but you might want to
* return other cells in a derived class.
*/
virtual cell_iterator first_cell ();
/**
* Return the next cell after @p cell which we want output for. If there
* are no more cells, <tt>#dofs->end()</tt> shall be returned.
*
* The default implementation returns the next active cell, but you might
* want to return other cells in a derived class. Note that the default
* implementation assumes that the given @p cell is active, which is
* guaranteed as long as first_cell() is also used from the default
* implementation. Overloading only one of the two functions might not be a
* good idea.
*/
virtual cell_iterator next_cell (const cell_iterator &cell);
private:
/**
* Return the first cell produced by the first_cell()/next_cell() function
* pair that is locally owned. If this object operates on a non-distributed
* triangulation, the result equals what first_cell() returns.
*/
virtual cell_iterator first_locally_owned_cell ();
/**
* Return the next cell produced by the next_cell() function that is locally
* owned. If this object operates on a non-distributed triangulation, the
* result equals what first_cell() returns.
*/
virtual cell_iterator next_locally_owned_cell (const cell_iterator &cell);
/**
* Build one patch. This function is called in a WorkStream context.
*
* The first argument here is the iterator, the second the scratch data
* object. All following are tied to particular values when calling
* WorkStream::run(). The function does not take a CopyData object but
* rather allocates one on its own stack for memory access efficiency
* reasons.
*/
void build_one_patch
(const std::pair<cell_iterator, unsigned int> *cell_and_index,
internal::DataOut::ParallelData<DoFHandlerType::dimension, DoFHandlerType::space_dimension> &scratch_data,
const unsigned int n_subdivisions,
const CurvedCellRegion curved_cell_region);
};
DEAL_II_NAMESPACE_CLOSE
#endif
|