Flow123d  release_3.0.0-506-g34af125
Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
DOFHandlerMultiDim Class Reference

Provides the numbering of the finite element degrees of freedom on the computational mesh. More...

#include <dofhandler.hh>

Inheritance diagram for DOFHandlerMultiDim:
Inheritance graph
[legend]
Collaboration diagram for DOFHandlerMultiDim:
Collaboration graph
[legend]

Public Member Functions

 DOFHandlerMultiDim (Mesh &_mesh, bool make_elem_part=true)
 Constructor. More...
 
void distribute_dofs (std::shared_ptr< DiscreteSpace > ds, bool sequential=false)
 Distributes degrees of freedom on the mesh needed for the given discrete space. More...
 
unsigned int get_dof_indices (const ElementAccessor< 3 > &cell, std::vector< LongIdx > &indices) const override
 Returns the global indices of dofs associated to the cell. More...
 
unsigned int get_loc_dof_indices (const ElementAccessor< 3 > &cell, std::vector< LongIdx > &indices) const override
 Returns the indices of dofs associated to the cell on the local process. More...
 
int el_index (int loc_el) const
 Returns the global index of local element. More...
 
LongIdx edge_index (int loc_edg) const
 Returns the global index of local edge. More...
 
LongIdx nb_index (int loc_nb) const
 Returns the global index of local neighbour. More...
 
unsigned int n_dofs (ElementAccessor< 3 > cell) const
 Return number of dofs on given cell. More...
 
unsigned int n_loc_edges () const
 Returns number of local edges. More...
 
unsigned int n_loc_nb () const
 Returns number of local neighbours. More...
 
bool el_is_local (int index) const
 
template<unsigned int dim>
FiniteElement< dim > * fe (const ElementAccessor< 3 > &cell) const
 Returns finite element object for given space dimension. More...
 
const Dofcell_dof (ElementAccessor< 3 > cell, unsigned int idof) const
 Return dof on a given cell. More...
 
std::size_t hash () const override
 
 ~DOFHandlerMultiDim () override
 Destructor. More...
 
- Public Member Functions inherited from DOFHandlerBase
 DOFHandlerBase (Mesh &_mesh)
 Constructor. More...
 
const unsigned int n_global_dofs () const
 Getter for the number of all mesh dofs required by the given finite element. More...
 
const unsigned int lsize () const
 Returns the number of dofs on the current process. More...
 
const unsigned int loffset () const
 Returns the offset of the local part of dofs. More...
 
const unsigned int max_elem_dofs () const
 Returns max. number of dofs on one element. More...
 
Distributiondistr () const
 
Meshmesh () const
 Returns the mesh. More...
 
virtual ~DOFHandlerBase ()
 Destructor. More...
 

Private Member Functions

void make_elem_partitioning ()
 Prepare parallel distribution of elements, edges and neighbours. More...
 
void init_cell_starts ()
 Initialize vector of starting indices for elements. More...
 
void init_node_dof_starts (std::vector< LongIdx > &node_dof_starts)
 Initialize auxiliary vector of starting indices of nodal dofs. More...
 
void init_node_status (std::vector< short int > &node_status)
 Initialize node_status. More...
 
void receive_ghost_dofs (unsigned int proc, std::vector< LongIdx > &dofs)
 Obtain dof numbers on ghost elements from other processor. More...
 
void send_ghost_dofs (unsigned int proc)
 Send dof numbers to other processor. More...
 
void update_local_dofs (unsigned int proc, const std::vector< bool > &update_cells, const std::vector< LongIdx > &dofs, const std::vector< LongIdx > &node_dof_starts, std::vector< LongIdx > &node_dofs)
 Update dofs on local elements from ghost element dofs. More...
 
void create_sequential ()
 Communicate local dof indices to all processors. More...
 

Private Attributes

std::shared_ptr< DiscreteSpaceds_
 Pointer to the discrete space for which the handler distributes dofs. More...
 
std::vector< LongIdxcell_starts
 Starting indices for element dofs (parallel version). More...
 
std::vector< LongIdxdof_indices
 Dof numbers on local and ghost elements (parallel version). More...
 
std::vector< LongIdxcell_starts_seq
 Starting indices for element dofs (sequential version). More...
 
std::vector< LongIdxdof_indices_seq
 Dof numbers on mesh elements (sequential version). More...
 
LongIdxrow_4_el
 Global element index -> index according to partitioning. More...
 
LongIdxel_4_loc
 Local element index -> global element index. More...
 
Distributionel_ds_
 Distribution of elements. More...
 
vector< LongIdxedg_4_loc
 Local edge index -> global edge index. More...
 
vector< LongIdxnb_4_loc
 Local neighbour index -> global neighbour index. More...
 
vector< LongIdxnode_4_loc
 Indices of local nodes in mesh tree. More...
 
vector< LongIdxghost_4_loc
 Indices of ghost cells (neighbouring with local elements). More...
 
set< unsigned int > ghost_proc
 Processors of ghost elements. More...
 
map< unsigned int, vector< LongIdx > > ghost_proc_el
 Arrays of ghost cells for each neighbouring processor. More...
 

Static Private Attributes

static const int INVALID_NODE = 1
 
static const int VALID_NODE = 2
 
static const int ASSIGNED_NODE = 3
 
static const int INVALID_DOF = -1
 

Additional Inherited Members

- Protected Attributes inherited from DOFHandlerBase
unsigned int n_global_dofs_
 Number of global dofs assigned by the handler. More...
 
unsigned int lsize_
 Number of dofs associated to local process. More...
 
unsigned int loffset_
 Index of the first dof on the local process. More...
 
unsigned int max_elem_dofs_
 Max. number of dofs per element. More...
 
Meshmesh_
 Pointer to the mesh to which the dof handler is associated. More...
 
Distributiondof_ds_
 Distribution of dofs associated to local process. More...
 

Detailed Description

Provides the numbering of the finite element degrees of freedom on the computational mesh.

Class DOFHandlerMultiDim distributes the degrees of freedom (dof) for a particular triplet of 1d, 2d and 3d finite elements on the computational mesh and provides mappings between local and global dofs. The template parameter dim denotes the spatial dimension of the reference finite element.

Currently the functionality is restricted to finite elements with internal and nodal dofs, i.e. the neighboring elements can share only dofs on nodes.

Definition at line 153 of file dofhandler.hh.

Constructor & Destructor Documentation

DOFHandlerMultiDim::DOFHandlerMultiDim ( Mesh _mesh,
bool  make_elem_part = true 
)

Constructor.

Parameters
_meshThe mesh.
make_elem_partAllow switch off make_element_partitioning, necessary for boundary DOF handler.

Definition at line 238 of file dofhandler.cc.

DOFHandlerMultiDim::~DOFHandlerMultiDim ( )
override

Destructor.

Definition at line 646 of file dofhandler.cc.

Member Function Documentation

const Dof & DOFHandlerMultiDim::cell_dof ( ElementAccessor< 3 >  cell,
unsigned int  idof 
) const

Return dof on a given cell.

Parameters
cellMesh cell.
idofNumber of dof on the cell.

Definition at line 262 of file dofhandler.cc.

Here is the caller graph for this function:

void DOFHandlerMultiDim::create_sequential ( )
private

Communicate local dof indices to all processors.

Definition at line 550 of file dofhandler.cc.

Here is the caller graph for this function:

void DOFHandlerMultiDim::distribute_dofs ( std::shared_ptr< DiscreteSpace ds,
bool  sequential = false 
)

Distributes degrees of freedom on the mesh needed for the given discrete space.

By default, the dof handler is parallel, meaning that each processor has access to dofs on the local elements and on one layer of ghost elements (owned by neighbouring elements). This can be changed by setting sequential to true.

Parameters
dsThe discrete space consisting of finite elements for each mesh element.
sequentialIf true then each processor will have information about all dofs.

Definition at line 451 of file dofhandler.cc.

LongIdx DOFHandlerMultiDim::edge_index ( int  loc_edg) const
inline

Returns the global index of local edge.

Parameters
loc_edgLocal index of edge.

Definition at line 209 of file dofhandler.hh.

int DOFHandlerMultiDim::el_index ( int  loc_el) const
inline

Returns the global index of local element.

Parameters
loc_elLocal index of element.

Definition at line 202 of file dofhandler.hh.

Here is the caller graph for this function:

bool DOFHandlerMultiDim::el_is_local ( int  index) const

Returns true if element is on local process.

Parameters
indexGlobal element index.

Definition at line 719 of file dofhandler.cc.

Here is the caller graph for this function:

template<unsigned int dim>
FiniteElement<dim>* DOFHandlerMultiDim::fe ( const ElementAccessor< 3 > &  cell) const
inline

Returns finite element object for given space dimension.

Parameters
cellCell accessor.

Definition at line 247 of file dofhandler.hh.

unsigned int DOFHandlerMultiDim::get_dof_indices ( const ElementAccessor< 3 > &  cell,
std::vector< LongIdx > &  indices 
) const
overridevirtual

Returns the global indices of dofs associated to the cell.

Parameters
cellThe cell.
indicesArray of dof indices on the cell.

Implements DOFHandlerBase.

Definition at line 605 of file dofhandler.cc.

unsigned int DOFHandlerMultiDim::get_loc_dof_indices ( const ElementAccessor< 3 > &  cell,
std::vector< LongIdx > &  indices 
) const
overridevirtual

Returns the indices of dofs associated to the cell on the local process.

Parameters
cellThe cell.
indicesArray of dof indices on the cell.

Implements DOFHandlerBase.

Definition at line 626 of file dofhandler.cc.

std::size_t DOFHandlerMultiDim::hash ( ) const
overridevirtual

Implements DOFHandlerBase::hash.

Implements DOFHandlerBase.

Definition at line 725 of file dofhandler.cc.

void DOFHandlerMultiDim::init_cell_starts ( )
private

Initialize vector of starting indices for elements.

Definition at line 279 of file dofhandler.cc.

Here is the caller graph for this function:

void DOFHandlerMultiDim::init_node_dof_starts ( std::vector< LongIdx > &  node_dof_starts)
private

Initialize auxiliary vector of starting indices of nodal dofs.

Parameters
node_dof_startsVector of starting indices (output).

Definition at line 300 of file dofhandler.cc.

Here is the caller graph for this function:

void DOFHandlerMultiDim::init_node_status ( std::vector< short int > &  node_status)
private

Initialize node_status.

Set VALID_NODE for nodes owned by local elements and INVALID_NODE for nodes owned by ghost elements.

Parameters
node_statusVector of nodal status (output).

Definition at line 317 of file dofhandler.cc.

Here is the caller graph for this function:

void DOFHandlerMultiDim::make_elem_partitioning ( )
private

Prepare parallel distribution of elements, edges and neighbours.

Definition at line 651 of file dofhandler.cc.

Here is the caller graph for this function:

unsigned int DOFHandlerMultiDim::n_dofs ( ElementAccessor< 3 >  cell) const

Return number of dofs on given cell.

Parameters
cellCell accessor.

Definition at line 246 of file dofhandler.cc.

Here is the caller graph for this function:

unsigned int DOFHandlerMultiDim::n_loc_edges ( ) const
inline

Returns number of local edges.

Definition at line 228 of file dofhandler.hh.

unsigned int DOFHandlerMultiDim::n_loc_nb ( ) const
inline

Returns number of local neighbours.

Definition at line 233 of file dofhandler.hh.

LongIdx DOFHandlerMultiDim::nb_index ( int  loc_nb) const
inline

Returns the global index of local neighbour.

Parameters
loc_nbLocal index of neighbour.

Definition at line 216 of file dofhandler.hh.

void DOFHandlerMultiDim::receive_ghost_dofs ( unsigned int  proc,
std::vector< LongIdx > &  dofs 
)
private

Obtain dof numbers on ghost elements from other processor.

Parameters
procNeighbouring processor.
dofsArray where dofs are stored (output).

Definition at line 347 of file dofhandler.cc.

Here is the caller graph for this function:

void DOFHandlerMultiDim::send_ghost_dofs ( unsigned int  proc)
private

Send dof numbers to other processor.

Parameters
procNeighbouring processor.

Definition at line 368 of file dofhandler.cc.

Here is the caller graph for this function:

void DOFHandlerMultiDim::update_local_dofs ( unsigned int  proc,
const std::vector< bool > &  update_cells,
const std::vector< LongIdx > &  dofs,
const std::vector< LongIdx > &  node_dof_starts,
std::vector< LongIdx > &  node_dofs 
)
private

Update dofs on local elements from ghost element dofs.

Parameters
procNeighbouring processor.
update_cellsVector of global indices of elements which need to be updated from ghost elements.
dofsVector of dof indices on ghost elements from processor proc.
node_dof_startsVector of starting indices of nodal dofs.
node_dofsVector of nodal dof indices (output).

Definition at line 391 of file dofhandler.cc.

Here is the caller graph for this function:

Member Data Documentation

const int DOFHandlerMultiDim::ASSIGNED_NODE = 3
staticprivate

Definition at line 339 of file dofhandler.hh.

std::vector<LongIdx> DOFHandlerMultiDim::cell_starts
private

Starting indices for element dofs (parallel version).

E.g. dof_indices[cell_starts[idx]] = dof number for first dof on the cell with index idx within the paralle structure. To use with element accessor use the following:

ElementAccessor<3> cell; ... // i-th dof number on the cell dof_indices[cell_starts[row_4_el[cell.idx()]]+i] = ...

Only local and ghost elements are stored, but the vector has size mesh_->n_elements()+1.

Definition at line 360 of file dofhandler.hh.

std::vector<LongIdx> DOFHandlerMultiDim::cell_starts_seq
private

Starting indices for element dofs (sequential version).

This vector stores information about all mesh elements. See cell_starts for paralle version.

Definition at line 376 of file dofhandler.hh.

std::vector<LongIdx> DOFHandlerMultiDim::dof_indices
private

Dof numbers on local and ghost elements (parallel version).

Dofs are ordered accordingly with cell_starts and local dof order given by the finite element. See cell_starts for more description.

Definition at line 368 of file dofhandler.hh.

std::vector<LongIdx> DOFHandlerMultiDim::dof_indices_seq
private

Dof numbers on mesh elements (sequential version).

This vector stores information about all mesh elements. See dof_indices for paralle version.

Definition at line 384 of file dofhandler.hh.

std::shared_ptr<DiscreteSpace> DOFHandlerMultiDim::ds_
private

Pointer to the discrete space for which the handler distributes dofs.

Definition at line 344 of file dofhandler.hh.

vector<LongIdx> DOFHandlerMultiDim::edg_4_loc
private

Local edge index -> global edge index.

Definition at line 397 of file dofhandler.hh.

LongIdx* DOFHandlerMultiDim::el_4_loc
private

Local element index -> global element index.

Definition at line 391 of file dofhandler.hh.

Distribution* DOFHandlerMultiDim::el_ds_
private

Distribution of elements.

Definition at line 394 of file dofhandler.hh.

vector<LongIdx> DOFHandlerMultiDim::ghost_4_loc
private

Indices of ghost cells (neighbouring with local elements).

Definition at line 406 of file dofhandler.hh.

set<unsigned int> DOFHandlerMultiDim::ghost_proc
private

Processors of ghost elements.

Definition at line 409 of file dofhandler.hh.

map<unsigned int, vector<LongIdx> > DOFHandlerMultiDim::ghost_proc_el
private

Arrays of ghost cells for each neighbouring processor.

Definition at line 412 of file dofhandler.hh.

const int DOFHandlerMultiDim::INVALID_DOF = -1
staticprivate

Definition at line 340 of file dofhandler.hh.

const int DOFHandlerMultiDim::INVALID_NODE = 1
staticprivate

Flags used during distribution of dofs to mark node and dof status.

Definition at line 337 of file dofhandler.hh.

vector<LongIdx> DOFHandlerMultiDim::nb_4_loc
private

Local neighbour index -> global neighbour index.

Definition at line 400 of file dofhandler.hh.

vector<LongIdx> DOFHandlerMultiDim::node_4_loc
private

Indices of local nodes in mesh tree.

Definition at line 403 of file dofhandler.hh.

LongIdx* DOFHandlerMultiDim::row_4_el
private

Global element index -> index according to partitioning.

Definition at line 388 of file dofhandler.hh.

const int DOFHandlerMultiDim::VALID_NODE = 2
staticprivate

Definition at line 338 of file dofhandler.hh.


The documentation for this class was generated from the following files: