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, D> 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, D> finalAlphas{};
100
101#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
102 {
103 VectorType<T, D> localAlphas{};
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 = 1. / 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::ConstSparseStarIterator<
300 typename Domain<T, D>::DomainType, 1>
301 it(domain, startVector);
302 it.getIndices() < endVector; ++it) {
303
304 // if the center is an active grid point
305 // <1.0 since it could have been change by 0.5 max
306 if (std::abs(it.getCenter().getValue()) <= 1.0) {
307
308 int k = 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))
312 break;
313
314 // if there is at least one neighbor of opposite sign
315 if (k != 2 * D) {
316 if (it.getCenter().getDefinedValue() > 0.5) {
317 int j = 0;
318 for (; j < 2 * D; j++) {
319 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
320 if (it.getNeighbor(j).getDefinedValue() < -0.5)
321 break;
322 }
323 if (j == 2 * D) {
324 domainSegment.insertNextDefinedPoint(
325 it.getIndices(), it.getCenter().getDefinedValue());
326 if (updateData)
327 newDataSourceIds[p].push_back(it.getCenter().getPointId());
328 // if there is at least one active grid point, which is < -0.5
329 } else {
330 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
331 if (updateData)
332 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
333 }
334 } else if (it.getCenter().getDefinedValue() < -0.5) {
335 int j = 0;
336 for (; j < 2 * D; j++) {
337 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
338 if (it.getNeighbor(j).getDefinedValue() > 0.5)
339 break;
340 }
341
342 if (j == 2 * D) {
343 domainSegment.insertNextDefinedPoint(
344 it.getIndices(), it.getCenter().getDefinedValue());
345 if (updateData)
346 newDataSourceIds[p].push_back(it.getCenter().getPointId());
347 // if there is at least one active grid point, which is > 0.5
348 } else {
349 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
350 if (updateData)
351 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
352 }
353 } else {
354 domainSegment.insertNextDefinedPoint(
355 it.getIndices(), it.getCenter().getDefinedValue());
356 if (updateData)
357 newDataSourceIds[p].push_back(it.getCenter().getPointId());
358 }
359 } else {
360 domainSegment.insertNextUndefinedPoint(
361 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
364 }
365
366 } else { // if the center is not an active grid point
367 if (it.getCenter().getValue() >= 0) {
368 int usedNeighbor = -1;
369 T distance = Domain<T, D>::POS_VALUE;
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;
375 usedNeighbor = i;
376 }
377 }
378 }
379
380 if (distance <= cutoff) {
381 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
382 if (updateData)
383 newDataSourceIds[p].push_back(
384 it.getNeighbor(usedNeighbor).getPointId());
385 } else {
386 domainSegment.insertNextUndefinedPoint(it.getIndices(),
388 }
389
390 } else {
391 int usedNeighbor = -1;
392 T distance = Domain<T, D>::NEG_VALUE;
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) {
397 // distance = std::max(distance, value - T(1.0));
398 distance = value - 1.0;
399 usedNeighbor = i;
400 }
401 }
402 }
403
404 if (distance >= -cutoff) {
405 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
406 if (updateData)
407 newDataSourceIds[p].push_back(
408 it.getNeighbor(usedNeighbor).getPointId());
409 } else {
410 domainSegment.insertNextUndefinedPoint(it.getIndices(),
412 }
413 }
414 }
415 }
416 }
417
418 // now copy old data into new level set
419 if (updateData) {
420 auto &pointData = levelSets.back()->getPointData();
421 newlsDomain->getPointData().translateFromMultiData(pointData,
422 newDataSourceIds);
423 }
424
425 newDomain.finalize();
426 newDomain.segment();
427 levelSets.back()->deepCopy(newlsDomain);
428 levelSets.back()->finalize(finalWidth);
429 }
430
435 template <class DiscretizationSchemeType>
436 double integrateTime(DiscretizationSchemeType spatialScheme,
437 double maxTimeStep) {
438
439 auto &topDomain = levelSets.back()->getDomain();
440 auto &grid = levelSets.back()->getGrid();
441
442 typename PointData<T>::ScalarDataType *voidMarkerPointer;
443 if (ignoreVoids) {
444 MarkVoidPoints<T, D>(levelSets.back()).apply();
445 auto &pointData = levelSets.back()->getPointData();
446 voidMarkerPointer =
448 if (voidMarkerPointer == nullptr) {
449 VIENNACORE_LOG_WARNING("Advect: Cannot find void point markers. Not "
450 "ignoring void points.");
451 ignoreVoids = false;
452 }
453 }
454 const bool ignoreVoidPoints = ignoreVoids;
455 const bool useAdaptiveTimeStepping = adaptiveTimeStepping;
456 const auto adaptiveFactor = 1.0 / adaptiveTimeStepSubdivisions;
457
458 if (!storedRates.empty()) {
459 VIENNACORE_LOG_WARNING("Advect: Overwriting previously stored rates.");
460 }
461
462 storedRates.resize(topDomain.getNumberOfSegments());
463
464#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
465 {
466 int p = 0;
467#ifdef _OPENMP
468 p = omp_get_thread_num();
469#endif
470 viennahrle::Index<D> startVector =
471 (p == 0) ? grid.getMinGridPoint()
472 : topDomain.getSegmentation()[p - 1];
473
474 viennahrle::Index<D> endVector =
475 (p != static_cast<int>(topDomain.getNumberOfSegments() - 1))
476 ? topDomain.getSegmentation()[p]
477 : grid.incrementIndices(grid.getMaxGridPoint());
478
479 double tempMaxTimeStep = maxTimeStep;
480 // store the rates and value of underneath LS for this segment
481 auto &rates = storedRates[p];
482 rates.reserve(
483 topDomain.getNumberOfPoints() /
484 static_cast<double>((levelSets.back())->getNumberOfSegments()) +
485 10);
486
487 // an iterator for each level set
488 std::vector<ConstSparseIterator> iterators;
489 for (auto const &ls : levelSets) {
490 iterators.emplace_back(ls->getDomain());
491 }
492
493 DiscretizationSchemeType scheme(spatialScheme);
494
495 for (ConstSparseIterator it(topDomain, startVector);
496 it.getStartIndices() < endVector; it.next()) {
497
498 if (!it.isDefined() || std::abs(it.getValue()) > integrationCutoff)
499 continue;
500
501 T value = it.getValue();
502 double maxStepTime = 0;
503 double cfl = timeStepRatio;
504
505 for (int currentLevelSetId = levelSets.size() - 1;
506 currentLevelSetId >= 0; --currentLevelSetId) {
507
508 std::pair<T, T> gradNDissipation;
509
510 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
511 // check if there is any other levelset at the same point:
512 // if yes, take the velocity of the lowest levelset
513 for (unsigned lowerLevelSetId = 0;
514 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
515 // put iterator to same position as the top levelset
516 iterators[lowerLevelSetId].goToIndicesSequential(
517 it.getStartIndices());
518
519 // if the lower surface is actually outside, i.e. its LS value
520 // is lower or equal
521 if (iterators[lowerLevelSetId].getValue() <=
522 value + wrappingLayerEpsilon) {
523 gradNDissipation =
524 scheme(it.getStartIndices(), lowerLevelSetId);
525 break;
526 }
527 }
528 }
529
530 T velocity = gradNDissipation.first - gradNDissipation.second;
531 if (velocity > 0.) {
532 // Case 1: Growth / Deposition (Velocity > 0)
533 // Limit the time step based on the standard CFL condition.
534 maxStepTime += cfl / velocity;
535 rates.emplace_back(gradNDissipation,
536 -std::numeric_limits<T>::max());
537 break;
538 } else if (velocity == 0.) {
539 // Case 2: Static (Velocity == 0)
540 // No time step limit imposed by this point.
541 maxStepTime = std::numeric_limits<T>::max();
542 rates.emplace_back(gradNDissipation, std::numeric_limits<T>::max());
543 break;
544 } else {
545 // Case 3: Etching (Velocity < 0)
546 // Retrieve the interface location of the underlying material.
547 T valueBelow;
548 if (currentLevelSetId > 0) {
549 iterators[currentLevelSetId - 1].goToIndicesSequential(
550 it.getStartIndices());
551 valueBelow = iterators[currentLevelSetId - 1].getValue();
552 } else {
553 valueBelow = std::numeric_limits<T>::max();
554 }
555 // Calculate the top material thickness
556 T difference = std::abs(valueBelow - value);
557
558 if (difference >= cfl) {
559 // Sub-case 3a: Standard Advection
560 // Far from interface: Use full CFL time step.
561 maxStepTime -= cfl / velocity;
562 rates.emplace_back(gradNDissipation,
563 std::numeric_limits<T>::max());
564 break;
565 } else {
566 // Sub-case 3b: Interface Interaction
567 // Use adaptiveFactor as threshold.
568 if (useAdaptiveTimeStepping &&
569 difference > adaptiveFactor * cfl) {
570 // Adaptive Sub-stepping:
571 // Approaching boundary: Force small steps to gather
572 // flux statistics and prevent numerical overshoot ("Soft
573 // Landing").
574 maxStepTime -= adaptiveFactor * cfl / velocity;
575 rates.emplace_back(gradNDissipation,
576 std::numeric_limits<T>::max());
577 break;
578 } else {
579 // Terminal Step:
580 // Within tolerance: Snap to boundary, consume budget, and
581 // switch material.
582 cfl -= difference;
583 value = valueBelow;
584 maxStepTime -= difference / velocity;
585 rates.emplace_back(gradNDissipation, valueBelow);
586 }
587 }
588 }
589 }
590
591 if (maxStepTime < tempMaxTimeStep)
592 tempMaxTimeStep = maxStepTime;
593 }
594
595#pragma omp critical
596 {
597 // If a Lax Friedrichs scheme is selected the time step is
598 // reduced depending on the dissipation coefficients
599 // For Engquist Osher scheme this function is empty.
600 scheme.reduceTimeStepHamiltonJacobi(
601 tempMaxTimeStep, levelSets.back()->getGrid().getGridDelta());
602
603 // set global timestep maximum
604 if (tempMaxTimeStep < maxTimeStep)
605 maxTimeStep = tempMaxTimeStep;
606 }
607 } // end of parallel section
608
609 // maxTimeStep is now the maximum time step possible for all points
610 // and rates are stored in a vector
611 return maxTimeStep;
612 }
613
616 void computeRates(double maxTimeStep = std::numeric_limits<double>::max()) {
617 prepareLS();
618 if (spatialScheme == SpatialSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
619 auto is = lsInternal::EngquistOsher<T, D, 1>(levelSets.back(), velocities,
620 calculateNormalVectors);
621 currentTimeStep = integrateTime(is, maxTimeStep);
622 } else if (spatialScheme == SpatialSchemeEnum::ENGQUIST_OSHER_2ND_ORDER) {
623 auto is = lsInternal::EngquistOsher<T, D, 2>(levelSets.back(), velocities,
624 calculateNormalVectors);
625 currentTimeStep = integrateTime(is, maxTimeStep);
626 } else if (spatialScheme == SpatialSchemeEnum::LAX_FRIEDRICHS_1ST_ORDER) {
627 auto alphas = findGlobalAlphas();
628 auto is = lsInternal::LaxFriedrichs<T, D, 1>(levelSets.back(), velocities,
629 dissipationAlpha, alphas,
630 calculateNormalVectors);
631 currentTimeStep = integrateTime(is, maxTimeStep);
632 } else if (spatialScheme == SpatialSchemeEnum::LAX_FRIEDRICHS_2ND_ORDER) {
633 auto alphas = findGlobalAlphas();
634 auto is = lsInternal::LaxFriedrichs<T, D, 2>(levelSets.back(), velocities,
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);
668 } else if (spatialScheme == SpatialSchemeEnum::WENO_3RD_ORDER) {
669 // Instantiate WENO with order 3
670 auto is = lsInternal::WENO<T, D, 3>(levelSets.back(), velocities,
671 dissipationAlpha);
672 currentTimeStep = integrateTime(is, maxTimeStep);
673 } else if (spatialScheme == SpatialSchemeEnum::WENO_5TH_ORDER) {
674 // Instantiate WENO with order 5
675 auto is = lsInternal::WENO<T, D, 5>(levelSets.back(), velocities,
676 dissipationAlpha);
677 currentTimeStep = integrateTime(is, maxTimeStep);
678 } else {
679 VIENNACORE_LOG_ERROR("Advect: Discretization scheme not found.");
680 currentTimeStep = -1.;
681 }
682 }
683
684 // Level Sets below are also considered in order to adjust the advection
685 // depth accordingly if there would be a material change.
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!");
691 }
692
693 auto &topDomain = levelSets.back()->getDomain();
694
695 assert(dt >= 0. && "No time step set!");
696 assert(storedRates.size() == topDomain.getNumberOfSegments());
697
698 // reduce to one layer thickness and apply new values directly to the
699 // domain segments --> DO NOT CHANGE SEGMENTATION HERE (true parameter)
700 Reduce<T, D>(levelSets.back(), 1, true).apply();
701
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());
707
708 const bool checkDiss = checkDissipation;
709
710#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
711 {
712 int p = 0;
713#ifdef _OPENMP
714 p = omp_get_thread_num();
715#endif
716 auto itRS = storedRates[p].cbegin();
717 auto &segment = topDomain.getDomainSegment(p);
718 const unsigned maxId = segment.getNumberOfPoints();
719
720 if (saveVelocities) {
721 velocityVectors[p].resize(maxId);
722 dissipationVectors[p].resize(maxId);
723 }
724
725 for (unsigned localId = 0; localId < maxId; ++localId) {
726 T &value = segment.definedValues[localId];
727
728 // Skip points that were not part of computeRates (outer layers)
729 if (std::abs(value) > integrationCutoff)
730 continue;
731
732 double time = dt;
733
734 // if there is a change in materials during one time step, deduct
735 // the time taken to advect up to the end of the top material and
736 // set the LS value to the one below
737 auto const [gradient, dissipation] = itRS->first;
738 T velocity = gradient - dissipation;
739 // check if dissipation is too high and would cause a change in
740 // direction of the velocity
741 if (checkDiss && (gradient < 0 && velocity > 0) ||
742 (gradient > 0 && velocity < 0)) {
743 velocity = 0;
744 }
745
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;
750 ++itRS; // advance the TempStopRates iterator by one
751
752 // recalculate velocity and rate
753 velocity = itRS->first.first - itRS->first.second;
754 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
755 (itRS->first.first > 0 && velocity < 0)) {
756 velocity = 0;
757 }
758 rate = time * velocity;
759 }
760
761 // now deduct the velocity times the time step we take
762 value -= rate;
763
764 if (saveVelocities) {
765 velocityVectors[p][localId] = rate;
766 dissipationVectors[p][localId] = itRS->first.second;
767 }
768
769 // this is run when two materials are close but the velocity is too slow
770 // to actually reach the second material, to get rid of the extra
771 // entry in the TempRatesStop
772 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
773 ++itRS;
774
775 // advance the TempStopRates iterator by one
776 ++itRS;
777 }
778 } // end of parallel section
779
780 if (saveVelocities) {
781 auto &pointData = levelSets.back()->getPointData();
782
783 typename PointData<T>::ScalarDataType vels;
784 typename PointData<T>::ScalarDataType diss;
785
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()));
793 }
794 pointData.insertReplaceScalarData(std::move(vels), velocityLabel);
795 pointData.insertReplaceScalarData(std::move(diss), dissipationLabel);
796 }
797
798 // clear the stored rates since surface has changed
799 storedRates.clear();
800 }
801
802 void adjustLowerLayers() {
803 // Adjust all level sets below the advected one
804 if (spatialScheme !=
806 for (unsigned i = 0; i < levelSets.size() - 1; ++i) {
808 levelSets[i], levelSets.back(),
810 .apply();
811 }
812 }
813 }
814
817 double advect(double maxTimeStep) {
818 switch (temporalScheme) {
821 *this, maxTimeStep);
824 *this, maxTimeStep);
826 default:
828 *this, maxTimeStep);
829 }
830 }
831
832public:
833 static constexpr char velocityLabel[] = "AdvectionVelocities";
834 static constexpr char dissipationLabel[] = "Dissipation";
835
836 Advect() = default;
837
838 explicit Advect(SmartPointer<Domain<T, D>> passedlsDomain) {
839 levelSets.push_back(passedlsDomain);
840 }
841
842 Advect(SmartPointer<Domain<T, D>> passedlsDomain,
843 SmartPointer<VelocityField<T>> passedVelocities) {
844 levelSets.push_back(passedlsDomain);
845 velocities = passedVelocities;
846 }
847
848 Advect(std::vector<SmartPointer<Domain<T, D>>> passedlsDomains,
849 SmartPointer<VelocityField<T>> passedVelocities)
850 : levelSets(passedlsDomains) {
851 velocities = passedVelocities;
852 }
853
856 void insertNextLevelSet(SmartPointer<Domain<T, D>> passedlsDomain) {
857 levelSets.push_back(passedlsDomain);
858 }
859
860 void clearLevelSets() { levelSets.clear(); }
861
864 void setVelocityField(SmartPointer<VelocityField<T>> passedVelocities) {
865 velocities = passedVelocities;
866 }
867
873 void setAdvectionTime(double time) { advectionTime = time; }
874
879 void setSingleStep(bool singleStep) { performOnlySingleStep = singleStep; }
880
885 void setTimeStepRatio(const double &cfl) { timeStepRatio = cfl; }
886
891 void setCalculateNormalVectors(bool cnv) { calculateNormalVectors = cnv; }
892
899 void setIgnoreVoids(bool iV) { ignoreVoids = iV; }
900
904 void setAdaptiveTimeStepping(bool aTS = true, unsigned subdivisions = 20) {
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.");
909 subdivisions = 1;
910 }
911 adaptiveTimeStepSubdivisions = subdivisions;
912 }
913
916 void setSaveAdvectionVelocities(bool sAV) { saveAdvectionVelocities = sAV; }
917
920 double getAdvectedTime() const { return advectedTime; }
921
923 double getCurrentTimeStep() const { return currentTimeStep; }
924
926 unsigned getNumberOfTimeSteps() const { return numberOfTimeSteps; }
927
929 double getTimeStepRatio() const { return timeStepRatio; }
930
932 bool getCalculateNormalVectors() const { return calculateNormalVectors; }
933
936 void setSpatialScheme(SpatialSchemeEnum scheme) { spatialScheme = scheme; }
937
938 // Deprecated and will be removed in future versions:
939 // use setSpatialScheme instead
940 [[deprecated("Use setSpatialScheme instead")]] void
941 setIntegrationScheme(IntegrationSchemeEnum scheme) {
942 VIENNACORE_LOG_WARNING(
943 "Advect::setIntegrationScheme is deprecated and will be removed in "
944 "future versions. Use setSpatialScheme instead.");
945 spatialScheme = scheme;
946 }
947
949 void setTemporalScheme(TemporalSchemeEnum scheme) { temporalScheme = scheme; }
950
955 void setDissipationAlpha(const double &a) { dissipationAlpha = a; }
956
957 // Sets the velocity to 0 if the dissipation is too high
958 void setCheckDissipation(bool check) { checkDissipation = check; }
959
962 void setUpdatePointData(bool update) { updatePointData = update; }
963
967 std::function<bool(SmartPointer<Domain<T, D>>)> callback) {
968 velocityUpdateCallback = callback;
969 }
970
971 // Prepare the levelset for advection, based on the provided spatial
972 // discretization scheme.
973 void prepareLS() {
974 // check whether a level set and velocities have been given
975 if (levelSets.empty()) {
976 VIENNACORE_LOG_ERROR("No level sets passed to Advect.");
977 return;
978 }
979
980 if (spatialScheme == SpatialSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
982 } else if (spatialScheme == SpatialSchemeEnum::ENGQUIST_OSHER_2ND_ORDER) {
984 } else if (spatialScheme == SpatialSchemeEnum::LAX_FRIEDRICHS_1ST_ORDER) {
986 } else if (spatialScheme == SpatialSchemeEnum::LAX_FRIEDRICHS_2ND_ORDER) {
988 } else if (spatialScheme ==
991 levelSets.back());
992 } else if (spatialScheme ==
995 } else if (spatialScheme ==
998 } else if (spatialScheme ==
1001 } else if (spatialScheme ==
1004 } else if (spatialScheme ==
1007 levelSets.back());
1008 } else if (spatialScheme == SpatialSchemeEnum::WENO_3RD_ORDER) {
1009 lsInternal::WENO<T, D, 3>::prepareLS(levelSets.back());
1010 } else if (spatialScheme == SpatialSchemeEnum::WENO_5TH_ORDER) {
1011 lsInternal::WENO<T, D, 5>::prepareLS(levelSets.back());
1012 } else {
1013 VIENNACORE_LOG_ERROR("Advect: Discretization scheme not found.");
1014 }
1015 }
1016
1017 void apply() {
1018 // check whether a level set and velocities have been given
1019 if (levelSets.empty()) {
1020 VIENNACORE_LOG_ERROR("No level sets passed to Advect. Not advecting.");
1021 return;
1022 }
1023 if (velocities == nullptr) {
1024 VIENNACORE_LOG_ERROR(
1025 "No velocity field passed to Advect. Not advecting.");
1026 return;
1027 }
1028
1029 if (advectionTime == 0.) {
1030 advectedTime = advect(std::numeric_limits<double>::max());
1031 numberOfTimeSteps = 1;
1032 } else {
1033 double currentTime = 0.0;
1034 numberOfTimeSteps = 0;
1035 while (currentTime < advectionTime) {
1036 currentTime += advect(advectionTime - currentTime);
1037 ++numberOfTimeSteps;
1038 if (performOnlySingleStep)
1039 break;
1040 }
1041 advectedTime = currentTime;
1042 }
1043 }
1044};
1045
1046// add all template specializations for this class
1048
1049} // namespace viennals
constexpr int D
Definition Epitaxy.cpp:12
double T
Definition Epitaxy.cpp:13
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:29
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:33
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:49
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:51
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:30
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:223
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:864
static constexpr char velocityLabel[]
Definition lsAdvect.hpp:833
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:856
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:966
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:962
void setSpatialScheme(SpatialSchemeEnum scheme)
Set which spatial discretization scheme should be used out of the ones specified in SpatialSchemeEnum...
Definition lsAdvect.hpp:936
bool getCalculateNormalVectors() const
Get whether normal vectors were calculated.
Definition lsAdvect.hpp:932
void setIgnoreVoids(bool iV)
Set whether level set values, which are not part of the "top" geometrically connected part of values,...
Definition lsAdvect.hpp:899
Advect()=default
void clearLevelSets()
Definition lsAdvect.hpp:860
void setCalculateNormalVectors(bool cnv)
Set whether normal vectors should be calculated at each level set point. Defaults to true....
Definition lsAdvect.hpp:891
void apply()
Definition lsAdvect.hpp:1017
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:879
void prepareLS()
Definition lsAdvect.hpp:973
double getAdvectedTime() const
Get by how much the physical time was advanced during the last apply() call.
Definition lsAdvect.hpp:920
void setTemporalScheme(TemporalSchemeEnum scheme)
Set which time integration scheme should be used.
Definition lsAdvect.hpp:949
void setDissipationAlpha(const double &a)
Set the alpha dissipation coefficient. For lsLaxFriedrichs, this is used as the alpha value....
Definition lsAdvect.hpp:955
void setCheckDissipation(bool check)
Definition lsAdvect.hpp:958
Advect(std::vector< SmartPointer< Domain< T, D > > > passedlsDomains, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:848
void setIntegrationScheme(IntegrationSchemeEnum scheme)
Definition lsAdvect.hpp:941
double getCurrentTimeStep() const
Return the last applied time step.
Definition lsAdvect.hpp:923
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:873
unsigned getNumberOfTimeSteps() const
Get how many advection steps were performed during the last apply() call.
Definition lsAdvect.hpp:926
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:904
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:885
Advect(SmartPointer< Domain< T, D > > passedlsDomain)
Definition lsAdvect.hpp:838
Advect(SmartPointer< Domain< T, D > > passedlsDomain, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:842
double getTimeStepRatio() const
Get the value of the CFL number.
Definition lsAdvect.hpp:929
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:916
static constexpr char dissipationLabel[]
Definition lsAdvect.hpp:834
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