![]() |
Reference documentation for deal.II version 8.1.0
|
#include <tria_accessor.h>
Public Types | |
typedef TriaAccessor< dim, dim, spacedim >::AccessorData | AccessorData |
typedef Triangulation< dim, spacedim > | Container |
![]() | |
typedef TriaAccessorBase < structdim, dim, spacedim > ::AccessorData | AccessorData |
![]() | |
typedef void * | LocalData |
Public Member Functions | |
DeclException0 (ExcRefineCellNotActive) | |
DeclException0 (ExcCellFlaggedForRefinement) | |
DeclException0 (ExcCellFlaggedForCoarsening) | |
template<> | |
internal::SubfaceCase< 1 > | subface_case (const unsigned int) const |
template<> | |
internal::SubfaceCase< 1 > | subface_case (const unsigned int) const |
template<> | |
internal::SubfaceCase< 1 > | subface_case (const unsigned int) const |
template<> | |
internal::SubfaceCase< 2 > | subface_case (const unsigned int face_no) const |
template<> | |
internal::SubfaceCase< 2 > | subface_case (const unsigned int face_no) const |
template<> | |
internal::SubfaceCase< 3 > | subface_case (const unsigned int face_no) const |
Constructors | |
CellAccessor (const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *local_data=0) | |
CellAccessor (const TriaAccessor< dim, dim, spacedim > &cell_accessor) | |
template<int structdim2, int dim2, int spacedim2> | |
CellAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
template<int structdim2, int dim2, int spacedim2> | |
CellAccessor (const TriaAccessor< structdim2, dim2, spacedim2 > &) | |
Accessing sub-objects and neighbors | |
TriaIterator< CellAccessor < dim, spacedim > > | child (const unsigned int i) const |
TriaIterator< TriaAccessor < dim-1, dim, spacedim > > | face (const unsigned int i) const |
unsigned int | face_index (const unsigned int i) const |
TriaIterator< CellAccessor < dim, spacedim > > | neighbor_child_on_subface (const unsigned int face_no, const unsigned int subface_no) const |
TriaIterator< CellAccessor < dim, spacedim > > | neighbor (const unsigned int i) const |
int | neighbor_index (const unsigned int i) const |
int | neighbor_level (const unsigned int i) const |
unsigned int | neighbor_of_neighbor (const unsigned int neighbor) const |
bool | neighbor_is_coarser (const unsigned int neighbor) const |
std::pair< unsigned int, unsigned int > | neighbor_of_coarser_neighbor (const unsigned int neighbor) const |
unsigned int | neighbor_face_no (const unsigned int neighbor) const |
Dealing with boundary indicators | |
bool | at_boundary (const unsigned int i) const |
bool | at_boundary () const |
bool | has_boundary_lines () const |
Dealing with refinement indicators | |
RefinementCase< dim > | refine_flag_set () const |
void | set_refine_flag (const RefinementCase< dim > ref_case=RefinementCase< dim >::isotropic_refinement) const |
void | clear_refine_flag () const |
bool | flag_for_face_refinement (const unsigned int face_no, const RefinementCase< dim-1 > &face_refinement_case=RefinementCase< dim-1 >::isotropic_refinement) const |
bool | flag_for_line_refinement (const unsigned int line_no) const |
::internal::SubfaceCase< dim > | subface_case (const unsigned int face_no) const |
bool | coarsen_flag_set () const |
void | set_coarsen_flag () const |
void | clear_coarsen_flag () const |
Dealing with material indicators | |
types::material_id | material_id () const |
void | set_material_id (const types::material_id new_material_id) const |
void | recursively_set_material_id (const types::material_id new_material_id) const |
Dealing with subdomain indicators | |
types::subdomain_id | subdomain_id () const |
void | set_subdomain_id (const types::subdomain_id new_subdomain_id) const |
types::subdomain_id | level_subdomain_id () const |
void | set_level_subdomain_id (const types::subdomain_id new_level_subdomain_id) const |
void | recursively_set_subdomain_id (const types::subdomain_id new_subdomain_id) const |
Dealing with codim 1 cell orientation | |
bool | direction_flag () const |
TriaIterator< CellAccessor < dim, spacedim > > | parent () const |
Other functions | |
bool | active () const |
bool | operator< (const CellAccessor< dim, spacedim > &other) const |
bool | is_locally_owned () const |
bool | is_locally_owned_on_level () const |
bool | is_ghost () const |
bool | is_artificial () const |
bool | point_inside (const Point< spacedim > &p) const |
void | set_neighbor (const unsigned int i, const TriaIterator< CellAccessor< dim, spacedim > > &pointer) const |
CellId | id () const |
![]() | |
TriaAccessor (const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *local_data=0) | |
TriaAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &) | |
TriaAccessor (const TriaAccessor< structdim2, dim2, spacedim2 > &) | |
bool | used () const |
int | parent_index () const |
unsigned int | vertex_index (const unsigned int i) const |
Point< spacedim > & | vertex (const unsigned int i) const |
typename::internal::Triangulation::Iterators < dim, spacedim > ::line_iterator | line (const unsigned int i) const |
unsigned int | line_index (const unsigned int i) const |
typename::internal::Triangulation::Iterators < dim, spacedim > ::quad_iterator | quad (const unsigned int i) const |
unsigned int | quad_index (const unsigned int i) const |
bool | face_orientation (const unsigned int face) const |
bool | face_flip (const unsigned int face) const |
bool | face_rotation (const unsigned int face) const |
bool | line_orientation (const unsigned int line) const |
bool | has_children () const |
unsigned int | n_children () const |
unsigned int | number_of_children () const |
unsigned int | max_refinement_depth () const |
TriaIterator< TriaAccessor < structdim, dim, spacedim > > | child (const unsigned int i) const |
TriaIterator< TriaAccessor < structdim, dim, spacedim > > | isotropic_child (const unsigned int i) const |
RefinementCase< structdim > | refinement_case () const |
int | child_index (const unsigned int i) const |
int | isotropic_child_index (const unsigned int i) const |
types::boundary_id | boundary_indicator () const |
void | set_boundary_indicator (const types::boundary_id) const |
void | set_all_boundary_indicators (const types::boundary_id) const |
bool | at_boundary () const |
const Boundary< dim, spacedim > & | get_boundary () const |
bool | user_flag_set () const |
void | set_user_flag () const |
void | clear_user_flag () const |
void | recursively_set_user_flag () const |
void | recursively_clear_user_flag () const |
void | clear_user_data () const |
void | set_user_pointer (void *p) const |
void | clear_user_pointer () const |
void * | user_pointer () const |
void | recursively_set_user_pointer (void *p) const |
void | recursively_clear_user_pointer () const |
void | set_user_index (const unsigned int p) const |
void | clear_user_index () const |
unsigned int | user_index () const |
void | recursively_set_user_index (const unsigned int p) const |
void | recursively_clear_user_index () const |
double | diameter () const |
double | extent_in_direction (const unsigned int axis) const |
double | minimum_vertex_distance () const |
Point< spacedim > | center () const |
Point< spacedim > | barycenter () const |
double | measure () const |
bool | is_translation_of (const TriaIterator< TriaAccessor< structdim, dim, spacedim > > &o) const |
![]() | |
int | level () const |
int | index () const |
IteratorState::IteratorStates | state () const |
const Triangulation< dim, spacedim > & | get_triangulation () const |
Protected Member Functions | |
unsigned int | neighbor_of_neighbor_internal (const unsigned int neighbor) const |
template<int dim_, int spacedim_> | |
bool | point_inside_codim (const Point< spacedim_ > &p) const |
![]() | |
TriaAccessorBase (const Triangulation< dim, spacedim > *parent=0, const int level=-1, const int index=-1, const AccessorData *=0) | |
TriaAccessorBase (const TriaAccessorBase &) | |
void | copy_from (const TriaAccessorBase &) |
TriaAccessorBase & | operator= (const TriaAccessorBase &) |
bool | operator< (const TriaAccessorBase &other) const |
void | operator= (const TriaAccessorBase *) |
bool | operator== (const TriaAccessorBase &) const |
bool | operator!= (const TriaAccessorBase &) const |
::internal::Triangulation::TriaObjects <::internal::Triangulation::TriaObject < structdim > > & | objects () const |
void | operator++ () |
void | operator-- () |
Private Member Functions | |
void | set_direction_flag (const bool new_direction_flag) const |
void | operator= (const CellAccessor< dim, spacedim > &) |
Friends | |
template<int , int > | |
class | Triangulation |
struct | ::internal::Triangulation::Implementation |
Additional Inherited Members | |
![]() | |
static const unsigned int | space_dimension = spacedim |
static const unsigned int | dimension = dim |
static const unsigned int | structure_dimension = structdim |
![]() | |
typedef void | AccessorData |
![]() | |
typename::internal::TriaAccessor::PresentLevelType < structdim, dim >::type | present_level |
int | present_index |
const Triangulation< dim, spacedim > * | tria |
This class allows access to a cell: a line in one dimension, a quad in two dimension, etc.
The following refers to any dimension:
This class allows access to a cell
, which is a line in 1D and a quad in 2D. Cells have more functionality than lines or quads by themselves, for example they can be flagged for refinement, they have neighbors, they have the possibility to check whether they are at the boundary etc. This class offers access to all this data.
Definition at line 2222 of file tria_accessor.h.
typedef TriaAccessor<dim,dim,spacedim>::AccessorData CellAccessor< dim, spacedim >::AccessorData |
Propagate the AccessorData type into the present class.
Definition at line 2229 of file tria_accessor.h.
typedef Triangulation<dim, spacedim> CellAccessor< dim, spacedim >::Container |
Define the type of the container this is part of.
Definition at line 2235 of file tria_accessor.h.
|
inline |
Constructor.
Definition at line 2469 of file tria_accessor.templates.h.
|
inline |
Copy constructor.
Definition at line 2481 of file tria_accessor.templates.h.
CellAccessor< dim, spacedim >::CellAccessor | ( | const InvalidAccessor< structdim2, dim2, spacedim2 > & | ) |
Conversion constructor. This constructor exists to make certain constructs simpler to write in dimension independent code. For example, it allows assigning a face iterator to a line iterator, an operation that is useful in 2d but doesn't make any sense in 3d. The constructor here exists for the purpose of making the code conform to C++ but it will unconditionally abort; in other words, assigning a face iterator to a line iterator is better put into an if-statement that checks that the dimension is two, and assign to a quad iterator in 3d (an operator that, without this constructor would be illegal if we happen to compile for 2d).
Definition at line 3169 of file tria_accessor.h.
CellAccessor< dim, spacedim >::CellAccessor | ( | const TriaAccessor< structdim2, dim2, spacedim2 > & | ) |
Another conversion operator between objects that don't make sense, just like the previous one.
Definition at line 3199 of file tria_accessor.h.
|
inline |
Return a pointer to the ith
child. Overloaded version which returns a more reasonable iterator class.
Definition at line 2887 of file tria_accessor.templates.h.
|
inline |
Return an iterator to the ith
face of this cell.
Definition at line 2540 of file tria_accessor.templates.h.
|
inline |
Return the (global) index of the ith
face of this cell.
Definition at line 2550 of file tria_accessor.templates.h.
TriaIterator<CellAccessor<dim, spacedim> > CellAccessor< dim, spacedim >::neighbor_child_on_subface | ( | const unsigned int | face_no, |
const unsigned int | subface_no | ||
) | const |
Return an iterator to that cell that neighbors the present cell on the given face and subface number.
To succeed, the present cell must not be further refined, and the neighbor on the given face must be further refined exactly once; the returned cell is then a child of that neighbor.
The function may not be called in 1d, since there we have no subfaces. The implementation of this function is rather straightforward in 2d, by first determining which face of the neighbor cell the present cell is bordering on (this is what the neighbor_of_neighbor
function does), and then asking GeometryInfo::child_cell_on_subface
for the index of the child.
However, the situation is more complicated in 3d, since there faces may have more than one orientation, and we have to use face_orientation
, face_flip
and face_rotation
for both this and the neighbor cell to figure out which cell we want to have.
This can lead to surprising results: if we are sitting on a cell and are asking for a cell behind subface sf
, then this means that we are considering the subface for the face in the natural direction for the present cell. However, if the face as seen from this cell has face_orientation()==false
, then the child of the face that separates the present cell from the neighboring cell's child is not necessarily the sf-th
child of the face of this cell. This is so because the subface_no
on a cell corresponds to the subface with respect to the intrinsic ordering of the present cell, whereas children of face iterators are computed with respect to the intrinsic ordering of faces; these two orderings are only identical if the face orientation is true
, and reversed otherwise.
Similarly, effects of face_flip()==true
and face_rotation()==true()
, both of which indicate a non-standard face have to be considered.
Fortunately, this is only very rarely of concern, since usually one simply wishes to loop over all finer neighbors at a given face of an active cell. Only in the process of refinement of a Triangulation we want to set neighbor information for both our child cells and the neighbor's children. Since we can respect orientation of faces from our current cell in that case, we do NOT respect face_orientation, face_flip and face_rotation of the present cell within this function, i.e. the returned neighbor's child is behind subface subface
concerning the intrinsic ordering of the given face.
|
inline |
Return a pointer to the ith
neighbor. If the neighbor does not exist, an invalid iterator is returned.
Note (cf. TriaLevel<0>): The neighbor of a cell has at most the same level as this cell, i.e. it may or may not be refined.
Definition at line 2871 of file tria_accessor.templates.h.
|
inline |
Return the index of the ith
neighbor. If the neighbor does not exist, its index is -1.
Definition at line 2575 of file tria_accessor.templates.h.
|
inline |
Return the level of the ith
neighbor. If the neighbor does not exist, its level is -1.
Definition at line 2587 of file tria_accessor.templates.h.
unsigned int CellAccessor< dim, spacedim >::neighbor_of_neighbor | ( | const unsigned int | neighbor | ) | const |
Return the how-many'th neighbor this cell is of cell->neighbor(neighbor)
, i.e. return the face_no
such that cell->neighbor(neighbor)->neighbor(face_no)==cell
. This function is the right one if you want to know how to get back from a neighbor to the present cell.
Note that this operation is only useful if the neighbor is not coarser than the present cell. If the neighbor is coarser this function throws an exception. Use the neighbor_of_coarser_neighbor
function in that case.
bool CellAccessor< dim, spacedim >::neighbor_is_coarser | ( | const unsigned int | neighbor | ) | const |
Return, whether the neighbor is coarser then the present cell. This is important in case of ansiotropic refinement where this information does not depend on the levels of the cells.
Note, that in an anisotropic setting, a cell can only be coarser than another one at a given face, not on a general basis. The face of the finer cell is contained in the corresponding face of the coarser cell, the finer face is either a child or a grandchild of the coarser face.
std::pair<unsigned int, unsigned int> CellAccessor< dim, spacedim >::neighbor_of_coarser_neighbor | ( | const unsigned int | neighbor | ) | const |
This function is a generalization of the neighbor_of_neighbor
function for the case of a coarser neighbor. It returns a pair of numbers, face_no and subface_no, with the following property, if the neighbor is not refined: cell->neighbor(neighbor)->neighbor_child_on_subface(face_no,subface_no)==cell
. In 3D, a coarser neighbor can still be refined. In that case subface_no denotes the child index of the neighbors face that relates to our face: cell->neighbor(neighbor)->face(face_no)->child(subface_no)==cell->face(neighbor)
. This case in 3d and how it can happen is discussed in the introduction of the step-30 tutorial program.
This function is impossible for dim==1
.
|
inline |
This function is a generalization of the neighbor_of_neighbor
and the neighbor_of_coarser_neighbor
functions. It checks whether the neighbor is coarser or not and calls the respective function. In both cases, only the face_no is returned.
Definition at line 2994 of file tria_accessor.templates.h.
bool CellAccessor< dim, spacedim >::at_boundary | ( | const unsigned int | i | ) | const |
Return whether the ith
vertex or face (depending on the dimension) is part of the boundary. This is true, if the ith
neighbor does not exist.
bool CellAccessor< dim, spacedim >::at_boundary | ( | ) | const |
Return whether the cell is at the boundary. Being at the boundary is defined by one face being on the boundary. Note that this does not catch cases where only one vertex of a quad or of a hex is at the boundary, or where only one line of a hex is at the boundary while the interiors of all faces are in the interior of the domain. For the latter case, the has_boundary_lines
function is the right one to ask.
bool CellAccessor< dim, spacedim >::has_boundary_lines | ( | ) | const |
This is a slight variation to the at_boundary
function: for 1 and 2 dimensions, it is equivalent, for three dimensions it returns whether at least one of the 12 lines of the hexahedron is at a boundary. This, of course, includes the case where a whole face is at the boundary, but also some other cases.
|
inline |
Return the RefinementCase<dim>
this cell was flagged to be refined with.
Definition at line 2599 of file tria_accessor.templates.h.
|
inline |
Flag the cell pointed to for refinement. This function is only allowed for active cells.
Definition at line 2617 of file tria_accessor.templates.h.
|
inline |
Clear the refinement flag.
Definition at line 2631 of file tria_accessor.templates.h.
|
inline |
Modify the refinement flag of the cell to ensure (at least) the given refinement case face_refinement_case
at face face_no
, taking into account orientation, flip and rotation of the face. Return, whether the refinement flag had to be modified. This function is only allowed for active cells.
Definition at line 2643 of file tria_accessor.templates.h.
|
inline |
Modify the refinement flag of the cell to ensure that line face_no
will be refined. Return, whether the refinement flag had to be modified. This function is only allowed for active cells.
Definition at line 2675 of file tria_accessor.templates.h.
::internal::SubfaceCase<dim> CellAccessor< dim, spacedim >::subface_case | ( | const unsigned int | face_no | ) | const |
Return the SubfaceCase of face face_no
. Note that this is not identical to asking cell->face(face_no)->refinement_case()
since the latter returns a RefinementCase<dim-1> and thus only considers one (anisotropic) refinement, whereas this function considers the complete refinement situation including possible refinement of the face's children. This function may only be called for active cells in 2d and 3d.
|
inline |
Return whether the coarsen flag is set or not.
Definition at line 2829 of file tria_accessor.templates.h.
|
inline |
Flag the cell pointed to for coarsening. This function is only allowed for active cells.
Definition at line 2847 of file tria_accessor.templates.h.
|
inline |
Clear the coarsen flag.
Definition at line 2860 of file tria_accessor.templates.h.
types::material_id CellAccessor< dim, spacedim >::material_id | ( | ) | const |
void CellAccessor< dim, spacedim >::set_material_id | ( | const types::material_id | new_material_id | ) | const |
void CellAccessor< dim, spacedim >::recursively_set_material_id | ( | const types::material_id | new_material_id | ) | const |
Set the material id of this cell and all its children (and grand-children, and so on) to the given value.
See the glossary for more information.
types::subdomain_id CellAccessor< dim, spacedim >::subdomain_id | ( | ) | const |
Return the subdomain id of this cell.
See the glossary for more information.
void CellAccessor< dim, spacedim >::set_subdomain_id | ( | const types::subdomain_id | new_subdomain_id | ) | const |
Set the subdomain id of this cell.
See the glossary for more information. This function should not be called if you use a parallel::distributed::Triangulation object.
types::subdomain_id CellAccessor< dim, spacedim >::level_subdomain_id | ( | ) | const |
Get the level subdomain id of this cell. This is used for parallel multigrid.
void CellAccessor< dim, spacedim >::set_level_subdomain_id | ( | const types::subdomain_id | new_level_subdomain_id | ) | const |
Set the level subdomain id of this cell. This is used for parallel multigrid.
void CellAccessor< dim, spacedim >::recursively_set_subdomain_id | ( | const types::subdomain_id | new_subdomain_id | ) | const |
Set the subdomain id of this cell (if it is active) or all its terminal children (and grand-children, and so on, as long as they have no children of their own) to the given value. Since the subdomain id is a concept that is only defined for cells that are active (i.e., have no children of their own), this function only sets the subdomain ids for all children and grand children of this cell that are actually active, skipping intermediate child cells.
See the glossary for more information. This function should not be called if you use a parallel::distributed::Triangulation object since there the subdomain id is implicitly defined by which processor you're on.
bool CellAccessor< dim, spacedim >::direction_flag | ( | ) | const |
Return the orientation of this cell.
For the meaning of this flag, see GlossDirectionFlag .
TriaIterator<CellAccessor<dim,spacedim> > CellAccessor< dim, spacedim >::parent | ( | ) | const |
Return an iterator to the parent. Throws an exception if this cell has no parent, i.e. has level 0.
|
inline |
Test whether the cell has children (this is the criterion for activity of a cell).
See the glossary for more information.
Definition at line 2903 of file tria_accessor.templates.h.
|
inline |
Ordering of accessors. This function implements a total ordering of cells even on a parallel::distributed::Triangulation. This function first compares level_subdomain_id(). If these are equal, and both cells are active, it compares subdomain_id(). If this is inconclusive, TriaAccessorBase::operator < () is called.
Definition at line 3010 of file tria_accessor.templates.h.
|
inline |
Return whether this cell is owned by the current processor or is owned by another processor. The function always returns true if applied to an object of type Triangulation, but may yield false if the triangulation is of type parallel::distributed::Triangulation.
See the glossary and the Parallel computing with multiple processors using distributed memory module for more information.
!is_ghost() && !is_artificial()
.Definition at line 2913 of file tria_accessor.templates.h.
|
inline |
Return true if either the Triangulation is not distributed or if level_subdomain_id() is equal to the id of the current processor.
Definition at line 2936 of file tria_accessor.templates.h.
|
inline |
Return whether this cell exists in the global mesh but (i) is owned by another processor, i.e. has a subdomain_id different from the one the current processor owns and (ii) is adjacent to a cell owned by the current processor.
This function only makes sense if the triangulation used is of kind parallel::distributed::Triangulation. In all other cases, the returned value is always false.
See the glossary and the Parallel computing with multiple processors using distributed memory module for more information.
!is_locally_owned() && !is_artificial()
.Definition at line 2955 of file tria_accessor.templates.h.
|
inline |
Return whether this cell is artificial, i.e. it isn't one of the cells owned by the current processor, and it also doesn't border on one. As a consequence, it exists in the mesh to ensure that each processor has all coarse mesh cells and that the 2:1 ratio of neighboring cells is maintained, but it is not one of the cells we should work on on the current processor. In particular, there is no guarantee that this cell isn't, in fact, further refined on one of the other processors.
This function only makes sense if the triangulation used is of kind parallel::distributed::Triangulation. In all other cases, the returned value is always false.
See the glossary and the Parallel computing with multiple processors using distributed memory module for more information.
!is_ghost() && !is_artificial()
.Definition at line 2979 of file tria_accessor.templates.h.
bool CellAccessor< dim, spacedim >::point_inside | ( | const Point< spacedim > & | p | ) | const |
Test whether the point p
is inside this cell. Points on the boundary are counted as being inside the cell.
Note that this function assumes that the mapping between unit cell and real cell is (bi-, tri-)linear, i.e. that faces in 2d and edges in 3d are straight lines. If you have higher order transformations, results may be different as to whether a point is in- or outside the cell in real space.
In case of codim>0, the point is first projected to the manifold where the cell is embedded and then check if this projection is inside the cell.
void CellAccessor< dim, spacedim >::set_neighbor | ( | const unsigned int | i, |
const TriaIterator< CellAccessor< dim, spacedim > > & | pointer | ||
) | const |
Set the neighbor i
of this cell to the cell pointed to by pointer
.
This function shouldn't really be public (but needs to for various reasons in order not to make a long list of functions friends): it modifies internal data structures and may leave things. Do not use it from application codes.
CellId CellAccessor< dim, spacedim >::id | ( | ) | const |
Return a unique ID for the current cell. This ID is constructed from the path in the hierarchy from the coarse father cell and works correctly in parallel computations.
Note: This operation takes O(log(level)) time.
Definition at line 3211 of file tria_accessor.h.
|
protected |
This function assumes that the neighbor is not coarser than the current cell. In this case it returns the neighbor_of_neighbor() value. If, however, the neighbor is coarser this function returns an invalid_unsigned_int
.
This function is not for public use. Use the function neighbor_of_neighbor() instead which throws an exception if called for a coarser neighbor. If neighbor is indeed coarser (you get to know this by e.g. the neighbor_is_coarser() function) then the neighbor_of_coarser_neighbor() function should be call. If you'd like to know only the face_no
which is required to get back from the neighbor to the present cell then simply use the neighbor_face_no() function which can be used for coarser as well as noncoarser neighbors.
|
protected |
As for any codim>0 we can use a similar code and c++ does not allow partial templates. we use this auxiliary function that is then called from point_inside.
|
private |
Set the orientation of this cell.
For the meaning of this flag, see GlossDirectionFlag .
|
private |
Copy operator. This is normally used in a context like iterator a,b; *a=*b;
. Since the meaning is to copy the object pointed to by b
to the object pointed to by a
and since accessors are not real but virtual objects, this operation is not useful for iterators on triangulations. We declare this function here private, thus it may not be used from outside. Furthermore it is not implemented and will give a linker error if used anyway.