66template <
class T,
int D>
class Advect {
67 using ConstSparseIterator =
68 viennahrle::ConstSparseIterator<typename Domain<T, D>::DomainType>;
69 using hrleIndexType = viennahrle::IndexType;
71 std::vector<SmartPointer<Domain<T, D>>> levelSets;
72 SmartPointer<VelocityField<T>> velocities =
nullptr;
75 double timeStepRatio = 0.4999;
76 double dissipationAlpha = 1.0;
77 bool calculateNormalVectors =
true;
78 bool ignoreVoids =
false;
79 double advectionTime = 0.;
80 bool performOnlySingleStep =
false;
81 double advectedTime = 0.;
82 unsigned numberOfTimeSteps = 0;
83 bool saveAdvectionVelocities =
false;
84 bool updatePointData =
true;
85 bool checkDissipation =
true;
86 bool adaptiveTimeStepping =
false;
87 double adaptiveTimeStepThreshold = 0.05;
88 static constexpr double wrappingLayerEpsilon = 1e-4;
92 std::vector<std::vector<std::pair<std::pair<T, T>,
T>>> storedRates;
93 double currentTimeStep = -1.;
95 VectorType<T, 3> findGlobalAlphas()
const {
97 auto &topDomain = levelSets.back()->getDomain();
98 auto &grid = levelSets.back()->getGrid();
100 const T gridDelta = grid.getGridDelta();
101 const T deltaPos = gridDelta;
102 const T deltaNeg = -gridDelta;
104 VectorType<T, 3> finalAlphas = {0., 0., 0.};
106#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
108 VectorType<T, 3> localAlphas = {0., 0., 0.};
111 p = omp_get_thread_num();
113 viennahrle::Index<D> startVector =
114 (p == 0) ? grid.getMinGridPoint()
115 : topDomain.getSegmentation()[p - 1];
116 viennahrle::Index<D> endVector =
117 (p !=
static_cast<int>(topDomain.getNumberOfSegments() - 1))
118 ? topDomain.getSegmentation()[p]
119 : grid.incrementIndices(grid.getMaxGridPoint());
122 std::vector<ConstSparseIterator> iterators;
123 for (
auto const &
ls : levelSets) {
124 iterators.emplace_back(
ls->getDomain());
128 viennahrle::ConstSparseStarIterator<typename Domain<T, D>::DomainType, 1>
129 neighborIterator(topDomain);
131 for (ConstSparseIterator it(topDomain, startVector);
132 it.getStartIndices() < endVector; ++it) {
134 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
137 const T value = it.getValue();
138 const auto indices = it.getStartIndices();
142 for (
unsigned lowerLevelSetId = 0; lowerLevelSetId < levelSets.size();
145 iterators[lowerLevelSetId].goToIndicesSequential(indices);
149 if (iterators[lowerLevelSetId].getValue() <=
150 value + wrappingLayerEpsilon) {
153 neighborIterator.goToIndicesSequential(indices);
156 for (
unsigned i = 0; i <
D; ++i) {
157 coords[i] = indices[i] * gridDelta;
160 Vec3D<T> normal = {};
161 T normalModulus = 0.;
162 for (
unsigned i = 0; i <
D; ++i) {
163 const T phiPos = neighborIterator.getNeighbor(i).getValue();
164 const T phiNeg = neighborIterator.getNeighbor(i +
D).getValue();
166 T diffPos = (phiPos - value) / deltaPos;
167 T diffNeg = (phiNeg - value) / deltaNeg;
169 normal[i] = (diffNeg + diffPos) * 0.5;
170 normalModulus += normal[i] * normal[i];
172 normalModulus = std::sqrt(normalModulus);
173 for (
unsigned i = 0; i <
D; ++i)
174 normal[i] /= normalModulus;
176 T scaVel = velocities->getScalarVelocity(
177 coords, lowerLevelSetId, normal,
178 neighborIterator.getCenter().getPointId());
179 auto vecVel = velocities->getVectorVelocity(
180 coords, lowerLevelSetId, normal,
181 neighborIterator.getCenter().getPointId());
183 for (
unsigned i = 0; i <
D; ++i) {
184 T tempAlpha = std::abs((scaVel + vecVel[i]) * normal[i]);
185 localAlphas[i] = std::max(localAlphas[i], tempAlpha);
196 for (
unsigned i = 0; i <
D; ++i) {
197 finalAlphas[i] = std::max(finalAlphas[i], localAlphas[i]);
210 auto &grid = levelSets.back()->getGrid();
211 auto newlsDomain = SmartPointer<Domain<T, D>>::New(grid);
212 auto &newDomain = newlsDomain->getDomain();
213 auto &domain = levelSets.back()->getDomain();
215 newDomain.initialize(domain.getNewSegmentation(),
216 domain.getAllocation() *
217 (2.0 / levelSets.back()->getLevelSetWidth()));
219 const bool updateData = updatePointData;
222 std::vector<std::vector<unsigned>> newDataSourceIds;
224 newDataSourceIds.resize(newDomain.getNumberOfSegments());
226#ifdef DEBUG_LS_ADVECT_HPP
228 auto mesh = SmartPointer<Mesh<T>>::New();
234#pragma omp parallel num_threads(newDomain.getNumberOfSegments())
238 p = omp_get_thread_num();
240 auto &domainSegment = newDomain.getDomainSegment(p);
242 viennahrle::Index<D> startVector =
243 (p == 0) ? grid.getMinGridPoint()
244 : newDomain.getSegmentation()[p - 1];
246 viennahrle::Index<D> endVector =
247 (p !=
static_cast<int>(newDomain.getNumberOfSegments() - 1))
248 ? newDomain.getSegmentation()[p]
249 : grid.incrementIndices(grid.getMaxGridPoint());
254 newDataSourceIds[p].reserve(2.5 * domainSegment.getNumberOfPoints());
257 it(domain, startVector);
258 it.getIndices() < endVector; ++it) {
262 if (std::abs(it.getCenter().getValue()) <= 1.0) {
265 for (; k < 2 *
D; k++)
266 if (std::signbit(it.getNeighbor(k).getValue() - 1e-7) !=
267 std::signbit(it.getCenter().getValue() + 1e-7))
272 if (it.getCenter().getDefinedValue() > 0.5) {
274 for (; j < 2 *
D; j++) {
275 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
276 if (it.getNeighbor(j).getDefinedValue() < -0.5)
280 domainSegment.insertNextDefinedPoint(
281 it.getIndices(), it.getCenter().getDefinedValue());
283 newDataSourceIds[p].push_back(it.getCenter().getPointId());
286 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
288 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
290 }
else if (it.getCenter().getDefinedValue() < -0.5) {
292 for (; j < 2 *
D; j++) {
293 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
294 if (it.getNeighbor(j).getDefinedValue() > 0.5)
299 domainSegment.insertNextDefinedPoint(
300 it.getIndices(), it.getCenter().getDefinedValue());
302 newDataSourceIds[p].push_back(it.getCenter().getPointId());
305 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
307 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
310 domainSegment.insertNextDefinedPoint(
311 it.getIndices(), it.getCenter().getDefinedValue());
313 newDataSourceIds[p].push_back(it.getCenter().getPointId());
316 domainSegment.insertNextUndefinedPoint(
317 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
323 if (it.getCenter().getValue() >= 0) {
324 int usedNeighbor = -1;
326 for (
int i = 0; i < 2 *
D; i++) {
327 T value = it.getNeighbor(i).getValue();
328 if (std::abs(value) <= 1.0 && (value < 0.)) {
329 if (distance > value + 1.0) {
330 distance = value + 1.0;
336 if (distance <= 1.) {
337 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
339 newDataSourceIds[p].push_back(
340 it.getNeighbor(usedNeighbor).getPointId());
342 domainSegment.insertNextUndefinedPoint(it.getIndices(),
347 int usedNeighbor = -1;
349 for (
int i = 0; i < 2 *
D; i++) {
350 T value = it.getNeighbor(i).getValue();
351 if (std::abs(value) <= 1.0 && (value > 0)) {
352 if (distance < value - 1.0) {
354 distance = value - 1.0;
360 if (distance >= -1.) {
361 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
363 newDataSourceIds[p].push_back(
364 it.getNeighbor(usedNeighbor).getPointId());
366 domainSegment.insertNextUndefinedPoint(it.getIndices(),
376 auto &pointData = levelSets.back()->getPointData();
377 newlsDomain->getPointData().translateFromMultiData(pointData,
381 newDomain.finalize();
383 levelSets.back()->deepCopy(newlsDomain);
384 levelSets.back()->finalize(2);
389 double advect(
double maxTimeStep) {
391 if (levelSets.empty()) {
392 VIENNACORE_LOG_ERROR(
"No level sets passed to Advect. Not advecting.");
393 return std::numeric_limits<double>::max();
395 if (velocities ==
nullptr) {
396 VIENNACORE_LOG_ERROR(
397 "No velocity field passed to Advect. Not advecting.");
398 return std::numeric_limits<double>::max();
401 if (currentTimeStep < 0. || storedRates.empty())
413 if (integrationScheme !=
415 for (
unsigned i = 0; i < levelSets.size() - 1; ++i) {
428 template <
class IntegrationSchemeType>
429 double integrateTime(IntegrationSchemeType IntegrationScheme,
430 double maxTimeStep) {
432 auto &topDomain = levelSets.back()->getDomain();
433 auto &grid = levelSets.back()->getGrid();
438 auto &pointData = levelSets.back()->getPointData();
441 if (voidMarkerPointer ==
nullptr) {
442 VIENNACORE_LOG_WARNING(
"Advect: Cannot find void point markers. Not "
443 "ignoring void points.");
447 const bool ignoreVoidPoints = ignoreVoids;
448 const bool useAdaptiveTimeStepping = adaptiveTimeStepping;
449 const double atsThreshold = adaptiveTimeStepThreshold;
451 if (!storedRates.empty()) {
452 VIENNACORE_LOG_WARNING(
"Advect: Overwriting previously stored rates.");
455 storedRates.resize(topDomain.getNumberOfSegments());
457#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
461 p = omp_get_thread_num();
463 viennahrle::Index<D> startVector =
464 (p == 0) ? grid.getMinGridPoint()
465 : topDomain.getSegmentation()[p - 1];
467 viennahrle::Index<D> endVector =
468 (p !=
static_cast<int>(topDomain.getNumberOfSegments() - 1))
469 ? topDomain.getSegmentation()[p]
470 : grid.incrementIndices(grid.getMaxGridPoint());
472 double tempMaxTimeStep = maxTimeStep;
473 auto &tempRates = storedRates[p];
475 topDomain.getNumberOfPoints() /
476 static_cast<double>((levelSets.back())->getNumberOfSegments()) +
480 std::vector<ConstSparseIterator> iterators;
481 for (
auto const &
ls : levelSets) {
482 iterators.emplace_back(
ls->getDomain());
485 IntegrationSchemeType scheme(IntegrationScheme);
487 for (ConstSparseIterator it(topDomain, startVector);
488 it.getStartIndices() < endVector; ++it) {
490 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
493 T value = it.getValue();
494 double maxStepTime = 0;
495 double cfl = timeStepRatio;
497 for (
int currentLevelSetId = levelSets.size() - 1;
498 currentLevelSetId >= 0; --currentLevelSetId) {
500 std::pair<T, T> gradNDissipation;
502 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
505 for (
unsigned lowerLevelSetId = 0;
506 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
508 iterators[lowerLevelSetId].goToIndicesSequential(
509 it.getStartIndices());
513 if (iterators[lowerLevelSetId].getValue() <=
514 value + wrappingLayerEpsilon) {
516 scheme(it.getStartIndices(), lowerLevelSetId);
522 T velocity = gradNDissipation.first - gradNDissipation.second;
526 maxStepTime += cfl / velocity;
527 tempRates.push_back(std::make_pair(gradNDissipation,
528 -std::numeric_limits<T>::max()));
530 }
else if (velocity == 0.) {
533 maxStepTime = std::numeric_limits<T>::max();
534 tempRates.push_back(std::make_pair(gradNDissipation,
535 std::numeric_limits<T>::max()));
541 if (currentLevelSetId > 0) {
542 iterators[currentLevelSetId - 1].goToIndicesSequential(
543 it.getStartIndices());
544 valueBelow = iterators[currentLevelSetId - 1].getValue();
546 valueBelow = std::numeric_limits<T>::max();
549 T difference = std::abs(valueBelow - value);
551 if (difference >= cfl) {
554 maxStepTime -= cfl / velocity;
555 tempRates.push_back(std::make_pair(
556 gradNDissipation, std::numeric_limits<T>::max()));
561 if (useAdaptiveTimeStepping && difference > atsThreshold * cfl) {
565 maxStepTime -= atsThreshold * cfl / velocity;
566 tempRates.push_back(std::make_pair(
567 gradNDissipation, std::numeric_limits<T>::min()));
573 std::make_pair(gradNDissipation, valueBelow));
576 maxStepTime -= difference / velocity;
582 if (maxStepTime < tempMaxTimeStep)
583 tempMaxTimeStep = maxStepTime;
591 scheme.reduceTimeStepHamiltonJacobi(
592 tempMaxTimeStep, levelSets.back()->getGrid().getGridDelta());
595 if (tempMaxTimeStep < maxTimeStep)
596 maxTimeStep = tempMaxTimeStep;
608 if (timeStepRatio >= 0.5) {
609 VIENNACORE_LOG_WARNING(
610 "Integration time step ratio should be smaller than 0.5. "
611 "Advection might fail!");
614 auto &topDomain = levelSets.back()->getDomain();
615 auto &grid = levelSets.back()->getGrid();
617 assert(currentTimeStep >= 0. &&
"No time step set!");
618 assert(storedRates.size() == topDomain.getNumberOfSegments());
624 const bool saveVelocities = saveAdvectionVelocities;
625 std::vector<std::vector<double>> dissipationVectors(
626 levelSets.back()->getNumberOfSegments());
627 std::vector<std::vector<double>> velocityVectors(
628 levelSets.back()->getNumberOfSegments());
630 const bool checkDiss = checkDissipation;
632#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
636 p = omp_get_thread_num();
638 auto itRS = storedRates[p].cbegin();
639 auto &segment = topDomain.getDomainSegment(p);
640 const unsigned maxId = segment.getNumberOfPoints();
642 if (saveVelocities) {
643 velocityVectors[p].resize(maxId);
644 dissipationVectors[p].resize(maxId);
647 for (
unsigned localId = 0; localId < maxId; ++localId) {
648 T &value = segment.definedValues[localId];
649 double time = currentTimeStep;
654 auto const [gradient, dissipation] = itRS->first;
655 T velocity = gradient - dissipation;
657 if (checkDiss && (gradient < 0 && velocity > 0) ||
658 (gradient > 0 && velocity < 0)) {
662 T rate = time * velocity;
663 while (std::abs(itRS->second - value) < std::abs(rate)) {
664 time -= std::abs((itRS->second - value) / velocity);
665 value = itRS->second;
669 velocity = itRS->first.first - itRS->first.second;
670 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
671 (itRS->first.first > 0 && velocity < 0)) {
674 rate = time * velocity;
680 if (saveVelocities) {
681 velocityVectors[p][localId] = rate;
682 dissipationVectors[p][localId] = itRS->first.second;
688 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
696 if (saveVelocities) {
697 auto &pointData = levelSets.back()->getPointData();
702 for (
unsigned i = 0; i < velocityVectors.size(); ++i) {
703 vels.insert(vels.end(),
704 std::make_move_iterator(velocityVectors[i].begin()),
705 std::make_move_iterator(velocityVectors[i].end()));
706 diss.insert(diss.end(),
707 std::make_move_iterator(dissipationVectors[i].begin()),
708 std::make_move_iterator(dissipationVectors[i].end()));
710 pointData.insertReplaceScalarData(std::move(vels),
velocityLabel);
725 levelSets.push_back(passedlsDomain);
730 levelSets.push_back(passedlsDomain);
731 velocities = passedVelocities;
736 : levelSets(passedlsDomains) {
737 velocities = passedVelocities;
743 levelSets.push_back(passedlsDomain);
751 velocities = passedVelocities;
796 adaptiveTimeStepThreshold = threshold;
822 integrationScheme = scheme;
842 if (levelSets.empty()) {
843 VIENNACORE_LOG_ERROR(
"No level sets passed to Advect.");
849 }
else if (integrationScheme ==
852 }
else if (integrationScheme ==
855 }
else if (integrationScheme ==
858 }
else if (integrationScheme ==
863 }
else if (integrationScheme ==
866 }
else if (integrationScheme ==
869 }
else if (integrationScheme ==
872 }
else if (integrationScheme ==
875 }
else if (integrationScheme ==
883 VIENNACORE_LOG_ERROR(
"Advect: Integration scheme not found.");
889 if (advectionTime == 0.) {
890 advectedTime = advect(std::numeric_limits<double>::max());
891 numberOfTimeSteps = 1;
893 double currentTime = 0.0;
894 numberOfTimeSteps = 0;
895 while (currentTime < advectionTime) {
896 currentTime += advect(advectionTime - currentTime);
898 if (performOnlySingleStep)
901 advectedTime = currentTime;
912 calculateNormalVectors);
913 currentTimeStep = integrateTime(is, maxTimeStep);
914 }
else if (integrationScheme ==
917 calculateNormalVectors);
918 currentTimeStep = integrateTime(is, maxTimeStep);
919 }
else if (integrationScheme ==
921 auto alphas = findGlobalAlphas();
923 dissipationAlpha, alphas,
924 calculateNormalVectors);
925 currentTimeStep = integrateTime(is, maxTimeStep);
926 }
else if (integrationScheme ==
928 auto alphas = findGlobalAlphas();
930 dissipationAlpha, alphas,
931 calculateNormalVectors);
932 currentTimeStep = integrateTime(is, maxTimeStep);
933 }
else if (integrationScheme ==
937 levelSets.back(), velocities);
938 currentTimeStep = integrateTime(is, maxTimeStep);
939 }
else if (integrationScheme ==
942 levelSets.back(), velocities, dissipationAlpha);
943 currentTimeStep = integrateTime(is, maxTimeStep);
944 }
else if (integrationScheme ==
947 levelSets.back(), velocities, dissipationAlpha);
948 currentTimeStep = integrateTime(is, maxTimeStep);
949 }
else if (integrationScheme ==
952 levelSets.back(), velocities, dissipationAlpha);
953 currentTimeStep = integrateTime(is, maxTimeStep);
954 }
else if (integrationScheme ==
957 levelSets.back(), velocities, dissipationAlpha);
958 currentTimeStep = integrateTime(is, maxTimeStep);
959 }
else if (integrationScheme ==
962 levelSets.back(), velocities, dissipationAlpha);
963 currentTimeStep = integrateTime(is, maxTimeStep);
967 calculateNormalVectors);
968 currentTimeStep = integrateTime(is, maxTimeStep);
970 VIENNACORE_LOG_ERROR(
"Advect: Integration scheme not found.");
971 currentTimeStep = -1.;