/usr/include/deal.II/fe/fe_base.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 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 | //---------------------------------------------------------------------------
// $Id: fe_base.h 18744 2009-04-26 21:37:54Z bangerth $
// Version: $Name$
//
// Copyright (C) 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 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__fe_base_h
#define __deal2__fe_base_h
#include <base/config.h>
#include <base/exceptions.h>
#include <base/subscriptor.h>
#include <base/point.h>
#include <base/tensor.h>
#include <base/table.h>
#include <base/vector_slice.h>
#include <base/geometry_info.h>
#include <lac/full_matrix.h>
#include <fe/fe_update_flags.h>
#include <fe/mapping.h>
#include <string>
#include <vector>
DEAL_II_NAMESPACE_OPEN
template<int dim, int spacedim> class FESystem;
/**
* A namespace solely for the purpose of defining the Domination enum as well
* as associated operators.
*/
namespace FiniteElementDomination
{
/**
* An enum that describes the
* outcome of comparing two elements for
* mutual domination. If one element
* dominates another, then the
* restriction of the space described by
* the dominated element to a face of the
* cell is strictly larger than that of
* the dominating element. For example,
* in 2-d Q(2) elements dominate Q(4)
* elements, because the traces of Q(4)
* elements are quartic polynomials which
* is a space strictly larger than the
* quadratic polynomials (the restriction
* of the Q(2) element). In general, Q(k)
* dominates Q(k') if $k\le k'$.
*
* This enum is used in the
* FiniteElement::compare_for_face_domination()
* function that is used in the context
* of hp finite element methods when
* determining what to do at faces where
* two different finite elements meet
* (see the @ref hp_paper "hp paper" for a more detailed
* description of the following). In that
* case, the degrees of freedom of one
* side need to be constrained to those
* on the other side. The determination
* which side is which is based on the
* outcome of a comparison for mutual
* domination: the dominated side is
* constrained to the dominating one.
*
* A similar situation happens in 3d, where
* we have to consider different elements
* meeting at only an edge, not an entire
* face. Such comparisons are then
* implemented in the
* FiniteElement::compare_for_line_domination()
* function.
*
* Note that there are situations where
* neither side dominates. The @ref hp_paper "hp paper"
* lists two case, with the simpler one
* being that a $Q_2\times Q_1$
* vector-valued element (i.e. a
* <code>FESystem(FE_Q(2),1,FE_Q(1),1)</code>)
* meets a $Q_1\times Q_2$ element: here,
* for each of the two vector-components,
* we can define a domination
* relationship, but it is different for
* the two components.
*
* It is clear that the concept of
* domination doesn't matter for
* discontinuous elements. However,
* discontinuous elements may be part of
* vector-valued elements and may
* therefore be compared against each
* other for domination. They should
* return
* <code>either_element_can_dominate</code>
* in that case. Likewise, when comparing
* two identical finite elements, they
* should return this code; the reason is
* that we can not decide which element
* will dominate at the time we look at
* the first component of, for example,
* two $Q_2\times Q_1$ and $Q_2\times
* Q_2$ elements, and have to keep our
* options open until we get to the
* second base element.
*
* More details on domination can be found
* in the @ref hp_paper "hp paper".
*/
enum Domination
{
this_element_dominates,
other_element_dominates,
neither_element_dominates,
either_element_can_dominate
};
/**
* A generalization of the binary
* <code>and</code> operator to a comparison
* relationship. The way this works is
* pretty much as when you would want to
* define a comparison relationship for
* vectors: either all elements of the
* first vector are smaller, equal, or
* larger than those of the second vector,
* or some are and some are not.
*
* This operator is pretty much the same:
* if both arguments are
* <code>this_element_dominates</code> or
* <code>other_element_dominates</code>,
* then the returned value is that
* value. On the other hand, if one of the
* values is
* <code>either_element_can_dominate</code>,
* then the returned value is that of the
* other argument. If either argument is
* <code>neither_element_dominates</code>,
* or if the two arguments are
* <code>this_element_dominates</code> and
* <code>other_element_dominates</code>,
* then the returned value is
* <code>neither_element_dominates</code>.
*/
inline Domination operator & (const Domination d1,
const Domination d2);
}
/**
* Dimension independent data for finite elements. See the derived
* class FiniteElement class for information on its use. All
* its data are available to the implementation in a concrete finite
* element class.
*
* @ingroup febase
* @author Wolfgang Bangerth, Guido Kanschat, 1998, 1999, 2000, 2001, 2003, 2005
*/
template <int dim>
class FiniteElementData
{
public:
/**
* Enumerator for the different
* types of continuity a finite
* element may have. Continuity
* is measured by the Sobolev
* space containing the
* constructed finite element
* space and is also called this
* way.
*
* Note that certain continuities
* may imply others. For
* instance, a function in
* <i>H<sup>1</sup></i> is in
* <i>H<sup>curl</sup></i> and
* <i>H<sup>div</sup></i> as
* well.
*
* If you are interested in
* continuity in the classical
* sense, then the following
* relations hold:
*
* <ol>
*
* <li> <i>H<sup>1</sup></i>
* implies that the function is
* continuous over cell
* boundaries.
*
* <li> <i>H<sup>2</sup></i>
* implies that the function is
* continuously differentiable
* over cell boundaries.
*
* <li> <i>L<sup>2</sup></i>
* indicates that the element is
* discontinuous. Since
* discontinuous elements have no
* topological couplings between
* grid cells and code may
* actually depend on this
* property, <i>L<sup>2</sup></i>
* conformity is handled in a
* special way in the sense that
* it is <b>not</b> implied by
* any higher conformity.
* </ol>
*
* In order to test if a finite
* element conforms to a certain
* space, use
* FiniteElementData<dim>::conforms().
*/
enum Conformity
{
/**
* Indicates incompatible
* continuities of a
* system.
*/
unknown = 0x00,
/**
* Discontinuous
* elements. See above!
*/
L2 = 0x01,
/**
* Conformity with the
* space
* <i>H<sup>curl</sup></i>
* (continuous tangential
* component of a vector
* field)
*/
Hcurl = 0x02,
/**
* Conformity with the
* space
* <i>H<sup>div</sup></i>
* (continuous normal
* component of a vector
* field)
*/
Hdiv = 0x04,
/**
* Conformity with the
* space
* <i>H<sup>1</sup></i>
* (continuous)
*/
H1 = Hcurl | Hdiv,
/**
* Conformity with the
* space
* <i>H<sup>2</sup></i>
* (continuously
* differentiable)
*/
H2 = 0x0e
};
/**
* The dimension of the finite
* element, which is the template
* parameter <tt>dim</tt>
*/
static const unsigned int dimension = dim;
/**
* Number of degrees of freedom on
* a vertex.
*/
const unsigned int dofs_per_vertex;
/** Number of degrees of freedom
* in a line; not including the
* degrees of freedom on the
* vertices of the line.
*/
const unsigned int dofs_per_line;
/** Number of degrees of freedom
* in a quadrilateral; not
* including the degrees of
* freedom on the lines and
* vertices of the
* quadrilateral.
*/
const unsigned int dofs_per_quad;
/** Number of degrees of freedom
* in a hexahedron; not
* including the degrees of
* freedom on the
* quadrilaterals, lines and
* vertices of the hecahedron.
*/
const unsigned int dofs_per_hex;
/**
* First index of dof on a line.
*/
const unsigned int first_line_index;
/**
* First index of dof on a quad.
*/
const unsigned int first_quad_index;
/**
* First index of dof on a hexahedron.
*/
const unsigned int first_hex_index;
/**
* First index of dof on a line for face data.
*/
const unsigned int first_face_line_index;
/**
* First index of dof on a quad for face data.
*/
const unsigned int first_face_quad_index;
/**
* Number of degrees of freedom
* on a face. This is the
* accumulated number of degrees
* of freedom on all the objects
* of dimension up to
* <tt>dim-1</tt> constituting a
* face.
*/
const unsigned int dofs_per_face;
/**
* Total number of degrees of freedom
* on a cell. This is the
* accumulated number of degrees
* of freedom on all the objects
* of dimension up to
* <tt>dim</tt> constituting a
* cell.
*/
const unsigned int dofs_per_cell;
/**
* Number of vector components of
* this finite element, and
* dimension of the image
* space. For vector-valued
* finite elements (i.e. when
* this number is greater than
* one), the number of vector
* components is in many cases
* equal to the number of base
* elements glued together with
* the help of the FESystem
* class. However, for elements
* like the Nedelec element, the
* number is greater than one
* even though we only have one
* base element.
*/
const unsigned int components;
/**
* The number of vector blocks a
* BlockVector for this element
* should contain. For primitive
* elements, this is equal to
* #components, but for vector
* valued base elements it may
* differ. Actually, in systems
* this is the sum of the base
* element multiplicities.
*/
const unsigned int blocks;
/**
* Maximal polynomial degree of a
* shape function in a single
* coordinate direction.
*/
const unsigned int degree;
/**
* Indicate the space this element conforms to.
*/
const Conformity conforming_space;
/**
* Default
* constructor. Constructs an
* element with no dofs. Checking
* n_dofs_per_cell() is therefore
* a good way to check if
* something went wrong.
*/
FiniteElementData ();
/**
* Constructor, computing all
* necessary values from the
* distribution of dofs to
* geometrcal objects.
*
* @param dofs_per_object Number
* of dofs on geometrical objects
* for each dimension. In this
* vector, entry 0 refers to dofs
* on vertices, entry 1 on lines
* and so on. Its length must be
* <i>dim+1</i>.
* @param n_components Number of
* vector components of the
* element.
* @param degree
* Maximal polynomial degree in a
* single direction.
* @param conformity The finite
* element space has continuity
* of this Sobolev space.
* @param n_blocks Number of
* vector blocks.
*/
FiniteElementData (const std::vector<unsigned int> &dofs_per_object,
const unsigned int n_components,
const unsigned int degree,
const Conformity conformity = unknown,
const unsigned int n_blocks = numbers::invalid_unsigned_int);
/**
* Number of dofs per vertex.
*/
unsigned int n_dofs_per_vertex () const;
/**
* Number of dofs per line. Not
* including dofs on lower
* dimensional objects.
*/
unsigned int n_dofs_per_line () const;
/**
* Number of dofs per quad. Not
* including dofs on lower
* dimensional objects.
*/
unsigned int n_dofs_per_quad () const;
/**
* Number of dofs per hex. Not
* including dofs on lower
* dimensional objects.
*/
unsigned int n_dofs_per_hex () const;
/**
* Number of dofs per face,
* accumulating degrees of
* freedom of all lower
* dimensional objects.
*/
unsigned int n_dofs_per_face () const;
/**
* Number of dofs per cell,
* accumulating degrees of
* freedom of all lower
* dimensional objects.
*/
unsigned int n_dofs_per_cell () const;
/**
* Return the number of degrees
* per structdim-dimensional
* object. For structdim==0, the
* function therefore returns
* dofs_per_vertex, for
* structdim==1 dofs_per_line,
* etc. This function is mostly
* used to allow some template
* trickery for functions that
* should work on all sorts of
* objects without wanting to use
* the different names (vertex,
* line, ...) associated with
* these objects.
*/
template <int structdim>
unsigned int n_dofs_per_object () const;
/**
* Number of components.
*/
unsigned int n_components () const;
/**
* Number of blocks.
*/
unsigned int n_blocks () const;
/**
* Maximal polynomial degree of a
* shape function in a single
* coordinate direction.
*
* This function can be used to
* determine the optimal
* quadrature rule.
*/
unsigned int tensor_degree () const;
/**
* Test whether a finite element
* space conforms to a certain
* Sobolev space.
*
* @note This function will
* return a true value even if
* the finite element space has
* higher regularity than asked
* for.
*/
bool conforms (const Conformity) const;
/**
* Given an index in the natural
* ordering of indices on a face,
* return the index of the same
* degree of freedom on the cell.
*
*/
unsigned int face_to_cell_index (const unsigned int face_index,
const unsigned int face,
const bool face_orientation = true,
const bool face_flip = false,
const bool face_rotation = false) const;
/**
* @deprecated This function is
* just a special version of
* face_to_cell_index for the face
* zero. It is therefore of
* limited use in aplications and
* most of the time, the other
* function is what is required.
*
* Given an index in the natural
* ordering of indices on a face,
* return the index of an
* equivalent degree of freedom
* on the cell.
*
* To explain the concept,
* consider the case where we
* would like to know whether a
* degree of freedom on a face is
* primitive. Unfortunately, the
* is_primitive() function in the
* FiniteElement class takes a
* cell index, so we would need
* to find the cell index of the
* shape function that
* corresponds to the present
* face index. This function does
* that.
*
* Code implementing this would
* then look like this:
* @code
* for (i=0; i<dofs_per_face; ++i)
* if (fe.is_primitive(fe.face_to_equivalent_cell_index(i)))
* ... do whatever
* @endcode
*
*/
unsigned int face_to_equivalent_cell_index (const unsigned int index) const;
/**
* Comparison operator.
*/
bool operator == (const FiniteElementData &) const;
};
// --------- inline and template functions ---------------
#ifndef DOXYGEN
namespace FiniteElementDomination
{
inline
Domination operator & (const Domination d1,
const Domination d2)
{
// go through the entire list of
// possibilities. note that if we were
// into speed, obfuscation and cared
// enough, we could implement this
// operator by doing a bitwise & (and) if
// we gave these values to the enum
// values: neither_element_dominates=0,
// this_element_dominates=1,
// other_element_dominates=2,
// either_element_can_dominate=3
// =this_element_dominates|other_element_dominates
switch (d1)
{
case this_element_dominates:
if ((d2 == this_element_dominates) ||
(d2 == either_element_can_dominate))
return this_element_dominates;
else
return neither_element_dominates;
case other_element_dominates:
if ((d2 == other_element_dominates) ||
(d2 == either_element_can_dominate))
return other_element_dominates;
else
return neither_element_dominates;
case neither_element_dominates:
return neither_element_dominates;
case either_element_can_dominate:
return d2;
default:
// shouldn't get here
Assert (false, ExcInternalError());
}
return neither_element_dominates;
}
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::n_dofs_per_vertex () const
{
return dofs_per_vertex;
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::n_dofs_per_line () const
{
return dofs_per_line;
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::n_dofs_per_quad () const
{
return dofs_per_quad;
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::n_dofs_per_hex () const
{
return dofs_per_hex;
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::n_dofs_per_face () const
{
return dofs_per_face;
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::n_dofs_per_cell () const
{
return dofs_per_cell;
}
template <int dim>
template <int structdim>
inline
unsigned int
FiniteElementData<dim>::n_dofs_per_object () const
{
switch (structdim)
{
case 0:
return dofs_per_vertex;
case 1:
return dofs_per_line;
case 2:
return dofs_per_quad;
case 3:
return dofs_per_hex;
default:
Assert (false, ExcInternalError());
}
return numbers::invalid_unsigned_int;
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::n_components () const
{
return components;
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::n_blocks () const
{
return blocks;
}
template <int dim>
inline
unsigned int
FiniteElementData<dim>::tensor_degree () const
{
return degree;
}
template <int dim>
inline
bool
FiniteElementData<dim>::conforms (const Conformity space) const
{
return ((space & conforming_space) == space);
}
template <>
inline
unsigned int
FiniteElementData<1>::
face_to_equivalent_cell_index (const unsigned int index) const
{
Assert (index < dofs_per_face,
ExcIndexRange (index, 0, dofs_per_face));
return index;
}
template <>
inline
unsigned int
FiniteElementData<2>::
face_to_equivalent_cell_index (const unsigned int index) const
{
Assert (index < dofs_per_face,
ExcIndexRange (index, 0, dofs_per_face));
// on face 0, the vertices are 0 and 2
if (index < this->dofs_per_vertex)
return index;
else if (index < 2*this->dofs_per_vertex)
return index + this->dofs_per_vertex;
else
// this is a dof on line 0, so on the
// cell there are now 4 vertices instead
// of only 2 ahead of this one
return index + 2*this->dofs_per_vertex;
}
template <>
inline
unsigned int
FiniteElementData<3>::
face_to_equivalent_cell_index (const unsigned int index) const
{
// this case is just way too
// complicated. fall back to
// face_to_cell_index
return face_to_cell_index(index, 0, true);
}
#endif // DOXYGEN
DEAL_II_NAMESPACE_CLOSE
#endif
|