Reference documentation for deal.II version GIT relicensing-249-g48dc7357c7 2024-03-29 12:30:02+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
Classes | Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | List of all members
internal::TriangulationImplementation::TriaObjects Class Reference

#include <deal.II/grid/tria_objects.h>

Classes

struct  BoundaryOrMaterialId
 
struct  UserData
 

Public Types

enum  UserDataType { data_unknown , data_pointer , data_index }
 

Public Member Functions

 TriaObjects ()
 
 TriaObjects (const unsigned int structdim)
 
unsigned int n_objects () const
 
ArrayView< intget_bounding_object_indices (const unsigned int index)
 
template<int structdim, int dim, int spacedim>
::TriaRawIterator<::TriaAccessor< structdim, dim, spacedim > > next_free_single_object (const Triangulation< dim, spacedim > &tria)
 
template<int structdim, int dim, int spacedim>
::TriaRawIterator<::TriaAccessor< structdim, dim, spacedim > > next_free_pair_object (const Triangulation< dim, spacedim > &tria)
 
template<int dim, int spacedim>
Triangulation< dim, spacedim >::raw_hex_iterator next_free_hex (const Triangulation< dim, spacedim > &tria, const unsigned int level)
 
void *& user_pointer (const unsigned int i)
 
const void * user_pointer (const unsigned int i) const
 
unsigned intuser_index (const unsigned int i)
 
unsigned int user_index (const unsigned int i) const
 
void clear_user_data (const unsigned int i)
 
void clear_user_data ()
 
void clear_user_flags ()
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
template<int structdim_, int dim, int spacedim>
::TriaRawIterator<::TriaAccessor< structdim_, dim, spacedim > > next_free_single_object (const Triangulation< dim, spacedim > &tria)
 
template<int structdim_, int dim, int spacedim>
::TriaRawIterator<::TriaAccessor< structdim_, dim, spacedim > > next_free_pair_object (const Triangulation< dim, spacedim > &tria)
 
template<int dim, int spacedim>
typename::Triangulation< dim, spacedim >::raw_hex_iterator next_free_hex (const ::Triangulation< dim, spacedim > &tria, const unsigned int level)
 

Static Public Member Functions

static ::ExceptionBaseExcPointerIndexClash ()
 

Public Attributes

unsigned int structdim
 
std::vector< intcells
 
std::vector< intchildren
 
std::vector< std::uint8_t > refinement_cases
 
std::vector< boolused
 
std::vector< booluser_flags
 
std::vector< BoundaryOrMaterialIdboundary_or_material_id
 
std::vector< types::manifold_idmanifold_id
 
unsigned int next_free_single
 
unsigned int next_free_pair
 
bool reverse_order_next_free_single
 
std::vector< UserDatauser_data
 
UserDataType user_data_type
 

Detailed Description

General template for information belonging to the geometrical objects of a triangulation, i.e. lines, quads, hexahedra... Apart from the vector of objects additional information is included, namely vectors indicating the children, the used-status, user-flags, material-ids..

Objects of these classes are included in the TriaLevel and TriaFaces classes.

Definition at line 52 of file tria_objects.h.

Member Enumeration Documentation

◆ UserDataType

Enum describing the possible types of userdata.

Enumerator
data_unknown 

No userdata used yet.

data_pointer 

UserData contains pointers.

data_index 

UserData contains indices.

Definition at line 345 of file tria_objects.h.

Constructor & Destructor Documentation

◆ TriaObjects() [1/2]

internal::TriangulationImplementation::TriaObjects::TriaObjects ( )
inline

Constructor resetting some data.

Definition at line 474 of file tria_objects.h.

◆ TriaObjects() [2/2]

internal::TriangulationImplementation::TriaObjects::TriaObjects ( const unsigned int  structdim)
inline

Constructor for a specific dimension.

Definition at line 483 of file tria_objects.h.

Member Function Documentation

◆ n_objects()

unsigned int internal::TriangulationImplementation::TriaObjects::n_objects ( ) const
inline

Return number of geometric objects stored by this class.

Definition at line 374 of file tria_objects.h.

◆ get_bounding_object_indices()

ArrayView< int > internal::TriangulationImplementation::TriaObjects::get_bounding_object_indices ( const unsigned int  index)
inline

