diff --git a/CMakeLists.txt b/CMakeLists.txt index ae99d7a5..79008e4c 100644 --- a/CMakeLists.txt +++ b/CMakeLists.txt @@ -2,7 +2,7 @@ cmake_minimum_required(VERSION 3.20 FATAL_ERROR) project( ViennaLS LANGUAGES CXX - VERSION 4.1.2) + VERSION 4.1.3) # -------------------------------------------------------------------------------------------------------- # Library options diff --git a/README.md b/README.md index 87229791..919a23ae 100644 --- a/README.md +++ b/README.md @@ -138,7 +138,7 @@ We recommend using [CPM.cmake](https://github.com/cpm-cmake/CPM.cmake) to consum * Installation with CPM ```cmake - CPMAddPackage("gh:viennatools/viennals@4.1.2") + CPMAddPackage("gh:viennatools/viennals@4.1.3") ``` * With a local installation diff --git a/include/viennals/lsAdvect.hpp b/include/viennals/lsAdvect.hpp index 7e669129..a9c32d03 100644 --- a/include/viennals/lsAdvect.hpp +++ b/include/viennals/lsAdvect.hpp @@ -77,8 +77,125 @@ template class Advect { unsigned numberOfTimeSteps = 0; bool saveAdvectionVelocities = false; bool updatePointData = true; + bool checkDissipation = true; static constexpr double wrappingLayerEpsilon = 1e-4; + hrleVectorType findGlobalAlphas() const { + + auto &topDomain = levelSets.back()->getDomain(); + auto &grid = levelSets.back()->getGrid(); + + const T gridDelta = grid.getGridDelta(); + const T deltaPos = gridDelta; + const T deltaNeg = -gridDelta; + + hrleVectorType finalAlphas(0., 0., 0.); + +#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments()) + { + int p = 0; +#ifdef _OPENMP + p = omp_get_thread_num(); +#endif + + hrleVectorType localAlphas(0., 0., 0.); + + hrleVectorType startVector = + (p == 0) ? grid.getMinGridPoint() + : topDomain.getSegmentation()[p - 1]; + + hrleVectorType endVector = + (p != static_cast(topDomain.getNumberOfSegments() - 1)) + ? topDomain.getSegmentation()[p] + : grid.incrementIndices(grid.getMaxGridPoint()); + + // an iterator for each level set + std::vector::DomainType>> + iterators; + for (auto it = levelSets.begin(); it != levelSets.end(); ++it) { + iterators.push_back( + hrleConstSparseIterator::DomainType>( + (*it)->getDomain())); + } + + // neighborIterator for the top level set + hrleConstSparseStarIterator, 1> neighborIterator( + topDomain); + + for (hrleSparseIterator::DomainType> it( + topDomain, startVector); + it.getStartIndices() < endVector; ++it) { + + if (!it.isDefined() || std::abs(it.getValue()) > 0.5) + continue; + + const T value = it.getValue(); + const auto indices = it.getStartIndices(); + + // check if there is any other levelset at the same point: + // if yes, take the velocity of the lowest levelset + for (unsigned lowerLevelSetId = 0; lowerLevelSetId < levelSets.size(); + ++lowerLevelSetId) { + // put iterator to same position as the top levelset + iterators[lowerLevelSetId].goToIndicesSequential(indices); + + // if the lower surface is actually outside, i.e. its LS value + // is lower or equal + if (iterators[lowerLevelSetId].getValue() <= + value + wrappingLayerEpsilon) { + + // move neighborIterator to current position + neighborIterator.goToIndicesSequential(indices); + + Vec3D coords; + for (unsigned i = 0; i < D; ++i) { + coords[i] = indices[i] * gridDelta; + } + + Vec3D normal = {}; + T normalModulus = 0.; + for (unsigned i = 0; i < D; ++i) { + const T phiPos = neighborIterator.getNeighbor(i).getValue(); + const T phiNeg = neighborIterator.getNeighbor(i + D).getValue(); + + T diffPos = (phiPos - value) / deltaPos; + T diffNeg = (phiNeg - value) / deltaNeg; + + normal[i] = (diffNeg + diffPos) * 0.5; + normalModulus += normal[i] * normal[i]; + } + normalModulus = std::sqrt(normalModulus); + for (unsigned i = 0; i < D; ++i) + normal[i] /= normalModulus; + + T scaVel = velocities->getScalarVelocity( + coords, lowerLevelSetId, normal, + neighborIterator.getCenter().getPointId()); + auto vecVel = velocities->getVectorVelocity( + coords, lowerLevelSetId, normal, + neighborIterator.getCenter().getPointId()); + + for (unsigned i = 0; i < D; ++i) { + T tempAlpha = std::abs((scaVel + vecVel[i]) * normal[i]); + localAlphas[i] = std::max(localAlphas[i], tempAlpha); + } + + break; + } + } + } + +#pragma omp critical + { + for (unsigned i = 0; i < D; ++i) { + finalAlphas[i] = std::max(finalAlphas[i], localAlphas[i]); + } + } + } // end of parallel section + + return finalAlphas; + } + void rebuildLS() { // TODO: this function uses manhatten distances for renormalisation, // since this is the quickest. For visualisation applications, better @@ -293,14 +410,16 @@ template class Advect { currentTime = integrateTime(is, maxTimeStep); } else if (integrationScheme == IntegrationSchemeEnum::LAX_FRIEDRICHS_1ST_ORDER) { + auto alphas = findGlobalAlphas(); auto is = lsInternal::LaxFriedrichs(levelSets.back(), velocities, - dissipationAlpha, + dissipationAlpha, alphas, calculateNormalVectors); currentTime = integrateTime(is, maxTimeStep); } else if (integrationScheme == IntegrationSchemeEnum::LAX_FRIEDRICHS_2ND_ORDER) { + auto alphas = findGlobalAlphas(); auto is = lsInternal::LaxFriedrichs(levelSets.back(), velocities, - dissipationAlpha, + dissipationAlpha, alphas, calculateNormalVectors); currentTime = integrateTime(is, maxTimeStep); } else if (integrationScheme == @@ -381,7 +500,7 @@ template class Advect { auto &topDomain = levelSets.back()->getDomain(); auto &grid = levelSets.back()->getGrid(); - std::vector>> totalTempRates; + std::vector, T>>> totalTempRates; totalTempRates.resize((levelSets.back())->getNumberOfSegments()); typename PointData::ScalarDataType *voidMarkerPointer; @@ -399,16 +518,6 @@ template class Advect { } } - // The global alpha values are stored in the integration scheme. - if (integrationScheme == IntegrationSchemeEnum::LAX_FRIEDRICHS_1ST_ORDER) { - lsInternal::advect::findGlobalAlpha( - IntegrationScheme, levelSets, velocities); - } else if (integrationScheme == - IntegrationSchemeEnum::LAX_FRIEDRICHS_2ND_ORDER) { - lsInternal::advect::findGlobalAlpha( - IntegrationScheme, levelSets, velocities); - } - const bool ignoreVoidPoints = ignoreVoids; #pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments()) @@ -458,7 +567,7 @@ template class Advect { for (int currentLevelSetId = levelSets.size() - 1; currentLevelSetId >= 0; --currentLevelSetId) { - T velocity = 0; + std::pair gradNDissipation; if (!(ignoreVoidPoints && (*voidMarkerPointer)[it.getPointId()])) { // check if there is any other levelset at the same point: @@ -473,7 +582,8 @@ template class Advect { // is lower or equal if (iterators[lowerLevelSetId].getValue() <= value + wrappingLayerEpsilon) { - velocity = scheme(it.getStartIndices(), lowerLevelSetId); + gradNDissipation = + scheme(it.getStartIndices(), lowerLevelSetId); break; } } @@ -491,16 +601,17 @@ template class Advect { // if velocity is positive, set maximum time step possible without // violating the cfl condition + T velocity = gradNDissipation.first - gradNDissipation.second; if (velocity > 0.) { maxStepTime += cfl / velocity; - tempRates.push_back( - std::make_pair(velocity, -std::numeric_limits::max())); + tempRates.push_back(std::make_pair(gradNDissipation, + -std::numeric_limits::max())); break; // if velocity is 0, maximum time step is infinite } else if (velocity == 0.) { maxStepTime = std::numeric_limits::max(); - tempRates.push_back( - std::make_pair(velocity, std::numeric_limits::max())); + tempRates.push_back(std::make_pair(gradNDissipation, + std::numeric_limits::max())); break; // if the velocity is negative apply the velocity for as long as // possible without infringing on material below @@ -509,14 +620,14 @@ template class Advect { if (difference >= cfl) { maxStepTime -= cfl / velocity; - tempRates.push_back( - std::make_pair(velocity, std::numeric_limits::max())); + tempRates.push_back(std::make_pair( + gradNDissipation, std::numeric_limits::max())); break; } else { maxStepTime -= difference / velocity; // the second part of the pair indicates how far we can move // in this time step until the end of the material is reached - tempRates.push_back(std::make_pair(velocity, valueBelow)); + tempRates.push_back(std::make_pair(gradNDissipation, valueBelow)); cfl -= difference; // use new LS value for next calculations value = valueBelow; @@ -540,16 +651,20 @@ template class Advect { if (tempMaxTimeStep < maxTimeStep) maxTimeStep = tempMaxTimeStep; } - } + } // end of parallel section // reduce to one layer thickness and apply new values directly to the // domain segments --> DO NOT CHANGE SEGMENTATION HERE (true parameter) Reduce(levelSets.back(), 1, true).apply(); const bool saveVelocities = saveAdvectionVelocities; + std::vector> dissipationVectors( + levelSets.back()->getNumberOfSegments()); std::vector> velocityVectors( levelSets.back()->getNumberOfSegments()); + const bool checkDiss = checkDissipation; + #pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments()) { int p = 0; @@ -557,13 +672,13 @@ template class Advect { p = omp_get_thread_num(); #endif - typename std::vector>::const_iterator itRS = - totalTempRates[p].begin(); + auto itRS = totalTempRates[p].cbegin(); auto &segment = topDomain.getDomainSegment(p); unsigned maxId = segment.getNumberOfPoints(); if (saveVelocities) { velocityVectors[p].resize(maxId); + dissipationVectors[p].resize(maxId); } for (unsigned localId = 0; localId < maxId; ++localId) { @@ -573,17 +688,24 @@ template class Advect { // if there is a change in materials during one time step, deduct // the time taken to advect up to the end of the top material and // set the LS value to the one below - while (std::abs(itRS->second - value) < std::abs(time * itRS->first)) { - time -= std::abs((itRS->second - value) / itRS->first); + T velocity = itRS->first.first - itRS->first.second; + if (checkDiss && (itRS->first.first < 0 && velocity > 0) || + (itRS->first.first > 0 && velocity < 0)) { + velocity = 0; + } + T rate = time * velocity; + while (std::abs(itRS->second - value) < std::abs(rate)) { + time -= std::abs((itRS->second - value) / velocity); value = itRS->second; ++itRS; } // now deduct the velocity times the time step we take - value -= time * itRS->first; + value -= rate; if (saveVelocities) { - velocityVectors[p][localId] = time * itRS->first; + velocityVectors[p][localId] = rate; + dissipationVectors[p][localId] = itRS->first.second; } // this @@ -596,23 +718,24 @@ template class Advect { // advance the TempStopRates iterator by one ++itRS; } - } + } // end of parallel section if (saveVelocities) { auto &pointData = levelSets.back()->getPointData(); - // delete if already exists - if (int i = pointData.getScalarDataIndex(velocityLabel); i != -1) { - pointData.eraseScalarData(i); - } typename PointData::ScalarDataType vels; + typename PointData::ScalarDataType diss; for (unsigned i = 0; i < velocityVectors.size(); ++i) { vels.insert(vels.end(), std::make_move_iterator(velocityVectors[i].begin()), std::make_move_iterator(velocityVectors[i].end())); + diss.insert(diss.end(), + std::make_move_iterator(dissipationVectors[i].begin()), + std::make_move_iterator(dissipationVectors[i].end())); } - pointData.insertNextScalarData(std::move(vels), velocityLabel); + pointData.insertReplaceScalarData(std::move(vels), velocityLabel); + pointData.insertReplaceScalarData(std::move(diss), dissipationLabel); } return maxTimeStep; @@ -620,6 +743,7 @@ template class Advect { public: static constexpr char velocityLabel[] = "AdvectionVelocities"; + static constexpr char dissipationLabel[] = "Dissipation"; Advect() {} @@ -713,6 +837,9 @@ template class Advect { /// scaling factor for the calculated alpha values. void setDissipationAlpha(const double &a) { dissipationAlpha = a; } + // Sets the velocity to 0 if the dissipation is too high + void setCheckDissipation(bool check) { checkDissipation = check; } + /// Set whether the point data in the old LS should /// be translated to the advected LS. Defaults to true. void setUpdatePointData(bool update) { updatePointData = update; } diff --git a/include/viennals/lsEnquistOsher.hpp b/include/viennals/lsEnquistOsher.hpp index 93878384..8c2391f9 100644 --- a/include/viennals/lsEnquistOsher.hpp +++ b/include/viennals/lsEnquistOsher.hpp @@ -38,7 +38,8 @@ template class EnquistOsher { levelSet->getDomain())), calculateNormalVectors(calcNormal) {} - T operator()(const hrleVectorType &indices, int material) { + std::pair operator()(const hrleVectorType &indices, + int material) { auto &grid = levelSet->getGrid(); double gridDelta = grid.getGridDelta(); @@ -160,11 +161,9 @@ template class EnquistOsher { } } - return vel_grad; + return {vel_grad, 0.}; } - void setFinalAlphas(const hrleVectorType &) {} - void reduceTimeStepHamiltonJacobi(double &MaxTimeStep, hrleCoordType gridDelta) {} }; diff --git a/include/viennals/lsLaxFriedrichs.hpp b/include/viennals/lsLaxFriedrichs.hpp index 1ba4614c..2dcf2aaf 100644 --- a/include/viennals/lsLaxFriedrichs.hpp +++ b/include/viennals/lsLaxFriedrichs.hpp @@ -20,9 +20,9 @@ template class LaxFriedrichs { SmartPointer> levelSet; SmartPointer> velocities; hrleSparseStarIterator, order> neighborIterator; - bool calculateNormalVectors = true; const double alphaFactor = 1.0; - hrleVectorType finalAlphas; + hrleVectorType const finalAlphas; + const bool calculateNormalVectors = true; static T pow2(const T &value) { return value * value; } @@ -34,22 +34,16 @@ template class LaxFriedrichs { } LaxFriedrichs(SmartPointer> passedlsDomain, - SmartPointer> vel, double a = 1.0, - bool calcNormal = true) + SmartPointer> vel, double alpha, + hrleVectorType &alphas, bool calcNormal) : levelSet(passedlsDomain), velocities(vel), neighborIterator(hrleSparseStarIterator, order>( levelSet->getDomain())), - calculateNormalVectors(calcNormal), alphaFactor(a) { - for (int i = 0; i < 3; ++i) { - finalAlphas[i] = 1.0; - } - } - - void setFinalAlphas(const hrleVectorType &alphas) { - finalAlphas = alphas; - } + alphaFactor(alpha), finalAlphas(alphas), + calculateNormalVectors(calcNormal) {} - T operator()(const hrleVectorType &indices, int material) { + std::pair operator()(const hrleVectorType &indices, + int material) { auto &grid = levelSet->getGrid(); double gridDelta = grid.getGridDelta(); @@ -70,7 +64,6 @@ template class LaxFriedrichs { Vec3D normalVector = {}; T normalModulus = 0; - const bool calcNormals = calculateNormalVectors; for (int i = 0; i < D; i++) { // iterate over dimensions @@ -126,16 +119,16 @@ template class LaxFriedrichs { gradPos[i] = diffNeg; gradNeg[i] = diffPos; - if (calcNormals) { + if (calculateNormalVectors) { normalVector[i] = (diffNeg + diffPos) * 0.5; normalModulus += normalVector[i] * normalVector[i]; } grad += pow2((diffNeg + diffPos) * 0.5); - dissipation += finalAlphas[i] * (diffPos - diffNeg) * 0.5; + dissipation += alphaFactor * finalAlphas[i] * (diffPos - diffNeg) * 0.5; } - if (calcNormals) { + if (calculateNormalVectors) { normalModulus = std::sqrt(normalModulus); for (unsigned i = 0; i < D; ++i) { normalVector[i] /= normalModulus; @@ -165,7 +158,7 @@ template class LaxFriedrichs { } } - return totalGrad - ((totalGrad != 0.) ? alphaFactor * dissipation : 0); + return {totalGrad, ((totalGrad != 0.) ? dissipation : 0)}; } void reduceTimeStepHamiltonJacobi(double &MaxTimeStep, @@ -181,171 +174,5 @@ template class LaxFriedrichs { timeStep = alpha_maxCFL / timeStep; MaxTimeStep = std::min(timeStep, MaxTimeStep); } - - hrleVectorType - calcAlpha(const hrleVectorType &indices, int material) { - const T gridDelta = levelSet->getGrid().getGridDelta(); - - // move neighborIterator to current position - neighborIterator.goToIndicesSequential(indices); - - Vec3D coords; - for (unsigned i = 0; i < D; ++i) { - coords[i] = indices[i] * gridDelta; - } - - const T deltaPos = gridDelta; - const T deltaNeg = -gridDelta; - - Vec3D normal = {}; - T normalModulus = 0.; - const T phi0 = neighborIterator.getCenter().getValue(); - for (unsigned i = 0; i < D; ++i) { - const T phiPos = neighborIterator.getNeighbor(i).getValue(); - const T phiNeg = neighborIterator.getNeighbor(i + D).getValue(); - - T diffPos = (phiPos - phi0) / deltaPos; - T diffNeg = (phiNeg - phi0) / deltaNeg; - - normal[i] = (diffNeg + diffPos) * 0.5; - normalModulus += normal[i] * normal[i]; - } - normalModulus = std::sqrt(normalModulus); - // normalise normal vector - for (unsigned i = 0; i < D; ++i) - normal[i] /= normalModulus; - - T scaVel = velocities->getScalarVelocity( - coords, material, normal, neighborIterator.getCenter().getPointId()); - auto vecVel = velocities->getVectorVelocity( - coords, material, normal, neighborIterator.getCenter().getPointId()); - - hrleVectorType alpha(0., 0., 0.); - for (unsigned i = 0; i < D; ++i) { - T tempAlpha = std::abs((scaVel + vecVel[i]) * normal[i]); - alpha[i] = std::max(alpha[i], tempAlpha); - } - - return alpha; - } }; - -namespace advect { -using namespace viennals; -template -void findGlobalAlpha(IntegrationSchemeType &integrationScheme, - std::vector>> &levelSets, - SmartPointer> velocities) { - - auto &topDomain = levelSets.back()->getDomain(); - auto &grid = levelSets.back()->getGrid(); - - const T gridDelta = grid.getGridDelta(); - const T deltaPos = gridDelta; - const T deltaNeg = -gridDelta; - - hrleVectorType finalAlphas(0., 0., 0.); - -#pragma omp parallel num_threads((levelSets.back())->getNumberOfSegments()) - { - int p = 0; -#ifdef _OPENMP - p = omp_get_thread_num(); -#endif - - hrleVectorType localAlphas(0., 0., 0.); - - hrleVectorType startVector = - (p == 0) ? grid.getMinGridPoint() : topDomain.getSegmentation()[p - 1]; - - hrleVectorType endVector = - (p != static_cast(topDomain.getNumberOfSegments() - 1)) - ? topDomain.getSegmentation()[p] - : grid.incrementIndices(grid.getMaxGridPoint()); - - // an iterator for each level set - std::vector::DomainType>> - iterators; - for (auto it = levelSets.begin(); it != levelSets.end(); ++it) { - iterators.push_back(hrleSparseIterator::DomainType>( - (*it)->getDomain())); - } - - // neighborIterator for the top level set - hrleSparseStarIterator, order> neighborIterator(topDomain); - - for (hrleSparseIterator::DomainType> it(topDomain, - startVector); - it.getStartIndices() < endVector; ++it) { - - if (!it.isDefined() || std::abs(it.getValue()) > 0.5) - continue; - - T value = it.getValue(); - - // check if there is any other levelset at the same point: - // if yes, take the velocity of the lowest levelset - for (unsigned lowerLevelSetId = 0; lowerLevelSetId < levelSets.size(); - ++lowerLevelSetId) { - // put iterator to same position as the top levelset - auto indices = it.getStartIndices(); - iterators[lowerLevelSetId].goToIndicesSequential(indices); - - // if the lower surface is actually outside, i.e. its LS value - // is lower or equal - if (iterators[lowerLevelSetId].getValue() <= value + 1e-4) { - - // move neighborIterator to current position - neighborIterator.goToIndicesSequential(indices); - - Vec3D coords; - for (unsigned i = 0; i < D; ++i) { - coords[i] = indices[i] * gridDelta; - } - - Vec3D normal = {}; - T normalModulus = 0.; - const T phi0 = neighborIterator.getCenter().getValue(); - for (unsigned i = 0; i < D; ++i) { - const T phiPos = neighborIterator.getNeighbor(i).getValue(); - const T phiNeg = neighborIterator.getNeighbor(i + D).getValue(); - - T diffPos = (phiPos - phi0) / deltaPos; - T diffNeg = (phiNeg - phi0) / deltaNeg; - - normal[i] = (diffNeg + diffPos) * 0.5; - normalModulus += normal[i] * normal[i]; - } - normalModulus = std::sqrt(normalModulus); - for (unsigned i = 0; i < D; ++i) - normal[i] /= normalModulus; - - T scaVel = velocities->getScalarVelocity( - coords, lowerLevelSetId, normal, - neighborIterator.getCenter().getPointId()); - auto vecVel = velocities->getVectorVelocity( - coords, lowerLevelSetId, normal, - neighborIterator.getCenter().getPointId()); - - for (unsigned i = 0; i < D; ++i) { - T tempAlpha = std::abs((scaVel + vecVel[i]) * normal[i]); - localAlphas[i] = std::max(localAlphas[i], tempAlpha); - } - - break; - } - } - } - -#pragma omp critical - { - for (unsigned i = 0; i < D; ++i) { - finalAlphas[i] = std::max(finalAlphas[i], localAlphas[i]); - } - } - } // end of parallel section - - integrationScheme.setFinalAlphas(finalAlphas); -} -} // namespace advect } // namespace lsInternal diff --git a/include/viennals/lsLocalLaxFriedrichs.hpp b/include/viennals/lsLocalLaxFriedrichs.hpp index d855cdd1..65ff67da 100644 --- a/include/viennals/lsLocalLaxFriedrichs.hpp +++ b/include/viennals/lsLocalLaxFriedrichs.hpp @@ -63,11 +63,8 @@ template class LocalLaxFriedrichs { } } - void setFinalAlphas(const hrleVectorType &alphas) { - finalAlphas = alphas; - } - - T operator()(const hrleVectorType &indices, int material) { + std::pair operator()(const hrleVectorType &indices, + int material) { auto &grid = levelSet->getGrid(); double gridDelta = grid.getGridDelta(); @@ -241,9 +238,7 @@ template class LocalLaxFriedrichs { dissipation += alphaFactor * alpha[i] * (gradNeg[i] - gradPos[i]) * 0.5; } - // std::cout << neighborIterator.getCenter().getPointId() << " dissipation: - // " << dissipation << std::endl; - return totalGrad - ((totalGrad != 0.) ? dissipation : 0); + return {totalGrad, ((totalGrad != 0.) ? dissipation : 0)}; } void reduceTimeStepHamiltonJacobi(double &MaxTimeStep, diff --git a/include/viennals/lsLocalLaxFriedrichsAnalytical.hpp b/include/viennals/lsLocalLaxFriedrichsAnalytical.hpp index 5b10d8c9..c6e7891d 100644 --- a/include/viennals/lsLocalLaxFriedrichsAnalytical.hpp +++ b/include/viennals/lsLocalLaxFriedrichsAnalytical.hpp @@ -62,11 +62,8 @@ template class LocalLaxFriedrichsAnalytical { } } - void setFinalAlphas(const hrleVectorType &alphas) { - finalAlphas = alphas; - } - - T operator()(const hrleVectorType &indices, int material) { + std::pair operator()(const hrleVectorType &indices, + int material) { auto &grid = levelSet->getGrid(); double gridDelta = grid.getGridDelta(); @@ -224,9 +221,7 @@ template class LocalLaxFriedrichsAnalytical { dissipation += alpha[i] * (gradNeg[i] - gradPos[i]) * 0.5; } - // std::cout << neighborIterator.getCenter().getPointId() << " dissipation: - // " << dissipation << std::endl; - return totalGrad - ((totalGrad != 0.) ? dissipation : 0); + return {totalGrad, ((totalGrad != 0.) ? dissipation : 0)}; } void reduceTimeStepHamiltonJacobi(double &MaxTimeStep, diff --git a/include/viennals/lsLocalLocalLaxFriedrichs.hpp b/include/viennals/lsLocalLocalLaxFriedrichs.hpp index 0d7b9a3b..f512c737 100644 --- a/include/viennals/lsLocalLocalLaxFriedrichs.hpp +++ b/include/viennals/lsLocalLocalLaxFriedrichs.hpp @@ -43,11 +43,8 @@ template class LocalLocalLaxFriedrichs { } } - void setFinalAlphas(const hrleVectorType &alphas) { - finalAlphas = alphas; - } - - T operator()(const hrleVectorType &indices, int material) { + std::pair operator()(const hrleVectorType &indices, + int material) { auto &grid = levelSet->getGrid(); double gridDelta = grid.getGridDelta(); @@ -169,7 +166,7 @@ template class LocalLocalLaxFriedrichs { dissipation += alphaFactor * alpha * (gradNeg[i] - gradPos[i]) * 0.5; } - return totalGrad - ((totalGrad != 0.) ? dissipation : 0); + return {totalGrad, ((totalGrad != 0.) ? dissipation : 0)}; } void reduceTimeStepHamiltonJacobi(double &MaxTimeStep, diff --git a/include/viennals/lsPointData.hpp b/include/viennals/lsPointData.hpp index 188044e6..ab34e1b8 100644 --- a/include/viennals/lsPointData.hpp +++ b/include/viennals/lsPointData.hpp @@ -81,6 +81,46 @@ class PointData { vectorDataLabels.push_back(label); } + /// insert or replace scalar data array + void insertReplaceScalarData(const ScalarDataType &scalars, + std::string label = "Scalars") { + if (int i = getScalarDataIndex(label); i != -1) { + scalarData[i] = scalars; + } else { + insertNextScalarData(scalars, label); + } + } + + /// insert or replace scalar data array + void insertReplaceScalarData(ScalarDataType &&scalars, + std::string label = "Scalars") { + if (int i = getScalarDataIndex(label); i != -1) { + scalarData[i] = std::move(scalars); + } else { + insertNextScalarData(std::move(scalars), label); + } + } + + /// insert or replace vector data array + void insertReplaceVectorData(const VectorDataType &vectors, + std::string label = "Vectors") { + if (int i = getVectorDataIndex(label); i != -1) { + vectorData[i] = vectors; + } else { + insertNextVectorData(vectors, label); + } + } + + /// insert new vector data array + void insertReplaceVectorData(VectorDataType &&vectors, + std::string label = "Vectors") { + if (int i = getVectorDataIndex(label); i != -1) { + vectorData[i] = std::move(vectors); + } else { + insertNextVectorData(std::move(vectors), label); + } + } + /// get the number of different scalar data arrays saved unsigned getScalarDataSize() const { return scalarData.size(); } diff --git a/include/viennals/lsStencilLocalLaxFriedrichsScalar.hpp b/include/viennals/lsStencilLocalLaxFriedrichsScalar.hpp index 90105d89..38ae3f4d 100644 --- a/include/viennals/lsStencilLocalLaxFriedrichsScalar.hpp +++ b/include/viennals/lsStencilLocalLaxFriedrichsScalar.hpp @@ -143,17 +143,13 @@ template class StencilLocalLaxFriedrichsScalar { neighborIterator(hrleSparseBoxIterator>( levelSet->getDomain(), static_cast(scheme) + 1 + order)), alphaFactor(a), numStencilPoints(std::pow(2 * order + 1, D)) { - for (int i = 0; i < 3; ++i) { finalAlphas[i] = 0; } } - void setFinalAlphas(const hrleVectorType &alphas) { - finalAlphas = alphas; - } - - T operator()(const hrleVectorType &indices, int material) { + std::pair operator()(const hrleVectorType &indices, + int material) { auto &grid = levelSet->getGrid(); double gridDelta = grid.getGridDelta(); @@ -204,7 +200,7 @@ template class StencilLocalLaxFriedrichsScalar { } if (scalarVelocity == T(0)) { - return 0; + return {0, 0}; } T hamiltonian = @@ -326,7 +322,7 @@ template class StencilLocalLaxFriedrichsScalar { dissipation += maxAlpha * gradientDiff[d]; } - return hamiltonian - dissipation; + return {hamiltonian, dissipation}; } } diff --git a/include/viennals/lsToSurfaceMesh.hpp b/include/viennals/lsToSurfaceMesh.hpp index f3b6bf4e..a89f2125 100644 --- a/include/viennals/lsToSurfaceMesh.hpp +++ b/include/viennals/lsToSurfaceMesh.hpp @@ -190,7 +190,8 @@ template class ToSurfaceMesh { } } - mesh->insertNextElement(nod_numbers); // insert new surface element + if (!triangleMisformed(nod_numbers)) + mesh->insertNextElement(nod_numbers); // insert new surface element } } @@ -200,6 +201,18 @@ template class ToSurfaceMesh { newDataSourceIds); } } + +private: + static bool inline triangleMisformed( + const std::array &nodeNumbers) { + if constexpr (D == 3) { + return nodeNumbers[0] == nodeNumbers[1] || + nodeNumbers[0] == nodeNumbers[2] || + nodeNumbers[1] == nodeNumbers[2]; + } else { + return nodeNumbers[0] == nodeNumbers[1]; + } + } }; // add all template specialisations for this class diff --git a/include/viennals/lsWriteVisualizationMesh.hpp b/include/viennals/lsWriteVisualizationMesh.hpp index 5d06fa0a..456d0154 100644 --- a/include/viennals/lsWriteVisualizationMesh.hpp +++ b/include/viennals/lsWriteVisualizationMesh.hpp @@ -50,7 +50,7 @@ template class WriteVisualizationMesh { bool extractVolumeMesh = true; bool extractHullMesh = false; bool bottomRemoved = false; - static constexpr double LSEpsilon = 1e-2; + double LSEpsilon = 1e-2; /// This function removes duplicate points and agjusts the pointIDs in the /// cells @@ -226,8 +226,7 @@ template class WriteVisualizationMesh { // boundary conditions template vtkSmartPointer - LS2RectiLinearGrid(SmartPointer> levelSet, - const double LSOffset = 0., + LS2RectiLinearGrid(SmartPointer> levelSet, const double LSOffset, int infiniteMinimum = std::numeric_limits::max(), int infiniteMaximum = -std::numeric_limits::max()) { @@ -437,6 +436,8 @@ template class WriteVisualizationMesh { materialMap = passedMaterialMap; } + void setWrappingLayerEpsilon(double epsilon) { LSEpsilon = epsilon; } + void apply() { // check if level sets have enough layers for (unsigned i = 0; i < levelSets.size(); ++i) { diff --git a/pyproject.toml b/pyproject.toml index ae63c85f..06dc742f 100644 --- a/pyproject.toml +++ b/pyproject.toml @@ -3,7 +3,7 @@ requires = ["scikit-build-core", "pybind11"] build-backend = "scikit_build_core.build" [project] -version = "4.1.2" +version = "4.1.3" name = "ViennaLS" readme = "README.md" license = {file = "LICENSE"} diff --git a/tests/Advection/Advection.cpp b/tests/Advection/Advection.cpp index a555aaa1..ff69423c 100644 --- a/tests/Advection/Advection.cpp +++ b/tests/Advection/Advection.cpp @@ -45,61 +45,72 @@ class velocityField : public ls::VelocityField { int main() { constexpr int D = 3; - omp_set_num_threads(4); - - double gridDelta = 0.4999999; - - auto sphere1 = ls::SmartPointer>::New(gridDelta); - - double origin[3] = {5., 0., 0.}; - double radius = 7.3; - - ls::MakeGeometry( - sphere1, ls::SmartPointer>::New(origin, radius)) - .apply(); - - // { - // std::cout << "Extracting..." << std::endl; - // auto mesh = ls::SmartPointer>::New(); - // ls::ToSurfaceMesh(sphere1, mesh).apply(); - // ls::VTKWriter(mesh, "before.vtk").apply(); - // } - - // instantiate velocities - auto velocities = ls::SmartPointer::New(); - - // std::cout << "Advecting" << std::endl; - - // Fill point data with original point IDs to see how LS changed - { - typename ls::PointData::ScalarDataType pointIDs( - sphere1->getNumberOfPoints()); - std::iota(std::begin(pointIDs), std::end(pointIDs), 0); - sphere1->getPointData().insertNextScalarData(pointIDs, "originalIDs"); - } - - ls::Advect advectionKernel; - advectionKernel.insertNextLevelSet(sphere1); - advectionKernel.setVelocityField(velocities); - advectionKernel.setIntegrationScheme( - ls::IntegrationSchemeEnum::ENGQUIST_OSHER_1ST_ORDER); - advectionKernel.setSaveAdvectionVelocities(true); - - double time = 0.; - for (unsigned i = 0; time < 2.0 && i < 1e3; ++i) { - advectionKernel.apply(); - time += advectionKernel.getAdvectedTime(); - - // std::string fileName = std::to_string(i) + ".vtp"; - // auto mesh = ls::SmartPointer>::New(); - // ls::ToMesh(sphere1, mesh).apply(); - // ls::VTKWriter(mesh, "points_" + fileName).apply(); - // ls::ToSurfaceMesh(sphere1, mesh).apply(); - // ls::VTKWriter(mesh, "surface_" + fileName).apply(); + omp_set_num_threads(1); + + double gridDelta = 0.6999999; + + std::vector integrationSchemes = { + // ls::IntegrationSchemeEnum::ENGQUIST_OSHER_1ST_ORDER, + // ls::IntegrationSchemeEnum::ENGQUIST_OSHER_2ND_ORDER, + ls::IntegrationSchemeEnum::LAX_FRIEDRICHS_1ST_ORDER}; + // ls::IntegrationSchemeEnum::LAX_FRIEDRICHS_2ND_ORDER, + // ls::IntegrationSchemeEnum::LOCAL_LOCAL_LAX_FRIEDRICHS_1ST_ORDER, + // ls::IntegrationSchemeEnum::LOCAL_LOCAL_LAX_FRIEDRICHS_2ND_ORDER, + // ls::IntegrationSchemeEnum::LOCAL_LAX_FRIEDRICHS_1ST_ORDER, + // ls::IntegrationSchemeEnum::LOCAL_LAX_FRIEDRICHS_2ND_ORDER}; + + for (auto integrationScheme : integrationSchemes) { + auto sphere1 = ls::SmartPointer>::New(gridDelta); + + double origin[3] = {5., 0., 0.}; + double radius = 7.3; + + ls::MakeGeometry( + sphere1, ls::SmartPointer>::New(origin, radius)) + .apply(); + + // { + // std::cout << "Extracting..." << std::endl; + // auto mesh = ls::SmartPointer>::New(); + // ls::ToSurfaceMesh(sphere1, mesh).apply(); + // ls::VTKWriter(mesh, "before.vtk").apply(); + // } + + // instantiate velocities + auto velocities = ls::SmartPointer::New(); + + // std::cout << "Advecting" << std::endl; + + // Fill point data with original point IDs to see how LS changed + { + typename ls::PointData::ScalarDataType pointIDs( + sphere1->getNumberOfPoints()); + std::iota(std::begin(pointIDs), std::end(pointIDs), 0); + sphere1->getPointData().insertNextScalarData(pointIDs, "originalIDs"); + } + + ls::Advect advectionKernel; + advectionKernel.insertNextLevelSet(sphere1); + advectionKernel.setVelocityField(velocities); + advectionKernel.setIntegrationScheme(integrationScheme); + advectionKernel.setSaveAdvectionVelocities(true); + + double time = 0.; + for (unsigned i = 0; time < 1.0 && i < 1e2; ++i) { + advectionKernel.apply(); + time += advectionKernel.getAdvectedTime(); + + std::string fileName = std::to_string(i) + ".vtp"; + auto mesh = ls::SmartPointer>::New(); + ls::ToMesh(sphere1, mesh).apply(); + ls::VTKWriter(mesh, "points_" + fileName).apply(); + ls::ToSurfaceMesh(sphere1, mesh).apply(); + ls::VTKWriter(mesh, "surface_" + fileName).apply(); + } + + LSTEST_ASSERT_VALID_LS(sphere1, double, D) } - LSTEST_ASSERT_VALID_LS(sphere1, double, D) - // std::cout << sphere1->getNumberOfPoints() << std::endl; // Prune(sphere1).apply();