Flow123d  release_3.0.0-973-g92f55e826
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 "fem/dofhandler.hh"
20 #include "fem/finite_element.hh"
21 #include "fem/dh_cell_accessor.hh"
22 #include "mesh/mesh.h"
23 #include "mesh/duplicate_nodes.h"
24 #include "mesh/partitioning.hh"
25 #include "mesh/long_idx.hh"
26 #include "mesh/accessors.hh"
27 #include "mesh/range_wrapper.hh"
28 #include "mesh/neighbours.h"
29 #include "la/distribution.hh"
30 
31 
32 
33 
35 {}
36 
37 
38 
39 
40 
41 
42 DOFHandlerMultiDim::DOFHandlerMultiDim(Mesh& _mesh, bool make_elem_part)
43  : DOFHandlerBase(_mesh),
44  ds_(nullptr),
45  is_parallel_(true),
46  dh_seq_(nullptr),
47  scatter_to_seq_(nullptr),
48  el_ds_(nullptr)
49 {
50  if (make_elem_part) make_elem_partitioning();
51 }
52 
53 
54 std::shared_ptr<DOFHandlerMultiDim> DOFHandlerMultiDim::sequential()
55 {
57  return dh_seq_;
58 }
59 
60 
61 std::shared_ptr<VecScatter> DOFHandlerMultiDim::sequential_scatter()
62 {
64  return scatter_to_seq_;
65 }
66 
67 
69 {
70  // get number of dofs per element and then set up cell_starts
72  for (auto cell : this->local_range())
73  {
74  cell_starts[cell.local_idx()+1] = cell.n_dofs();
75  max_elem_dofs_ = max( (int)max_elem_dofs_, (int)cell.n_dofs() );
76  }
77  for (unsigned int i=0; i<cell_starts.size()-1; ++i)
78  cell_starts[i+1] += cell_starts[i];
79 }
80 
81 
83 {
84  // initialize dofs on nodes
85  // We must separate dofs for dimensions because FE functions
86  // may be discontinuous on nodes shared by different
87  // dimensions.
88  unsigned int n_node_dofs = 0;
89 
90  for (unsigned int nid=0; nid<mesh_->tree->n_nodes(); nid++)
91  {
92  node_dof_starts.push_back(n_node_dofs);
93  n_node_dofs += ds_->n_node_dofs(nid);
94  }
95  node_dof_starts.push_back(n_node_dofs);
96 }
97 
98 
100 {
101  // mark local dofs
102  for (auto cell : this->own_range())
103  {
104  for (unsigned int n=0; n<cell.dim()+1; n++)
105  {
106  unsigned int nid = mesh_->tree->objects(cell.dim())[mesh_->tree->obj_4_el()[cell.elm_idx()]].nodes[n];
107  node_status[nid] = VALID_NODE;
108  }
109  }
110 
111  // unmark dofs on ghost cells from lower procs
112  for (auto cell : this->ghost_range())
113  {
114  if (cell.elm().proc() < el_ds_->myp())
115  {
116  for (unsigned int n=0; n<cell.dim()+1; n++)
117  {
118  unsigned int nid = mesh_->tree->objects(cell.dim())[mesh_->tree->obj_4_el()[cell.elm_idx()]].nodes[n];
119  node_status[nid] = INVALID_NODE;
120  }
121  }
122  }
123 }
124 
125 
127 {
128  // send number of elements required from the other processor
129  unsigned int n_elems = ghost_proc_el[proc].size();
130  MPI_Send(&n_elems, 1, MPI_UNSIGNED, proc, 0, MPI_COMM_WORLD);
131 
132  // send indices of elements required
133  MPI_Send(&(ghost_proc_el[proc][0]), n_elems, MPI_LONG_IDX, proc, 1, MPI_COMM_WORLD);
134 
135  // receive numbers of dofs on required elements
136  vector<unsigned int> n_dofs(n_elems);
137  MPI_Recv(&(n_dofs[0]), n_elems, MPI_UNSIGNED, proc, 2, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
138 
139  // receive dofs on required elements
140  unsigned int n_dofs_sum = 0;
141  for (auto nd : n_dofs) n_dofs_sum += nd;
142  dofs.resize(n_dofs_sum);
143  MPI_Recv(&(dofs[0]), n_dofs_sum, MPI_LONG_IDX, proc, 3, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
144 }
145 
146 
147 void DOFHandlerMultiDim::send_ghost_dofs(unsigned int proc)
148 {
149  // receive number of elements required by the other processor
150  unsigned int n_elems;
151  MPI_Recv(&n_elems, 1, MPI_UNSIGNED, proc, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
152 
153  // receive indices of elements required
154  vector<LongIdx> elems(n_elems);
155  MPI_Recv(&(elems[0]), n_elems, MPI_LONG_IDX, proc, 1, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
156 
157  // send numbers of dofs on required elements
158  vector<unsigned int> n_dofs;
159  for (LongIdx el : elems)
160  {
161  auto cell = this->cell_accessor_from_element(el);
162  n_dofs.push_back(cell_starts[cell.local_idx()+1] - cell_starts[cell.local_idx()]);
163  }
164  MPI_Send(&(n_dofs[0]), n_elems, MPI_UNSIGNED, proc, 2, MPI_COMM_WORLD);
165 
166  // send dofs on the required elements
167  vector<LongIdx> dofs;
168  for (LongIdx el : elems)
169  {
170  auto cell = this->cell_accessor_from_element(el);
171  dofs.insert(dofs.end(), dof_indices.begin()+cell_starts[cell.local_idx()], dof_indices.begin()+cell_starts[cell.local_idx()+1]);
172  }
173  MPI_Send(&(dofs[0]), dofs.size(), MPI_LONG_IDX, proc, 3, MPI_COMM_WORLD);
174 }
175 
176 
178  const std::vector<bool> &update_cells,
179  const std::vector<LongIdx> &dofs,
180  const std::vector<LongIdx> &node_dof_starts,
181  std::vector<LongIdx> &node_dofs
182  )
183 {
184  // update dof_indices on ghost cells
185  unsigned int dof_offset=0;
186  for (unsigned int gid=0; gid<ghost_proc_el[proc].size(); gid++)
187  {
188  DHCellAccessor dh_cell = this->cell_accessor_from_element( ghost_proc_el[proc][gid] );
189 
190  for (unsigned dof=0; dof<dh_cell.n_dofs(); dof++)
191  dof_indices[cell_starts[dh_cell.local_idx()]+dof] = dofs[dof_offset+dof];
192 
193  vector<unsigned int> loc_node_dof_count(dh_cell.elm()->n_nodes(), 0);
194  for (unsigned int idof = 0; idof<dh_cell.n_dofs(); ++idof)
195  {
196  if (dh_cell.cell_dof(idof).dim == 0)
197  { // update nodal dof
198  unsigned int dof_nface_idx = dh_cell.cell_dof(idof).n_face_idx;
199  unsigned int nid = mesh_->tree->objects(dh_cell.dim())[mesh_->tree->obj_4_el()[dh_cell.elm_idx()]].nodes[dof_nface_idx];
200 
201  if (node_dofs[node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx]] == INVALID_DOF)
202  {
203  node_dofs[node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx]] = dofs[dof_offset+idof];
204  local_to_global_dof_idx_.push_back(dofs[dof_offset+idof]);
205  }
206 
207  loc_node_dof_count[dof_nface_idx]++;
208  }
209  else if (dh_cell.cell_dof(idof).dim == dh_cell.dim())
210  {
211  local_to_global_dof_idx_.push_back(dofs[dof_offset+idof]);
212  }
213  }
214 
215  dof_offset += dh_cell.n_dofs();
216  }
217 
218  // update dof_indices on local elements
219  for (auto cell : this->own_range())
220  {
221  if (!update_cells[cell.local_idx()]) continue;
222 
223  // loop over element dofs
224  vector<unsigned int> loc_node_dof_count(cell.elm()->n_nodes(), 0);
225  for (unsigned int idof = 0; idof<cell.n_dofs(); ++idof)
226  {
227  if (cell.cell_dof(idof).dim == 0)
228  {
229  unsigned int dof_nface_idx = cell.cell_dof(idof).n_face_idx;
230  if (dof_indices[cell_starts[cell.local_idx()]+idof] == INVALID_DOF)
231  { // update nodal dof
232  unsigned int nid = mesh_->tree->objects(cell.dim())[mesh_->tree->obj_4_el()[cell.elm_idx()]].nodes[dof_nface_idx];
233  dof_indices[cell_starts[cell.local_idx()]+idof] = node_dofs[node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx]];
234  }
235  loc_node_dof_count[dof_nface_idx]++;
236  }
237  }
238  }
239 }
240 
241 
242 void DOFHandlerMultiDim::distribute_dofs(std::shared_ptr<DiscreteSpace> ds)
243 {
244  // First check if dofs are already distributed.
245  OLD_ASSERT(ds_ == nullptr, "Attempt to distribute DOFs multiple times!");
246 
247  ds_ = ds;
248 
249  std::vector<LongIdx> node_dofs, node_dof_starts;
251  std::vector<bool> update_cells(el_ds_->lsize(), false);
252  unsigned int next_free_dof = 0;
253 
255  init_node_dof_starts(node_dof_starts);
256  node_dofs.resize(node_dof_starts[node_dof_starts.size()-1], (LongIdx)INVALID_DOF);
257  init_node_status(node_status);
258 
259  // Distribute dofs on local elements.
260  dof_indices.resize(cell_starts[cell_starts.size()-1]);
261  local_to_global_dof_idx_.reserve(dof_indices.size());
262  for (auto cell : this->own_range())
263  {
264 
265  // loop over element dofs
266  vector<unsigned int> loc_node_dof_count(cell.elm()->n_nodes(), 0);
267  for (unsigned int idof = 0; idof<cell.n_dofs(); ++idof)
268  {
269  unsigned int dof_dim = cell.cell_dof(idof).dim;
270  unsigned int dof_nface_idx = cell.cell_dof(idof).n_face_idx;
271 
272  if (dof_dim == 0)
273  { // add dofs shared by nodes
274  unsigned int nid = mesh_->tree->objects(cell.dim())[mesh_->tree->obj_4_el()[cell.elm_idx()]].nodes[dof_nface_idx];
275  unsigned int node_dof_idx = node_dof_starts[nid]+loc_node_dof_count[dof_nface_idx];
276 
277  switch (node_status[nid])
278  {
279  case VALID_NODE:
280  for (int i=0; i<node_dof_starts[nid+1] - node_dof_starts[nid]; i++)
281  {
282  local_to_global_dof_idx_.push_back(next_free_dof);
283  node_dofs[node_dof_starts[nid]+i] = next_free_dof++;
284  }
285  node_status[nid] = ASSIGNED_NODE;
286  break;
287  case INVALID_NODE:
288  node_dofs[node_dof_idx] = INVALID_DOF;
289  update_cells[cell.local_idx()] = true;
290  break;
291  }
292  dof_indices[cell_starts[cell.local_idx()]+idof] = node_dofs[node_dof_idx];
293  loc_node_dof_count[dof_nface_idx]++;
294  }
295  else if (dof_dim == cell.dim())
296  {
297  // add dofs owned only by the element
298  local_to_global_dof_idx_.push_back(next_free_dof);
299  dof_indices[cell_starts[cell.local_idx()]+idof] = next_free_dof++;
300  }
301  else
302  ASSERT(false).error("Unsupported dof n_face.");
303  }
304  }
305  node_status.clear();
306 
307  lsize_ = next_free_dof;
308 
309  // communicate n_dofs across all processes
310  dof_ds_ = std::make_shared<Distribution>(lsize_, PETSC_COMM_WORLD);
311  n_global_dofs_ = dof_ds_->size();
312 
313  // shift dof indices
314  loffset_ = dof_ds_->get_starts_array()[dof_ds_->myp()];
315  if (loffset_ > 0)
316  {
317  for (unsigned int i=0; i<dof_indices.size(); i++)
318  if (dof_indices[i] != INVALID_DOF)
319  dof_indices[i] += loffset_;
320  for (auto &i : local_to_global_dof_idx_)
321  i += loffset_;
322  }
323 
324  // communicate dofs from ghost cells
325  // first propagate from lower procs to higher procs and then vice versa
326  for (unsigned int from_higher = 0; from_higher < 2; from_higher++)
327  {
328  for (unsigned int proc : ghost_proc)
329  {
330  if ((proc > el_ds_->myp()) == from_higher)
331  { // receive dofs from master processor
332  vector<LongIdx> dofs;
333  receive_ghost_dofs(proc, dofs);
334 
335  // update dof_indices and node_dofs on ghost elements
336  update_local_dofs(proc, update_cells, dofs, node_dof_starts, node_dofs);
337 
338  }
339  else
340  send_ghost_dofs(proc);
341  }
342  }
343  update_cells.clear();
344  node_dofs.clear();
345  node_dof_starts.clear();
346 }
347 
348 
350 {
351  if (dh_seq_ != nullptr) return;
352 
353  if ( !is_parallel_ )
354  {
355  dh_seq_ = std::make_shared<DOFHandlerMultiDim>(*this);
356  return;
357  }
358 
359  dh_seq_ = std::make_shared<DOFHandlerMultiDim>(*mesh_);
360 
361  dh_seq_->n_global_dofs_ = n_global_dofs_;
362  dh_seq_->lsize_ = n_global_dofs_;
363  dh_seq_->loffset_ = 0;
364  dh_seq_->mesh_ = mesh_;
365  dh_seq_->dof_ds_ = dof_ds_; // should be sequential distribution
366  dh_seq_->ds_ = ds_;
367  dh_seq_->is_parallel_ = false;
368 
371  &dh_seq_->max_elem_dofs_,
372  1,
373  MPI_UNSIGNED,
374  MPI_MAX,
376 
377  for (unsigned int i=0; i<mesh_->n_elements(); i++) dh_seq_->global_to_local_el_idx_[i] = mesh_->get_row_4_el()[i];
378 
379  // Auxiliary vectors cell_starts_loc and dof_indices_loc contain
380  // only local element data (without ghost elements).
381  // Then it is possible to create sequential vectors by simple reduce/gather operation.
382  vector<LongIdx> cell_starts_loc(mesh_->n_elements()+1, 0);
383  vector<LongIdx> dof_indices_loc;
384 
385  // construct cell_starts_loc
386  for (auto cell : this->own_range())
387  {
388  cell_starts_loc[mesh_->get_row_4_el()[cell.elm_idx()]+1] = cell.n_dofs();
389  }
390  for (unsigned int i=0; i<mesh_->n_elements(); ++i)
391  cell_starts_loc[i+1] += cell_starts_loc[i];
392 
393  // construct dof_indices_loc
394  dof_indices_loc.resize(cell_starts_loc[mesh_->n_elements()]);
395  for (auto cell : this->own_range())
396  {
397  for (unsigned int idof=0; idof<cell.n_dofs(); idof++)
398  dof_indices_loc[cell_starts_loc[mesh_->get_row_4_el()[cell.elm_idx()]]+idof] = dof_indices[cell_starts[cell.local_idx()]+idof];
399  }
400 
401  Distribution distr(dof_indices_loc.size(), PETSC_COMM_WORLD);
402  dh_seq_->cell_starts.resize(mesh_->n_elements()+1);
403  dh_seq_->dof_indices.resize(distr.size());
404 
405  MPI_Allreduce( cell_starts_loc.data(),
406  dh_seq_->cell_starts.data(),
407  cell_starts_loc.size(),
408  MPI_LONG_IDX,
409  MPI_SUM,
410  MPI_COMM_WORLD );
411 
412  MPI_Allgatherv( dof_indices_loc.data(),
413  dof_indices_loc.size(),
414  MPI_LONG_IDX,
415  dh_seq_->dof_indices.data(),
416  (const int *)distr.get_lsizes_array(),
417  (const int *)distr.get_starts_array(),
418  MPI_LONG_IDX,
419  MPI_COMM_WORLD );
420 
421  // create scatter from parallel to sequential vector
422  Vec v_from;
423  VecCreateMPI(PETSC_COMM_WORLD, lsize_, PETSC_DETERMINE, &v_from);
424  scatter_to_seq_ = std::make_shared<VecScatter>();
425  VecScatterCreateToAll(v_from, scatter_to_seq_.get(), NULL);
426  VecDestroy(&v_from);
427 
428  // create scatter for sequential dof handler (Warning: not tested)
429  Vec v_seq;
430  VecCreateSeq(PETSC_COMM_SELF, n_global_dofs_, &v_seq);
431  dh_seq_->scatter_to_seq_ = std::make_shared<VecScatter>();
432  VecScatterCreateToAll(v_seq, dh_seq_->scatter_to_seq_.get(), NULL);
433  VecDestroy(&v_seq);
434 }
435 
436 
438 {
439  if (is_parallel_ && el_ds_->np() > 1)
440  { // for parallel DH create vector with ghost values
442  VectorMPI vec(lsize_, ghost_dofs);
443  return vec;
444  } else {
445  VectorMPI vec(lsize_, PETSC_COMM_SELF);
446  return vec;
447  }
448 }
449 
450 
452 {
453  unsigned int ndofs = 0;
454  LongIdx cell_idx = global_to_local_el_idx_.at(cell.idx());
455  ndofs = cell_starts[cell_idx+1]-cell_starts[cell_idx];
456  for (unsigned int k=0; k<ndofs; k++)
457  indices[k] = dof_indices[cell_starts[cell_idx]+k];
458 
459  return ndofs;
460 }
461 
462 
463 
465 {
466  unsigned int ndofs = 0;
467  LongIdx cell_idx = global_to_local_el_idx_.at(cell.idx());
468  ndofs = cell_starts[cell_idx+1]-cell_starts[cell_idx];
469  for (unsigned int k=0; k<ndofs; k++)
470  indices[k] = cell_starts[cell_idx]+k;
471 
472  return ndofs;
473 }
474 
475 
477 {}
478 
479 
480 
482 {
483  // create local arrays of elements
484  el_ds_ = mesh_->get_el_ds();
485 
486  // create local array of edges
487  for (unsigned int iedg=0; iedg<mesh_->n_edges(); iedg++)
488  {
489  bool is_edge_local = false;
490  Edge *edg = &mesh_->edges[iedg];
491  for (int sid=0; sid<edg->n_sides; sid++)
492  if ( el_is_local(edg->side(sid)->element().idx()) )
493  {
494  is_edge_local = true;
495  break;
496  }
497  if (is_edge_local)
498  edg_4_loc.push_back(iedg);
499  }
500 
501  // create local array of neighbours
502  for (unsigned int inb=0; inb<mesh_->n_vb_neighbours(); inb++)
503  {
504  Neighbour *nb = &mesh_->vb_neighbours_[inb];
505  if ( el_is_local(nb->element().idx())
506  || el_is_local(nb->side()->element().idx()) )
507  nb_4_loc.push_back(inb);
508  }
509 
510  // create array of local nodes
511  std::vector<bool> node_is_local(mesh_->tree->n_nodes(), false);
512  for (auto cell : this->own_range())
513  {
514  unsigned int obj_idx = mesh_->tree->obj_4_el()[cell.elm_idx()];
515  for (unsigned int nid=0; nid<cell.elm()->n_nodes(); nid++)
516  node_is_local[mesh_->tree->objects(cell.dim())[obj_idx].nodes[nid]] = true;
517  }
518 
519  // init global to local element map with locally owned elements (later add ghost elements)
520  for ( unsigned int iel = 0; iel < el_ds_->lsize(); iel++ )
522 
523  // create array of local ghost cells
524  for ( auto cell : mesh_->elements_range() )
525  {
526  if (cell.proc() != el_ds_->myp())
527  {
528  bool has_local_node = false;
529  unsigned int obj_idx = mesh_->tree->obj_4_el()[cell.idx()];
530  for (unsigned int nid=0; nid<cell->n_nodes(); nid++)
531  if (node_is_local[mesh_->tree->objects(cell->dim())[obj_idx].nodes[nid]])
532  {
533  has_local_node = true;
534  break;
535  }
536  if (has_local_node)
537  {
538  ghost_4_loc.push_back(cell.idx());
539  ghost_proc.insert(cell.proc());
540  ghost_proc_el[cell.proc()].push_back(cell.idx());
541  global_to_local_el_idx_[cell.idx()] = el_ds_->lsize() - 1 + ghost_4_loc.size();
542  }
543  }
544  }
545  for (auto nb : nb_4_loc)
546  {
547  auto cell = mesh_->vb_neighbours_[nb].element();
548  if (!el_is_local(cell.idx()) && find(ghost_4_loc.begin(), ghost_4_loc.end(), cell.idx()) == ghost_4_loc.end())
549  {
550  ghost_4_loc.push_back(cell.idx());
551  ghost_proc.insert(cell.proc());
552  ghost_proc_el[cell.proc()].push_back(cell.idx());
553  global_to_local_el_idx_[cell.idx()] = el_ds_->lsize() - 1 + ghost_4_loc.size();
554  }
555  cell = mesh_->vb_neighbours_[nb].side()->element();
556  if (!el_is_local(cell.idx()) && find(ghost_4_loc.begin(), ghost_4_loc.end(), cell.idx()) == ghost_4_loc.end())
557  {
558  ghost_4_loc.push_back(cell.idx());
559  ghost_proc.insert(cell.proc());
560  ghost_proc_el[cell.proc()].push_back(cell.idx());
561  global_to_local_el_idx_[cell.idx()] = el_ds_->lsize() - 1 + ghost_4_loc.size();
562  }
563  }
564 }
565 
566 
567 bool DOFHandlerMultiDim::el_is_local(int index) const
568 {
569  return el_ds_->is_local(mesh_->get_row_4_el()[index]);
570 }
571 
572 
573 std::size_t DOFHandlerMultiDim::hash() const {
574  return this->n_global_dofs_;
575 }
576 
577 
579  auto bgn_it = make_iter<DHCellAccessor>( DHCellAccessor(this, 0) );
580  auto end_it = make_iter<DHCellAccessor>( DHCellAccessor(this, el_ds_->lsize()) );
581  return Range<DHCellAccessor>(bgn_it, end_it);
582 }
583 
584 
586  auto bgn_it = make_iter<DHCellAccessor>( DHCellAccessor(this, 0) );
587  auto end_it = make_iter<DHCellAccessor>( DHCellAccessor(this, el_ds_->lsize()+ghost_4_loc.size()) );
588  return Range<DHCellAccessor>(bgn_it, end_it);
589 }
590 
591 
593  auto bgn_it = make_iter<DHCellAccessor>( DHCellAccessor(this, el_ds_->lsize()) );
594  auto end_it = make_iter<DHCellAccessor>( DHCellAccessor(this, el_ds_->lsize()+ghost_4_loc.size()) );
595  return Range<DHCellAccessor>(bgn_it, end_it);
596 }
597 
598 
600  auto map_it = global_to_local_el_idx_.find((LongIdx)elm_idx); // find in global to local map
601  ASSERT( map_it != global_to_local_el_idx_.end() )(elm_idx).error("DH accessor can be create only for own or ghost elements!\n");
602  return DHCellAccessor(this, map_it->second);
603 }
604 
605 
607  stringstream s;
609 
610  s << "DOFHandlerMultiDim structure:" << endl;
611  s << "- is parallel: " << (is_parallel_?"true":"false") << endl;
612  s << "- proc id: " << el_ds_->myp() << endl;
613  s << "- global number of dofs: " << n_global_dofs_ << endl;
614  s << "- number of locally owned cells: " << el_ds_->lsize() << endl;
615  s << "- number of ghost cells: " << ghost_4_loc.size() << endl;
616  s << "- dofs on locally owned cells:" << endl;
617 
618  for (auto cell : own_range())
619  {
620  auto ndofs = cell.get_dof_indices(dofs);
621  s << "-- cell " << cell.elm().index() << ": ";
622  for (unsigned int idof=0; idof<ndofs; idof++) s << dofs[idof] << " "; s << endl;
623  }
624  s << "- dofs on ghost cells:" << endl;
625  for (auto cell : ghost_range())
626  {
627  auto ndofs = cell.get_dof_indices(dofs);
628  s << "-- cell " << cell.elm().index() << ": ";
629  for (unsigned int idof=0; idof<ndofs; idof++) s << dofs[idof] << " "; s << endl;
630  }
631  s << "- locally owned dofs (" << lsize_ << "): ";
632  for (unsigned int i=0; i<lsize_; i++) s << local_to_global_dof_idx_[i] << " "; s << endl;
633  s << "- ghost dofs (" << local_to_global_dof_idx_.size() - lsize_ << "): ";
634  for (unsigned int i=lsize_; i<local_to_global_dof_idx_.size(); i++) s << local_to_global_dof_idx_[i] << " "; s << endl;
635  s << "- global-to-local-cell map:" << endl;
636  for (auto cell : global_to_local_el_idx_) s << "-- " << cell.first << " -> " << cell.second << " " << endl;
637  s << endl;
638 
639  printf("%s", s.str().c_str());
640 }
641 
642 
643 
644 
645 
646 
647 
648 
649 
650 
651 
652 
653 
654 
655 
656 
657 
DOFHandlerMultiDim::ghost_range
Range< DHCellAccessor > ghost_range() const
Returns range over ghosts DOF handler cells.
Definition: dofhandler.cc:592
DOFHandlerMultiDim::dof_indices
std::vector< LongIdx > dof_indices
Dof numbers on local and ghost elements.
Definition: dofhandler.hh:388
Mesh::get_el_ds
Distribution * get_el_ds() const
Definition: mesh.h:168
DOFHandlerMultiDim::local_to_global_dof_idx_
std::vector< LongIdx > local_to_global_dof_idx_
Maps local and ghost dof indices to global ones.
Definition: dofhandler.hh:396
Distribution::np
unsigned int np() const
get num of processors
Definition: distribution.hh:105
DOFHandlerMultiDim::INVALID_DOF
static const int INVALID_DOF
Definition: dofhandler.hh:350
Dof::n_face_idx
unsigned int n_face_idx
Index of n-face to which the dof is associated.
Definition: finite_element.hh:99
DOFHandlerMultiDim::global_to_local_el_idx_
std::unordered_map< LongIdx, LongIdx > global_to_local_el_idx_
Maps global element index into local/ghost index (obsolete).
Definition: dofhandler.hh:401
DOFHandlerMultiDim::ghost_proc
set< unsigned int > ghost_proc
Processors of ghost elements.
Definition: dofhandler.hh:416
Distribution::lsize
unsigned int lsize(int proc) const
get local size
Definition: distribution.hh:115
neighbours.h
DOFHandlerMultiDim::ds
std::shared_ptr< DiscreteSpace > ds() const
Return pointer to discrete space for which the handler distributes dofs.
Definition: dofhandler.hh:265
Distribution::myp
unsigned int myp() const
get my processor
Definition: distribution.hh:107
ASSERT
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
Definition: asserts.hh:346
LongIdx
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
Definition: long_idx.hh:22
Mesh::elements_range
virtual Range< ElementAccessor< 3 > > elements_range() const
Returns range of bulk elements.
Definition: mesh.cc:1116
distribution.hh
Support classes for parallel programing.
DOFHandlerMultiDim::DOFHandlerMultiDim
DOFHandlerMultiDim(Mesh &_mesh, bool make_elem_part=true)
Constructor.
Definition: dofhandler.cc:42
Mesh::tree
DuplicateNodes * tree
Definition: mesh.h:277
DuplicateNodes::obj_4_el
const std::vector< unsigned int > & obj_4_el() const
Definition: duplicate_nodes.h:110
Neighbour::side
SideIter side()
Definition: neighbours.h:147
MPI_MAX
#define MPI_MAX
Definition: mpi.h:197
DOFHandlerMultiDim::ghost_4_loc
vector< LongIdx > ghost_4_loc
Indices of ghost cells (neighbouring with local elements).
Definition: dofhandler.hh:413
std::vector< LongIdx >
Mesh::n_edges
unsigned int n_edges() const
Definition: mesh.h:141
ElementAccessor< 3 >
MPI_Send
#define MPI_Send(buf, count, datatype, dest, tag, comm)
Definition: mpi.h:263
DOFHandlerBase::mesh_
Mesh * mesh_
Pointer to the mesh to which the dof handler is associated.
Definition: dofhandler.hh:127
DOFHandlerBase::~DOFHandlerBase
virtual ~DOFHandlerBase()
Destructor.
Definition: dofhandler.cc:34
Mesh::get_row_4_el
LongIdx * get_row_4_el() const
Definition: mesh.h:171
DHCellAccessor::dim
unsigned int dim() const
Return dimension of element appropriate to cell.
Definition: dh_cell_accessor.hh:94
DOFHandlerMultiDim::el_ds_
Distribution * el_ds_
Distribution of elements.
Definition: dofhandler.hh:404
dofhandler.hh
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
DOFHandlerMultiDim::local_range
Range< DHCellAccessor > local_range() const
Returns range over own and ghost cells of DOF handler.
Definition: dofhandler.cc:585
DOFHandlerMultiDim::ghost_proc_el
map< unsigned int, vector< LongIdx > > ghost_proc_el
Arrays of ghost cells for each neighbouring processor.
Definition: dofhandler.hh:419
duplicate_nodes.h
Neighbour
Definition: neighbours.h:119
Mesh::n_vb_neighbours
unsigned int n_vb_neighbours() const
Definition: mesh.cc:236
DOFHandlerMultiDim::init_node_dof_starts
void init_node_dof_starts(std::vector< LongIdx > &node_dof_starts)
Initialize auxiliary vector of starting indices of nodal dofs.
Definition: dofhandler.cc:82
Dof::dim
unsigned int dim
Association to n-face of given dimension (point, line, triangle, tetrahedron.
Definition: finite_element.hh:96
DOFHandlerMultiDim::get_loc_dof_indices
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.
Definition: dofhandler.cc:464
DOFHandlerMultiDim::VALID_NODE
static const int VALID_NODE
Definition: dofhandler.hh:348
Distribution
Definition: distribution.hh:50
MPI_Recv
#define MPI_Recv(buf, count, datatype, source, tag, comm, status)
Definition: mpi.h:271
dh_cell_accessor.hh
accessors.hh
DOFHandlerMultiDim::dh_seq_
std::shared_ptr< DOFHandlerMultiDim > dh_seq_
Sequential dof handler associated to the current (parallel) one.
Definition: dofhandler.hh:360
Edge::n_sides
int n_sides
Definition: edges.h:36
DOFHandlerMultiDim::el_is_local
bool el_is_local(int index) const
Definition: dofhandler.cc:567
Edge::side
SideIter side(const unsigned int i) const
Definition: edges.h:31
DOFHandlerMultiDim::update_local_dofs
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.
Definition: dofhandler.cc:177
DHCellAccessor::local_idx
unsigned int local_idx() const
Return local index to element (index of DOF handler).
Definition: dh_cell_accessor.hh:53
DOFHandlerMultiDim::own_range
Range< DHCellAccessor > own_range() const
Returns range of DOF handler cells (only range of own without ghost cells)
Definition: dofhandler.cc:578
finite_element.hh
Abstract class for description of finite elements.
DOFHandlerMultiDim::distribute_dofs
void distribute_dofs(std::shared_ptr< DiscreteSpace > ds)
Distributes degrees of freedom on the mesh needed for the given discrete space.
Definition: dofhandler.cc:242
DOFHandlerMultiDim::receive_ghost_dofs
void receive_ghost_dofs(unsigned int proc, std::vector< LongIdx > &dofs)
Obtain dof numbers on ghost elements from other processor.
Definition: dofhandler.cc:126
DHCellAccessor::elm
const ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Definition: dh_cell_accessor.hh:66
DuplicateNodes::n_nodes
unsigned int n_nodes() const
Definition: duplicate_nodes.h:104
Neighbour::element
ElementAccessor< 3 > element()
Definition: neighbours.h:163
DOFHandlerMultiDim::create_vector
VectorMPI create_vector()
Allocates PETSc vector according to the dof distribution.
Definition: dofhandler.cc:437
DOFHandlerBase::n_global_dofs_
unsigned int n_global_dofs_
Number of global dofs assigned by the handler.
Definition: dofhandler.hh:109
DOFHandlerBase::loffset_
unsigned int loffset_
Index of the first dof on the local process.
Definition: dofhandler.hh:119
DOFHandlerMultiDim::sequential
std::shared_ptr< DOFHandlerMultiDim > sequential()
Returns sequential version of the current dof handler.
Definition: dofhandler.cc:54
Side::element
ElementAccessor< 3 > element() const
Definition: side_impl.hh:53
DOFHandlerMultiDim::INVALID_NODE
static const int INVALID_NODE
Definition: dofhandler.hh:347
mesh.h
DHCellAccessor::n_dofs
unsigned int n_dofs() const
Return number of dofs on given cell.
Definition: dh_cell_accessor.hh:332
DOFHandlerBase::distr
std::shared_ptr< Distribution > distr() const
Definition: dofhandler.hh:73
DHCellAccessor::elm_idx
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
Definition: dh_cell_accessor.hh:59
DOFHandlerMultiDim::DHCellAccessor
friend class DHCellAccessor
Definition: dofhandler.hh:273
partitioning.hh
MPI_SUM
#define MPI_SUM
Definition: mpi.h:196
Mesh::edges
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:261
DOFHandlerMultiDim::cell_accessor_from_element
const DHCellAccessor cell_accessor_from_element(unsigned int elm_idx) const
Return DHCellAccessor appropriate to ElementAccessor of given idx.
Definition: dofhandler.cc:599
Mesh
Definition: mesh.h:80
DOFHandlerMultiDim::ASSIGNED_NODE
static const int ASSIGNED_NODE
Definition: dofhandler.hh:349
MPI_STATUS_IGNORE
#define MPI_STATUS_IGNORE
Definition: mpi.h:203
DOFHandlerMultiDim::~DOFHandlerMultiDim
~DOFHandlerMultiDim() override
Destructor.
Definition: dofhandler.cc:476
DOFHandlerMultiDim::nb_4_loc
vector< LongIdx > nb_4_loc
Local neighbour index -> global neighbour index.
Definition: dofhandler.hh:410
Range
Range helper class.
Definition: mesh.h:52
DOFHandlerMultiDim::ds_
std::shared_ptr< DiscreteSpace > ds_
Pointer to the discrete space for which the handler distributes dofs.
Definition: dofhandler.hh:354
DOFHandlerMultiDim::get_dof_indices
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.
Definition: dofhandler.cc:451
DOFHandlerBase
Definition: dofhandler.hh:47
OLD_ASSERT
#define OLD_ASSERT(...)
Definition: global_defs.h:131
DHCellAccessor
Definition: dh_cell_accessor.hh:36
MPI_Allreduce
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
Definition: mpi.h:612
fmt::printf
void printf(BasicWriter< Char > &w, BasicCStringRef< Char > format, ArgList args)
Definition: printf.h:444
Edge
Definition: edges.h:26
DOFHandlerMultiDim::cell_starts
std::vector< LongIdx > cell_starts
Starting indices for element dofs.
Definition: dofhandler.hh:380
DOFHandlerMultiDim::print
void print() const
Output structure of dof handler.
Definition: dofhandler.cc:606
DuplicateNodes::objects
const std::vector< MeshObject > & objects(unsigned int dim) const
Definition: duplicate_nodes.h:107
DOFHandlerMultiDim::make_elem_partitioning
void make_elem_partitioning()
Prepare parallel distribution of elements, edges and neighbours.
Definition: dofhandler.cc:481
DOFHandlerBase::max_elem_dofs_
unsigned int max_elem_dofs_
Max. number of dofs per element.
Definition: dofhandler.hh:122
DOFHandlerMultiDim::create_sequential
void create_sequential()
Communicate local dof indices to all processors and create new sequential dof handler.
Definition: dofhandler.cc:349
DOFHandlerBase::lsize_
unsigned int lsize_
Number of dofs associated to local process.
Definition: dofhandler.hh:114
MPI_COMM_WORLD
#define MPI_COMM_WORLD
Definition: mpi.h:123
std::vector::data
T data
Definition: doxy_dummy_defs.hh:7
DOFHandlerMultiDim::scatter_to_seq_
std::shared_ptr< VecScatter > scatter_to_seq_
Scatter context for parallel to sequential vectors.
Definition: dofhandler.hh:363
ElementAccessor::idx
unsigned int idx() const
Return local idx of element in boundary / bulk part of element vector.
Definition: accessors.hh:111
DOFHandlerMultiDim::sequential_scatter
std::shared_ptr< VecScatter > sequential_scatter()
Returns scatter context from parallel to sequential vectors.
Definition: dofhandler.cc:61
MPI_Allgatherv
#define MPI_Allgatherv(sendbuf, sendcount, sendtype, recvbuf, recvcounts, displs, recvtype, comm)
Definition: mpi.h:587
long_idx.hh
DOFHandlerMultiDim::is_parallel_
bool is_parallel_
Indicator for parallel/sequential dof handler.
Definition: dofhandler.hh:357
MPI_LONG_IDX
#define MPI_LONG_IDX
Definition: long_idx.hh:24
VectorMPI
Definition: vector_mpi.hh:42
Mesh::n_elements
virtual unsigned int n_elements(bool boundary=false) const
Returns count of boundary or bulk elements.
Definition: mesh.h:355
MPI_UNSIGNED
#define MPI_UNSIGNED
Definition: mpi.h:165
DOFHandlerMultiDim::hash
std::size_t hash() const override
Definition: dofhandler.cc:573
DOFHandlerMultiDim::send_ghost_dofs
void send_ghost_dofs(unsigned int proc)
Send dof numbers to other processor.
Definition: dofhandler.cc:147
DHCellAccessor::cell_dof
const Dof & cell_dof(unsigned int idof) const
Return dof on a given cell.
Definition: dh_cell_accessor.hh:349
DOFHandlerMultiDim::init_node_status
void init_node_status(std::vector< short int > &node_status)
Initialize node_status.
Definition: dofhandler.cc:99
Distribution::is_local
bool is_local(unsigned int idx) const
identify local index
Definition: distribution.hh:120
DOFHandlerBase::dof_ds_
std::shared_ptr< Distribution > dof_ds_
Distribution of dofs associated to local process.
Definition: dofhandler.hh:132
Mesh::vb_neighbours_
vector< Neighbour > vb_neighbours_
Definition: mesh.h:282
DOFHandlerMultiDim::init_cell_starts
void init_cell_starts()
Initialize vector of starting indices for elements.
Definition: dofhandler.cc:68
Element::n_nodes
unsigned int n_nodes() const
Definition: elements.h:129
range_wrapper.hh
Implementation of range helper class.
DOFHandlerMultiDim::edg_4_loc
vector< LongIdx > edg_4_loc
Local edge index -> global edge index.
Definition: dofhandler.hh:407
Mesh::get_el_4_loc
LongIdx * get_el_4_loc() const
Definition: mesh.h:174