Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
dof_distribution.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3 * Author: Andreas Vogel
4 *
5 * This file is part of UG4.
6 *
7 * UG4 is free software: you can redistribute it and/or modify it under the
8 * terms of the GNU Lesser General Public License version 3 (as published by the
9 * Free Software Foundation) with the following additional attribution
10 * requirements (according to LGPL/GPL v3 §7):
11 *
12 * (1) The following notice must be displayed in the Appropriate Legal Notices
13 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14 *
15 * (2) The following notice must be displayed at a prominent place in the
16 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17 *
18 * (3) The following bibliography is recommended for citation and must be
19 * preserved in all covered files:
20 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21 * parallel geometric multigrid solver on hierarchically distributed grids.
22 * Computing and visualization in science 16, 4 (2013), 151-164"
23 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24 * flexible software system for simulating pde based models on high performance
25 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26 *
27 * This program is distributed in the hope that it will be useful,
28 * but WITHOUT ANY WARRANTY; without even the implied warranty of
29 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 * GNU Lesser General Public License for more details.
31 */
32
33#ifndef __H__UG__LIB_DISC__DOF_MANAGER__DOF_DISTRIBUTION__
34#define __H__UG__LIB_DISC__DOF_MANAGER__DOF_DISTRIBUTION__
35
39#include "dof_index_storage.h"
40#include "dof_count.h"
41
42#ifdef UG_PARALLEL
44#endif
45
46namespace ug{
47
48class IGridFunction;
49
51{
52 public:
57 SmartPtr<SurfaceView> spSurfView,
58 const GridLevel& level, bool bGrouped,
59 SmartPtr<DoFIndexStorage> spDoFIndexStorage = SPNULL);
60
63
66
72
74 const GridLevel& grid_level() const {return m_gridLevel;}
75
76 public:
77 template <typename TElem>
84
85 template <int dim>
92
95 template <typename TElem>
98
99 template <typename TElem>
102 {return m_spSurfView->begin<TElem>(m_gridLevel, validStates);}
103
104 template <typename TElem>
107
108 template <typename TElem>
111 {return m_spSurfView->end<TElem>(m_gridLevel, validStates);}
112
113 template <typename TElem>
115 {return m_spSurfView->begin<TElem>(m_gridLevel, defaultValidSurfState());}
116
117 template <typename TElem>
120 {return m_spSurfView->begin<TElem>(m_gridLevel, validStates);}
121
122 template <typename TElem>
125
126 template <typename TElem>
129 {return m_spSurfView->end<TElem>(m_gridLevel, validStates);}
131
134 template <typename TElem>
136 {return m_spSurfView->begin<TElem>(si, m_gridLevel, defaultValidSurfState());}
137
138 template <typename TElem>
141 {return m_spSurfView->begin<TElem>(si, m_gridLevel, validStates);}
142
143 template <typename TElem>
145 {return m_spSurfView->end<TElem>(si, m_gridLevel, defaultValidSurfState());}
146
147 template <typename TElem>
150 {return m_spSurfView->end<TElem>(si, m_gridLevel, validStates);}
151
152 template <typename TElem>
154 {return m_spSurfView->begin<TElem>(si, m_gridLevel, defaultValidSurfState());}
155
156 template <typename TElem>
158 begin(int si, SurfaceView::SurfaceConstants validStates) const
159 {return m_spSurfView->begin<TElem>(si, m_gridLevel, validStates);}
160
161 template <typename TElem>
162 typename traits<TElem>::const_iterator end(int si) const
163 {return m_spSurfView->end<TElem>(si, m_gridLevel, defaultValidSurfState());}
164
165 template <typename TElem>
167 end(int si, SurfaceView::SurfaceConstants validStates) const
168 {return m_spSurfView->end<TElem>(si, m_gridLevel, validStates);}
170
173
175 template <typename TElem, typename TBaseElem>
176 void collect_associated(std::vector<TBaseElem*>& vAssElem,
177 TElem* elem, bool clearContainer = true) const{
178// if(grid_level().is_level())
179// CollectAssociated(vAssElem, *m_pMG, elem, clearContainer);
180// else{
181 m_spSurfView->collect_associated(vAssElem, elem, grid_level(), clearContainer);
182// }
183 }
184
185 template <typename TBaseElem>
187 (
188 std::vector<TBaseElem*>& vAssElem,
189 GridObject* elem,
190 bool clearContainer = true
191 ) const
192 {
193 if (dynamic_cast<Vertex*>(elem))
194 collect_associated<Vertex, TBaseElem>(vAssElem, dynamic_cast<Vertex*>(elem), clearContainer);
195 else if (dynamic_cast<Edge*>(elem))
196 collect_associated<Edge, TBaseElem>(vAssElem, dynamic_cast<Edge*>(elem), clearContainer);
197 else if (dynamic_cast<Face*>(elem))
198 collect_associated<Face, TBaseElem>(vAssElem, dynamic_cast<Face*>(elem), clearContainer);
199 else if (dynamic_cast<Volume*>(elem))
200 collect_associated<Volume, TBaseElem>(vAssElem, dynamic_cast<Volume*>(elem), clearContainer);
201 else
202 UG_THROW("Element is neither Vertex nor Edge, Face or Volume. "
203 "Other elements not implemented.");
204 }
205
207 template <class TGeomObj>
208 bool is_contained(TGeomObj* obj) const{
209 return m_spSurfView->is_contained(obj, grid_level());
210 }
211
212
213 public:
215 bool grouped() const {return m_bGrouped;}
216
218 size_t num_indices() const {return m_numIndex;}
219
221 size_t num_indices(int si) const {return m_vNumIndexOnSubset[si];}
222
223 public:
225
239 void indices(GridObject* elem, LocalIndices& ind, bool bHang = false) const;
240 void indices(Vertex* elem, LocalIndices& ind, bool bHang = false) const;
241 void indices(Edge* elem, LocalIndices& ind, bool bHang = false) const;
242 void indices(Face* elem, LocalIndices& ind, bool bHang = false) const;
243 void indices(Volume* elem, LocalIndices& ind, bool bHang = false) const;
245
247
264 size_t dof_indices(GridObject* elem, size_t fct, std::vector<DoFIndex>& ind,
265 bool bHang = false, bool bClear = true) const;
266 size_t dof_indices(Vertex* elem, size_t fct, std::vector<DoFIndex>& ind,
267 bool bHang = false, bool bClear = true) const;
268 size_t dof_indices(Edge* elem, size_t fct, std::vector<DoFIndex>& ind,
269 bool bHang = false, bool bClear = true) const;
270 size_t dof_indices(Face* elem, size_t fct, std::vector<DoFIndex>& ind,
271 bool bHang = false, bool bClear = true) const;
272 size_t dof_indices(Volume* elem, size_t fct, std::vector<DoFIndex>& ind,
273 bool bHang = false, bool bClear = true) const;
275
277
290 size_t inner_dof_indices(GridObject* elem, size_t fct, std::vector<DoFIndex>& ind,
291 bool bClear = true) const;
292 size_t inner_dof_indices(Vertex* elem, size_t fct, std::vector<DoFIndex>& ind,
293 bool bClear = true) const;
294 size_t inner_dof_indices(Edge* elem, size_t fct, std::vector<DoFIndex>& ind,
295 bool bClear = true) const;
296 size_t inner_dof_indices(Face* elem, size_t fct, std::vector<DoFIndex>& ind,
297 bool bClear = true) const;
298 size_t inner_dof_indices(Volume* elem, size_t fct, std::vector<DoFIndex>& ind,
299 bool bClear = true) const;
301
303
314 size_t algebra_indices(GridObject* elem, std::vector<size_t>& ind,
315 bool bClear = true) const;
316 size_t algebra_indices(Vertex* elem, std::vector<size_t>& ind,
317 bool bClear = true) const;
318 size_t algebra_indices(Edge* elem, std::vector<size_t>& ind,
319 bool bClear = true) const;
320 size_t algebra_indices(Face* elem, std::vector<size_t>& ind,
321 bool bClear = true) const;
322 size_t algebra_indices(Volume* elem, std::vector<size_t>& ind,
323 bool bClear = true) const;
325
326 size_t inner_algebra_indices_for_fct(GridObject* elem, std::vector<size_t>& ind,
327 bool bClear, int fct) const;
328
330
341 size_t inner_algebra_indices(GridObject* elem, std::vector<size_t>& ind,
342 bool bClear = true) const;
343 size_t inner_algebra_indices(Vertex* elem, std::vector<size_t>& ind,
344 bool bClear = true) const;
345 size_t inner_algebra_indices(Edge* elem, std::vector<size_t>& ind,
346 bool bClear = true) const;
347 size_t inner_algebra_indices(Face* elem, std::vector<size_t>& ind,
348 bool bClear = true) const;
349 size_t inner_algebra_indices(Volume* elem, std::vector<size_t>& ind,
350 bool bClear = true) const;
352
353 protected:
354 template <typename TBaseElem>
355 void _indices(TBaseElem* elem, LocalIndices& ind, bool bHang = false) const;
356
357 template<typename TBaseElem>
358 size_t _dof_indices(TBaseElem* elem, size_t fct,
359 std::vector<DoFIndex>& ind,
360 bool bHang = false, bool bClear = true) const;
361
362 template<typename TBaseElem>
363 size_t _inner_dof_indices(TBaseElem* elem, size_t fct,
364 std::vector<DoFIndex>& ind,
365 bool bClear = true) const;
366
367 template<typename TBaseElem>
368 size_t _algebra_indices(TBaseElem* elem, std::vector<size_t>& ind,
369 bool bClear = true) const;
370
371 template<typename TBaseElem>
372 size_t _inner_algebra_indices(TBaseElem* elem, std::vector<size_t>& ind,
373 bool bClear = true) const;
374
376 template <typename TBaseElem>
377 void changable_indices(std::vector<size_t>& vIndex,
378 const std::vector<TBaseElem*>& vElem) const;
379
381 template<typename TBaseElem>
382 void indices_on_vertex(TBaseElem* elem, const ReferenceObjectID roid,
383 LocalIndices& ind,
384 const Grid::SecureVertexContainer& vElem) const;
385
387 template <typename TConstraining, typename TConstrained, typename TBaseElem>
389 const typename Grid::traits<TBaseElem>::secure_container& vSubElem) const;
390
391 template <typename TBaseElem,typename TConstraining, typename TConstrained, typename TSubElem>
392 void constrained_edge_indices(TBaseElem* elem,LocalIndices& ind,
393 const typename Grid::traits<TSubElem>::secure_container& vSubElem) const;
394
395 template <typename TBaseElem,typename TConstraining, typename TConstrained, typename TSubElem>
396 void constrained_face_indices(TBaseElem* elem,LocalIndices& ind,
397 const typename Grid::traits<TSubElem>::secure_container& vSubElem) const;
398
399 // sorts indices on constrained edges
400 template <typename TBaseElem,typename TConstraining, typename TConstrained>
401 void sort_constrained_edges(std::vector<size_t>& sortedInd,TBaseElem* elem,TConstraining* constrainingObj,size_t objIndex) const;
402
403 // sorts indices on constrained faces
404 template <typename TBaseElem,typename TConstraining, typename TConstrained>
405 void sort_constrained_faces(std::vector<size_t>& sortedInd,TBaseElem* elem,TConstraining* constrainingObj,size_t objIndex) const;
406
408 template<typename TBaseElem, typename TSubBaseElem>
409 void indices(TBaseElem* elem, const ReferenceObjectID roid,
410 LocalIndices& ind,
411 const typename Grid::traits<TSubBaseElem>::secure_container& vElem) const;
412
414 template<typename TBaseElem, typename TSubBaseElem>
415 void dof_indices(TBaseElem* elem, const ReferenceObjectID roid,
416 size_t fct, std::vector<DoFIndex>& ind,
417 const typename Grid::traits<TSubBaseElem>::secure_container& vElem) const;
418
420 template <typename TConstraining, typename TConstrained, typename TBaseElem>
421 void constrained_vertex_dof_indices(size_t fct,std::vector<DoFIndex>& ind,
422 const typename Grid::traits<TBaseElem>::secure_container& vSubElem) const;
423
425 template <typename TBaseElem,typename TConstraining, typename TConstrained, typename TSubElem>
426 void constrained_edge_dof_indices(TBaseElem* elem,size_t fct,std::vector<DoFIndex>& ind,
427 const typename Grid::traits<TSubElem>::secure_container& vSubElem) const;
428
430 template <typename TBaseElem,typename TConstraining, typename TConstrained, typename TSubElem>
431 void constrained_face_dof_indices(TBaseElem* elem,size_t fct,std::vector<DoFIndex>& ind,
432 const typename Grid::traits<TSubElem>::secure_container& vSubElem) const;
433
435 template <typename TBaseElem>
436 size_t extract_inner_algebra_indices(TBaseElem* elem,
437 std::vector<size_t>& ind) const;
438
440 template<typename TBaseElem>
442 std::vector<size_t>& ind) const;
443
444 protected:
446 template <typename TBaseObject>
447 void add(TBaseObject* obj, const ReferenceObjectID roid, const int si);
448
450 void check_subsets();
451
454 template <typename TElem>
455 inline size_t& obj_index(TElem* obj) {return m_spDoFIndexStorage->obj_index(obj);}
456
457 template <typename TElem>
458 inline const size_t& obj_index(TElem* obj) const {return m_spDoFIndexStorage->obj_index(obj);}
460
461 protected:
464
468
471
474
477
480
481 protected:
484
486 std::vector<size_t> m_vNumIndexOnSubset;
487
488 public:
490 void get_connections(std::vector<std::vector<size_t> >& vvConnection) const;
491
493 void permute_indices(const std::vector<size_t>& vIndNew);
494
496 void reinit();
497
498 protected:
500 template <typename TBaseElem>
501 void reinit();
502
503 template <typename TBaseElem>
504 void permute_indices(const std::vector<size_t>& vNewInd);
505
506 public:
507 template <typename TBaseElem>
508 void get_connections(std::vector<std::vector<size_t> >& vvConnection) const;
509
510 public:
512 void manage_grid_function(IGridFunction& gridFct);
513
516
517 public:
519 void permute_values(const std::vector<size_t>& vIndNew);
520
522 void copy_values(const std::vector<std::pair<size_t, size_t> >& vIndexMap,
523 bool bDisjunct);
524
526 void resize_values(size_t newSize);
527
528 protected:
530
538 std::vector<IGridFunction*> m_vpGridFunction;
539
540#ifdef UG_PARALLEL
541 public:
547
548 protected:
551
553
554 void reinit_index_layout(IndexLayout& layout, int keyType);
555
556 template <typename TBaseElem>
557 void add_indices_from_layouts(IndexLayout& indexLayout, int keyType);
558#endif
559
560 public:
561 DoFCount dof_count() const;
562
563 protected:
564 template <typename TBaseElem>
565 void sum_dof_count(DoFCount& cnt) const;
566
567};
568
569} // end namespace ug
570
571#endif /* __H__UG__LIB_DISC__DOF_MANAGER__DOF_DISTRIBUTION__ */
Definition smart_pointer.h:296
Definition smart_pointer.h:108
Definition dof_count.h:43
Definition dof_distribution.h:51
const size_t & obj_index(TElem *obj) const
Definition dof_distribution.h:458
void constrained_face_indices(TBaseElem *elem, LocalIndices &ind, const typename Grid::traits< TSubElem >::secure_container &vSubElem) const
Definition dof_distribution.cpp:1063
std::vector< size_t > m_vNumIndexOnSubset
number of distributed indices on each subset
Definition dof_distribution.h:486
size_t _algebra_indices(TBaseElem *elem, std::vector< size_t > &ind, bool bClear=true) const
Definition dof_distribution.cpp:300
traits< TElem >::iterator begin()
Definition dof_distribution.h:96
SmartPtr< AlgebraLayouts > m_spAlgebraLayouts
algebra layouts
Definition dof_distribution.h:550
traits< TElem >::const_iterator begin() const
Definition dof_distribution.h:114
ConstSmartPtr< AlgebraLayouts > layouts() const
Definition dof_distribution.h:544
MultiGrid * m_pMG
Definition dof_distribution.h:467
void constrained_vertex_indices(LocalIndices &ind, const typename Grid::traits< TBaseElem >::secure_container &vSubElem) const
extract dofs on constrained objects
Definition dof_distribution.cpp:853
void _indices(TBaseElem *elem, LocalIndices &ind, bool bHang=false) const
Definition dof_distribution.cpp:1118
size_t _inner_dof_indices(TBaseElem *elem, size_t fct, std::vector< DoFIndex > &ind, bool bClear=true) const
Definition dof_distribution.cpp:410
void resize_values(size_t newSize)
changes values in managed functions, number of indices changed
Definition dof_distribution.cpp:1326
void constrained_edge_indices(TBaseElem *elem, LocalIndices &ind, const typename Grid::traits< TSubElem >::secure_container &vSubElem) const
Definition dof_distribution.cpp:1008
traits< TElem >::iterator end(SurfaceView::SurfaceConstants validStates)
Definition dof_distribution.h:110
void reinit()
initializes the indices
Definition dof_distribution.cpp:1430
traits< TElem >::iterator end()
Definition dof_distribution.h:105
traits< TElem >::iterator begin(int si, SurfaceView::SurfaceConstants validStates)
Definition dof_distribution.h:140
const GridLevel & grid_level() const
returns grid level
Definition dof_distribution.h:74
size_t dof_indices(GridObject *elem, size_t fct, std::vector< DoFIndex > &ind, bool bHang=false, bool bClear=true) const
extracts all multiindices for a function (sorted)
Definition dof_distribution.cpp:1227
traits< TElem >::iterator begin(int si)
Definition dof_distribution.h:135
void indices(GridObject *elem, LocalIndices &ind, bool bHang=false) const
extracts all indices of the element (sorted)
Definition dof_distribution.cpp:1211
SmartPtr< MultiGrid > m_spMG
Multi Grid.
Definition dof_distribution.h:466
void check_subsets()
checks that subset assignment is ok
Definition dof_distribution.cpp:94
ConstSmartPtr< SurfaceView > surface_view() const
returns the surface view
Definition dof_distribution.h:65
traits< TElem >::iterator end(int si, SurfaceView::SurfaceConstants validStates)
Definition dof_distribution.h:149
size_t algebra_indices(GridObject *elem, std::vector< size_t > &ind, bool bClear=true) const
extracts all algebra indices of an element (not sorted)
Definition dof_distribution.cpp:1260
SmartPtr< MultiGrid > multi_grid()
Definition dof_distribution.h:69
void add_indices_from_layouts(IndexLayout &indexLayout, int keyType)
Definition dof_distribution.cpp:1501
size_t _inner_algebra_indices(TBaseElem *elem, std::vector< size_t > &ind, bool bClear=true) const
Definition dof_distribution.cpp:243
void sort_constrained_faces(std::vector< size_t > &sortedInd, TBaseElem *elem, TConstraining *constrainingObj, size_t objIndex) const
Definition dof_distribution.cpp:957
traits< TElem >::const_iterator begin(SurfaceView::SurfaceConstants validStates) const
Definition dof_distribution.h:119
void collect_associated(std::vector< TBaseElem * > &vAssElem, TElem *elem, bool clearContainer=true) const
returns the adjacent elements
Definition dof_distribution.h:176
void copy_values(const std::vector< std::pair< size_t, size_t > > &vIndexMap, bool bDisjunct)
swaps values in managed functions, if indices swapped
Definition dof_distribution.cpp:1318
void unmanage_grid_function(IGridFunction &gridFct)
unregisters a grid function for adaptation management
Definition dof_distribution.cpp:1304
size_t _dof_indices(TBaseElem *elem, size_t fct, std::vector< DoFIndex > &ind, bool bHang=false, bool bClear=true) const
Definition dof_distribution.cpp:665
void get_connections(std::vector< std::vector< size_t > > &vvConnection) const
returns the connections
Definition dof_distribution.cpp:1697
traits< TElem >::const_iterator begin(int si, SurfaceView::SurfaceConstants validStates) const
Definition dof_distribution.h:158
size_t m_numIndex
number of distributed indices on whole domain
Definition dof_distribution.h:483
SurfaceView::SurfaceConstants defaultValidSurfState() const
returns the default valid surface state
Definition dof_distribution.cpp:135
void reinit_layouts_and_communicator()
Definition dof_distribution.cpp:1448
traits< TElem >::const_iterator end(int si) const
Definition dof_distribution.h:162
void constrained_face_dof_indices(TBaseElem *elem, size_t fct, std::vector< DoFIndex > &ind, const typename Grid::traits< TSubElem >::secure_container &vSubElem) const
multi indices on constrained faces
Definition dof_distribution.cpp:590
size_t inner_algebra_indices_for_fct(GridObject *elem, std::vector< size_t > &ind, bool bClear, int fct) const
Definition dof_distribution.cpp:256
size_t inner_algebra_indices(GridObject *elem, std::vector< size_t > &ind, bool bClear=true) const
extracts all algebra indices in the inner of the element (not sorted)
Definition dof_distribution.cpp:1276
void indices_on_vertex(TBaseElem *elem, const ReferenceObjectID roid, LocalIndices &ind, const Grid::SecureVertexContainer &vElem) const
extracts the indices of the vertices
Definition dof_distribution.cpp:729
traits< TElem >::iterator begin(SurfaceView::SurfaceConstants validStates)
Definition dof_distribution.h:101
void reinit_index_layout(IndexLayout &layout, int keyType)
Definition dof_distribution.cpp:1483
traits< TElem >::const_iterator end(SurfaceView::SurfaceConstants validStates) const
Definition dof_distribution.h:128
size_t inner_dof_indices(GridObject *elem, size_t fct, std::vector< DoFIndex > &ind, bool bClear=true) const
extracts all multiindices of a function in the inner (sorted)
Definition dof_distribution.cpp:1243
size_t num_indices() const
return the number of dofs distributed
Definition dof_distribution.h:218
void constrained_vertex_dof_indices(size_t fct, std::vector< DoFIndex > &ind, const typename Grid::traits< TBaseElem >::secure_container &vSubElem) const
multi indices on constrained vertices
Definition dof_distribution.cpp:457
size_t & obj_index(TElem *obj)
Definition dof_distribution.h:455
traits< TElem >::const_iterator end(int si, SurfaceView::SurfaceConstants validStates) const
Definition dof_distribution.h:167
bool m_bGrouped
grouping
Definition dof_distribution.h:463
SmartPtr< AlgebraLayouts > layouts()
Definition dof_distribution.h:545
size_t extract_inner_algebra_indices(TBaseElem *elem, std::vector< size_t > &ind) const
adds all algebra indices of an geom object to the LocalIndices
Definition dof_distribution.cpp:190
traits< TElem >::iterator end(int si)
Definition dof_distribution.h:144
void permute_indices(const std::vector< size_t > &vIndNew)
renames the indices
Definition dof_distribution.cpp:1754
void manage_grid_function(IGridFunction &gridFct)
registers a grid function for adaptation management
Definition dof_distribution.cpp:1293
void add(TBaseObject *obj, const ReferenceObjectID roid, const int si)
adds indices to a geometric object
Definition dof_distribution.cpp:147
void collect_associated(std::vector< TBaseElem * > &vAssElem, GridObject *elem, bool clearContainer=true) const
Definition dof_distribution.h:187
bool is_contained(TGeomObj *obj) const
returns if the grid object is part of the dof distribution
Definition dof_distribution.h:208
bool grouped() const
returns if dofs are grouped
Definition dof_distribution.h:215
SmartPtr< SurfaceView > m_spSurfView
MultiGrid Subset Handler.
Definition dof_distribution.h:473
void changable_indices(std::vector< size_t > &vIndex, const std::vector< TBaseElem * > &vElem) const
returns indices, that can be changed on the element
Definition dof_distribution.cpp:1183
ConstSmartPtr< MultiGrid > multi_grid() const
Definition dof_distribution.h:70
GridLevel m_gridLevel
Grid level.
Definition dof_distribution.h:476
traits< TElem >::const_iterator begin(int si) const
Definition dof_distribution.h:153
void sort_constrained_edges(std::vector< size_t > &sortedInd, TBaseElem *elem, TConstraining *constrainingObj, size_t objIndex) const
Definition dof_distribution.cpp:906
SmartPtr< MGSubsetHandler > m_spMGSH
Subset Handler.
Definition dof_distribution.h:470
traits< TElem >::const_iterator end() const
Definition dof_distribution.h:123
DoFCount dof_count() const
Definition dof_distribution.cpp:1609
void constrained_edge_dof_indices(TBaseElem *elem, size_t fct, std::vector< DoFIndex > &ind, const typename Grid::traits< TSubElem >::secure_container &vSubElem) const
multi indices on constrained edges
Definition dof_distribution.cpp:515
~DoFDistribution()
destructor
Definition dof_distribution.cpp:91
void permute_values(const std::vector< size_t > &vIndNew)
permutes values in managed functions, if indices permuted
Definition dof_distribution.cpp:1311
size_t num_indices(int si) const
return the number of dofs distributed on subset si
Definition dof_distribution.h:221
std::vector< IGridFunction * > m_vpGridFunction
managed grid functions
Definition dof_distribution.h:538
void sum_dof_count(DoFCount &cnt) const
Definition dof_distribution.cpp:1571
SmartPtr< DoFIndexStorage > m_spDoFIndexStorage
DoF-Index Memory Storage.
Definition dof_distribution.h:479
Definition dof_distribution_info.h:152
Base-class for edges.
Definition grid_base_objects.h:397
Faces are 2-dimensional objects.
Definition grid_base_objects.h:510
traits< Vertex >::secure_container SecureVertexContainer
Container to store associated vertices.
Definition grid.h:162
Definition grid_level.h:42
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
Base class for all Grid Functions.
Definition grid_function.h:56
Definition local_algebra.h:50
Definition multi_grid.h:72
SurfaceConstants
Definition surface_view.h:72
Base-class for all vertex-types.
Definition grid_base_objects.h:231
Volumes are 3-dimensional objects.
Definition grid_base_objects.h:754
const NullSmartPtr SPNULL
The equivalent to NULL for smart pointers.
Definition smart_pointer.h:90
#define UG_THROW(msg)
Definition error.h:57
the ug namespace
ReferenceObjectID
these ids are used to identify the shape of a geometric object.
Definition grid_base_objects.h:74
Definition dof_distribution.h:87
domain_traits< dim >::grid_base_object grid_base_object
Definition dof_distribution.h:88
SurfaceView::traits< grid_base_object >::iterator iterator
Definition dof_distribution.h:89
SurfaceView::traits< grid_base_object >::const_iterator const_iterator
Definition dof_distribution.h:90
Definition dof_distribution.h:79
TElem grid_object
Definition dof_distribution.h:80
SurfaceView::traits< TElem >::const_iterator const_iterator
Definition dof_distribution.h:82
SurfaceView::traits< TElem >::iterator iterator
Definition dof_distribution.h:81
PointerConstArray< TElem * > secure_container
Definition grid.h:146
ConstSurfaceViewElementIterator< TElem > const_iterator
Definition surface_view.h:300
SurfaceViewElementIterator< TElem > iterator
Definition surface_view.h:299
Definition domain_traits.h:53