/usr/include/deal.II/grid/tria_objects.h is in libdeal.ii-dev 6.3.1-1.1.
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 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 | //---------------------------------------------------------------------------
// $Id: tria_objects.h 18092 2009-01-06 04:36:35Z bangerth $
// Version: $Name$
//
// Copyright (C) 2006, 2007, 2008, 2009 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
// to the file deal.II/doc/license.html for the text and
// further information on this license.
//
//---------------------------------------------------------------------------
#ifndef __deal2__tria_objects_h
#define __deal2__tria_objects_h
#include <base/config.h>
#include <base/exceptions.h>
#include <base/geometry_info.h>
#include <grid/tria_object.h>
#include <vector>
DEAL_II_NAMESPACE_OPEN
//TODO: See if we can unify the class hierarchy a bit using partial
// specialization of classes here, e.g. declare a general class
// TriaObjects<dim,G>.
// To consider for this: in that case we would have to duplicate quite a few
// things, e.g. TriaObjectsHex is derived from TriaObjects<Hexahedron> and
// declares mainly additional data. This would have to be changed in case of a
// partial specialization.
//TODO: The TriaObjects class contains a std::vector<G>. This is only an
//efficient storage scheme if G is relatively well packed, i.e. it's not a
//bool and then an integer and then a double, etc. Verify that this is
//actually the case.
template <int dim, int spacedim> class Triangulation;
namespace internal
{
namespace Triangulation
{
/**
* 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.
*
* @ingroup grid
* @author Tobias Leicht, Guido Kanschat, 2006, 2007
*/
template <typename G>
class TriaObjects
{
public:
/**
* Constructor resetting some data.
*/
TriaObjects();
/**
* Vector of the objects belonging to
* this level. The index of the object
* equals the index in this container.
*/
std::vector<G> cells;
/**
* 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.
*/
std::vector<int> children;
/**
* 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.
*/
std::vector<RefinementCase<G::dimension> > refinement_cases;
/**
* Vector storing whether an object is
* used in the @p cells vector.
*
* Since it is difficult to delete
* elements in a @p vector, when an
* element is not needed any more
* (e.g. after derefinement), it is
* not deleted from the list, but
* rather the according @p used flag
* is set to @p false.
*/
std::vector<bool> used;
/**
* 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().
*/
std::vector<bool> user_flags;
/**
* 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 254. 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 254;
* the id 255 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.
*/
std::vector<unsigned char> material_id;
/**
* Assert that enough space
* is allocated to
* accomodate
* <code>new_objs_in_pairs</code>
* new objects, stored in
* pairs, plus
* <code>new_obj_single</code>
* stored individually.
* This function does not
* only call
* <code>vector::reserve()</code>,
* but does really append
* the needed elements.
*
* In 2D e.g. refined lines have to be
* stored in pairs, whereas new lines in the
* interior of refined cells can be stored as
* single lines.
*/
void reserve_space (const unsigned int new_objs_in_pairs,
const unsigned int new_objs_single = 0);
/**
* Return an iterator to the
* next free slot for a
* single line. Only
* implemented for
* <code>G=TriaObject<1>
* </code>.
*/
template <int dim, int spacedim>
typename dealii::Triangulation<dim,spacedim>::raw_line_iterator next_free_single_line (const dealii::Triangulation<dim,spacedim> &tria);
/**
* Return an iterator to the
* next free slot for a pair
* of lines. Only implemented
* for <code>G=TriaObject<1>
* </code>.
*/
template <int dim, int spacedim>
typename dealii::Triangulation<dim,spacedim>::raw_line_iterator next_free_pair_line (const dealii::Triangulation<dim,spacedim> &tria);
/**
* Return an iterator to the
* next free slot for a
* single quad. Only
* implemented for
* <code>G=TriaObject@<2@>
* </code>.
*/
template <int dim, int spacedim>
typename dealii::Triangulation<dim,spacedim>::raw_quad_iterator next_free_single_quad (const dealii::Triangulation<dim,spacedim> &tria);
/**
* Return an iterator to the
* next free slot for a pair
* of quads. Only implemented
* for <code>G=TriaObject@<2@>
* </code>.
*/
template <int dim, int spacedim>
typename dealii::Triangulation<dim,spacedim>::raw_quad_iterator next_free_pair_quad (const dealii::Triangulation<dim,spacedim> &tria);
/**
* Return an iterator to the
* next free slot for a pair
* of hexes. Only implemented
* for
* <code>G=Hexahedron</code>.
*/
template <int dim, int spacedim>
typename dealii::Triangulation<dim,spacedim>::raw_hex_iterator next_free_hex (const dealii::Triangulation<dim,spacedim> &tria,
const unsigned int level);
/**
* Clear all the data contained in this object.
*/
void clear();
/**
* The orientation of the
* face number <code>face</code>
* of the cell with number
* <code>cell</code>. The return
* value is <code>true</code>, if
* the normal vector points
* the usual way
* (GeometryInfo::unit_normal_orientation)
* and <code>false</code> else.
*
* The result is always
* <code>true</code> in this
* class, but derived classes
* will reimplement this.
*
* @warning There is a bug in
* the class hierarchy right
* now. Avoid ever calling
* this function through a
* reference, since you might
* end up with the base class
* function instead of the
* derived class. Still, we
* do not want to make it
* virtual for efficiency
* reasons.
*/
bool face_orientation(const unsigned int cell, const unsigned int face) const;
/**
* Access to user pointers.
*/
void*& user_pointer(const unsigned int i);
/**
* Read-only access to user pointers.
*/
const void* user_pointer(const unsigned int i) const;
/**
* Access to user indices.
*/
unsigned int& user_index(const unsigned int i);
/**
* Read-only access to user pointers.
*/
unsigned int user_index(const unsigned int i) const;
/**
* Reset user data to zero.
*/
void clear_user_data(const unsigned int i);
/**
* Clear all user pointers or
* indices and reset their
* type, such that the next
* access may be aither or.
*/
void clear_user_data();
/**
* Clear all user flags.
*/
void clear_user_flags();
/**
* Check the memory consistency of the
* different containers. Should only be
* called with the prepro flag @p DEBUG
* set. The function should be called from
* the functions of the higher
* TriaLevel classes.
*/
void monitor_memory (const unsigned int true_dimension) const;
/**
* Determine an estimate for the
* memory consumption (in bytes)
* of this object.
*/
unsigned int memory_consumption () const;
/**
* Exception
*/
DeclException3 (ExcMemoryWasted,
char*, int, int,
<< "The container " << arg1 << " contains "
<< arg2 << " elements, but it`s capacity is "
<< arg3 << ".");
/**
* Exception
* @ingroup Exceptions
*/
DeclException2 (ExcMemoryInexact,
int, int,
<< "The containers have sizes " << arg1 << " and "
<< arg2 << ", which is not as expected.");
/**
* Exception
*/
DeclException2 (ExcWrongIterator,
char*, char*,
<< "You asked for the next free " << arg1 << "_iterator, "
"but you can only ask for " << arg2 <<"_iterators.");
/**
* Triangulation objacts can
* either access a user
* pointer or a user
* index. What you tried to
* do is trying to access one
* of those after using the
* other.
*
* @ingroup Exceptions
*/
DeclException0 (ExcPointerIndexClash);
protected:
/**
* Counter for next_free_single_* functions
*/
unsigned int next_free_single;
/**
* Counter for next_free_pair_* functions
*/
unsigned int next_free_pair;
/**
* Bool flag for next_free_single_* functions
*/
bool reverse_order_next_free_single;
/**
* The data type storing user
* pointers or user indices.
*/
union UserData
{
/// The entry used as user pointer.
void* p;
/// The entry used as user index.
unsigned int i;
/// Default constructor
UserData()
{
p = 0;
}
};
/**
* Enum descibing the
* possible types of
* userdata.
*/
enum UserDataType
{
/// No userdata used yet.
data_unknown,
/// UserData contains pointers.
data_pointer,
/// UserData contains indices.
data_index
};
/**
* 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.
*/
std::vector<UserData> user_data;
/**
* 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.
*/
mutable UserDataType user_data_type;
};
/**
* For hexahedrons the data of TriaObjects needs to be extended, as we can obtain faces
* (quads) in non-standard-orientation, therefore we declare a class TriaObjectsHex, which
* additionaly contains a bool-vector of the face-orientations.
* @ingroup grid
*/
class TriaObjectsHex: public TriaObjects<TriaObject<3> >
{
public:
/**
* The orientation of the
* face number <code>face</code>
* of the cell with number
* <code>cell</code>. The return
* value is <code>true</code>, if
* the normal vector points
* the usual way
* (GeometryInfo::unit_normal_orientation)
* and <code>false</code> if they
* point in opposite
* direction.
*/
bool face_orientation(const unsigned int cell, const unsigned int face) const;
/**
* For edges, we enforce a
* standard convention that
* opposite edges should be
* parallel. Now, that's
* enforcable in most cases,
* and we have code that
* makes sure that if a mesh
* allows this to happen,
* that we have this
* convention. We also know
* that it is always possible
* to have opposite faces
* have parallel normal
* vectors. (For both things,
* see the Agelek, Anderson,
* Bangerth, Barth paper
* mentioned in the
* publications list.)
*
* The problem is that we
* originally had another
* condition, namely that
* faces 0, 2 and 6 have
* normals that point into
* the cell, while the other
* faces have normals that
* point outward. It turns
* out that this is not
* always possible. In
* effect, we have to store
* whether the normal vector
* of each face of each cell
* follows this convention or
* not. If this is so, then
* this variable stores a
* @p true value, otherwise
* a @p false value.
*
* In effect, this field has
* <code>6*n_cells</code> elements,
* being the number of cells
* times the six faces each
* has.
*/
std::vector<bool> face_orientations;
/**
* flip = rotation by 180 degrees
*/
std::vector<bool> face_flips;
/**
* rotation by 90 degrees
*/
std::vector<bool> face_rotations;
/**
* Assert that enough space is
* allocated to accomodate
* <code>new_objs</code> new objects.
* This function does not only call
* <code>vector::reserve()</code>, but
* does really append the needed
* elements.
*/
void reserve_space (const unsigned int new_objs);
/**
* Clear all the data contained in this object.
*/
void clear();
/**
* Check the memory consistency of the
* different containers. Should only be
* called with the prepro flag @p DEBUG
* set. The function should be called from
* the functions of the higher
* TriaLevel classes.
*/
void monitor_memory (const unsigned int true_dimension) const;
/**
* Determine an estimate for the
* memory consumption (in bytes)
* of this object.
*/
unsigned int memory_consumption () const;
};
/**
* For quadrilaterals in 3D the data of TriaObjects needs to be extended, as we
* can obtain faces (quads) with lines in non-standard-orientation, therefore we
* declare a class TriaObjectsQuad3D, which additionaly contains a bool-vector
* of the line-orientations.
* @ingroup grid
*/
class TriaObjectsQuad3D: public TriaObjects<TriaObject<2> >
{
public:
/**
* The orientation of the
* face number <code>face</code>
* of the cell with number
* <code>cell</code>. The return
* value is <code>true</code>, if
* the normal vector points
* the usual way
* (GeometryInfo::unit_normal_orientation)
* and <code>false</code> if they
* point in opposite
* direction.
*/
bool face_orientation(const unsigned int cell, const unsigned int face) const;
/**
* In effect, this field has
* <code>4*n_quads</code> elements,
* being the number of quads
* times the four lines each
* has.
*/
std::vector<bool> line_orientations;
/**
* Assert that enough space
* is allocated to
* accomodate
* <code>new_quads_in_pairs</code>
* new quads, stored in
* pairs, plus
* <code>new_quads_single</code>
* stored individually.
* This function does not
* only call
* <code>vector::reserve()</code>,
* but does really append
* the needed elements.
*/
void reserve_space (const unsigned int new_quads_in_pairs,
const unsigned int new_quads_single = 0);
/**
* Clear all the data contained in this object.
*/
void clear();
/**
* Check the memory consistency of the
* different containers. Should only be
* called with the prepro flag @p DEBUG
* set. The function should be called from
* the functions of the higher
* TriaLevel classes.
*/
void monitor_memory (const unsigned int true_dimension) const;
/**
* Determine an estimate for the
* memory consumption (in bytes)
* of this object.
*/
unsigned int memory_consumption () const;
};
//----------------------------------------------------------------------//
template<typename G>
inline
bool
TriaObjects<G>::face_orientation(const unsigned int, const unsigned int) const
{
return true;
}
template<typename G>
inline
void*&
TriaObjects<G>::user_pointer (const unsigned int i)
{
#ifdef DEBUG
Assert(user_data_type == data_unknown || user_data_type == data_pointer,
ExcPointerIndexClash());
user_data_type = data_pointer;
#endif
Assert(i<user_data.size(), ExcIndexRange(i,0,user_data.size()));
return user_data[i].p;
}
template<typename G>
inline
const void*
TriaObjects<G>::user_pointer (const unsigned int i) const
{
#ifdef DEBUG
Assert(user_data_type == data_unknown || user_data_type == data_pointer,
ExcPointerIndexClash());
user_data_type = data_pointer;
#endif
Assert(i<user_data.size(), ExcIndexRange(i,0,user_data.size()));
return user_data[i].p;
}
template<typename G>
inline
unsigned int&
TriaObjects<G>::user_index (const unsigned int i)
{
#ifdef DEBUG
Assert(user_data_type == data_unknown || user_data_type == data_index,
ExcPointerIndexClash());
user_data_type = data_index;
#endif
Assert(i<user_data.size(), ExcIndexRange(i,0,user_data.size()));
return user_data[i].i;
}
template<typename G>
inline
void
TriaObjects<G>::clear_user_data (const unsigned int i)
{
Assert(i<user_data.size(), ExcIndexRange(i,0,user_data.size()));
user_data[i].i = 0;
}
template <typename G>
inline
TriaObjects<G>::TriaObjects()
:
user_data_type(data_unknown)
{}
template<typename G>
inline
unsigned int TriaObjects<G>::user_index (const unsigned int i) const
{
#ifdef DEBUG
Assert(user_data_type == data_unknown || user_data_type == data_index,
ExcPointerIndexClash());
user_data_type = data_index;
#endif
Assert(i<user_data.size(), ExcIndexRange(i,0,user_data.size()));
return user_data[i].i;
}
template<typename G>
inline
void TriaObjects<G>::clear_user_data ()
{
user_data_type = data_unknown;
for (unsigned int i=0;i<user_data.size();++i)
user_data[i].p = 0;
}
template<typename G>
inline
void TriaObjects<G>::clear_user_flags ()
{
user_flags.assign(user_flags.size(),false);
}
//----------------------------------------------------------------------//
inline
bool
TriaObjectsHex::face_orientation(const unsigned int cell, const unsigned int face) const
{
Assert (cell < face_orientations.size() / GeometryInfo<3>::faces_per_cell,
ExcIndexRange(0, cell, face_orientations.size() / GeometryInfo<3>::faces_per_cell));
Assert (face < GeometryInfo<3>::faces_per_cell,
ExcIndexRange(0, face, GeometryInfo<3>::faces_per_cell));
return face_orientations[cell * GeometryInfo<3>::faces_per_cell
+ face];
}
//----------------------------------------------------------------------//
inline
bool
TriaObjectsQuad3D::face_orientation(const unsigned int cell, const unsigned int face) const
{
return line_orientations[cell * GeometryInfo<2>::faces_per_cell
+ face];
}
// declaration of explicit specializations
template<>
void
TriaObjects<TriaObject<1> >::reserve_space (const unsigned int new_lines_in_pairs,
const unsigned int new_lines_single);
template<>
void
TriaObjects<TriaObject<2> >::reserve_space (const unsigned int new_quads_in_pairs,
const unsigned int new_quads_single);
template<>
void
TriaObjects<TriaObject<2> >::monitor_memory (const unsigned int) const;
}
}
DEAL_II_NAMESPACE_CLOSE
#endif
|