Reference documentation for deal.II version 8.4.2
fe_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__fe_base_h
17 #define dealii__fe_base_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/base/exceptions.h>
21 #include <deal.II/base/subscriptor.h>
22 #include <deal.II/base/point.h>
23 #include <deal.II/base/tensor.h>
24 #include <deal.II/base/table.h>
25 #include <deal.II/base/vector_slice.h>
26 #include <deal.II/base/geometry_info.h>
27 #include <deal.II/lac/full_matrix.h>
28 #include <deal.II/lac/block_indices.h>
29 #include <deal.II/fe/fe_update_flags.h>
30 
31 #include <string>
32 #include <vector>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 template<int dim, int spacedim> class FESystem;
37 
38 
44 {
99  {
100  this_element_dominates,
101  other_element_dominates,
102  neither_element_dominates,
103  either_element_can_dominate,
104  no_requirements
105  };
106 
107 
126  inline Domination operator & (const Domination d1,
127  const Domination d2);
128 }
129 
130 
147 template <int dim>
149 {
150 public:
182  {
186  unknown = 0x00,
187 
191  L2 = 0x01,
192 
197  Hcurl = 0x02,
198 
203  Hdiv = 0x04,
204 
208  H1 = Hcurl | Hdiv,
209 
214  H2 = 0x0e
215  };
216 
221  static const unsigned int dimension = dim;
222 
226  const unsigned int dofs_per_vertex;
227 
232  const unsigned int dofs_per_line;
233 
238  const unsigned int dofs_per_quad;
239 
244  const unsigned int dofs_per_hex;
245 
249  const unsigned int first_line_index;
250 
254  const unsigned int first_quad_index;
255 
259  const unsigned int first_hex_index;
260 
264  const unsigned int first_face_line_index;
265 
269  const unsigned int first_face_quad_index;
270 
276  const unsigned int dofs_per_face;
277 
283  const unsigned int dofs_per_cell;
284 
293  const unsigned int components;
294 
299  const unsigned int degree;
300 
305 
312 
351  FiniteElementData (const std::vector<unsigned int> &dofs_per_object,
352  const unsigned int n_components,
353  const unsigned int degree,
354  const Conformity conformity = unknown,
355  const BlockIndices &block_indices = BlockIndices());
356 
360  unsigned int n_dofs_per_vertex () const;
361 
365  unsigned int n_dofs_per_line () const;
366 
370  unsigned int n_dofs_per_quad () const;
371 
375  unsigned int n_dofs_per_hex () const;
376 
381  unsigned int n_dofs_per_face () const;
382 
387  unsigned int n_dofs_per_cell () const;
388 
397  template <int structdim>
398  unsigned int n_dofs_per_object () const;
399 
405  unsigned int n_components () const;
406 
412  unsigned int n_blocks () const;
413 
417  const BlockIndices &block_indices() const;
418 
427  bool is_primitive () const;
428 
435  unsigned int tensor_degree () const;
436 
443  bool conforms (const Conformity) const;
444 
448  bool operator == (const FiniteElementData &) const;
449 
450 protected:
451 
458  void set_primitivity(const bool value);
459 
460 private:
467 };
468 
469 
470 
471 // --------- inline and template functions ---------------
472 
473 
474 #ifndef DOXYGEN
475 
476 namespace FiniteElementDomination
477 {
478  inline
480  const Domination d2)
481  {
482  // go through the entire list of possibilities. note that if we were into
483  // speed, obfuscation and cared enough, we could implement this operator
484  // by doing a bitwise & (and) if we gave these values to the enum values:
485  // neither_element_dominates=0, this_element_dominates=1,
486  // other_element_dominates=2, either_element_can_dominate=3
487  // =this_element_dominates|other_element_dominates
488  switch (d1)
489  {
490  case this_element_dominates:
491  if ((d2 == this_element_dominates) ||
492  (d2 == either_element_can_dominate) ||
493  (d2 == no_requirements))
494  return this_element_dominates;
495  else
496  return neither_element_dominates;
497 
498  case other_element_dominates:
499  if ((d2 == other_element_dominates) ||
500  (d2 == either_element_can_dominate) ||
501  (d2 == no_requirements))
502  return other_element_dominates;
503  else
504  return neither_element_dominates;
505 
506  case neither_element_dominates:
507  return neither_element_dominates;
508 
509  case either_element_can_dominate:
510  if (d2 == no_requirements)
511  return either_element_can_dominate;
512  else
513  return d2;
514 
515  case no_requirements:
516  return d2;
517 
518  default:
519  // shouldn't get here
520  Assert (false, ExcInternalError());
521  }
522 
523  return neither_element_dominates;
524  }
525 }
526 
527 
528 template <int dim>
529 inline
530 unsigned int
532 {
533  return dofs_per_vertex;
534 }
535 
536 
537 
538 template <int dim>
539 inline
540 unsigned int
542 {
543  return dofs_per_line;
544 }
545 
546 
547 
548 template <int dim>
549 inline
550 unsigned int
552 {
553  return dofs_per_quad;
554 }
555 
556 
557 
558 template <int dim>
559 inline
560 unsigned int
562 {
563  return dofs_per_hex;
564 }
565 
566 
567 
568 template <int dim>
569 inline
570 unsigned int
572 {
573  return dofs_per_face;
574 }
575 
576 
577 
578 template <int dim>
579 inline
580 unsigned int
582 {
583  return dofs_per_cell;
584 }
585 
586 
587 
588 template <int dim>
589 template <int structdim>
590 inline
591 unsigned int
593 {
594  switch (structdim)
595  {
596  case 0:
597  return dofs_per_vertex;
598  case 1:
599  return dofs_per_line;
600  case 2:
601  return dofs_per_quad;
602  case 3:
603  return dofs_per_hex;
604  default:
605  Assert (false, ExcInternalError());
606  }
608 }
609 
610 
611 
612 template <int dim>
613 inline
614 unsigned int
616 {
617  return components;
618 }
619 
620 
621 
622 template <int dim>
623 inline
624 bool
626 {
627  return cached_primitivity;
628 }
629 
630 
631 template <int dim>
632 inline
633 void
635 {
636  cached_primitivity = value;
637 }
638 
639 
640 template <int dim>
641 inline
642 const BlockIndices &
644 {
645  return block_indices_data;
646 }
647 
648 
649 
650 template <int dim>
651 inline
652 unsigned int
654 {
655  return block_indices_data.size();
656 }
657 
658 
659 
660 template <int dim>
661 inline
662 unsigned int
664 {
665  return degree;
666 }
667 
668 
669 template <int dim>
670 inline
671 bool
673 {
674  return ((space & conforming_space) == space);
675 }
676 
677 
678 #endif // DOXYGEN
679 
680 
681 DEAL_II_NAMESPACE_CLOSE
682 
683 #endif
const unsigned int first_hex_index
Definition: fe_base.h:259
bool is_primitive() const
static const unsigned int invalid_unsigned_int
Definition: types.h:164
const unsigned int components
Definition: fe_base.h:293
unsigned int n_dofs_per_object() const
Auxiliary class aiding in the handling of block structures like in BlockVector or FESystem...
Definition: block_indices.h:54
const unsigned int dofs_per_quad
Definition: fe_base.h:238
const unsigned int degree
Definition: fe_base.h:299
unsigned int n_dofs_per_vertex() const
Domination operator&(const Domination d1, const Domination d2)
unsigned int n_dofs_per_hex() const
const unsigned int dofs_per_line
Definition: fe_base.h:232
unsigned int n_dofs_per_face() const
unsigned int n_blocks() const
const unsigned int first_face_line_index
Definition: fe_base.h:264
const unsigned int first_quad_index
Definition: fe_base.h:254
unsigned int n_dofs_per_line() const
const unsigned int dofs_per_hex
Definition: fe_base.h:244
#define Assert(cond, exc)
Definition: exceptions.h:294
Definition: fe.h:34
bool cached_primitivity
Definition: fe_base.h:466
const unsigned int dofs_per_cell
Definition: fe_base.h:283
const BlockIndices block_indices_data
Definition: fe_base.h:311
unsigned int n_dofs_per_quad() const
unsigned int n_components() const
unsigned int n_dofs_per_cell() const
const unsigned int first_face_quad_index
Definition: fe_base.h:269
bool conforms(const Conformity) const
const Conformity conforming_space
Definition: fe_base.h:304
const unsigned int dofs_per_face
Definition: fe_base.h:276
const unsigned int first_line_index
Definition: fe_base.h:249
const unsigned int dofs_per_vertex
Definition: fe_base.h:226
void set_primitivity(const bool value)
unsigned int size() const
const BlockIndices & block_indices() const
unsigned int tensor_degree() const