Skip to content
Snippets Groups Projects
Commit 56276b30 authored by Sebastian Eibl's avatar Sebastian Eibl
Browse files

[API] added reset function to ReduceProperty

parent 9bd7bff3
Branches
Tags
No related merge requests found
...@@ -89,15 +89,15 @@ int main( int argc, char ** argv ) ...@@ -89,15 +89,15 @@ int main( int argc, char ** argv )
WALBERLA_LOG_DEVEL(timestep); WALBERLA_LOG_DEVEL(timestep);
linkedCells.clear(); linkedCells.clear();
storage->forEachParticle(true, kernel::SelectAll(), ac, ipilc, ac, linkedCells); storage->forEachParticle(true, kernel::SelectAll(), ac, ipilc, ac, linkedCells);
storage->forEachParticle(true, kernel::SelectAll(), ac, vvPreForce, ac); storage->forEachParticle(true, kernel::SelectLocal(), ac, vvPreForce, ac);
linkedCells.forEachParticlePairHalf(true, kernel::SelectAll(), ac, lj, ac); linkedCells.forEachParticlePairHalf(true, kernel::SelectAll(), ac, lj, ac);
const real_t coeff = real_t(0.2); const real_t coeff = real_t(0.2);
storage->forEachParticle(true, storage->forEachParticle(true,
kernel::SelectAll(), kernel::SelectLocal(),
ac, ac,
[coeff](const size_t idx, auto& access){ access.setForce(idx, -coeff*access.getPosition(idx) + access.getForce(idx)); }, [coeff](const size_t idx, auto& access){ access.setForce(idx, -coeff*access.getPosition(idx) + access.getForce(idx)); },
ac); ac);
storage->forEachParticle(true, kernel::SelectAll(), ac, vvPostForce, ac); storage->forEachParticle(true, kernel::SelectLocal(), ac, vvPostForce, ac);
vtkWriter->write(); vtkWriter->write();
} }
......
...@@ -29,6 +29,7 @@ ...@@ -29,6 +29,7 @@
#include <mesa_pd/data/DataTypes.h> #include <mesa_pd/data/DataTypes.h>
#include <mesa_pd/data/Flags.h> #include <mesa_pd/data/Flags.h>
#include <mesa_pd/data/ParticleStorage.h> #include <mesa_pd/data/ParticleStorage.h>
#include <mesa_pd/mpi/notifications/reset.h>
#include <core/mpi/BufferSystem.h> #include <core/mpi/BufferSystem.h>
#include <core/logging/Logging.h> #include <core/logging/Logging.h>
...@@ -100,6 +101,7 @@ void ReduceProperty::operator()(data::ParticleStorage& ps) const ...@@ -100,6 +101,7 @@ void ReduceProperty::operator()(data::ParticleStorage& ps) const
} }
sb << Notification( p ); sb << Notification( p );
reset<Notification>( p );
} else } else
{ {
//local particles should receive the property and sum it up //local particles should receive the property and sum it up
......
...@@ -28,6 +28,7 @@ ...@@ -28,6 +28,7 @@
#include <mesa_pd/data/ContactHistory.h> #include <mesa_pd/data/ContactHistory.h>
#include <mesa_pd/data/DataTypes.h> #include <mesa_pd/data/DataTypes.h>
#include <mesa_pd/mpi/notifications/reset.h>
#include <core/mpi/BufferDataTypeExtensions.h> #include <core/mpi/BufferDataTypeExtensions.h>
#include <core/mpi/Datatype.h> #include <core/mpi/Datatype.h>
...@@ -56,6 +57,12 @@ public: ...@@ -56,6 +57,12 @@ public:
const data::Particle& p_; const data::Particle& p_;
}; };
template <>
void reset<ContactHistoryNotification>(data::Particle& p)
{
p.setNewContactHistory(std::map<walberla::id_t, walberla::mesa_pd::data::ContactHistory>());
}
void reduce(data::Particle&& p, const ContactHistoryNotification::Parameters& objparam) void reduce(data::Particle&& p, const ContactHistoryNotification::Parameters& objparam)
{ {
auto& ch = p.getNewContactHistoryRef(); auto& ch = p.getNewContactHistoryRef();
......
...@@ -29,6 +29,7 @@ ...@@ -29,6 +29,7 @@
#include <mesa_pd/data/DataTypes.h> #include <mesa_pd/data/DataTypes.h>
#include <mesa_pd/data/ParticleStorage.h> #include <mesa_pd/data/ParticleStorage.h>
#include <mesa_pd/mpi/notifications/NotificationType.h> #include <mesa_pd/mpi/notifications/NotificationType.h>
#include <mesa_pd/mpi/notifications/reset.h>
#include <core/mpi/Datatype.h> #include <core/mpi/Datatype.h>
#include <core/mpi/RecvBuffer.h> #include <core/mpi/RecvBuffer.h>
...@@ -55,6 +56,13 @@ public: ...@@ -55,6 +56,13 @@ public:
const data::Particle& p_; const data::Particle& p_;
}; };
template <>
void reset<ForceTorqueNotification>(data::Particle& p)
{
p.setForce( Vec3(real_t(0)) );
p.setTorque( Vec3(real_t(0)) );
}
void reduce(data::Particle&& p, const ForceTorqueNotification::Parameters& objparam) void reduce(data::Particle&& p, const ForceTorqueNotification::Parameters& objparam)
{ {
p.getForceRef() += objparam.force_; p.getForceRef() += objparam.force_;
......
...@@ -29,6 +29,7 @@ ...@@ -29,6 +29,7 @@
#include <mesa_pd/data/DataTypes.h> #include <mesa_pd/data/DataTypes.h>
#include <mesa_pd/data/Flags.h> #include <mesa_pd/data/Flags.h>
#include <mesa_pd/data/ParticleStorage.h> #include <mesa_pd/data/ParticleStorage.h>
#include <mesa_pd/mpi/notifications/reset.h>
#include <core/mpi/BufferSystem.h> #include <core/mpi/BufferSystem.h>
#include <core/logging/Logging.h> #include <core/logging/Logging.h>
...@@ -100,6 +101,7 @@ void ReduceProperty::operator()(data::ParticleStorage& ps) const ...@@ -100,6 +101,7 @@ void ReduceProperty::operator()(data::ParticleStorage& ps) const
} }
sb << Notification( p ); sb << Notification( p );
reset<Notification>( p );
} else } else
{ {
//local particles should receive the property and sum it up //local particles should receive the property and sum it up
......
...@@ -28,6 +28,7 @@ ...@@ -28,6 +28,7 @@
#include <mesa_pd/data/ContactHistory.h> #include <mesa_pd/data/ContactHistory.h>
#include <mesa_pd/data/DataTypes.h> #include <mesa_pd/data/DataTypes.h>
#include <mesa_pd/mpi/notifications/reset.h>
#include <core/mpi/BufferDataTypeExtensions.h> #include <core/mpi/BufferDataTypeExtensions.h>
#include <core/mpi/Datatype.h> #include <core/mpi/Datatype.h>
...@@ -56,6 +57,12 @@ public: ...@@ -56,6 +57,12 @@ public:
const data::Particle& p_; const data::Particle& p_;
}; };
template <>
void reset<ContactHistoryNotification>(data::Particle& p)
{
p.setNewContactHistory(std::map<walberla::id_t, walberla::mesa_pd::data::ContactHistory>());
}
void reduce(data::Particle&& p, const ContactHistoryNotification::Parameters& objparam) void reduce(data::Particle&& p, const ContactHistoryNotification::Parameters& objparam)
{ {
auto& ch = p.getNewContactHistoryRef(); auto& ch = p.getNewContactHistoryRef();
......
...@@ -29,6 +29,7 @@ ...@@ -29,6 +29,7 @@
#include <mesa_pd/data/DataTypes.h> #include <mesa_pd/data/DataTypes.h>
#include <mesa_pd/data/ParticleStorage.h> #include <mesa_pd/data/ParticleStorage.h>
#include <mesa_pd/mpi/notifications/NotificationType.h> #include <mesa_pd/mpi/notifications/NotificationType.h>
#include <mesa_pd/mpi/notifications/reset.h>
#include <core/mpi/Datatype.h> #include <core/mpi/Datatype.h>
#include <core/mpi/RecvBuffer.h> #include <core/mpi/RecvBuffer.h>
...@@ -55,6 +56,13 @@ public: ...@@ -55,6 +56,13 @@ public:
const data::Particle& p_; const data::Particle& p_;
}; };
template <>
void reset<ForceTorqueNotification>(data::Particle& p)
{
p.setForce( Vec3(real_t(0)) );
p.setTorque( Vec3(real_t(0)) );
}
void reduce(data::Particle&& p, const ForceTorqueNotification::Parameters& objparam) void reduce(data::Particle&& p, const ForceTorqueNotification::Parameters& objparam)
{ {
p.getForceRef() += objparam.force_; p.getForceRef() += objparam.force_;
......
//======================================================================================================================
//
// 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 reset.h
//! \author Sebastian Eibl <sebastian.eibl@fau.de>
//
//======================================================================================================================
#pragma once
#include <mesa_pd/data/ParticleStorage.h>
namespace walberla {
namespace mesa_pd {
template <class Notification>
void reset(data::Particle& /*p*/)
{
WALBERLA_ABORT("not implemented!");
}
} // mpi
} // walberla
...@@ -95,7 +95,7 @@ void main( int argc, char ** argv ) ...@@ -95,7 +95,7 @@ void main( int argc, char ** argv )
WALBERLA_CHECK_FLOAT_EQUAL( pIt->getForce(), Vec3(real_t(28)) ); WALBERLA_CHECK_FLOAT_EQUAL( pIt->getForce(), Vec3(real_t(28)) );
} else } else
{ {
WALBERLA_CHECK_FLOAT_EQUAL( pIt->getForce(), Vec3(real_t(walberla::mpi::MPIManager::instance()->rank())) ); WALBERLA_CHECK_FLOAT_EQUAL( pIt->getForce(), Vec3(0) );
} }
} }
......
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