53template <
class T,
int D>
class Advect {
54 using ConstSparseIterator =
55 viennahrle::ConstSparseIterator<typename Domain<T, D>::DomainType>;
56 using hrleIndexType = viennahrle::IndexType;
61 std::vector<SmartPointer<Domain<T, D>>> levelSets;
62 SmartPointer<VelocityField<T>> velocities =
nullptr;
65 double timeStepRatio = 0.4999;
66 double dissipationAlpha = 1.0;
67 bool calculateNormalVectors =
true;
68 bool ignoreVoids =
false;
69 double advectionTime = 0.;
70 bool performOnlySingleStep =
false;
71 double advectedTime = 0.;
72 unsigned numberOfTimeSteps = 0;
73 bool saveAdvectionVelocities =
false;
74 bool updatePointData =
true;
75 bool checkDissipation =
true;
76 double integrationCutoff = 0.5;
77 bool adaptiveTimeStepping =
false;
78 unsigned adaptiveTimeStepSubdivisions = 20;
79 static constexpr double wrappingLayerEpsilon = 1e-4;
80 SmartPointer<Domain<T, D>> originalLevelSet =
nullptr;
84 std::vector<std::vector<std::pair<std::pair<T, T>,
T>>> storedRates;
85 double currentTimeStep = -1.;
87 VectorType<T, 3> findGlobalAlphas()
const {
89 auto &topDomain = levelSets.back()->getDomain();
90 auto &grid = levelSets.back()->getGrid();
92 const T gridDelta = grid.getGridDelta();
93 const T deltaPos = gridDelta;
94 const T deltaNeg = -gridDelta;
96 VectorType<T, 3> finalAlphas = {0., 0., 0.};
98#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
100 VectorType<T, 3> localAlphas = {0., 0., 0.};
103 p = omp_get_thread_num();
105 viennahrle::Index<D> startVector =
106 (p == 0) ? grid.getMinGridPoint()
107 : topDomain.getSegmentation()[p - 1];
108 viennahrle::Index<D> endVector =
109 (p !=
static_cast<int>(topDomain.getNumberOfSegments() - 1))
110 ? topDomain.getSegmentation()[p]
111 : grid.incrementIndices(grid.getMaxGridPoint());
114 std::vector<ConstSparseIterator> iterators;
115 for (
auto const &
ls : levelSets) {
116 iterators.emplace_back(
ls->getDomain());
120 viennahrle::ConstSparseStarIterator<typename Domain<T, D>::DomainType, 1>
121 neighborIterator(topDomain);
123 for (ConstSparseIterator it(topDomain, startVector);
124 it.getStartIndices() < endVector; ++it) {
126 if (!it.isDefined() || std::abs(it.getValue()) > integrationCutoff)
129 const T value = it.getValue();
130 const auto indices = it.getStartIndices();
134 for (
unsigned lowerLevelSetId = 0; lowerLevelSetId < levelSets.size();
137 iterators[lowerLevelSetId].goToIndicesSequential(indices);
141 if (iterators[lowerLevelSetId].getValue() <=
142 value + wrappingLayerEpsilon) {
145 neighborIterator.goToIndicesSequential(indices);
148 for (
unsigned i = 0; i <
D; ++i) {
149 coords[i] = indices[i] * gridDelta;
152 Vec3D<T> normal = {};
153 T normalModulus = 0.;
154 for (
unsigned i = 0; i <
D; ++i) {
155 const T phiPos = neighborIterator.getNeighbor(i).getValue();
156 const T phiNeg = neighborIterator.getNeighbor(i +
D).getValue();
158 T diffPos = (phiPos - value) / deltaPos;
159 T diffNeg = (phiNeg - value) / deltaNeg;
161 normal[i] = (diffNeg + diffPos) * 0.5;
162 normalModulus += normal[i] * normal[i];
164 normalModulus = 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 void combineLevelSets(
double wTarget,
double wSource) {
201 auto &domainDest = levelSets.back()->getDomain();
202 auto &grid = levelSets.back()->getGrid();
204#pragma omp parallel num_threads(domainDest.getNumberOfSegments())
208 p = omp_get_thread_num();
210 auto &segDest = domainDest.getDomainSegment(p);
212 viennahrle::Index<D> start = (p == 0)
213 ? grid.getMinGridPoint()
214 : domainDest.getSegmentation()[p - 1];
215 viennahrle::Index<D> end =
216 (p !=
static_cast<int>(domainDest.getNumberOfSegments()) - 1)
217 ? domainDest.getSegmentation()[p]
218 : grid.incrementIndices(grid.getMaxGridPoint());
220 ConstSparseIterator itDest(domainDest, start);
221 ConstSparseIterator itTarget(originalLevelSet->getDomain(), start);
223 unsigned definedValueIndex = 0;
224 for (; itDest.getStartIndices() < end; ++itDest) {
225 if (itDest.isDefined()) {
226 itTarget.goToIndicesSequential(itDest.getStartIndices());
227 T valSource = itDest.getValue();
228 T valTarget = itTarget.getValue();
229 segDest.definedValues[definedValueIndex++] =
230 wTarget * valTarget + wSource * valSource;
241 auto &grid = levelSets.back()->getGrid();
242 auto newlsDomain = SmartPointer<Domain<T, D>>::New(grid);
243 auto &newDomain = newlsDomain->getDomain();
244 auto &domain = levelSets.back()->getDomain();
256 newDomain.initialize(domain.getNewSegmentation(),
257 domain.getAllocation() *
258 (2.0 / levelSets.back()->getLevelSetWidth()));
260 const bool updateData = updatePointData;
263 std::vector<std::vector<unsigned>> newDataSourceIds;
265 newDataSourceIds.resize(newDomain.getNumberOfSegments());
267#ifdef DEBUG_LS_ADVECT_HPP
269 auto mesh = SmartPointer<Mesh<T>>::New();
275#pragma omp parallel num_threads(newDomain.getNumberOfSegments())
279 p = omp_get_thread_num();
281 auto &domainSegment = newDomain.getDomainSegment(p);
283 viennahrle::Index<D> startVector =
284 (p == 0) ? grid.getMinGridPoint()
285 : newDomain.getSegmentation()[p - 1];
287 viennahrle::Index<D> endVector =
288 (p !=
static_cast<int>(newDomain.getNumberOfSegments() - 1))
289 ? newDomain.getSegmentation()[p]
290 : grid.incrementIndices(grid.getMaxGridPoint());
295 newDataSourceIds[p].reserve(2.5 * domainSegment.getNumberOfPoints());
298 it(domain, startVector);
299 it.getIndices() < endVector; ++it) {
303 if (std::abs(it.getCenter().getValue()) <= 1.0) {
306 for (; k < 2 *
D; k++)
307 if (std::signbit(it.getNeighbor(k).getValue() - 1e-7) !=
308 std::signbit(it.getCenter().getValue() + 1e-7))
313 if (it.getCenter().getDefinedValue() > 0.5) {
315 for (; j < 2 *
D; j++) {
316 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
317 if (it.getNeighbor(j).getDefinedValue() < -0.5)
321 domainSegment.insertNextDefinedPoint(
322 it.getIndices(), it.getCenter().getDefinedValue());
324 newDataSourceIds[p].push_back(it.getCenter().getPointId());
327 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
329 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
331 }
else if (it.getCenter().getDefinedValue() < -0.5) {
333 for (; j < 2 *
D; j++) {
334 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
335 if (it.getNeighbor(j).getDefinedValue() > 0.5)
340 domainSegment.insertNextDefinedPoint(
341 it.getIndices(), it.getCenter().getDefinedValue());
343 newDataSourceIds[p].push_back(it.getCenter().getPointId());
346 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
348 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
351 domainSegment.insertNextDefinedPoint(
352 it.getIndices(), it.getCenter().getDefinedValue());
354 newDataSourceIds[p].push_back(it.getCenter().getPointId());
357 domainSegment.insertNextUndefinedPoint(
358 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
364 if (it.getCenter().getValue() >= 0) {
365 int usedNeighbor = -1;
367 for (
int i = 0; i < 2 *
D; i++) {
368 T value = it.getNeighbor(i).getValue();
369 if (std::abs(value) <= 1.0 && (value < 0.)) {
370 if (distance > value + 1.0) {
371 distance = value + 1.0;
377 if (distance <= cutoff) {
378 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
380 newDataSourceIds[p].push_back(
381 it.getNeighbor(usedNeighbor).getPointId());
383 domainSegment.insertNextUndefinedPoint(it.getIndices(),
388 int usedNeighbor = -1;
390 for (
int i = 0; i < 2 *
D; i++) {
391 T value = it.getNeighbor(i).getValue();
392 if (std::abs(value) <= 1.0 && (value > 0)) {
393 if (distance < value - 1.0) {
395 distance = value - 1.0;
401 if (distance >= -cutoff) {
402 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
404 newDataSourceIds[p].push_back(
405 it.getNeighbor(usedNeighbor).getPointId());
407 domainSegment.insertNextUndefinedPoint(it.getIndices(),
417 auto &pointData = levelSets.back()->getPointData();
418 newlsDomain->getPointData().translateFromMultiData(pointData,
422 newDomain.finalize();
424 levelSets.back()->deepCopy(newlsDomain);
425 levelSets.back()->finalize(finalWidth);
432 template <
class DiscretizationSchemeType>
433 double integrateTime(DiscretizationSchemeType spatialScheme,
434 double maxTimeStep) {
436 auto &topDomain = levelSets.back()->getDomain();
437 auto &grid = levelSets.back()->getGrid();
442 auto &pointData = levelSets.back()->getPointData();
445 if (voidMarkerPointer ==
nullptr) {
446 VIENNACORE_LOG_WARNING(
"Advect: Cannot find void point markers. Not "
447 "ignoring void points.");
451 const bool ignoreVoidPoints = ignoreVoids;
452 const bool useAdaptiveTimeStepping = adaptiveTimeStepping;
454 if (!storedRates.empty()) {
455 VIENNACORE_LOG_WARNING(
"Advect: Overwriting previously stored rates.");
458 storedRates.resize(topDomain.getNumberOfSegments());
460#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
464 p = omp_get_thread_num();
466 viennahrle::Index<D> startVector =
467 (p == 0) ? grid.getMinGridPoint()
468 : topDomain.getSegmentation()[p - 1];
470 viennahrle::Index<D> endVector =
471 (p !=
static_cast<int>(topDomain.getNumberOfSegments() - 1))
472 ? topDomain.getSegmentation()[p]
473 : grid.incrementIndices(grid.getMaxGridPoint());
475 double tempMaxTimeStep = maxTimeStep;
476 auto &tempRates = storedRates[p];
478 topDomain.getNumberOfPoints() /
479 static_cast<double>((levelSets.back())->getNumberOfSegments()) +
483 std::vector<ConstSparseIterator> iterators;
484 for (
auto const &
ls : levelSets) {
485 iterators.emplace_back(
ls->getDomain());
488 DiscretizationSchemeType scheme(spatialScheme);
490 for (ConstSparseIterator it(topDomain, startVector);
491 it.getStartIndices() < endVector; ++it) {
493 if (!it.isDefined() || std::abs(it.getValue()) > integrationCutoff)
496 T value = it.getValue();
497 double maxStepTime = 0;
498 double cfl = timeStepRatio;
500 for (
int currentLevelSetId = levelSets.size() - 1;
501 currentLevelSetId >= 0; --currentLevelSetId) {
503 std::pair<T, T> gradNDissipation;
505 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
508 for (
unsigned lowerLevelSetId = 0;
509 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
511 iterators[lowerLevelSetId].goToIndicesSequential(
512 it.getStartIndices());
516 if (iterators[lowerLevelSetId].getValue() <=
517 value + wrappingLayerEpsilon) {
519 scheme(it.getStartIndices(), lowerLevelSetId);
525 T velocity = gradNDissipation.first - gradNDissipation.second;
529 maxStepTime += cfl / velocity;
530 tempRates.push_back(std::make_pair(gradNDissipation,
531 -std::numeric_limits<T>::max()));
533 }
else if (velocity == 0.) {
536 maxStepTime = std::numeric_limits<T>::max();
537 tempRates.push_back(std::make_pair(gradNDissipation,
538 std::numeric_limits<T>::max()));
544 if (currentLevelSetId > 0) {
545 iterators[currentLevelSetId - 1].goToIndicesSequential(
546 it.getStartIndices());
547 valueBelow = iterators[currentLevelSetId - 1].getValue();
549 valueBelow = std::numeric_limits<T>::max();
552 T difference = std::abs(valueBelow - value);
554 if (difference >= cfl) {
557 maxStepTime -= cfl / velocity;
558 tempRates.push_back(std::make_pair(
559 gradNDissipation, std::numeric_limits<T>::max()));
564 auto adaptiveFactor = 1.0 / adaptiveTimeStepSubdivisions;
565 if (useAdaptiveTimeStepping && difference > 0.2 * cfl) {
570 maxStepTime -= adaptiveFactor * cfl / velocity;
571 tempRates.push_back(std::make_pair(
572 gradNDissipation, std::numeric_limits<T>::min()));
578 std::make_pair(gradNDissipation, valueBelow));
581 maxStepTime -= difference / velocity;
587 if (maxStepTime < tempMaxTimeStep)
588 tempMaxTimeStep = maxStepTime;
596 scheme.reduceTimeStepHamiltonJacobi(
597 tempMaxTimeStep, levelSets.back()->getGrid().getGridDelta());
600 if (tempMaxTimeStep < maxTimeStep)
601 maxTimeStep = tempMaxTimeStep;
612 void computeRates(
double maxTimeStep = std::numeric_limits<double>::max()) {
616 calculateNormalVectors);
617 currentTimeStep = integrateTime(is, maxTimeStep);
620 calculateNormalVectors);
621 currentTimeStep = integrateTime(is, maxTimeStep);
623 auto alphas = findGlobalAlphas();
625 dissipationAlpha, alphas,
626 calculateNormalVectors);
627 currentTimeStep = integrateTime(is, maxTimeStep);
629 auto alphas = findGlobalAlphas();
631 dissipationAlpha, alphas,
632 calculateNormalVectors);
633 currentTimeStep = integrateTime(is, maxTimeStep);
634 }
else if (spatialScheme ==
637 levelSets.back(), velocities);
638 currentTimeStep = integrateTime(is, maxTimeStep);
639 }
else if (spatialScheme ==
642 levelSets.back(), velocities, dissipationAlpha);
643 currentTimeStep = integrateTime(is, maxTimeStep);
644 }
else if (spatialScheme ==
647 levelSets.back(), velocities, dissipationAlpha);
648 currentTimeStep = integrateTime(is, maxTimeStep);
649 }
else if (spatialScheme ==
652 levelSets.back(), velocities, dissipationAlpha);
653 currentTimeStep = integrateTime(is, maxTimeStep);
654 }
else if (spatialScheme ==
657 levelSets.back(), velocities, dissipationAlpha);
658 currentTimeStep = integrateTime(is, maxTimeStep);
659 }
else if (spatialScheme ==
662 levelSets.back(), velocities, dissipationAlpha);
663 currentTimeStep = integrateTime(is, maxTimeStep);
668 currentTimeStep = integrateTime(is, maxTimeStep);
670 VIENNACORE_LOG_ERROR(
"Advect: Discretization scheme not found.");
671 currentTimeStep = -1.;
677 void updateLevelSet(
double dt) {
678 if (timeStepRatio >= 0.5) {
679 VIENNACORE_LOG_WARNING(
680 "Integration time step ratio should be smaller than 0.5. "
681 "Advection might fail!");
684 auto &topDomain = levelSets.back()->getDomain();
686 assert(dt >= 0. &&
"No time step set!");
687 assert(storedRates.size() == topDomain.getNumberOfSegments());
693 const bool saveVelocities = saveAdvectionVelocities;
694 std::vector<std::vector<double>> dissipationVectors(
695 levelSets.back()->getNumberOfSegments());
696 std::vector<std::vector<double>> velocityVectors(
697 levelSets.back()->getNumberOfSegments());
699 const bool checkDiss = checkDissipation;
701#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
705 p = omp_get_thread_num();
707 auto itRS = storedRates[p].cbegin();
708 auto &segment = topDomain.getDomainSegment(p);
709 const unsigned maxId = segment.getNumberOfPoints();
711 if (saveVelocities) {
712 velocityVectors[p].resize(maxId);
713 dissipationVectors[p].resize(maxId);
716 for (
unsigned localId = 0; localId < maxId; ++localId) {
717 T &value = segment.definedValues[localId];
720 if (std::abs(value) > integrationCutoff)
728 auto const [gradient, dissipation] = itRS->first;
729 T velocity = gradient - dissipation;
731 if (checkDiss && (gradient < 0 && velocity > 0) ||
732 (gradient > 0 && velocity < 0)) {
736 T rate = time * velocity;
737 while (std::abs(itRS->second - value) < std::abs(rate)) {
738 time -= std::abs((itRS->second - value) / velocity);
739 value = itRS->second;
743 velocity = itRS->first.first - itRS->first.second;
744 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
745 (itRS->first.first > 0 && velocity < 0)) {
748 rate = time * velocity;
754 if (saveVelocities) {
755 velocityVectors[p][localId] = rate;
756 dissipationVectors[p][localId] = itRS->first.second;
762 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
770 if (saveVelocities) {
771 auto &pointData = levelSets.back()->getPointData();
776 for (
unsigned i = 0; i < velocityVectors.size(); ++i) {
777 vels.insert(vels.end(),
778 std::make_move_iterator(velocityVectors[i].begin()),
779 std::make_move_iterator(velocityVectors[i].end()));
780 diss.insert(diss.end(),
781 std::make_move_iterator(dissipationVectors[i].begin()),
782 std::make_move_iterator(dissipationVectors[i].end()));
784 pointData.insertReplaceScalarData(std::move(vels),
velocityLabel);
792 void adjustLowerLayers() {
796 for (
unsigned i = 0; i < levelSets.size() - 1; ++i) {
798 levelSets[i], levelSets.back(),
807 double advect(
double maxTimeStep) {
808 switch (temporalScheme) {
829 levelSets.push_back(passedlsDomain);
834 levelSets.push_back(passedlsDomain);
835 velocities = passedVelocities;
840 : levelSets(passedlsDomains) {
841 velocities = passedVelocities;
847 levelSets.push_back(passedlsDomain);
855 velocities = passedVelocities;
895 adaptiveTimeStepping = aTS;
896 if (subdivisions < 1) {
897 VIENNACORE_LOG_WARNING(
"Advect: Adaptive time stepping subdivisions must "
898 "be at least 1. Setting to 1.");
901 adaptiveTimeStepSubdivisions = subdivisions;
930 [[deprecated(
"Use setSpatialScheme instead")]]
void
932 VIENNACORE_LOG_WARNING(
933 "Advect::setIntegrationScheme is deprecated and will be removed in "
934 "future versions. Use setSpatialScheme instead.");
935 spatialScheme = scheme;
958 if (levelSets.empty()) {
959 VIENNACORE_LOG_ERROR(
"No level sets passed to Advect.");
971 }
else if (spatialScheme ==
975 }
else if (spatialScheme ==
978 }
else if (spatialScheme ==
981 }
else if (spatialScheme ==
984 }
else if (spatialScheme ==
987 }
else if (spatialScheme ==
995 VIENNACORE_LOG_ERROR(
"Advect: Discretization scheme not found.");
1001 if (levelSets.empty()) {
1002 VIENNACORE_LOG_ERROR(
"No level sets passed to Advect. Not advecting.");
1005 if (velocities ==
nullptr) {
1006 VIENNACORE_LOG_ERROR(
1007 "No velocity field passed to Advect. Not advecting.");
1011 if (advectionTime == 0.) {
1012 advectedTime = advect(std::numeric_limits<double>::max());
1013 numberOfTimeSteps = 1;
1015 double currentTime = 0.0;
1016 numberOfTimeSteps = 0;
1017 while (currentTime < advectionTime) {
1018 currentTime += advect(advectionTime - currentTime);
1019 ++numberOfTimeSteps;
1020 if (performOnlySingleStep)
1023 advectedTime = currentTime;