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// Integration schemes
20#include <lsEngquistOsher.hpp>
21#include <lsLaxFriedrichs.hpp>
26#include <lsWENO5.hpp>
27
28// Velocity accessor
29#include <lsVelocityField.hpp>
30
31// #define DEBUG_LS_ADVECT_HPP
32#ifdef DEBUG_LS_ADVECT_HPP
33#include <lsToMesh.hpp>
34#include <lsVTKWriter.hpp>
35#endif
36
37namespace viennals {
38
39using namespace viennacore;
40
56
66template <class T, int D> class Advect {
67 using ConstSparseIterator =
68 viennahrle::ConstSparseIterator<typename Domain<T, D>::DomainType>;
69 using hrleIndexType = viennahrle::IndexType;
70
71 std::vector<SmartPointer<Domain<T, D>>> levelSets;
72 SmartPointer<VelocityField<T>> velocities = nullptr;
73 IntegrationSchemeEnum integrationScheme =
75 double timeStepRatio = 0.4999;
76 double dissipationAlpha = 1.0;
77 bool calculateNormalVectors = true;
78 bool ignoreVoids = false;
79 double advectionTime = 0.;
80 bool performOnlySingleStep = false;
81 double advectedTime = 0.;
82 unsigned numberOfTimeSteps = 0;
83 bool saveAdvectionVelocities = false;
84 bool updatePointData = true;
85 bool checkDissipation = true;
86 bool adaptiveTimeStepping = false;
87 double adaptiveTimeStepThreshold = 0.05;
88 static constexpr double wrappingLayerEpsilon = 1e-4;
89
90 // this vector will hold the maximum time step for each point and the
91 // corresponding velocity
92 std::vector<std::vector<std::pair<std::pair<T, T>, T>>> storedRates;
93 double currentTimeStep = -1.;
94
95 VectorType<T, 3> findGlobalAlphas() const {
96
97 auto &topDomain = levelSets.back()->getDomain();
98 auto &grid = levelSets.back()->getGrid();
99
100 const T gridDelta = grid.getGridDelta();
101 const T deltaPos = gridDelta;
102 const T deltaNeg = -gridDelta;
103
104 VectorType<T, 3> finalAlphas = {0., 0., 0.};
105
106#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
107 {
108 VectorType<T, 3> localAlphas = {0., 0., 0.};
109 int p = 0;
110#ifdef _OPENMP
111 p = omp_get_thread_num();
112#endif
113 viennahrle::Index<D> startVector =
114 (p == 0) ? grid.getMinGridPoint()
115 : topDomain.getSegmentation()[p - 1];
116 viennahrle::Index<D> endVector =
117 (p != static_cast<int>(topDomain.getNumberOfSegments() - 1))
118 ? topDomain.getSegmentation()[p]
119 : grid.incrementIndices(grid.getMaxGridPoint());
120
121 // an iterator for each level set
122 std::vector<ConstSparseIterator> iterators;
123 for (auto const &ls : levelSets) {
124 iterators.emplace_back(ls->getDomain());
125 }
126
127 // neighborIterator for the top level set
128 viennahrle::ConstSparseStarIterator<typename Domain<T, D>::DomainType, 1>
129 neighborIterator(topDomain);
130
131 for (ConstSparseIterator it(topDomain, startVector);
132 it.getStartIndices() < endVector; ++it) {
133
134 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
135 continue;
136
137 const T value = it.getValue();
138 const auto indices = it.getStartIndices();
139
140 // check if there is any other levelset at the same point:
141 // if yes, take the velocity of the lowest levelset
142 for (unsigned lowerLevelSetId = 0; lowerLevelSetId < levelSets.size();
143 ++lowerLevelSetId) {
144 // put iterator to same position as the top levelset
145 iterators[lowerLevelSetId].goToIndicesSequential(indices);
146
147 // if the lower surface is actually outside, i.e. its LS value
148 // is lower or equal
149 if (iterators[lowerLevelSetId].getValue() <=
150 value + wrappingLayerEpsilon) {
151
152 // move neighborIterator to current position
153 neighborIterator.goToIndicesSequential(indices);
154
155 Vec3D<T> coords;
156 for (unsigned i = 0; i < D; ++i) {
157 coords[i] = indices[i] * gridDelta;
158 }
159
160 Vec3D<T> normal = {};
161 T normalModulus = 0.;
162 for (unsigned i = 0; i < D; ++i) {
163 const T phiPos = neighborIterator.getNeighbor(i).getValue();
164 const T phiNeg = neighborIterator.getNeighbor(i + D).getValue();
165
166 T diffPos = (phiPos - value) / deltaPos;
167 T diffNeg = (phiNeg - value) / deltaNeg;
168
169 normal[i] = (diffNeg + diffPos) * 0.5;
170 normalModulus += normal[i] * normal[i];
171 }
172 normalModulus = std::sqrt(normalModulus);
173 for (unsigned i = 0; i < D; ++i)
174 normal[i] /= normalModulus;
175
176 T scaVel = velocities->getScalarVelocity(
177 coords, lowerLevelSetId, normal,
178 neighborIterator.getCenter().getPointId());
179 auto vecVel = velocities->getVectorVelocity(
180 coords, lowerLevelSetId, normal,
181 neighborIterator.getCenter().getPointId());
182
183 for (unsigned i = 0; i < D; ++i) {
184 T tempAlpha = std::abs((scaVel + vecVel[i]) * normal[i]);
185 localAlphas[i] = std::max(localAlphas[i], tempAlpha);
186 }
187
188 // exit material loop
189 break;
190 }
191 }
192 }
193
194#pragma omp critical
195 {
196 for (unsigned i = 0; i < D; ++i) {
197 finalAlphas[i] = std::max(finalAlphas[i], localAlphas[i]);
198 }
199 }
200 } // end of parallel section
201
202 return finalAlphas;
203 }
204
205 void rebuildLS() {
206 // TODO: this function uses Manhattan distances for renormalisation,
207 // since this is the quickest. For visualisation applications, better
208 // renormalisation is needed, so it might be good to implement
209 // Euler distance renormalisation as an option
210 auto &grid = levelSets.back()->getGrid();
211 auto newlsDomain = SmartPointer<Domain<T, D>>::New(grid);
212 auto &newDomain = newlsDomain->getDomain();
213 auto &domain = levelSets.back()->getDomain();
214
215 newDomain.initialize(domain.getNewSegmentation(),
216 domain.getAllocation() *
217 (2.0 / levelSets.back()->getLevelSetWidth()));
218
219 const bool updateData = updatePointData;
220 // save how data should be transferred to new level set
221 // list of indices into the old pointData vector
222 std::vector<std::vector<unsigned>> newDataSourceIds;
223 if (updateData)
224 newDataSourceIds.resize(newDomain.getNumberOfSegments());
225
226#ifdef DEBUG_LS_ADVECT_HPP
227 {
228 auto mesh = SmartPointer<Mesh<T>>::New();
229 ToMesh<T, D>(levelSets.back(), mesh).apply();
230 VTKWriter<T>(mesh, "Advect_beforeRebuild.vtk").apply();
231 }
232#endif
233
234#pragma omp parallel num_threads(newDomain.getNumberOfSegments())
235 {
236 int p = 0;
237#ifdef _OPENMP
238 p = omp_get_thread_num();
239#endif
240 auto &domainSegment = newDomain.getDomainSegment(p);
241
242 viennahrle::Index<D> startVector =
243 (p == 0) ? grid.getMinGridPoint()
244 : newDomain.getSegmentation()[p - 1];
245
246 viennahrle::Index<D> endVector =
247 (p != static_cast<int>(newDomain.getNumberOfSegments() - 1))
248 ? newDomain.getSegmentation()[p]
249 : grid.incrementIndices(grid.getMaxGridPoint());
250
251 // reserve a bit more to avoid reallocation
252 // would expect number of points to roughly double
253 if (updateData)
254 newDataSourceIds[p].reserve(2.5 * domainSegment.getNumberOfPoints());
255
256 for (viennahrle::SparseStarIterator<typename Domain<T, D>::DomainType, 1>
257 it(domain, startVector);
258 it.getIndices() < endVector; ++it) {
259
260 // if the center is an active grid point
261 // <1.0 since it could have been change by 0.5 max
262 if (std::abs(it.getCenter().getValue()) <= 1.0) {
263
264 int k = 0;
265 for (; k < 2 * D; k++)
266 if (std::signbit(it.getNeighbor(k).getValue() - 1e-7) !=
267 std::signbit(it.getCenter().getValue() + 1e-7))
268 break;
269
270 // if there is at least one neighbor of opposite sign
271 if (k != 2 * D) {
272 if (it.getCenter().getDefinedValue() > 0.5) {
273 int j = 0;
274 for (; j < 2 * D; j++) {
275 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
276 if (it.getNeighbor(j).getDefinedValue() < -0.5)
277 break;
278 }
279 if (j == 2 * D) {
280 domainSegment.insertNextDefinedPoint(
281 it.getIndices(), it.getCenter().getDefinedValue());
282 if (updateData)
283 newDataSourceIds[p].push_back(it.getCenter().getPointId());
284 // if there is at least one active grid point, which is < -0.5
285 } else {
286 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
287 if (updateData)
288 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
289 }
290 } else if (it.getCenter().getDefinedValue() < -0.5) {
291 int j = 0;
292 for (; j < 2 * D; j++) {
293 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
294 if (it.getNeighbor(j).getDefinedValue() > 0.5)
295 break;
296 }
297
298 if (j == 2 * D) {
299 domainSegment.insertNextDefinedPoint(
300 it.getIndices(), it.getCenter().getDefinedValue());
301 if (updateData)
302 newDataSourceIds[p].push_back(it.getCenter().getPointId());
303 // if there is at least one active grid point, which is > 0.5
304 } else {
305 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
306 if (updateData)
307 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
308 }
309 } else {
310 domainSegment.insertNextDefinedPoint(
311 it.getIndices(), it.getCenter().getDefinedValue());
312 if (updateData)
313 newDataSourceIds[p].push_back(it.getCenter().getPointId());
314 }
315 } else {
316 domainSegment.insertNextUndefinedPoint(
317 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
320 }
321
322 } else { // if the center is not an active grid point
323 if (it.getCenter().getValue() >= 0) {
324 int usedNeighbor = -1;
325 T distance = Domain<T, D>::POS_VALUE;
326 for (int i = 0; i < 2 * D; i++) {
327 T value = it.getNeighbor(i).getValue();
328 if (std::abs(value) <= 1.0 && (value < 0.)) {
329 if (distance > value + 1.0) {
330 distance = value + 1.0;
331 usedNeighbor = i;
332 }
333 }
334 }
335
336 if (distance <= 1.) {
337 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
338 if (updateData)
339 newDataSourceIds[p].push_back(
340 it.getNeighbor(usedNeighbor).getPointId());
341 } else {
342 domainSegment.insertNextUndefinedPoint(it.getIndices(),
344 }
345
346 } else {
347 int usedNeighbor = -1;
348 T distance = Domain<T, D>::NEG_VALUE;
349 for (int i = 0; i < 2 * D; i++) {
350 T value = it.getNeighbor(i).getValue();
351 if (std::abs(value) <= 1.0 && (value > 0)) {
352 if (distance < value - 1.0) {
353 // distance = std::max(distance, value - T(1.0));
354 distance = value - 1.0;
355 usedNeighbor = i;
356 }
357 }
358 }
359
360 if (distance >= -1.) {
361 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
362 if (updateData)
363 newDataSourceIds[p].push_back(
364 it.getNeighbor(usedNeighbor).getPointId());
365 } else {
366 domainSegment.insertNextUndefinedPoint(it.getIndices(),
368 }
369 }
370 }
371 }
372 }
373
374 // now copy old data into new level set
375 if (updateData) {
376 auto &pointData = levelSets.back()->getPointData();
377 newlsDomain->getPointData().translateFromMultiData(pointData,
378 newDataSourceIds);
379 }
380
381 newDomain.finalize();
382 newDomain.segment();
383 levelSets.back()->deepCopy(newlsDomain);
384 levelSets.back()->finalize(2);
385 }
386
389 double advect(double maxTimeStep) {
390 // check whether a level set and velocities have been given
391 if (levelSets.empty()) {
392 VIENNACORE_LOG_ERROR("No level sets passed to Advect. Not advecting.");
393 return std::numeric_limits<double>::max();
394 }
395 if (velocities == nullptr) {
396 VIENNACORE_LOG_ERROR(
397 "No velocity field passed to Advect. Not advecting.");
398 return std::numeric_limits<double>::max();
399 }
400
401 if (currentTimeStep < 0. || storedRates.empty())
402 applyIntegration(maxTimeStep);
403
404 moveSurface();
405
406 rebuildLS();
407
408 // Adjust all level sets below the advected one
409 // This means, that when the top level-set and one below
410 // are etched, the lower one is moved with the top level-set
411 // TODO: Adjust lower layers also when they have grown,
412 // to allow for two different growth rates of materials
413 if (integrationScheme !=
415 for (unsigned i = 0; i < levelSets.size() - 1; ++i) {
416 BooleanOperation<T, D>(levelSets[i], levelSets.back(),
418 .apply();
419 }
420 }
421
422 return getCurrentTimeStep();
423 }
424
428 template <class IntegrationSchemeType>
429 double integrateTime(IntegrationSchemeType IntegrationScheme,
430 double maxTimeStep) {
431
432 auto &topDomain = levelSets.back()->getDomain();
433 auto &grid = levelSets.back()->getGrid();
434
435 typename PointData<T>::ScalarDataType *voidMarkerPointer;
436 if (ignoreVoids) {
437 MarkVoidPoints<T, D>(levelSets.back()).apply();
438 auto &pointData = levelSets.back()->getPointData();
439 voidMarkerPointer =
441 if (voidMarkerPointer == nullptr) {
442 VIENNACORE_LOG_WARNING("Advect: Cannot find void point markers. Not "
443 "ignoring void points.");
444 ignoreVoids = false;
445 }
446 }
447 const bool ignoreVoidPoints = ignoreVoids;
448 const bool useAdaptiveTimeStepping = adaptiveTimeStepping;
449 const double atsThreshold = adaptiveTimeStepThreshold;
450
451 if (!storedRates.empty()) {
452 VIENNACORE_LOG_WARNING("Advect: Overwriting previously stored rates.");
453 }
454
455 storedRates.resize(topDomain.getNumberOfSegments());
456
457#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
458 {
459 int p = 0;
460#ifdef _OPENMP
461 p = omp_get_thread_num();
462#endif
463 viennahrle::Index<D> startVector =
464 (p == 0) ? grid.getMinGridPoint()
465 : topDomain.getSegmentation()[p - 1];
466
467 viennahrle::Index<D> endVector =
468 (p != static_cast<int>(topDomain.getNumberOfSegments() - 1))
469 ? topDomain.getSegmentation()[p]
470 : grid.incrementIndices(grid.getMaxGridPoint());
471
472 double tempMaxTimeStep = maxTimeStep;
473 auto &tempRates = storedRates[p];
474 tempRates.reserve(
475 topDomain.getNumberOfPoints() /
476 static_cast<double>((levelSets.back())->getNumberOfSegments()) +
477 10);
478
479 // an iterator for each level set
480 std::vector<ConstSparseIterator> iterators;
481 for (auto const &ls : levelSets) {
482 iterators.emplace_back(ls->getDomain());
483 }
484
485 IntegrationSchemeType scheme(IntegrationScheme);
486
487 for (ConstSparseIterator it(topDomain, startVector);
488 it.getStartIndices() < endVector; ++it) {
489
490 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
491 continue;
492
493 T value = it.getValue();
494 double maxStepTime = 0;
495 double cfl = timeStepRatio;
496
497 for (int currentLevelSetId = levelSets.size() - 1;
498 currentLevelSetId >= 0; --currentLevelSetId) {
499
500 std::pair<T, T> gradNDissipation;
501
502 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
503 // check if there is any other levelset at the same point:
504 // if yes, take the velocity of the lowest levelset
505 for (unsigned lowerLevelSetId = 0;
506 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
507 // put iterator to same position as the top levelset
508 iterators[lowerLevelSetId].goToIndicesSequential(
509 it.getStartIndices());
510
511 // if the lower surface is actually outside, i.e. its LS value
512 // is lower or equal
513 if (iterators[lowerLevelSetId].getValue() <=
514 value + wrappingLayerEpsilon) {
515 gradNDissipation =
516 scheme(it.getStartIndices(), lowerLevelSetId);
517 break;
518 }
519 }
520 }
521
522 T velocity = gradNDissipation.first - gradNDissipation.second;
523 if (velocity > 0.) {
524 // Case 1: Growth / Deposition (Velocity > 0)
525 // Limit the time step based on the standard CFL condition.
526 maxStepTime += cfl / velocity;
527 tempRates.push_back(std::make_pair(gradNDissipation,
528 -std::numeric_limits<T>::max()));
529 break;
530 } else if (velocity == 0.) {
531 // Case 2: Static (Velocity == 0)
532 // No time step limit imposed by this point.
533 maxStepTime = std::numeric_limits<T>::max();
534 tempRates.push_back(std::make_pair(gradNDissipation,
535 std::numeric_limits<T>::max()));
536 break;
537 } else {
538 // Case 3: Etching (Velocity < 0)
539 // Retrieve the interface location of the underlying material.
540 T valueBelow;
541 if (currentLevelSetId > 0) {
542 iterators[currentLevelSetId - 1].goToIndicesSequential(
543 it.getStartIndices());
544 valueBelow = iterators[currentLevelSetId - 1].getValue();
545 } else {
546 valueBelow = std::numeric_limits<T>::max();
547 }
548 // Calculate the top material thickness
549 T difference = std::abs(valueBelow - value);
550
551 if (difference >= cfl) {
552 // Sub-case 3a: Standard Advection
553 // Far from interface: Use full CFL time step.
554 maxStepTime -= cfl / velocity;
555 tempRates.push_back(std::make_pair(
556 gradNDissipation, std::numeric_limits<T>::max()));
557 break;
558
559 } else {
560 // Sub-case 3b: Interface Interaction
561 if (useAdaptiveTimeStepping && difference > atsThreshold * cfl) {
562 // Adaptive Sub-stepping:
563 // Approaching boundary: Force small steps to gather flux
564 // statistics and prevent numerical overshoot ("Soft Landing").
565 maxStepTime -= atsThreshold * cfl / velocity;
566 tempRates.push_back(std::make_pair(
567 gradNDissipation, std::numeric_limits<T>::min()));
568 } else {
569 // Terminal Step:
570 // Within tolerance: Snap to boundary, consume budget, and
571 // switch material.
572 tempRates.push_back(
573 std::make_pair(gradNDissipation, valueBelow));
574 cfl -= difference;
575 value = valueBelow;
576 maxStepTime -= difference / velocity;
577 }
578 }
579 }
580 }
581
582 if (maxStepTime < tempMaxTimeStep)
583 tempMaxTimeStep = maxStepTime;
584 }
585
586#pragma omp critical
587 {
588 // If a Lax Friedrichs scheme is selected the time step is
589 // reduced depending on the dissipation coefficients
590 // For Engquist Osher scheme this function is empty.
591 scheme.reduceTimeStepHamiltonJacobi(
592 tempMaxTimeStep, levelSets.back()->getGrid().getGridDelta());
593
594 // set global timestep maximum
595 if (tempMaxTimeStep < maxTimeStep)
596 maxTimeStep = tempMaxTimeStep;
597 }
598 } // end of parallel section
599
600 // maxTimeStep is now the maximum time step possible for all points
601 // and rates are stored in a vector
602 return maxTimeStep;
603 }
604
605 // Level Sets below are also considered in order to adjust the advection
606 // depth accordingly if there would be a material change.
607 void moveSurface() {
608 if (timeStepRatio >= 0.5) {
609 VIENNACORE_LOG_WARNING(
610 "Integration time step ratio should be smaller than 0.5. "
611 "Advection might fail!");
612 }
613
614 auto &topDomain = levelSets.back()->getDomain();
615 auto &grid = levelSets.back()->getGrid();
616
617 assert(currentTimeStep >= 0. && "No time step set!");
618 assert(storedRates.size() == topDomain.getNumberOfSegments());
619
620 // reduce to one layer thickness and apply new values directly to the
621 // domain segments --> DO NOT CHANGE SEGMENTATION HERE (true parameter)
622 Reduce<T, D>(levelSets.back(), 1, true).apply();
623
624 const bool saveVelocities = saveAdvectionVelocities;
625 std::vector<std::vector<double>> dissipationVectors(
626 levelSets.back()->getNumberOfSegments());
627 std::vector<std::vector<double>> velocityVectors(
628 levelSets.back()->getNumberOfSegments());
629
630 const bool checkDiss = checkDissipation;
631
632#pragma omp parallel num_threads(topDomain.getNumberOfSegments())
633 {
634 int p = 0;
635#ifdef _OPENMP
636 p = omp_get_thread_num();
637#endif
638 auto itRS = storedRates[p].cbegin();
639 auto &segment = topDomain.getDomainSegment(p);
640 const unsigned maxId = segment.getNumberOfPoints();
641
642 if (saveVelocities) {
643 velocityVectors[p].resize(maxId);
644 dissipationVectors[p].resize(maxId);
645 }
646
647 for (unsigned localId = 0; localId < maxId; ++localId) {
648 T &value = segment.definedValues[localId];
649 double time = currentTimeStep;
650
651 // if there is a change in materials during one time step, deduct
652 // the time taken to advect up to the end of the top material and
653 // set the LS value to the one below
654 auto const [gradient, dissipation] = itRS->first;
655 T velocity = gradient - dissipation;
656 // check if dissipation is too high
657 if (checkDiss && (gradient < 0 && velocity > 0) ||
658 (gradient > 0 && velocity < 0)) {
659 velocity = 0;
660 }
661
662 T rate = time * velocity;
663 while (std::abs(itRS->second - value) < std::abs(rate)) {
664 time -= std::abs((itRS->second - value) / velocity);
665 value = itRS->second;
666 ++itRS; // advance the TempStopRates iterator by one
667
668 // recalculate velocity and rate
669 velocity = itRS->first.first - itRS->first.second;
670 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
671 (itRS->first.first > 0 && velocity < 0)) {
672 velocity = 0;
673 }
674 rate = time * velocity;
675 }
676
677 // now deduct the velocity times the time step we take
678 value -= rate;
679
680 if (saveVelocities) {
681 velocityVectors[p][localId] = rate;
682 dissipationVectors[p][localId] = itRS->first.second;
683 }
684
685 // this is run when two materials are close but the velocity is too slow
686 // to actually reach the second material, to get rid of the extra
687 // entry in the TempRatesStop
688 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
689 ++itRS;
690
691 // advance the TempStopRates iterator by one
692 ++itRS;
693 }
694 } // end of parallel section
695
696 if (saveVelocities) {
697 auto &pointData = levelSets.back()->getPointData();
698
699 typename PointData<T>::ScalarDataType vels;
700 typename PointData<T>::ScalarDataType diss;
701
702 for (unsigned i = 0; i < velocityVectors.size(); ++i) {
703 vels.insert(vels.end(),
704 std::make_move_iterator(velocityVectors[i].begin()),
705 std::make_move_iterator(velocityVectors[i].end()));
706 diss.insert(diss.end(),
707 std::make_move_iterator(dissipationVectors[i].begin()),
708 std::make_move_iterator(dissipationVectors[i].end()));
709 }
710 pointData.insertReplaceScalarData(std::move(vels), velocityLabel);
711 pointData.insertReplaceScalarData(std::move(diss), dissipationLabel);
712 }
713
714 // clear the stored rates since surface has changed
715 storedRates.clear();
716 }
717
718public:
719 static constexpr char velocityLabel[] = "AdvectionVelocities";
720 static constexpr char dissipationLabel[] = "Dissipation";
721
722 Advect() = default;
723
724 explicit Advect(SmartPointer<Domain<T, D>> passedlsDomain) {
725 levelSets.push_back(passedlsDomain);
726 }
727
728 Advect(SmartPointer<Domain<T, D>> passedlsDomain,
729 SmartPointer<VelocityField<T>> passedVelocities) {
730 levelSets.push_back(passedlsDomain);
731 velocities = passedVelocities;
732 }
733
734 Advect(std::vector<SmartPointer<Domain<T, D>>> passedlsDomains,
735 SmartPointer<VelocityField<T>> passedVelocities)
736 : levelSets(passedlsDomains) {
737 velocities = passedVelocities;
738 }
739
742 void insertNextLevelSet(SmartPointer<Domain<T, D>> passedlsDomain) {
743 levelSets.push_back(passedlsDomain);
744 }
745
746 void clearLevelSets() { levelSets.clear(); }
747
750 void setVelocityField(SmartPointer<VelocityField<T>> passedVelocities) {
751 velocities = passedVelocities;
752 }
753
759 void setAdvectionTime(double time) { advectionTime = time; }
760
765 void setSingleStep(bool singleStep) { performOnlySingleStep = singleStep; }
766
771 void setTimeStepRatio(const double &cfl) { timeStepRatio = cfl; }
772
777 void setCalculateNormalVectors(bool cnv) { calculateNormalVectors = cnv; }
778
785 void setIgnoreVoids(bool iV) { ignoreVoids = iV; }
786
790 void setAdaptiveTimeStepping(bool aTS) { adaptiveTimeStepping = aTS; }
791
795 void setAdaptiveTimeStepThreshold(double threshold) {
796 adaptiveTimeStepThreshold = threshold;
797 }
798
801 void setSaveAdvectionVelocities(bool sAV) { saveAdvectionVelocities = sAV; }
802
805 double getAdvectedTime() const { return advectedTime; }
806
808 double getCurrentTimeStep() const { return currentTimeStep; }
809
811 unsigned getNumberOfTimeSteps() const { return numberOfTimeSteps; }
812
814 double getTimeStepRatio() const { return timeStepRatio; }
815
817 bool getCalculateNormalVectors() const { return calculateNormalVectors; }
818
822 integrationScheme = scheme;
823 }
824
829 void setDissipationAlpha(const double &a) { dissipationAlpha = a; }
830
831 // Sets the velocity to 0 if the dissipation is too high
832 void setCheckDissipation(bool check) { checkDissipation = check; }
833
836 void setUpdatePointData(bool update) { updatePointData = update; }
837
838 // Prepare the levelset for advection, based on the provided integration
839 // scheme.
840 void prepareLS() {
841 // check whether a level set and velocities have been given
842 if (levelSets.empty()) {
843 VIENNACORE_LOG_ERROR("No level sets passed to Advect.");
844 return;
845 }
846
847 if (integrationScheme == IntegrationSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
849 } else if (integrationScheme ==
852 } else if (integrationScheme ==
855 } else if (integrationScheme ==
858 } else if (integrationScheme ==
862 levelSets.back());
863 } else if (integrationScheme ==
866 } else if (integrationScheme ==
869 } else if (integrationScheme ==
872 } else if (integrationScheme ==
875 } else if (integrationScheme ==
878 levelSets.back());
879 } else if (integrationScheme == IntegrationSchemeEnum::WENO_5TH_ORDER) {
880 // WENO5 requires a stencil radius of 3 (template parameter 3)
882 } else {
883 VIENNACORE_LOG_ERROR("Advect: Integration scheme not found.");
884 }
885 }
886
888 void apply() {
889 if (advectionTime == 0.) {
890 advectedTime = advect(std::numeric_limits<double>::max());
891 numberOfTimeSteps = 1;
892 } else {
893 double currentTime = 0.0;
894 numberOfTimeSteps = 0;
895 while (currentTime < advectionTime) {
896 currentTime += advect(advectionTime - currentTime);
897 ++numberOfTimeSteps;
898 if (performOnlySingleStep)
899 break;
900 }
901 advectedTime = currentTime;
902 }
903 }
904
907 void
908 applyIntegration(double maxTimeStep = std::numeric_limits<double>::max()) {
909 prepareLS();
910 if (integrationScheme == IntegrationSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
911 auto is = lsInternal::EngquistOsher<T, D, 1>(levelSets.back(), velocities,
912 calculateNormalVectors);
913 currentTimeStep = integrateTime(is, maxTimeStep);
914 } else if (integrationScheme ==
916 auto is = lsInternal::EngquistOsher<T, D, 2>(levelSets.back(), velocities,
917 calculateNormalVectors);
918 currentTimeStep = integrateTime(is, maxTimeStep);
919 } else if (integrationScheme ==
921 auto alphas = findGlobalAlphas();
922 auto is = lsInternal::LaxFriedrichs<T, D, 1>(levelSets.back(), velocities,
923 dissipationAlpha, alphas,
924 calculateNormalVectors);
925 currentTimeStep = integrateTime(is, maxTimeStep);
926 } else if (integrationScheme ==
928 auto alphas = findGlobalAlphas();
929 auto is = lsInternal::LaxFriedrichs<T, D, 2>(levelSets.back(), velocities,
930 dissipationAlpha, alphas,
931 calculateNormalVectors);
932 currentTimeStep = integrateTime(is, maxTimeStep);
933 } else if (integrationScheme ==
937 levelSets.back(), velocities);
938 currentTimeStep = integrateTime(is, maxTimeStep);
939 } else if (integrationScheme ==
942 levelSets.back(), velocities, dissipationAlpha);
943 currentTimeStep = integrateTime(is, maxTimeStep);
944 } else if (integrationScheme ==
947 levelSets.back(), velocities, dissipationAlpha);
948 currentTimeStep = integrateTime(is, maxTimeStep);
949 } else if (integrationScheme ==
952 levelSets.back(), velocities, dissipationAlpha);
953 currentTimeStep = integrateTime(is, maxTimeStep);
954 } else if (integrationScheme ==
957 levelSets.back(), velocities, dissipationAlpha);
958 currentTimeStep = integrateTime(is, maxTimeStep);
959 } else if (integrationScheme ==
962 levelSets.back(), velocities, dissipationAlpha);
963 currentTimeStep = integrateTime(is, maxTimeStep);
964 } else if (integrationScheme == IntegrationSchemeEnum::WENO_5TH_ORDER) {
965 // Instantiate WENO5 with order 3 (neighbors +/- 3)
966 auto is = lsInternal::WENO5<T, D, 3>(levelSets.back(), velocities,
967 calculateNormalVectors);
968 currentTimeStep = integrateTime(is, maxTimeStep);
969 } else {
970 VIENNACORE_LOG_ERROR("Advect: Integration scheme not found.");
971 currentTimeStep = -1.;
972 }
973 }
974};
975
976// add all template specializations for this class
978
979} // namespace viennals
constexpr int D
Definition Epitaxy.cpp:11
double T
Definition Epitaxy.cpp:12
Engquist-Osher integration scheme based on the upwind integration scheme. Offers high performance but...
Definition lsEngquistOsher.hpp:18
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsEngquistOsher.hpp:28
Lax Friedrichs integration scheme with constant alpha value for dissipation. This alpha value should ...
Definition lsLaxFriedrichs.hpp:19
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLaxFriedrichs.hpp:32
Lax Friedrichs integration scheme, which uses alpha values provided by the user in getDissipationAlph...
Definition lsLocalLaxFriedrichsAnalytical.hpp:20
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLaxFriedrichsAnalytical.hpp:48
Lax Friedrichs integration scheme, which uses a first neighbour stencil to calculate the alpha values...
Definition lsLocalLaxFriedrichs.hpp:20
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLaxFriedrichs.hpp:49
Lax Friedrichs integration scheme, which considers only the current point for alpha calculation....
Definition lsLocalLocalLaxFriedrichs.hpp:18
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLocalLaxFriedrichs.hpp:29
Stencil Local Lax Friedrichs Integration 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:750
static constexpr char velocityLabel[]
Definition lsAdvect.hpp:719
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:742
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:836
bool getCalculateNormalVectors() const
Get whether normal vectors were calculated.
Definition lsAdvect.hpp:817
void setIgnoreVoids(bool iV)
Set whether level set values, which are not part of the "top" geometrically connected part of values,...
Definition lsAdvect.hpp:785
Advect()=default
void clearLevelSets()
Definition lsAdvect.hpp:746
void setCalculateNormalVectors(bool cnv)
Set whether normal vectors should be calculated at each level set point. Defaults to true....
Definition lsAdvect.hpp:777
void apply()
Perform the advection.
Definition lsAdvect.hpp:888
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:765
void setAdaptiveTimeStepThreshold(double threshold)
Set the threshold (in fraction of the CFL condition) below which adaptive time stepping is applied....
Definition lsAdvect.hpp:795
void prepareLS()
Definition lsAdvect.hpp:840
double getAdvectedTime() const
Get by how much the physical time was advanced during the last apply() call.
Definition lsAdvect.hpp:805
void setAdaptiveTimeStepping(bool aTS)
Set whether adaptive time stepping should be used when approaching material boundaries during etching...
Definition lsAdvect.hpp:790
void setDissipationAlpha(const double &a)
Set the alpha dissipation coefficient. For lsLaxFriedrichs, this is used as the alpha value....
Definition lsAdvect.hpp:829
void setCheckDissipation(bool check)
Definition lsAdvect.hpp:832
Advect(std::vector< SmartPointer< Domain< T, D > > > passedlsDomains, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:734
void setIntegrationScheme(IntegrationSchemeEnum scheme)
Set which integration scheme should be used out of the ones specified in IntegrationSchemeEnum.
Definition lsAdvect.hpp:821
double getCurrentTimeStep() const
Return the last applied time step.
Definition lsAdvect.hpp:808
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:759
unsigned getNumberOfTimeSteps() const
Get how many advection steps were performed during the last apply() call.
Definition lsAdvect.hpp:811
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:771
void applyIntegration(double maxTimeStep=std::numeric_limits< double >::max())
This function applies the integration scheme and calculates the rates and the maximum time step,...
Definition lsAdvect.hpp:908
Advect(SmartPointer< Domain< T, D > > passedlsDomain)
Definition lsAdvect.hpp:724
Advect(SmartPointer< Domain< T, D > > passedlsDomain, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:728
double getTimeStepRatio() const
Get the value of the CFL number.
Definition lsAdvect.hpp:814
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:801
static constexpr char dissipationLabel[]
Definition lsAdvect.hpp:720
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:37
IntegrationSchemeEnum
Enumeration for the different Integration schemes used by the advection kernel.
Definition lsAdvect.hpp:43
@ LOCAL_LOCAL_LAX_FRIEDRICHS_2ND_ORDER
Definition lsAdvect.hpp:50
@ STENCIL_LOCAL_LAX_FRIEDRICHS_1ST_ORDER
Definition lsAdvect.hpp:53
@ LOCAL_LOCAL_LAX_FRIEDRICHS_1ST_ORDER
Definition lsAdvect.hpp:49
@ LAX_FRIEDRICHS_2ND_ORDER
Definition lsAdvect.hpp:47
@ LOCAL_LAX_FRIEDRICHS_1ST_ORDER
Definition lsAdvect.hpp:51
@ ENGQUIST_OSHER_2ND_ORDER
Definition lsAdvect.hpp:45
@ LAX_FRIEDRICHS_1ST_ORDER
Definition lsAdvect.hpp:46
@ LOCAL_LAX_FRIEDRICHS_2ND_ORDER
Definition lsAdvect.hpp:52
@ ENGQUIST_OSHER_1ST_ORDER
Definition lsAdvect.hpp:44
@ LOCAL_LAX_FRIEDRICHS_ANALYTICAL_1ST_ORDER
Definition lsAdvect.hpp:48
@ WENO_5TH_ORDER
Definition lsAdvect.hpp:54
@ INTERSECT
Definition lsBooleanOperation.hpp:26