Flow123d  jenkins-Flow123d-linux-release-multijob-282
dofhandler.cc
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4  *
5  * Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6  * especially for academic research:
7  * Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8  *
9  * This program is free software; you can redistribute it and/or modify it under the terms
10  * of the GNU General Public License version 3 as published by the Free Software Foundation.
11  *
12  * This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13  * without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14  * See the GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License along with this program; if not,
17  * write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18  *
19  *
20  * $Id$
21  * $Revision$
22  * $LastChangedBy$
23  * $LastChangedDate$
24  *
25  * @file
26  * @brief Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between local and global dofs.
27  * @author Jan Stebel
28  */
29 
30 
31 #include "fem/dofhandler.hh"
32 #include "fem/finite_element.hh"
33 #include "mesh/mesh.h"
34 #include "mesh/partitioning.hh"
35 #include "la/distribution.hh"
36 
37 
38 
39 
40 
41 
42 
43 //template<unsigned int dim, unsigned int spacedim> inline
44 //DOFHandler<dim,spacedim>::DOFHandler(Mesh & _mesh)
45 //: DOFHandlerBase(_mesh),
46 // finite_element(0)
47 //{
48 // object_dofs = new int**[mesh->n_elements()];
49 // for (int i=0; i<mesh->n_elements(); i++)
50 // object_dofs[i] = NULL;
51 //}
52 //
53 //
54 //template<unsigned int dim, unsigned int spacedim> inline
55 //DOFHandler<dim,spacedim>::~DOFHandler()
56 //{
57 // for (int i=0; i<mesh->n_elements(); i++)
58 // if (object_dofs[i] != NULL)
59 // {
60 // for (int j=0; j<mesh->element[i].dim(); j++)
61 // if (object_dofs[i][j] != NULL)
62 // delete[] object_dofs[i][j];
63 //
64 // delete[] object_dofs[i];
65 // }
66 // delete[] object_dofs;
67 //}
68 //
69 //
70 //template<unsigned int dim, unsigned int spacedim> inline
71 //void DOFHandler<dim,spacedim>::distribute_dofs(FiniteElement<dim,spacedim> & fe, const unsigned int offset)
72 //{
73 // // First check if dofs are already distributed.
74 // ASSERT(finite_element == 0, "Attempt to distribute DOFs multiple times!");
75 //
76 // unsigned int next_free_dof = offset;
77 // unsigned int n_obj_dofs[dim+1];
78 //
79 // finite_element = &fe;
80 // global_dof_offset = offset;
81 //
82 // for (unsigned int dm=0; dm <= dim; dm++)
83 // {
84 // n_obj_dofs[dm] = 0;
85 // for (unsigned int m=0; m<dof_multiplicities.size(); m++)
86 // n_obj_dofs[dm] += fe.n_object_dofs(dm, dof_multiplicities[m])*dof_multiplicities[m];
87 // }
88 //
89 // // Broadcast partition of elements to all processes.
90 // int *loc_part;
91 // int myp = mesh->get_part()->get_init_distr()->myp();
92 // if (myp == 0)
93 // {
94 // loc_part = (int*)mesh->get_part()->get_loc_part();
95 // }
96 // else
97 // {
98 // loc_part = new int[mesh->n_elements()];
99 // }
100 // MPI_Bcast(loc_part, mesh->n_elements(), MPI_INT, 0, mesh->get_part()->get_init_distr()->get_comm());
101 //
102 // // Distribute element dofs.
103 // // First we distribute dofs on elements associated to process 0 and so on.
104 // for (int proc=0; proc<mesh->get_part()->get_init_distr()->np(); proc++)
105 // {
106 // if (proc == myp)
107 // loffset_ = next_free_dof;
108 //
109 // FOR_ELEMENTS(mesh, cell)
110 // {
111 // if (cell->dim() != dim) continue;
112 //
113 // if (loc_part[cell.index()] != proc) continue;
114 //
115 // // distribute dofs
116 // // TODO: For the moment we distribute only dofs associated to the cell
117 // // In the future we want to distribute dofs on vertices, lines,
118 // // and triangles as well.
119 // object_dofs[cell.index()] = new int*[dim+1];
120 // for (int i=0; i<dim+1; i++)
121 // object_dofs[cell.index()][i] = NULL;
122 // object_dofs[cell.index()][dim] = new int[n_obj_dofs[dim]];
123 //
124 // for (unsigned int i=0; i<n_obj_dofs[dim]; i++)
125 // object_dofs[cell.index()][dim][i] = next_free_dof++;
126 // }
127 //
128 // if (proc == myp)
129 // lsize_ = next_free_dof - loffset_;
130 // }
131 //
132 // // Finally we free the unused array loc_part.
133 // if (mesh->get_part()->get_init_distr()->myp() != 0)
134 // delete[] loc_part;
135 //
136 //// FOR_ELEMENTS(mesh,cell)
137 //// {
138 //// // skip cells of different dimension
139 //// if (cell->dim() != dim) continue;
140 ////
141 //// // distribute dofs
142 //// // TODO: For the moment we distribute only dofs associated to the cell
143 //// // In the future we want to distribute dofs on vertices, lines,
144 //// // and triangles as well.
145 //// object_dofs[dim][cell] = new int[n_obj_dofs[dim]];
146 //// for (unsigned int i=0; i<n_obj_dofs[dim]; i++)
147 //// object_dofs[dim][cell][i] = next_free_dof++;
148 //// }
149 //
150 //
151 //// FOR_ELEMENTS(mesh,cell)
152 //// {
153 //// // skip cells of different dimension
154 //// if (cell->dim != dim) continue;
155 ////
156 //// // distribute dofs
157 //// for (int dm=0; dm<=dim; dm++)
158 //// {
159 //// for (int i=0; i<cell->n_sides_by_dim(dm); i++)
160 //// {
161 //// void *side = cell->side_by_dim(dm, i);
162 //// // check if node has already assigned dofs, otherwise
163 //// // distribute
164 //// if (object_dofs[dm].find(side) == object_dofs[dm].end())
165 //// {
166 //// object_dofs[dm][side] = new int[n_obj_dofs[dm]];
167 //// for (int i=0; i<n_obj_dofs[dm]; i++)
168 //// object_dofs[dm][side][i] = next_free_dof++;
169 //// }
170 //// }
171 //// }
172 //// }
173 //
174 // n_dofs = next_free_dof - offset;
175 //}
176 //
177 //
178 //
179 //template<unsigned int dim, unsigned int spacedim> inline
180 //const unsigned int DOFHandler<dim,spacedim>::n_local_dofs()
181 //{
182 // return finite_element->n_dofs();
183 //}
184 //
185 //
186 //
187 //template<unsigned int dim, unsigned int spacedim>
188 //void DOFHandler<dim,spacedim>::get_dof_indices(const CellIterator &cell, unsigned int indices[])
189 //{
190 //// void *side;
191 //// unsigned int offset, pid;
192 //
193 // for (unsigned int k=0; k<finite_element->n_object_dofs(dim,DOF_SINGLE); k++)
194 // indices[k] = object_dofs[cell.index()][dim][k];
195 //
196 //// indices.clear();
197 ////
198 //// get_object_dof_indices<0>(cell, indices);
199 //// get_object_dof_indices<1>(cell, indices);
200 //// get_object_dof_indices<2>(cell, indices);
201 //// get_object_dof_indices<3>(cell, indices);
202 //}
203 //
204 ////template<unsigned int dim> template<unsigned int obj_dim> inline void DOFHandler<dim>::get_object_dof_indices(const CellIterator &cell, unsigned int indices[])
205 ////{
206 // // TODO: implement for lower dimensional objects
207 //
208 //// void *side;
209 //// unsigned int offset, pid;
210 ////
211 //// // loop over cell points/lines/triangles/tetrahedra
212 //// for (int i=0; i<n_simplex_objects<dim>(obj_dim); i++)
213 //// {
214 //// side = cell->side_by_dim(obj_dim,i);
215 //// pid = permutation_id<dim,obj_dim>(cell,i);
216 //// offset = 0;
217 //// // loop over dof multiplicities (single dofs, pairs, triples, sextuples)
218 //// for (vector<unsigned int>::iterator m=dof_multiplicities.begin(); m!=dof_multiplicities.end(); m++)
219 //// {
220 //// // loop over particular single dofs/dof pairs/triples/sextuples
221 //// for (int j=0; j<finite_element.n_object_dofs(obj_dim,*m); j++)
222 //// {
223 //// // loop over particular dofs (the single dof/2 dofs in the pair etc.)
224 //// for (int k=0; k<*m; k++)
225 //// indices.push_back(object_dofs[obj_dim][side][offset+Simplex<obj_dim>::pair_permutations[pid][k]]);
226 ////
227 //// offset += *m;
228 //// }
229 //// }
230 //// }
231 ////}
232 //
233 //template<unsigned int dim, unsigned int spacedim> inline
234 //void DOFHandler<dim,spacedim>::get_dof_values(const CellIterator &cell, const Vec &values, double local_values[])
235 //{
236 // unsigned int indices[finite_element->n_dofs()];
237 //
238 // get_dof_indices(cell, indices);
239 // VecGetValues(values, finite_element->n_dofs(), (PetscInt *)indices, local_values);
240 //}
241 
242 
243 
245  : DOFHandlerBase(_mesh),
246  fe1d_(0),
247  fe2d_(0),
248  fe3d_(0)
249 {
250  object_dofs = new int**[mesh_->n_elements()];
251  for (unsigned int i=0; i<mesh_->n_elements(); i++)
252  object_dofs[i] = NULL;
253 
255 }
256 
259  const unsigned int offset)
260 {
261  // First check if dofs are already distributed.
262  ASSERT((fe1d_ == 0) && (fe2d_ == 0) && (fe3d_ == 0), "Attempt to distribute DOFs multiple times!");
263 
264  unsigned int next_free_dof = offset;
265  // n_obj_dofs[element_dim][general_face_dim]
266  // for every dim of ref. element, number of DoFs on its generalized face (face, edge, vertex)
267  unsigned int n_obj_dofs[4][4];
268 
269 
270  fe1d_ = &fe1d;
271  fe2d_ = &fe2d;
272  fe3d_ = &fe3d;
274 
275  for (unsigned int dm=0; dm <= 1; dm++)
276  {
277  n_obj_dofs[1][dm] = 0;
278  for (unsigned int m=0; m<dof_multiplicities.size(); m++)
279  n_obj_dofs[1][dm] += fe1d.n_object_dofs(dm, dof_multiplicities[m])*dof_multiplicities[m];
280  }
281  for (unsigned int dm=0; dm <= 2; dm++)
282  {
283  n_obj_dofs[2][dm] = 0;
284  for (unsigned int m=0; m<dof_multiplicities.size(); m++)
285  n_obj_dofs[2][dm] += fe2d.n_object_dofs(dm, dof_multiplicities[m])*dof_multiplicities[m];
286  }
287  for (unsigned int dm=0; dm <= 3; dm++)
288  {
289  n_obj_dofs[3][dm] = 0;
290  for (unsigned int m=0; m<dof_multiplicities.size(); m++)
291  n_obj_dofs[3][dm] += fe3d.n_object_dofs(dm, dof_multiplicities[m])*dof_multiplicities[m];
292  }
293 
294  // Broadcast partition of elements to all processes.
295  int *loc_part;
296  unsigned int myp = mesh_->get_part()->get_init_distr()->myp();
297  if (myp == 0)
298  {
299  loc_part = (int*)mesh_->get_part()->get_loc_part();
300  }
301  else
302  {
303  loc_part = new int[mesh_->n_elements()];
304  }
306 
307  // Distribute element dofs.
308  // First we distribute dofs on elements associated to process 0 and so on.
309  for (unsigned int proc=0; proc<mesh_->get_part()->get_init_distr()->np(); proc++)
310  {
311  if (proc == myp)
312  loffset_ = next_free_dof;
313 
314  FOR_ELEMENTS(mesh_, cell)
315  {
316  if (loc_part[cell.index()] != (int)proc) continue;
317 
318  unsigned int dim = cell->dim();
319 
320  // distribute dofs
321  // TODO: For the moment we distribute only dofs associated to the cell
322  // In the future we want to distribute dofs on vertices, lines,
323  // and triangles as well.
324  object_dofs[cell.index()] = new int*[dim+1];
325  for (unsigned int i=0; i<dim+1; i++)
326  object_dofs[cell.index()][i] = NULL;
327  object_dofs[cell.index()][dim] = new int[n_obj_dofs[dim][dim]];
328 
329  for (unsigned int i=0; i<n_obj_dofs[dim][dim]; i++)
330  object_dofs[cell.index()][dim][i] = next_free_dof++;
331  }
332 
333  if (proc == myp) {
334  lsize_ = next_free_dof - loffset_;
335  ds_ = new Distribution(lsize_, PETSC_COMM_WORLD);
336  }
337  }
338 
339  // Finally we free the unused array loc_part.
340  if (mesh_->get_part()->get_init_distr()->myp() != 0)
341  delete[] loc_part;
342 
343  n_dofs = next_free_dof - offset;
344 }
345 
346 void DOFHandlerMultiDim::get_dof_indices(const CellIterator &cell, unsigned int indices[]) const
347 {
348  unsigned int dim = cell->dim();
349  switch (dim)
350  {
351  case 1:
352  for (unsigned int k=0; k<fe1d_->n_object_dofs(dim,DOF_SINGLE); k++)
353  indices[k] = object_dofs[cell.index()][dim][k];
354  break;
355  case 2:
356  for (unsigned int k=0; k<fe2d_->n_object_dofs(dim,DOF_SINGLE); k++)
357  indices[k] = object_dofs[cell.index()][dim][k];
358  break;
359  case 3:
360  for (unsigned int k=0; k<fe3d_->n_object_dofs(dim,DOF_SINGLE); k++)
361  indices[k] = object_dofs[cell.index()][dim][k];
362  break;
363  }
364 }
365 
366 void DOFHandlerMultiDim::get_dof_values(const CellIterator &cell, const Vec &values, double local_values[]) const
367 {
368  int ndofs=0;
369 
370  switch (cell->dim())
371  {
372  case 1:
373  ndofs = fe1d_->n_dofs();
374  break;
375  case 2:
376  ndofs = fe2d_->n_dofs();
377  break;
378  case 3:
379  ndofs = fe3d_->n_dofs();
380  break;
381  }
382 
383  unsigned int indices[ndofs];
384 
385  get_dof_indices(cell, indices);
386  VecGetValues(values, ndofs, (PetscInt *)indices, local_values);
387 }
388 
390 {
391  for (ElementFullIter elem=mesh_->element.begin(); elem!=mesh_->element.end(); ++elem)
392  if (object_dofs[elem.index()] != NULL)
393  {
394  for (unsigned int j=0; j<elem->dim(); j++)
395  if (object_dofs[elem.index()][j] != NULL)
396  delete[] object_dofs[elem.index()][j];
397 
398  delete[] object_dofs[elem.index()];
399  }
400  delete[] object_dofs;
401 }
402 
403 
404 
406 {
407  // create local arrays of elements
408  int *id_4_old = new int[mesh_->n_elements()];
409  int i = 0;
410  FOR_ELEMENTS(mesh_, ele) id_4_old[i++] = ele.index();
412  delete[] id_4_old;
413 
414  // create local array of edges
415  for (unsigned int iedg=0; iedg<mesh_->edges.size(); iedg++)
416  {
417  bool is_edge_local = false;
418  Edge *edg = &mesh_->edges[iedg];
419  for (int sid=0; sid<edg->n_sides; sid++)
420  if (el_is_local(edg->side(sid)->element().index()))
421  {
422  is_edge_local = true;
423  break;
424  }
425  if (is_edge_local)
426  edg_4_loc.push_back(iedg);
427  }
428 
429  // create local array of neighbours
430  for (unsigned int inb=0; inb<mesh_->vb_neighbours_.size(); inb++)
431  {
432  Neighbour *nb = &mesh_->vb_neighbours_[inb];
433  if (el_is_local(mesh_->element.index(nb->element()))
434  || el_is_local(nb->side()->element().index()))
435  nb_4_loc.push_back(inb);
436  }
437 }
438 
439 
440 bool DOFHandlerMultiDim::el_is_local(int index) const
441 {
442  return el_ds_->is_local(row_4_el[index]);
443 }
444 
445 
446 
447 
448 
449 template<> FiniteElement<1,3> *DOFHandlerMultiDim::fe<1>() const { return fe1d_; }
450 template<> FiniteElement<2,3> *DOFHandlerMultiDim::fe<2>() const { return fe2d_; }
451 template<> FiniteElement<3,3> *DOFHandlerMultiDim::fe<3>() const { return fe3d_; }
452 
453 
454 //template class DOFHandler<0,3>;
455 //template class DOFHandler<1,3>;
456 //template class DOFHandler<2,3>;
457 //template class DOFHandler<3,3>;
Declaration of class which handles the ordering of degrees of freedom (dof) and mappings between loca...
vector< int > edg_4_loc
Local edge index -&gt; global edge index.
Definition: dofhandler.hh:386
DOFHandlerMultiDim(Mesh &_mesh)
Constructor.
Definition: dofhandler.cc:244
void id_maps(int n_ids, int *id_4_old, Distribution *&new_ds, int *&id_4_loc, int *&new_4_id)
FiniteElement< 3, 3 > * fe3d_
Definition: dofhandler.hh:366
#define FOR_ELEMENTS(_mesh_, __i)
Definition: mesh.h:364
FiniteElement< 1, 3 > * fe1d_
Pointer to the finite element class for which the handler distributes dofs.
Definition: dofhandler.hh:364
Distribution * el_ds_
Distribution of elements.
Definition: dofhandler.hh:383
???
Definition: mesh.h:109
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:257
int index() const
Definition: sys_vector.hh:88
int n_sides
Definition: edges.h:48
unsigned int global_dof_offset
Index of first global dof.
Definition: dofhandler.hh:113
Definition: edges.h:38
Partitioning * get_part()
Definition: mesh.cc:191
const std::vector< DofMultiplicity > dof_multiplicities
unsigned int lsize_
Number of dofs associated to local process.
Definition: dofhandler.hh:134
unsigned int n_elements() const
Definition: mesh.h:141
const unsigned int offset() const
Returns the number of the first global dof handled by this DOFHandler.
Definition: dofhandler.hh:78
#define ASSERT(...)
Definition: global_defs.h:121
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:144
SideIter side()
unsigned int index(const T *pointer) const
Definition: sys_vector.hh:383
FullIter begin()
Definition: sys_vector.hh:393
unsigned int np() const
get num of processors
void make_elem_partitioning()
Prepare parallel distribution of elements, edges and neighbours.
Definition: dofhandler.cc:405
Distribution * ds_
Distribution of dofs associated to local process.
Definition: dofhandler.hh:149
ElementIter element()
int *** object_dofs
Number of dofs associated to geometrical entities.
Definition: dofhandler.hh:375
unsigned int myp() const
get my processor
Support classes for parallel programing.
FiniteElement< 2, 3 > * fe2d_
Definition: dofhandler.hh:365
vector< Neighbour > vb_neighbours_
Definition: mesh.h:235
unsigned int loffset_
Index of the first dof on the local process.
Definition: dofhandler.hh:139
const unsigned int n_object_dofs(unsigned int object_dim, DofMultiplicity multiplicity)
Returns the number of single dofs/dof pairs/triples/sextuples that lie on a single geometric entity o...
void get_dof_indices(const CellIterator &cell, unsigned int indices[]) const override
Returns the global indices of dofs associated to the cell.
Definition: dofhandler.cc:346
ElementFullIter element() const
Definition: side_impl.hh:41
#define MPI_INT
Definition: mpi.h:160
~DOFHandlerMultiDim() override
Destructor.
Definition: dofhandler.cc:389
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:216
unsigned int n_dofs
Number of global dofs assigned by the handler.
Definition: dofhandler.hh:129
const Distribution * get_init_distr() const
Definition: partitioning.cc:58
const int * get_loc_part() const
Definition: partitioning.cc:65
#define MPI_Bcast(buffer, count, datatype, root, comm)
Definition: mpi.h:534
int * row_4_el
Global element index -&gt; index according to partitioning.
Definition: dofhandler.hh:379
const unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
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:366
bool el_is_local(int index) const
Definition: dofhandler.cc:440
vector< int > nb_4_loc
Local neighbour index -&gt; global neighbour index.
Definition: dofhandler.hh:389
FullIter end()
Returns FullFullIterer of the fictions past the end element.
Definition: sys_vector.hh:397
int * el_4_loc
Local element index -&gt; global element index.
Definition: dofhandler.hh:381
SideIter side(const unsigned int i) const
Definition: edges.h:43
ElementVector element
Vector of elements of the mesh.
Definition: mesh.h:205