From f4b9a8d6ddeb254629db50eb508f6d32c237faa6 Mon Sep 17 00:00:00 2001 From: Sebastian Eibl <sebastian.eibl@fau.de> Date: Wed, 13 Jun 2018 13:21:15 +0200 Subject: [PATCH] made marshalling functions use smart pointers --- src/mesh/pe/Types.h | 3 ++ src/mesh/pe/communication/ConvexPolyhedron.h | 13 +++++---- src/pe/communication/DynamicMarshalling.h | 6 ++-- src/pe/communication/Instantiate.h | 2 +- src/pe/communication/ParseMessage.h | 9 +++--- src/pe/communication/rigidbody/Box.h | 13 +++++---- src/pe/communication/rigidbody/Capsule.h | 13 +++++---- src/pe/communication/rigidbody/Ellipsoid.h | 13 +++++---- src/pe/communication/rigidbody/Sphere.h | 13 +++++---- src/pe/communication/rigidbody/Squirmer.h | 13 +++++---- src/pe/communication/rigidbody/Union.h | 29 ++++++++++--------- src/pe/rigidbody/StorageDataHandling.h | 6 ++-- tests/pe/Marshalling.cpp | 30 ++++++++++++++++---- 13 files changed, 96 insertions(+), 67 deletions(-) diff --git a/src/mesh/pe/Types.h b/src/mesh/pe/Types.h index 842ba1132..6d94fad9c 100644 --- a/src/mesh/pe/Types.h +++ b/src/mesh/pe/Types.h @@ -19,6 +19,8 @@ // //====================================================================================================================== +#include <memory> + namespace walberla { namespace mesh { namespace pe { @@ -28,6 +30,7 @@ class ConvexPolyhedron; typedef ConvexPolyhedron ConvexPolyhedronType; //!< Type of the convex polyhedron geometric primitive. typedef ConvexPolyhedron* ConvexPolyhedronID; //!< Handle for a convexpolyhedron primitive. typedef const ConvexPolyhedron* ConstConvexPolyhedronID; //!< Handle for a constant convex polyhedron primitive. +using ConvexPolyhedronPtr = std::unique_ptr<ConvexPolyhedron>; } // namespace pe } // namespace mesh diff --git a/src/mesh/pe/communication/ConvexPolyhedron.h b/src/mesh/pe/communication/ConvexPolyhedron.h index ef9bbb21c..df14bb176 100644 --- a/src/mesh/pe/communication/ConvexPolyhedron.h +++ b/src/mesh/pe/communication/ConvexPolyhedron.h @@ -69,16 +69,17 @@ namespace pe { namespace communication { template<> -inline mesh::pe::ConvexPolyhedronID instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, mesh::pe::ConvexPolyhedronID& newBody ) +inline mesh::pe::ConvexPolyhedronPtr instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, mesh::pe::ConvexPolyhedronID& newBody ) { mesh::pe::ConvexPolyhedronParameters subobjparam; unmarshal( buffer, subobjparam ); correctBodyPosition(domain, block.center(), subobjparam.gpos_); - newBody = new mesh::pe::ConvexPolyhedron( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.mesh_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); - newBody->setLinearVel( subobjparam.v_ ); - newBody->setAngularVel( subobjparam.w_ ); - newBody->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); - return newBody; + auto cp = std::make_unique<mesh::pe::ConvexPolyhedron>( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.mesh_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); + cp->setLinearVel( subobjparam.v_ ); + cp->setAngularVel( subobjparam.w_ ); + cp->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); + newBody = static_cast<mesh::pe::ConvexPolyhedronID>(cp.get()); + return cp; } } // namespace communication diff --git a/src/pe/communication/DynamicMarshalling.h b/src/pe/communication/DynamicMarshalling.h index e2d2c4a6c..032708bff 100644 --- a/src/pe/communication/DynamicMarshalling.h +++ b/src/pe/communication/DynamicMarshalling.h @@ -91,7 +91,7 @@ private: , block_(block) {} template< typename BodyType > - BodyID operator()( BodyType* bd) { return instantiate( buffer_, domain_, block_, bd ); } + BodyPtr operator()( BodyType* bd) { return instantiate( buffer_, domain_, block_, bd ); } }; public: @@ -105,10 +105,10 @@ public: * The rigid body is casted dynamically to its original type and then marshalled. For recognition * an identifying tag is prepended. */ - static BodyID execute(mpi::RecvBuffer& buffer, const id_t typeID, const math::AABB& domain, const math::AABB& block) + static BodyPtr execute(mpi::RecvBuffer& buffer, const id_t typeID, const math::AABB& domain, const math::AABB& block) { UnmarshalFunctor func(buffer, domain, block); - return SingleCast<BodyTypeTuple, UnmarshalFunctor, BodyID>::execute (typeID, func); + return SingleCast<BodyTypeTuple, UnmarshalFunctor, BodyPtr>::execute (typeID, func); } }; diff --git a/src/pe/communication/Instantiate.h b/src/pe/communication/Instantiate.h index 6096c6552..f10a27ddd 100644 --- a/src/pe/communication/Instantiate.h +++ b/src/pe/communication/Instantiate.h @@ -55,7 +55,7 @@ void correctBodyPosition(const math::AABB& domain, const Vec3& center, Vec3& pos } template < class BodyT > -BodyT* instantiate( mpi::RecvBuffer& /*buffer*/, const math::AABB& /*domain*/, const math::AABB& /*block*/, BodyT*& /*newBody*/ ) +std::unique_ptr<BodyT> instantiate( mpi::RecvBuffer& /*buffer*/, const math::AABB& /*domain*/, const math::AABB& /*block*/, BodyT*& /*newBody*/ ) { WALBERLA_ABORT( "Body instantiation not implemented! (" << demangle(typeid(BodyT).name()) << ")" ); } diff --git a/src/pe/communication/ParseMessage.h b/src/pe/communication/ParseMessage.h index 96b422cd7..deea7db96 100644 --- a/src/pe/communication/ParseMessage.h +++ b/src/pe/communication/ParseMessage.h @@ -64,17 +64,18 @@ void parseMessage(Owner sender, mpi::RecvBuffer& rb, const domain_decomposition: WALBERLA_LOG_DETAIL( "Received " << objparam.geomType_ << " copy notification from neighboring process with rank " << sender ); - BodyID obj = UnmarshalDynamically<BodyTypeTuple>::execute(rb, objparam.geomType_, blockStorage.getDomain(), block.getAABB()); + BodyPtr obj = UnmarshalDynamically<BodyTypeTuple>::execute(rb, objparam.geomType_, blockStorage.getDomain(), block.getAABB()); obj->setRemote( true ); + obj->MPITrait.setBlockState( sender.blockID_ ); + if (shadowStorage.find( obj->getSystemID() ) == shadowStorage.end()) { - WALBERLA_LOG_DETAIL( "Adding new shadow copy with id " << obj->getSystemID() << " to domain " << block.getAABB() << ".\n" << obj); - shadowStorage.add( obj ); + WALBERLA_LOG_DETAIL( "Adding new shadow copy with id " << obj->getSystemID() << " to domain " << block.getAABB() << ".\n" << *obj); + shadowStorage.add( std::move(obj) ); } else { WALBERLA_LOG_DETAIL( "Shadow copy with id " << obj->getSystemID() << " already existend."); } - obj->MPITrait.setBlockState( sender.blockID_ ); WALBERLA_LOG_DETAIL( "Processed " << objparam.geomType_ << " copy notification." ); diff --git a/src/pe/communication/rigidbody/Box.h b/src/pe/communication/rigidbody/Box.h index 3a4c35a20..77f543cbc 100644 --- a/src/pe/communication/rigidbody/Box.h +++ b/src/pe/communication/rigidbody/Box.h @@ -59,16 +59,17 @@ void marshal( mpi::SendBuffer& buffer, const Box& obj ); */ void unmarshal( mpi::RecvBuffer& buffer, BoxParameters& objparam ); //************************************************************************************************* -inline BoxID instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, BoxID& newBody ) +inline BoxPtr instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, BoxID& newBody ) { BoxParameters subobjparam; unmarshal( buffer, subobjparam ); correctBodyPosition(domain, block.center(), subobjparam.gpos_); - newBody = new Box( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.lengths_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); - newBody->setLinearVel( subobjparam.v_ ); - newBody->setAngularVel( subobjparam.w_ ); - newBody->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); - return newBody; + auto bx = std::make_unique<Box>( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.lengths_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); + bx->setLinearVel( subobjparam.v_ ); + bx->setAngularVel( subobjparam.w_ ); + bx->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); + newBody = bx.get(); + return bx; } } // namespace communication diff --git a/src/pe/communication/rigidbody/Capsule.h b/src/pe/communication/rigidbody/Capsule.h index 7f99b1d1d..3ef02f07d 100644 --- a/src/pe/communication/rigidbody/Capsule.h +++ b/src/pe/communication/rigidbody/Capsule.h @@ -61,16 +61,17 @@ void unmarshal( mpi::RecvBuffer& buffer, CapsuleParameters& objparam ); //************************************************************************************************* -inline CapsuleID instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, CapsuleID& newBody ) +inline CapsulePtr instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, CapsuleID& newBody ) { CapsuleParameters subobjparam; unmarshal( buffer, subobjparam ); correctBodyPosition(domain, block.center(), subobjparam.gpos_); - newBody = new Capsule( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.radius_, subobjparam.length_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); - newBody->setLinearVel( subobjparam.v_ ); - newBody->setAngularVel( subobjparam.w_ ); - newBody->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); - return newBody; + auto cp = std::make_unique<Capsule>( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.radius_, subobjparam.length_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); + cp->setLinearVel( subobjparam.v_ ); + cp->setAngularVel( subobjparam.w_ ); + cp->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); + newBody = cp.get(); + return cp; } } // namespace communication diff --git a/src/pe/communication/rigidbody/Ellipsoid.h b/src/pe/communication/rigidbody/Ellipsoid.h index 57916d268..a292c5b1d 100644 --- a/src/pe/communication/rigidbody/Ellipsoid.h +++ b/src/pe/communication/rigidbody/Ellipsoid.h @@ -61,16 +61,17 @@ void unmarshal( mpi::RecvBuffer& buffer, EllipsoidParameters& objparam ); //************************************************************************************************* -inline EllipsoidID instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, EllipsoidID& newBody ) +inline EllipsoidPtr instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, EllipsoidID& newBody ) { EllipsoidParameters subobjparam; unmarshal( buffer, subobjparam ); correctBodyPosition(domain, block.center(), subobjparam.gpos_); - newBody = new Ellipsoid( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.semiAxes_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); - newBody->setLinearVel( subobjparam.v_ ); - newBody->setAngularVel( subobjparam.w_ ); - newBody->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); - return newBody; + auto el = std::make_unique<Ellipsoid>( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.semiAxes_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); + el->setLinearVel( subobjparam.v_ ); + el->setAngularVel( subobjparam.w_ ); + el->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); + newBody = el.get(); + return el; } } // namespace communication diff --git a/src/pe/communication/rigidbody/Sphere.h b/src/pe/communication/rigidbody/Sphere.h index 55881f65a..3e7da5601 100644 --- a/src/pe/communication/rigidbody/Sphere.h +++ b/src/pe/communication/rigidbody/Sphere.h @@ -61,16 +61,17 @@ void unmarshal( mpi::RecvBuffer& buffer, SphereParameters& objparam ); //************************************************************************************************* -inline SphereID instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, SphereID& newBody ) +inline SpherePtr instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, SphereID& newBody ) { SphereParameters subobjparam; unmarshal( buffer, subobjparam ); correctBodyPosition(domain, block.center(), subobjparam.gpos_); - newBody = new Sphere( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.radius_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); - newBody->setLinearVel( subobjparam.v_ ); - newBody->setAngularVel( subobjparam.w_ ); - newBody->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); - return newBody; + auto sp = std::make_unique<Sphere>( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.radius_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); + sp->setLinearVel( subobjparam.v_ ); + sp->setAngularVel( subobjparam.w_ ); + sp->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); + newBody = sp.get(); + return sp; } } // namespace communication diff --git a/src/pe/communication/rigidbody/Squirmer.h b/src/pe/communication/rigidbody/Squirmer.h index 9029c9358..3c9b3fbe0 100644 --- a/src/pe/communication/rigidbody/Squirmer.h +++ b/src/pe/communication/rigidbody/Squirmer.h @@ -61,16 +61,17 @@ void unmarshal( mpi::RecvBuffer& buffer, SquirmerParameters& objparam ); //************************************************************************************************* -inline SquirmerID instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, SquirmerID& newBody ) +inline SquirmerPtr instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, SquirmerID& newBody ) { SquirmerParameters subobjparam; unmarshal( buffer, subobjparam ); correctBodyPosition(domain, block.center(), subobjparam.gpos_); - newBody = new Squirmer( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.radius_, subobjparam.squirmerVelocity_, subobjparam.squirmerBeta_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); - newBody->setLinearVel( subobjparam.v_ ); - newBody->setAngularVel( subobjparam.w_ ); - newBody->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); - return newBody; + auto sq = std::make_unique<Squirmer>( subobjparam.sid_, subobjparam.uid_, subobjparam.gpos_, subobjparam.rpos_, subobjparam.q_, subobjparam.radius_, subobjparam.squirmerVelocity_, subobjparam.squirmerBeta_, subobjparam.material_, false, subobjparam.communicating_, subobjparam.infiniteMass_ ); + sq->setLinearVel( subobjparam.v_ ); + sq->setAngularVel( subobjparam.w_ ); + sq->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); + newBody = sq.get(); + return sq; } } // namespace communication diff --git a/src/pe/communication/rigidbody/Union.h b/src/pe/communication/rigidbody/Union.h index 596f54118..544c131ef 100644 --- a/src/pe/communication/rigidbody/Union.h +++ b/src/pe/communication/rigidbody/Union.h @@ -110,22 +110,22 @@ void unmarshal( mpi::RecvBuffer& buffer, UnionParameters& objparam ) template <typename BodyTypeTuple> -inline Union<BodyTypeTuple>* instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, Union<BodyTypeTuple>*& newBody ) +inline std::unique_ptr<Union<BodyTypeTuple>> instantiate( mpi::RecvBuffer& buffer, const math::AABB& domain, const math::AABB& block, Union<BodyTypeTuple>*& newBody ) { UnionParameters subobjparam; unmarshal( buffer, subobjparam ); correctBodyPosition(domain, block.center(), subobjparam.gpos_); - newBody = new Union<BodyTypeTuple>( subobjparam.sid_, - subobjparam.uid_, - subobjparam.gpos_, - subobjparam.rpos_, - subobjparam.q_, - false, - subobjparam.communicating_, - subobjparam.infiniteMass_ ); - newBody->setLinearVel( subobjparam.v_ ); - newBody->setAngularVel( subobjparam.w_ ); - newBody->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); + auto un = std::make_unique<Union<BodyTypeTuple>>( subobjparam.sid_, + subobjparam.uid_, + subobjparam.gpos_, + subobjparam.rpos_, + subobjparam.q_, + false, + subobjparam.communicating_, + subobjparam.infiniteMass_ ); + un->setLinearVel( subobjparam.v_ ); + un->setAngularVel( subobjparam.w_ ); + un->MPITrait.setOwner( subobjparam.mpiTrait_.owner_ ); // Decoding the contained primitives for( size_t i = 0; i < subobjparam.size_; ++i ) @@ -134,10 +134,11 @@ inline Union<BodyTypeTuple>* instantiate( mpi::RecvBuffer& buffer, const math::A buffer >> type; BodyPtr obj( UnmarshalDynamically<BodyTypeTuple>::execute(buffer, type, domain, block) ); obj->setRemote( true ); - newBody->add(std::move(obj)); + un->add(std::move(obj)); } - return newBody; + newBody = un.get(); + return un; } } // namespace communication diff --git a/src/pe/rigidbody/StorageDataHandling.h b/src/pe/rigidbody/StorageDataHandling.h index ae407845d..2d49a98c8 100644 --- a/src/pe/rigidbody/StorageDataHandling.h +++ b/src/pe/rigidbody/StorageDataHandling.h @@ -214,16 +214,16 @@ void StorageDataHandling<BodyTuple>::deserializeImpl( IBlock * const block, cons typename RigidBodyCopyNotification::Parameters objparam; unmarshal( buffer, objparam ); - BodyID bd = UnmarshalDynamically<BodyTuple>::execute(buffer, objparam.geomType_, block->getBlockStorage().getDomain(), block->getAABB()); + auto bd = UnmarshalDynamically<BodyTuple>::execute(buffer, objparam.geomType_, block->getBlockStorage().getDomain(), block->getAABB()); bd->setRemote( false ); bd->MPITrait.setOwner(Owner(MPIManager::instance()->rank(), block->getId().getID())); if ( !block->getAABB().contains( bd->getPosition()) ) { - WALBERLA_ABORT("Loaded body not contained within block!\n" << "aabb: " << block->getAABB() << "\nparticle:" << bd ); + WALBERLA_ABORT("Loaded body not contained within block!\n" << "aabb: " << block->getAABB() << "\nparticle:" << *bd ); } WALBERLA_ASSERT_EQUAL(localBodyStorage.find( bd->getSystemID() ), localBodyStorage.end()); - localBodyStorage.add(bd); + localBodyStorage.add(std::move(bd)); --numBodies; } diff --git a/tests/pe/Marshalling.cpp b/tests/pe/Marshalling.cpp index b815a5d96..985269bd6 100644 --- a/tests/pe/Marshalling.cpp +++ b/tests/pe/Marshalling.cpp @@ -41,11 +41,14 @@ using namespace walberla::pe::communication; typedef boost::tuple<Sphere> UnionTypeTuple; typedef Union< UnionTypeTuple > UnionT; typedef UnionT* UnionID; +typedef std::unique_ptr<UnionT> UnionPtr; typedef boost::tuple<Box, Capsule, Sphere, Squirmer, UnionT, Ellipsoid> BodyTuple ; void testBox() { + WALBERLA_LOG_INFO_ON_ROOT("*** testBox ***"); + MaterialID iron = Material::find("iron"); Box b1(759846, 1234794, Vec3(real_c(1), real_c(2), real_c(3)), Vec3(0,0,0), Quat(), Vec3(1,2,3), iron, false, true, false); @@ -56,7 +59,8 @@ void testBox() MarshalDynamically<BodyTuple>::execute(sb, b1); mpi::RecvBuffer rb(sb); - BoxID b2 = static_cast<BoxID> (UnmarshalDynamically<BodyTuple>::execute(rb, Box::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)))); + auto bPtr = UnmarshalDynamically<BodyTuple>::execute(rb, Box::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100))); + BoxID b2 = static_cast<BoxID>(bPtr.get()); WALBERLA_CHECK_FLOAT_EQUAL(b1.getPosition(), b2->getPosition()); WALBERLA_CHECK_FLOAT_EQUAL(b1.getLinearVel(), b2->getLinearVel()); @@ -68,6 +72,8 @@ void testBox() void testCapsule() { + WALBERLA_LOG_INFO_ON_ROOT("*** testCapsule ***"); + MaterialID iron = Material::find("iron"); Capsule c1(759846, 1234794, Vec3(real_c(1), real_c(2), real_c(3)), Vec3(0,0,0), Quat(), 5, 7, iron, false, false, false); @@ -78,7 +84,8 @@ void testCapsule() MarshalDynamically<BodyTuple>::execute(sb, c1); mpi::RecvBuffer rb(sb); - CapsuleID c2 = static_cast<CapsuleID> (UnmarshalDynamically<BodyTuple>::execute(rb, Capsule::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)))); + auto cPtr = UnmarshalDynamically<BodyTuple>::execute(rb, Capsule::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100))); + CapsuleID c2 = static_cast<CapsuleID> (cPtr.get()); WALBERLA_CHECK_FLOAT_EQUAL(c1.getPosition(), c2->getPosition()); WALBERLA_CHECK_FLOAT_EQUAL(c1.getLinearVel(), c2->getLinearVel()); @@ -91,6 +98,8 @@ void testCapsule() void testSphere() { + WALBERLA_LOG_INFO_ON_ROOT("*** testSphere ***"); + MaterialID iron = Material::find("iron"); Sphere s1(759846, 1234794, Vec3(real_c(1), real_c(2), real_c(3)), Vec3(0,0,0), Quat(), 5, iron, false, false, false); @@ -101,7 +110,8 @@ void testSphere() MarshalDynamically<BodyTuple>::execute(sb, s1); mpi::RecvBuffer rb(sb); - SphereID s2 = static_cast<SphereID> (UnmarshalDynamically<BodyTuple>::execute(rb, Sphere::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)))); + auto sPtr = UnmarshalDynamically<BodyTuple>::execute(rb, Sphere::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100))); + SphereID s2 = static_cast<SphereID> (sPtr.get()); WALBERLA_CHECK_FLOAT_EQUAL(s1.getPosition(), s2->getPosition()); WALBERLA_CHECK_FLOAT_EQUAL(s1.getLinearVel(), s2->getLinearVel()); @@ -113,6 +123,8 @@ void testSphere() void testSquirmer() { + WALBERLA_LOG_INFO_ON_ROOT("*** testSquirmer ***"); + MaterialID iron = Material::find("iron"); Squirmer s1(759846, 1234794, Vec3(real_c(1), real_c(2), real_c(3)), Vec3(0,0,0), Quat(), real_c(5), real_c(0.1), real_c(4.93), iron, false, false, false); @@ -123,7 +135,8 @@ void testSquirmer() MarshalDynamically<BodyTuple>::execute(sb, s1); mpi::RecvBuffer rb(sb); - SquirmerID s2 = static_cast<SquirmerID> (UnmarshalDynamically<BodyTuple>::execute(rb, Squirmer::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)))); + auto sPtr = UnmarshalDynamically<BodyTuple>::execute(rb, Squirmer::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100))); + SquirmerID s2 = static_cast<SquirmerID> (sPtr.get()); WALBERLA_CHECK_FLOAT_EQUAL(s1.getSquirmerVelocity(), s2->getSquirmerVelocity()); WALBERLA_CHECK_FLOAT_EQUAL(s1.getSquirmerBeta(), s2->getSquirmerBeta()); @@ -131,6 +144,8 @@ void testSquirmer() void testEllipsoid() { + WALBERLA_LOG_INFO_ON_ROOT("*** testEllipsoid ***"); + MaterialID iron = Material::find("iron"); Ellipsoid e1(759847, 1234795, Vec3(real_c(1), real_c(2), real_c(3)), Vec3(0,0,0), Quat(), Vec3(3,1,5), iron, false, false, false); @@ -141,7 +156,8 @@ void testEllipsoid() MarshalDynamically<BodyTuple>::execute(sb, e1); mpi::RecvBuffer rb(sb); - EllipsoidID e2 = static_cast<EllipsoidID> (UnmarshalDynamically<BodyTuple>::execute(rb, Ellipsoid::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)))); + auto ePtr = UnmarshalDynamically<BodyTuple>::execute(rb, Ellipsoid::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100))); + EllipsoidID e2 = static_cast<EllipsoidID>(ePtr.get()); WALBERLA_CHECK_FLOAT_EQUAL(e1.getPosition(), e2->getPosition()); WALBERLA_CHECK_FLOAT_EQUAL(e1.getLinearVel(), e2->getLinearVel()); @@ -153,6 +169,7 @@ void testEllipsoid() void testUnion() { + WALBERLA_LOG_INFO_ON_ROOT("*** testUnion ***"); UnionT u1(159, 423, Vec3(real_c(1), real_c(2), real_c(3)), Vec3(0,0,0), Quat(), false, false, false); SphereID s11 = createSphere< UnionTypeTuple >(&u1, 1234794, Vec3(real_c(1), real_c(2), real_c(3)), 2); SphereID s21 = createSphere< UnionTypeTuple >(&u1, 4567789, Vec3(real_c(3), real_c(2), real_c(3)), real_c(1.5)); @@ -163,7 +180,8 @@ void testUnion() MarshalDynamically<BodyTuple>::execute(sb, u1); mpi::RecvBuffer rb(sb); - UnionID u2 = static_cast<UnionID> (UnmarshalDynamically<BodyTuple>::execute(rb, UnionT::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)))); + auto uPtr = UnmarshalDynamically<BodyTuple>::execute(rb, UnionT::getStaticTypeID(), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100)), math::AABB(Vec3(-100,-100,-100), Vec3(100,100,100))); + UnionID u2 = static_cast<UnionID>(uPtr.get()); WALBERLA_CHECK_NOT_NULLPTR( u2 ); WALBERLA_CHECK_EQUAL(u1.size(), 2); -- GitLab