16 #ifndef dealii_tria_objects_h 17 #define dealii_tria_objects_h 44 template <
int dim,
int spacedim>
46 template <
class Accessor>
48 template <
int,
int,
int>
54 namespace TriangulationImplementation
153 template <
class Archive>
155 serialize(Archive &ar,
const unsigned int version);
192 const unsigned int new_objs_single = 0);
205 template <
int dim,
int spacedim>
220 template <
int dim,
int spacedim>
228 template <
int dim,
int spacedim>
231 const unsigned int level);
319 template <
class Archive>
321 serialize(Archive &ar,
const unsigned int version);
330 <<
"The containers have sizes " << arg1 <<
" and " << arg2
331 <<
", which is not as expected.");
385 template <
class Archive>
387 serialize(Archive &ar,
const unsigned int version);
506 template <
class Archive>
508 serialize(Archive &ar,
const unsigned int version);
547 const unsigned int new_quads_single = 0);
574 template <
class Archive>
576 serialize(Archive &ar,
const unsigned int version);
582 template <
typename G>
590 template <
typename G>
599 template <
typename G>
600 template <
class Archive>
621 template <
typename G>
624 const unsigned int)
const 630 template <
typename G>
643 template <
typename G>
656 template <
typename G>
657 inline unsigned int &
669 template <
typename G>
678 template <
typename G>
687 template <
typename G>
691 Assert(user_data_type == data_unknown || user_data_type == data_index,
692 ExcPointerIndexClash());
693 user_data_type = data_index;
696 return user_data[i].i;
700 template <
typename G>
704 user_data_type = data_unknown;
705 for (
unsigned int i = 0; i < user_data.size(); ++i)
706 user_data[i].p =
nullptr;
710 template <
typename G>
714 user_flags.assign(user_flags.size(),
false);
718 template <
typename G>
719 template <
class Archive>
729 template <
typename G>
730 template <
class Archive>
735 ar & refinement_cases;
738 ar & boundary_or_material_id;
740 ar &next_free_single &next_free_pair &reverse_order_next_free_single;
741 ar &user_data &user_data_type;
745 template <
class Archive>
747 TriaObjectsHex::serialize(Archive &ar,
const unsigned int version)
751 ar &face_orientations &face_flips &face_rotations;
755 template <
class Archive>
757 TriaObjectsQuad3D::serialize(Archive &ar,
const unsigned int version)
761 ar &line_orientations;
768 TriaObjectsHex::face_orientation(
const unsigned int cell,
769 const unsigned int face)
const 772 face_orientations.size() /
782 TriaObjectsQuad3D::face_orientation(
const unsigned int cell,
783 const unsigned int face)
const 792 template <
int dim,
int spacedim>
800 int pos = next_free_single, last = used.size() - 1;
801 if (!reverse_order_next_free_single)
805 for (; pos < last; ++pos)
815 reverse_order_next_free_single =
true;
816 next_free_single = used.size() - 1;
817 pos = used.size() - 1;
820 next_free_single = pos + 1;
823 if (reverse_order_next_free_single)
827 for (; pos >= 0; --pos)
831 next_free_single = pos - 1;
834 return ::TriaRawIterator<
838 return ::TriaRawIterator<
845 template <
int dim,
int spacedim>
853 int pos = next_free_pair, last = used.size() - 1;
854 for (; pos < last; ++pos)
864 return ::TriaRawIterator<
867 next_free_pair = pos + 2;
869 return ::TriaRawIterator<
static const unsigned int invalid_unsigned_int
std::vector< bool > line_orientations
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcMemoryInexact(int arg1, int arg2)
unsigned int next_free_single
unsigned int next_free_pair
void serialize(Archive &ar, const unsigned int version)
#define AssertIndexRange(index, range)
void monitor_memory(const unsigned int true_dimension) const
std::vector< bool > face_rotations
types::material_id material_id
bool reverse_order_next_free_single
std::vector< int > children
bool face_orientation(const unsigned int cell, const unsigned int face) const
unsigned int & user_index(const unsigned int i)
std::vector< BoundaryOrMaterialId > boundary_or_material_id
std::vector< bool > user_flags
std::vector< UserData > user_data
::TriaRawIterator<::TriaAccessor< G::dimension, dim, spacedim > > next_free_single_object(const Triangulation< dim, spacedim > &tria)
#define Assert(cond, exc)
static ::ExceptionBase & ExcPointerIndexClash()
std::vector< bool > face_orientations
#define DeclException0(Exception0)
#define DEAL_II_NAMESPACE_CLOSE
types::boundary_id boundary_id
void reserve_space(const unsigned int new_objs_in_pairs, const unsigned int new_objs_single=0)
std::vector< RefinementCase< G::dimension > > refinement_cases
Triangulation< dim, spacedim >::raw_hex_iterator next_free_hex(const Triangulation< dim, spacedim > &tria, const unsigned int level)
std::vector< bool > face_flips
UserData contains indices.
#define DEAL_II_NAMESPACE_OPEN
typename IteratorSelector::raw_hex_iterator raw_hex_iterator
std::vector< types::manifold_id > manifold_id
UserData contains pointers.
::TriaRawIterator<::TriaAccessor< G::dimension, dim, spacedim > > next_free_pair_object(const Triangulation< dim, spacedim > &tria)
static std::size_t memory_consumption()
void *& user_pointer(const unsigned int i)
const types::material_id invalid_material_id
UserDataType user_data_type