Reference documentation for deal.II version 9.1.0-pre
fe_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2017 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_fe_base_h
17 #define dealii_fe_base_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/geometry_info.h>
23 #include <deal.II/base/point.h>
24 #include <deal.II/base/subscriptor.h>
25 #include <deal.II/base/table.h>
26 #include <deal.II/base/tensor.h>
27 #include <deal.II/base/vector_slice.h>
28 
29 #include <deal.II/fe/fe_update_flags.h>
30 
31 #include <deal.II/lac/block_indices.h>
32 #include <deal.II/lac/full_matrix.h>
33 
34 #include <string>
35 #include <vector>
36 
37 DEAL_II_NAMESPACE_OPEN
38 
44 {
99  {
120  };
121 
122 
141  inline Domination operator&(const Domination d1, const Domination d2);
142 } // namespace FiniteElementDomination
143 
144 
161 template <int dim>
163 {
164 public:
196  {
200  unknown = 0x00,
201 
205  L2 = 0x01,
206 
211  Hcurl = 0x02,
212 
217  Hdiv = 0x04,
218 
222  H1 = Hcurl | Hdiv,
223 
228  H2 = 0x0e
229  };
230 
235  static const unsigned int dimension = dim;
236 
240  const unsigned int dofs_per_vertex;
241 
246  const unsigned int dofs_per_line;
247 
252  const unsigned int dofs_per_quad;
253 
258  const unsigned int dofs_per_hex;
259 
263  const unsigned int first_line_index;
264 
268  const unsigned int first_quad_index;
269 
273  const unsigned int first_hex_index;
274 
278  const unsigned int first_face_line_index;
279 
283  const unsigned int first_face_quad_index;
284 
290  const unsigned int dofs_per_face;
291 
297  const unsigned int dofs_per_cell;
298 
307  const unsigned int components;
308 
313  const unsigned int degree;
314 
319 
326 
365  FiniteElementData(const std::vector<unsigned int> &dofs_per_object,
366  const unsigned int n_components,
367  const unsigned int degree,
368  const Conformity conformity = unknown,
369  const BlockIndices &block_indices = BlockIndices());
370 
374  unsigned int
375  n_dofs_per_vertex() const;
376 
380  unsigned int
381  n_dofs_per_line() const;
382 
386  unsigned int
387  n_dofs_per_quad() const;
388 
392  unsigned int
393  n_dofs_per_hex() const;
394 
399  unsigned int
400  n_dofs_per_face() const;
401 
406  unsigned int
407  n_dofs_per_cell() const;
408 
417  template <int structdim>
418  unsigned int
419  n_dofs_per_object() const;
420 
426  unsigned int
427  n_components() const;
428 
434  unsigned int
435  n_blocks() const;
436 
440  const BlockIndices &
441  block_indices() const;
442 
449  unsigned int
450  tensor_degree() const;
451 
458  bool
459  conforms(const Conformity) const;
460 
464  bool
465  operator==(const FiniteElementData &) const;
466 };
467 
468 
469 
470 // --------- inline and template functions ---------------
471 
472 
473 #ifndef DOXYGEN
474 
475 namespace FiniteElementDomination
476 {
477  inline Domination operator&(const Domination d1, const Domination d2)
478  {
479  // go through the entire list of possibilities. note that if we were into
480  // speed, obfuscation and cared enough, we could implement this operator
481  // by doing a bitwise & (and) if we gave these values to the enum values:
482  // neither_element_dominates=0, this_element_dominates=1,
483  // other_element_dominates=2, either_element_can_dominate=3
484  // =this_element_dominates|other_element_dominates
485  switch (d1)
486  {
488  if ((d2 == this_element_dominates) ||
489  (d2 == either_element_can_dominate) || (d2 == no_requirements))
490  return this_element_dominates;
491  else
493 
495  if ((d2 == other_element_dominates) ||
496  (d2 == either_element_can_dominate) || (d2 == no_requirements))
498  else
500 
503 
505  if (d2 == no_requirements)
507  else
508  return d2;
509 
510  case no_requirements:
511  return d2;
512 
513  default:
514  // shouldn't get here
515  Assert(false, ExcInternalError());
516  }
517 
519  }
520 } // namespace FiniteElementDomination
521 
522 
523 template <int dim>
524 inline unsigned int
526 {
527  return dofs_per_vertex;
528 }
529 
530 
531 
532 template <int dim>
533 inline unsigned int
535 {
536  return dofs_per_line;
537 }
538 
539 
540 
541 template <int dim>
542 inline unsigned int
544 {
545  return dofs_per_quad;
546 }
547 
548 
549 
550 template <int dim>
551 inline unsigned int
553 {
554  return dofs_per_hex;
555 }
556 
557 
558 
559 template <int dim>
560 inline unsigned int
562 {
563  return dofs_per_face;
564 }
565 
566 
567 
568 template <int dim>
569 inline unsigned int
571 {
572  return dofs_per_cell;
573 }
574 
575 
576 
577 template <int dim>
578 template <int structdim>
579 inline unsigned int
581 {
582  switch (structdim)
583  {
584  case 0:
585  return dofs_per_vertex;
586  case 1:
587  return dofs_per_line;
588  case 2:
589  return dofs_per_quad;
590  case 3:
591  return dofs_per_hex;
592  default:
593  Assert(false, ExcInternalError());
594  }
596 }
597 
598 
599 
600 template <int dim>
601 inline unsigned int
603 {
604  return components;
605 }
606 
607 
608 
609 template <int dim>
610 inline const BlockIndices &
612 {
613  return block_indices_data;
614 }
615 
616 
617 
618 template <int dim>
619 inline unsigned int
621 {
622  return block_indices_data.size();
623 }
624 
625 
626 
627 template <int dim>
628 inline unsigned int
630 {
631  return degree;
632 }
633 
634 
635 template <int dim>
636 inline bool
638 {
639  return ((space & conforming_space) == space);
640 }
641 
642 
643 #endif // DOXYGEN
644 
645 
646 DEAL_II_NAMESPACE_CLOSE
647 
648 #endif
const unsigned int first_hex_index
Definition: fe_base.h:273
static const unsigned int invalid_unsigned_int
Definition: types.h:173
const unsigned int components
Definition: fe_base.h:307
const unsigned int dofs_per_quad
Definition: fe_base.h:252
const unsigned int degree
Definition: fe_base.h:313
const unsigned int dofs_per_line
Definition: fe_base.h:246
const unsigned int first_face_line_index
Definition: fe_base.h:278
const BlockIndices & block_indices() const
unsigned int tensor_degree() const
unsigned int n_dofs_per_face() const
const unsigned int first_quad_index
Definition: fe_base.h:268
const unsigned int dofs_per_hex
Definition: fe_base.h:258
#define Assert(cond, exc)
Definition: exceptions.h:1227
unsigned int n_components() const
unsigned int n_dofs_per_vertex() const
const unsigned int dofs_per_cell
Definition: fe_base.h:297
const BlockIndices block_indices_data
Definition: fe_base.h:325
Domination operator&(const Domination d1, const Domination d2)
bool conforms(const Conformity) const
unsigned int n_dofs_per_object() const
const unsigned int first_face_quad_index
Definition: fe_base.h:283
unsigned int n_dofs_per_line() const
const Conformity conforming_space
Definition: fe_base.h:318
unsigned int n_blocks() const
const unsigned int dofs_per_face
Definition: fe_base.h:290
const unsigned int first_line_index
Definition: fe_base.h:263
unsigned int n_dofs_per_cell() const
unsigned int size() const
const unsigned int dofs_per_vertex
Definition: fe_base.h:240
unsigned int n_dofs_per_hex() const
unsigned int n_dofs_per_quad() const
static::ExceptionBase & ExcInternalError()