Reference documentation for deal.II version 8.5.1
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 
41 {
96  {
117  };
118 
119 
138  inline Domination operator & (const Domination d1,
139  const Domination d2);
140 }
141 
142 
159 template <int dim>
161 {
162 public:
194  {
198  unknown = 0x00,
199 
203  L2 = 0x01,
204 
209  Hcurl = 0x02,
210 
215  Hdiv = 0x04,
216 
220  H1 = Hcurl | Hdiv,
221 
226  H2 = 0x0e
227  };
228 
233  static const unsigned int dimension = dim;
234 
238  const unsigned int dofs_per_vertex;
239 
244  const unsigned int dofs_per_line;
245 
250  const unsigned int dofs_per_quad;
251 
256  const unsigned int dofs_per_hex;
257 
261  const unsigned int first_line_index;
262 
266  const unsigned int first_quad_index;
267 
271  const unsigned int first_hex_index;
272 
276  const unsigned int first_face_line_index;
277 
281  const unsigned int first_face_quad_index;
282 
288  const unsigned int dofs_per_face;
289 
295  const unsigned int dofs_per_cell;
296 
305  const unsigned int components;
306 
311  const unsigned int degree;
312 
317 
324 
363  FiniteElementData (const std::vector<unsigned int> &dofs_per_object,
364  const unsigned int n_components,
365  const unsigned int degree,
366  const Conformity conformity = unknown,
368 
372  unsigned int n_dofs_per_vertex () const;
373 
377  unsigned int n_dofs_per_line () const;
378 
382  unsigned int n_dofs_per_quad () const;
383 
387  unsigned int n_dofs_per_hex () const;
388 
393  unsigned int n_dofs_per_face () const;
394 
399  unsigned int n_dofs_per_cell () const;
400 
409  template <int structdim>
410  unsigned int n_dofs_per_object () const;
411 
417  unsigned int n_components () const;
418 
424  unsigned int n_blocks () const;
425 
429  const BlockIndices &block_indices() const;
430 
437  unsigned int tensor_degree () const;
438 
445  bool conforms (const Conformity) const;
446 
450  bool operator == (const FiniteElementData &) const;
451 };
452 
453 
454 
455 // --------- inline and template functions ---------------
456 
457 
458 #ifndef DOXYGEN
459 
460 namespace FiniteElementDomination
461 {
462  inline
464  const Domination d2)
465  {
466  // go through the entire list of possibilities. note that if we were into
467  // speed, obfuscation and cared enough, we could implement this operator
468  // by doing a bitwise & (and) if we gave these values to the enum values:
469  // neither_element_dominates=0, this_element_dominates=1,
470  // other_element_dominates=2, either_element_can_dominate=3
471  // =this_element_dominates|other_element_dominates
472  switch (d1)
473  {
475  if ((d2 == this_element_dominates) ||
476  (d2 == either_element_can_dominate) ||
477  (d2 == no_requirements))
478  return this_element_dominates;
479  else
481 
483  if ((d2 == other_element_dominates) ||
484  (d2 == either_element_can_dominate) ||
485  (d2 == no_requirements))
487  else
489 
492 
494  if (d2 == no_requirements)
496  else
497  return d2;
498 
499  case no_requirements:
500  return d2;
501 
502  default:
503  // shouldn't get here
504  Assert (false, ExcInternalError());
505  }
506 
508  }
509 }
510 
511 
512 template <int dim>
513 inline
514 unsigned int
516 {
517  return dofs_per_vertex;
518 }
519 
520 
521 
522 template <int dim>
523 inline
524 unsigned int
526 {
527  return dofs_per_line;
528 }
529 
530 
531 
532 template <int dim>
533 inline
534 unsigned int
536 {
537  return dofs_per_quad;
538 }
539 
540 
541 
542 template <int dim>
543 inline
544 unsigned int
546 {
547  return dofs_per_hex;
548 }
549 
550 
551 
552 template <int dim>
553 inline
554 unsigned int
556 {
557  return dofs_per_face;
558 }
559 
560 
561 
562 template <int dim>
563 inline
564 unsigned int
566 {
567  return dofs_per_cell;
568 }
569 
570 
571 
572 template <int dim>
573 template <int structdim>
574 inline
575 unsigned int
577 {
578  switch (structdim)
579  {
580  case 0:
581  return dofs_per_vertex;
582  case 1:
583  return dofs_per_line;
584  case 2:
585  return dofs_per_quad;
586  case 3:
587  return dofs_per_hex;
588  default:
589  Assert (false, ExcInternalError());
590  }
592 }
593 
594 
595 
596 template <int dim>
597 inline
598 unsigned int
600 {
601  return components;
602 }
603 
604 
605 
606 template <int dim>
607 inline
608 const BlockIndices &
610 {
611  return block_indices_data;
612 }
613 
614 
615 
616 template <int dim>
617 inline
618 unsigned int
620 {
621  return block_indices_data.size();
622 }
623 
624 
625 
626 template <int dim>
627 inline
628 unsigned int
630 {
631  return degree;
632 }
633 
634 
635 template <int dim>
636 inline
637 bool
638 FiniteElementData<dim>::conforms (const Conformity space) const
639 {
640  return ((space & conforming_space) == space);
641 }
642 
643 
644 #endif // DOXYGEN
645 
646 
647 DEAL_II_NAMESPACE_CLOSE
648 
649 #endif
const unsigned int first_hex_index
Definition: fe_base.h:271
static const unsigned int invalid_unsigned_int
Definition: types.h:170
const unsigned int components
Definition: fe_base.h:305
unsigned int n_dofs_per_object() const
const unsigned int dofs_per_quad
Definition: fe_base.h:250
const unsigned int degree
Definition: fe_base.h:311
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:244
unsigned int n_dofs_per_face() const
unsigned int n_blocks() const
const unsigned int first_face_line_index
Definition: fe_base.h:276
bool operator==(const FiniteElementData &) const
Definition: fe_data.cc:72
const unsigned int first_quad_index
Definition: fe_base.h:266
unsigned int n_dofs_per_line() const
const unsigned int dofs_per_hex
Definition: fe_base.h:256
#define Assert(cond, exc)
Definition: exceptions.h:313
const unsigned int dofs_per_cell
Definition: fe_base.h:295
static const unsigned int dimension
Definition: fe_base.h:233
const BlockIndices block_indices_data
Definition: fe_base.h:323
FiniteElementData(const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices())
Definition: fe_data.cc:23
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:281
bool conforms(const Conformity) const
const Conformity conforming_space
Definition: fe_base.h:316
const unsigned int dofs_per_face
Definition: fe_base.h:288
const unsigned int first_line_index
Definition: fe_base.h:261
const unsigned int dofs_per_vertex
Definition: fe_base.h:238
unsigned int size() const
const BlockIndices & block_indices() const
unsigned int tensor_degree() const
static ::ExceptionBase & ExcInternalError()