Return a view on the indices of the objects that bound the index-th object stored by the current object. For example, if the current object stores cells, then this function returns the equivalent of an array containing the indices of the faces that bound the index-th cell.

Definition at line 385 of file tria_objects.h.

◆ next_free_single_object() [1/2]

template<int structdim, int dim, int spacedim>
::TriaRawIterator<::TriaAccessor< structdim, dim, spacedim > > internal::TriangulationImplementation::TriaObjects::next_free_single_object ( const Triangulation< dim, spacedim > &  tria)

Return an iterator to the next free slot for a single object. This function is only used by Triangulation::execute_refinement() in 3d.

Warning
Interestingly, this function is not used for 1d or 2d triangulations, where it seems the authors of the refinement function insist on reimplementing its contents.
Todo:
This function is not instantiated for the codim-one case

◆ next_free_pair_object() [1/2]

template<int structdim, int dim, int spacedim>
::TriaRawIterator<::TriaAccessor< structdim, dim, spacedim > > internal::TriangulationImplementation::TriaObjects::next_free_pair_object ( const Triangulation< dim, spacedim > &  tria)

Return an iterator to the next free slot for a pair of objects. This function is only used by Triangulation::execute_refinement() in 3d.

Warning
Interestingly, this function is not used for 1d or 2d triangulations, where it seems the authors of the refinement function insist on reimplementing its contents.
Todo:
This function is not instantiated for the codim-one case

◆ next_free_hex() [1/2]

template<int dim, int spacedim>
Triangulation< dim, spacedim >::raw_hex_iterator internal::TriangulationImplementation::TriaObjects::next_free_hex ( const Triangulation< dim, spacedim > &  tria,
const unsigned int  level 
)

Return an iterator to the next free slot for a pair of hexes. Only implemented for G=Hexahedron.

◆ user_pointer() [1/2]

void *& internal::TriangulationImplementation::TriaObjects::user_pointer ( const unsigned int  i)
inline

Access to user pointers.

Definition at line 431 of file tria_objects.h.

◆ user_pointer() [2/2]

const void * internal::TriangulationImplementation::TriaObjects::user_pointer ( const unsigned int  i) const
inline

Read-only access to user pointers.

Definition at line 443 of file tria_objects.h.

◆ user_index() [1/2]

unsigned int & internal::TriangulationImplementation::TriaObjects::user_index ( const unsigned int  i)
inline

Access to user indices.

Definition at line 455 of file tria_objects.h.

◆ user_index() [2/2]

unsigned int internal::TriangulationImplementation::TriaObjects::user_index ( const unsigned int  i) const
inline

Read-only access to user pointers.

Definition at line 493 of file tria_objects.h.

◆ clear_user_data() [1/2]

void internal::TriangulationImplementation::TriaObjects::clear_user_data ( const unsigned int  i)
inline

Reset user data to zero.

Definition at line 467 of file tria_objects.h.

◆ clear_user_data() [2/2]

void internal::TriangulationImplementation::TriaObjects::clear_user_data ( )
inline

Clear all user pointers or indices and reset their type, such that the next access may be either or.

Definition at line 505 of file tria_objects.h.

◆ clear_user_flags()

void internal::TriangulationImplementation::TriaObjects::clear_user_flags ( )
inline

Clear all user flags.

Definition at line 514 of file tria_objects.h.

◆ memory_consumption()

std::size_t internal::TriangulationImplementation::TriaObjects::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 65 of file tria_objects.cc.

◆ serialize()

template<class Archive >
void internal::TriangulationImplementation::TriaObjects::serialize ( Archive &  ar,
const unsigned int  version 
)

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

Definition at line 532 of file tria_objects.h.

◆ next_free_single_object() [2/2]

template<int structdim_, int dim, int spacedim>
::TriaRawIterator<::TriaAccessor< structdim_, dim, spacedim > > internal::TriangulationImplementation::TriaObjects::next_free_single_object ( const Triangulation< dim, spacedim > &  tria)

Definition at line 550 of file tria_objects.h.

◆ next_free_pair_object() [2/2]

template<int structdim_, int dim, int spacedim>
::TriaRawIterator<::TriaAccessor< structdim_, dim, spacedim > > internal::TriangulationImplementation::TriaObjects::next_free_pair_object ( const Triangulation< dim, spacedim > &  tria)

Definition at line 604 of file tria_objects.h.

◆ next_free_hex() [2/2]

