Flow123d  release_2.2.0-914-gf1a3a4f
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 "mesh/partitioning.hh"
22 #include "la/distribution.hh"
23 
24 
25 
26 
27 
28 
29 
30 //template<unsigned int dim, unsigned int spacedim> inline
31 //DOFHandler<dim,spacedim>::DOFHandler(Mesh & _mesh)
32 //: DOFHandlerBase(_mesh),
33 // finite_element(0)
34 //{
35 // object_dofs = new int**[mesh->n_elements()];
36 // for (int i=0; i<mesh->n_elements(); i++)
37 // object_dofs[i] = NULL;
38 //}
39 //
40 //
41 //template<unsigned int dim, unsigned int spacedim> inline
42 //DOFHandler<dim,spacedim>::~DOFHandler()
43 //{
44 // for (int i=0; i<mesh->n_elements(); i++)
45 // if (object_dofs[i] != NULL)
46 // {
47 // for (int j=0; j<mesh->element[i].dim(); j++)
48 // if (object_dofs[i][j] != NULL)
49 // delete[] object_dofs[i][j];
50 //
51 // delete[] object_dofs[i];
52 // }
53 // delete[] object_dofs;
54 //}
55 //
56 //
57 //template<unsigned int dim, unsigned int spacedim> inline
58 //void DOFHandler<dim,spacedim>::distribute_dofs(FiniteElement<dim,spacedim> & fe, const unsigned int offset)
59 //{
60 // // First check if dofs are already distributed.
61 // OLD_ASSERT(finite_element == 0, "Attempt to distribute DOFs multiple times!");
62 //
63 // unsigned int next_free_dof = offset;
64 // unsigned int n_obj_dofs[dim+1];
65 //
66 // finite_element = &fe;
67 // global_dof_offset = offset;
68 //
69 // for (unsigned int dm=0; dm <= dim; dm++)
70 // {
71 // n_obj_dofs[dm] = 0;
72 // for (unsigned int m=0; m<dof_multiplicities.size(); m++)
73 // n_obj_dofs[dm] += fe.n_object_dofs(dm, dof_multiplicities[m])*dof_multiplicities[m];
74 // }
75 //
76 // // Broadcast partition of elements to all processes.
77 // IdxInt *loc_part;
78 // int myp = mesh->get_part()->get_init_distr()->myp();
79 // if (myp == 0)
80 // {
81 // loc_part = (int*)mesh->get_part()->get_loc_part();
82 // }
83 // else
84 // {
85 // loc_part = new IdxInt[mesh->n_elements()];
86 // }
87 // MPI_Bcast(loc_part, mesh->n_elements(), MPI_INT, 0, mesh->get_part()->get_init_distr()->get_comm());
88 //
89 // // Distribute element dofs.
90 // // First we distribute dofs on elements associated to process 0 and so on.
91 // for (int proc=0; proc<mesh->get_part()->get_init_distr()->np(); proc++)
92 // {
93 // if (proc == myp)
94 // loffset_ = next_free_dof;
95 //
96 // FOR_ELEMENTS(mesh, cell)
97 // {
98 // if (cell->dim() != dim) continue;
99 //
100 // if (loc_part[cell.index()] != proc) continue;
101 //
102 // // distribute dofs
103 // // TODO: For the moment we distribute only dofs associated to the cell
104 // // In the future we want to distribute dofs on vertices, lines,
105 // // and triangles as well.
106 // object_dofs[cell.index()] = new int*[dim+1];
107 // for (int i=0; i<dim+1; i++)
108 // object_dofs[cell.index()][i] = NULL;
109 // object_dofs[cell.index()][dim] = new int[n_obj_dofs[dim]];
110 //
111 // for (unsigned int i=0; i<n_obj_dofs[dim]; i++)
112 // object_dofs[cell.index()][dim][i] = next_free_dof++;
113 // }
114 //
115 // if (proc == myp)
116 // lsize_ = next_free_dof - loffset_;
117 // }
118 //
119 // // Finally we free the unused array loc_part.
120 // if (mesh->get_part()->get_init_distr()->myp() != 0)
121 // delete[] loc_part;
122 //
123 //// FOR_ELEMENTS(mesh,cell)
124 //// {
125 //// // skip cells of different dimension
126 //// if (cell->dim() != dim) continue;
127 ////
128 //// // distribute dofs
129 //// // TODO: For the moment we distribute only dofs associated to the cell
130 //// // In the future we want to distribute dofs on vertices, lines,
131 //// // and triangles as well.
132 //// object_dofs[dim][cell] = new int[n_obj_dofs[dim]];
133 //// for (unsigned int i=0; i<n_obj_dofs[dim]; i++)
134 //// object_dofs[dim][cell][i] = next_free_dof++;
135 //// }
136 //
137 //
138 //// FOR_ELEMENTS(mesh,cell)
139 //// {
140 //// // skip cells of different dimension
141 //// if (cell->dim != dim) continue;
142 ////
143 //// // distribute dofs
144 //// for (int dm=0; dm<=dim; dm++)
145 //// {
146 //// for (int i=0; i<cell->n_sides_by_dim(dm); i++)
147 //// {
148 //// void *side = cell->side_by_dim(dm, i);
149 //// // check if node has already assigned dofs, otherwise
150 //// // distribute
151 //// if (object_dofs[dm].find(side) == object_dofs[dm].end())
152 //// {
153 //// object_dofs[dm][side] = new int[n_obj_dofs[dm]];
154 //// for (int i=0; i<n_obj_dofs[dm]; i++)
155 //// object_dofs[dm][side][i] = next_free_dof++;
156 //// }
157 //// }
158 //// }
159 //// }
160 //
161 // n_dofs = next_free_dof - offset;
162 //}
163 //
164 //
165 //
166 //template<unsigned int dim, unsigned int spacedim> inline
167 //const unsigned int DOFHandler<dim,spacedim>::n_local_dofs()
168 //{
169 // return finite_element->n_dofs();
170 //}
171 //
172 //
173 //
174 //template<unsigned int dim, unsigned int spacedim>
175 //void DOFHandler<dim,spacedim>::get_dof_indices(const CellIterator &cell, unsigned int indices[])
176 //{
177 //// void *side;
178 //// unsigned int offset, pid;
179 //
180 // for (unsigned int k=0; k<finite_element->n_object_dofs(dim,DOF_SINGLE); k++)
181 // indices[k] = object_dofs[cell.index()][dim][k];
182 //
183 //// indices.clear();
184 ////
185 //// get_object_dof_indices<0>(cell, indices);
186 //// get_object_dof_indices<1>(cell, indices);
187 //// get_object_dof_indices<2>(cell, indices);
188 //// get_object_dof_indices<3>(cell, indices);
189 //}
190 //
191 ////template<unsigned int dim> template<unsigned int obj_dim> inline void DOFHandler<dim>::get_object_dof_indices(const CellIterator &cell, unsigned int indices[])
192 ////{
193 // // TODO: implement for lower dimensional objects
194 //
195 //// void *side;
196 //// unsigned int offset, pid;
197 ////
198 //// // loop over cell points/lines/triangles/tetrahedra
199 //// for (int i=0; i<n_simplex_objects<dim>(obj_dim); i++)
200 //// {
201 //// side = cell->side_by_dim(obj_dim,i);
202 //// pid = permutation_id<dim,obj_dim>(cell,i);
203 //// offset = 0;
204 //// // loop over dof multiplicities (single dofs, pairs, triples, sextuples)
205 //// for (vector<unsigned int>::iterator m=dof_multiplicities.begin(); m!=dof_multiplicities.end(); m++)
206 //// {
207 //// // loop over particular single dofs/dof pairs/triples/sextuples
208 //// for (int j=0; j<finite_element.n_object_dofs(obj_dim,*m); j++)
209 //// {
210 //// // loop over particular dofs (the single dof/2 dofs in the pair etc.)
211 //// for (int k=0; k<*m; k++)
212 //// indices.push_back(object_dofs[obj_dim][side][offset+Simplex<obj_dim>::pair_permutations[pid][k]]);
213 ////
214 //// offset += *m;
215 //// }
216 //// }
217 //// }
218 ////}
219 //
220 //template<unsigned int dim, unsigned int spacedim> inline
221 //void DOFHandler<dim,spacedim>::get_dof_values(const CellIterator &cell, const Vec &values, double local_values[])
222 //{
223 // unsigned int indices[finite_element->n_dofs()];
224 //
225 // get_dof_indices(cell, indices);
226 // VecGetValues(values, finite_element->n_dofs(), (PetscInt *)indices, local_values);
227 //}
228 
229 
230 
232  : DOFHandlerBase(_mesh),
233  fe1d_(0),
234  fe2d_(0),
235  fe3d_(0)
236 {
237  object_dofs = new IdxInt**[mesh_->n_elements()];
238  for (unsigned int i=0; i<mesh_->n_elements(); i++)
239  object_dofs[i] = NULL;
240 
242 }
243 
246  const unsigned int offset)
247 {
248  // First check if dofs are already distributed.
249  OLD_ASSERT((fe1d_ == 0) && (fe2d_ == 0) && (fe3d_ == 0), "Attempt to distribute DOFs multiple times!");
250 
251  unsigned int next_free_dof = offset;
252  // n_obj_dofs[element_dim][general_face_dim]
253  // for every dim of ref. element, number of DoFs on its generalized face (face, edge, vertex)
254  unsigned int n_obj_dofs[4][4];
255 
256 
257  fe1d_ = &fe1d;
258  fe2d_ = &fe2d;
259  fe3d_ = &fe3d;
261  max_elem_dofs_ = max(fe1d_->n_dofs(), max(fe2d_->n_dofs(), fe3d_->n_dofs()));
262 
263  for (unsigned int dm=0; dm <= 3; dm++)
264  for (unsigned int dd=1; dd <= 3; dd++)
265  n_obj_dofs[dd][dm] = 0;
266  for (unsigned int d=0; d<fe1d.n_dofs(); d++)
267  n_obj_dofs[1][fe1d.dof(d).dim]++;
268  for (unsigned int d=0; d<fe2d.n_dofs(); d++)
269  n_obj_dofs[2][fe2d.dof(d).dim]++;
270  for (unsigned int d=0; d<fe3d.n_dofs(); d++)
271  n_obj_dofs[3][fe3d.dof(d).dim]++;
272 
273  // Broadcast partition of elements to all processes.
274  IdxInt *loc_part;
275  unsigned int myp = mesh_->get_part()->get_init_distr()->myp();
276  if (myp == 0)
277  {
278  loc_part = (IdxInt*)mesh_->get_part()->get_loc_part();
279  }
280  else
281  {
282  loc_part = new int[mesh_->n_elements()];
283  }
285 
286  // Distribute element dofs.
287  // First we distribute dofs on elements associated to process 0 and so on.
288  for (unsigned int proc=0; proc<mesh_->get_part()->get_init_distr()->np(); proc++)
289  {
290  if (proc == myp)
291  loffset_ = next_free_dof;
292 
293  FOR_ELEMENTS(mesh_, cell)
294  {
295  if (loc_part[cell.index()] != (int)proc) continue;
296 
297  unsigned int dim = cell->dim();
298 
299  // distribute dofs
300  // TODO: For the moment we distribute only dofs associated to the cell
301  // In the future we want to distribute dofs on vertices, lines,
302  // and triangles as well.
303  object_dofs[cell.index()] = new IdxInt*[dim+1];
304  for (unsigned int i=0; i<dim+1; i++)
305  object_dofs[cell.index()][i] = NULL;
306  object_dofs[cell.index()][dim] = new int[n_obj_dofs[dim][dim]];
307 
308  for (unsigned int i=0; i<n_obj_dofs[dim][dim]; i++)
309  object_dofs[cell.index()][dim][i] = next_free_dof++;
310  }
311 
312  if (proc == myp) {
313  lsize_ = next_free_dof - loffset_;
314  ds_ = new Distribution(lsize_, PETSC_COMM_WORLD);
315  }
316  }
317 
318  // Finally we free the unused array loc_part.
319  if (mesh_->get_part()->get_init_distr()->myp() != 0)
320  delete[] loc_part;
321 
322  n_dofs = next_free_dof - offset;
323 }
324 
326 {
327  unsigned int dim = cell->dim();
328  unsigned int n_objects_dofs;
329 
330  switch (dim)
331  {
332  case 1:
333  n_objects_dofs = fe1d_->n_dofs();
334  break;
335  case 2:
336  n_objects_dofs = fe2d_->n_dofs();
337  break;
338  case 3:
339  n_objects_dofs = fe3d_->n_dofs();
340  break;
341  }
342 
343  ASSERT_LE(n_objects_dofs, indices.size()).error();
344 
345  for (unsigned int k = 0; k < n_objects_dofs; k++)
346  indices[k] = object_dofs[cell.index()][dim][k];
347 
348  return n_objects_dofs;
349 }
350 
352 {
353  /*
354  * TODO: This method is currently wrong (the shift by loc_offset_ is not enough for
355  * local <-> global mapping. Only local to global mapping should be used in future.
356  */
357  unsigned int dim = cell->dim();
358  unsigned int n_objects_dofs;
359 
360  switch (dim)
361  {
362  case 1:
363  n_objects_dofs = fe1d_->n_dofs();
364  break;
365  case 2:
366  n_objects_dofs = fe2d_->n_dofs();
367  break;
368  case 3:
369  n_objects_dofs = fe3d_->n_dofs();
370  break;
371  }
372 
373  ASSERT_LE(n_objects_dofs, indices.size()).error();
374 
375  for (unsigned int k = 0; k < n_objects_dofs; k++)
376  indices[k] = object_dofs[cell.index()][dim][k] - loffset_;
377 
378  return n_objects_dofs;
379 }
380 
381 void DOFHandlerMultiDim::get_dof_values(const CellIterator &cell, const Vec &values, double local_values[]) const
382 {
383  int ndofs=0;
384 
385  switch (cell->dim())
386  {
387  case 1:
388  ndofs = fe1d_->n_dofs();
389  break;
390  case 2:
391  ndofs = fe2d_->n_dofs();
392  break;
393  case 3:
394  ndofs = fe3d_->n_dofs();
395  break;
396  }
397 
398  std::vector<IdxInt> indices(ndofs);
399 
400  get_dof_indices(cell, indices);
401  VecGetValues(values, ndofs, (PetscInt *)indices[0], local_values);
402 }
403 
405 {
406  for (ElementFullIter elem=mesh_->element.begin(); elem!=mesh_->element.end(); ++elem)
407  if (object_dofs[elem.index()] != NULL)
408  {
409  for (unsigned int j=0; j<=elem->dim(); j++)
410  if (object_dofs[elem.index()][j] != NULL)
411  delete[] object_dofs[elem.index()][j];
412 
413  delete[] object_dofs[elem.index()];
414  }
415  delete[] object_dofs;
416 }
417 
418 
419 
421 {
422  // create local arrays of elements
423  el_ds_ = mesh_->get_el_ds();
426 
427  // create local array of edges
428  for (unsigned int iedg=0; iedg<mesh_->n_edges(); iedg++)
429  {
430  bool is_edge_local = false;
431  Edge *edg = &mesh_->edges[iedg];
432  for (int sid=0; sid<edg->n_sides; sid++)
433  if (el_is_local(edg->side(sid)->element().index()))
434  {
435  is_edge_local = true;
436  break;
437  }
438  if (is_edge_local)
439  edg_4_loc.push_back(iedg);
440  }
441 
442  // create local array of neighbours
443  for (unsigned int inb=0; inb<mesh_->n_vb_neighbours(); inb++)
444  {
445  Neighbour *nb = &mesh_->vb_neighbours_[inb];
446  if (el_is_local(mesh_->element.index(nb->element()))
447  || el_is_local(nb->side()->element().index()))
448  nb_4_loc.push_back(inb);
449  }
450 }
451 
452 
453 bool DOFHandlerMultiDim::el_is_local(int index) const
454 {
455  return el_ds_->is_local(row_4_el[index]);
456 }
457 
458 
459 std::size_t DOFHandlerMultiDim::hash() const {
460  return this->n_dofs;
461 }
462 
463 
464 
465 
466 
467 template<> FiniteElement<1,3> *DOFHandlerMultiDim::fe<1>() const { return fe1d_; }
468 template<> FiniteElement<2,3> *DOFHandlerMultiDim::fe<2>() const { return fe2d_; }
469 template<> FiniteElement<3,3> *DOFHandlerMultiDim::fe<3>() const { return fe3d_; }
470 
471 
472 //template class DOFHandler<0,3>;
473 //template class DOFHandler<1,3>;
474 //template class DOFHandler<2,3>;
475 //template class DOFHandler<3,3>;
int IdxInt
Define integers that are indices into large arrays (elements, nodes, dofs etc.)
Definition: mesh.h:85
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
vector< IdxInt > nb_4_loc
Local neighbour index -> global neighbour index.
Definition: dofhandler.hh:413
DOFHandlerMultiDim(Mesh &_mesh)
Constructor.
Definition: dofhandler.cc:231
IdxInt *** object_dofs
Number of dofs associated to geometrical entities.
Definition: dofhandler.hh:399
FiniteElement< 3, 3 > * fe3d_
Definition: dofhandler.hh:390
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:480
FiniteElement< 1, 3 > * fe1d_
Pointer to the finite element class for which the handler distributes dofs.
Definition: dofhandler.hh:388
IdxInt * get_el_4_loc() const
Definition: mesh.h:195
Distribution * el_ds_
Distribution of elements.
Definition: dofhandler.hh:407
unsigned int dim
Association to n-face of given dimension (point, line, triangle, tetrahedron.
Definition: mesh.h:99
void distribute_dofs(FiniteElement< 1, 3 > &fe1d, FiniteElement< 2, 3 > &fe2d, FiniteElement< 3, 3 > &fe3d, const unsigned int offset=0)
Distributes degrees of freedom on the mesh needed for the given finite elements.
Definition: dofhandler.cc:244
vector< IdxInt > edg_4_loc
Local edge index -> global edge index.
Definition: dofhandler.hh:410
#define ASSERT_LE(a, b)
Definition of comparative assert macro (Less or Equal)
Definition: asserts.hh:303
int index() const
Definition: sys_vector.hh:78
Dof dof(unsigned int i) const
Returns i -th degree of freedom.
int n_sides
Definition: edges.h:36
unsigned int global_dof_offset
Index of first global dof.
Definition: dofhandler.hh:127
Definition: edges.h:26
unsigned int get_dof_indices(const CellIterator &cell, std::vector< IdxInt > &indices) const override
Returns the global indices of dofs associated to the cell.
Definition: dofhandler.cc:325
Partitioning * get_part()
Definition: mesh.cc:218
unsigned int lsize_
Number of dofs associated to local process.
Definition: dofhandler.hh:148
#define OLD_ASSERT(...)
Definition: global_defs.h:131
unsigned int n_vb_neighbours() const
Definition: mesh.h:208
unsigned int get_loc_dof_indices(const CellIterator &cell, std::vector< IdxInt > &indices) const override
Returns the indices of dofs associated to the cell on the local process.
Definition: dofhandler.cc:351
unsigned int n_elements() const
Definition: mesh.h:156
const unsigned int offset() const
Returns the number of the first global dof handled by this DOFHandler.
Definition: dofhandler.hh:68
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:161
SideIter side()
unsigned int index(const T *pointer) const
Definition: sys_vector.hh:373
unsigned int max_elem_dofs_
Max. number of dofs per element.
Definition: dofhandler.hh:156
FullIter begin()
Definition: sys_vector.hh:383
unsigned int np() const
get num of processors
Distribution * get_el_ds() const
Definition: mesh.h:189
void make_elem_partitioning()
Prepare parallel distribution of elements, edges and neighbours.
Definition: dofhandler.cc:420
Distribution * ds_
Distribution of dofs associated to local process.
Definition: dofhandler.hh:166
ElementIter element()
unsigned int myp() const
get my processor
Support classes for parallel programing.
FiniteElement< 2, 3 > * fe2d_
Definition: dofhandler.hh:389
vector< Neighbour > vb_neighbours_
Definition: mesh.h:290
unsigned int loffset_
Index of the first dof on the local process.
Definition: dofhandler.hh:153
ElementFullIter element() const
Definition: side_impl.hh:52
#define MPI_INT
Definition: mpi.h:160
IdxInt * row_4_el
Global element index -> index according to partitioning.
Definition: dofhandler.hh:403
~DOFHandlerMultiDim() override
Destructor.
Definition: dofhandler.cc:404
MPI_Comm get_comm() const
Returns communicator.
std::vector< Edge > edges
Vector of MH edges, this should not be part of the geometrical mesh.
Definition: mesh.h:271
unsigned int n_dofs
Number of global dofs assigned by the handler.
Definition: dofhandler.hh:143
const Distribution * get_init_distr() const
Definition: partitioning.cc:74
#define MPI_Bcast(buffer, count, datatype, root, comm)
Definition: mpi.h:534
const unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
const IdxInt * get_loc_part() const
Definition: partitioning.cc:81
Abstract class for description of finite elements.
void get_dof_values(const CellIterator &cell, const Vec &values, double local_values[]) const override
Returns the dof values associated to the cell.
Definition: dofhandler.cc:381
unsigned int n_edges() const
Definition: mesh.h:164
bool el_is_local(int index) const
Definition: dofhandler.cc:453
FullIter end()
Returns FullFullIterer of the fictions past the end element.
Definition: sys_vector.hh:387
IdxInt * el_4_loc
Local element index -> global element index.
Definition: dofhandler.hh:405
SideIter side(const unsigned int i) const
Definition: edges.h:31
IdxInt * get_row_4_el() const
Definition: mesh.h:192
std::size_t hash() const override
Definition: dofhandler.cc:459
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:260