Reference documentation for deal.II version 8.5.1
|
#include <deal.II/dofs/dof_handler.h>
Classes | |
class | MGVertexDoFs |
Public Types | |
typedef ActiveSelector::active_cell_iterator | active_cell_iterator |
typedef ActiveSelector::cell_iterator | cell_iterator |
Public Member Functions | |
DoFHandler () | |
DoFHandler (const Triangulation< dim, spacedim > &tria) | |
virtual | ~DoFHandler () |
void | initialize (const Triangulation< dim, spacedim > &tria, const FiniteElement< dim, spacedim > &fe) |
virtual void | distribute_dofs (const FiniteElement< dim, spacedim > &fe) |
virtual void | distribute_mg_dofs (const FiniteElement< dim, spacedim > &fe) |
bool | has_level_dofs () const |
bool | has_active_dofs () const |
void | initialize_local_block_info () |
virtual void | clear () |
void | renumber_dofs (const std::vector< types::global_dof_index > &new_numbers) |
void | renumber_dofs (const unsigned int level, const std::vector< types::global_dof_index > &new_numbers) |
unsigned int | max_couplings_between_dofs () const |
unsigned int | max_couplings_between_boundary_dofs () const |
Cell iterator functions | |
cell_iterator | begin (const unsigned int level=0) const |
active_cell_iterator | begin_active (const unsigned int level=0) const |
cell_iterator | end () const |
cell_iterator | end (const unsigned int level) const |
active_cell_iterator | end_active (const unsigned int level) const |
level_cell_iterator | begin_mg (const unsigned int level=0) const |
level_cell_iterator | end_mg (const unsigned int level) const |
level_cell_iterator | end_mg () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Attributes | |
static const unsigned int | dimension = dim |
static const unsigned int | space_dimension = spacedim |
static const types::global_dof_index | invalid_dof_index = numbers::invalid_dof_index |
static const unsigned int | default_fe_index = 0 |
Related Functions | |
(Note that these are not member functions.) | |
Functions to support code that generically uses both DoFHandler and hp::DoFHandler | |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_cell (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_face (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | max_dofs_per_vertex (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
unsigned int | n_components (const DoFHandler< dim, spacedim > &dh) |
template<int dim, int spacedim> | |
bool | fe_is_primitive (const DoFHandler< dim, spacedim > &dh) |
Cell iterator functions returning ranges of iterators | |
BlockInfo | block_info_object |
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > | tria |
SmartPointer< const FiniteElement< dim, spacedim >, DoFHandler< dim, spacedim > > | selected_fe |
std_cxx11::shared_ptr<::internal::DoFHandler::Policy::PolicyBase< dim, spacedim > > | policy |
::internal::DoFHandler::NumberCache | number_cache |
std::vector<::internal::DoFHandler::NumberCache > | mg_number_cache |
std::vector< types::global_dof_index > | vertex_dofs |
std::vector< MGVertexDoFs > | mg_vertex_dofs |
std::vector<::internal::DoFHandler::DoFLevel< dim > * > | levels |
std::vector<::internal::DoFHandler::DoFLevel< dim > * > | mg_levels |
::internal::DoFHandler::DoFFaces< dim > * | faces |
::internal::DoFHandler::DoFFaces< dim > * | mg_faces |
template<int , class , bool > | |
class | DoFAccessor |
template<class , bool > | |
class | DoFCellAccessor |
struct | ::internal::DoFAccessor::Implementation |
struct | ::internal::DoFCellAccessor::Implementation |
struct | ::internal::DoFHandler::Implementation |
struct | ::internal::DoFHandler::Policy::Implementation |
IteratorRange< cell_iterator > | cell_iterators () const |
IteratorRange< active_cell_iterator > | active_cell_iterators () const |
IteratorRange< level_cell_iterator > | mg_cell_iterators () const |
IteratorRange< cell_iterator > | cell_iterators_on_level (const unsigned int level) const |
IteratorRange< active_cell_iterator > | active_cell_iterators_on_level (const unsigned int level) const |
IteratorRange< level_cell_iterator > | mg_cell_iterators_on_level (const unsigned int level) const |
types::global_dof_index | n_dofs () const |
types::global_dof_index | n_dofs (const unsigned int level) const |
types::global_dof_index | n_boundary_dofs () const |
template<typename number > | |
types::global_dof_index | n_boundary_dofs (const std::map< types::boundary_id, const Function< spacedim, number > *> &boundary_ids) const |
types::global_dof_index | n_boundary_dofs (const std::set< types::boundary_id > &boundary_ids) const |
const BlockInfo & | block_info () const |
unsigned int | n_locally_owned_dofs () const |
const IndexSet & | locally_owned_dofs () const |
const IndexSet & | locally_owned_mg_dofs (const unsigned int level) const |
const std::vector< IndexSet > & | locally_owned_dofs_per_processor () const |
const std::vector< IndexSet > & | locally_owned_mg_dofs_per_processor (const unsigned int level) const |
const std::vector< types::global_dof_index > & | n_locally_owned_dofs_per_processor () const |
const FiniteElement< dim, spacedim > & | get_fe () const |
const Triangulation< dim, spacedim > & | get_tria () const 1 |
const Triangulation< dim, spacedim > & | get_triangulation () const |
virtual std::size_t | memory_consumption () const |
template<class Archive > | |
void | save (Archive &ar, const unsigned int version) const |
template<class Archive > | |
void | load (Archive &ar, const unsigned int version) |
static ::ExceptionBase & | ExcRenumberingIncomplete () |
static ::ExceptionBase & | ExcGridsDoNotMatch () |
static ::ExceptionBase & | ExcInvalidBoundaryIndicator () |
static ::ExceptionBase & | ExcNewNumbersNotConsecutive (types::global_dof_index arg1) |
static ::ExceptionBase & | ExcInvalidLevel (int arg1) |
static ::ExceptionBase & | ExcFacesHaveNoLevel () |
static ::ExceptionBase & | ExcEmptyLevel (int arg1) |
DoFHandler (const DoFHandler &) | |
DoFHandler & | operator= (const DoFHandler &) |
void | clear_mg_space () |
void | clear_space () |
void | reserve_space () |
template<int structdim> | |
types::global_dof_index | get_dof_index (const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index) const |
template<int structdim> | |
void | set_dof_index (const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index) const |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Given a triangulation and a description of a finite element, this class enumerates degrees of freedom on all vertices, edges, faces, and cells of the triangulation. As a result, it also provides a basis for a discrete space \(V_h\) whose elements are finite element functions defined on each cell by a FiniteElement object. This class satisfies the MeshType concept requirements.
It is first used in the step-2 tutorial program.
For each vertex, line, quad, etc, this class stores a list of the indices of degrees of freedom living on this object. These indices refer to the unconstrained degrees of freedom, i.e. constrained degrees of freedom are numbered in the same way as unconstrained ones, and are only later eliminated. This leads to the fact that indices in global vectors and matrices also refer to all degrees of freedom and some kind of condensation is needed to restrict the systems of equations to the unconstrained degrees of freedom only. The actual layout of storage of the indices is described in the internal::DoFHandler::DoFLevel class documentation.
The class offers iterators to traverse all cells, in much the same way as the Triangulation class does. Using the begin() and end() functions (and companions, like begin_active()), one can obtain iterators to walk over cells, and query the degree of freedom structures as well as the triangulation data. These iterators are built on top of those of the Triangulation class, but offer the additional information on degrees of freedom functionality compared to pure triangulation iterators. The order in which dof iterators are presented by the ++
and –
operators is the same as that for the corresponding iterators traversing the triangulation on which this DoFHandler is constructed.
The spacedim
parameter has to be used if one wants to solve problems on surfaces. If not specified, this parameter takes the default value =dim
implying that we want to solve problems in a domain whose dimension equals the dimension of the space in which it is embedded.
The degrees of freedom (‘dofs’) are distributed on the given triangulation by the function distribute_dofs(). It gets passed a finite element object describing how many degrees of freedom are located on vertices, lines, etc. It traverses the triangulation cell by cell and numbers the dofs of that cell if not yet numbered. For non-multigrid algorithms, only active cells are considered. Active cells are defined to be those cells which have no children, i.e. they are the most refined ones.
Since the triangulation is traversed starting with the cells of the coarsest active level and going to more refined levels, the lowest numbers for dofs are given to the largest cells as well as their bounding lines and vertices, with the dofs of more refined cells getting higher numbers.
This numbering implies very large bandwidths of the resulting matrices and is thus vastly suboptimal for some solution algorithms. For this reason, the DoFRenumbering class offers several algorithms to reorder the dof numbering according. See there for a discussion of the implemented algorithms.
Upon construction, this class takes a reference to a triangulation object. In most cases, this will be a reference to an object of type Triangulation, i.e. the class that represents triangulations that entirely reside on a single processor. However, it can also be of type parallel::distributed::Triangulation (see, for example, step-32, step-40 and in particular the Parallel computing with multiple processors using distributed memory module) in which case the DoFHandler object will proceed to only manage degrees of freedom on locally owned and ghost cells. This process is entirely transparent to the used.
The DoFRenumbering class offers a number of renumbering schemes like the Cuthill-McKee scheme. Basically, the function sets up an array in which for each degree of freedom we store the new index this DoF should have after renumbering. Using this array, the renumber_dofs() function of the present class is called, which actually performs the change from old DoF indices to the ones given in the array. In some cases, however, a user may want to compute her own renumbering order; in this case, one can allocate an array with one element per degree of freedom and fill it with the number that the respective degree of freedom shall be assigned. This number may, for example, be obtained by sorting the support points of the degrees of freedom in downwind direction. Then call the renumber_dofs(vector<types::global_dof_index>)
function with the array, which converts old into new degree of freedom indices.
Like many other classes in deal.II, the DoFHandler class can stream its contents to an archive using BOOST's serialization facilities. The data so stored can later be retrieved again from the archive to restore the contents of this object. This facility is frequently used to save the state of a program to disk for possible later resurrection, often in the context of checkpoint/restart strategies for long running computations or on computers that aren't very reliable (e.g. on very large clusters where individual nodes occasionally fail and then bring down an entire MPI job).
The model for doing so is similar for the DoFHandler class as it is for the Triangulation class (see the section in the general documentation of that class). In particular, the load() function does not exactly restore the same state as was stored previously using the save() function. Rather, the function assumes that you load data into a DoFHandler object that is already associated with a triangulation that has a content that matches the one that was used when the data was saved. Likewise, the load() function assumes that the current object is already associated with a finite element object that matches the one that was associated with it when data was saved; the latter can be achieved by calling DoFHandler::distribute_dofs() using the same kind of finite element before re-loading data from the serialization archive.
Definition at line 29 of file block_info.h.
DoFHandler< dim, spacedim >::DoFHandler | ( | ) |
Standard constructor, not initializing any data. After constructing an object with this constructor, use initialize() to make a valid DoFHandler.
Definition at line 678 of file dof_handler.cc.
DoFHandler< dim, spacedim >::DoFHandler | ( | const Triangulation< dim, spacedim > & | tria | ) |
Constructor. Take tria
as the triangulation to work on.
Definition at line 654 of file dof_handler.cc.
|
virtual |
Destructor.
Definition at line 688 of file dof_handler.cc.
|
private |
Copy constructor. I can see no reason why someone might want to use it, so I don't provide it. Since this class has pointer members, making it private prevents the compiler to provide it's own, incorrect one if anyone chose to copy such an object.
void DoFHandler< dim, spacedim >::initialize | ( | const Triangulation< dim, spacedim > & | tria, |
const FiniteElement< dim, spacedim > & | fe | ||
) |
Assign a Triangulation and a FiniteElement to the DoFHandler and compute the distribution of degrees of freedom over the mesh.
Definition at line 697 of file dof_handler.cc.
|
virtual |
Go through the triangulation and "distribute" the degrees of freedoms needed for the given finite element. "Distributing" degrees of freedom involved allocating memory to store the information that describes it (e.g., whether it is located on a vertex, edge, face, etc) and to sequentially enumerate all degrees of freedom. In other words, while the mesh and the finite element object by themselves simply define a finite element space \(V_h\), the process of distributing degrees of freedom makes sure that there is a basis for this space and that the shape functions of this basis are enumerated in an indexable, predictable way.
The purpose of this function is first discussed in the introduction to the step-2 tutorial program.
clear
member function which also releases the lock of this object to the finite element. Definition at line 1079 of file dof_handler.cc.
|
virtual |
Distribute level degrees of freedom on each level for geometric multigrid. The active DoFs need to be distributed using distribute_dofs() before calling this function and the fe
needs to be identical to the finite element passed to distribute_dofs().
Definition at line 1115 of file dof_handler.cc.
bool DoFHandler< dim, spacedim >::has_level_dofs | ( | ) | const |
This function returns whether this DoFHandler has DoFs distributed on each multigrid level or in other words if distribute_mg_dofs() has been called.
bool DoFHandler< dim, spacedim >::has_active_dofs | ( | ) | const |
This function returns whether this DoFHandler has active DoFs. This is equivalent to asking whether (i) distribute_dofs() has been called and (ii) the finite element for which degrees of freedom have been distributed actually has degrees of freedom (which is not the case for FE_Nothing, for example).
If this object is based on a parallel::distributed::Triangulation, then the current function returns true if any partition of the parallel DoFHandler object has any degrees of freedom. In other words, the function returns true even if the Triangulation does not own any active cells on the current MPI process, but at least one process owns cells and at least this one process has any degrees of freedom associated with it.
void DoFHandler< dim, spacedim >::initialize_local_block_info | ( | ) |
After distribute_dofs() with an FESystem element, the block structure of global and level vectors is stored in a BlockInfo object accessible with block_info(). This function initializes the local block structure on each cell in the same object.
Definition at line 1172 of file dof_handler.cc.
|
virtual |
Clear all data of this object and especially delete the lock this object has to the finite element used the last time when distribute_dofs
was called.
Definition at line 1180 of file dof_handler.cc.
void DoFHandler< dim, spacedim >::renumber_dofs | ( | const std::vector< types::global_dof_index > & | new_numbers | ) |
Renumber degrees of freedom based on a list of new dof numbers for all the dofs.
This function is called by the functions in DoFRenumbering function after computing the ordering of the degrees of freedom. This function is called, for example, by the functions in the DoFRenumbering namespace, but it can of course also be called from user code.
i
is currently locally owned, then new_numbers[locally_owned_dofs().index_within_set(i)]
returns the new global DoF index of i
. Since the IndexSet of locally_owned_dofs() is complete in the sequential case, the latter convention for the content of the array reduces to the former in the case that only one processor participates in the mesh. Definition at line 1194 of file dof_handler.cc.
void DoFHandler< dim, spacedim >::renumber_dofs | ( | const unsigned int | level, |
const std::vector< types::global_dof_index > & | new_numbers | ||
) |
The same function as above, but renumber the degrees of freedom of a single level of a multigrid hierarchy.
Definition at line 1232 of file dof_handler.cc.
unsigned int DoFHandler< dim, spacedim >::max_couplings_between_dofs | ( | ) | const |
Return the maximum number of degrees of freedom a degree of freedom in the given triangulation with the given finite element may couple with. This is the maximum number of entries per line in the system matrix; this information can therefore be used upon construction of the SparsityPattern object.
The returned number is not really the maximum number but an estimate based on the finite element and the maximum number of cells meeting at a vertex. The number holds for the constrained matrix as well.
The determination of the number of couplings can be done by simple picture drawing. An example can be found in the implementation of this function.
Definition at line 1431 of file dof_handler.cc.
unsigned int DoFHandler< dim, spacedim >::max_couplings_between_boundary_dofs | ( | ) | const |
Return the number of degrees of freedom located on the boundary another dof on the boundary can couple with.
The number is the same as for max_couplings_between_dofs() in one dimension less.
Definition at line 1440 of file dof_handler.cc.
DoFHandler< dim, spacedim >::cell_iterator DoFHandler< dim, spacedim >::begin | ( | const unsigned int | level = 0 | ) | const |
Iterator to the first used cell on level level
.
Definition at line 728 of file dof_handler.cc.
DoFHandler< dim, spacedim >::active_cell_iterator DoFHandler< dim, spacedim >::begin_active | ( | const unsigned int | level = 0 | ) | const |
Iterator to the first active cell on level level
. If the given level does not contain any active cells (i.e., all cells on this level are further refined, then this function returns end_active(level)
so that loops of the kind
have zero iterations, as may be expected if there are no active cells on this level.
Definition at line 740 of file dof_handler.cc.
DoFHandler< dim, spacedim >::cell_iterator DoFHandler< dim, spacedim >::end | ( | ) | const |
Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.
Definition at line 756 of file dof_handler.cc.
DoFHandler< dim, spacedim >::cell_iterator DoFHandler< dim, spacedim >::end | ( | const unsigned int | level | ) | const |
Return an iterator which is the first iterator not on the given level. If level
is the last level, then this returns end()
.
Definition at line 767 of file dof_handler.cc.
DoFHandler< dim, spacedim >::active_cell_iterator DoFHandler< dim, spacedim >::end_active | ( | const unsigned int | level | ) | const |
Return an active iterator which is the first active iterator not on the given level. If level
is the last level, then this returns end()
.
Definition at line 778 of file dof_handler.cc.
DoFHandler< dim, spacedim >::level_cell_iterator DoFHandler< dim, spacedim >::begin_mg | ( | const unsigned int | level = 0 | ) | const |
Iterator to the first used cell on level level
. This returns a level_cell_iterator that returns level dofs when dof_indices() is called.
Definition at line 790 of file dof_handler.cc.
DoFHandler< dim, spacedim >::level_cell_iterator DoFHandler< dim, spacedim >::end_mg | ( | const unsigned int | level | ) | const |
Iterator past the last cell on level level
. This returns a level_cell_iterator that returns level dofs when dof_indices() is called.
Definition at line 803 of file dof_handler.cc.
DoFHandler< dim, spacedim >::level_cell_iterator DoFHandler< dim, spacedim >::end_mg | ( | ) | const |
Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.
Definition at line 816 of file dof_handler.cc.
types::global_dof_index DoFHandler< dim, spacedim >::n_dofs | ( | ) | const |
Return the global number of degrees of freedom. If the current object handles all degrees of freedom itself (even if you may intend to solve your linear system in parallel, such as in step-17 or step-18), then this number equals the number of locally owned degrees of freedom since this object doesn't know anything about what you want to do with it and believes that it owns every degree of freedom it knows about.
On the other hand, if this object operates on a parallel::distributed::Triangulation object, then this function returns the global number of degrees of freedom, accumulated over all processors.
In either case, included in the returned number are those DoFs which are constrained by hanging nodes, see Constraints on degrees of freedom.
types::global_dof_index DoFHandler< dim, spacedim >::n_dofs | ( | const unsigned int | level | ) | const |
The (global) number of multilevel degrees of freedom on a given level.
If no level degrees of freedom have been assigned to this level, returns numbers::invalid_dof_index. Else returns the number of degrees of freedom on this level.
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs | ( | ) | const |
Return the number of degrees of freedom located on the boundary.
Definition at line 975 of file dof_handler.cc.
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs | ( | const std::map< types::boundary_id, const Function< spacedim, number > *> & | boundary_ids | ) | const |
Return the number of degrees of freedom located on those parts of the boundary which have a boundary indicator listed in the given set. The reason that a map
rather than a set
is used is the same as described in the section on the make_boundary_sparsity_pattern
function.
The type of boundary_ids equals typename FunctionMap<spacedim,number>::type .
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs | ( | const std::set< types::boundary_id > & | boundary_ids | ) | const |
Same function, but with different data type of the argument, which is here simply a list of the boundary indicators under consideration.
Definition at line 1014 of file dof_handler.cc.
const BlockInfo& DoFHandler< dim, spacedim >::block_info | ( | ) | const |
Access to an object informing of the block structure of the dof handler.
If an FESystem is used in distribute_dofs(), degrees of freedom naturally split into several blocks. For each base element as many blocks appear as its multiplicity.
At the end of distribute_dofs(), the number of degrees of freedom in each block is counted, and stored in a BlockInfo object, which can be accessed here. If you have previously called distribute_mg_dofs(), the same is done on each level of the multigrid hierarchy. Additionally, the block structure on each cell can be generated in this object by calling initialize_local_block_info().
unsigned int DoFHandler< dim, spacedim >::n_locally_owned_dofs | ( | ) | const |
Return the number of degrees of freedom that belong to this process.
If this is a sequential job, then the result equals that produced by n_dofs(). On the other hand, if we are operating on a parallel::distributed::Triangulation, then it includes only the degrees of freedom that the current processor owns. Note that in this case this does not include all degrees of freedom that have been distributed on the current processor's image of the mesh: in particular, some of the degrees of freedom on the interface between the cells owned by this processor and cells owned by other processors may be theirs, and degrees of freedom on ghost cells are also not necessarily included.
const IndexSet& DoFHandler< dim, spacedim >::locally_owned_dofs | ( | ) | const |
Return an IndexSet describing the set of locally owned DoFs as a subset of 0..n_dofs(). The number of elements of this set equals n_locally_owned_dofs().
const IndexSet& DoFHandler< dim, spacedim >::locally_owned_mg_dofs | ( | const unsigned int | level | ) | const |
Return an IndexSet describing the set of locally owned DoFs used for the given multigrid level as a subset of 0..n_dofs(level).
const std::vector<IndexSet>& DoFHandler< dim, spacedim >::locally_owned_dofs_per_processor | ( | ) | const |
Return a vector that stores the locally owned DoFs of each processor. If you are only interested in the number of elements each processor owns then n_locally_owned_dofs_per_processor() is a better choice.
If this is a sequential job, then the vector has a single element that equals the IndexSet representing the entire range [0,n_dofs()].
const std::vector<types::global_dof_index>& DoFHandler< dim, spacedim >::n_locally_owned_dofs_per_processor | ( | ) | const |
Return a vector that stores the number of degrees of freedom each processor that participates in this triangulation owns locally. The sum of all these numbers equals the number of degrees of freedom that exist globally, i.e. what n_dofs() returns.
Each element of the vector returned by this function equals the number of elements of the corresponding sets returned by global_dof_indices().
If this is a sequential job, then the vector has a single element equal to n_dofs().
const FiniteElement<dim,spacedim>& DoFHandler< dim, spacedim >::get_fe | ( | ) | const |
Return a constant reference to the selected finite element object.
const Triangulation<dim,spacedim>& DoFHandler< dim, spacedim >::get_tria | ( | ) | const |
Return a constant reference to the triangulation underlying this object.
const Triangulation<dim,spacedim>& DoFHandler< dim, spacedim >::get_triangulation | ( | ) | const |
Return a constant reference to the triangulation underlying this object.
|
virtual |
Determine an estimate for the memory consumption (in bytes) of this object.
This function is made virtual, since a dof handler object might be accessed through a pointers to this base class, although the actual object might be a derived class.
Definition at line 1050 of file dof_handler.cc.
void DoFHandler< dim, spacedim >::save | ( | Archive & | ar, |
const unsigned int | version | ||
) | const |
Write the data of this object to a stream for the purpose of serialization.
void DoFHandler< dim, spacedim >::load | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Read the data of this object from a stream for the purpose of serialization.
|
private |
Copy operator. I can see no reason why someone might want to use it, so I don't provide it. Since this class has pointer members, making it private prevents the compiler to provide it's own, incorrect one if anyone chose to copy such an object.
|
private |
Free all used memory.
Definition at line 1478 of file dof_handler.cc.
Make accessor objects friends.
Definition at line 1089 of file dof_handler.h.
|
related |
Maximal number of degrees of freedom on a cell.
|
related |
Maximal number of degrees of freedom on a face.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
|
related |
Maximal number of degrees of freedom on a vertex.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
|
related |
Number of vector components in the finite element object used by this DoFHandler.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
|
related |
Find out whether the FiniteElement used by this DoFHandler is primitive or not.
This function exists for both non-hp and hp DoFHandlers, to allow for a uniform interface to query this property.
|
static |
Make the dimension available in function templates.
Definition at line 271 of file dof_handler.h.
|
static |
Make the space dimension available in function templates.
Definition at line 276 of file dof_handler.h.
|
static |
When the arrays holding the DoF indices are set up, but before they are filled with actual values, they are set to an invalid value, in order to monitor possible problems. This invalid value is the constant defined here.
Please note that you should not rely on it having a certain value, but rather take its symbolic name.
Definition at line 287 of file dof_handler.h.
|
static |
The default index of the finite element to be used on a given cell. Since the present class only supports the same finite element to be used on all cells, the index of the finite element needs to be the same on all cells anyway, and by convention we pick zero for this value. The situation is different for hp objects (i.e. the hp::DoFHandler class) where different finite element indices may be used on different cells, and the default index there corresponds to an invalid value.
Definition at line 298 of file dof_handler.h.
|
private |
An object containing information on the block structure.
Definition at line 913 of file dof_handler.h.
|
private |
Address of the triangulation to work on.
Definition at line 919 of file dof_handler.h.
|
private |
Store a pointer to the finite element given latest for the distribution of dofs. In order to avoid destruction of the object before the lifetime of the DoF handler, we subscribe to the finite element object. To unlock the FE before the end of the lifetime of this DoF handler, use the clear()
function (this clears all data of this object as well, though).
Definition at line 930 of file dof_handler.h.
|
private |
An object that describes how degrees of freedom should be distributed and renumbered.
Definition at line 936 of file dof_handler.h.
|
private |
A structure that contains all sorts of numbers that characterize the degrees of freedom this object works on.
For most members of this structure, there is an accessor function in this class that returns its value.
Definition at line 945 of file dof_handler.h.
|
private |
Data structure like number_cache, but for each multigrid level.
Definition at line 950 of file dof_handler.h.
|
private |
Array to store the indices for degrees of freedom located at vertices.
Definition at line 1061 of file dof_handler.h.
|
private |
An array to store the indices for level degrees of freedom located at vertices.
Definition at line 1067 of file dof_handler.h.
|
private |
Space to store the DoF numbers for the different levels. Analogous to the levels[]
tree of the Triangulation objects.
Definition at line 1073 of file dof_handler.h.
|
private |
Space to store DoF numbers of faces. They are not stored in levels
since faces are not organized hierarchically, but in a flat array.
Definition at line 1082 of file dof_handler.h.