From 33a5cb698d0553c75ab99deb6891b904baa90a42 Mon Sep 17 00:00:00 2001 From: Lukas Werner <lks.werner@fau.de> Date: Mon, 12 Mar 2018 15:31:14 +0100 Subject: [PATCH] Hopefully satisfy CI --- apps/tutorials/pe/01_ConfinedGas.cpp | 2 +- src/pe/ccd/HashGrids.h | 40 +++++----- src/pe/raytracing/Intersects.h | 4 + src/pe/raytracing/Ray.h | 4 +- src/pe/raytracing/Raytracer.cpp | 2 +- src/pe/raytracing/ShadingFunctions.h | 10 +++ tests/pe/Raytracing.cpp | 105 ++++++++++++--------------- 7 files changed, 82 insertions(+), 85 deletions(-) diff --git a/apps/tutorials/pe/01_ConfinedGas.cpp b/apps/tutorials/pe/01_ConfinedGas.cpp index 20e41ea70..52795b4d3 100644 --- a/apps/tutorials/pe/01_ConfinedGas.cpp +++ b/apps/tutorials/pe/01_ConfinedGas.cpp @@ -138,7 +138,7 @@ int main( int argc, char ** argv ) if (cfg == NULL) { WALBERLA_ABORT("raytracer needs a working config"); } - Raytracer raytracer(forest, storageID, globalBodyStorage, cfg->getBlock("raytracing")); + Raytracer raytracer(forest, storageID, globalBodyStorage, ccdID, cfg->getBlock("raytracing")); WALBERLA_LOG_INFO_ON_ROOT("*** INTEGRATOR ***"); //! [Integrator] diff --git a/src/pe/ccd/HashGrids.h b/src/pe/ccd/HashGrids.h index 305a25906..c13b1b4ba 100644 --- a/src/pe/ccd/HashGrids.h +++ b/src/pe/ccd/HashGrids.h @@ -189,11 +189,11 @@ private: BodyID getRayIntersectingBody(const raytracing::Ray& ray, const AABB& blockAABB, real_t& t, Vec3& n) const; template<typename BodyTuple> - const BodyID getBodyIntersectionForBlockCell(const Vector3<int32_t>& blockCell, - const int8_t cellNormalAxis, const int8_t cellNormalDir, - const AABB& blockAABB, - const raytracing::Ray& ray, - real_t& t_closest, Vec3& n_closest) const; + BodyID getBodyIntersectionForBlockCell(const Vector3<int32_t>& blockCell, + const int8_t cellNormalAxis, const int8_t cellNormalDir, + const AABB& blockAABB, + const raytracing::Ray& ray, + real_t& t_closest, Vec3& n_closest) const; void clear(); //@} @@ -512,11 +512,11 @@ void HashGrids::HashGrid::processBodies( BodyID* bodies, size_t bodyCount, Conta * cell into the intersected cell (and thus, possibly intersect with the ray as well). */ template<typename BodyTuple> -const BodyID HashGrids::HashGrid::getBodyIntersectionForBlockCell(const Vector3<int32_t>& blockCell, - const int8_t cellNormalAxis, const int8_t cellNormalDir, - const AABB& blockAABB, - const raytracing::Ray& ray, - real_t& t_closest, Vec3& n_closest) const { +BodyID HashGrids::HashGrid::getBodyIntersectionForBlockCell(const Vector3<int32_t>& blockCell, + const int8_t cellNormalAxis, const int8_t cellNormalDir, + const AABB& blockAABB, + const raytracing::Ray& ray, + real_t& t_closest, Vec3& n_closest) const { real_t t_local; Vec3 n_local; BodyID body = NULL; @@ -647,14 +647,14 @@ BodyID HashGrids::HashGrid::getRayIntersectingBody(const raytracing::Ray& ray, c (stepZ >= 0) ? (firstCell[2]+1)*cellSpan_-firstPoint[2] : firstPoint[2]-firstCell[2]*cellSpan_); // tMax: distance along the ray to the next cell change in the axis direction - real_t tMaxX = (ray.xDir() != 0) ? std::abs(near[0]*ray.xInvDir()) : inf; - real_t tMaxY = (ray.yDir() != 0) ? std::abs(near[1]*ray.yInvDir()) : inf; - real_t tMaxZ = (ray.zDir() != 0) ? std::abs(near[2]*ray.zInvDir()) : inf; + real_t tMaxX = (!realIsEqual(ray.xDir(), 0)) ? std::abs(near[0]*ray.xInvDir()) : inf; + real_t tMaxY = (!realIsEqual(ray.yDir(), 0)) ? std::abs(near[1]*ray.yInvDir()) : inf; + real_t tMaxZ = (!realIsEqual(ray.zDir(), 0)) ? std::abs(near[2]*ray.zInvDir()) : inf; // tDelta: how far along the ray must be moved to encounter a new cell in the specified axis direction - real_t tDeltaX = (ray.xDir() != 0) ? std::abs(cellSpan_*ray.xInvDir()) : inf; - real_t tDeltaY = (ray.yDir() != 0) ? std::abs(cellSpan_*ray.yInvDir()) : inf; - real_t tDeltaZ = (ray.zDir() != 0) ? std::abs(cellSpan_*ray.zInvDir()) : inf; + real_t tDeltaX = (!realIsEqual(ray.xDir(), 0)) ? std::abs(cellSpan_*ray.xInvDir()) : inf; + real_t tDeltaY = (!realIsEqual(ray.yDir(), 0)) ? std::abs(cellSpan_*ray.yInvDir()) : inf; + real_t tDeltaZ = (!realIsEqual(ray.zDir(), 0)) ? std::abs(cellSpan_*ray.zInvDir()) : inf; Vector3<int32_t> currentCell = firstCell; @@ -675,14 +675,11 @@ BodyID HashGrids::HashGrid::getRayIntersectingBody(const raytracing::Ray& ray, c int8_t blockCellNormalAxis; int8_t blockCellNormalDir; - static int earlyCutoffs = 0; - while (true) { if (tMaxX < tMaxY) { if (tMaxX < tMaxZ) { #if !defined(HASHGRIDS_DISABLE_EARLY_CUTOFF) if (tRayOriginToGrid+tMaxX-tDeltaX > t_closest) { - earlyCutoffs++; break; } #endif @@ -696,7 +693,6 @@ BodyID HashGrids::HashGrid::getRayIntersectingBody(const raytracing::Ray& ray, c } else { #if !defined(HASHGRIDS_DISABLE_EARLY_CUTOFF) if (tRayOriginToGrid+tMaxZ-tDeltaZ > t_closest) { - earlyCutoffs++; break; } #endif @@ -712,7 +708,6 @@ BodyID HashGrids::HashGrid::getRayIntersectingBody(const raytracing::Ray& ray, c if (tMaxY < tMaxZ) { #if !defined(HASHGRIDS_DISABLE_EARLY_CUTOFF) if (tRayOriginToGrid+tMaxY-tDeltaY > t_closest) { - earlyCutoffs++; break; } #endif @@ -726,7 +721,6 @@ BodyID HashGrids::HashGrid::getRayIntersectingBody(const raytracing::Ray& ray, c } else { #if !defined(HASHGRIDS_DISABLE_EARLY_CUTOFF) if (tRayOriginToGrid+tMaxZ-tDeltaZ > t_closest) { - earlyCutoffs++; break; } #endif @@ -760,7 +754,7 @@ BodyID HashGrids::HashGrid::getRayIntersectingBody(const raytracing::Ray& ray, c * \return BodyID of closest body, NULL if none found. */ template<typename BodyTuple> -BodyID HashGrids::getClosestBodyIntersectingWithRay(const raytracing::Ray& ray, const AABB& blockAABB, +BodyID HashGrids::getClosestBodyIntersectingWithRay(const raytracing::Ray& ray, const AABB& blockAABB, real_t& t, Vec3& n) { real_t inf = std::numeric_limits<real_t>::max(); diff --git a/src/pe/raytracing/Intersects.h b/src/pe/raytracing/Intersects.h index 5e9df8397..20fb957b6 100644 --- a/src/pe/raytracing/Intersects.h +++ b/src/pe/raytracing/Intersects.h @@ -290,6 +290,10 @@ inline bool intersects(const CapsuleID capsule, const Ray& ray, real_t& t, Vec3& } inline bool intersects(const BodyID body, const Ray& ray, real_t& t, Vec3& n) { + WALBERLA_UNUSED(body); + WALBERLA_UNUSED(ray); + WALBERLA_UNUSED(t); + WALBERLA_UNUSED(n); WALBERLA_ABORT("This ray - body intersection test is not implemented yet!"); return false; } diff --git a/src/pe/raytracing/Ray.h b/src/pe/raytracing/Ray.h index 82d659dcc..a95948053 100644 --- a/src/pe/raytracing/Ray.h +++ b/src/pe/raytracing/Ray.h @@ -78,7 +78,7 @@ public: /*!\brief Returns the normalized direction vector of the ray for a given axis. */ inline real_t getDirection (size_t axis) const { - WALBERLA_ASSERT(axis >= 0 && axis <= 2, "No valid axis index passed."); + WALBERLA_ASSERT(axis <= 2, "No valid axis index passed."); return direction_[axis]; } @@ -108,7 +108,7 @@ public: /*!\brief Returns the inverse of the direction vector of the ray for a given axis. */ - inline const real_t getInvDirection (size_t axis) const { + inline real_t getInvDirection (size_t axis) const { WALBERLA_ASSERT(axis >= 0 && axis <= 2, "No valid axis index passed."); return inv_direction_[axis]; } diff --git a/src/pe/raytracing/Raytracer.cpp b/src/pe/raytracing/Raytracer.cpp index 992734640..5b9be64b6 100644 --- a/src/pe/raytracing/Raytracer.cpp +++ b/src/pe/raytracing/Raytracer.cpp @@ -199,7 +199,7 @@ void Raytracer::setupView_() { */ void Raytracer::setupFilenameRankWidth_() { int numProcesses = mpi::MPIManager::instance()->numProcesses(); - filenameRankWidth_ = uint8_c(log10(numProcesses))+1; + filenameRankWidth_ = uint8_c(log10(numProcesses)+1); } void Raytracer::setupMPI_() { diff --git a/src/pe/raytracing/ShadingFunctions.h b/src/pe/raytracing/ShadingFunctions.h index 78fe48cff..e8b440f24 100644 --- a/src/pe/raytracing/ShadingFunctions.h +++ b/src/pe/raytracing/ShadingFunctions.h @@ -61,6 +61,7 @@ inline ShadingParameters defaultShadingParams (const BodyID body) { } inline ShadingParameters whiteShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(1, 1, 1), Color(0.9, 0.9, 0.9), Color(0, 0, 0), @@ -69,6 +70,7 @@ inline ShadingParameters whiteShadingParams (const BodyID body) { } inline ShadingParameters blackShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(0, 0, 0), Color(0, 0, 0), Color(0.1, 0.1, 0.1), @@ -77,6 +79,7 @@ inline ShadingParameters blackShadingParams (const BodyID body) { } inline ShadingParameters lightGreyShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(0.82, 0.82, 0.82), Color(0.5, 0.5, 0.5), Color(0, 0, 0), @@ -85,6 +88,7 @@ inline ShadingParameters lightGreyShadingParams (const BodyID body) { } inline ShadingParameters greyShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(0.5, 0.5, 0.5), Color(0.4, 0.4, 0.4), Color(0.1, 0.1, 0.1), @@ -93,6 +97,7 @@ inline ShadingParameters greyShadingParams (const BodyID body) { } inline ShadingParameters darkGreyShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(0.2, 0.2, 0.2), Color(0.06, 0.06, 0.06), Color(0.1, 0.1, 0.1), @@ -101,6 +106,7 @@ inline ShadingParameters darkGreyShadingParams (const BodyID body) { } inline ShadingParameters redShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(1, 0, 0), Color(0.5, 0, 0), Color(0.1, 0.1, 0.1), @@ -109,6 +115,7 @@ inline ShadingParameters redShadingParams (const BodyID body) { } inline ShadingParameters greenShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(0, 0.72, 0), Color(0, 0.41, 0), Color(0.1, 0.1, 0.1), @@ -117,6 +124,7 @@ inline ShadingParameters greenShadingParams (const BodyID body) { } inline ShadingParameters blueShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(0.15, 0.44, 0.91), Color(0, 0, 0.4), Color(0.1, 0.1, 0.1), @@ -125,6 +133,7 @@ inline ShadingParameters blueShadingParams (const BodyID body) { } inline ShadingParameters yellowShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(1, 0.96, 0), Color(0.5, 0.48, 0), Color(0, 0, 0), @@ -133,6 +142,7 @@ inline ShadingParameters yellowShadingParams (const BodyID body) { } inline ShadingParameters violetShadingParams (const BodyID body) { + WALBERLA_UNUSED(body); ShadingParameters s(Color(0.6, 0, 0.9), Color(0.5, 0, 0.8), Color(0, 0, 0), diff --git a/tests/pe/Raytracing.cpp b/tests/pe/Raytracing.cpp index dd40130f8..09c616b37 100644 --- a/tests/pe/Raytracing.cpp +++ b/tests/pe/Raytracing.cpp @@ -233,7 +233,7 @@ ShadingParameters customBodyToShadingParams(const BodyID body) { } } -void RaytracerTest(Raytracer::Algorithm raytracingAlgorithm = Raytracer::RAYTRACE_HASHGRIDS) { +void RaytracerTest(Raytracer::Algorithm raytracingAlgorithm = Raytracer::RAYTRACE_HASHGRIDS, const std::string& outputFolder = ".") { WALBERLA_LOG_INFO("Raytracer"); shared_ptr<BodyStorage> globalBodyStorage = make_shared<BodyStorage>(); shared_ptr<BlockForest> forest = createBlockForest(AABB(0,0,0,10,10,10), Vec3(1,1,1), Vec3(false, false, false)); @@ -287,7 +287,7 @@ void RaytracerTest(Raytracer::Algorithm raytracingAlgorithm = Raytracer::RAYTRAC //raytracer.setTBufferOutputDirectory("tbuffer"); //raytracer.setTBufferOutputEnabled(true); - raytracer.setImageOutputDirectory("image"); + raytracer.setImageOutputDirectory(outputFolder); raytracer.setImageOutputEnabled(true); //raytracer.setLocalImageOutputEnabled(true); @@ -336,7 +336,7 @@ ShadingParameters customSpheresBodyToShadingParams(const BodyID body) { } } -void RaytracerSpheresTest(Raytracer::Algorithm raytracingAlgorithm = Raytracer::RAYTRACE_HASHGRIDS) { +void RaytracerSpheresTest(Raytracer::Algorithm raytracingAlgorithm = Raytracer::RAYTRACE_HASHGRIDS, const std::string& outputFolder = ".") { WALBERLA_LOG_INFO("Raytracer Spheres Scene"); shared_ptr<BodyStorage> globalBodyStorage = make_shared<BodyStorage>(); shared_ptr<BlockForest> forest = createBlockForest(AABB(0,0,0,10,10,10), Vec3(1,1,1), Vec3(false, false, false)); @@ -373,7 +373,7 @@ void RaytracerSpheresTest(Raytracer::Algorithm raytracingAlgorithm = Raytracer:: } } - raytracer.setImageOutputDirectory("image"); + raytracer.setImageOutputDirectory(outputFolder); raytracer.setImageOutputEnabled(true); raytracer.setRaytracingAlgorithm(raytracingAlgorithm); @@ -389,7 +389,7 @@ ShadingParameters customHashGridsBodyToShadingParams(const BodyID body) { } -void HashGridsTest(Raytracer::Algorithm raytracingAlgorithm, +void HashGridsTest(Raytracer::Algorithm raytracingAlgorithm, const std::string& outputFolder, size_t boxes, size_t capsules, size_t spheres, size_t numberOfViews = 1, real_t boxLenMin = 0.1, real_t boxLenMax = 0.2, bool boxRotation = false, real_t capRadiusMin = 0.1, real_t capRadiusMax = 0.2, real_t capLenMin = 0.1, real_t capLenMax = 0.3, @@ -539,17 +539,7 @@ void HashGridsTest(Raytracer::Algorithm raytracingAlgorithm, Color(0.2,0.2,0.2), real_t(2), customHashGridsBodyToShadingParams); - switch (raytracingAlgorithm) { - case Raytracer::RAYTRACE_NAIVE: - raytracer.setImageOutputDirectory("image/naive"); - break; - case Raytracer::RAYTRACE_HASHGRIDS: - raytracer.setImageOutputDirectory("image/hashgrids"); - break; - case Raytracer::RAYTRACE_COMPARE_BOTH: - raytracer.setImageOutputDirectory("image/comparison"); - break; - } + raytracer.setImageOutputDirectory(outputFolder); raytracer.setImageOutputEnabled(true); raytracer.setFilenameTimestepWidth(12); WALBERLA_LOG_INFO("output #" << i << " to: " << (boxes*100000000 + capsules*10000 + spheres) << " in " << raytracer.getImageOutputDirectory()); @@ -571,7 +561,7 @@ ShadingParameters customArtifactsBodyToShadingParams(const BodyID body) { return defaultShadingParams(body); } -void raytraceArtifactsForest(Raytracer::Algorithm raytracingAlgorithm, +void raytraceArtifactsForest(Raytracer::Algorithm raytracingAlgorithm, const std::string& outputFolder, const shared_ptr<BlockStorage> forest, const BlockDataID storageID, const shared_ptr<BodyStorage> globalBodyStorage, const BlockDataID ccdID, @@ -594,17 +584,7 @@ void raytraceArtifactsForest(Raytracer::Algorithm raytracingAlgorithm, Color(0.2,0.2,0.2), real_t(2), customArtifactsBodyToShadingParams); - switch (raytracingAlgorithm) { - case Raytracer::RAYTRACE_NAIVE: - raytracer.setImageOutputDirectory("image/naive"); - break; - case Raytracer::RAYTRACE_HASHGRIDS: - raytracer.setImageOutputDirectory("image/hashgrids"); - break; - case Raytracer::RAYTRACE_COMPARE_BOTH: - raytracer.setImageOutputDirectory("image/comparison"); - break; - } + raytracer.setImageOutputDirectory(outputFolder); raytracer.setImageOutputEnabled(true); raytracer.setFilenameTimestepWidth(timestepWidth); raytracer.setRaytracingAlgorithm(raytracingAlgorithm); @@ -617,7 +597,8 @@ void raytraceArtifactsForest(Raytracer::Algorithm raytracingAlgorithm, } } -void HashGridsArtifactsTest(Raytracer::Algorithm raytracingAlgorithm, size_t boxes, real_t boxLenMin = 0.1, real_t boxLenMax = 0.2) { +void HashGridsArtifactsTest(Raytracer::Algorithm raytracingAlgorithm, const std::string& outputFolder, + size_t boxes, real_t boxLenMin = 0.1, real_t boxLenMax = 0.2) { WALBERLA_LOG_INFO_ON_ROOT("HashGrids Artifacts Test - In negative Z direction"); WALBERLA_LOG_INFO(" Generating " << boxes << " boxes"); @@ -651,13 +632,14 @@ void HashGridsArtifactsTest(Raytracer::Algorithm raytracingAlgorithm, size_t box WALBERLA_CHECK(box_ != NULL); } - raytraceArtifactsForest(raytracingAlgorithm, + raytraceArtifactsForest(raytracingAlgorithm, outputFolder, forest, storageID, globalBodyStorage, ccdID, Vec3(2, 2, 7), Vec3(2, 2, 4), Vec3(0,1,0), boxes, 3); } -void HashGridsFromNegativeArtifactsTest(Raytracer::Algorithm raytracingAlgorithm, size_t boxes, real_t boxLenMin = 0.1, real_t boxLenMax = 0.2) { +void HashGridsFromNegativeArtifactsTest(Raytracer::Algorithm raytracingAlgorithm, const std::string& outputFolder, + size_t boxes, real_t boxLenMin = 0.1, real_t boxLenMax = 0.2) { WALBERLA_LOG_INFO_ON_ROOT("HashGrids Artifacts Test - In positive Z direction"); WALBERLA_LOG_INFO_ON_ROOT(" Generating " << boxes << " boxes"); @@ -695,13 +677,14 @@ void HashGridsFromNegativeArtifactsTest(Raytracer::Algorithm raytracingAlgorithm WALBERLA_CHECK(box_ != NULL); } - raytraceArtifactsForest(raytracingAlgorithm, + raytraceArtifactsForest(raytracingAlgorithm, outputFolder, forest, storageID, globalBodyStorage, ccdID, Vec3(2, 2, -3), Vec3(2, 2, 0), Vec3(0,1,0), boxes, 4); } -void HashGridsFromNegativeXArtifactsTest(Raytracer::Algorithm raytracingAlgorithm, size_t boxes, real_t boxLenMin = 0.1, real_t boxLenMax = 0.2) { +void HashGridsFromNegativeXArtifactsTest(Raytracer::Algorithm raytracingAlgorithm, const std::string& outputFolder, + size_t boxes, real_t boxLenMin = 0.1, real_t boxLenMax = 0.2) { WALBERLA_LOG_INFO_ON_ROOT("HashGrids Artifacts Test - In positive X direction"); WALBERLA_LOG_INFO_ON_ROOT(" Generating " << boxes << " boxes"); @@ -737,7 +720,7 @@ void HashGridsFromNegativeXArtifactsTest(Raytracer::Algorithm raytracingAlgorith WALBERLA_CHECK(box_ != NULL); } - raytraceArtifactsForest(raytracingAlgorithm, + raytraceArtifactsForest(raytracingAlgorithm, outputFolder, forest, storageID, globalBodyStorage, ccdID, Vec3(-3, 2, 2), Vec3(0, 2, 2), Vec3(0,0,1), boxes, 6); @@ -748,7 +731,7 @@ Vec3 minCornerToGpos(const Vec3& minCorner, real_t lengths) { return minCorner + Vec3(lengths/2, lengths/2, lengths/2); } -void HashGridsTestScene(Raytracer::Algorithm raytracingAlgorithm = Raytracer::RAYTRACE_HASHGRIDS) { +void HashGridsTestScene(Raytracer::Algorithm raytracingAlgorithm = Raytracer::RAYTRACE_HASHGRIDS, const std::string& outputFolder = ".") { WALBERLA_LOG_INFO_ON_ROOT("HashGrids Test Scene"); shared_ptr<BodyStorage> globalBodyStorage = make_shared<BodyStorage>(); @@ -842,18 +825,7 @@ void HashGridsTestScene(Raytracer::Algorithm raytracingAlgorithm = Raytracer::RA Color(0.2,0.2,0.2), real_t(2)); - switch (raytracingAlgorithm) { - case Raytracer::RAYTRACE_NAIVE: - raytracer.setImageOutputDirectory("image/naive"); - break; - case Raytracer::RAYTRACE_HASHGRIDS: - raytracer.setImageOutputDirectory("image/hashgrids"); - break; - case Raytracer::RAYTRACE_COMPARE_BOTH: - raytracer.setImageOutputDirectory("image/comparison"); - break; - } - + raytracer.setImageOutputDirectory(outputFolder); raytracer.setRaytracingAlgorithm(raytracingAlgorithm); raytracer.setImageOutputEnabled(true); raytracer.setFilenameTimestepWidth(1); @@ -878,44 +850,61 @@ int main( int argc, char** argv ) math::seedRandomGenerator( static_cast<unsigned int>(1337 * mpi::MPIManager::instance()->worldRank()) ); Raytracer::Algorithm algorithm = Raytracer::RAYTRACE_COMPARE_BOTH; + bool outputToFoldersEnabled = false; + + std::string outputFolder = "."; + if (outputToFoldersEnabled) { + switch (algorithm) { + case Raytracer::RAYTRACE_NAIVE: + outputFolder= "image/naive"; + break; + case Raytracer::RAYTRACE_HASHGRIDS: + outputFolder = "image/hashgrids"; + break; + case Raytracer::RAYTRACE_COMPARE_BOTH: + outputFolder = "image/comparison"; + break; + } + } //SphereIntersectsTest(); //PlaneIntersectsTest(); //BoxIntersectsTest(); //AABBIntersectsTest(); //CapsuleIntersectsTest(); - //RaytracerTest(algorithm); - //RaytracerSpheresTest(algorithm); + //RaytracerTest(algorithm, outputFolder); + //RaytracerSpheresTest(algorithm, outputFolder); - HashGridsTestScene(algorithm); + HashGridsTestScene(algorithm, outputFolder); std::vector<size_t> boxes = {127, 70, 20, 150}; std::vector<size_t> capsules = {127, 60, 140, 100}; std::vector<size_t> spheres = {0, 50, 40, 120}; for (size_t i = 0; i < boxes.size(); ++i) { - HashGridsTest(algorithm, boxes[i], capsules[i], spheres[i], 1); + HashGridsTest(algorithm, outputFolder, boxes[i], capsules[i], spheres[i], 1); } - //HashGridsTest(algorithm, + //HashGridsTest(algorithm, outputFolder, // 60, 60, 3, 1, // 0.1, 0.3, true, // 0.1, 0.2, 0.1, 0.2, // 0.5, 0.6); - //HashGridsTest(algorithm, + //HashGridsTest(algorithm, outputFolder, // 750, 0, 0, 1, // 0.2, 0.3); - //HashGridsTest(algorithm, + //HashGridsTest(algorithm, outputFolder, // 400, 0, 0, 1, // 0.1, 0.3); - HashGridsArtifactsTest(algorithm, 750, 0.2, 0.3); - HashGridsFromNegativeArtifactsTest(algorithm, 750, 0.2, 0.3); - HashGridsFromNegativeXArtifactsTest(algorithm, 750, 0.2, 0.3); + HashGridsArtifactsTest(algorithm, outputFolder, 750, 0.2, 0.3); + HashGridsFromNegativeArtifactsTest(algorithm, outputFolder, 750, 0.2, 0.3); + HashGridsFromNegativeXArtifactsTest(algorithm, outputFolder, 750, 0.2, 0.3); - //HashGridsTest(algorithm, 9999, 0, 4000, 1, + //HashGridsTest(algorithm, outputFolder, + // 9999, 0, 4000, 1, // 0.1, 0.2); return EXIT_SUCCESS; -- GitLab