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 <lsConcepts.hpp>
21#include <lsEnquistOsher.hpp>
22#include <lsLaxFriedrichs.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
55
65template <class T, int D> class Advect {
66 std::vector<SmartPointer<Domain<T, D>>> levelSets;
67 SmartPointer<VelocityField<T>> velocities = nullptr;
68 IntegrationSchemeEnum integrationScheme =
70 double timeStepRatio = 0.4999;
71 double dissipationAlpha = 1.0;
72 bool calculateNormalVectors = true;
73 bool ignoreVoids = false;
74 double advectionTime = 0.;
75 bool performOnlySingleStep = false;
76 double advectedTime = 0.;
77 unsigned numberOfTimeSteps = 0;
78 bool saveAdvectionVelocities = false;
79 bool updatePointData = true;
80 bool checkDissipation = true;
81 static constexpr double wrappingLayerEpsilon = 1e-4;
82
83 hrleVectorType<T, 3> findGlobalAlphas() const {
84
85 auto &topDomain = levelSets.back()->getDomain();
86 auto &grid = levelSets.back()->getGrid();
87
88 const T gridDelta = grid.getGridDelta();
89 const T deltaPos = gridDelta;
90 const T deltaNeg = -gridDelta;
91
92 hrleVectorType<T, 3> finalAlphas(0., 0., 0.);
93
94#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
95 {
96 int p = 0;
97#ifdef _OPENMP
98 p = omp_get_thread_num();
99#endif
100
101 hrleVectorType<T, 3> localAlphas(0., 0., 0.);
102
103 hrleVectorType<hrleIndexType, D> startVector =
104 (p == 0) ? grid.getMinGridPoint()
105 : topDomain.getSegmentation()[p - 1];
106
107 hrleVectorType<hrleIndexType, D> endVector =
108 (p != static_cast<int>(topDomain.getNumberOfSegments() - 1))
109 ? topDomain.getSegmentation()[p]
110 : grid.incrementIndices(grid.getMaxGridPoint());
111
112 // an iterator for each level set
113 std::vector<hrleConstSparseIterator<typename Domain<T, D>::DomainType>>
114 iterators;
115 for (auto it = levelSets.begin(); it != levelSets.end(); ++it) {
116 iterators.push_back(
117 hrleConstSparseIterator<typename Domain<T, D>::DomainType>(
118 (*it)->getDomain()));
119 }
120
121 // neighborIterator for the top level set
122 hrleConstSparseStarIterator<hrleDomain<T, D>, 1> neighborIterator(
123 topDomain);
124
125 for (hrleSparseIterator<typename Domain<T, D>::DomainType> it(
126 topDomain, startVector);
127 it.getStartIndices() < endVector; ++it) {
128
129 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
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 T diffPos = (phiPos - value) / deltaPos;
162 T diffNeg = (phiNeg - value) / deltaNeg;
163
164 normal[i] = (diffNeg + diffPos) * 0.5;
165 normalModulus += normal[i] * normal[i];
166 }
167 normalModulus = std::sqrt(normalModulus);
168 for (unsigned i = 0; i < D; ++i)
169 normal[i] /= normalModulus;
170
171 T scaVel = velocities->getScalarVelocity(
172 coords, lowerLevelSetId, normal,
173 neighborIterator.getCenter().getPointId());
174 auto vecVel = velocities->getVectorVelocity(
175 coords, lowerLevelSetId, normal,
176 neighborIterator.getCenter().getPointId());
177
178 for (unsigned i = 0; i < D; ++i) {
179 T tempAlpha = std::abs((scaVel + vecVel[i]) * normal[i]);
180 localAlphas[i] = std::max(localAlphas[i], tempAlpha);
181 }
182
183 break;
184 }
185 }
186 }
187
188#pragma omp critical
189 {
190 for (unsigned i = 0; i < D; ++i) {
191 finalAlphas[i] = std::max(finalAlphas[i], localAlphas[i]);
192 }
193 }
194 } // end of parallel section
195
196 return finalAlphas;
197 }
198
199 void rebuildLS() {
200 // TODO: this function uses manhatten distances for renormalisation,
201 // since this is the quickest. For visualisation applications, better
202 // renormalisation is needed, so it might be good to implement
203 // Euler distance renormalisation as an option
204 auto &grid = levelSets.back()->getGrid();
205 auto newlsDomain = SmartPointer<Domain<T, D>>::New(grid);
206 typename Domain<T, D>::DomainType &newDomain = newlsDomain->getDomain();
207 typename Domain<T, D>::DomainType &domain = levelSets.back()->getDomain();
208
209 newDomain.initialize(domain.getNewSegmentation(),
210 domain.getAllocation() *
211 (2.0 / levelSets.back()->getLevelSetWidth()));
212
213 const bool updateData = updatePointData;
214 // save how data should be transferred to new level set
215 // list of indices into the old pointData vector
216 std::vector<std::vector<unsigned>> newDataSourceIds;
217 if (updateData)
218 newDataSourceIds.resize(newDomain.getNumberOfSegments());
219
220#ifdef DEBUG_LS_ADVECT_HPP
221 {
222 auto mesh = SmartPointer<Mesh<T>>::New();
223 ToMesh<T, D>(levelSets.back(), mesh).apply();
224 VTKWriter<T>(mesh, "Advect_beforeRebuild.vtk").apply();
225 }
226#endif
227
228#pragma omp parallel num_threads(newDomain.getNumberOfSegments())
229 {
230 int p = 0;
231#ifdef _OPENMP
232 p = omp_get_thread_num();
233#endif
234
235 auto &domainSegment = newDomain.getDomainSegment(p);
236
237 hrleVectorType<hrleIndexType, D> startVector =
238 (p == 0) ? grid.getMinGridPoint()
239 : newDomain.getSegmentation()[p - 1];
240
241 hrleVectorType<hrleIndexType, D> endVector =
242 (p != static_cast<int>(newDomain.getNumberOfSegments() - 1))
243 ? newDomain.getSegmentation()[p]
244 : grid.incrementIndices(grid.getMaxGridPoint());
245
246 // reserve a bit more to avoid reallocation
247 // would expect number of points to roughly double
248 if (updateData)
249 newDataSourceIds[p].reserve(2.5 * domainSegment.getNumberOfPoints());
250
251 for (hrleSparseStarIterator<typename Domain<T, D>::DomainType, 1> it(
252 domain, startVector);
253 it.getIndices() < endVector; ++it) {
254
255 // if the center is an active grid point
256 // <1.0 since it could have been change by 0.5 max
257 if (std::abs(it.getCenter().getValue()) <= 1.0) {
258
259 int k = 0;
260 for (; k < 2 * D; k++)
261 if (std::signbit(it.getNeighbor(k).getValue() - 1e-7) !=
262 std::signbit(it.getCenter().getValue() + 1e-7))
263 break;
264
265 // if there is at least one neighbor of opposite sign
266 if (k != 2 * D) {
267 if (it.getCenter().getDefinedValue() > 0.5) {
268 int j = 0;
269 for (; j < 2 * D; j++) {
270 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
271 if (it.getNeighbor(j).getDefinedValue() < -0.5)
272 break;
273 }
274 if (j == 2 * D) {
275 domainSegment.insertNextDefinedPoint(
276 it.getIndices(), it.getCenter().getDefinedValue());
277 if (updateData)
278 newDataSourceIds[p].push_back(it.getCenter().getPointId());
279 // if there is at least one active grid point, which is < -0.5
280 } else {
281 domainSegment.insertNextDefinedPoint(it.getIndices(), 0.5);
282 if (updateData)
283 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
284 }
285 } else if (it.getCenter().getDefinedValue() < -0.5) {
286 int j = 0;
287 for (; j < 2 * D; j++) {
288 if (std::abs(it.getNeighbor(j).getValue()) <= 1.0)
289 if (it.getNeighbor(j).getDefinedValue() > 0.5)
290 break;
291 }
292
293 if (j == 2 * D) {
294 domainSegment.insertNextDefinedPoint(
295 it.getIndices(), it.getCenter().getDefinedValue());
296 if (updateData)
297 newDataSourceIds[p].push_back(it.getCenter().getPointId());
298 // if there is at least one active grid point, which is > 0.5
299 } else {
300 domainSegment.insertNextDefinedPoint(it.getIndices(), -0.5);
301 if (updateData)
302 newDataSourceIds[p].push_back(it.getNeighbor(j).getPointId());
303 }
304 } else {
305 domainSegment.insertNextDefinedPoint(
306 it.getIndices(), it.getCenter().getDefinedValue());
307 if (updateData)
308 newDataSourceIds[p].push_back(it.getCenter().getPointId());
309 }
310 } else {
311 domainSegment.insertNextUndefinedPoint(
312 it.getIndices(), (it.getCenter().getDefinedValue() < 0)
315 }
316
317 } else { // if the center is not an active grid point
318 if (it.getCenter().getValue() >= 0) {
319 int usedNeighbor = -1;
320 T distance = Domain<T, D>::POS_VALUE;
321 for (int i = 0; i < 2 * D; i++) {
322 T value = it.getNeighbor(i).getValue();
323 if (std::abs(value) <= 1.0 && (value < 0.)) {
324 if (distance > value + 1.0) {
325 distance = value + 1.0;
326 usedNeighbor = i;
327 }
328 }
329 }
330
331 if (distance <= 1.) {
332 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
333 if (updateData)
334 newDataSourceIds[p].push_back(
335 it.getNeighbor(usedNeighbor).getPointId());
336 } else {
337 domainSegment.insertNextUndefinedPoint(it.getIndices(),
339 }
340
341 } else {
342 int usedNeighbor = -1;
343 T distance = Domain<T, D>::NEG_VALUE;
344 for (int i = 0; i < 2 * D; i++) {
345 T value = it.getNeighbor(i).getValue();
346 if (std::abs(value) <= 1.0 && (value > 0)) {
347 if (distance < value - 1.0) {
348 // distance = std::max(distance, value - T(1.0));
349 distance = value - 1.0;
350 usedNeighbor = i;
351 }
352 }
353 }
354
355 if (distance >= -1.) {
356 domainSegment.insertNextDefinedPoint(it.getIndices(), distance);
357 if (updateData)
358 newDataSourceIds[p].push_back(
359 it.getNeighbor(usedNeighbor).getPointId());
360 } else {
361 domainSegment.insertNextUndefinedPoint(it.getIndices(),
363 }
364 }
365 }
366 }
367 }
368
369 // now copy old data into new level set
370 if (updateData) {
371 auto &pointData = levelSets.back()->getPointData();
372 newlsDomain->getPointData().translateFromMultiData(pointData,
373 newDataSourceIds);
374 }
375
376 newDomain.finalize();
377 newDomain.segment();
378 levelSets.back()->deepCopy(newlsDomain);
379 levelSets.back()->finalize(2);
380 }
381
384 double advect(double maxTimeStep = std::numeric_limits<double>::max()) {
385 // check whether a level set and velocites have been given
386 if (levelSets.size() < 1) {
387 Logger::getInstance()
388 .addWarning("No level sets passed to Advect. Not advecting.")
389 .print();
390 return std::numeric_limits<double>::max();
391 }
392 if (velocities == nullptr) {
393 Logger::getInstance()
394 .addWarning("No velocity field passed to Advect. Not advecting.")
395 .print();
396 return std::numeric_limits<double>::max();
397 }
398
399 prepareLS();
400
401 double currentTime = 0.;
402 if (integrationScheme == IntegrationSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
403 auto is = lsInternal::EnquistOsher<T, D, 1>(levelSets.back(), velocities,
404 calculateNormalVectors);
405 currentTime = integrateTime(is, maxTimeStep);
406 } else if (integrationScheme ==
408 auto is = lsInternal::EnquistOsher<T, D, 2>(levelSets.back(), velocities,
409 calculateNormalVectors);
410 currentTime = integrateTime(is, maxTimeStep);
411 } else if (integrationScheme ==
413 auto alphas = findGlobalAlphas();
414 auto is = lsInternal::LaxFriedrichs<T, D, 1>(levelSets.back(), velocities,
415 dissipationAlpha, alphas,
416 calculateNormalVectors);
417 currentTime = integrateTime(is, maxTimeStep);
418 } else if (integrationScheme ==
420 auto alphas = findGlobalAlphas();
421 auto is = lsInternal::LaxFriedrichs<T, D, 2>(levelSets.back(), velocities,
422 dissipationAlpha, alphas,
423 calculateNormalVectors);
424 currentTime = integrateTime(is, maxTimeStep);
425 } else if (integrationScheme ==
429 levelSets.back(), velocities);
430 currentTime = integrateTime(is, maxTimeStep);
431 } else if (integrationScheme ==
434 levelSets.back(), velocities, dissipationAlpha);
435 currentTime = integrateTime(is, maxTimeStep);
436 } else if (integrationScheme ==
439 levelSets.back(), velocities, dissipationAlpha);
440 currentTime = integrateTime(is, maxTimeStep);
441 } else if (integrationScheme ==
444 levelSets.back(), velocities, dissipationAlpha);
445 currentTime = integrateTime(is, maxTimeStep);
446 } else if (integrationScheme ==
449 levelSets.back(), velocities, dissipationAlpha);
450 currentTime = integrateTime(is, maxTimeStep);
451 } else if (integrationScheme ==
454 levelSets.back(), velocities, dissipationAlpha);
455 currentTime = integrateTime(is, maxTimeStep);
456 } else {
457 Logger::getInstance()
458 .addWarning("Advect: Integration scheme not found. Not advecting.")
459 .print();
460
461 // if no correct scheme was found return infinity
462 // to stop advection
463 return std::numeric_limits<double>::max();
464 }
465
466 rebuildLS();
467
468 // Adjust all level sets below the advected one
469 // This means, that when the top level-set and one below
470 // are etched, the lower one is moved with the top level-set
471 // TODO: Adjust lower layers also when they have grown,
472 // to allow for two different growth rates of materials
473 if (integrationScheme !=
475 for (unsigned i = 0; i < levelSets.size() - 1; ++i) {
476 BooleanOperation<T, D>(levelSets[i], levelSets.back(),
478 .apply();
479 }
480 }
481
482 return currentTime;
483 }
484
489 template <class IntegrationSchemeType>
490 double
491 integrateTime(IntegrationSchemeType IntegrationScheme,
492 double maxTimeStep = std::numeric_limits<double>::max()) {
493 if (timeStepRatio >= 0.5) {
494 Logger::getInstance()
495 .addWarning("Integration time step ratio should be smaller than 0.5. "
496 "Advection might fail!")
497 .print();
498 }
499
500 auto &topDomain = levelSets.back()->getDomain();
501 auto &grid = levelSets.back()->getGrid();
502
503 std::vector<std::vector<std::pair<std::pair<T, T>, T>>> totalTempRates;
504 totalTempRates.resize((levelSets.back())->getNumberOfSegments());
505
506 typename PointData<T>::ScalarDataType *voidMarkerPointer;
507 if (ignoreVoids) {
508 MarkVoidPoints<T, D>(levelSets.back()).apply();
509 auto &pointData = levelSets.back()->getPointData();
510 voidMarkerPointer =
512 if (voidMarkerPointer == nullptr) {
513 Logger::getInstance()
514 .addWarning("Advect: Cannot find void point markers. Not "
515 "ignoring void points.")
516 .print();
517 ignoreVoids = false;
518 }
519 }
520
521 const bool ignoreVoidPoints = ignoreVoids;
522
523#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
524 {
525 int p = 0;
526#ifdef _OPENMP
527 p = omp_get_thread_num();
528#endif
529
530 hrleVectorType<hrleIndexType, D> startVector =
531 (p == 0) ? grid.getMinGridPoint()
532 : topDomain.getSegmentation()[p - 1];
533
534 hrleVectorType<hrleIndexType, D> endVector =
535 (p != static_cast<int>(topDomain.getNumberOfSegments() - 1))
536 ? topDomain.getSegmentation()[p]
537 : grid.incrementIndices(grid.getMaxGridPoint());
538
539 double tempMaxTimeStep = maxTimeStep;
540 auto &tempRates = totalTempRates[p];
541 tempRates.reserve(topDomain.getNumberOfPoints() /
542 double((levelSets.back())->getNumberOfSegments()) +
543 10);
544
545 // an iterator for each level set
546 std::vector<hrleSparseIterator<typename Domain<T, D>::DomainType>>
547 iterators;
548 for (auto it = levelSets.begin(); it != levelSets.end(); ++it) {
549 iterators.push_back(
550 hrleSparseIterator<typename Domain<T, D>::DomainType>(
551 (*it)->getDomain()));
552 }
553
554 IntegrationSchemeType scheme(IntegrationScheme);
555
556 for (hrleSparseIterator<typename Domain<T, D>::DomainType> it(
557 topDomain, startVector);
558 it.getStartIndices() < endVector; ++it) {
559
560 if (!it.isDefined() || std::abs(it.getValue()) > 0.5)
561 continue;
562
563 T value = it.getValue();
564 double maxStepTime = 0;
565 double cfl = timeStepRatio;
566
567 for (int currentLevelSetId = levelSets.size() - 1;
568 currentLevelSetId >= 0; --currentLevelSetId) {
569
570 std::pair<T, T> gradNDissipation;
571
572 if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) {
573 // check if there is any other levelset at the same point:
574 // if yes, take the velocity of the lowest levelset
575 for (unsigned lowerLevelSetId = 0;
576 lowerLevelSetId < levelSets.size(); ++lowerLevelSetId) {
577 // put iterator to same position as the top levelset
578 iterators[lowerLevelSetId].goToIndicesSequential(
579 it.getStartIndices());
580
581 // if the lower surface is actually outside, i.e. its LS value
582 // is lower or equal
583 if (iterators[lowerLevelSetId].getValue() <=
584 value + wrappingLayerEpsilon) {
585 gradNDissipation =
586 scheme(it.getStartIndices(), lowerLevelSetId);
587 break;
588 }
589 }
590 }
591
592 T valueBelow;
593 // get value of material below (before in levelSets list)
594 if (currentLevelSetId > 0) {
595 iterators[currentLevelSetId - 1].goToIndicesSequential(
596 it.getStartIndices());
597 valueBelow = iterators[currentLevelSetId - 1].getValue();
598 } else {
599 valueBelow = std::numeric_limits<T>::max();
600 }
601
602 // if velocity is positive, set maximum time step possible without
603 // violating the cfl condition
604 T velocity = gradNDissipation.first - gradNDissipation.second;
605 if (velocity > 0.) {
606 maxStepTime += cfl / velocity;
607 tempRates.push_back(std::make_pair(gradNDissipation,
608 -std::numeric_limits<T>::max()));
609 break;
610 // if velocity is 0, maximum time step is infinite
611 } else if (velocity == 0.) {
612 maxStepTime = std::numeric_limits<T>::max();
613 tempRates.push_back(std::make_pair(gradNDissipation,
614 std::numeric_limits<T>::max()));
615 break;
616 // if the velocity is negative apply the velocity for as long as
617 // possible without infringing on material below
618 } else {
619 T difference = std::abs(valueBelow - value);
620
621 if (difference >= cfl) {
622 maxStepTime -= cfl / velocity;
623 tempRates.push_back(std::make_pair(
624 gradNDissipation, std::numeric_limits<T>::max()));
625 break;
626 } else {
627 maxStepTime -= difference / velocity;
628 // the second part of the pair indicates how far we can move
629 // in this time step until the end of the material is reached
630 tempRates.push_back(std::make_pair(gradNDissipation, valueBelow));
631 cfl -= difference;
632 // use new LS value for next calculations
633 value = valueBelow;
634 }
635 }
636 }
637
638 if (maxStepTime < tempMaxTimeStep)
639 tempMaxTimeStep = maxStepTime;
640 }
641
642#pragma omp critical
643 {
644 // If a Lax Friedrichs scheme is selected the time step is
645 // reduced depending on the dissipation coefficients
646 // For Engquist Osher scheme this function is empty.
647 scheme.reduceTimeStepHamiltonJacobi(
648 tempMaxTimeStep, levelSets.back()->getGrid().getGridDelta());
649
650 // set global timestep maximum
651 if (tempMaxTimeStep < maxTimeStep)
652 maxTimeStep = tempMaxTimeStep;
653 }
654 } // end of parallel section
655
656 // reduce to one layer thickness and apply new values directly to the
657 // domain segments --> DO NOT CHANGE SEGMENTATION HERE (true parameter)
658 Reduce<T, D>(levelSets.back(), 1, true).apply();
659
660 const bool saveVelocities = saveAdvectionVelocities;
661 std::vector<std::vector<double>> dissipationVectors(
662 levelSets.back()->getNumberOfSegments());
663 std::vector<std::vector<double>> velocityVectors(
664 levelSets.back()->getNumberOfSegments());
665
666 const bool checkDiss = checkDissipation;
667
668#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments())
669 {
670 int p = 0;
671#ifdef _OPENMP
672 p = omp_get_thread_num();
673#endif
674
675 auto itRS = totalTempRates[p].cbegin();
676 auto &segment = topDomain.getDomainSegment(p);
677 unsigned maxId = segment.getNumberOfPoints();
678
679 if (saveVelocities) {
680 velocityVectors[p].resize(maxId);
681 dissipationVectors[p].resize(maxId);
682 }
683
684 for (unsigned localId = 0; localId < maxId; ++localId) {
685 T &value = segment.definedValues[localId];
686 double time = maxTimeStep;
687
688 // if there is a change in materials during one time step, deduct
689 // the time taken to advect up to the end of the top material and
690 // set the LS value to the one below
691 T velocity = itRS->first.first - itRS->first.second;
692 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
693 (itRS->first.first > 0 && velocity < 0)) {
694 velocity = 0;
695 }
696 T rate = time * velocity;
697 while (std::abs(itRS->second - value) < std::abs(rate)) {
698 time -= std::abs((itRS->second - value) / velocity);
699 value = itRS->second;
700 ++itRS;
701 // recalculate velocity and rate
702 velocity = itRS->first.first - itRS->first.second;
703 if (checkDiss && (itRS->first.first < 0 && velocity > 0) ||
704 (itRS->first.first > 0 && velocity < 0)) {
705 velocity = 0;
706 }
707 rate = time * velocity;
708 }
709
710 // now deduct the velocity times the time step we take
711 value -= rate;
712
713 if (saveVelocities) {
714 velocityVectors[p][localId] = rate;
715 dissipationVectors[p][localId] = itRS->first.second;
716 }
717
718 // this
719 // is run when two materials are close but the velocity is too slow
720 // to actually reach the second material, to get rid of the extra
721 // entry in the TempRatesStop
722 while (std::abs(itRS->second) != std::numeric_limits<T>::max())
723 ++itRS;
724
725 // advance the TempStopRates iterator by one
726 ++itRS;
727 }
728 } // end of parallel section
729
730 if (saveVelocities) {
731 auto &pointData = levelSets.back()->getPointData();
732
733 typename PointData<T>::ScalarDataType vels;
734 typename PointData<T>::ScalarDataType diss;
735
736 for (unsigned i = 0; i < velocityVectors.size(); ++i) {
737 vels.insert(vels.end(),
738 std::make_move_iterator(velocityVectors[i].begin()),
739 std::make_move_iterator(velocityVectors[i].end()));
740 diss.insert(diss.end(),
741 std::make_move_iterator(dissipationVectors[i].begin()),
742 std::make_move_iterator(dissipationVectors[i].end()));
743 }
744 pointData.insertReplaceScalarData(std::move(vels), velocityLabel);
745 pointData.insertReplaceScalarData(std::move(diss), dissipationLabel);
746 }
747
748 return maxTimeStep;
749 }
750
751public:
752 static constexpr char velocityLabel[] = "AdvectionVelocities";
753 static constexpr char dissipationLabel[] = "Dissipation";
754
756
757 Advect(SmartPointer<Domain<T, D>> passedlsDomain) {
758 levelSets.push_back(passedlsDomain);
759 }
760
761 Advect(SmartPointer<Domain<T, D>> passedlsDomain,
762 SmartPointer<VelocityField<T>> passedVelocities) {
763 levelSets.push_back(passedlsDomain);
764 velocities = passedVelocities;
765 }
766
767 Advect(std::vector<SmartPointer<Domain<T, D>>> passedlsDomains,
768 SmartPointer<VelocityField<T>> passedVelocities)
769 : levelSets(passedlsDomains) {
770 velocities = passedVelocities;
771 }
772
775 void insertNextLevelSet(SmartPointer<Domain<T, D>> passedlsDomain) {
776 levelSets.push_back(passedlsDomain);
777 }
778
781 void setVelocityField(SmartPointer<VelocityField<T>> passedVelocities) {
782 velocities = passedVelocities;
783 }
784
790 void setAdvectionTime(double time) { advectionTime = time; }
791
796 void setSingleStep(bool singleStep) { performOnlySingleStep = singleStep; }
797
802 void setTimeStepRatio(const double &cfl) { timeStepRatio = cfl; }
803
808 void setCalculateNormalVectors(bool cnv) { calculateNormalVectors = cnv; }
809
816 void setIgnoreVoids(bool iV) { ignoreVoids = iV; }
817
820 void setSaveAdvectionVelocities(bool sAV) { saveAdvectionVelocities = sAV; }
821
824 double getAdvectedTime() { return advectedTime; }
825
827 unsigned getNumberOfTimeSteps() { return numberOfTimeSteps; }
828
830 double getTimeStepRatio() { return timeStepRatio; }
831
833 bool getCalculateNormalVectors() { return calculateNormalVectors; }
834
838 integrationScheme = scheme;
839 }
840
845 void setDissipationAlpha(const double &a) { dissipationAlpha = a; }
846
847 // Sets the velocity to 0 if the dissipation is too high
848 void setCheckDissipation(bool check) { checkDissipation = check; }
849
852 void setUpdatePointData(bool update) { updatePointData = update; }
853
854 // Prepare the levelset for advection, based on the provided integration
855 // scheme.
856 void prepareLS() {
857 // check whether a level set and velocities have been given
858 if (levelSets.size() < 1) {
859 Logger::getInstance()
860 .addWarning("No level sets passed to Advect.")
861 .print();
862 return;
863 }
864
865 if (integrationScheme == IntegrationSchemeEnum::ENGQUIST_OSHER_1ST_ORDER) {
867 } else if (integrationScheme ==
870 } else if (integrationScheme ==
873 } else if (integrationScheme ==
876 } else if (integrationScheme ==
880 levelSets.back());
881 } else if (integrationScheme ==
884 } else if (integrationScheme ==
887 } else if (integrationScheme ==
890 } else if (integrationScheme ==
893 } else if (integrationScheme ==
896 levelSets.back());
897 } else {
898 Logger::getInstance()
899 .addWarning("Advect: Integration scheme not found. Not advecting.")
900 .print();
901 }
902 }
903
905 void apply() {
906 if (advectionTime == 0.) {
907 advectedTime = advect();
908 numberOfTimeSteps = 1;
909 } else {
910 double currentTime = 0.0;
911 numberOfTimeSteps = 0;
912 while (currentTime < advectionTime) {
913 currentTime += advect(advectionTime - currentTime);
914 ++numberOfTimeSteps;
915 if (performOnlySingleStep)
916 break;
917 }
918 advectedTime = currentTime;
919 }
920 }
921};
922
923// add all template specializations for this class
925
926} // namespace viennals
Engquist osher integration scheme based on the upwind integration scheme. Offers high performance but...
Definition lsEnquistOsher.hpp:19
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsEnquistOsher.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:31
Lax Friedrichs integration scheme, which uses alpha values provided by the user in getDissipationAlph...
Definition lsLocalLaxFriedrichsAnalytical.hpp:21
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLaxFriedrichsAnalytical.hpp:47
Lax Friedrichs integration scheme, which uses a first neighbour stencil to calculate the alpha values...
Definition lsLocalLaxFriedrichs.hpp:21
static void prepareLS(SmartPointer< viennals::Domain< T, D > > passedlsDomain)
Definition lsLocalLaxFriedrichs.hpp:48
Lax Friedrichs integration scheme, which considers only the current point for alpha calculation....
Definition lsLocalLocalLaxFriedrichs.hpp:19
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:21
static void prepareLS(LevelSetType passedlsDomain)
Definition lsStencilLocalLaxFriedrichsScalar.hpp:131
void setVelocityField(SmartPointer< VelocityField< T > > passedVelocities)
Set the velocity field used for advection. This should be a concrete implementation of lsVelocityFiel...
Definition lsAdvect.hpp:781
static constexpr char velocityLabel[]
Definition lsAdvect.hpp:752
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:775
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:852
bool getCalculateNormalVectors()
Get whether normal vectors were caluclated.
Definition lsAdvect.hpp:833
void setIgnoreVoids(bool iV)
Set whether level set values, which are not part of the "top" geometrically connected part of values,...
Definition lsAdvect.hpp:816
void setCalculateNormalVectors(bool cnv)
Set whether normal vectors should be calculated at each level set point. Defaults to true....
Definition lsAdvect.hpp:808
void apply()
Perform the advection.
Definition lsAdvect.hpp:905
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:796
unsigned getNumberOfTimeSteps()
Get how many advection steps were performed during the last apply() call.
Definition lsAdvect.hpp:827
double getTimeStepRatio()
Get the value of the CFL number.
Definition lsAdvect.hpp:830
void prepareLS()
Definition lsAdvect.hpp:856
void setDissipationAlpha(const double &a)
Set the alpha dissipation coefficient. For lsLaxFriedrichs, this is used as the alpha value....
Definition lsAdvect.hpp:845
void setCheckDissipation(bool check)
Definition lsAdvect.hpp:848
Advect(std::vector< SmartPointer< Domain< T, D > > > passedlsDomains, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:767
void setIntegrationScheme(IntegrationSchemeEnum scheme)
Set which integration scheme should be used out of the ones specified in IntegrationSchemeEnum.
Definition lsAdvect.hpp:837
double getAdvectedTime()
Get by how much the physical time was advanced during the last apply() call.
Definition lsAdvect.hpp:824
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:790
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:802
Advect()
Definition lsAdvect.hpp:755
Advect(SmartPointer< Domain< T, D > > passedlsDomain)
Definition lsAdvect.hpp:757
Advect(SmartPointer< Domain< T, D > > passedlsDomain, SmartPointer< VelocityField< T > > passedVelocities)
Definition lsAdvect.hpp:761
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:820
static constexpr char dissipationLabel[]
Definition lsAdvect.hpp:753
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:312
Class containing all information about the level set, including the dimensions of the domain,...
Definition lsDomain.hpp:27
static constexpr T NEG_VALUE
Definition lsDomain.hpp:54
void deepCopy(const SmartPointer< Domain< T, D > > passedDomain)
copy all values of "passedDomain" to this Domain
Definition lsDomain.hpp:120
unsigned getNumberOfSegments() const
returns the number of segments, the levelset is split into. This is useful for algorithm parallelisat...
Definition lsDomain.hpp:149
void finalize(int newWidth)
this function sets a new levelset width and finalizes the levelset, so it is ready for use by other a...
Definition lsDomain.hpp:113
static constexpr T POS_VALUE
Definition lsDomain.hpp:53
hrleDomain< T, D > DomainType
Definition lsDomain.hpp:32
This class is used to mark points of the level set which are enclosed in a void.
Definition lsMarkVoidPoints.hpp:30
void apply()
Definition lsMarkVoidPoints.hpp:154
static constexpr char voidPointLabel[]
Definition lsMarkVoidPoints.hpp:86
std::vector< T > ScalarDataType
Definition lsPointData.hpp:22
ScalarDataType * getScalarData(int index)
Definition lsPointData.hpp:130
Reduce the level set size to the specified width. This means all level set points with value <= 0....
Definition lsReduce.hpp:16
void apply()
Reduces the leveleSet to the specified number of layers. The largest value in the levelset is thus wi...
Definition lsReduce.hpp:57
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:31
void apply()
Definition lsVTKWriter.hpp:88
Abstract class defining the interface for the velocity field used during advection using lsAdvect.
Definition lsVelocityField.hpp:13
#define PRECOMPILE_PRECISION_DIMENSION(className)
Definition lsPreCompileMacros.hpp:24
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
@ INTERSECT
Definition lsBooleanOperation.hpp:26
constexpr int D
Definition pyWrap.cpp:66
double T
Definition pyWrap.cpp:64