ViennaLS
Loading...
Searching...
No Matches
lsAdvect.hpp
Go to the documentation of this file.
1#pragma once
2
4
5#include <limits>
6#include <vector>
7
8#include <hrleSparseIterator.hpp>
9#include <hrleSparseStarIterator.hpp>
10
11#include <vcLogger.hpp>
12#include <vcSmartPointer.hpp>
13
15#include <lsDomain.hpp>
16#include <lsMarkVoidPoints.hpp>
17#include <lsReduce.hpp>
18
19// Integration schemes
20#include <lsConcepts.hpp>
21#include <lsEnquistOsher.hpp>
22#include <lsLaxFriedrichs.hpp>
27
28// Velocity accessor
29#include <lsVelocityField.hpp>
30
31// #define DEBUG_LS_ADVECT_HPP
32#ifdef DEBUG_LS_ADVECT_HPP
33#include <lsToMesh.hpp>
34#include <lsVTKWriter.hpp>
35#endif
36
38template <class IntegrationSchemeType, class T, int D,
39 lsConcepts::IsNotSame<IntegrationSchemeType,
42void reduceTimeStepHamiltonJacobi(IntegrationSchemeType &, double &,
43 hrleCoordType) {}
44} // namespace lsInternal::advect
45
46namespace viennals {
47
48using namespace viennacore;
49
64
74template <class T, int D> class Advect {
75 std::vector<SmartPointer<Domain<T, D>>> levelSets;
76 SmartPointer<VelocityField<T>> velocities = nullptr;
77 IntegrationSchemeEnum integrationScheme =
79 double timeStepRatio = 0.4999;
80 double dissipationAlpha = 1.0;
81 bool calculateNormalVectors = true;
82 bool ignoreVoids = false;
83 double advectionTime = 0.;
84 bool performOnlySingleStep = false;
85 double advectedTime = 0.;
86 unsigned numberOfTimeSteps = 0;
87 bool saveAdvectionVelocities = false;
88 bool updatePointData = true;
89 static constexpr double wrappingLayerEpsilon = 1e-4;
90
91 void rebuildLS() {
92 // TODO: this function uses manhatten distances for renormalisation,
93 // since this is the quickest. For visualisation applications, better
94 // renormalisation is needed, so it might be good to implement
95 // Euler distance renormalisation as an option
96 auto &grid = levelSets.back()->getGrid();
97 auto newlsDomain = SmartPointer<Domain<T, D>>::New(grid);
98 typename Domain<T, D>::DomainType &newDomain = newlsDomain->getDomain();
99 typename Domain<T, D>::DomainType &domain = levelSets.back()->getDomain();
100
101 newDomain.initialize(domain.getNewSegmentation(),
102 domain.getAllocation() *
103 (2.0 / levelSets.back()->getLevelSetWidth()));
104
105 const bool updateData = updatePointData;
106 // save how data should be transferred to new level set
107 // list of indices into the old pointData vector
108 std::vector<std::vector<unsigned>> newDataSourceIds;
109 if (updateData)
110 newDataSourceIds.resize(newDomain.getNumberOfSegments());
111
112#ifdef DEBUG_LS_ADVECT_HPP
113 {
114 auto mesh = SmartPointer<Mesh<T>>::New();
115 ToMesh<T, D>(levelSets.back(), mesh).apply();
116 VTKWriter<T>(mesh, "Advect_beforeRebuild.vtk").apply();
117 }
118#endif
119
120#pragma omp parallel num_threads(newDomain.getNumberOfSegments())
121 {
122 int p = 0;
123#ifdef _OPENMP
124 p = omp_get_thread_num();
125#endif
126
127 auto &domainSegment = newDomain.getDomainSegment(p);
128
129 hrleVectorType<hrleIndexType, D> startVector =
130 (p == 0) ? grid.getMinGridPoint()
131 : newDomain.getSegmentation()[p - 1];
132
133 hrleVectorType<hrleIndexType, D> endVector =
134 (p != static_cast<int>(newDomain.getNumberOfSegments() - 1))
135 ? newDomain.getSegmentation()[p]
136 : grid.incrementIndices(grid.getMaxGridPoint());
137
138 // reserve a bit more to avoid reallocation
139 // would expect number of points to roughly double
140 if (updateData)
141 newDataSourceIds[p].reserve(2.5 * domainSegment.getNumberOfPoints());
142
143 for (hrleSparseStarIterator<typename Domain<T, D>::DomainType, 1> it(
144 domain, startVector);
145 it.getIndices() < endVector; ++it) {
146
147 // if the center is an active grid point
148 // <1.0 since it could have been change by 0.5 max
149 if (std::abs(it.getCenter().getValue()) <= 1.0) {
150
151 int k = 0;
152 for (; k < 2 * D; k++)
153 if (std::signbit(it.getNeighbor(k).getValue() - 1e-7) !=
154 std::signbit(it.getCenter().getValue() + 1e-7))
155 break;
156
157 // if there is at least one neighbor of opposite sign
158 if (k != 2 * D) {
159 if (it.getCenter().getDefinedValue() > 0.5) {
160 int j = 0;
161 for (; j < 2 * D; j++) {
162 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
163 if (it.getNeighbor(j).getDefinedValue() < -0.5)
164 break;
165 }
166 if (j == 2 * D) {
167 domainSegment.insertNextDefinedPoint(
168 it.getIndices(), it.getCenter().getDefinedValue());
169 if (updateData)
170 newDataSourceIds[p].push_back(it.getCenter().getPointId());
171 // if there is at least one active grid point, which is < -0.5
172 } else {
173 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
174 if (updateData)
175 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
176 }
177 } else if (it.getCenter().getDefinedValue() < -0.5) {
178 int j = 0;
179 for (; j < 2 * D; j++) {
180 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
181 if (it.getNeighbor(j).getDefinedValue() > 0.5)
182 break;
183 }
184
185 if (j == 2 * D) {
186 domainSegment.insertNextDefinedPoint(
187 it.getIndices(), it.getCenter().getDefinedValue());
188 if (updateData)
189 newDataSourceIds[p].push_back(it.getCenter().getPointId());
190 // if there is at least one active grid point, which is > 0.5
191 } else {
192 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
193 if (updateData)
194 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
195 }
196 } else {
197 domainSegment.insertNextDefinedPoint(
198 it.getIndices(), it.getCenter().getDefinedValue());
199 if (updateData)
200 newDataSourceIds[p].push_back(it.getCenter().getPointId());
201 }
202 } else {
203 domainSegment.insertNextUndefinedPoint(
204 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
207 }
208
209 } else { // if the center is not an active grid point
210 if (it.getCenter().getValue() >= 0) {
211 int usedNeighbor = -1;
212 T distance = Domain<T, D>::POS_VALUE;
213 for (int i = 0; i < 2 * D; i++) {
214 T value = it.getNeighbor(i).getValue();
215 if (std::abs(value) <= 1.0 && (value < 0.)) {
216 if (distance > value + 1.0) {
217 distance = value + 1.0;
218 usedNeighbor = i;
219 }
220 }
221 }
222
223 if (distance <= 1.) {
224 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
225 if (updateData)
226 newDataSourceIds[p].push_back(
227 it.getNeighbor(usedNeighbor).getPointId());
228 } else {
229 domainSegment.insertNextUndefinedPoint(it.getIndices(),
231 }
232
233 } else {
234 int usedNeighbor = -1;
235 T distance = Domain<T, D>::NEG_VALUE;
236 for (int i = 0; i < 2 * D; i++) {
237 T value = it.getNeighbor(i).getValue();
238 if (std::abs(value) <= 1.0 && (value > 0)) {
239 if (distance < value - 1.0) {
240 // distance = std::max(distance, value - T(1.0));
241 distance = value - 1.0;
242 usedNeighbor = i;
243 }
244 }
245 }
246
247 if (distance >= -1.) {
248 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
249 if (updateData)
250 newDataSourceIds[p].push_back(
251 it.getNeighbor(usedNeighbor).getPointId());
252 } else {
253 domainSegment.insertNextUndefinedPoint(it.getIndices(),
255 }
256 }
257 }
258 }
259 }
260
261 // now copy old data into new level set
262 if (updateData) {
263 auto &pointData = levelSets.back()->getPointData();
264 newlsDomain->getPointData().translateFromMultiData(pointData,
265 newDataSourceIds);
266 }
267
268 newDomain.finalize();
269 newDomain.segment();
270 levelSets.back()->deepCopy(newlsDomain);
271 levelSets.back()->finalize(2);
272 }
273
276 double advect(double maxTimeStep = std::numeric_limits<double>::max()) {
277 // check whether a level set and velocites have been given
278 if (levelSets.size() < 1) {
279 Logger::getInstance()
280 .addWarning("No level sets passed to Advect. Not advecting.")
281 .print();
282 return std::numeric_limits<double>::max();
283 }
284 if (velocities == nullptr) {
285 Logger::getInstance()
286 .addWarning("No velocity field passed to Advect. Not advecting.")
287 .print();
288 return std::numeric_limits<double>::max();
289 }
290
291 prepareLS();
292
293 double currentTime = 0.;
294 if (integrationScheme == IntegrationSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
295 auto is = lsInternal::EnquistOsher<T, D, 1>(levelSets.back(), velocities,
296 calculateNormalVectors);
297 currentTime = integrateTime(is, maxTimeStep);
298 } else if (integrationScheme ==
300 auto is = lsInternal::EnquistOsher<T, D, 2>(levelSets.back(), velocities,
301 calculateNormalVectors);
302 currentTime = integrateTime(is, maxTimeStep);
303 } else if (integrationScheme ==
305 auto is = lsInternal::LaxFriedrichs<T, D, 1>(levelSets.back(), velocities,
306 dissipationAlpha,
307 calculateNormalVectors);
308 currentTime = integrateTime(is, maxTimeStep);
309 } else if (integrationScheme ==
311 auto is = lsInternal::LaxFriedrichs<T, D, 2>(levelSets.back(), velocities,
312 dissipationAlpha,
313 calculateNormalVectors);
314 currentTime = integrateTime(is, maxTimeStep);
315 } else if (integrationScheme ==
319 levelSets.back(), velocities);
320 currentTime = integrateTime(is, maxTimeStep);
321 } else if (integrationScheme ==
324 levelSets.back(), velocities, dissipationAlpha);
325 currentTime = integrateTime(is, maxTimeStep);
326 } else if (integrationScheme ==
329 levelSets.back(), velocities, dissipationAlpha);
330 currentTime = integrateTime(is, maxTimeStep);
331 } else if (integrationScheme ==
334 levelSets.back(), velocities, dissipationAlpha);
335 currentTime = integrateTime(is, maxTimeStep);
336 } else if (integrationScheme ==
339 levelSets.back(), velocities, dissipationAlpha);
340 currentTime = integrateTime(is, maxTimeStep);
341 } else if (integrationScheme ==
344 levelSets.back(), velocities, dissipationAlpha);
345 currentTime = integrateTime(is, maxTimeStep);
346 } else {
347 Logger::getInstance()
348 .addWarning("Advect: Integration scheme not found. Not advecting.")
349 .print();
350
351 // if no correct scheme was found return infinity
352 // to stop advection
353 return std::numeric_limits<double>::max();
354 }
355
356 rebuildLS();
357
358 // Adjust all level sets below the advected one
359 // This means, that when the top level-set and one below
360 // are etched, the lower one is moved with the top level-set
361 // TODO: Adjust lower layers also when they have grown,
362 // to allow for two different growth rates of materials
363 if (integrationScheme !=
365 for (unsigned i = 0; i < levelSets.size() - 1; ++i) {
366 BooleanOperation<T, D>(levelSets[i], levelSets.back(),
368 .apply();
369 }
370 }
371
372 return currentTime;
373 }
374
379 template <class IntegrationSchemeType>
380 double
381 integrateTime(IntegrationSchemeType IntegrationScheme,
382 double maxTimeStep = std::numeric_limits<double>::max()) {
383 if (timeStepRatio >= 0.5) {
384 Logger::getInstance()
385 .addWarning("Integration time step ratio should be smaller than 0.5. "
386 "Advection might fail!")
387 .print();
388 }
389
390 auto &topDomain = levelSets.back()->getDomain();
391 auto &grid = levelSets.back()->getGrid();
392
393 std::vector<std::vector<std::pair<T, T>>> totalTempRates;
394 totalTempRates.resize((levelSets.back())->getNumberOfSegments());
395
396 typename PointData<T>::ScalarDataType *voidMarkerPointer;
397 if (ignoreVoids) {
398 MarkVoidPoints<T, D>(levelSets.back()).apply();
399 auto &pointData = levelSets.back()->getPointData();
400 voidMarkerPointer =
402 if (voidMarkerPointer == nullptr) {
403 Logger::getInstance()
404 .addWarning("Advect: Cannot find void point markers. Not "
405 "ignoring void points.")
406 .print();
407 ignoreVoids = false;
408 }
409 }
410
411 const bool ignoreVoidPoints = ignoreVoids;
412
413#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
414 {
415 int p = 0;
416#ifdef _OPENMP
417 p = omp_get_thread_num();
418#endif
419
420 hrleVectorType<hrleIndexType, D> startVector =
421 (p == 0) ? grid.getMinGridPoint()
422 : topDomain.getSegmentation()[p - 1];
423
424 hrleVectorType<hrleIndexType, D> endVector =
425 (p != static_cast<int>(topDomain.getNumberOfSegments() - 1))
426 ? topDomain.getSegmentation()[p]
427 : grid.incrementIndices(grid.getMaxGridPoint());
428
429 double tempMaxTimeStep = maxTimeStep;
430 auto &tempRates = totalTempRates[p];
431 tempRates.reserve(topDomain.getNumberOfPoints() /
432 double((levelSets.back())->getNumberOfSegments()) +
433 10);
434
435 // an iterator for each level set
436 std::vector<hrleSparseIterator<typename Domain<T, D>::DomainType>>
437 iterators;
438 for (auto it = levelSets.begin(); it != levelSets.end(); ++it) {
439 iterators.push_back(
440 hrleSparseIterator<typename Domain<T, D>::DomainType>(
441 (*it)->getDomain()));
442 }
443
444 IntegrationSchemeType scheme(IntegrationScheme);
445
446 for (hrleSparseIterator<typename Domain<T, D>::DomainType> it(
447 topDomain, startVector);
448 it.getStartIndices() < endVector; ++it) {
449
450 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
451 continue;
452
453 T value = it.getValue();
454 double maxStepTime = 0;
455 double cfl = timeStepRatio;
456
457 for (int currentLevelSetId = levelSets.size() - 1;
458 currentLevelSetId >= 0; --currentLevelSetId) {
459
460 T velocity = 0;
461
462 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
463 // check if there is any other levelset at the same point:
464 // if yes, take the velocity of the lowest levelset
465 for (unsigned lowerLevelSetId = 0;
466 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
467 // put iterator to same position as the top levelset
468 iterators[lowerLevelSetId].goToIndicesSequential(
469 it.getStartIndices());
470
471 // if the lower surface is actually outside, i.e. its LS value
472 // is lower or equal
473 if (iterators[lowerLevelSetId].getValue() <=
474 value + wrappingLayerEpsilon) {
475 velocity = scheme(it.getStartIndices(), lowerLevelSetId);
476 break;
477 }
478 }
479 }
480
481 T valueBelow;
482 // get value of material below (before in levelSets list)
483 if (currentLevelSetId > 0) {
484 iterators[currentLevelSetId - 1].goToIndicesSequential(
485 it.getStartIndices());
486 valueBelow = iterators[currentLevelSetId - 1].getValue();
487 } else {
488 valueBelow = std::numeric_limits<T>::max();
489 }
490
491 // if velocity is positive, set maximum time step possible without
492 // violating the cfl condition
493 if (velocity > 0.) {
494 maxStepTime += cfl / velocity;
495 tempRates.push_back(
496 std::make_pair(velocity, -std::numeric_limits<T>::max()));
497 break;
498 // if velocity is 0, maximum time step is infinite
499 } else if (velocity == 0.) {
500 maxStepTime = std::numeric_limits<T>::max();
501 tempRates.push_back(
502 std::make_pair(velocity, std::numeric_limits<T>::max()));
503 break;
504 // if the velocity is negative apply the velocity for as long as
505 // possible without infringing on material below
506 } else {
507 T difference = std::abs(valueBelow - value);
508
509 if (difference >= cfl) {
510 maxStepTime -= cfl / velocity;
511 tempRates.push_back(
512 std::make_pair(velocity, std::numeric_limits<T>::max()));
513 break;
514 } else {
515 maxStepTime -= difference / velocity;
516 // the second part of the pair indicates how far we can move
517 // in this time step until the end of the material is reached
518 tempRates.push_back(std::make_pair(velocity, valueBelow));
519 cfl -= difference;
520 // use new LS value for next calculations
521 value = valueBelow;
522 }
523 }
524 }
525
526 if (maxStepTime < tempMaxTimeStep)
527 tempMaxTimeStep = maxStepTime;
528 }
529
530#pragma omp critical
531 {
532 // If scheme is STENCIL_LOCAL_LAX_FRIEDRICHS the time step is
533 // reduced depending on the dissipation coefficients For all
534 // remaining schemes this function is empty.
536 T, D>(
537 scheme, tempMaxTimeStep,
538 levelSets.back()->getGrid().getGridDelta());
539
540 // set global timestep maximum
541 if (tempMaxTimeStep < maxTimeStep)
542 maxTimeStep = tempMaxTimeStep;
543 }
544 }
545
546 // reduce to one layer thickness and apply new values directly to the
547 // domain segments --> DO NOT CHANGE SEGMENTATION HERE (true parameter)
548 Reduce<T, D>(levelSets.back(), 1, true).apply();
549
550 const bool saveVelocities = saveAdvectionVelocities;
551 std::vector<std::vector<double>> velocityVectors(
552 levelSets.back()->getNumberOfSegments());
553
554#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
555 {
556 int p = 0;
557#ifdef _OPENMP
558 p = omp_get_thread_num();
559#endif
560
561 typename std::vector<std::pair<T, T>>::const_iterator itRS =
562 totalTempRates[p].begin();
563 auto &segment = topDomain.getDomainSegment(p);
564 unsigned maxId = segment.getNumberOfPoints();
565
566 if (saveVelocities) {
567 velocityVectors[p].resize(maxId);
568 }
569
570 for (unsigned localId = 0; localId < maxId; ++localId) {
571 T &value = segment.definedValues[localId];
572 double time = maxTimeStep;
573
574 // if there is a change in materials during one time step, deduct
575 // the time taken to advect up to the end of the top material and
576 // set the LS value to the one below
577 while (std::abs(itRS->second - value) < std::abs(time * itRS->first)) {
578 time -= std::abs((itRS->second - value) / itRS->first);
579 value = itRS->second;
580 ++itRS;
581 }
582
583 // now deduct the velocity times the time step we take
584 value -= time * itRS->first;
585
586 if (saveVelocities) {
587 velocityVectors[p][localId] = time * itRS->first;
588 }
589
590 // this
591 // is run when two materials are close but the velocity is too slow
592 // to actually reach the second material, to get rid of the extra
593 // entry in the TempRatesStop
594 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
595 ++itRS;
596
597 // advance the TempStopRates iterator by one
598 ++itRS;
599 }
600 }
601
602 if (saveVelocities) {
603 auto &pointData = levelSets.back()->getPointData();
604 // delete if already exists
605 if (int i = pointData.getScalarDataIndex(velocityLabel); i != -1) {
606 pointData.eraseScalarData(i);
607 }
608
609 typename PointData<T>::ScalarDataType vels;
610
611 for (unsigned i = 0; i < velocityVectors.size(); ++i) {
612 vels.insert(vels.end(),
613 std::make_move_iterator(velocityVectors[i].begin()),
614 std::make_move_iterator(velocityVectors[i].end()));
615 }
616 pointData.insertNextScalarData(std::move(vels), velocityLabel);
617 }
618
619 return maxTimeStep;
620 }
621
622public:
623 static constexpr char velocityLabel[] = "AdvectionVelocities";
624
626
627 Advect(SmartPointer<Domain<T, D>> passedlsDomain) {
628 levelSets.push_back(passedlsDomain);
629 }
630
631 Advect(SmartPointer<Domain<T, D>> passedlsDomain,
632 SmartPointer<VelocityField<T>> passedVelocities) {
633 levelSets.push_back(passedlsDomain);
634 velocities = passedVelocities;
635 }
636
637 Advect(std::vector<SmartPointer<Domain<T, D>>> passedlsDomains,
638 SmartPointer<VelocityField<T>> passedVelocities)
639 : levelSets(passedlsDomains) {
640 velocities = passedVelocities;
641 }
642
645 void insertNextLevelSet(SmartPointer<Domain<T, D>> passedlsDomain) {
646 levelSets.push_back(passedlsDomain);
647 }
648
651 void setVelocityField(SmartPointer<VelocityField<T>> passedVelocities) {
652 velocities = passedVelocities;
653 }
654
660 void setAdvectionTime(double time) { advectionTime = time; }
661
666 void setSingleStep(bool singleStep) { performOnlySingleStep = singleStep; }
667
672 void setTimeStepRatio(const double &cfl) { timeStepRatio = cfl; }
673
678 void setCalculateNormalVectors(bool cnv) { calculateNormalVectors = cnv; }
679
686 void setIgnoreVoids(bool iV) { ignoreVoids = iV; }
687
690 void setSaveAdvectionVelocities(bool sAV) { saveAdvectionVelocities = sAV; }
691
694 double getAdvectedTime() { return advectedTime; }
695
697 unsigned getNumberOfTimeSteps() { return numberOfTimeSteps; }
698
700 double getTimeStepRatio() { return timeStepRatio; }
701
703 bool getCalculateNormalVectors() { return calculateNormalVectors; }
704
708 integrationScheme = scheme;
709 }
710
715 void setDissipationAlpha(const double &a) { dissipationAlpha = a; }
716
719 void setUpdatePointData(bool update) { updatePointData = update; }
720
721 // Prepare the levelset for advection, based on the provided integration
722 // scheme.
723 void prepareLS() {
724 // check whether a level set and velocities have been given
725 if (levelSets.size() < 1) {
726 Logger::getInstance()
727 .addWarning("No level sets passed to Advect.")
728 .print();
729 return;
730 }
731
732 if (integrationScheme == IntegrationSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
734 } else if (integrationScheme ==
737 } else if (integrationScheme ==
740 } else if (integrationScheme ==
743 } else if (integrationScheme ==
747 levelSets.back());
748 } else if (integrationScheme ==
751 } else if (integrationScheme ==
754 } else if (integrationScheme ==
757 } else if (integrationScheme ==
760 } else if (integrationScheme ==
763 levelSets.back());
764 } else {
765 Logger::getInstance()
766 .addWarning("Advect: Integration scheme not found. Not advecting.")
767 .print();
768 }
769 }
770
772 void apply() {
773 if (advectionTime == 0.) {
774 advectedTime = advect();
775 numberOfTimeSteps = 1;
776 } else {
777 double currentTime = 0.0;
778 numberOfTimeSteps = 0;
779 while (currentTime < advectionTime) {
780 currentTime += advect(advectionTime - currentTime);
781 ++numberOfTimeSteps;
782 if (performOnlySingleStep)
783 break;
784 }
785 advectedTime = currentTime;
786 }
787 }
788};
789
790// add all template specializations for this class
792
793} // namespace viennals
Engquist osher integration scheme based on the upwind integration scheme. Offers high performance but...
Definition lsEnquistOsher.hpp:19
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsEnquistOsher.hpp:28
Lax Friedrichs integration scheme with constant alpha value for dissipation. This alpha value should ...
Definition lsLaxFriedrichs.hpp:19
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLaxFriedrichs.hpp:29
Lax Friedrichs integration scheme, which uses alpha values provided by the user in getDissipationAlph...
Definition lsLocalLaxFriedrichsAnalytical.hpp:21
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLaxFriedrichsAnalytical.hpp:46
Lax Friedrichs integration scheme, which uses a first neighbour stencil to calculate the alpha values...
Definition lsLocalLaxFriedrichs.hpp:21
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLaxFriedrichs.hpp:47
Lax Friedrichs integration scheme, which considers only the current point for alpha calculation....
Definition lsLocalLocalLaxFriedrichs.hpp:19
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLocalLaxFriedrichs.hpp:28
Stencil Local Lax Friedrichs Integration Scheme. It uses a stencil of order around active points,...
Definition lsStencilLocalLaxFriedrichsScalar.hpp:21
static void prepareLS(LevelSetType passedlsDomain)
Definition lsStencilLocalLaxFriedrichsScalar.hpp:133
This class is used to advance level sets over time. Level sets are passed to the constructor in an st...
Definition lsAdvect.hpp:74
void setVelocityField(SmartPointer< VelocityField< T > > passedVelocities)
Set the velocity field used for advection. This should be a concrete implementation of lsVelocityFiel...
Definition lsAdvect.hpp:651
static constexpr char velocityLabel[]
Definition lsAdvect.hpp:623
void insertNextLevelSet(SmartPointer< Domain< T, D > > passedlsDomain)
Pushes the passed level set to the back of the list of level sets used for advection.
Definition lsAdvect.hpp:645
void setUpdatePointData(bool update)
Set whether the point data in the old LS should be translated to the advected LS. Defaults to true.
Definition lsAdvect.hpp:719
bool getCalculateNormalVectors()
Get whether normal vectors were caluclated.
Definition lsAdvect.hpp:703
void setIgnoreVoids(bool iV)
Set whether level set values, which are not part of the "top" geometrically connected part of values,...
Definition lsAdvect.hpp:686
void setCalculateNormalVectors(bool cnv)
Set whether normal vectors should be calculated at each level set point. Defaults to true....
Definition lsAdvect.hpp:678
void apply()
Perform the advection.
Definition lsAdvect.hpp:772
void setSingleStep(bool singleStep)
If set to true, only a single advection step will be performed, even if the advection time set with s...
Definition lsAdvect.hpp:666
unsigned getNumberOfTimeSteps()
Get how many advection steps were performed during the last apply() call.
Definition lsAdvect.hpp:697
double getTimeStepRatio()
Get the value of the CFL number.
Definition lsAdvect.hpp:700
void prepareLS()
Definition lsAdvect.hpp:723
void setDissipationAlpha(const double &a)
Set the alpha dissipation coefficient. For lsLaxFriedrichs, this is used as the alpha value....
Definition lsAdvect.hpp:715
Advect(std::vector< SmartPointer< Domain< T, D > > > passedlsDomains, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:637
void setIntegrationScheme(IntegrationSchemeEnum scheme)
Set which integration scheme should be used out of the ones specified in IntegrationSchemeEnum.
Definition lsAdvect.hpp:707
double getAdvectedTime()
Get by how much the physical time was advanced during the last apply() call.
Definition lsAdvect.hpp:694
void setAdvectionTime(double time)
Set the time until when the level set should be advected. If this takes more than one advection step,...
Definition lsAdvect.hpp:660
void setTimeStepRatio(const double &cfl)
Set the CFL condition to use during advection. The CFL condition sets the maximum distance a surface ...
Definition lsAdvect.hpp:672
Advect()
Definition lsAdvect.hpp:625
Advect(SmartPointer< Domain< T, D > > passedlsDomain)
Definition lsAdvect.hpp:627
Advect(SmartPointer< Domain< T, D > > passedlsDomain, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:631
void setSaveAdvectionVelocities(bool sAV)
Set whether the velocities applied to each point should be saved in the level set for debug purposes.
Definition lsAdvect.hpp:690
This class is used to perform boolean operations on two level sets and write the resulting level set ...
Definition lsBooleanOperation.hpp:43
void apply()
Perform operation.
Definition lsBooleanOperation.hpp:312
Class containing all information about the level set, including the dimensions of the domain,...
Definition lsDomain.hpp:28
void deepCopy(const SmartPointer< Domain< T, D > > passedDomain)
copy all values of "passedDomain" to this Domain
Definition lsDomain.hpp:121
unsigned getNumberOfSegments() const
returns the number of segments, the levelset is split into. This is useful for algorithm parallelisat...
Definition lsDomain.hpp:150
hrleDomain< T, D > DomainType
Definition lsDomain.hpp:33
void finalize(int newWidth)
this function sets a new levelset width and finalizes the levelset, so it is ready for use by other a...
Definition lsDomain.hpp:114
This class is used to mark points of the level set which are enclosed in a void.
Definition lsMarkVoidPoints.hpp:30
void apply()
Definition lsMarkVoidPoints.hpp:154
std::vector< T > ScalarDataType
Definition lsPointData.hpp:22
ScalarDataType * getScalarData(int index)
Definition lsPointData.hpp:90
Reduce the level set size to the specified width. This means all level set points with value <= 0....
Definition lsReduce.hpp:16
void apply()
Reduces the leveleSet to the specified number of layers. The largest value in the levelset is thus wi...
Definition lsReduce.hpp:57
Extract the regular grid, on which the level set values are defined, to an explicit Mesh<>....
Definition lsToMesh.hpp:19
void apply()
Definition lsToMesh.hpp:49
Class handling the output of an Mesh<> to VTK file types.
Definition lsVTKWriter.hpp:31
void apply()
Definition lsVTKWriter.hpp:88
Abstract class defining the interface for the velocity field used during advection using lsAdvect.
Definition lsVelocityField.hpp:13
#define PRECOMPILE_PRECISION_DIMENSION(className)
Definition lsPreCompileMacros.hpp:24
std::enable_if_t<!std::is_same< A, B >::value, AssignType > IsNotSame
Definition lsConcepts.hpp:20
constexpr AssignType assignable
Definition lsConcepts.hpp:10
Definition lsAdvect.hpp:37
void reduceTimeStepHamiltonJacobi(IntegrationSchemeType &, double &, hrleCoordType)
Definition lsAdvect.hpp:42
Definition lsAdvect.hpp:46
IntegrationSchemeEnum
Enumeration for the different Integration schemes used by the advection kernel.
Definition lsAdvect.hpp:52
constexpr int D
Definition pyWrap.cpp:65
double T
Definition pyWrap.cpp:63