Line data Source code
1 : /*
2 : * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Sebastian Reiter
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 : #include <cassert>
34 : #include "subset_handler_multi_grid.h"
35 :
36 : using namespace std;
37 :
38 : namespace ug
39 : {
40 : ////////////////////////////////////////////////////////////////////////
41 : // MultiGridSubsetHandler implementation
42 0 : MultiGridSubsetHandler::
43 0 : MultiGridSubsetHandler(uint supportedElements) :
44 : ISubsetHandler(supportedElements),
45 : m_aSharedEntryVRT("MGSubsetHandler_SharedListEntryVRT", false),
46 : m_aSharedEntryEDGE("MGSubsetHandler_SharedListEntryEDGE", false),
47 : m_aSharedEntryFACE("MGSubsetHandler_SharedListEntryFACE", false),
48 0 : m_aSharedEntryVOL("MGSubsetHandler_SharedListEntryVOL", false)
49 : {
50 0 : m_numSubsets = 0;
51 0 : m_pMG = NULL;
52 0 : }
53 :
54 0 : MultiGridSubsetHandler::
55 0 : MultiGridSubsetHandler(MultiGrid& mg, uint supportedElements) :
56 : ISubsetHandler(supportedElements),
57 : m_aSharedEntryVRT("MGSubsetHandler_SharedListEntryVRT", false),
58 : m_aSharedEntryEDGE("MGSubsetHandler_SharedListEntryEDGE", false),
59 : m_aSharedEntryFACE("MGSubsetHandler_SharedListEntryFACE", false),
60 0 : m_aSharedEntryVOL("MGSubsetHandler_SharedListEntryVOL", false)
61 : {
62 0 : m_numSubsets = 0;
63 0 : m_pMG = NULL;
64 0 : assign_grid(mg);
65 0 : }
66 :
67 0 : MultiGridSubsetHandler::
68 0 : MultiGridSubsetHandler(const MultiGridSubsetHandler& sh) :
69 0 : ISubsetHandler(sh.m_supportedElements),
70 : m_aSharedEntryVRT("MGSubsetHandler_SharedListEntryVRT", false),
71 : m_aSharedEntryEDGE("MGSubsetHandler_SharedListEntryEDGE", false),
72 : m_aSharedEntryFACE("MGSubsetHandler_SharedListEntryFACE", false),
73 0 : m_aSharedEntryVOL("MGSubsetHandler_SharedListEntryVOL", false)
74 : {
75 0 : MultiGrid* pGrid = sh.m_pMG;
76 :
77 0 : if(pGrid){
78 : //TODO: remove virtual function calls from constructor
79 0 : assign_grid(*pGrid);
80 0 : assign_subset_handler(sh);
81 : }
82 0 : }
83 :
84 0 : MultiGridSubsetHandler::~MultiGridSubsetHandler()
85 : {
86 0 : if(m_pGrid != NULL){
87 0 : erase_subset_lists_impl();
88 0 : detach_data();
89 : }
90 0 : }
91 :
92 0 : void MultiGridSubsetHandler::grid_to_be_destroyed(Grid* grid)
93 : {
94 : assert((m_pGrid == grid) && "ERROR in MultiGridSubsetHandler::grid_to_be_destroyed(...): Grids do not match.");
95 0 : cleanup();
96 0 : ISubsetHandler::grid_to_be_destroyed(grid);
97 0 : }
98 :
99 0 : void MultiGridSubsetHandler::cleanup()
100 : {
101 0 : erase_subset_lists_impl();
102 : m_levels.clear();
103 :
104 0 : if(m_pGrid){
105 0 : detach_data();
106 :
107 : // unregister the previously registered callback
108 0 : m_callbackId = MessageHub::SPCallbackId(NULL);
109 :
110 0 : m_pMG = NULL;
111 : }
112 :
113 : //ISubsetHandler::set_grid(NULL);
114 0 : }
115 :
116 0 : void MultiGridSubsetHandler::assign_grid(MultiGrid& mg)
117 : {
118 0 : if(m_pMG == &mg)
119 : return;
120 :
121 0 : if(m_pMG)
122 0 : cleanup();
123 :
124 0 : m_pMG = &mg;
125 :
126 0 : ISubsetHandler::set_grid(&mg);
127 :
128 : // attach shared entries
129 0 : if(elements_are_supported(SHE_VERTEX))
130 0 : m_pGrid->attach_to_vertices(m_aSharedEntryVRT);
131 0 : if(elements_are_supported(SHE_EDGE))
132 0 : m_pGrid->attach_to_edges(m_aSharedEntryEDGE);
133 0 : if(elements_are_supported(SHE_FACE))
134 0 : m_pGrid->attach_to_faces(m_aSharedEntryFACE);
135 0 : if(elements_are_supported(SHE_VOLUME))
136 0 : m_pGrid->attach_to_volumes(m_aSharedEntryVOL);
137 :
138 : // register the callback
139 0 : m_callbackId = m_pMG->message_hub()->register_class_callback(
140 0 : this, &ug::MultiGridSubsetHandler::multigrid_changed);
141 :
142 0 : level_required(m_pMG->num_levels());
143 : }
144 :
145 0 : void MultiGridSubsetHandler::detach_data()
146 : {
147 0 : if(elements_are_supported(SHE_VERTEX))
148 0 : m_pGrid->detach_from_vertices(m_aSharedEntryVRT);
149 0 : if(elements_are_supported(SHE_EDGE))
150 0 : m_pGrid->detach_from_edges(m_aSharedEntryEDGE);
151 0 : if(elements_are_supported(SHE_FACE))
152 0 : m_pGrid->detach_from_faces(m_aSharedEntryFACE);
153 0 : if(elements_are_supported(SHE_VOLUME))
154 0 : m_pGrid->detach_from_volumes(m_aSharedEntryVOL);
155 0 : }
156 :
157 0 : void MultiGridSubsetHandler::erase_subset_lists()
158 : {
159 0 : erase_subset_lists_impl();
160 0 : }
161 :
162 0 : void MultiGridSubsetHandler::erase_subset_lists_impl()
163 : {
164 0 : for(size_t level = 0; level < m_levels.size(); ++level)
165 : {
166 0 : for(size_t i = 0; i < m_levels[level].size(); ++i)
167 0 : delete m_levels[level][i];
168 :
169 : m_levels[level].clear();
170 : }
171 :
172 0 : m_numSubsets = 0;
173 0 : }
174 :
175 0 : void MultiGridSubsetHandler::clear_subset_lists(int index)
176 : {
177 0 : if(m_pGrid)
178 : {
179 0 : for(size_t level = 0; level < m_levels.size(); ++level){
180 0 : section_container<Vertex>(index, level).clear();
181 0 : section_container<Edge>(index, level).clear();
182 0 : section_container<Face>(index, level).clear();
183 0 : section_container<Volume>(index, level).clear();
184 : }
185 : }
186 0 : }
187 :
188 : template<class TElem>
189 : void
190 0 : MultiGridSubsetHandler::
191 : assign_subset_impl(TElem* elem, int subsetIndex)
192 : {
193 : assert((m_pGrid != NULL) && "ERROR in SubsetHandler::assign_subset(): No grid assigned to SubsetHandler.");
194 :
195 0 : int level = m_pMG->get_level(elem);
196 :
197 : // check if we have to remove elem from a subset.
198 0 : int oldIndex = get_subset_index(elem);
199 :
200 0 : if(oldIndex != -1)
201 0 : section_container<TElem>(oldIndex, level).erase(get_list_iterator(elem), elem->container_section());
202 :
203 : // add the element to the subset.
204 0 : if(subsetIndex != -1)
205 : {
206 : subset_required(subsetIndex);
207 : level_required(level);
208 0 : section_container<TElem>(subsetIndex, level).insert(elem, elem->container_section());
209 : subset_assigned(elem, subsetIndex);
210 : }
211 : else {
212 : //TODO: iterator is useless!
213 : subset_assigned(elem, -1);
214 : }
215 0 : }
216 :
217 0 : void MultiGridSubsetHandler::assign_subset(Vertex* elem, int subsetIndex)
218 : {
219 0 : if(elements_are_supported(SHE_VERTEX))
220 0 : assign_subset_impl(elem, subsetIndex);
221 0 : }
222 :
223 0 : void MultiGridSubsetHandler::assign_subset(Edge* elem, int subsetIndex)
224 : {
225 0 : if(elements_are_supported(SHE_EDGE))
226 0 : assign_subset_impl(elem, subsetIndex);
227 0 : }
228 :
229 0 : void MultiGridSubsetHandler::assign_subset(Face* elem, int subsetIndex)
230 : {
231 0 : if(elements_are_supported(SHE_FACE))
232 0 : assign_subset_impl(elem, subsetIndex);
233 0 : }
234 :
235 0 : void MultiGridSubsetHandler::assign_subset(Volume* elem, int subsetIndex)
236 : {
237 0 : if(elements_are_supported(SHE_VOLUME))
238 0 : assign_subset_impl(elem, subsetIndex);
239 0 : }
240 :
241 0 : void MultiGridSubsetHandler::
242 : change_subset_indices(int indOld, int indNew)
243 : {
244 0 : if(m_pGrid)
245 : {
246 0 : if(elements_are_supported(SHE_VERTEX))
247 0 : change_elem_subset_indices<Vertex>(indOld, indNew);
248 0 : if(elements_are_supported(SHE_EDGE))
249 0 : change_elem_subset_indices<Edge>(indOld, indNew);
250 0 : if(elements_are_supported(SHE_FACE))
251 0 : change_elem_subset_indices<Face>(indOld, indNew);
252 0 : if(elements_are_supported(SHE_VOLUME))
253 0 : change_elem_subset_indices<Volume>(indOld, indNew);
254 : }
255 0 : }
256 :
257 0 : void MultiGridSubsetHandler::add_required_subset_lists(int maxIndex)
258 : {
259 0 : while((int)num_subsets_in_list() <= maxIndex)
260 0 : add_subset_to_all_levels();
261 0 : }
262 :
263 0 : void MultiGridSubsetHandler::erase_subset_lists(int index)
264 : {
265 0 : for(size_t level = 0; level < m_levels.size(); ++level)
266 : {
267 : SubsetVec& subsets = m_levels[level];
268 0 : delete subsets[index];
269 0 : for(uint i = index + 1; i < num_subsets_in_list(); ++i)
270 0 : subsets[i-1] = subsets[i];
271 0 : subsets.resize(subsets.size() - 1);
272 : }
273 0 : m_numSubsets--;
274 0 : }
275 :
276 0 : void MultiGridSubsetHandler::swap_subset_lists(int ind1, int ind2)
277 : {
278 0 : for(size_t level = 0; level < m_levels.size(); ++level)
279 : {
280 : SubsetVec& subsets = m_levels[level];
281 0 : Subset* pTmp = subsets[ind1];
282 0 : subsets[ind1] = subsets[ind2];
283 0 : subsets[ind2] = pTmp;
284 : }
285 0 : }
286 :
287 0 : void MultiGridSubsetHandler::move_subset_lists(int indexFrom, int indexTo)
288 : {
289 : int moveDir = 0;
290 : // we have to distinguish two cases
291 0 : if(indexFrom < indexTo)
292 : moveDir = 1;
293 0 : else if(indexTo < indexFrom)
294 : moveDir = -1;
295 :
296 : if(moveDir != 0)
297 : {
298 0 : for(size_t level = 0; level < m_levels.size(); ++level)
299 : {
300 : SubsetVec& subsets = m_levels[level];
301 :
302 : // store pointer to the from-subset
303 0 : Subset* pFrom = subsets[indexFrom];
304 :
305 0 : for(int i = indexFrom; i != indexTo; i+= moveDir)
306 : {
307 0 : int iNext = i+moveDir;
308 :
309 : // move pointer
310 0 : subsets[i] = subsets[iNext];
311 : }
312 :
313 : // assign stored pointer
314 0 : subsets[indexTo] = pFrom;
315 : }
316 : }
317 0 : }
318 :
319 0 : void MultiGridSubsetHandler::join_subset_lists(int target, int src1, int src2)
320 : {
321 0 : for(size_t level = 0; level < m_levels.size(); ++level){
322 : SubsetVec& subsets = m_levels[level];
323 :
324 : // store pointer to the from-subset
325 0 : Subset& t = *subsets[target];
326 0 : Subset& s1 = *subsets[src1];
327 0 : Subset& s2 = *subsets[src2];
328 :
329 0 : if(target != src1){
330 0 : t.m_vertices.append(s1.m_vertices);
331 0 : t.m_edges.append(s1.m_edges);
332 0 : t.m_faces.append(s1.m_faces);
333 0 : t.m_volumes.append(s1.m_volumes);
334 : }
335 0 : if(target != src2){
336 0 : t.m_vertices.append(s2.m_vertices);
337 0 : t.m_edges.append(s2.m_edges);
338 0 : t.m_faces.append(s2.m_faces);
339 0 : t.m_volumes.append(s2.m_volumes);
340 : }
341 : }
342 :
343 0 : if(target != src1)
344 0 : clear_subset_lists(src1);
345 0 : if(target != src2)
346 0 : clear_subset_lists(src2);
347 0 : }
348 :
349 : /*
350 : void MultiGridSubsetHandler::
351 : register_subset_elements_at_pipe()
352 : {
353 : for(size_t l = 0; l < num_levels(); ++l)
354 : {
355 : for(size_t i = 0; i < num_subsets_in_list(); ++i)
356 : {
357 : // register vertices
358 : for(VertexIterator iter = begin<Vertex>(i, l);
359 : iter != end<Vertex>(i, l); ++iter)
360 : register_at_pipe(*iter);
361 :
362 : // register edges
363 : for(EdgeIterator iter = begin<Edge>(i, l);
364 : iter != end<Edge>(i, l); ++iter)
365 : register_at_pipe(*iter);
366 :
367 : // register faces
368 : for(FaceIterator iter = begin<Face>(i, l);
369 : iter != end<Face>(i, l); ++iter)
370 : register_at_pipe(*iter);
371 :
372 : // register volumes
373 : for(VolumeIterator iter = begin<Volume>(i, l);
374 : iter != end<Volume>(i, l); ++iter)
375 : register_at_pipe(*iter);
376 : }
377 : }
378 : }
379 : */
380 :
381 : GridObjectCollection
382 0 : MultiGridSubsetHandler::
383 : get_grid_objects(int subsetIndex, int level) const
384 : {
385 0 : subset_required(subsetIndex);
386 0 : level_required(level);
387 :
388 : return GridObjectCollection(&m_levels[level][subsetIndex]->m_vertices,
389 : &m_levels[level][subsetIndex]->m_edges,
390 : &m_levels[level][subsetIndex]->m_faces,
391 0 : &m_levels[level][subsetIndex]->m_volumes);
392 : }
393 :
394 : GridObjectCollection
395 0 : MultiGridSubsetHandler::
396 : get_grid_objects_in_subset(int subsetIndex) const
397 : {
398 0 : subset_required(subsetIndex);
399 0 : GridObjectCollection goc(m_levels.size());
400 0 : for(size_t i = 0; i < m_levels.size(); ++i)
401 : {
402 0 : goc.add_level( &m_levels[i][subsetIndex]->m_vertices,
403 : &m_levels[i][subsetIndex]->m_edges,
404 : &m_levels[i][subsetIndex]->m_faces,
405 0 : &m_levels[i][subsetIndex]->m_volumes);
406 : }
407 :
408 0 : return goc;
409 : }
410 :
411 : GridObjectCollection
412 0 : MultiGridSubsetHandler::
413 : get_grid_objects_in_level(int level) const
414 : {
415 0 : level_required(level);
416 : uint numSubsets = num_subsets_in_list();
417 0 : GridObjectCollection goc(numSubsets);
418 0 : for(uint i = 0; i < numSubsets; ++i)
419 : {
420 0 : goc.add_level( &m_levels[level][i]->m_vertices,
421 : &m_levels[level][i]->m_edges,
422 : &m_levels[level][i]->m_faces,
423 0 : &m_levels[level][i]->m_volumes);
424 : }
425 :
426 0 : return goc;
427 : }
428 :
429 0 : MultiGridSubsetHandler::Subset* MultiGridSubsetHandler::new_subset()
430 : {
431 :
432 0 : Subset* sub = new Subset;
433 0 : if(elements_are_supported(SHE_VERTEX))
434 0 : sub->m_vertices.get_container().set_pipe(
435 0 : &m_pGrid->get_attachment_pipe<Vertex>(), m_aSharedEntryVRT);
436 0 : if(elements_are_supported(SHE_EDGE))
437 0 : sub->m_edges.get_container().set_pipe(
438 0 : &m_pGrid->get_attachment_pipe<Edge>(), m_aSharedEntryEDGE);
439 0 : if(elements_are_supported(SHE_FACE))
440 0 : sub->m_faces.get_container().set_pipe(
441 0 : &m_pGrid->get_attachment_pipe<Face>(), m_aSharedEntryFACE);
442 0 : if(elements_are_supported(SHE_VOLUME))
443 0 : sub->m_volumes.get_container().set_pipe(
444 0 : &m_pGrid->get_attachment_pipe<Volume>(), m_aSharedEntryVOL);
445 0 : return sub;
446 : }
447 :
448 0 : void MultiGridSubsetHandler::add_level()
449 : {
450 0 : m_levels.push_back(SubsetVec());
451 0 : int topLevel = m_levels.size() - 1;
452 0 : for(uint i = 0; i < num_subsets_in_list(); ++i)
453 0 : m_levels[topLevel].push_back(new_subset());
454 0 : }
455 :
456 0 : void MultiGridSubsetHandler::add_subset_to_all_levels()
457 : {
458 0 : for(size_t i = 0; i < m_levels.size(); ++i)
459 0 : m_levels[i].push_back(new_subset());
460 :
461 0 : m_numSubsets++;
462 0 : }
463 :
464 0 : void MultiGridSubsetHandler::
465 : multigrid_changed(const GridMessage_MultiGridChanged& gm)
466 : {
467 0 : if(gm.message_type() == GMMGCT_LEVEL_ADDED)
468 : level_required(gm.num_levels_in_grid() - 1);
469 0 : }
470 :
471 : /*
472 : size_t MultiGridSubsetHandler::
473 : collect_subset_elements(std::vector<Vertex*>& vrtsOut, int subsetIndex) const
474 : {
475 : return collect_subset_elements_impl(vrtsOut, subsetIndex);
476 : }
477 :
478 : size_t MultiGridSubsetHandler::
479 : collect_subset_elements(std::vector<Edge*>& edgesOut, int subsetIndex) const
480 : {
481 : return collect_subset_elements_impl(edgesOut, subsetIndex);
482 : }
483 :
484 : size_t MultiGridSubsetHandler::
485 : collect_subset_elements(std::vector<Face*>& facesOut, int subsetIndex) const
486 : {
487 : return collect_subset_elements_impl(facesOut, subsetIndex);
488 : }
489 :
490 : size_t MultiGridSubsetHandler::
491 : collect_subset_elements(std::vector<Volume*>& volsOut, int subsetIndex) const
492 : {
493 : return collect_subset_elements_impl(volsOut, subsetIndex);
494 : }
495 :
496 : template <class TElem>
497 : size_t MultiGridSubsetHandler::
498 : collect_subset_elements_impl(std::vector<TElem*>& elemsOut, int subsetIndex) const
499 : {
500 : typedef typename geometry_traits<TElem>::iterator ElemIter;
501 : typedef typename geometry_traits<TElem>::const_iterator ConstElemIter;
502 :
503 : elemsOut.clear();
504 : if(!m_pGrid)
505 : return 0;
506 :
507 : if(subsetIndex < 0){
508 : // iterate over all elements of the underlying grid and compare indices
509 : Grid& grid = *m_pGrid;
510 :
511 : for(ElemIter iter = grid.begin<TElem>(); iter != grid.end<TElem>(); ++iter)
512 : {
513 : if(get_subset_index(*iter) == -1)
514 : elemsOut.push_back(*iter);
515 : }
516 : }
517 : else{
518 : elemsOut.reserve(num<TElem>(subsetIndex));
519 : for(size_t i = 0; i < num_levels(); ++i){
520 : for(ConstElemIter iter = begin<TElem>(subsetIndex, i);
521 : iter != end<TElem>(subsetIndex, i); ++iter)
522 : {
523 : elemsOut.push_back(*iter);
524 : }
525 : }
526 : }
527 :
528 : return elemsOut.size();
529 : }
530 : */
531 : }// end of namespace
|