Skip to content
Snippets Groups Projects
Commit 73e67f6e authored by Christoph Rettinger's avatar Christoph Rettinger
Browse files

Merge branch 'code_coverage' into 'master'

Increase Code Coverage of MESA-PD

See merge request walberla/walberla!360
parents db9b8919 43b8bacf
No related merge requests found
......@@ -20,22 +20,34 @@
#pragma once
#include <core/mpi/MPIManager.h>
#include <mesa_pd/domain/IDomain.h>
namespace walberla {
namespace mesa_pd {
namespace domain {
/**
* Every process assumes the whole simulation space belongs to its subdomain.
*/
class InfiniteDomain : public IDomain
{
public:
bool isContainedInProcessSubdomain(const uint_t /*rank*/, const Vec3& /*pt*/) const override {return true;}
int findContainingProcessRank(const Vec3& /*pt*/) const override {return walberla::mpi::MPIManager::instance()->rank();}
///Everything belongs to the calling process.
bool isContainedInProcessSubdomain(const uint_t rank, const Vec3& /*pt*/) const override {return rank==rank_;}
///Everything belongs to the calling process.
int findContainingProcessRank(const Vec3& /*pt*/) const override {return static_cast<int>(rank_);}
///Nothing to do here since domain is infinite.
void periodicallyMapToDomain(Vec3& /*pt*/) const override {}
///If I own everything I do not have neighbors.
std::vector<uint_t> getNeighborProcesses() const override {return {};}
///Everything belongs to my subdomain.
bool intersectsWithProcessSubdomain(const uint_t rank, const Vec3& /*pt*/, const real_t& /*radius*/) const override
{ return int_c(rank)==walberla::mpi::MPIManager::instance()->rank() ? true : false;}
{ return rank==rank_;}
///Nothing to do here.
void correctParticlePosition(Vec3& /*pt*/) const override {}
private:
const uint_t rank_ = static_cast<uint_t>(MPIManager::instance()->rank());
};
} //namespace domain
......
......@@ -37,6 +37,9 @@ waLBerla_execute_test( NAME MESA_PD_COMMON_IntersectionRatio )
waLBerla_compile_test( NAME MESA_PD_ContactDetection FILES ContactDetection.cpp DEPENDS blockforest core pe)
waLBerla_execute_test( NAME MESA_PD_ContactDetection PROCESSES 8 )
waLBerla_compile_test( NAME MESA_PD_Data_ContactHistory FILES data/ContactHistory.cpp DEPENDS core )
waLBerla_execute_test( NAME MESA_PD_Data_ContactHistory )
waLBerla_compile_test( NAME MESA_PD_Data_Flags FILES data/Flags.cpp DEPENDS core )
waLBerla_execute_test( NAME MESA_PD_Data_Flags )
......@@ -66,6 +69,9 @@ waLBerla_execute_test( NAME MESA_PD_Domain_DistanceCalculation )
waLBerla_compile_test( NAME MESA_PD_Domain_DynamicRefinement FILES domain/DynamicRefinement.cpp DEPENDS blockforest core pe )
waLBerla_execute_test( NAME MESA_PD_Domain_DynamicRefinement PROCESSES 8)
waLBerla_compile_test( NAME MESA_PD_Domain_InfiniteDomain FILES domain/InfiniteDomain.cpp DEPENDS core )
waLBerla_execute_test( NAME MESA_PD_Domain_InfiniteDomain )
waLBerla_compile_test( NAME MESA_PD_Domain_SerializeDeserialize FILES domain/SerializeDeserialize.cpp DEPENDS blockforest core pe)
waLBerla_execute_test( NAME MESA_PD_Domain_SerializeDeserialize PROCESSES 8 )
......
//======================================================================================================================
//
// This file is part of waLBerla. waLBerla is free software: you can
// redistribute it and/or modify it under the terms of the GNU General Public
// License as published by the Free Software Foundation, either version 3 of
// the License, or (at your option) any later version.
//
// waLBerla is distributed in the hope that it will be useful, but WITHOUT
// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
// FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
// for more details.
//
// You should have received a copy of the GNU General Public License along
// with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>.
//
//! \file ContactHistory.cpp
//! \author Sebastian Eibl <sebastian.eibl@fau.de>
//
//======================================================================================================================
#include <mesa_pd/data/ContactHistory.h>
#include <core/Environment.h>
#include <core/logging/Logging.h>
#include <algorithm>
#include <iostream>
namespace walberla {
using namespace walberla::mesa_pd;
void basic_test()
{
//init data structures
data::ContactHistory cs;
cs.setImpactVelocityMagnitude(1.23456_r);
cs.setIsSticking(true);
cs.setTangentialSpringDisplacement(Vec3(1.23_r,2.345_r,3.56_r));
mpi::SendBuffer sb;
sb << cs;
mpi::RecvBuffer rb(sb);
data::ContactHistory cs_recv;
rb >> cs_recv;
WALBERLA_CHECK_IDENTICAL(cs.getImpactVelocityMagnitude(), cs_recv.getImpactVelocityMagnitude());
WALBERLA_CHECK_IDENTICAL(cs.getIsSticking(), cs_recv.getIsSticking());
WALBERLA_CHECK_IDENTICAL(cs.getTangentialSpringDisplacement(), cs_recv.getTangentialSpringDisplacement());
WALBERLA_LOG_DEVEL( cs );
}
int main( int argc, char ** argv )
{
Environment env(argc, argv);
WALBERLA_UNUSED(env);
mpi::MPIManager::instance()->useWorldComm();
basic_test();
return EXIT_SUCCESS;
}
} //namespace walberla
int main( int argc, char ** argv )
{
return walberla::main(argc, argv);
}
//======================================================================================================================
//
// This file is part of waLBerla. waLBerla is free software: you can
// redistribute it and/or modify it under the terms of the GNU General Public
// License as published by the Free Software Foundation, either version 3 of
// the License, or (at your option) any later version.
//
// waLBerla is distributed in the hope that it will be useful, but WITHOUT
// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
// FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
// for more details.
//
// You should have received a copy of the GNU General Public License along
// with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>.
//
//! \file InfiniteDomain.cpp
//! \author Sebastian Eibl <sebastian.eibl@fau.de>
//
//======================================================================================================================
#include <mesa_pd/domain/InfiniteDomain.h>
#include <core/Environment.h>
#include <core/logging/Logging.h>
#include <iostream>
namespace walberla {
using namespace walberla::mesa_pd;
void main( int argc, char ** argv )
{
Environment env(argc, argv);
WALBERLA_UNUSED(env);
mpi::MPIManager::instance()->useWorldComm();
auto rank = mpi::MPIManager::instance()->rank();
auto randomPoint = Vec3(1.23_r,2.34_r,3.45_r);
domain::InfiniteDomain domain;
WALBERLA_CHECK(domain.isContainedInProcessSubdomain(static_cast<uint_t>(rank), randomPoint));
WALBERLA_CHECK(!domain.isContainedInProcessSubdomain(static_cast<uint_t>(rank + 1), randomPoint));
WALBERLA_CHECK_EQUAL(domain.findContainingProcessRank(randomPoint), rank);
auto pt = randomPoint;
domain.periodicallyMapToDomain(pt);
WALBERLA_CHECK_IDENTICAL(pt, randomPoint);
WALBERLA_CHECK_EQUAL(domain.getNeighborProcesses().size(), 0);
WALBERLA_CHECK(domain.intersectsWithProcessSubdomain(static_cast<uint_t>(rank), randomPoint, 1_r));
WALBERLA_CHECK(!domain.intersectsWithProcessSubdomain(static_cast<uint_t>(rank + 1), randomPoint, 1_r));
domain.correctParticlePosition(pt);
WALBERLA_CHECK_IDENTICAL(pt, randomPoint);
}
}
int main( int argc, char ** argv )
{
walberla::main(argc, argv);
return EXIT_SUCCESS;
}
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment