ViennaLS
Loading...
Searching...
No Matches
lsAdvect.hpp
Go to the documentation of this file.
1#pragma once
2
4
5#include <functional>
6#include <limits>
7#include <vector>
8
9#include <hrleSparseIterator.hpp>
10#include <hrleSparseStarIterator.hpp>
11
12#include <vcLogger.hpp>
13#include <vcSmartPointer.hpp>
14
16#include <lsDomain.hpp>
17#include <lsMarkVoidPoints.hpp>
18#include <lsReduce.hpp>
19
20// Spatial discretization schemes
21#include <lsEngquistOsher.hpp>
22#include <lsLaxFriedrichs.hpp>
27#include <lsWENO.hpp>
28
29// Include implementation of time integration schemes
31
32// Velocity accessor
33#include <lsVelocityField.hpp>
34
35// #define DEBUG_LS_ADVECT_HPP
36#ifdef DEBUG_LS_ADVECT_HPP
37#include <lsToMesh.hpp>
38#include <lsVTKWriter.hpp>
39#endif
40
41namespace viennals {
42
43using namespace viennacore;
44
54template <class T, int D> class Advect {
55 using ConstSparseIterator =
56 viennahrle::ConstSparseIterator<typename Domain<T, D>::DomainType>;
57 using hrleIndexType = viennahrle::IndexType;
58
59 // Allow the time integration struct to access private members
61
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 =
83 nullptr;
84
85 // this vector will hold the maximum time step for each point and the
86 // corresponding velocity
87 std::vector<std::vector<std::pair<std::pair<T, T>, T>>> storedRates;
88 double currentTimeStep = -1.;
89
90 VectorType<T, 3> findGlobalAlphas() const {
91
92 auto &topDomain = levelSets.back()->getDomain();
93 auto &grid = levelSets.back()->getGrid();
94
95 const T gridDelta = grid.getGridDelta();
96 const T deltaPos = gridDelta;
97 const T deltaNeg = -gridDelta;
98
99 VectorType<T, 3> finalAlphas = {0., 0., 0.};
100
101#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
102 {
103 VectorType<T, 3> localAlphas = {0., 0., 0.};
104 int p = 0;
105#ifdef _OPENMP
106 p = omp_get_thread_num();
107#endif
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());
115
116 // an iterator for each level set
117 std::vector<ConstSparseIterator> iterators;
118 for (auto const &ls : levelSets) {
119 iterators.emplace_back(ls->getDomain());
120 }
121
122 // neighborIterator for the top level set
123 viennahrle::ConstSparseStarIterator<typename Domain<T, D>::DomainType, 1>
124 neighborIterator(topDomain);
125
126 for (ConstSparseIterator it(topDomain, startVector);
127 it.getStartIndices() < endVector; ++it) {
128
129 if (!it.isDefined() || std::abs(it.getValue()) > integrationCutoff)
130 continue;
131
132 const T value = it.getValue();
133 const auto indices = it.getStartIndices();
134
135 // check if there is any other levelset at the same point:
136 // if yes, take the velocity of the lowest levelset
137 for (unsigned lowerLevelSetId = 0; lowerLevelSetId < levelSets.size();
138 ++lowerLevelSetId) {
139 // put iterator to same position as the top levelset
140 iterators[lowerLevelSetId].goToIndicesSequential(indices);
141
142 // if the lower surface is actually outside, i.e. its LS value
143 // is lower or equal
144 if (iterators[lowerLevelSetId].getValue() <=
145 value + wrappingLayerEpsilon) {
146
147 // move neighborIterator to current position
148 neighborIterator.goToIndicesSequential(indices);
149
150 Vec3D<T> coords{};
151 for (unsigned i = 0; i < D; ++i) {
152 coords[i] = indices[i] * gridDelta;
153 }
154
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();
160
161 normal[i] = phiPos - phiNeg;
162 normalModulus += normal[i] * normal[i];
163 }
164 normalModulus = std::sqrt(normalModulus);
165 for (unsigned i = 0; i < D; ++i)
166 normal[i] /= normalModulus;
167
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());
174
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);
178 }
179
180 // exit material loop
181 break;
182 }
183 }
184 }
185
186#pragma omp critical
187 {
188 for (unsigned i = 0; i < D; ++i) {
189 finalAlphas[i] = std::max(finalAlphas[i], localAlphas[i]);
190 }
191 }
192 } // end of parallel section
193
194 return finalAlphas;
195 }
196
197 // Helper function for linear combination:
198 // target = wTarget * target + wSource * source
199 bool combineLevelSets(T wTarget, T wSource) {
200 // Calculate required expansion width based on CFL and RK steps
201 int steps = 1;
202 if (temporalScheme == TemporalSchemeEnum::RUNGE_KUTTA_2ND_ORDER) {
203 steps = 2;
204 } else if (temporalScheme == TemporalSchemeEnum::RUNGE_KUTTA_3RD_ORDER) {
205 steps = 3;
206 }
207
208 // Expand both level sets to ensure sufficient overlap
209 int expansionWidth = std::ceil(2.0 * steps * timeStepRatio + 1);
210 viennals::Expand<T, D>(levelSets.back(), expansionWidth).apply();
211 viennals::Expand<T, D>(initialLevelSets.back(), expansionWidth).apply();
212
213 bool movedDown = false;
214
215 viennals::BooleanOperation<T, D> op(levelSets.back(),
216 initialLevelSets.back(),
219 [wTarget, wSource, &movedDown](const T &a, const T &b) {
222 T res = wSource * a + wTarget * b;
223 if (res > b)
224 movedDown = true;
225 return std::make_pair(res, true);
226 }
228 return std::make_pair(a, false);
229 return std::make_pair(b, false);
230 });
231 op.apply();
232
233 rebuildLS();
234
235 return movedDown;
236 }
237
238 void rebuildLS() {
239 // TODO: this function uses Manhattan distances for renormalisation,
240 // since this is the quickest. For visualisation applications, better
241 // renormalisation is needed, so it might be good to implement
242 // Euler distance renormalisation as an option
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();
247
248 // Determine cutoff and width based on discretization scheme to avoid
249 // immediate re-expansion
250 T cutoff = 1.0;
251 int finalWidth = 2;
252 if (spatialScheme ==
254 cutoff = 1.5;
255 finalWidth = 3;
256 }
257
258 newDomain.initialize(domain.getNewSegmentation(),
259 domain.getAllocation() *
260 (2.0 / levelSets.back()->getLevelSetWidth()));
261
262 const bool updateData = updatePointData;
263 // save how data should be transferred to new level set
264 // list of indices into the old pointData vector
265 std::vector<std::vector<unsigned>> newDataSourceIds;
266 if (updateData)
267 newDataSourceIds.resize(newDomain.getNumberOfSegments());
268
269#ifdef DEBUG_LS_ADVECT_HPP
270 {
271 auto mesh = SmartPointer<Mesh<T>>::New();
272 ToMesh<T, D>(levelSets.back(), mesh).apply();
273 VTKWriter<T>(mesh, "Advect_beforeRebuild.vtk").apply();
274 }
275#endif
276
277#pragma omp parallel num_threads(newDomain.getNumberOfSegments())
278 {
279 int p = 0;
280#ifdef _OPENMP
281 p = omp_get_thread_num();
282#endif
283 auto &domainSegment = newDomain.getDomainSegment(p);
284
285 viennahrle::Index<D> startVector =
286 (p == 0) ? grid.getMinGridPoint()
287 : newDomain.getSegmentation()[p - 1];
288
289 viennahrle::Index<D> endVector =
290 (p != static_cast<int>(newDomain.getNumberOfSegments() - 1))
291 ? newDomain.getSegmentation()[p]
292 : grid.incrementIndices(grid.getMaxGridPoint());
293
294 // reserve a bit more to avoid reallocation
295 // would expect number of points to roughly double
296 if (updateData)
297 newDataSourceIds[p].reserve(2.5 * domainSegment.getNumberOfPoints());
298
299 for (viennahrle::SparseStarIterator<typename Domain<T, D>::DomainType, 1>
300 it(domain, startVector);
301 it.getIndices() < endVector; ++it) {
302
303 // if the center is an active grid point
304 // <1.0 since it could have been change by 0.5 max
305 if (std::abs(it.getCenter().getValue()) <= 1.0) {
306
307 int k = 0;
308 for (; k < 2 * D; k++)
309 if (std::signbit(it.getNeighbor(k).getValue() - 1e-7) !=
310 std::signbit(it.getCenter().getValue() + 1e-7))
311 break;
312
313 // if there is at least one neighbor of opposite sign
314 if (k != 2 * D) {
315 if (it.getCenter().getDefinedValue() > 0.5) {
316 int j = 0;
317 for (; j < 2 * D; j++) {
318 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
319 if (it.getNeighbor(j).getDefinedValue() < -0.5)
320 break;
321 }
322 if (j == 2 * D) {
323 domainSegment.insertNextDefinedPoint(
324 it.getIndices(), it.getCenter().getDefinedValue());
325 if (updateData)
326 newDataSourceIds[p].push_back(it.getCenter().getPointId());
327 // if there is at least one active grid point, which is < -0.5
328 } else {
329 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
330 if (updateData)
331 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
332 }
333 } else if (it.getCenter().getDefinedValue() < -0.5) {
334 int j = 0;
335 for (; j < 2 * D; j++) {
336 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
337 if (it.getNeighbor(j).getDefinedValue() > 0.5)
338 break;
339 }
340
341 if (j == 2 * D) {
342 domainSegment.insertNextDefinedPoint(
343 it.getIndices(), it.getCenter().getDefinedValue());
344 if (updateData)
345 newDataSourceIds[p].push_back(it.getCenter().getPointId());
346 // if there is at least one active grid point, which is > 0.5
347 } else {
348 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
349 if (updateData)
350 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
351 }
352 } else {
353 domainSegment.insertNextDefinedPoint(
354 it.getIndices(), it.getCenter().getDefinedValue());
355 if (updateData)
356 newDataSourceIds[p].push_back(it.getCenter().getPointId());
357 }
358 } else {
359 domainSegment.insertNextUndefinedPoint(
360 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
363 }
364
365 } else { // if the center is not an active grid point
366 if (it.getCenter().getValue() >= 0) {
367 int usedNeighbor = -1;
368 T distance = Domain<T, D>::POS_VALUE;
369 for (int i = 0; i < 2 * D; i++) {
370 T value = it.getNeighbor(i).getValue();
371 if (std::abs(value) <= 1.0 && (value < 0.)) {
372 if (distance > value + 1.0) {
373 distance = value + 1.0;
374 usedNeighbor = i;
375 }
376 }
377 }
378
379 if (distance <= cutoff) {
380 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
381 if (updateData)
382 newDataSourceIds[p].push_back(
383 it.getNeighbor(usedNeighbor).getPointId());
384 } else {
385 domainSegment.insertNextUndefinedPoint(it.getIndices(),
387 }
388
389 } else {
390 int usedNeighbor = -1;
391 T distance = Domain<T, D>::NEG_VALUE;
392 for (int i = 0; i < 2 * D; i++) {
393 T value = it.getNeighbor(i).getValue();
394 if (std::abs(value) <= 1.0 && (value > 0)) {
395 if (distance < value - 1.0) {
396 // distance = std::max(distance, value - T(1.0));
397 distance = value - 1.0;
398 usedNeighbor = i;
399 }
400 }
401 }
402
403 if (distance >= -cutoff) {
404 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
405 if (updateData)
406 newDataSourceIds[p].push_back(
407 it.getNeighbor(usedNeighbor).getPointId());
408 } else {
409 domainSegment.insertNextUndefinedPoint(it.getIndices(),
411 }
412 }
413 }
414 }
415 }
416
417 // now copy old data into new level set
418 if (updateData) {
419 auto &pointData = levelSets.back()->getPointData();
420 newlsDomain->getPointData().translateFromMultiData(pointData,
421 newDataSourceIds);
422 }
423
424 newDomain.finalize();
425 newDomain.segment();
426 levelSets.back()->deepCopy(newlsDomain);
427 levelSets.back()->finalize(finalWidth);
428 }
429
434 template <class DiscretizationSchemeType>
435 double integrateTime(DiscretizationSchemeType spatialScheme,
436 double maxTimeStep) {
437
438 auto &topDomain = levelSets.back()->getDomain();
439 auto &grid = levelSets.back()->getGrid();
440
441 typename PointData<T>::ScalarDataType *voidMarkerPointer;
442 if (ignoreVoids) {
443 MarkVoidPoints<T, D>(levelSets.back()).apply();
444 auto &pointData = levelSets.back()->getPointData();
445 voidMarkerPointer =
447 if (voidMarkerPointer == nullptr) {
448 VIENNACORE_LOG_WARNING("Advect: Cannot find void point markers. Not "
449 "ignoring void points.");
450 ignoreVoids = false;
451 }
452 }
453 const bool ignoreVoidPoints = ignoreVoids;
454 const bool useAdaptiveTimeStepping = adaptiveTimeStepping;
455 const auto adaptiveFactor = 1.0 / adaptiveTimeStepSubdivisions;
456
457 if (!storedRates.empty()) {
458 VIENNACORE_LOG_WARNING("Advect: Overwriting previously stored rates.");
459 }
460
461 storedRates.resize(topDomain.getNumberOfSegments());
462
463#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
464 {
465 int p = 0;
466#ifdef _OPENMP
467 p = omp_get_thread_num();
468#endif
469 viennahrle::Index<D> startVector =
470 (p == 0) ? grid.getMinGridPoint()
471 : topDomain.getSegmentation()[p - 1];
472
473 viennahrle::Index<D> endVector =
474 (p != static_cast<int>(topDomain.getNumberOfSegments() - 1))
475 ? topDomain.getSegmentation()[p]
476 : grid.incrementIndices(grid.getMaxGridPoint());
477
478 double tempMaxTimeStep = maxTimeStep;
479 // store the rates and value of underneath LS for this segment
480 auto &rates = storedRates[p];
481 rates.reserve(
482 topDomain.getNumberOfPoints() /
483 static_cast<double>((levelSets.back())->getNumberOfSegments()) +
484 10);
485
486 // an iterator for each level set
487 std::vector<ConstSparseIterator> iterators;
488 for (auto const &ls : levelSets) {
489 iterators.emplace_back(ls->getDomain());
490 }
491
492 DiscretizationSchemeType scheme(spatialScheme);
493
494 for (ConstSparseIterator it(topDomain, startVector);
495 it.getStartIndices() < endVector; ++it) {
496
497 if (!it.isDefined() || std::abs(it.getValue()) > integrationCutoff)
498 continue;
499
500 T value = it.getValue();
501 double maxStepTime = 0;
502 double cfl = timeStepRatio;
503
504 for (int currentLevelSetId = levelSets.size() - 1;
505 currentLevelSetId >= 0; --currentLevelSetId) {
506
507 std::pair<T, T> gradNDissipation;
508
509 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
510 // check if there is any other levelset at the same point:
511 // if yes, take the velocity of the lowest levelset
512 for (unsigned lowerLevelSetId = 0;
513 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
514 // put iterator to same position as the top levelset
515 iterators[lowerLevelSetId].goToIndicesSequential(
516 it.getStartIndices());
517
518 // if the lower surface is actually outside, i.e. its LS value
519 // is lower or equal
520 if (iterators[lowerLevelSetId].getValue() <=
521 value + wrappingLayerEpsilon) {
522 gradNDissipation =
523 scheme(it.getStartIndices(), lowerLevelSetId);
524 break;
525 }
526 }
527 }
528
529 T velocity = gradNDissipation.first - gradNDissipation.second;
530 if (velocity > 0.) {
531 // Case 1: Growth / Deposition (Velocity > 0)
532 // Limit the time step based on the standard CFL condition.
533 maxStepTime += cfl / velocity;
534 rates.emplace_back(gradNDissipation,
535 -std::numeric_limits<T>::max());
536 break;
537 } else if (velocity == 0.) {
538 // Case 2: Static (Velocity == 0)
539 // No time step limit imposed by this point.
540 maxStepTime = std::numeric_limits<T>::max();
541 rates.emplace_back(gradNDissipation, std::numeric_limits<T>::max());
542 break;
543 } else {
544 // Case 3: Etching (Velocity < 0)
545 // Retrieve the interface location of the underlying material.
546 T valueBelow;
547 if (currentLevelSetId > 0) {
548 iterators[currentLevelSetId - 1].goToIndicesSequential(
549 it.getStartIndices());
550 valueBelow = iterators[currentLevelSetId - 1].getValue();
551 } else {
552 valueBelow = std::numeric_limits<T>::max();
553 }
554 // Calculate the top material thickness
555 T difference = std::abs(valueBelow - value);
556
557 if (difference >= cfl) {
558 // Sub-case 3a: Standard Advection
559 // Far from interface: Use full CFL time step.
560 maxStepTime -= cfl / velocity;
561 rates.emplace_back(gradNDissipation,
562 std::numeric_limits<T>::max());
563 break;
564 } else {
565 // Sub-case 3b: Interface Interaction
566 // Use adaptiveFactor as threshold.
567 if (useAdaptiveTimeStepping &&
568 difference > adaptiveFactor * cfl) {
569 // Adaptive Sub-stepping:
570 // Approaching boundary: Force small steps to gather
571 // flux statistics and prevent numerical overshoot ("Soft
572 // Landing").
573 maxStepTime -= adaptiveFactor * cfl / velocity;
574 rates.emplace_back(gradNDissipation,
575 std::numeric_limits<T>::max());
576 break;
577 } else {
578 // Terminal Step:
579 // Within tolerance: Snap to boundary, consume budget, and
580 // switch material.
581 cfl -= difference;
582 value = valueBelow;
583 maxStepTime -= difference / velocity;
584 rates.emplace_back(gradNDissipation, valueBelow);
585 }
586 }
587 }
588 }
589
590 if (maxStepTime < tempMaxTimeStep)
591 tempMaxTimeStep = maxStepTime;
592 }
593
594#pragma omp critical
595 {
596 // If a Lax Friedrichs scheme is selected the time step is
597 // reduced depending on the dissipation coefficients
598 // For Engquist Osher scheme this function is empty.
599 scheme.reduceTimeStepHamiltonJacobi(
600 tempMaxTimeStep, levelSets.back()->getGrid().getGridDelta());
601
602 // set global timestep maximum
603 if (tempMaxTimeStep < maxTimeStep)
604 maxTimeStep = tempMaxTimeStep;
605 }
606 } // end of parallel section
607
608 // maxTimeStep is now the maximum time step possible for all points
609 // and rates are stored in a vector
610 return maxTimeStep;
611 }
612
615 void computeRates(double maxTimeStep = std::numeric_limits<double>::max()) {
616 prepareLS();
617 if (spatialScheme == SpatialSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
618 auto is = lsInternal::EngquistOsher<T, D, 1>(levelSets.back(), velocities,
619 calculateNormalVectors);
620 currentTimeStep = integrateTime(is, maxTimeStep);
621 } else if (spatialScheme == SpatialSchemeEnum::ENGQUIST_OSHER_2ND_ORDER) {
622 auto is = lsInternal::EngquistOsher<T, D, 2>(levelSets.back(), velocities,
623 calculateNormalVectors);
624 currentTimeStep = integrateTime(is, maxTimeStep);
625 } else if (spatialScheme == SpatialSchemeEnum::LAX_FRIEDRICHS_1ST_ORDER) {
626 auto alphas = findGlobalAlphas();
627 auto is = lsInternal::LaxFriedrichs<T, D, 1>(levelSets.back(), velocities,
628 dissipationAlpha, alphas,
629 calculateNormalVectors);
630 currentTimeStep = integrateTime(is, maxTimeStep);
631 } else if (spatialScheme == SpatialSchemeEnum::LAX_FRIEDRICHS_2ND_ORDER) {
632 auto alphas = findGlobalAlphas();
633 auto is = lsInternal::LaxFriedrichs<T, D, 2>(levelSets.back(), velocities,
634 dissipationAlpha, alphas,
635 calculateNormalVectors);
636 currentTimeStep = integrateTime(is, maxTimeStep);
637 } else if (spatialScheme ==
640 levelSets.back(), velocities);
641 currentTimeStep = integrateTime(is, maxTimeStep);
642 } else if (spatialScheme ==
645 levelSets.back(), velocities, dissipationAlpha);
646 currentTimeStep = integrateTime(is, maxTimeStep);
647 } else if (spatialScheme ==
650 levelSets.back(), velocities, dissipationAlpha);
651 currentTimeStep = integrateTime(is, maxTimeStep);
652 } else if (spatialScheme ==
655 levelSets.back(), velocities, dissipationAlpha);
656 currentTimeStep = integrateTime(is, maxTimeStep);
657 } else if (spatialScheme ==
660 levelSets.back(), velocities, dissipationAlpha);
661 currentTimeStep = integrateTime(is, maxTimeStep);
662 } else if (spatialScheme ==
665 levelSets.back(), velocities, dissipationAlpha);
666 currentTimeStep = integrateTime(is, maxTimeStep);
667 } else if (spatialScheme == SpatialSchemeEnum::WENO_3RD_ORDER) {
668 // Instantiate WENO with order 3
669 auto is = lsInternal::WENO<T, D, 3>(levelSets.back(), velocities,
670 dissipationAlpha);
671 currentTimeStep = integrateTime(is, maxTimeStep);
672 } else if (spatialScheme == SpatialSchemeEnum::WENO_5TH_ORDER) {
673 // Instantiate WENO with order 5
674 auto is = lsInternal::WENO<T, D, 5>(levelSets.back(), velocities,
675 dissipationAlpha);
676 currentTimeStep = integrateTime(is, maxTimeStep);
677 } else {
678 VIENNACORE_LOG_ERROR("Advect: Discretization scheme not found.");
679 currentTimeStep = -1.;
680 }
681 }
682
683 // Level Sets below are also considered in order to adjust the advection
684 // depth accordingly if there would be a material change.
685 void updateLevelSet(double dt) {
686 if (timeStepRatio >= 0.5) {
687 VIENNACORE_LOG_WARNING(
688 "Integration time step ratio should be smaller than 0.5. "
689 "Advection might fail!");
690 }
691
692 auto &topDomain = levelSets.back()->getDomain();
693
694 assert(dt >= 0. && "No time step set!");
695 assert(storedRates.size() == topDomain.getNumberOfSegments());
696
697 // reduce to one layer thickness and apply new values directly to the
698 // domain segments --> DO NOT CHANGE SEGMENTATION HERE (true parameter)
699 Reduce<T, D>(levelSets.back(), 1, true).apply();
700
701 const bool saveVelocities = saveAdvectionVelocities;
702 std::vector<std::vector<double>> dissipationVectors(
703 levelSets.back()->getNumberOfSegments());
704 std::vector<std::vector<double>> velocityVectors(
705 levelSets.back()->getNumberOfSegments());
706
707 const bool checkDiss = checkDissipation;
708
709#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
710 {
711 int p = 0;
712#ifdef _OPENMP
713 p = omp_get_thread_num();
714#endif
715 auto itRS = storedRates[p].cbegin();
716 auto &segment = topDomain.getDomainSegment(p);
717 const unsigned maxId = segment.getNumberOfPoints();
718
719 if (saveVelocities) {
720 velocityVectors[p].resize(maxId);
721 dissipationVectors[p].resize(maxId);
722 }
723
724 for (unsigned localId = 0; localId < maxId; ++localId) {
725 T &value = segment.definedValues[localId];
726
727 // Skip points that were not part of computeRates (outer layers)
728 if (std::abs(value) > integrationCutoff)
729 continue;
730
731 double time = dt;
732
733 // if there is a change in materials during one time step, deduct
734 // the time taken to advect up to the end of the top material and
735 // set the LS value to the one below
736 auto const [gradient, dissipation] = itRS->first;
737 T velocity = gradient - dissipation;
738 // check if dissipation is too high and would cause a change in
739 // direction of the velocity
740 if (checkDiss && (gradient < 0 && velocity > 0) ||
741 (gradient > 0 && velocity < 0)) {
742 velocity = 0;
743 }
744
745 T rate = time * velocity;
746 while (std::abs(itRS->second - value) < std::abs(rate)) {
747 time -= std::abs((itRS->second - value) / velocity);
748 value = itRS->second;
749 ++itRS; // advance the TempStopRates iterator by one
750
751 // recalculate velocity and rate
752 velocity = itRS->first.first - itRS->first.second;
753 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
754 (itRS->first.first > 0 && velocity < 0)) {
755 velocity = 0;
756 }
757 rate = time * velocity;
758 }
759
760 // now deduct the velocity times the time step we take
761 value -= rate;
762
763 if (saveVelocities) {
764 velocityVectors[p][localId] = rate;
765 dissipationVectors[p][localId] = itRS->first.second;
766 }
767
768 // this is run when two materials are close but the velocity is too slow
769 // to actually reach the second material, to get rid of the extra
770 // entry in the TempRatesStop
771 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
772 ++itRS;
773
774 // advance the TempStopRates iterator by one
775 ++itRS;
776 }
777 } // end of parallel section
778
779 if (saveVelocities) {
780 auto &pointData = levelSets.back()->getPointData();
781
782 typename PointData<T>::ScalarDataType vels;
783 typename PointData<T>::ScalarDataType diss;
784
785 for (unsigned i = 0; i < velocityVectors.size(); ++i) {
786 vels.insert(vels.end(),
787 std::make_move_iterator(velocityVectors[i].begin()),
788 std::make_move_iterator(velocityVectors[i].end()));
789 diss.insert(diss.end(),
790 std::make_move_iterator(dissipationVectors[i].begin()),
791 std::make_move_iterator(dissipationVectors[i].end()));
792 }
793 pointData.insertReplaceScalarData(std::move(vels), velocityLabel);
794 pointData.insertReplaceScalarData(std::move(diss), dissipationLabel);
795 }
796
797 // clear the stored rates since surface has changed
798 storedRates.clear();
799 }
800
801 void adjustLowerLayers() {
802 // Adjust all level sets below the advected one
803 if (spatialScheme !=
805 for (unsigned i = 0; i < levelSets.size() - 1; ++i) {
807 levelSets[i], levelSets.back(),
809 .apply();
810 }
811 }
812 }
813
816 double advect(double maxTimeStep) {
817 switch (temporalScheme) {
820 *this, maxTimeStep);
823 *this, maxTimeStep);
825 default:
827 *this, maxTimeStep);
828 }
829 }
830
831public:
832 static constexpr char velocityLabel[] = "AdvectionVelocities";
833 static constexpr char dissipationLabel[] = "Dissipation";
834
835 Advect() = default;
836
837 explicit Advect(SmartPointer<Domain<T, D>> passedlsDomain) {
838 levelSets.push_back(passedlsDomain);
839 }
840
841 Advect(SmartPointer<Domain<T, D>> passedlsDomain,
842 SmartPointer<VelocityField<T>> passedVelocities) {
843 levelSets.push_back(passedlsDomain);
844 velocities = passedVelocities;
845 }
846
847 Advect(std::vector<SmartPointer<Domain<T, D>>> passedlsDomains,
848 SmartPointer<VelocityField<T>> passedVelocities)
849 : levelSets(passedlsDomains) {
850 velocities = passedVelocities;
851 }
852
855 void insertNextLevelSet(SmartPointer<Domain<T, D>> passedlsDomain) {
856 levelSets.push_back(passedlsDomain);
857 }
858
859 void clearLevelSets() { levelSets.clear(); }
860
863 void setVelocityField(SmartPointer<VelocityField<T>> passedVelocities) {
864 velocities = passedVelocities;
865 }
866
872 void setAdvectionTime(double time) { advectionTime = time; }
873
878 void setSingleStep(bool singleStep) { performOnlySingleStep = singleStep; }
879
884 void setTimeStepRatio(const double &cfl) { timeStepRatio = cfl; }
885
890 void setCalculateNormalVectors(bool cnv) { calculateNormalVectors = cnv; }
891
898 void setIgnoreVoids(bool iV) { ignoreVoids = iV; }
899
903 void setAdaptiveTimeStepping(bool aTS = true, unsigned subdivisions = 20) {
904 adaptiveTimeStepping = aTS;
905 if (subdivisions < 1) {
906 VIENNACORE_LOG_WARNING("Advect: Adaptive time stepping subdivisions must "
907 "be at least 1. Setting to 1.");
908 subdivisions = 1;
909 }
910 adaptiveTimeStepSubdivisions = subdivisions;
911 }
912
915 void setSaveAdvectionVelocities(bool sAV) { saveAdvectionVelocities = sAV; }
916
919 double getAdvectedTime() const { return advectedTime; }
920
922 double getCurrentTimeStep() const { return currentTimeStep; }
923
925 unsigned getNumberOfTimeSteps() const { return numberOfTimeSteps; }
926
928 double getTimeStepRatio() const { return timeStepRatio; }
929
931 bool getCalculateNormalVectors() const { return calculateNormalVectors; }
932
935 void setSpatialScheme(SpatialSchemeEnum scheme) { spatialScheme = scheme; }
936
937 // Deprecated and will be removed in future versions:
938 // use setSpatialScheme instead
939 [[deprecated("Use setSpatialScheme instead")]] void
940 setIntegrationScheme(IntegrationSchemeEnum scheme) {
941 VIENNACORE_LOG_WARNING(
942 "Advect::setIntegrationScheme is deprecated and will be removed in "
943 "future versions. Use setSpatialScheme instead.");
944 spatialScheme = scheme;
945 }
946
948 void setTemporalScheme(TemporalSchemeEnum scheme) { temporalScheme = scheme; }
949
954 void setDissipationAlpha(const double &a) { dissipationAlpha = a; }
955
956 // Sets the velocity to 0 if the dissipation is too high
957 void setCheckDissipation(bool check) { checkDissipation = check; }
958
961 void setUpdatePointData(bool update) { updatePointData = update; }
962
966 std::function<bool(SmartPointer<Domain<T, D>>)> callback) {
967 velocityUpdateCallback = callback;
968 }
969
970 // Prepare the levelset for advection, based on the provided spatial
971 // discretization scheme.
972 void prepareLS() {
973 // check whether a level set and velocities have been given
974 if (levelSets.empty()) {
975 VIENNACORE_LOG_ERROR("No level sets passed to Advect.");
976 return;
977 }
978
979 if (spatialScheme == SpatialSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
981 } else if (spatialScheme == SpatialSchemeEnum::ENGQUIST_OSHER_2ND_ORDER) {
983 } else if (spatialScheme == SpatialSchemeEnum::LAX_FRIEDRICHS_1ST_ORDER) {
985 } else if (spatialScheme == SpatialSchemeEnum::LAX_FRIEDRICHS_2ND_ORDER) {
987 } else if (spatialScheme ==
990 levelSets.back());
991 } else if (spatialScheme ==
994 } else if (spatialScheme ==
997 } else if (spatialScheme ==
1000 } else if (spatialScheme ==
1003 } else if (spatialScheme ==
1006 levelSets.back());
1007 } else if (spatialScheme == SpatialSchemeEnum::WENO_3RD_ORDER) {
1008 lsInternal::WENO<T, D, 3>::prepareLS(levelSets.back());
1009 } else if (spatialScheme == SpatialSchemeEnum::WENO_5TH_ORDER) {
1010 lsInternal::WENO<T, D, 5>::prepareLS(levelSets.back());
1011 } else {
1012 VIENNACORE_LOG_ERROR("Advect: Discretization scheme not found.");
1013 }
1014 }
1015
1016 void apply() {
1017 // check whether a level set and velocities have been given
1018 if (levelSets.empty()) {
1019 VIENNACORE_LOG_ERROR("No level sets passed to Advect. Not advecting.");
1020 return;
1021 }
1022 if (velocities == nullptr) {
1023 VIENNACORE_LOG_ERROR(
1024 "No velocity field passed to Advect. Not advecting.");
1025 return;
1026 }
1027
1028 if (advectionTime == 0.) {
1029 advectedTime = advect(std::numeric_limits<double>::max());
1030 numberOfTimeSteps = 1;
1031 } else {
1032 double currentTime = 0.0;
1033 numberOfTimeSteps = 0;
1034 while (currentTime < advectionTime) {
1035 currentTime += advect(advectionTime - currentTime);
1036 ++numberOfTimeSteps;
1037 if (performOnlySingleStep)
1038 break;
1039 }
1040 advectedTime = currentTime;
1041 }
1042 }
1043};
1044
1045// add all template specializations for this class
1047
1048} // namespace viennals
constexpr int D
Definition Epitaxy.cpp:11
double T
Definition Epitaxy.cpp:12
Engquist-Osher spatial discretization scheme based on the upwind spatial discretization scheme....
Definition lsEngquistOsher.hpp:18
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsEngquistOsher.hpp:28
Lax Friedrichs spatial discretization scheme with constant alpha value for dissipation....
Definition lsLaxFriedrichs.hpp:19
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLaxFriedrichs.hpp:32
Lax Friedrichs spatial discretization scheme, which uses alpha values provided by the user in getDiss...
Definition lsLocalLaxFriedrichsAnalytical.hpp:20
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLaxFriedrichsAnalytical.hpp:48
Lax Friedrichs spatial discretization scheme, which uses a first neighbour stencil to calculate the a...
Definition lsLocalLaxFriedrichs.hpp:20
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLaxFriedrichs.hpp:49
Lax Friedrichs spatial discretization scheme, which considers only the current point for alpha calcul...
Definition lsLocalLocalLaxFriedrichs.hpp:18
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLocalLaxFriedrichs.hpp:29
Stencil Local Lax Friedrichs Discretization Scheme. It uses a stencil of order around active points,...
Definition lsStencilLocalLaxFriedrichsScalar.hpp:33
static void prepareLS(LevelSetType passedlsDomain)
Definition lsStencilLocalLaxFriedrichsScalar.hpp:197
Weighted Essentially Non-Oscillatory (WENO) scheme. This kernel acts as the grid-interface for the ma...
Definition lsWENO.hpp:20
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsWENO.hpp:42
void setVelocityField(SmartPointer< VelocityField< T > > passedVelocities)
Set the velocity field used for advection. This should be a concrete implementation of lsVelocityFiel...
Definition lsAdvect.hpp:863
static constexpr char velocityLabel[]
Definition lsAdvect.hpp:832
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:855
void setVelocityUpdateCallback(std::function< bool(SmartPointer< Domain< T, D > >)> callback)
Set a callback function that is called after the level set has been updated during intermediate time ...
Definition lsAdvect.hpp:965
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:961
void setSpatialScheme(SpatialSchemeEnum scheme)
Set which spatial discretization scheme should be used out of the ones specified in SpatialSchemeEnum...
Definition lsAdvect.hpp:935
bool getCalculateNormalVectors() const
Get whether normal vectors were calculated.
Definition lsAdvect.hpp:931
void setIgnoreVoids(bool iV)
Set whether level set values, which are not part of the "top" geometrically connected part of values,...
Definition lsAdvect.hpp:898
Advect()=default
void clearLevelSets()
Definition lsAdvect.hpp:859
void setCalculateNormalVectors(bool cnv)
Set whether normal vectors should be calculated at each level set point. Defaults to true....
Definition lsAdvect.hpp:890
void apply()
Definition lsAdvect.hpp:1016
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:878
void prepareLS()
Definition lsAdvect.hpp:972
double getAdvectedTime() const
Get by how much the physical time was advanced during the last apply() call.
Definition lsAdvect.hpp:919
void setTemporalScheme(TemporalSchemeEnum scheme)
Set which time integration scheme should be used.
Definition lsAdvect.hpp:948
void setDissipationAlpha(const double &a)
Set the alpha dissipation coefficient. For lsLaxFriedrichs, this is used as the alpha value....
Definition lsAdvect.hpp:954
void setCheckDissipation(bool check)
Definition lsAdvect.hpp:957
Advect(std::vector< SmartPointer< Domain< T, D > > > passedlsDomains, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:847
void setIntegrationScheme(IntegrationSchemeEnum scheme)
Definition lsAdvect.hpp:940
double getCurrentTimeStep() const
Return the last applied time step.
Definition lsAdvect.hpp:922
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:872
unsigned getNumberOfTimeSteps() const
Get how many advection steps were performed during the last apply() call.
Definition lsAdvect.hpp:925
void setAdaptiveTimeStepping(bool aTS=true, unsigned subdivisions=20)
Set whether adaptive time stepping should be used when approaching material boundaries during etching...
Definition lsAdvect.hpp:903
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:884
Advect(SmartPointer< Domain< T, D > > passedlsDomain)
Definition lsAdvect.hpp:837
Advect(SmartPointer< Domain< T, D > > passedlsDomain, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:841
double getTimeStepRatio() const
Get the value of the CFL number.
Definition lsAdvect.hpp:928
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:915
static constexpr char dissipationLabel[]
Definition lsAdvect.hpp:833
This class is used to perform boolean operations on two level sets and write the resulting level set ...
Definition lsBooleanOperation.hpp:45
void setBooleanOperationComparator(ComparatorType passedOperationComp)
Set the comparator to be used when the BooleanOperation is set to CUSTOM.
Definition lsBooleanOperation.hpp:307
void apply()
Perform operation.
Definition lsBooleanOperation.hpp:319
Class containing all information about the level set, including the dimensions of the domain,...
Definition lsDomain.hpp:28
viennahrle::Domain< T, D > DomainType
Definition lsDomain.hpp:33
static constexpr T NEG_VALUE
Definition lsDomain.hpp:53
static constexpr T POS_VALUE
Definition lsDomain.hpp:52
Expands the levelSet to the specified number of layers. The largest value in the levelset is thus wid...
Definition lsExpand.hpp:17
void apply()
Apply the expansion to the specified width.
Definition lsExpand.hpp:44
This class is used to mark points of the level set which are enclosed in a void.
Definition lsMarkVoidPoints.hpp:28
void apply()
Definition lsMarkVoidPoints.hpp:154
static constexpr char voidPointLabel[]
Definition lsMarkVoidPoints.hpp:86
std::vector< T > ScalarDataType
Definition lsPointData.hpp:24
ScalarDataType * getScalarData(int index)
Definition lsPointData.hpp:134
Reduce the level set size to the specified width. This means all level set points with value <= 0....
Definition lsReduce.hpp:14
void apply()
Reduces the leveleSet to the specified number of layers. The largest value in the levelset is thus wi...
Definition lsReduce.hpp:55
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:33
void apply()
Definition lsVTKWriter.hpp:136
Abstract class defining the interface for the velocity field used during advection using lsAdvect.
Definition lsVelocityField.hpp:11
#define PRECOMPILE_PRECISION_DIMENSION(className)
Definition lsPreCompileMacros.hpp:24
Definition lsAdvect.hpp:41
SpatialSchemeEnum
Enumeration for the different spatial discretization schemes used by the advection kernel.
Definition lsAdvectIntegrationSchemes.hpp:10
@ LOCAL_LOCAL_LAX_FRIEDRICHS_2ND_ORDER
Definition lsAdvectIntegrationSchemes.hpp:17
@ STENCIL_LOCAL_LAX_FRIEDRICHS_1ST_ORDER
Definition lsAdvectIntegrationSchemes.hpp:20
@ LOCAL_LOCAL_LAX_FRIEDRICHS_1ST_ORDER
Definition lsAdvectIntegrationSchemes.hpp:16
@ LAX_FRIEDRICHS_2ND_ORDER
Definition lsAdvectIntegrationSchemes.hpp:14
@ LOCAL_LAX_FRIEDRICHS_1ST_ORDER
Definition lsAdvectIntegrationSchemes.hpp:18
@ ENGQUIST_OSHER_2ND_ORDER
Definition lsAdvectIntegrationSchemes.hpp:12
@ LAX_FRIEDRICHS_1ST_ORDER
Definition lsAdvectIntegrationSchemes.hpp:13
@ LOCAL_LAX_FRIEDRICHS_2ND_ORDER
Definition lsAdvectIntegrationSchemes.hpp:19
@ WENO_3RD_ORDER
Definition lsAdvectIntegrationSchemes.hpp:21
@ ENGQUIST_OSHER_1ST_ORDER
Definition lsAdvectIntegrationSchemes.hpp:11
@ LOCAL_LAX_FRIEDRICHS_ANALYTICAL_1ST_ORDER
Definition lsAdvectIntegrationSchemes.hpp:15
@ WENO_5TH_ORDER
Definition lsAdvectIntegrationSchemes.hpp:22
TemporalSchemeEnum
Enumeration for the different time integration schemes used to select the advection kernel.
Definition lsAdvectIntegrationSchemes.hpp:31
@ RUNGE_KUTTA_2ND_ORDER
Definition lsAdvectIntegrationSchemes.hpp:33
@ RUNGE_KUTTA_3RD_ORDER
Definition lsAdvectIntegrationSchemes.hpp:34
@ FORWARD_EULER
Definition lsAdvectIntegrationSchemes.hpp:32
@ INTERSECT
Definition lsBooleanOperation.hpp:28
@ CUSTOM
Definition lsBooleanOperation.hpp:32
Definition lsAdvectIntegrationSchemes.hpp:43
static double evolveRungeKutta3(AdvectType &kernel, double maxTimeStep)
Definition lsAdvectIntegrationSchemes.hpp:110
static double evolveForwardEuler(AdvectType &kernel, double maxTimeStep, bool updateLowerLayers=true)
Definition lsAdvectIntegrationSchemes.hpp:46
static double evolveRungeKutta2(AdvectType &kernel, double maxTimeStep)
Definition lsAdvectIntegrationSchemes.hpp:61