G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsAdaptiveRefUtils.h
Go to the documentation of this file.
1
14#pragma once
15
16
17#include <iostream>
18
19namespace gismo
20{
21
22enum MarkingStrategy
23{
24 GARU=1,
25 PUCA=2,
26 BULK=3,
27 PBULK=4
28};
29
30template <class T>
31void gsMarkThreshold( const std::vector<T> & elError, T refParameter,
32 std::vector<bool> & elMarked)
33{
34 // First, conduct a brutal search for the maximum local error
35 const T maxErr = *std::max_element(elError.begin(), elError.end() );
36
37 // Compute the threshold:
38 const T Thr = refParameter * maxErr;
39
40 elMarked.resize( elError.size() );
41 // Now just check for each element, whether the local error
42 // is above the computed threshold or not, and mark accordingly.
43
44 typename std::vector<T>::const_iterator err = elError.begin();
45 for(std::vector<bool>::iterator i = elMarked.begin(); i!= elMarked.end(); ++i, ++err)
46 *i = ( *err >= Thr );
47}
48
49template <class T>
50void gsMarkPercentage( const std::vector<T> & elError, T refParameter,
51 std::vector<bool> & elMarked)
52{
53 T Thr = (T)(0);
54
55 // Total number of elements:
56 size_t NE = elError.size();
57 // The vector of local errors will need to be sorted,
58 // which will be done on a copy:
59 std::vector<T> elErrCopy = elError;
60
61 // Compute the index from which the refinement should start,
62 // once the vector is sorted.
63 size_t idxRefineStart = cast<T,size_t>( math::floor( refParameter * (T)(NE) ) );
64 // ...and just to be sure we are in range:
65 if( idxRefineStart == elErrCopy.size() )
66 {
67 GISMO_ASSERT(idxRefineStart >= 1, "idxRefineStart can't get negative");
68 idxRefineStart -= 1;
69 }
70
71 // Sort the list using bubblesort.
72 // After each loop, the largest elements are at the end
73 // of the list. Since we are only interested in the largest elements,
74 // it is enough to run the sorting until enough "largest" elements
75 // have been found, i.e., until we have reached indexRefineStart
76 size_t lastSwapDone = elErrCopy.size() - 1;
77 size_t lastCheckIdx = lastSwapDone;
78
79 bool didSwap;
80 T tmp;
81 do{
82 didSwap = false;
83 lastCheckIdx = lastSwapDone;
84 for( size_t i=0; i < lastCheckIdx; i++)
85 if( elErrCopy[i] > elErrCopy[i+1] )
86 {
87 tmp = elErrCopy[i];
88 elErrCopy[i] = elErrCopy[i+1];
89 elErrCopy[i+1] = tmp;
90
91 didSwap = true;
92 lastSwapDone = i;
93 }
94 }while( didSwap && (lastSwapDone+1 >= idxRefineStart ) );
95
96 // Compute the threshold:
97 Thr = elErrCopy[ idxRefineStart ];
98 elMarked.resize( elError.size() );
99 // Now just check for each element, whether the local error
100 // is above the computed threshold or not, and mark accordingly.
101 for( size_t i=0; i < elError.size(); i++)
102 ( elError[i] >= Thr ? elMarked[i] = true : elMarked[i] = false );
103}
104
105
106template <class T>
107void gsMarkFraction( const std::vector<T> & elError, T refParameter,
108 std::vector<bool> & elMarked)
109{
110 T Thr = (T)(0);
111
112 // The vector of local errors will need to be sorted,
113 // which will be done on a copy:
114 std::vector<T> elErrCopy = elError;
115
116 // Compute the sum, i.e., the global/total error
117 T totalError = (T)(0);
118 for( size_t i = 0; i < elErrCopy.size(); ++i)
119 totalError += elErrCopy[i];
120
121 // We want to mark just enough cells such that their
122 // cummulated errors add up to a certain fraction
123 // of the total error.
124 T errorMarkSum = (1-refParameter) * totalError;
125 T cummulErrMarked = 0;
126
127 T tmp;
128 GISMO_ASSERT(elErrCopy.size() >= 1, "elErrCopy needs at least 1 element");
129 size_t lastSwapDone = elErrCopy.size() - 1;
130 do{
131 for( size_t i=0; i < lastSwapDone; i++)
132 if( elErrCopy[i] > elErrCopy[i+1] )
133 {
134 tmp = elErrCopy[i];
135 elErrCopy[i] = elErrCopy[i+1];
136 elErrCopy[i+1] = tmp;
137 }
138
139 cummulErrMarked += elErrCopy[ lastSwapDone ];
140 lastSwapDone -= 1;
141
142 }while( cummulErrMarked < errorMarkSum && lastSwapDone > 0 );
143
144 // Compute the threshold:
145 Thr = elErrCopy[ lastSwapDone + 1 ];
146 elMarked.resize( elError.size() );
147 // Now just check for each element, whether the local error
148 // is above the computed threshold or not, and mark accordingly.
149 for( size_t i=0; i < elError.size(); i++)
150 ( elError[i] >= Thr ? elMarked[i] = true : elMarked[i] = false );
151}
152
153
195template <class T>
196void gsMarkElementsForRef( const std::vector<T> & elError, int refCriterion,
197 T refParameter, std::vector<bool> & elMarked)
198{
199 switch (refCriterion)
200 {
201 case GARU:
202 gsMarkThreshold(elError,refParameter,elMarked);
203 break;
204 case PUCA:
205 gsMarkPercentage(elError,refParameter,elMarked);
206 break;
207 case BULK:
208 gsMarkFraction(elError,refParameter,elMarked);
209 break;
210 default:
211 GISMO_ERROR("unknown marking strategy");
212 }
213}
214
215template <class T>
216std::vector<gsMatrix<T>> gsGetRefinementBoxes( gsMultiBasis<T> & basis,
217 const std::vector<bool> & elMarked)
218{
219 GISMO_ASSERT(basis.totalElements()==elMarked.size(),
220 "Incorrect input, num_el="<<basis.totalElements()
221 <<", mark_vec="<<elMarked.size() );
222
223 std::vector<gsMatrix<T>> result(basis.nBases());
224 const short_t dim = basis.dim();
225
226 // numMarked: Number of marked cells on current patch, also currently marked cell
227 // poffset : offset index for the first element on a patch
228 // globalCount: counter for the current global element index
229 index_t numMarked, poffset = 0, globalCount = 0;
230
231 // refBoxes: contains marked boxes on a given patch
232 gsMatrix<T> refBoxes;
233
234 for (size_t pn=0; pn < basis.nBases(); ++pn )// for all patches
235 {
236 // Get number of elements to be refined on this patch
237 const size_t numEl = basis[pn].numElements();
238 numMarked = std::count_if(elMarked.begin() + poffset,
239 elMarked.begin() + poffset + numEl,
240 GS_BIND2ND(std::equal_to<bool>(), true) );
241 poffset += numEl;
242 refBoxes.resize(dim, 2*numMarked);
243 numMarked = 0;// counting current patch element to be refined
244
245 // for all elements in patch pn
246 typename gsBasis<T>::domainIter domIt = basis.basis(pn).makeDomainIterator();
247 for (; domIt->good(); domIt->next())
248 {
249 if( elMarked[ globalCount++ ] ) // refine this element ?
250 {
251 // Construct degenerate box by setting both
252 // corners equal to the center
253 refBoxes.col(2*numMarked ) = domIt->lowerCorner();
254 refBoxes.col(2*numMarked+1) = domIt->upperCorner();
255
256 // Advance marked cells counter
257 numMarked++;
258 }
259 }
260
261 result[pn] = refBoxes;
262 }
263 return result;
264}
265
293template <class T>
295 const std::vector<bool> & elMarked,
296 index_t refExtension = 0)
297{
298 GISMO_ASSERT(basis.totalElements()==elMarked.size(),
299 "Incorrect input, num_el="<<basis.totalElements()
300 <<", mark_vec="<<elMarked.size() );
301
302 std::vector<gsMatrix<T>> boxes = gsGetRefinementBoxes(basis,elMarked);
303 for (size_t pn=0; pn < basis.nBases(); ++pn )// for all patches
304 {
305 // Refine all of the found refBoxes in this patch
306 basis.refine( pn, boxes[pn], refExtension );
307 }
308}
309
310template <class T>
311void gsRefineWithBoxes( gsMultiBasis<T> & basis,
312 const std::vector<gsMatrix<T>> & boxes,
313 index_t refExtension = 0)
314{
315 for (size_t pn=0; pn < basis.nBases(); ++pn )// for all patches
316 {
317 // Refine all of the found refBoxes in this patch
318 basis.refine( pn, boxes[pn], refExtension );
319 }
320}
321
322template <class T>
323void gsRefineMarkedElements(gsMultiPatch<T> & mp,
324 const std::vector<bool> & elMarked,
325 int refExtension = 0)
326{
327 const int dim = mp.dim();
328
329 // numMarked: Number of marked cells on current patch, also currently marked cell
330 // poffset : offset index for the first element on a patch
331 // globalCount: counter for the current global element index
332 int numMarked, poffset = 0, globalCount = 0;
333
334 // refBoxes: contains marked boxes on a given patch
335 gsMatrix<T> refBoxes;
336
337 for (size_t pn=0; pn < mp.nPatches(); ++pn )// for all patches
338 {
339 // Get number of elements to be refined on this patch
340 const int numEl = mp[pn].basis().numElements();
341 numMarked = std::count_if(elMarked.begin() + poffset,
342 elMarked.begin() + poffset + numEl,
343 GS_BIND2ND(std::equal_to<bool>(), true) );
344
345 poffset += numEl;
346 refBoxes.resize(dim, 2*numMarked);
347 numMarked = 0;// counting current patch element to be refined
348
349 // for all elements in patch pn
350 typename gsBasis<T>::domainIter domIt = mp.patch(pn).basis().makeDomainIterator();
351 for (; domIt->good(); domIt->next())
352 {
353 if( elMarked[ globalCount++ ] ) // refine this element ?
354 {
355 // Construct degenerate box by setting both
356 // corners equal to the center
357 refBoxes.col(2*numMarked ) =
358 refBoxes.col(2*numMarked+1) = domIt->centerPoint();
359
360 // Advance marked cells counter
361 numMarked++;
362 }
363 }
364 // Refine all of the found refBoxes in this patch
365 std::vector<index_t> elements = mp.patch(pn).basis().asElements(refBoxes, refExtension);
366 mp.patch(pn).refineElements( elements );
367 }
368}
369
370template <class T>
371void gsRefineMarkedFunctions(gsMultiPatch<T> & mp,
372 const std::vector<bool> & funMarked,
373 int refExtension = 0)
374{
375 GISMO_ASSERT(funMarked.size()==gsMultiBasis<T>(mp).totalSize(),"Vector of marked functions must have the same size as the total number of basis functions (over all patches)!");
376
377 const int dim = mp.dim();
378
379 int numMarked, poffset = 0, globalCount = 0;
380
381 // refBoxes: contains marked boxes on a given patch
382 gsMatrix<T> refBoxes;
383
384 for (size_t pn=0; pn < mp.nPatches(); ++pn )// for all patches
385 {
386 // Get number of elements to be refined on this patch
387 const int numFun = mp[pn].basis().size();
388 numMarked = std::count_if(funMarked.begin() + poffset,
389 funMarked.begin() + poffset + numFun,
390 GS_BIND2ND(std::equal_to<bool>(), true) );
391
392 poffset += numFun;
393 refBoxes.resize(dim, 2*numMarked);
394 numMarked = 0;// counting current patch element to be refined
395 for (index_t i = 0; i != mp.basis(pn).size(); ++i)
396 if (funMarked[globalCount++])
397 {
398 refBoxes.block(0,2*numMarked,2,2) = mp.basis(pn).support(i);
399 numMarked++;
400 }
401
402 // Refine all of the found refBoxes in this patch
403 std::vector<index_t> elements = mp.patch(pn).basis().asElements(refBoxes, refExtension);
404
405 mp.patch(pn).refineElements( elements );
406 }
407}
408
436template <class T>
438 const std::vector<bool> & elMarked,
439 index_t refExtension = 0)
440{
441 GISMO_ASSERT(basis.numElements()==elMarked.size(),
442 "Incorrect input, num_el="<<basis.numElements()
443 <<", mark_vec="<<elMarked.size() );
444
445 const short_t dim = basis.dim();
446
447 // numMarked: Number of marked cells on current patch, also currently marked cell
448 // poffset : offset index for the first element on a patch
449 // globalCount: counter for the current global element index
450 index_t numMarked, poffset = 0, globalCount = 0;
451
452 // refBoxes: contains marked boxes on a given patch
453 gsMatrix<T> refBoxes;
454
455 for (size_t pn=0; pn < basis.nBases(); ++pn )// for all patches
456 {
457 // Get number of elements to be refined on this patch
458 const size_t numEl = basis[pn].numElements();
459 numMarked = std::count_if(elMarked.begin() + poffset,
460 elMarked.begin() + poffset + numEl,
461 GS_BIND2ND(std::equal_to<bool>(), true) );
462 poffset += numEl;
463 refBoxes.resize(dim, 2*numMarked);
464 numMarked = 0;// counting current patch element to be refined
465
466 // for all elements in patch pn
467 typename gsBasis<T>::domainIter domIt = basis.basis(pn).makeDomainIterator();
468 for (; domIt->good(); domIt->next())
469 {
470 if( elMarked[ globalCount++ ] ) // refine this element ?
471 {
472 // Construct degenerate box by setting both
473 // corners equal to the center
474 refBoxes.col(2*numMarked ) = domIt->lowerCorner();
475 refBoxes.col(2*numMarked+1) = domIt->upperCorner();
476
477 // Advance marked cells counter
478 numMarked++;
479 }
480 }
481
482 // Refine all of the found refBoxes in this patch
483 basis.unrefine( pn, refBoxes, refExtension );
484 }
485}
486
487template <class T>
488void gsProcessMarkedElements(gsMultiBasis<T> & basis,
489 const std::vector<bool> & elRefined,
490 const std::vector<bool> & elCoarsened,
491 index_t refExtension = 0,
492 index_t crsExtension = 0)
493{
494 const short_t dim = basis.dim();
495
496 // numMarked: Number of marked cells on current patch, also currently marked cell
497 // poffset : offset index for the first element on a patch
498 // globalCount: counter for the current global element index
499 index_t numRefined, numCoarsened, poffset = 0, globalCount = 0;
500
501 // refBoxes: contains marked boxes on a given patch
502 gsMatrix<T> refBoxes, crsBoxes;
503
504 for (size_t pn=0; pn < basis.nBases(); ++pn )// for all patches
505 {
506 // Get number of elements to be refined on this patch
507 const size_t numEl = basis[pn].numElements();
508 numRefined = std::count_if(elRefined.begin() + poffset,
509 elRefined.begin() + poffset + numEl,
510 GS_BIND2ND(std::equal_to<bool>(), true) );
511
512 numCoarsened = std::count_if(elCoarsened.begin() + poffset,
513 elCoarsened.begin() + poffset + numEl,
514 GS_BIND2ND(std::equal_to<bool>(), true) );
515
516
517 poffset += numEl;
518
519 refBoxes.resize(dim, 2*numRefined);
520 crsBoxes.resize(dim, 2*numCoarsened);
521 numRefined = 0;// counting current patch element to be refined
522 numCoarsened = 0;// counting current patch element to be refined
523
524 // for all elements in patch pn
525 typename gsBasis<T>::domainIter domIt = basis.basis(pn).makeDomainIterator();
526 for (; domIt->good(); domIt->next())
527 {
528 if( elRefined[ globalCount ] ) // refine this element ?
529 {
530 // Construct degenerate box by setting both
531 // corners equal to the center
532 refBoxes.col(2*numRefined ) = domIt->lowerCorner();
533 refBoxes.col(2*numRefined+1) = domIt->upperCorner();
534
535 // Advance marked cells counter
536 numRefined++;
537 }
538 if( elCoarsened[ globalCount ] ) // refine this element ?
539 {
540 // Construct degenerate box by setting both
541 // corners equal to the center
542 crsBoxes.col(2*numCoarsened ) = domIt->lowerCorner();
543 crsBoxes.col(2*numCoarsened+1) = domIt->upperCorner();
544
545 // Advance marked cells counter
546 numCoarsened++;
547 }
548
549 globalCount++;
550 }
551 // Refine all of the found refBoxes in this patch
552 basis.unrefine( pn, crsBoxes, crsExtension);
553 basis.refine( pn, refBoxes, refExtension );
554 }
555}
556
557template <class T>
558void gsProcessMarkedElements(gsMultiPatch<T> & mp,
559 const std::vector<bool> & elRefined,
560 const std::vector<bool> & elCoarsened,
561 index_t refExtension = 0,
562 index_t crsExtension = 0)
563{
564 const int dim = mp.dim();
565
566 // numMarked: Number of marked cells on current patch, also currently marked cell
567 // poffset : offset index for the first element on a patch
568 // globalCount: counter for the current global element index
569 index_t numRefined, numCoarsened, poffset = 0, globalCount = 0;
570
571 // refBoxes: contains marked boxes on a given patch
572 gsMatrix<T> refBoxes, crsBoxes;
573
574 for (size_t pn=0; pn < mp.nPatches(); ++pn )// for all patches
575 {
576 // Get number of elements to be refined on this patch
577 const size_t numEl = mp[pn].basis().numElements();
578 numRefined = std::count_if(elRefined.begin() + poffset,
579 elRefined.begin() + poffset + numEl,
580 GS_BIND2ND(std::equal_to<bool>(), true) );
581
582 numCoarsened = std::count_if(elCoarsened.begin() + poffset,
583 elCoarsened.begin() + poffset + numEl,
584 GS_BIND2ND(std::equal_to<bool>(), true) );
585
586
587 poffset += numEl;
588 refBoxes.resize(dim, 2*numRefined);
589 crsBoxes.resize(dim, 2*numCoarsened);
590 numRefined = 0;// counting current patch element to be refined
591 numCoarsened = 0;// counting current patch element to be refined
592
593 // for all elements in patch pn
594 typename gsBasis<T>::domainIter domIt = mp.patch(pn).basis().makeDomainIterator();
595 for (; domIt->good(); domIt->next())
596 {
597 if( elRefined[ globalCount ] ) // refine this element ?
598 {
599 // Construct degenerate box by setting both
600 // corners equal to the center
601 refBoxes.col(2*numRefined ) = domIt->lowerCorner();
602 refBoxes.col(2*numRefined+1) = domIt->upperCorner();
603
604 // Advance marked cells counter
605 numRefined++;
606 }
607 if( elCoarsened[ globalCount ] ) // refine this element ?
608 {
609 // Construct degenerate box by setting both
610 // corners equal to the center
611 crsBoxes.col(2*numCoarsened ) = domIt->lowerCorner();
612 crsBoxes.col(2*numCoarsened+1) = domIt->upperCorner();
613
614 // Advance marked cells counter
615 numCoarsened++;
616 }
617
618 globalCount++;
619 }
620
621 std::vector<index_t> elements;
622 // Unrefine all of the found refBoxes in this patch
623 elements = mp.patch(pn).basis().asElementsUnrefine(crsBoxes, crsExtension);
624 mp.patch(pn).unrefineElements( elements );
625
626 // Refine all of the found refBoxes in this patch
627 elements = mp.patch(pn).basis().asElements(refBoxes, refExtension);
628 mp.patch(pn).refineElements( elements );
629 }
630}
631
632
633template <class T>
634void gsUnrefineMarkedElements(gsMultiPatch<T> & mp,
635 const std::vector<bool> & elMarked,
636 int refExtension = 0)
637{
638 const int dim = mp.dim();
639
640 // numMarked: Number of marked cells on current patch, also currently marked cell
641 // poffset : offset index for the first element on a patch
642 // globalCount: counter for the current global element index
643 int numMarked, poffset = 0, globalCount = 0;
644
645 // refBoxes: contains marked boxes on a given patch
646 gsMatrix<T> refBoxes;
647
648 for (size_t pn=0; pn < mp.nPatches(); ++pn )// for all patches
649 {
650 // Get number of elements to be refined on this patch
651 const int numEl = mp[pn].basis().numElements();
652 numMarked = std::count_if(elMarked.begin() + poffset,
653 elMarked.begin() + poffset + numEl,
654 GS_BIND2ND(std::equal_to<bool>(), true) );
655
656 poffset += numEl;
657 refBoxes.resize(dim, 2*numMarked);
658 numMarked = 0;// counting current patch element to be refined
659
660 // for all elements in patch pn
661 typename gsBasis<T>::domainIter domIt = mp.patch(pn).basis().makeDomainIterator();
662 for (; domIt->good(); domIt->next())
663 {
664 if( elMarked[ globalCount++ ] ) // refine this element ?
665 {
666 // Construct degenerate box by setting both
667 // corners equal to the center
668 refBoxes.col(2*numMarked ) =
669 refBoxes.col(2*numMarked+1) = domIt->centerPoint();
670
671 // Advance marked cells counter
672 numMarked++;
673 }
674 }
675 // Refine all of the found refBoxes in this patch
676 std::vector<index_t> elements = mp.patch(pn).basis().asElementsUnrefine(refBoxes, refExtension);
677 mp.patch(pn).unrefineElements( elements );
678 }
679}
680
681template <class T>
682void gsUnrefineMarkedFunctions(gsMultiPatch<T> & mp,
683 const std::vector<bool> & funMarked,
684 int refExtension = 0)
685{
686 GISMO_ASSERT(funMarked.size()==gsMultiBasis<T>(mp).totalSize(),"Vector of marked functions must have the same size as the total number of basis functions (over all patches)!");
687
688 const int dim = mp.dim();
689
690 int numMarked, poffset = 0, globalCount = 0;
691
692 // refBoxes: contains marked boxes on a given patch
693 gsMatrix<T> refBoxes;
694
695 for (size_t pn=0; pn < mp.nPatches(); ++pn )// for all patches
696 {
697 // Get number of elements to be refined on this patch
698 const int numFun = mp[pn].basis().size();
699 numMarked = std::count_if(funMarked.begin() + poffset,
700 funMarked.begin() + poffset + numFun,
701 GS_BIND2ND(std::equal_to<bool>(), true) );
702
703 poffset += numFun;
704 refBoxes.resize(dim, 2*numMarked);
705 numMarked = 0;// counting current patch element to be refined
706 for (index_t i = 0; i != mp.basis(pn).size(); ++i)
707 if (funMarked[globalCount++])
708 {
709 refBoxes.block(0,2*numMarked,2,2) = mp.basis(pn).support(i);
710 numMarked++;
711 }
712
713 // Refine all of the found refBoxes in this patch
714 std::vector<index_t> elements = mp.patch(pn).basis().asElements(refBoxes, refExtension);
715
716 mp.patch(pn).unrefineElements( elements );
717 }
718}
719
720
721} // namespace gismo
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Holds a set of patch-wise bases and their topology information.
Definition gsMultiBasis.h:37
void gsMarkElementsForRef(const std::vector< T > &elError, int refCriterion, T refParameter, std::vector< bool > &elMarked)
Marks elements/cells for refinement.
Definition gsAdaptiveRefUtils.h:196
void gsRefineMarkedElements(gsMultiBasis< T > &basis, const std::vector< bool > &elMarked, index_t refExtension=0)
Refine a gsMultiBasis, based on a vector of element-markings.
Definition gsAdaptiveRefUtils.h:294
void gsUnrefineMarkedElements(gsMultiBasis< T > &basis, const std::vector< bool > &elMarked, index_t refExtension=0)
Unrefine a gsMultiBasis, based on a vector of element-markings.
Definition gsAdaptiveRefUtils.h:437
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_ERROR(message)
Definition gsDebug.h:118
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
The G+Smo namespace, containing all definitions for the library.