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