54template <
class T,
int D>
class Advect {
55 using ConstSparseIterator =
56 viennahrle::ConstSparseIterator<typename Domain<T, D>::DomainType>;
57 using hrleIndexType = viennahrle::IndexType;
62 std::vector<SmartPointer<Domain<T, D>>> levelSets;
63 SmartPointer<VelocityField<T>> velocities =
nullptr;
66 double timeStepRatio = 0.4999;
67 double dissipationAlpha = 1.0;
68 bool calculateNormalVectors =
true;
69 bool ignoreVoids =
false;
70 double advectionTime = 0.;
71 bool performOnlySingleStep =
false;
72 double advectedTime = 0.;
73 unsigned numberOfTimeSteps = 0;
74 bool saveAdvectionVelocities =
false;
75 bool updatePointData =
true;
76 bool checkDissipation =
true;
77 double integrationCutoff = 0.5;
78 bool adaptiveTimeStepping =
false;
79 unsigned adaptiveTimeStepSubdivisions = 20;
80 static constexpr double wrappingLayerEpsilon = 1e-4;
81 std::vector<SmartPointer<Domain<T, D>>> initialLevelSets;
82 std::function<bool(SmartPointer<
Domain<T, D>>)> velocityUpdateCallback =
87 std::vector<std::vector<std::pair<std::pair<T, T>,
T>>> storedRates;
88 double currentTimeStep = -1.;
90 VectorType<T, D> findGlobalAlphas()
const {
92 auto &topDomain = levelSets.back()->getDomain();
93 auto &grid = levelSets.back()->getGrid();
95 const T gridDelta = grid.getGridDelta();
96 const T deltaPos = gridDelta;
97 const T deltaNeg = -gridDelta;
99 VectorType<T, D> finalAlphas{};
101#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
103 VectorType<T, D> localAlphas{};
106 p = omp_get_thread_num();
108 viennahrle::Index<D> startVector =
109 (p == 0) ? grid.getMinGridPoint()
110 : topDomain.getSegmentation()[p - 1];
111 viennahrle::Index<D> endVector =
112 (p !=
static_cast<int>(topDomain.getNumberOfSegments() - 1))
113 ? topDomain.getSegmentation()[p]
114 : grid.incrementIndices(grid.getMaxGridPoint());
117 std::vector<ConstSparseIterator> iterators;
118 for (
auto const &
ls : levelSets) {
119 iterators.emplace_back(
ls->getDomain());
123 viennahrle::ConstSparseStarIterator<typename Domain<T, D>::DomainType, 1>
124 neighborIterator(topDomain);
126 for (ConstSparseIterator it(topDomain, startVector);
127 it.getStartIndices() < endVector; ++it) {
129 if (!it.isDefined() || std::abs(it.getValue()) > integrationCutoff)
132 const T value = it.getValue();
133 const auto indices = it.getStartIndices();
137 for (
unsigned lowerLevelSetId = 0; lowerLevelSetId < levelSets.size();
140 iterators[lowerLevelSetId].goToIndicesSequential(indices);
144 if (iterators[lowerLevelSetId].getValue() <=
145 value + wrappingLayerEpsilon) {
148 neighborIterator.goToIndicesSequential(indices);
151 for (
unsigned i = 0; i <
D; ++i) {
152 coords[i] = indices[i] * gridDelta;
156 T normalModulus = 0.;
157 for (
unsigned i = 0; i <
D; ++i) {
158 const T phiPos = neighborIterator.getNeighbor(i).getValue();
159 const T phiNeg = neighborIterator.getNeighbor(i +
D).getValue();
161 normal[i] = phiPos - phiNeg;
162 normalModulus += normal[i] * normal[i];
164 normalModulus = 1. / std::sqrt(normalModulus);
165 for (
unsigned i = 0; i <
D; ++i)
166 normal[i] *= normalModulus;
168 T scaVel = velocities->getScalarVelocity(
169 coords, lowerLevelSetId, normal,
170 neighborIterator.getCenter().getPointId());
171 auto vecVel = velocities->getVectorVelocity(
172 coords, lowerLevelSetId, normal,
173 neighborIterator.getCenter().getPointId());
175 for (
unsigned i = 0; i <
D; ++i) {
176 T tempAlpha = std::abs((scaVel + vecVel[i]) * normal[i]);
177 localAlphas[i] = std::max(localAlphas[i], tempAlpha);
188 for (
unsigned i = 0; i <
D; ++i) {
189 finalAlphas[i] = std::max(finalAlphas[i], localAlphas[i]);
199 bool combineLevelSets(
T wTarget,
T wSource) {
209 int expansionWidth = std::ceil(2.0 * steps * timeStepRatio + 1);
213 bool movedDown =
false;
216 initialLevelSets.back(),
219 [wTarget, wSource, &movedDown](
const T &a,
const T &b) {
222 T res = wSource * a + wTarget * b;
225 return std::make_pair(res,
true);
228 return std::make_pair(a,
false);
229 return std::make_pair(b,
false);
243 auto &grid = levelSets.back()->getGrid();
244 auto newlsDomain = SmartPointer<Domain<T, D>>::New(grid);
245 auto &newDomain = newlsDomain->getDomain();
246 auto &domain = levelSets.back()->getDomain();
258 newDomain.initialize(domain.getNewSegmentation(),
259 domain.getAllocation() *
260 (2.0 / levelSets.back()->getLevelSetWidth()));
262 const bool updateData = updatePointData;
265 std::vector<std::vector<unsigned>> newDataSourceIds;
267 newDataSourceIds.resize(newDomain.getNumberOfSegments());
269#ifdef DEBUG_LS_ADVECT_HPP
271 auto mesh = SmartPointer<Mesh<T>>::New();
277#pragma omp parallel num_threads(newDomain.getNumberOfSegments())
281 p = omp_get_thread_num();
283 auto &domainSegment = newDomain.getDomainSegment(p);
285 viennahrle::Index<D> startVector =
286 (p == 0) ? grid.getMinGridPoint()
287 : newDomain.getSegmentation()[p - 1];
289 viennahrle::Index<D> endVector =
290 (p !=
static_cast<int>(newDomain.getNumberOfSegments() - 1))
291 ? newDomain.getSegmentation()[p]
292 : grid.incrementIndices(grid.getMaxGridPoint());
297 newDataSourceIds[p].reserve(2.5 * domainSegment.getNumberOfPoints());
299 for (viennahrle::ConstSparseStarIterator<
301 it(domain, startVector);
302 it.getIndices() < endVector; ++it) {
306 if (std::abs(it.getCenter().getValue()) <= 1.0) {
309 for (; k < 2 *
D; k++)
310 if (std::signbit(it.getNeighbor(k).getValue() - 1e-7) !=
311 std::signbit(it.getCenter().getValue() + 1e-7))
316 if (it.getCenter().getDefinedValue() > 0.5) {
318 for (; j < 2 *
D; j++) {
319 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
320 if (it.getNeighbor(j).getDefinedValue() < -0.5)
324 domainSegment.insertNextDefinedPoint(
325 it.getIndices(), it.getCenter().getDefinedValue());
327 newDataSourceIds[p].push_back(it.getCenter().getPointId());
330 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
332 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
334 }
else if (it.getCenter().getDefinedValue() < -0.5) {
336 for (; j < 2 *
D; j++) {
337 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
338 if (it.getNeighbor(j).getDefinedValue() > 0.5)
343 domainSegment.insertNextDefinedPoint(
344 it.getIndices(), it.getCenter().getDefinedValue());
346 newDataSourceIds[p].push_back(it.getCenter().getPointId());
349 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
351 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
354 domainSegment.insertNextDefinedPoint(
355 it.getIndices(), it.getCenter().getDefinedValue());
357 newDataSourceIds[p].push_back(it.getCenter().getPointId());
360 domainSegment.insertNextUndefinedPoint(
361 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
367 if (it.getCenter().getValue() >= 0) {
368 int usedNeighbor = -1;
370 for (
int i = 0; i < 2 *
D; i++) {
371 T value = it.getNeighbor(i).getValue();
372 if (std::abs(value) <= 1.0 && (value < 0.)) {
373 if (distance > value + 1.0) {
374 distance = value + 1.0;
380 if (distance <= cutoff) {
381 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
383 newDataSourceIds[p].push_back(
384 it.getNeighbor(usedNeighbor).getPointId());
386 domainSegment.insertNextUndefinedPoint(it.getIndices(),
391 int usedNeighbor = -1;
393 for (
int i = 0; i < 2 *
D; i++) {
394 T value = it.getNeighbor(i).getValue();
395 if (std::abs(value) <= 1.0 && (value > 0)) {
396 if (distance < value - 1.0) {
398 distance = value - 1.0;
404 if (distance >= -cutoff) {
405 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
407 newDataSourceIds[p].push_back(
408 it.getNeighbor(usedNeighbor).getPointId());
410 domainSegment.insertNextUndefinedPoint(it.getIndices(),
420 auto &pointData = levelSets.back()->getPointData();
421 newlsDomain->getPointData().translateFromMultiData(pointData,
425 newDomain.finalize();
427 levelSets.back()->deepCopy(newlsDomain);
428 levelSets.back()->finalize(finalWidth);
435 template <
class DiscretizationSchemeType>
436 double integrateTime(DiscretizationSchemeType spatialScheme,
437 double maxTimeStep) {
439 auto &topDomain = levelSets.back()->getDomain();
440 auto &grid = levelSets.back()->getGrid();
445 auto &pointData = levelSets.back()->getPointData();
448 if (voidMarkerPointer ==
nullptr) {
449 VIENNACORE_LOG_WARNING(
"Advect: Cannot find void point markers. Not "
450 "ignoring void points.");
454 const bool ignoreVoidPoints = ignoreVoids;
455 const bool useAdaptiveTimeStepping = adaptiveTimeStepping;
456 const auto adaptiveFactor = 1.0 / adaptiveTimeStepSubdivisions;
458 if (!storedRates.empty()) {
459 VIENNACORE_LOG_WARNING(
"Advect: Overwriting previously stored rates.");
462 storedRates.resize(topDomain.getNumberOfSegments());
464#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
468 p = omp_get_thread_num();
470 viennahrle::Index<D> startVector =
471 (p == 0) ? grid.getMinGridPoint()
472 : topDomain.getSegmentation()[p - 1];
474 viennahrle::Index<D> endVector =
475 (p !=
static_cast<int>(topDomain.getNumberOfSegments() - 1))
476 ? topDomain.getSegmentation()[p]
477 : grid.incrementIndices(grid.getMaxGridPoint());
479 double tempMaxTimeStep = maxTimeStep;
481 auto &rates = storedRates[p];
483 topDomain.getNumberOfPoints() /
484 static_cast<double>((levelSets.back())->getNumberOfSegments()) +
488 std::vector<ConstSparseIterator> iterators;
489 for (
auto const &
ls : levelSets) {
490 iterators.emplace_back(
ls->getDomain());
493 DiscretizationSchemeType scheme(spatialScheme);
495 for (ConstSparseIterator it(topDomain, startVector);
496 it.getStartIndices() < endVector; it.next()) {
498 if (!it.isDefined() || std::abs(it.getValue()) > integrationCutoff)
501 T value = it.getValue();
502 double maxStepTime = 0;
503 double cfl = timeStepRatio;
505 for (
int currentLevelSetId = levelSets.size() - 1;
506 currentLevelSetId >= 0; --currentLevelSetId) {
508 std::pair<T, T> gradNDissipation;
510 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
513 for (
unsigned lowerLevelSetId = 0;
514 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
516 iterators[lowerLevelSetId].goToIndicesSequential(
517 it.getStartIndices());
521 if (iterators[lowerLevelSetId].getValue() <=
522 value + wrappingLayerEpsilon) {
524 scheme(it.getStartIndices(), lowerLevelSetId);
530 T velocity = gradNDissipation.first - gradNDissipation.second;
534 maxStepTime += cfl / velocity;
535 rates.emplace_back(gradNDissipation,
536 -std::numeric_limits<T>::max());
538 }
else if (velocity == 0.) {
541 maxStepTime = std::numeric_limits<T>::max();
542 rates.emplace_back(gradNDissipation, std::numeric_limits<T>::max());
548 if (currentLevelSetId > 0) {
549 iterators[currentLevelSetId - 1].goToIndicesSequential(
550 it.getStartIndices());
551 valueBelow = iterators[currentLevelSetId - 1].getValue();
553 valueBelow = std::numeric_limits<T>::max();
556 T difference = std::abs(valueBelow - value);
558 if (difference >= cfl) {
561 maxStepTime -= cfl / velocity;
562 rates.emplace_back(gradNDissipation,
563 std::numeric_limits<T>::max());
568 if (useAdaptiveTimeStepping &&
569 difference > adaptiveFactor * cfl) {
574 maxStepTime -= adaptiveFactor * cfl / velocity;
575 rates.emplace_back(gradNDissipation,
576 std::numeric_limits<T>::max());
584 maxStepTime -= difference / velocity;
585 rates.emplace_back(gradNDissipation, valueBelow);
591 if (maxStepTime < tempMaxTimeStep)
592 tempMaxTimeStep = maxStepTime;
600 scheme.reduceTimeStepHamiltonJacobi(
601 tempMaxTimeStep, levelSets.back()->getGrid().getGridDelta());
604 if (tempMaxTimeStep < maxTimeStep)
605 maxTimeStep = tempMaxTimeStep;
616 void computeRates(
double maxTimeStep = std::numeric_limits<double>::max()) {
620 calculateNormalVectors);
621 currentTimeStep = integrateTime(is, maxTimeStep);
624 calculateNormalVectors);
625 currentTimeStep = integrateTime(is, maxTimeStep);
627 auto alphas = findGlobalAlphas();
629 dissipationAlpha, alphas,
630 calculateNormalVectors);
631 currentTimeStep = integrateTime(is, maxTimeStep);
633 auto alphas = findGlobalAlphas();
635 dissipationAlpha, alphas,
636 calculateNormalVectors);
637 currentTimeStep = integrateTime(is, maxTimeStep);
638 }
else if (spatialScheme ==
641 levelSets.back(), velocities);
642 currentTimeStep = integrateTime(is, maxTimeStep);
643 }
else if (spatialScheme ==
646 levelSets.back(), velocities, dissipationAlpha);
647 currentTimeStep = integrateTime(is, maxTimeStep);
648 }
else if (spatialScheme ==
651 levelSets.back(), velocities, dissipationAlpha);
652 currentTimeStep = integrateTime(is, maxTimeStep);
653 }
else if (spatialScheme ==
656 levelSets.back(), velocities, dissipationAlpha);
657 currentTimeStep = integrateTime(is, maxTimeStep);
658 }
else if (spatialScheme ==
661 levelSets.back(), velocities, dissipationAlpha);
662 currentTimeStep = integrateTime(is, maxTimeStep);
663 }
else if (spatialScheme ==
666 levelSets.back(), velocities, dissipationAlpha);
667 currentTimeStep = integrateTime(is, maxTimeStep);
672 currentTimeStep = integrateTime(is, maxTimeStep);
677 currentTimeStep = integrateTime(is, maxTimeStep);
679 VIENNACORE_LOG_ERROR(
"Advect: Discretization scheme not found.");
680 currentTimeStep = -1.;
686 void updateLevelSet(
double dt) {
687 if (timeStepRatio >= 0.5) {
688 VIENNACORE_LOG_WARNING(
689 "Integration time step ratio should be smaller than 0.5. "
690 "Advection might fail!");
693 auto &topDomain = levelSets.back()->getDomain();
695 assert(dt >= 0. &&
"No time step set!");
696 assert(storedRates.size() == topDomain.getNumberOfSegments());
702 const bool saveVelocities = saveAdvectionVelocities;
703 std::vector<std::vector<double>> dissipationVectors(
704 levelSets.back()->getNumberOfSegments());
705 std::vector<std::vector<double>> velocityVectors(
706 levelSets.back()->getNumberOfSegments());
708 const bool checkDiss = checkDissipation;
710#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
714 p = omp_get_thread_num();
716 auto itRS = storedRates[p].cbegin();
717 auto &segment = topDomain.getDomainSegment(p);
718 const unsigned maxId = segment.getNumberOfPoints();
720 if (saveVelocities) {
721 velocityVectors[p].resize(maxId);
722 dissipationVectors[p].resize(maxId);
725 for (
unsigned localId = 0; localId < maxId; ++localId) {
726 T &value = segment.definedValues[localId];
729 if (std::abs(value) > integrationCutoff)
737 auto const [gradient, dissipation] = itRS->first;
738 T velocity = gradient - dissipation;
741 if (checkDiss && (gradient < 0 && velocity > 0) ||
742 (gradient > 0 && velocity < 0)) {
746 T rate = time * velocity;
747 while (std::abs(itRS->second - value) < std::abs(rate)) {
748 time -= std::abs((itRS->second - value) / velocity);
749 value = itRS->second;
753 velocity = itRS->first.first - itRS->first.second;
754 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
755 (itRS->first.first > 0 && velocity < 0)) {
758 rate = time * velocity;
764 if (saveVelocities) {
765 velocityVectors[p][localId] = rate;
766 dissipationVectors[p][localId] = itRS->first.second;
772 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
780 if (saveVelocities) {
781 auto &pointData = levelSets.back()->getPointData();
786 for (
unsigned i = 0; i < velocityVectors.size(); ++i) {
787 vels.insert(vels.end(),
788 std::make_move_iterator(velocityVectors[i].begin()),
789 std::make_move_iterator(velocityVectors[i].end()));
790 diss.insert(diss.end(),
791 std::make_move_iterator(dissipationVectors[i].begin()),
792 std::make_move_iterator(dissipationVectors[i].end()));
794 pointData.insertReplaceScalarData(std::move(vels),
velocityLabel);
802 void adjustLowerLayers() {
806 for (
unsigned i = 0; i < levelSets.size() - 1; ++i) {
808 levelSets[i], levelSets.back(),
817 double advect(
double maxTimeStep) {
818 switch (temporalScheme) {
839 levelSets.push_back(passedlsDomain);
844 levelSets.push_back(passedlsDomain);
845 velocities = passedVelocities;
850 : levelSets(passedlsDomains) {
851 velocities = passedVelocities;
857 levelSets.push_back(passedlsDomain);
865 velocities = passedVelocities;
905 adaptiveTimeStepping = aTS;
906 if (subdivisions < 1) {
907 VIENNACORE_LOG_WARNING(
"Advect: Adaptive time stepping subdivisions must "
908 "be at least 1. Setting to 1.");
911 adaptiveTimeStepSubdivisions = subdivisions;
940 [[deprecated(
"Use setSpatialScheme instead")]]
void
942 VIENNACORE_LOG_WARNING(
943 "Advect::setIntegrationScheme is deprecated and will be removed in "
944 "future versions. Use setSpatialScheme instead.");
945 spatialScheme = scheme;
967 std::function<
bool(SmartPointer<
Domain<T, D>>)> callback) {
968 velocityUpdateCallback = callback;
975 if (levelSets.empty()) {
976 VIENNACORE_LOG_ERROR(
"No level sets passed to Advect.");
988 }
else if (spatialScheme ==
992 }
else if (spatialScheme ==
995 }
else if (spatialScheme ==
998 }
else if (spatialScheme ==
1001 }
else if (spatialScheme ==
1004 }
else if (spatialScheme ==
1013 VIENNACORE_LOG_ERROR(
"Advect: Discretization scheme not found.");
1019 if (levelSets.empty()) {
1020 VIENNACORE_LOG_ERROR(
"No level sets passed to Advect. Not advecting.");
1023 if (velocities ==
nullptr) {
1024 VIENNACORE_LOG_ERROR(
1025 "No velocity field passed to Advect. Not advecting.");
1029 if (advectionTime == 0.) {
1030 advectedTime = advect(std::numeric_limits<double>::max());
1031 numberOfTimeSteps = 1;
1033 double currentTime = 0.0;
1034 numberOfTimeSteps = 0;
1035 while (currentTime < advectionTime) {
1036 currentTime += advect(advectionTime - currentTime);
1037 ++numberOfTimeSteps;
1038 if (performOnlySingleStep)
1041 advectedTime = currentTime;