65template <
class T,
int D>
class Advect {
66 std::vector<SmartPointer<Domain<T, D>>> levelSets;
67 SmartPointer<VelocityField<T>> velocities =
nullptr;
70 double timeStepRatio = 0.4999;
71 double dissipationAlpha = 1.0;
72 bool calculateNormalVectors =
true;
73 bool ignoreVoids =
false;
74 double advectionTime = 0.;
75 bool performOnlySingleStep =
false;
76 double advectedTime = 0.;
77 unsigned numberOfTimeSteps = 0;
78 bool saveAdvectionVelocities =
false;
79 bool updatePointData =
true;
80 bool checkDissipation =
true;
81 static constexpr double wrappingLayerEpsilon = 1e-4;
83 hrleVectorType<T, 3> findGlobalAlphas()
const {
85 auto &topDomain = levelSets.back()->getDomain();
86 auto &grid = levelSets.back()->getGrid();
88 const T gridDelta = grid.getGridDelta();
89 const T deltaPos = gridDelta;
90 const T deltaNeg = -gridDelta;
92 hrleVectorType<T, 3> finalAlphas(0., 0., 0.);
94#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
98 p = omp_get_thread_num();
101 hrleVectorType<T, 3> localAlphas(0., 0., 0.);
103 hrleVectorType<hrleIndexType, D> startVector =
104 (p == 0) ? grid.getMinGridPoint()
105 : topDomain.getSegmentation()[p - 1];
107 hrleVectorType<hrleIndexType, D> endVector =
108 (p !=
static_cast<int>(topDomain.getNumberOfSegments() - 1))
109 ? topDomain.getSegmentation()[p]
110 : grid.incrementIndices(grid.getMaxGridPoint());
113 std::vector<hrleConstSparseIterator<typename Domain<T, D>::DomainType>>
115 for (
auto it = levelSets.begin(); it != levelSets.end(); ++it) {
118 (*it)->getDomain()));
122 hrleConstSparseStarIterator<hrleDomain<T, D>, 1> neighborIterator(
126 topDomain, startVector);
127 it.getStartIndices() < endVector; ++it) {
129 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
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;
155 Vec3D<T> normal = {};
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 T diffPos = (phiPos - value) / deltaPos;
162 T diffNeg = (phiNeg - value) / deltaNeg;
164 normal[i] = (diffNeg + diffPos) * 0.5;
165 normalModulus += normal[i] * normal[i];
167 normalModulus = std::sqrt(normalModulus);
168 for (
unsigned i = 0; i <
D; ++i)
169 normal[i] /= normalModulus;
171 T scaVel = velocities->getScalarVelocity(
172 coords, lowerLevelSetId, normal,
173 neighborIterator.getCenter().getPointId());
174 auto vecVel = velocities->getVectorVelocity(
175 coords, lowerLevelSetId, normal,
176 neighborIterator.getCenter().getPointId());
178 for (
unsigned i = 0; i <
D; ++i) {
179 T tempAlpha = std::abs((scaVel + vecVel[i]) * normal[i]);
180 localAlphas[i] = std::max(localAlphas[i], tempAlpha);
190 for (
unsigned i = 0; i <
D; ++i) {
191 finalAlphas[i] = std::max(finalAlphas[i], localAlphas[i]);
204 auto &grid = levelSets.back()->getGrid();
205 auto newlsDomain = SmartPointer<Domain<T, D>>::New(grid);
209 newDomain.initialize(domain.getNewSegmentation(),
210 domain.getAllocation() *
211 (2.0 / levelSets.back()->getLevelSetWidth()));
213 const bool updateData = updatePointData;
216 std::vector<std::vector<unsigned>> newDataSourceIds;
220#ifdef DEBUG_LS_ADVECT_HPP
222 auto mesh = SmartPointer<Mesh<T>>::New();
228#pragma omp parallel num_threads(newDomain.getNumberOfSegments())
232 p = omp_get_thread_num();
235 auto &domainSegment = newDomain.getDomainSegment(p);
237 hrleVectorType<hrleIndexType, D> startVector =
238 (p == 0) ? grid.getMinGridPoint()
239 : newDomain.getSegmentation()[p - 1];
241 hrleVectorType<hrleIndexType, D> endVector =
243 ? newDomain.getSegmentation()[p]
244 : grid.incrementIndices(grid.getMaxGridPoint());
249 newDataSourceIds[p].reserve(2.5 * domainSegment.getNumberOfPoints());
252 domain, startVector);
253 it.getIndices() < endVector; ++it) {
257 if (std::abs(it.getCenter().getValue()) <= 1.0) {
260 for (; k < 2 *
D; k++)
261 if (std::signbit(it.getNeighbor(k).getValue() - 1e-7) !=
262 std::signbit(it.getCenter().getValue() + 1e-7))
267 if (it.getCenter().getDefinedValue() > 0.5) {
269 for (; j < 2 *
D; j++) {
270 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
271 if (it.getNeighbor(j).getDefinedValue() < -0.5)
275 domainSegment.insertNextDefinedPoint(
276 it.getIndices(), it.getCenter().getDefinedValue());
278 newDataSourceIds[p].push_back(it.getCenter().getPointId());
281 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
283 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
285 }
else if (it.getCenter().getDefinedValue() < -0.5) {
287 for (; j < 2 *
D; j++) {
288 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
289 if (it.getNeighbor(j).getDefinedValue() > 0.5)
294 domainSegment.insertNextDefinedPoint(
295 it.getIndices(), it.getCenter().getDefinedValue());
297 newDataSourceIds[p].push_back(it.getCenter().getPointId());
300 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
302 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
305 domainSegment.insertNextDefinedPoint(
306 it.getIndices(), it.getCenter().getDefinedValue());
308 newDataSourceIds[p].push_back(it.getCenter().getPointId());
311 domainSegment.insertNextUndefinedPoint(
312 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
318 if (it.getCenter().getValue() >= 0) {
319 int usedNeighbor = -1;
321 for (
int i = 0; i < 2 *
D; i++) {
322 T value = it.getNeighbor(i).getValue();
323 if (std::abs(value) <= 1.0 && (value < 0.)) {
324 if (distance > value + 1.0) {
325 distance = value + 1.0;
331 if (distance <= 1.) {
332 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
334 newDataSourceIds[p].push_back(
335 it.getNeighbor(usedNeighbor).getPointId());
337 domainSegment.insertNextUndefinedPoint(it.getIndices(),
342 int usedNeighbor = -1;
344 for (
int i = 0; i < 2 *
D; i++) {
345 T value = it.getNeighbor(i).getValue();
346 if (std::abs(value) <= 1.0 && (value > 0)) {
347 if (distance < value - 1.0) {
349 distance = value - 1.0;
355 if (distance >= -1.) {
356 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
358 newDataSourceIds[p].push_back(
359 it.getNeighbor(usedNeighbor).getPointId());
361 domainSegment.insertNextUndefinedPoint(it.getIndices(),
371 auto &pointData = levelSets.back()->getPointData();
372 newlsDomain->getPointData().translateFromMultiData(pointData,
378 levelSets.back()->
deepCopy(newlsDomain);
379 levelSets.back()->finalize(2);
384 double advect(
double maxTimeStep = std::numeric_limits<double>::max()) {
386 if (levelSets.size() < 1) {
387 Logger::getInstance()
388 .addWarning(
"No level sets passed to Advect. Not advecting.")
390 return std::numeric_limits<double>::max();
392 if (velocities ==
nullptr) {
393 Logger::getInstance()
394 .addWarning(
"No velocity field passed to Advect. Not advecting.")
396 return std::numeric_limits<double>::max();
401 double currentTime = 0.;
404 calculateNormalVectors);
405 currentTime = integrateTime(is, maxTimeStep);
406 }
else if (integrationScheme ==
409 calculateNormalVectors);
410 currentTime = integrateTime(is, maxTimeStep);
411 }
else if (integrationScheme ==
413 auto alphas = findGlobalAlphas();
415 dissipationAlpha, alphas,
416 calculateNormalVectors);
417 currentTime = integrateTime(is, maxTimeStep);
418 }
else if (integrationScheme ==
420 auto alphas = findGlobalAlphas();
422 dissipationAlpha, alphas,
423 calculateNormalVectors);
424 currentTime = integrateTime(is, maxTimeStep);
425 }
else if (integrationScheme ==
429 levelSets.back(), velocities);
430 currentTime = integrateTime(is, maxTimeStep);
431 }
else if (integrationScheme ==
434 levelSets.back(), velocities, dissipationAlpha);
435 currentTime = integrateTime(is, maxTimeStep);
436 }
else if (integrationScheme ==
439 levelSets.back(), velocities, dissipationAlpha);
440 currentTime = integrateTime(is, maxTimeStep);
441 }
else if (integrationScheme ==
444 levelSets.back(), velocities, dissipationAlpha);
445 currentTime = integrateTime(is, maxTimeStep);
446 }
else if (integrationScheme ==
449 levelSets.back(), velocities, dissipationAlpha);
450 currentTime = integrateTime(is, maxTimeStep);
451 }
else if (integrationScheme ==
454 levelSets.back(), velocities, dissipationAlpha);
455 currentTime = integrateTime(is, maxTimeStep);
457 Logger::getInstance()
458 .addWarning(
"Advect: Integration scheme not found. Not advecting.")
463 return std::numeric_limits<double>::max();
473 if (integrationScheme !=
475 for (
unsigned i = 0; i < levelSets.size() - 1; ++i) {
489 template <
class IntegrationSchemeType>
491 integrateTime(IntegrationSchemeType IntegrationScheme,
492 double maxTimeStep = std::numeric_limits<double>::max()) {
493 if (timeStepRatio >= 0.5) {
494 Logger::getInstance()
495 .addWarning(
"Integration time step ratio should be smaller than 0.5. "
496 "Advection might fail!")
500 auto &topDomain = levelSets.back()->getDomain();
501 auto &grid = levelSets.back()->getGrid();
503 std::vector<std::vector<std::pair<std::pair<T, T>,
T>>> totalTempRates;
504 totalTempRates.resize((levelSets.back())->getNumberOfSegments());
509 auto &pointData = levelSets.back()->getPointData();
512 if (voidMarkerPointer ==
nullptr) {
513 Logger::getInstance()
514 .addWarning(
"Advect: Cannot find void point markers. Not "
515 "ignoring void points.")
521 const bool ignoreVoidPoints = ignoreVoids;
523#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
527 p = omp_get_thread_num();
530 hrleVectorType<hrleIndexType, D> startVector =
531 (p == 0) ? grid.getMinGridPoint()
532 : topDomain.getSegmentation()[p - 1];
534 hrleVectorType<hrleIndexType, D> endVector =
535 (p !=
static_cast<int>(topDomain.getNumberOfSegments() - 1))
536 ? topDomain.getSegmentation()[p]
537 : grid.incrementIndices(grid.getMaxGridPoint());
539 double tempMaxTimeStep = maxTimeStep;
540 auto &tempRates = totalTempRates[p];
541 tempRates.reserve(topDomain.getNumberOfPoints() /
542 double((levelSets.back())->getNumberOfSegments()) +
546 std::vector<hrleSparseIterator<typename Domain<T, D>::DomainType>>
548 for (
auto it = levelSets.begin(); it != levelSets.end(); ++it) {
551 (*it)->getDomain()));
554 IntegrationSchemeType scheme(IntegrationScheme);
557 topDomain, startVector);
558 it.getStartIndices() < endVector; ++it) {
560 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
563 T value = it.getValue();
564 double maxStepTime = 0;
565 double cfl = timeStepRatio;
567 for (
int currentLevelSetId = levelSets.size() - 1;
568 currentLevelSetId >= 0; --currentLevelSetId) {
570 std::pair<T, T> gradNDissipation;
572 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
575 for (
unsigned lowerLevelSetId = 0;
576 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
578 iterators[lowerLevelSetId].goToIndicesSequential(
579 it.getStartIndices());
583 if (iterators[lowerLevelSetId].getValue() <=
584 value + wrappingLayerEpsilon) {
586 scheme(it.getStartIndices(), lowerLevelSetId);
594 if (currentLevelSetId > 0) {
595 iterators[currentLevelSetId - 1].goToIndicesSequential(
596 it.getStartIndices());
597 valueBelow = iterators[currentLevelSetId - 1].getValue();
599 valueBelow = std::numeric_limits<T>::max();
604 T velocity = gradNDissipation.first - gradNDissipation.second;
606 maxStepTime += cfl / velocity;
607 tempRates.push_back(std::make_pair(gradNDissipation,
608 -std::numeric_limits<T>::max()));
611 }
else if (velocity == 0.) {
612 maxStepTime = std::numeric_limits<T>::max();
613 tempRates.push_back(std::make_pair(gradNDissipation,
614 std::numeric_limits<T>::max()));
619 T difference = std::abs(valueBelow - value);
621 if (difference >= cfl) {
622 maxStepTime -= cfl / velocity;
623 tempRates.push_back(std::make_pair(
624 gradNDissipation, std::numeric_limits<T>::max()));
627 maxStepTime -= difference / velocity;
630 tempRates.push_back(std::make_pair(gradNDissipation, valueBelow));
638 if (maxStepTime < tempMaxTimeStep)
639 tempMaxTimeStep = maxStepTime;
647 scheme.reduceTimeStepHamiltonJacobi(
648 tempMaxTimeStep, levelSets.back()->getGrid().getGridDelta());
651 if (tempMaxTimeStep < maxTimeStep)
652 maxTimeStep = tempMaxTimeStep;
660 const bool saveVelocities = saveAdvectionVelocities;
661 std::vector<std::vector<double>> dissipationVectors(
662 levelSets.back()->getNumberOfSegments());
663 std::vector<std::vector<double>> velocityVectors(
664 levelSets.back()->getNumberOfSegments());
666 const bool checkDiss = checkDissipation;
668#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
672 p = omp_get_thread_num();
675 auto itRS = totalTempRates[p].cbegin();
676 auto &segment = topDomain.getDomainSegment(p);
677 unsigned maxId = segment.getNumberOfPoints();
679 if (saveVelocities) {
680 velocityVectors[p].resize(maxId);
681 dissipationVectors[p].resize(maxId);
684 for (
unsigned localId = 0; localId < maxId; ++localId) {
685 T &value = segment.definedValues[localId];
686 double time = maxTimeStep;
691 T velocity = itRS->first.first - itRS->first.second;
692 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
693 (itRS->first.first > 0 && velocity < 0)) {
696 T rate = time * velocity;
697 while (std::abs(itRS->second - value) < std::abs(rate)) {
698 time -= std::abs((itRS->second - value) / velocity);
699 value = itRS->second;
702 velocity = itRS->first.first - itRS->first.second;
703 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
704 (itRS->first.first > 0 && velocity < 0)) {
707 rate = time * velocity;
713 if (saveVelocities) {
714 velocityVectors[p][localId] = rate;
715 dissipationVectors[p][localId] = itRS->first.second;
722 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
730 if (saveVelocities) {
731 auto &pointData = levelSets.back()->getPointData();
736 for (
unsigned i = 0; i < velocityVectors.size(); ++i) {
737 vels.insert(vels.end(),
738 std::make_move_iterator(velocityVectors[i].begin()),
739 std::make_move_iterator(velocityVectors[i].end()));
740 diss.insert(diss.end(),
741 std::make_move_iterator(dissipationVectors[i].begin()),
742 std::make_move_iterator(dissipationVectors[i].end()));
744 pointData.insertReplaceScalarData(std::move(vels),
velocityLabel);
758 levelSets.push_back(passedlsDomain);
763 levelSets.push_back(passedlsDomain);
764 velocities = passedVelocities;
769 : levelSets(passedlsDomains) {
770 velocities = passedVelocities;
776 levelSets.push_back(passedlsDomain);
782 velocities = passedVelocities;
838 integrationScheme = scheme;
858 if (levelSets.size() < 1) {
859 Logger::getInstance()
860 .addWarning(
"No level sets passed to Advect.")
867 }
else if (integrationScheme ==
870 }
else if (integrationScheme ==
873 }
else if (integrationScheme ==
876 }
else if (integrationScheme ==
881 }
else if (integrationScheme ==
884 }
else if (integrationScheme ==
887 }
else if (integrationScheme ==
890 }
else if (integrationScheme ==
893 }
else if (integrationScheme ==
898 Logger::getInstance()
899 .addWarning(
"Advect: Integration scheme not found. Not advecting.")
906 if (advectionTime == 0.) {
907 advectedTime = advect();
908 numberOfTimeSteps = 1;
910 double currentTime = 0.0;
911 numberOfTimeSteps = 0;
912 while (currentTime < advectionTime) {
913 currentTime += advect(advectionTime - currentTime);
915 if (performOnlySingleStep)
918 advectedTime = currentTime;