Flow123d  JS_before_hm-1626-gde32303
dofhandler.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file dofhandler.cc
15  * @brief Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between local and global dofs.
16  * @author Jan Stebel
17  */
18 
19 #include "system/index_types.hh"
20 #include "fem/dofhandler.hh"
21 #include "fem/finite_element.hh"
22 #include "fem/fe_system.hh"
23 #include "fem/dh_cell_accessor.hh"
24 #include "mesh/mesh.h"
25 #include "mesh/duplicate_nodes.h"
26 #include "mesh/partitioning.hh"
27 #include "mesh/accessors.hh"
28 #include "mesh/range_wrapper.hh"
29 #include "mesh/neighbours.h"
30 #include "la/distribution.hh"
31 
32 
36 const int DOFHandlerMultiDim::INVALID_DOF = -1;
37 
38 
39 
40 
42 {}
43 
44 
45 
46 
47 
48 
49 DOFHandlerMultiDim::DOFHandlerMultiDim(Mesh& _mesh, bool make_elem_part)
50  : DOFHandlerBase(_mesh),
51  ds_(nullptr),
52  is_parallel_(true),
53  dh_seq_(nullptr),
54  scatter_to_seq_(nullptr),
55  el_ds_(nullptr)
56 {
57  if (make_elem_part) make_elem_partitioning();
58 }
59 
60 
61 std::shared_ptr<DOFHandlerMultiDim> DOFHandlerMultiDim::sequential()
62 {
64  return dh_seq_;
65 }
66 
67 
68 std::shared_ptr<VecScatter> DOFHandlerMultiDim::sequential_scatter()
69 {
71  return scatter_to_seq_;
72 }
73 
74 
76 {
77  // get number of dofs per element and then set up cell_starts
79  for (auto cell : this->local_range())
80  {
81  cell_starts[cell.local_idx()+1] = cell.n_dofs();
82  max_elem_dofs_ = max( (int)max_elem_dofs_, (int)cell.n_dofs() );
83  }
84  for (unsigned int i=0; i<cell_starts.size()-1; ++i)
85  cell_starts[i+1] += cell_starts[i];
86 }
87 
88 
90  std::vector<LongIdx> &node_dof_starts,
91  std::vector<LongIdx> &edge_dof_starts)
92 {
93  // initialize dofs on nodes
94  // We must separate dofs for dimensions because FE functions
95  // may be discontinuous on nodes shared by different
96  // dimensions.
97  unsigned int n_node_dofs = 0;
98  for (unsigned int nid=0; nid<mesh_->tree->n_nodes(); nid++)
99  {
100  node_dof_starts.push_back(n_node_dofs);
101  n_node_dofs += ds_->n_node_dofs(nid);
102  }
103  node_dof_starts.push_back(n_node_dofs);
104 
105  unsigned int n_edge_dofs = 0;
106  for (auto edge : mesh_->edge_range())
107  {
108  edge_dof_starts.push_back(n_edge_dofs);
109  n_edge_dofs += ds_->n_edge_dofs(edge);
110  }
111  edge_dof_starts.push_back(n_edge_dofs);
112 }
113 
114 
116  std::vector<short int> &node_status,
117  std::vector<short int> &edge_status)
118 {
119  // mark local dofs
120  for (auto cell : this->own_range())
121  {
122  for (unsigned int n=0; n<cell.dim()+1; n++)
123  {
124  unsigned int nid = mesh_->tree->objects(cell.dim())[mesh_->tree->obj_4_el()[cell.elm_idx()]].nodes[n];
125  node_status[nid] = VALID_NFACE;
126  }
127  }
128 
129  // unmark dofs on ghost cells from lower procs
130  for (auto cell : this->ghost_range())
131  {
132  if (cell.elm().proc() < el_ds_->myp())
133  {
134  for (unsigned int n=0; n<cell.dim()+1; n++)
135  {
136  unsigned int nid = mesh_->tree->objects(cell.dim())[mesh_->tree->obj_4_el()[cell.elm_idx()]].nodes[n];
137  node_status[nid] = INVALID_NFACE;
138  }
139  }
140  }
141 
142  // mark local edges
143  for (auto eid : edg_4_loc)
144  edge_status[eid] = VALID_NFACE;
145 
146  // unmark dofs on ghost cells from lower procs
147  for (auto cell : this->ghost_range())
148  {
149  if (cell.elm().proc() < el_ds_->myp())
150  {
151  for (unsigned int n=0; n<cell.dim()+1; n++)
152  {
153  unsigned int eid = cell.elm().side(n)->edge_idx();
154  edge_status[eid] = INVALID_NFACE;
155  }
156  }
157  }
158 }
159 
160 
162 {
163  // send number of elements required from the other processor
164  unsigned int n_elems = ghost_proc_el[proc].size();
165  MPI_Send(&n_elems, 1, MPI_UNSIGNED, proc, 0, MPI_COMM_WORLD);
166 
167  // send indices of elements required
168  MPI_Send(&(ghost_proc_el[proc][0]), n_elems, MPI_LONG_IDX, proc, 1, MPI_COMM_WORLD);
169 
170  // receive numbers of dofs on required elements
171  vector<unsigned int> n_dofs(n_elems);
172  MPI_Recv(&(n_dofs[0]), n_elems, MPI_UNSIGNED, proc, 2, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
173 
174  // receive dofs on required elements
175  unsigned int n_dofs_sum = 0;
176  for (auto nd : n_dofs) n_dofs_sum += nd;
177  dofs.resize(n_dofs_sum);
178  MPI_Recv(&(dofs[0]), n_dofs_sum, MPI_LONG_IDX, proc, 3, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
179 }
180 
181 
182 void DOFHandlerMultiDim::send_ghost_dofs(unsigned int proc)
183 {
184  // receive number of elements required by the other processor
185  unsigned int n_elems;
186  MPI_Recv(&n_elems, 1, MPI_UNSIGNED, proc, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
187 
188  // receive indices of elements required
189  vector<LongIdx> elems(n_elems);
190  MPI_Recv(&(elems[0]), n_elems, MPI_LONG_IDX, proc, 1, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
191 
192  // send numbers of dofs on required elements
193  vector<unsigned int> n_dofs;
194  for (LongIdx el : elems)
195  {
196  auto cell = this->cell_accessor_from_element(el);
197  n_dofs.push_back(cell_starts[cell.local_idx()+1] - cell_starts[cell.local_idx()]);
198  }
199  MPI_Send(&(n_dofs[0]), n_elems, MPI_UNSIGNED, proc, 2, MPI_COMM_WORLD);
200 
201  // send dofs on the required elements
202  vector<LongIdx> dofs;
203  for (LongIdx el : elems)
204  {
205  auto cell = this->cell_accessor_from_element(el);
206  for (LongIdx i=cell_starts[cell.local_idx()]; i<cell_starts[cell.local_idx()+1]; i++)
207  dofs.push_back(local_to_global_dof_idx_[dof_indices[i]]);
208  }
209  MPI_Send(&(dofs[0]), dofs.size(), MPI_LONG_IDX, proc, 3, MPI_COMM_WORLD);
210 }
211 
212 
214  const std::vector<bool> &update_cells,
215  const std::vector<LongIdx> &dofs,
216  const std::vector<LongIdx> &node_dof_starts,
217  std::vector<LongIdx> &node_dofs,
218  const std::vector<LongIdx> &edge_dof_starts,
219  std::vector<LongIdx> &edge_dofs)
220 {
221  // update dof_indices on ghost cells
222  unsigned int dof_offset=0;
223  for (unsigned int gid=0; gid<ghost_proc_el[proc].size(); gid++)
224  {
225  DHCellAccessor dh_cell = this->cell_accessor_from_element( ghost_proc_el[proc][gid] );
226 
227  vector<unsigned int> loc_node_dof_count(dh_cell.elm()->n_nodes(), 0);
228  vector<unsigned int> loc_edge_dof_count(dh_cell.elm()->n_sides(), 0);
229  for (unsigned int idof = 0; idof<dh_cell.n_dofs(); ++idof)
230  {
231  if (dh_cell.cell_dof(idof).dim == 0)
232  { // update nodal dof
233  unsigned int dof_nface_idx = dh_cell.cell_dof(idof).n_face_idx;
234  unsigned int nid = mesh_->tree->objects(dh_cell.dim())[mesh_->tree->obj_4_el()[dh_cell.elm_idx()]].nodes[dof_nface_idx];
235  unsigned int node_dof_idx = node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx];
236 
237  if (node_dofs[node_dof_idx] == INVALID_DOF)
238  {
239  node_dofs[node_dof_idx] = local_to_global_dof_idx_.size();
240  local_to_global_dof_idx_.push_back(dofs[dof_offset+idof]);
241  }
242  dof_indices[cell_starts[dh_cell.local_idx()]+idof] = node_dofs[node_dof_idx];
243 
244  loc_node_dof_count[dof_nface_idx]++;
245  }
246  else if (dh_cell.cell_dof(idof).dim == dh_cell.dim()-1)
247  { // update edge dof
248  unsigned int dof_nface_idx = dh_cell.cell_dof(idof).n_face_idx;
249  unsigned int eid = dh_cell.elm().side(dof_nface_idx)->edge_idx();
250  unsigned int edge_dof_idx = edge_dof_starts[eid]+loc_edge_dof_count[dof_nface_idx];
251 
252  if (edge_dofs[edge_dof_idx] == INVALID_DOF)
253  {
254  edge_dofs[edge_dof_idx] = local_to_global_dof_idx_.size();
255  local_to_global_dof_idx_.push_back(dofs[dof_offset+idof]);
256  }
257  dof_indices[cell_starts[dh_cell.local_idx()]+idof] = edge_dofs[edge_dof_idx];
258 
259  loc_edge_dof_count[dof_nface_idx]++;
260  } else if (dh_cell.cell_dof(idof).dim == dh_cell.dim())
261  {
262  dof_indices[cell_starts[dh_cell.local_idx()]+idof] = local_to_global_dof_idx_.size();
263  local_to_global_dof_idx_.push_back(dofs[dof_offset+idof]);
264  }
265  }
266 
267  dof_offset += dh_cell.n_dofs();
268  }
269 
270  // update dof_indices on local elements
271  for (auto cell : this->own_range())
272  {
273  if (!update_cells[cell.local_idx()]) continue;
274 
275  // loop over element dofs
276  vector<unsigned int> loc_node_dof_count(cell.elm()->n_nodes(), 0);
277  vector<unsigned int> loc_edge_dof_count(cell.elm()->n_sides(), 0);
278  for (unsigned int idof = 0; idof<cell.n_dofs(); ++idof)
279  {
280  unsigned int dof_nface_idx = cell.cell_dof(idof).n_face_idx;
281  if (cell.cell_dof(idof).dim == 0)
282  {
283  if (dof_indices[cell_starts[cell.local_idx()]+idof] == INVALID_DOF)
284  { // update nodal dof
285  unsigned int nid = mesh_->tree->objects(cell.dim())[mesh_->tree->obj_4_el()[cell.elm_idx()]].nodes[dof_nface_idx];
286  dof_indices[cell_starts[cell.local_idx()]+idof] = node_dofs[node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx]];
287  }
288  loc_node_dof_count[dof_nface_idx]++;
289  } else if (cell.cell_dof(idof).dim == cell.dim()-1)
290  {
291  if (dof_indices[cell_starts[cell.local_idx()]+idof] == INVALID_DOF)
292  { // update edge dof
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]];
295  }
296  loc_edge_dof_count[dof_nface_idx]++;
297  }
298  }
299  }
300 }
301 
302 
303 void DOFHandlerMultiDim::distribute_dofs(std::shared_ptr<DiscreteSpace> ds)
304 {
305  // First check if dofs are already distributed.
306  OLD_ASSERT(ds_ == nullptr, "Attempt to distribute DOFs multiple times!");
307 
308  ds_ = ds;
309 
310  std::vector<LongIdx> node_dofs, node_dof_starts, edge_dofs, edge_dof_starts;
312  edge_status(mesh_->n_edges(), INVALID_NFACE);
313  std::vector<bool> update_cells(el_ds_->lsize(), false);
314  unsigned int next_free_dof = 0;
315 
317  init_dof_starts(node_dof_starts, edge_dof_starts);
318  node_dofs.resize(node_dof_starts[node_dof_starts.size()-1], (LongIdx)INVALID_DOF);
319  edge_dofs.resize(edge_dof_starts[edge_dof_starts.size()-1], (LongIdx)INVALID_DOF);
320  init_status(node_status, edge_status);
321 
322  // Distribute dofs on local elements.
323  dof_indices.resize(cell_starts[cell_starts.size()-1]);
324  local_to_global_dof_idx_.reserve(dof_indices.size());
325  for (auto cell : this->own_range())
326  {
327 
328  // loop over element dofs
329  vector<unsigned int> loc_node_dof_count(cell.elm()->n_nodes(), 0);
330  vector<unsigned int> loc_edge_dof_count(cell.elm()->dim()+1, 0);
331  for (unsigned int idof = 0; idof<cell.n_dofs(); ++idof)
332  {
333  unsigned int dof_dim = cell.cell_dof(idof).dim;
334  unsigned int dof_nface_idx = cell.cell_dof(idof).n_face_idx;
335 
336  if (dof_dim == 0)
337  { // add dofs shared by nodes
338  unsigned int nid = mesh_->tree->objects(cell.dim())[mesh_->tree->obj_4_el()[cell.elm_idx()]].nodes[dof_nface_idx];
339  unsigned int node_dof_idx = node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx];
340 
341  switch (node_status[nid])
342  {
343  case VALID_NFACE:
344  for (int i=0; i<node_dof_starts[nid+1] - node_dof_starts[nid]; i++)
345  {
346  local_to_global_dof_idx_.push_back(next_free_dof);
347  node_dofs[node_dof_starts[nid]+i] = next_free_dof++;
348  }
349  node_status[nid] = ASSIGNED_NFACE;
350  break;
351  case INVALID_NFACE:
352  node_dofs[node_dof_idx] = INVALID_DOF;
353  update_cells[cell.local_idx()] = true;
354  break;
355  }
356  dof_indices[cell_starts[cell.local_idx()]+idof] = node_dofs[node_dof_idx];
357  loc_node_dof_count[dof_nface_idx]++;
358  }
359  else if (dof_dim == cell.dim()-1)
360  { // add dofs shared by edges
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])
364  {
365  case VALID_NFACE:
366  for (int i=0; i<edge_dof_starts[eid+1] - edge_dof_starts[eid]; i++)
367  {
368  local_to_global_dof_idx_.push_back(next_free_dof);
369  edge_dofs[edge_dof_starts[eid]+i] = next_free_dof++;
370  }
371  edge_status[eid] = ASSIGNED_NFACE;
372  break;
373  case INVALID_NFACE:
374  edge_dofs[edge_dof_idx] = INVALID_DOF;
375  update_cells[cell.local_idx()] = true;
376  break;
377  }
378  dof_indices[cell_starts[cell.local_idx()]+idof] = edge_dofs[edge_dof_idx];
379  loc_edge_dof_count[dof_nface_idx]++;
380  }
381  else if (dof_dim == cell.dim())
382  { // add dofs owned only by the element
383  local_to_global_dof_idx_.push_back(next_free_dof);
384  dof_indices[cell_starts[cell.local_idx()]+idof] = next_free_dof++;
385  }
386  else
387  ASSERT(false).error("Unsupported dof n_face.");
388  }
389  }
390  node_status.clear();
391  edge_status.clear();
392 
393  lsize_ = next_free_dof;
394 
395  // communicate n_dofs across all processes
396  dof_ds_ = std::make_shared<Distribution>(lsize_, PETSC_COMM_WORLD);
397  n_global_dofs_ = dof_ds_->size();
398 
399  // shift dof indices
400  loffset_ = dof_ds_->get_starts_array()[dof_ds_->myp()];
401  if (loffset_ > 0)
402  {
403  for (auto &i : local_to_global_dof_idx_)
404  i += loffset_;
405  }
406 
407  // communicate dofs from ghost cells
408  // first propagate from lower procs to higher procs and then vice versa
409  for (unsigned int from_higher = 0; from_higher < 2; from_higher++)
410  {
411  for (unsigned int proc : ghost_proc)
412  {
413  if ((proc > el_ds_->myp()) == from_higher)
414  { // receive dofs from master processor
415  vector<LongIdx> dofs;
416  receive_ghost_dofs(proc, dofs);
417 
418  // update dof_indices and node_dofs on ghost elements
419  update_local_dofs(proc,
420  update_cells,
421  dofs,
422  node_dof_starts,
423  node_dofs,
424  edge_dof_starts,
425  edge_dofs
426  );
427 
428  }
429  else
430  send_ghost_dofs(proc);
431  }
432  }
433  update_cells.clear();
434  node_dofs.clear();
435  node_dof_starts.clear();
436  edge_dofs.clear();
437  edge_dof_starts.clear();
438 }
439 
440 
442 {
443  if (dh_seq_ != nullptr) return;
444 
445  if ( !is_parallel_ )
446  {
447  dh_seq_ = std::make_shared<DOFHandlerMultiDim>(*this);
448  return;
449  }
450 
451  dh_seq_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
452 
453  dh_seq_->n_global_dofs_ = n_global_dofs_;
454  dh_seq_->lsize_ = n_global_dofs_;
455  dh_seq_->loffset_ = 0;
456  dh_seq_->mesh_ = mesh_;
457  dh_seq_->dof_ds_ = dof_ds_; // should be sequential distribution
458  dh_seq_->ds_ = ds_;
459  dh_seq_->is_parallel_ = false;
460  for (unsigned int i=0; i<n_global_dofs_; i++) dh_seq_->local_to_global_dof_idx_.push_back(i);
461 
464  &dh_seq_->max_elem_dofs_,
465  1,
466  MPI_UNSIGNED,
467  MPI_MAX,
469 
470  for (unsigned int i=0; i<mesh_->n_elements(); i++) dh_seq_->global_to_local_el_idx_[i] = mesh_->get_row_4_el()[i];
471 
472  // Auxiliary vectors cell_starts_loc and dof_indices_loc contain
473  // only local element data (without ghost elements).
474  // Then it is possible to create sequential vectors by simple reduce/gather operation.
475  vector<LongIdx> cell_starts_loc(mesh_->n_elements()+1, 0);
476  vector<LongIdx> dof_indices_loc;
477 
478  // construct cell_starts_loc
479  for (auto cell : this->own_range())
480  {
481  cell_starts_loc[mesh_->get_row_4_el()[cell.elm_idx()]+1] = cell.n_dofs();
482  }
483  for (unsigned int i=0; i<mesh_->n_elements(); ++i)
484  cell_starts_loc[i+1] += cell_starts_loc[i];
485 
486  // construct dof_indices_loc
487  dof_indices_loc.resize(cell_starts_loc[mesh_->n_elements()]);
488  for (auto cell : this->own_range())
489  {
490  for (unsigned int idof=0; idof<cell.n_dofs(); idof++)
491  dof_indices_loc[cell_starts_loc[mesh_->get_row_4_el()[cell.elm_idx()]]+idof] = local_to_global_dof_idx_[dof_indices[cell_starts[cell.local_idx()]+idof]];
492  }
493 
494  Distribution distr(dof_indices_loc.size(), PETSC_COMM_WORLD);
495  dh_seq_->cell_starts.resize(mesh_->n_elements()+1);
496  dh_seq_->dof_indices.resize(distr.size());
497 
498  MPI_Allreduce( cell_starts_loc.data(),
499  dh_seq_->cell_starts.data(),
500  cell_starts_loc.size(),
501  MPI_LONG_IDX,
502  MPI_SUM,
503  MPI_COMM_WORLD );
504 
505  MPI_Allgatherv( dof_indices_loc.data(),
506  dof_indices_loc.size(),
507  MPI_LONG_IDX,
508  dh_seq_->dof_indices.data(),
509  (const int *)distr.get_lsizes_array(),
510  (const int *)distr.get_starts_array(),
511  MPI_LONG_IDX,
512  MPI_COMM_WORLD );
513 
514  // create scatter from parallel to sequential vector
515  Vec v_from;
516  VecCreateMPI(PETSC_COMM_WORLD, lsize_, PETSC_DETERMINE, &v_from);
517  scatter_to_seq_ = std::make_shared<VecScatter>();
518  VecScatterCreateToAll(v_from, scatter_to_seq_.get(), NULL);
519  VecDestroy(&v_from);
520 
521  // create scatter for sequential dof handler (Warning: not tested)
522  Vec v_seq;
523  VecCreateSeq(PETSC_COMM_SELF, n_global_dofs_, &v_seq);
524  dh_seq_->scatter_to_seq_ = std::make_shared<VecScatter>();
525  VecScatterCreateToAll(v_seq, dh_seq_->scatter_to_seq_.get(), NULL);
526  VecDestroy(&v_seq);
527 }
528 
529 
531 {
532  if (is_parallel_ && el_ds_->np() > 1)
533  { // for parallel DH create vector with ghost values
535  VectorMPI vec(lsize_, ghost_dofs);
536  return vec;
537  } else {
538  VectorMPI vec(lsize_, PETSC_COMM_SELF);
539  return vec;
540  }
541 }
542 
543 
545 {
546  unsigned int ndofs = 0;
547  ndofs = cell_starts[cell.local_idx()+1]-cell_starts[cell.local_idx()];
548  indices.resize(ndofs);
549  for (unsigned int k=0; k<ndofs; k++)
550  indices[k] = local_to_global_dof_idx_[dof_indices[cell_starts[cell.local_idx()]+k]];
551 
552  return ndofs;
553 }
554 
555 
557 {}
558 
559 
560 
562 {
563  // create local arrays of elements
564  el_ds_ = mesh_->get_el_ds();
565 
566  // create local array of edges
567  for (auto edge : mesh_->edge_range())
568  {
569  bool is_edge_local = false;
570  for (uint sid=0; sid<edge.n_sides(); sid++)
571  if ( el_is_local(edge.side(sid)->element().idx()) )
572  {
573  is_edge_local = true;
574  break;
575  }
576  if (is_edge_local)
577  edg_4_loc.push_back(edge.idx());
578  }
579 
580  // create local array of neighbours
581  for (unsigned int inb=0; inb<mesh_->n_vb_neighbours(); inb++)
582  {
583  Neighbour *nb = &mesh_->vb_neighbours_[inb];
584  if ( el_is_local(nb->element().idx())
585  || el_is_local(nb->side()->element().idx()) )
586  nb_4_loc.push_back(inb);
587  }
588 
589  // init global to local element map with locally owned elements (later add ghost elements)
590  for ( unsigned int iel = 0; iel < el_ds_->lsize(); iel++ )
592 
593  // create array of local nodes
594  std::vector<bool> node_is_local(mesh_->tree->n_nodes(), false);
595  for (auto cell : this->own_range())
596  {
597  unsigned int obj_idx = mesh_->tree->obj_4_el()[cell.elm_idx()];
598  for (unsigned int nid=0; nid<cell.elm()->n_nodes(); nid++)
599  node_is_local[mesh_->tree->objects(cell.dim())[obj_idx].nodes[nid]] = true;
600  }
601 
602  // create array of local ghost cells
603  for ( auto cell : mesh_->elements_range() )
604  {
605  if (cell.proc() != el_ds_->myp())
606  {
607  bool has_local_node = false;
608  unsigned int obj_idx = mesh_->tree->obj_4_el()[cell.idx()];
609  for (unsigned int nid=0; nid<cell->n_nodes(); nid++)
610  if (node_is_local[mesh_->tree->objects(cell->dim())[obj_idx].nodes[nid]])
611  {
612  has_local_node = true;
613  break;
614  }
615  if (has_local_node)
616  {
617  ghost_4_loc.push_back(cell.idx());
618  ghost_proc.insert(cell.proc());
619  ghost_proc_el[cell.proc()].push_back(cell.idx());
620  global_to_local_el_idx_[cell.idx()] = el_ds_->lsize() - 1 + ghost_4_loc.size();
621  }
622  }
623  }
624  for (auto nb : nb_4_loc)
625  {
626  auto cell = mesh_->vb_neighbours_[nb].element();
627  if (!el_is_local(cell.idx()) && find(ghost_4_loc.begin(), ghost_4_loc.end(), cell.idx()) == ghost_4_loc.end())
628  {
629  ghost_4_loc.push_back(cell.idx());
630  ghost_proc.insert(cell.proc());
631  ghost_proc_el[cell.proc()].push_back(cell.idx());
632  global_to_local_el_idx_[cell.idx()] = el_ds_->lsize() - 1 + ghost_4_loc.size();
633  }
634  cell = mesh_->vb_neighbours_[nb].side()->element();
635  if (!el_is_local(cell.idx()) && find(ghost_4_loc.begin(), ghost_4_loc.end(), cell.idx()) == ghost_4_loc.end())
636  {
637  ghost_4_loc.push_back(cell.idx());
638  ghost_proc.insert(cell.proc());
639  ghost_proc_el[cell.proc()].push_back(cell.idx());
640  global_to_local_el_idx_[cell.idx()] = el_ds_->lsize() - 1 + ghost_4_loc.size();
641  }
642  }
643 }
644 
645 
646 bool DOFHandlerMultiDim::el_is_local(int index) const
647 {
648  return el_ds_->is_local(mesh_->get_row_4_el()[index]);
649 }
650 
651 
652 std::size_t DOFHandlerMultiDim::hash() const {
653  return this->n_global_dofs_;
654 }
655 
656 
658  auto bgn_it = make_iter<DHCellAccessor>( DHCellAccessor(this, 0) );
659  auto end_it = make_iter<DHCellAccessor>( DHCellAccessor(this, el_ds_->lsize()) );
660  return Range<DHCellAccessor>(bgn_it, end_it);
661 }
662 
663 
665  auto bgn_it = make_iter<DHCellAccessor>( DHCellAccessor(this, 0) );
666  auto end_it = make_iter<DHCellAccessor>( DHCellAccessor(this, el_ds_->lsize()+ghost_4_loc.size()) );
667  return Range<DHCellAccessor>(bgn_it, end_it);
668 }
669 
670 
672  auto bgn_it = make_iter<DHCellAccessor>( DHCellAccessor(this, el_ds_->lsize()) );
673  auto end_it = make_iter<DHCellAccessor>( DHCellAccessor(this, el_ds_->lsize()+ghost_4_loc.size()) );
674  return Range<DHCellAccessor>(bgn_it, end_it);
675 }
676 
677 
679  auto map_it = global_to_local_el_idx_.find((LongIdx)elm_idx); // find in global to local map
680  ASSERT( map_it != global_to_local_el_idx_.end() )(elm_idx).error("DH accessor can be create only for own or ghost elements!\n");
681  return DHCellAccessor(this, map_it->second);
682 }
683 
684 
686  stringstream s;
688 
689  s << "DOFHandlerMultiDim structure:" << endl;
690  s << "- is parallel: " << (is_parallel_?"true":"false") << endl;
691  s << "- proc id: " << el_ds_->myp() << endl;
692  s << "- global number of dofs: " << n_global_dofs_ << 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;
696 
697  for (auto cell : own_range())
698  {
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;
702  }
703  s << "- dofs on ghost cells:" << endl;
704  for (auto cell : ghost_range())
705  {
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;
709  }
710  s << "- locally owned dofs (" << lsize_ << "): ";
711  for (unsigned int i=0; i<lsize_; i++) s << local_to_global_dof_idx_[i] << " "; s << endl;
712  s << "- ghost dofs (" << local_to_global_dof_idx_.size() - lsize_ << "): ";
713  for (unsigned int i=lsize_; i<local_to_global_dof_idx_.size(); i++) s << local_to_global_dof_idx_[i] << " "; s << endl;
714  s << "- global-to-local-cell map:" << endl;
715  for (auto cell : global_to_local_el_idx_) s << "-- " << cell.first << " -> " << cell.second << " " << endl;
716  s << endl;
717 
718  printf("%s", s.str().c_str());
719 }
720 
721 
722 
723 
724 
725 
726 
727 
728 
729 SubDOFHandlerMultiDim::SubDOFHandlerMultiDim(std::shared_ptr<DOFHandlerMultiDim> dh, unsigned int component_idx)
730 : DOFHandlerMultiDim(*dh->mesh()),
731  parent_(dh),
732  fe_idx_(component_idx)
733 {
734  // create discrete space, we assume equal type of FE on each cell.
735  ASSERT_DBG( dynamic_cast<EqualOrderDiscreteSpace *>(dh->ds().get()) != nullptr )
736  .error("sub_handler can be used only with dof handler using EqualOrderDiscreteSpace!");
737  FESystem<0> *fe_sys0 = dynamic_cast<FESystem<0>*>( dh->ds()->fe()[0_d].get() );
738  FESystem<1> *fe_sys1 = dynamic_cast<FESystem<1>*>( dh->ds()->fe()[1_d].get() );
739  FESystem<2> *fe_sys2 = dynamic_cast<FESystem<2>*>( dh->ds()->fe()[2_d].get() );
740  FESystem<3> *fe_sys3 = dynamic_cast<FESystem<3>*>( dh->ds()->fe()[3_d].get() );
741  ASSERT_DBG( fe_sys0 != nullptr ).error("sub_handler assumes that dof handler uses FESystem<0>!");
742  ASSERT_DBG( fe_sys1 != nullptr ).error("sub_handler assumes that dof handler uses FESystem<1>!");
743  ASSERT_DBG( fe_sys2 != nullptr ).error("sub_handler assumes that dof handler uses FESystem<2>!");
744  ASSERT_DBG( fe_sys3 != nullptr ).error("sub_handler assumes that dof handler uses FESystem<3>!");
745  MixedPtr<FiniteElement> sub_mixed(fe_sys0->fe()[component_idx],
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);
750 
751  is_parallel_ = dh->is_parallel_;
752 
753  // create list of dofs for sub_dh on one cell
754  FESystemFunctionSpace *fs[4] = {
755  dynamic_cast<FESystemFunctionSpace*>( fe_sys0->function_space_.get() ),
756  dynamic_cast<FESystemFunctionSpace*>( fe_sys1->function_space_.get() ),
757  dynamic_cast<FESystemFunctionSpace*>( fe_sys2->function_space_.get() ),
758  dynamic_cast<FESystemFunctionSpace*>( fe_sys3->function_space_.get() ) };
759  for (unsigned int d=0; d<=3; d++)
760  ASSERT_DBG( fs[d] != nullptr ).error("Function space must be of type FESystemFunctionSpace!" );
761  vector<unsigned int> sub_fe_dofs[4];
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);
770 
772  dof_indices.resize(cell_starts[cell_starts.size()-1]);
773  // sub_local_indices maps local dofs of parent handler to local dofs of sub-handler
774  vector<LongIdx> sub_local_indices(dh->local_to_global_dof_idx_.size(), INVALID_DOF);
775  map<LongIdx,LongIdx> global_to_local_dof_idx;
776  // first add owned dofs to local_to_global_dof_idx_ and sub_local_indices
777  for (auto cell : dh->local_range())
778  {
779  LocDofVec cell_dof_indices = cell.get_loc_dof_indices();
780  for (auto sub_dof : sub_fe_dofs[cell.dim()])
781  {
782  if (cell_dof_indices[sub_dof] < static_cast<int>(dh->lsize_) &&
783  sub_local_indices[cell_dof_indices[sub_dof]] == INVALID_DOF)
784  {
785  sub_local_indices[cell_dof_indices[sub_dof]] = local_to_global_dof_idx_.size();
786  parent_dof_idx_.push_back(cell_dof_indices[sub_dof]);
787  global_to_local_dof_idx[parent_->local_to_global_dof_idx_[cell_dof_indices[sub_dof]]] = local_to_global_dof_idx_.size();
789  }
790  }
791  }
793  // then do the same for ghost dofs and set dof_indices
794  for (auto cell : dh->local_range())
795  {
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()])
799  {
800  if (sub_local_indices[cell_dof_indices[sub_dof]] == INVALID_DOF)
801  {
802  sub_local_indices[cell_dof_indices[sub_dof]] = local_to_global_dof_idx_.size();
803  parent_dof_idx_.push_back(cell_dof_indices[sub_dof]);
804  // temporarily we keep the global dof idx of parent dh, we replace it later from the owning processor
805  local_to_global_dof_idx_.push_back(parent_->local_to_global_dof_idx_[cell_dof_indices[sub_dof]]);
806  }
807  dof_indices[cell_starts[cell.local_idx()]+idof++] = sub_local_indices[cell_dof_indices[sub_dof]];
808  }
809  }
810 
811  dof_ds_ = std::make_shared<Distribution>(lsize_, PETSC_COMM_WORLD);
812  n_global_dofs_ = dof_ds_->size();
813  loffset_ = dof_ds_->get_starts_array()[dof_ds_->myp()];
814 
815  // shift dof indices
816  if (loffset_ > 0)
817  for (unsigned int i=0; i<lsize_; i++)
819 
820  // communicate ghost values
821  // first propagate from lower procs to higher procs and then vice versa
822  for (unsigned int from_higher = 0; from_higher < 2; from_higher++)
823  {
824  for (unsigned int proc : ghost_proc)
825  {
826  if ((proc > el_ds_->myp()) == from_higher)
827  { // receive dofs from master processor
828  vector<LongIdx> dofs;
829  receive_sub_ghost_dofs(proc, dofs);
830  }
831  else
832  send_sub_ghost_dofs(proc, global_to_local_dof_idx);
833  }
834  }
835 }
836 
837 
839 {
840  // send number of ghost dofs required from the other processor
842  for (unsigned int i=lsize_; i<local_to_global_dof_idx_.size(); i++)
843  if (parent_->dof_ds_->get_proc(parent_->local_to_global_dof_idx_[parent_dof_idx_[i]]) == proc)
844  dof_indices.push_back(parent_->local_to_global_dof_idx_[parent_dof_idx_[i]]);
845  unsigned int n_ghosts = dof_indices.size();
846  MPI_Send(&n_ghosts, 1, MPI_UNSIGNED, proc, 0, MPI_COMM_WORLD);
847 
848  // send indices of dofs required
849  MPI_Send(dof_indices.data(), n_ghosts, MPI_LONG_IDX, proc, 1, MPI_COMM_WORLD);
850 
851  // receive dofs
852  dofs.resize(n_ghosts);
853  MPI_Recv(dofs.data(), n_ghosts, MPI_LONG_IDX, proc, 2, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
854 
855  // update ghost dofs
856  unsigned int idof = 0;
857  for (unsigned int i=lsize_; i<local_to_global_dof_idx_.size(); i++)
858  if (parent_->dof_ds_->get_proc(parent_->local_to_global_dof_idx_[parent_dof_idx_[i]]) == proc)
859  local_to_global_dof_idx_[i] = dofs[idof++];
860 }
861 
862 
863 void SubDOFHandlerMultiDim::send_sub_ghost_dofs(unsigned int proc, const map<LongIdx,LongIdx> &global_to_local_dof_idx)
864 {
865  // receive number of dofs required by the other processor
866  unsigned int n_ghosts;
867  MPI_Recv(&n_ghosts, 1, MPI_UNSIGNED, proc, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
868 
869  // receive global indices of dofs required
870  vector<LongIdx> dof_indices(n_ghosts);
871  MPI_Recv(dof_indices.data(), n_ghosts, MPI_LONG_IDX, proc, 1, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
872 
873  // send global dof indices relative to the sub-handler
874  vector<LongIdx> dofs;
875  for (auto global_dof : dof_indices)
876  dofs.push_back(global_to_local_dof_idx.at(global_dof) + dof_ds_->begin());
877  MPI_Send(dofs.data(), dofs.size(), MPI_LONG_IDX, proc, 2, MPI_COMM_WORLD);
878 }
879 
880 
882 {
883  ASSERT_DBG( vec.size() == parent_->local_to_global_dof_idx_.size() ).error("Incompatible parent vector in update_subvector()!");
884  ASSERT_DBG( subvec.size() == local_to_global_dof_idx_.size() ).error("Incompatible subvector in update_subvector()!");
885 
886  for (unsigned int i=0; i<parent_dof_idx_.size(); i++)
887  subvec.set( i, vec.get(parent_dof_idx_[i]) );
888 }
889 
890 
892 {
893  ASSERT_DBG( vec.size() == parent_->local_to_global_dof_idx_.size() ).error("Incompatible parent vector in update_subvector()!");
894  ASSERT_DBG( subvec.size() == local_to_global_dof_idx_.size() ).error("Incompatible subvector in update_subvector()!");
895 
896  for (unsigned int i=0; i<parent_dof_idx_.size(); i++)
897  vec.set( parent_dof_idx_[i], subvec.get(i) );
898 }
899 
900 
901 
902 
903 
904 
905 
906 
907 
#define MPI_Recv(buf, count, datatype, source, tag, comm, status)
Definition: mpi.h:271
void receive_ghost_dofs(unsigned int proc, std::vector< LongIdx > &dofs)
Obtain dof numbers on ghost elements from other processor.
Definition: dofhandler.cc:161
std::shared_ptr< DOFHandlerMultiDim > sequential()
Returns sequential version of the current dof handler.
Definition: dofhandler.cc:61
unsigned int size() const
Return size of output data.
Definition: vector_mpi.hh:98
const Dof & cell_dof(unsigned int idof) const
Return dof on a given cell.
virtual ~DOFHandlerBase()
Destructor.
Definition: dofhandler.cc:41
void create_sequential()
Communicate local dof indices to all processors and create new sequential dof handler.
Definition: dofhandler.cc:441
void init_cell_starts()
Initialize vector of starting indices for elements.
Definition: dofhandler.cc:75
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
unsigned int n_nodes() const
Definition: elements.h:126
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
const DHCellAccessor cell_accessor_from_element(unsigned int elm_idx) const
Return DHCellAccessor appropriate to ElementAccessor of given idx.
Definition: dofhandler.cc:678
arma::Col< IntIdx > LocDofVec
Definition: index_types.hh:28
ArmaVec< double, N > vec
Definition: armor.hh:885
LongIdx * get_row_4_el() const
Definition: mesh.h:174
unsigned int uint
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.
Definition: dofhandler.cc:89
void update_parent_vector(VectorMPI &vec, const VectorMPI &subvec)
Update values in parent vector from values of subvector.
Definition: dofhandler.cc:891
const std::vector< MeshObject > & objects(unsigned int dim) const
unsigned int edge_idx() const
Returns global index of the edge connected to the side.
DuplicateNodes * tree
Definition: mesh.h:279
Range helper class.
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
Definition: fe_system.hh:147
Distribution * el_ds_
Distribution of elements.
Definition: dofhandler.hh:443
static const int INVALID_NFACE
Definition: dofhandler.hh:386
std::shared_ptr< DOFHandlerMultiDim > parent_
Parent dof handler.
Definition: dofhandler.hh:503
unsigned int dim() const
Return dimension of element appropriate to cell.
unsigned int dim
Association to n-face of given dimension (point, line, triangle, tetrahedron.
unsigned int local_idx() const
Return local index to element (index of DOF handler).
void distribute_dofs(std::shared_ptr< DiscreteSpace > ds)
Distributes degrees of freedom on the mesh needed for the given discrete space.
Definition: dofhandler.cc:303
#define MPI_SUM
Definition: mpi.h:196
Definition: mesh.h:77
Cell accessor allow iterate over DOF handler cells.
std::vector< IntIdx > dof_indices
Dof numbers on local and ghost elements.
Definition: dofhandler.hh:427
std::shared_ptr< VecScatter > sequential_scatter()
Returns scatter context from parallel to sequential vectors.
Definition: dofhandler.cc:68
Class FESystem for compound finite elements.
std::vector< LongIdx > local_to_global_dof_idx_
Maps local and ghost dof indices to global ones.
Definition: dofhandler.hh:435
SideIter side(const unsigned int loc_index)
#define MPI_Send(buf, count, datatype, dest, tag, comm)
Definition: mpi.h:263
SubDOFHandlerMultiDim(std::shared_ptr< DOFHandlerMultiDim > dh, unsigned int component_idx)
Creates a new dof handler for a component of FESystem.
Definition: dofhandler.cc:729
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
Definition: asserts.hh:347
void init_status(std::vector< short int > &node_status, std::vector< short int > &edge_status)
Initialize node_status and edge_status.
Definition: dofhandler.cc:115
std::shared_ptr< Distribution > dof_ds_
Distribution of dofs associated to local process.
Definition: dofhandler.hh:131
#define MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm)
Definition: mpi.h:587
Mesh * mesh() const
Returns the mesh.
Definition: dofhandler.hh:78
unsigned int n_dofs() const
Return number of dofs on given cell.
ElementAccessor< 3 > element() const
Returns iterator to the element of the side.
std::shared_ptr< DOFHandlerMultiDim > dh_seq_
Sequential dof handler associated to the current (parallel) one.
Definition: dofhandler.hh:399
unsigned int lsize_
Number of dofs associated to local process.
Definition: dofhandler.hh:113
unsigned int n_face_idx
Index of n-face to which the dof is associated.
map< unsigned int, vector< LongIdx > > ghost_proc_el
Arrays of ghost cells for each neighbouring processor.
Definition: dofhandler.hh:458
static const int ASSIGNED_NFACE
Definition: dofhandler.hh:388
void send_ghost_dofs(unsigned int proc)
Send dof numbers to other processor.
Definition: dofhandler.cc:182
std::vector< LongIdx > cell_starts
Starting indices for local (owned+ghost) element dofs.
Definition: dofhandler.hh:419
#define OLD_ASSERT(...)
Definition: global_defs.h:108
unsigned int n_vb_neighbours() const
Definition: mesh.cc:242
ElementAccessor< 3 > element()
Definition: neighbours.h:161
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
bool is_local(unsigned int idx) const
identify local index
Compound finite element on dim dimensional simplex.
Definition: fe_system.hh:101
Mesh * mesh_
Pointer to the mesh to which the dof handler is associated.
Definition: dofhandler.hh:126
SideIter side()
Definition: neighbours.h:145
unsigned int n_sides() const
Definition: elements.h:132
Provides the numbering of the finite element degrees of freedom on the computational mesh...
Definition: dofhandler.hh:151
static const int VALID_NFACE
Definition: dofhandler.hh:387
Range< DHCellAccessor > ghost_range() const
Returns range over ghosts DOF handler cells.
Definition: dofhandler.cc:671
vector< LongIdx > nb_4_loc
Local neighbour index -> global neighbour index.
Definition: dofhandler.hh:449
unsigned int max_elem_dofs_
Max. number of dofs per element.
Definition: dofhandler.hh:121
vector< LongIdx > ghost_4_loc
Indices of ghost cells (neighbouring with local elements).
Definition: dofhandler.hh:452
std::unordered_map< LongIdx, LongIdx > global_to_local_el_idx_
Maps global element index into local/ghost index (obsolete).
Definition: dofhandler.hh:440
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1241
unsigned int n_nodes() const
unsigned int np() const
get num of processors
Distribution * get_el_ds() const
Definition: mesh.h:171
void make_elem_partitioning()
Prepare parallel distribution of elements, edges and neighbours.
Definition: dofhandler.cc:561
#define MPI_LONG_IDX
Definition: index_types.hh:30
void print() const
Output structure of dof handler.
Definition: dofhandler.cc:685
std::shared_ptr< Distribution > distr() const
Definition: dofhandler.hh:73
std::vector< LongIdx > parent_dof_idx_
Local indices in the parent handler.
Definition: dofhandler.hh:509
unsigned int myp() const
get my processor
void set(unsigned int pos, double val)
Set value on given position.
Definition: vector_mpi.hh:111
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.
Definition: dofhandler.cc:863
Support classes for parallel programing.
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Definition: mpi.h:612
double get(unsigned int pos) const
Return value on given position.
Definition: vector_mpi.hh:105
vector< Neighbour > vb_neighbours_
Definition: mesh.h:284
#define MPI_STATUS_IGNORE
Definition: mpi.h:203
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: index_types.hh:24
vector< LongIdx > edg_4_loc
Local edge index -> global edge index.
Definition: dofhandler.hh:446
Range< DHCellAccessor > own_range() const
Returns range of DOF handler cells (only range of own without ghost cells)
Definition: dofhandler.cc:657
unsigned int loffset_
Index of the first dof on the local process.
Definition: dofhandler.hh:118
static const int INVALID_DOF
Definition: dofhandler.hh:389
std::shared_ptr< FunctionSpace > function_space_
Function space defining the FE.
set< unsigned int > ghost_proc
Processors of ghost elements.
Definition: dofhandler.hh:455
~DOFHandlerMultiDim() override
Destructor.
Definition: dofhandler.cc:556
#define ASSERT_DBG(expr)
#define MPI_COMM_WORLD
Definition: mpi.h:123
std::shared_ptr< DiscreteSpace > ds_
Pointer to the discrete space for which the handler distributes dofs.
Definition: dofhandler.hh:393
#define MPI_MAX
Definition: mpi.h:197
void receive_sub_ghost_dofs(unsigned int proc, vector< LongIdx > &dofs)
Get global dof indices of ghost dofs for sub-handler.
Definition: dofhandler.cc:838
void update_subvector(const VectorMPI &vec, VectorMPI &subvec)
Update values in subvector from parent vector.
Definition: dofhandler.cc:881
#define MPI_UNSIGNED
Definition: mpi.h:165
unsigned int get_dof_indices(const DHCellAccessor &cell, std::vector< LongIdx > &indices) const override
Returns the global indices of dofs associated to the cell.
Definition: dofhandler.cc:544
Abstract class for description of finite elements.
virtual VectorMPI create_vector()
Allocates PETSc vector according to the dof distribution.
Definition: dofhandler.cc:530
unsigned int n_global_dofs_
Number of global dofs assigned by the handler.
Definition: dofhandler.hh:108
unsigned int n_edges() const
Definition: mesh.h:151
bool el_is_local(int index) const
Definition: dofhandler.cc:646
DOFHandlerMultiDim(Mesh &_mesh, bool make_elem_part=true)
Constructor.
Definition: dofhandler.cc:49
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.
Definition: dofhandler.cc:213
Range< DHCellAccessor > local_range() const
Returns range over own and ghost cells of DOF handler.
Definition: dofhandler.cc:664
const std::vector< unsigned int > & obj_4_el() const
bool is_parallel_
Indicator for parallel/sequential dof handler.
Definition: dofhandler.hh:396
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:181
std::shared_ptr< VecScatter > scatter_to_seq_
Scatter context for parallel to sequential vectors.
Definition: dofhandler.hh:402
LongIdx * get_el_4_loc() const
Definition: mesh.h:177
Range< Edge > edge_range() const
Returns range of edges.
Definition: mesh.cc:1253
Implementation of range helper class.
friend class DHCellAccessor
Definition: dofhandler.hh:265
std::size_t hash() const override
Definition: dofhandler.cc:652
void printf(BasicWriter< Char > &w, BasicCStringRef< Char > format, ArgList args)
Definition: printf.h:444
std::shared_ptr< DiscreteSpace > ds() const
Return pointer to discrete space for which the handler distributes dofs.
Definition: dofhandler.hh:255
unsigned int lsize(int proc) const
get local size