Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Testsuite: coldet feature guard and test name #2344

Merged
merged 4 commits into from
Nov 7, 2018
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
64 changes: 18 additions & 46 deletions src/core/collision.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -19,6 +19,7 @@

#include <vector>

#include "boost/mpi/collectives.hpp"
#include "cells.hpp"
#include "collision.hpp"
#include "communication.hpp"
Expand All @@ -29,7 +30,9 @@
#include "particle_data.hpp"
#include "rotation.hpp"
#include "utils/mpi/all_compare.hpp"
#include "utils/mpi/gather_buffer.hpp"
#include "virtual_sites/VirtualSitesRelative.hpp"
#include <boost/serialization/serialization.hpp>

#ifdef COLLISION_DETECTION_DEBUG
#define TRACE(a) a
Expand All @@ -45,6 +48,16 @@ typedef struct {
int pp2; // 2nd particle id
} collision_struct;

namespace boost {
namespace serialization {
template <typename Archive>
void serialize(Archive &ar, collision_struct &c, const unsigned int) {
ar &c.pp1;
ar &c.pp2;
}
} // namespace serialization
} // namespace boost

// During force calculation, colliding particles are recorded in the queue
// The queue is processed after force calculation, when it is save to add
// particles
Expand Down Expand Up @@ -431,45 +444,9 @@ void glue_to_surface_bind_part_to_vs(const Particle *const p1,
#endif

std::vector<collision_struct> gather_global_collision_queue() {
std::vector<collision_struct> res;

std::vector<int> displacements(n_nodes); // offsets into collisions

// Initialize number of collisions gathered from all processors
std::vector<int> counts(n_nodes);
// Total number of collisions
int total_collisions;
int local_queue_size = local_collision_queue.size();
MPI_Allreduce(&local_queue_size, &total_collisions, 1, MPI_INT, MPI_SUM,
comm_cart);

if (total_collisions == 0)
return res;

// Gather number of collisions
MPI_Allgather(&local_queue_size, 1, MPI_INT, &(counts[0]), 1, MPI_INT,
comm_cart);

// initialize displacement information for all nodes
displacements[0] = 0;

// Find where to place collision information for each processor
std::vector<int> byte_counts(n_nodes);
for (int k = 1; k < n_nodes; k++)
displacements[k] =
displacements[k - 1] + (counts[k - 1]) * sizeof(collision_struct);

for (int k = 0; k < n_nodes; k++)
byte_counts[k] = counts[k] * sizeof(collision_struct);

// Allocate mem for the new collision info

res.resize(total_collisions);

// Gather collision information from all nodes and send it to all nodes
MPI_Allgatherv(&(local_collision_queue[0]), byte_counts[this_node], MPI_BYTE,
&(res[0]), &(byte_counts[0]), &(displacements[0]), MPI_BYTE,
comm_cart);
std::vector<collision_struct> res = local_collision_queue;
Utils::Mpi::gather_buffer(res, comm_cart);
boost::mpi::broadcast(comm_cart, res, 0);

return res;
}
Expand Down Expand Up @@ -621,15 +598,10 @@ void handle_collisions() {
} else { // We consider the pair because one particle
// is local to the node and the other is local or ghost

// Calculate initial position for new vs, which is in the local node's
// Use initial position for new vs, which is in the local node's
// domain
// Vs is moved afterwards and resorted after all collision s are handled
// Use position of non-ghost colliding particle.
Vector3d initial_pos;
if (p1->l.ghost)
initial_pos = p2->r.p;
else
initial_pos = p1->r.p;
const Vector3d initial_pos{my_left[0], my_left[1], my_left[2]};

// If we are in the two vs mode
// Virtual site related to first particle in the collision
Expand Down
4 changes: 2 additions & 2 deletions testsuite/python/collision_detection.py
Original file line number Diff line number Diff line change
Expand Up @@ -575,8 +575,8 @@ def test_glue_to_surface_random(self):
sigma=0,
cutoff=0)

#@ut.skipIf(not espressomd.has_features("AngleHarmonic"),"Tests skipped because AngleHarmonic not compiled in")
def test_AngleHarmonic(self):
@ut.skipIf(not espressomd.has_features("BOND_ANGLE"), "Tests skipped because AngleHarmonic not compiled in")
def test_bind_three_particles(self):
# Setup particles
self.s.part.clear()
dx = np.array((1, 0, 0))
Expand Down