template<int dim, int spacedim>
typename::Triangulation< dim, spacedim >::raw_hex_iterator internal::TriangulationImplementation::TriaObjects::next_free_hex ( const ::Triangulation< dim, spacedim > &  tria,
const unsigned int  level 
)

Definition at line 35 of file tria_objects.cc.

Member Data Documentation

◆ structdim

unsigned int internal::TriangulationImplementation::TriaObjects::structdim

Definition at line 65 of file tria_objects.h.

◆ cells

std::vector<int> internal::TriangulationImplementation::TriaObjects::cells

Vector of the objects belonging to this level. The index of the object equals the index in this container.

Definition at line 71 of file tria_objects.h.

◆ children

std::vector<int> internal::TriangulationImplementation::TriaObjects::children

Index of the even children of an object. Since when objects are refined, all children are created at the same time, they are appended to the list at least in pairs after each other. We therefore only store the index of the even children, the uneven follow immediately afterwards.

If an object has no children, -1 is stored in this list. An object is called active if it has no children. The function TriaAccessorBase::has_children() tests for this.

Definition at line 100 of file tria_objects.h.

◆ refinement_cases

std::vector<std::uint8_t> internal::TriangulationImplementation::TriaObjects::refinement_cases

Store the refinement case each of the cells is refined with. This vector might be replaced by vector<vector<bool> > (dim, vector<bool> (n_cells)) which is more memory efficient.

Definition at line 107 of file tria_objects.h.

◆ used

std::vector<bool> internal::TriangulationImplementation::TriaObjects::used

Vector storing whether an object is used in the cells vector.

Since it is difficult to delete elements in a vector, when an element is not needed any more (e.g. after derefinement), it is not deleted from the list, but rather the according used flag is set to false.

Definition at line 117 of file tria_objects.h.

◆ user_flags

std::vector<bool> internal::TriangulationImplementation::TriaObjects::user_flags

Make available a field for user data, one bit per object. This field is usually used when an operation runs over all cells and needs information whether another cell (e.g. a neighbor) has already been processed.

You can clear all used flags using Triangulation::clear_user_flags().

Definition at line 128 of file tria_objects.h.

◆ boundary_or_material_id

std::vector<BoundaryOrMaterialId> internal::TriangulationImplementation::TriaObjects::boundary_or_material_id

Store boundary and material data. For example, in one dimension, this field stores the material id of a line, which is a number between 0 and numbers::invalid_material_id-1. In more than one dimension, lines have no material id, but they may be at the boundary; then, we store the boundary indicator in this field, which denotes to which part of the boundary this line belongs and which boundary conditions hold on this part. The boundary indicator also is a number between zero and numbers::internal_face_boundary_id-1; the id numbers::internal_face_boundary_id is reserved for lines in the interior and may be used to check whether a line is at the boundary or not, which otherwise is not possible if you don't know which cell it belongs to.

Definition at line 179 of file tria_objects.h.

◆ manifold_id

std::vector<types::manifold_id> internal::TriangulationImplementation::TriaObjects::manifold_id

Store manifold ids. This field stores the manifold id of each object, which is a number between 0 and numbers::flat_manifold_id-1.

Definition at line 185 of file tria_objects.h.

◆ next_free_single

unsigned int internal::TriangulationImplementation::TriaObjects::next_free_single

Counter for next_free_single_* functions

Definition at line 297 of file tria_objects.h.

◆ next_free_pair

unsigned int internal::TriangulationImplementation::TriaObjects::next_free_pair

Counter for next_free_pair_* functions

Definition at line 302 of file tria_objects.h.

◆ reverse_order_next_free_single

bool internal::TriangulationImplementation::TriaObjects::reverse_order_next_free_single

Bool flag for next_free_single_* functions

Definition at line 307 of file tria_objects.h.

◆ user_data

std::vector<UserData> internal::TriangulationImplementation::TriaObjects::user_data

Pointer which is not used by the library but may be accessed and set by the user to handle data local to a line/quad/etc.

Definition at line 360 of file tria_objects.h.

◆ user_data_type

UserDataType internal::TriangulationImplementation::TriaObjects::user_data_type
mutable

In order to avoid confusion between user pointers and indices, this enum is set by the first function accessing either and subsequent access will not be allowed to change the type of data accessed.

Definition at line 367 of file tria_objects.h.


The documentation for this class was generated from the following files: