Flow123d
JS_before_hm-2207-gef9ee6d82
|
Go to the documentation of this file.
54 scatter_to_seq_(nullptr),
97 unsigned int n_node_dofs = 0;
100 node_dof_starts.push_back(n_node_dofs);
101 n_node_dofs +=
ds_->n_node_dofs(nid);
103 node_dof_starts.push_back(n_node_dofs);
105 unsigned int n_edge_dofs = 0;
108 edge_dof_starts.push_back(n_edge_dofs);
109 n_edge_dofs +=
ds_->n_edge_dofs(
edge);
111 edge_dof_starts.push_back(n_edge_dofs);
122 for (
unsigned int n=0; n<cell.dim()+1; n++)
134 for (
unsigned int n=0; n<cell.dim()+1; n++)
151 for (
unsigned int n=0; n<cell.dim()+1; n++)
153 unsigned int eid = cell.elm().side(n)->edge_idx();
175 unsigned int n_dofs_sum = 0;
176 for (
auto nd : n_dofs) n_dofs_sum += nd;
177 dofs.resize(n_dofs_sum);
185 unsigned int n_elems;
222 unsigned int dof_offset=0;
223 for (
unsigned int gid=0; gid<
ghost_proc_el[proc].size(); gid++)
229 for (
unsigned int idof = 0; idof<dh_cell.
n_dofs(); ++idof)
235 unsigned int node_dof_idx = node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx];
244 loc_node_dof_count[dof_nface_idx]++;
250 unsigned int edge_dof_idx = edge_dof_starts[eid]+loc_edge_dof_count[dof_nface_idx];
259 loc_edge_dof_count[dof_nface_idx]++;
267 dof_offset += dh_cell.
n_dofs();
273 if (!update_cells[cell.local_idx()])
continue;
278 for (
unsigned int idof = 0; idof<cell.n_dofs(); ++idof)
280 unsigned int dof_nface_idx = cell.cell_dof(idof).n_face_idx;
281 if (cell.cell_dof(idof).dim == 0)
286 dof_indices[
cell_starts[cell.local_idx()]+idof] = node_dofs[node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx]];
288 loc_node_dof_count[dof_nface_idx]++;
289 }
else if (cell.cell_dof(idof).dim == cell.dim()-1)
293 unsigned int eid = cell.elm().side(dof_nface_idx)->edge_idx();
294 dof_indices[
cell_starts[cell.local_idx()]+idof] = edge_dofs[edge_dof_starts[eid]+loc_edge_dof_count[dof_nface_idx]];
296 loc_edge_dof_count[dof_nface_idx]++;
314 unsigned int next_free_dof = 0;
331 for (
unsigned int idof = 0; idof<cell.n_dofs(); ++idof)
333 unsigned int dof_dim = cell.cell_dof(idof).dim;
334 unsigned int dof_nface_idx = cell.cell_dof(idof).n_face_idx;
339 unsigned int node_dof_idx = node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx];
341 switch (node_status[nid])
344 for (
int i=0; i<node_dof_starts[nid+1] - node_dof_starts[nid]; i++)
347 node_dofs[node_dof_starts[nid]+i] = next_free_dof++;
353 update_cells[cell.local_idx()] =
true;
357 loc_node_dof_count[dof_nface_idx]++;
359 else if (dof_dim == cell.dim()-1)
361 unsigned int eid = cell.elm().side(dof_nface_idx)->edge_idx();
362 unsigned int edge_dof_idx = edge_dof_starts[eid]+loc_edge_dof_count[dof_nface_idx];
363 switch (edge_status[eid])
366 for (
int i=0; i<edge_dof_starts[eid+1] - edge_dof_starts[eid]; i++)
369 edge_dofs[edge_dof_starts[eid]+i] = next_free_dof++;
375 update_cells[cell.local_idx()] =
true;
379 loc_edge_dof_count[dof_nface_idx]++;
381 else if (dof_dim == cell.dim())
396 dof_ds_ = std::make_shared<Distribution>(
lsize_, PETSC_COMM_WORLD);
409 for (
unsigned int from_higher = 0; from_higher < 2; from_higher++)
413 if ((proc >
el_ds_->
myp()) == from_higher)
433 update_cells.clear();
435 node_dof_starts.clear();
437 edge_dof_starts.clear();
443 if (
dh_seq_ !=
nullptr)
return;
447 dh_seq_ = std::make_shared<DOFHandlerMultiDim>(*
this);
484 cell_starts_loc[i+1] += cell_starts_loc[i];
490 for (
unsigned int idof=0; idof<cell.n_dofs(); idof++)
500 cell_starts_loc.size(),
506 dof_indices_loc.size(),
509 (
const int *)
distr.get_lsizes_array(),
510 (
const int *)
distr.get_starts_array(),
516 VecCreateMPI(PETSC_COMM_WORLD,
lsize_, PETSC_DETERMINE, &v_from);
524 dh_seq_->scatter_to_seq_ = std::make_shared<VecScatter>();
525 VecScatterCreateToAll(v_seq,
dh_seq_->scatter_to_seq_.get(), NULL);
546 unsigned int ndofs = 0;
548 indices.resize(ndofs);
549 for (
unsigned int k=0; k<ndofs; k++)
569 bool is_edge_local =
false;
570 for (
uint sid=0; sid<
edge.n_sides(); sid++)
573 is_edge_local =
true;
590 for (
unsigned int iel = 0; iel <
el_ds_->
lsize(); iel++ )
598 for (
unsigned int nid=0; nid<cell.elm()->n_nodes(); nid++)
607 bool has_local_node =
false;
609 for (
unsigned int nid=0; nid<cell->n_nodes(); nid++)
612 has_local_node =
true;
658 auto bgn_it = make_iter<DHCellAccessor>(
DHCellAccessor(
this, 0) );
665 auto bgn_it = make_iter<DHCellAccessor>(
DHCellAccessor(
this, 0) );
689 s <<
"DOFHandlerMultiDim structure:" << endl;
690 s <<
"- is parallel: " << (
is_parallel_?
"true":
"false") << endl;
691 s <<
"- proc id: " <<
el_ds_->
myp() << endl;
693 s <<
"- number of locally owned cells: " <<
el_ds_->
lsize() << endl;
694 s <<
"- number of ghost cells: " <<
ghost_4_loc.size() << endl;
695 s <<
"- dofs on locally owned cells:" << endl;
699 auto ndofs = cell.get_dof_indices(dofs);
700 s <<
"-- cell " << cell.elm().idx() <<
": ";
701 for (
unsigned int idof=0; idof<ndofs; idof++) s << dofs[idof] <<
" "; s << endl;
703 s <<
"- dofs on ghost cells:" << endl;
706 auto ndofs = cell.get_dof_indices(dofs);
707 s <<
"-- cell " << cell.elm().idx() <<
": ";
708 for (
unsigned int idof=0; idof<ndofs; idof++) s << dofs[idof] <<
" "; s << endl;
710 s <<
"- locally owned dofs (" <<
lsize_ <<
"): ";
714 s <<
"- global-to-local-cell map:" << endl;
718 printf(
"%s", s.str().c_str());
732 fe_idx_(component_idx)
736 .error(
"sub_handler can be used only with dof handler using EqualOrderDiscreteSpace!");
741 ASSERT( fe_sys0 !=
nullptr ).error(
"sub_handler assumes that dof handler uses FESystem<0>!");
742 ASSERT( fe_sys1 !=
nullptr ).error(
"sub_handler assumes that dof handler uses FESystem<1>!");
743 ASSERT( fe_sys2 !=
nullptr ).error(
"sub_handler assumes that dof handler uses FESystem<2>!");
744 ASSERT( fe_sys3 !=
nullptr ).error(
"sub_handler assumes that dof handler uses FESystem<3>!");
746 fe_sys1->
fe()[component_idx],
747 fe_sys2->
fe()[component_idx],
748 fe_sys3->
fe()[component_idx]);
749 ds_ = std::make_shared<EqualOrderDiscreteSpace>(
mesh_, sub_mixed);
759 for (
unsigned int d=0; d<=3; d++)
760 ASSERT( fs[d] !=
nullptr ).error(
"Function space must be of type FESystemFunctionSpace!" );
762 for (
unsigned int i=0; i<fe_sys0->
n_dofs(); i++)
763 if (fs[0]->
dof_indices()[i].fe_index == component_idx) sub_fe_dofs[0].push_back(i);
764 for (
unsigned int i=0; i<fe_sys1->
n_dofs(); i++)
765 if (fs[1]->
dof_indices()[i].fe_index == component_idx) sub_fe_dofs[1].push_back(i);
766 for (
unsigned int i=0; i<fe_sys2->
n_dofs(); i++)
767 if (fs[2]->
dof_indices()[i].fe_index == component_idx) sub_fe_dofs[2].push_back(i);
768 for (
unsigned int i=0; i<fe_sys3->
n_dofs(); i++)
769 if (fs[3]->
dof_indices()[i].fe_index == component_idx) sub_fe_dofs[3].push_back(i);
777 for (
auto cell : dh->local_range())
779 LocDofVec cell_dof_indices = cell.get_loc_dof_indices();
780 for (
auto sub_dof : sub_fe_dofs[cell.dim()])
782 if (cell_dof_indices[sub_dof] <
static_cast<int>(dh->lsize_) &&
783 sub_local_indices[cell_dof_indices[sub_dof]] ==
INVALID_DOF)
794 for (
auto cell : dh->local_range())
796 LocDofVec cell_dof_indices = cell.get_loc_dof_indices();
797 unsigned int idof = 0;
798 for (
auto sub_dof : sub_fe_dofs[cell.dim()])
800 if (sub_local_indices[cell_dof_indices[sub_dof]] ==
INVALID_DOF)
811 dof_ds_ = std::make_shared<Distribution>(
lsize_, PETSC_COMM_WORLD);
817 for (
unsigned int i=0; i<
lsize_; i++)
822 for (
unsigned int from_higher = 0; from_higher < 2; from_higher++)
826 if ((proc >
el_ds_->
myp()) == from_higher)
852 dofs.resize(n_ghosts);
856 unsigned int idof = 0;
866 unsigned int n_ghosts;
876 dofs.push_back(global_to_local_dof_idx.at(global_dof) +
dof_ds_->begin());
883 ASSERT(
vec.size() ==
parent_->local_to_global_dof_idx_.size() ).error(
"Incompatible parent vector in update_subvector()!");
893 ASSERT(
vec.size() ==
parent_->local_to_global_dof_idx_.size() ).error(
"Incompatible parent vector in update_subvector()!");
LongIdx * get_row_4_el() const
Range< DHCellAccessor > ghost_range() const
Returns range over ghosts DOF handler cells.
arma::Col< IntIdx > LocDofVec
std::vector< LongIdx > local_to_global_dof_idx_
Maps local and ghost dof indices to global ones.
unsigned int np() const
get num of processors
static const int INVALID_DOF
unsigned int n_face_idx
Index of n-face to which the dof is associated.
std::unordered_map< LongIdx, LongIdx > global_to_local_el_idx_
Maps global element index into local/ghost index (obsolete).
set< unsigned int > ghost_proc
Processors of ghost elements.
unsigned int lsize(int proc) const
get local size
std::shared_ptr< DiscreteSpace > ds() const
Return pointer to discrete space for which the handler distributes dofs.
unsigned int myp() const
get my processor
friend class SubDOFHandlerMultiDim
Distribution * get_el_ds() const
Support classes for parallel programing.
const std::vector< unsigned int > & obj_4_el() const
MeshBase * mesh_
Pointer to the mesh to which the dof handler is associated.
void update_subvector(const VectorMPI &vec, VectorMPI &subvec)
Update values in subvector from parent vector.
unsigned int n_edges() const
vector< LongIdx > ghost_4_loc
Indices of ghost cells (neighbouring with local elements).
#define MPI_Send(buf, count, datatype, dest, tag, comm)
virtual ~DOFHandlerBase()
Destructor.
unsigned int dim() const
Return dimension of element appropriate to cell.
Distribution * el_ds_
Distribution of elements.
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
Range< DHCellAccessor > local_range() const
Returns range over own and ghost cells of DOF handler.
map< unsigned int, vector< LongIdx > > ghost_proc_el
Arrays of ghost cells for each neighbouring processor.
Range< ElementAccessor< 3 > > elements_range() const
Returns range of mesh elements.
void receive_sub_ghost_dofs(unsigned int proc, vector< LongIdx > &dofs)
Get global dof indices of ghost dofs for sub-handler.
Class FESystem for compound finite elements.
unsigned int dim
Association to n-face of given dimension (point, line, triangle, tetrahedron.
ElementAccessor< 3 > element() const
const Neighbour & vb_neighbour(unsigned int nb) const
Return neighbour with given index.
void set(unsigned int pos, double val)
Set value on given position.
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
#define MPI_Recv(buf, count, datatype, source, tag, comm, status)
unsigned int edge_idx() const
Returns global index of the edge connected to the side.
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, const std::vector< LongIdx > &edge_dof_starts, std::vector< LongIdx > &edge_dofs)
Update dofs on local elements from ghost element dofs.
std::shared_ptr< DOFHandlerMultiDim > dh_seq_
Sequential dof handler associated to the current (parallel) one.
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
bool el_is_local(int index) const
Provides the numbering of the finite element degrees of freedom on the computational mesh.
void update_parent_vector(VectorMPI &vec, const VectorMPI &subvec)
Update values in parent vector from values of subvector.
unsigned int local_idx() const
Return local index to element (index of DOF handler).
static const int VALID_NFACE
Range< DHCellAccessor > own_range() const
Returns range of DOF handler cells (only range of own without ghost cells)
Abstract class for description of finite elements.
static const int INVALID_NFACE
void distribute_dofs(std::shared_ptr< DiscreteSpace > ds)
Distributes degrees of freedom on the mesh needed for the given discrete space.
void receive_ghost_dofs(unsigned int proc, std::vector< LongIdx > &dofs)
Obtain dof numbers on ghost elements from other processor.
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
unsigned int n_nodes() const
virtual VectorMPI create_vector()
Allocates PETSc vector according to the dof distribution.
unsigned int n_global_dofs_
Number of global dofs assigned by the handler.
unsigned int loffset_
Index of the first dof on the local process.
std::shared_ptr< DOFHandlerMultiDim > parent_
Parent dof handler.
std::shared_ptr< DOFHandlerMultiDim > sequential()
Returns sequential version of the current dof handler.
ElementAccessor< 3 > element() const
Returns iterator to the element of the side.
unsigned int n_dofs() const
Return number of dofs on given cell.
void init_status(std::vector< short int > &node_status, std::vector< short int > &edge_status)
Initialize node_status and edge_status.
std::shared_ptr< Distribution > distr() const
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
friend class DHCellAccessor
unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
void init_dof_starts(std::vector< LongIdx > &node_dof_starts, std::vector< LongIdx > &edge_dof_starts)
Initialize auxiliary vector of starting indices of nodal/edge dofs.
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
const DHCellAccessor cell_accessor_from_element(unsigned int elm_idx) const
Return DHCellAccessor appropriate to ElementAccessor of given idx.
unsigned int n_sides() const
#define MPI_STATUS_IGNORE
~DOFHandlerMultiDim() override
Destructor.
vector< LongIdx > nb_4_loc
Local neighbour index -> global neighbour index.
std::shared_ptr< DiscreteSpace > ds_
Pointer to the discrete space for which the handler distributes dofs.
Cell accessor allow iterate over DOF handler cells.
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
void printf(BasicWriter< Char > &w, BasicCStringRef< Char > format, ArgList args)
LongIdx * get_el_4_loc() const
std::vector< LongIdx > cell_starts
Starting indices for local (owned+ghost) element dofs.
void print() const
Output structure of dof handler.
const std::vector< MeshObject > & objects(unsigned int dim) const
void make_elem_partitioning()
Prepare parallel distribution of elements, edges and neighbours.
Base class for Mesh and BCMesh.
std::vector< IntIdx > dof_indices
Dof numbers on local and ghost elements.
unsigned int max_elem_dofs_
Max. number of dofs per element.
void create_sequential()
Communicate local dof indices to all processors and create new sequential dof handler.
unsigned int lsize_
Number of dofs associated to local process.
DOFHandlerMultiDim(MeshBase &_mesh, bool make_elem_part=true)
Constructor.
void send_sub_ghost_dofs(unsigned int proc, const map< LongIdx, LongIdx > &global_to_local_dof_idx)
Send global indices of dofs that are ghost on other processors.
unsigned int n_vb_neighbours() const
std::vector< LongIdx > parent_dof_idx_
Local indices in the parent handler.
std::shared_ptr< VecScatter > scatter_to_seq_
Scatter context for parallel to sequential vectors.
Range< Edge > edge_range() const
Return range of edges.
unsigned int idx() const
We need this method after replacing Region by RegionIdx, and movinf RegionDB instance into particular...
std::shared_ptr< FunctionSpace > function_space_
Function space defining the FE.
std::shared_ptr< VecScatter > sequential_scatter()
Returns scatter context from parallel to sequential vectors.
#define MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm)
double get(unsigned int pos) const
Return value on given position.
bool is_parallel_
Indicator for parallel/sequential dof handler.
unsigned int n_elements() const
static const int ASSIGNED_NFACE
std::size_t hash() const override
void send_ghost_dofs(unsigned int proc)
Send dof numbers to other processor.
SideIter side(const unsigned int loc_index)
const Dof & cell_dof(unsigned int idof) const
Return dof on a given cell.
Compound finite element on dim dimensional simplex.
bool is_local(unsigned int idx) const
identify local index
unsigned int get_dof_indices(const DHCellAccessor &cell, std::vector< LongIdx > &indices) const override
Returns the global indices of dofs associated to the cell.
const DuplicateNodes * duplicate_nodes() const
unsigned int size() const
Return size of output data.
std::shared_ptr< Distribution > dof_ds_
Distribution of dofs associated to local process.
void init_cell_starts()
Initialize vector of starting indices for elements.
unsigned int n_nodes() const
Implementation of range helper class.
vector< LongIdx > edg_4_loc
Local edge index -> global edge index.