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