Introduction
Introduction Statistics Contact Development Disclaimer Help
tinitialize constants to double precision, clean up after test - sphere - GPU-b…
git clone git://src.adamsgaard.dk/sphere
Log
Files
Refs
LICENSE
---
commit b528d4d01d9a1c17ed5e9a86ec8d7dbd76083063
parent 21ed15ef593a831b2124a48bbc662b9a227f85ac
Author: Anders Damsgaard <[email protected]>
Date: Mon, 11 Aug 2014 14:47:58 +0200
initialize constants to double precision, clean up after test
Diffstat:
M src/contactsearch.cuh | 47 ++++++++++++++++++++---------…
M tests/fluid_particle_interaction.py | 2 +-
2 files changed, 32 insertions(+), 17 deletions(-)
---
diff --git a/src/contactsearch.cuh b/src/contactsearch.cuh
t@@ -466,15 +466,15 @@ __global__ void interact(
//uint4 bonds = dev_bonds_sorted[idx_a];
// Initiate shear friction loss rate at 0.0
- Float es_dot = 0.0f;
- Float ev_dot = 0.0f;
+ Float es_dot = 0.0;
+ Float ev_dot = 0.0;
// Initiate pressure on particle at 0.0
- Float p = 0.0f;
+ Float p = 0.0;
// Allocate memory for temporal force/torque vector values
- Float3 F = MAKE_FLOAT3(0.0f, 0.0f, 0.0f);
- Float3 T = MAKE_FLOAT3(0.0f, 0.0f, 0.0f);
+ Float3 F = MAKE_FLOAT3(0.0, 0.0, 0.0);
+ Float3 T = MAKE_FLOAT3(0.0, 0.0, 0.0);
// Apply linear elastic, frictional contact model to registered contac…
if (devC_params.contactmodel == 2 || devC_params.contactmodel == 3) {
t@@ -484,7 +484,8 @@ __global__ void interact(
Float4 x_b, distmod;
Float4 vel_a = dev_vel_sorted[idx_a];
Float4 angvel4_a = dev_angvel_sorted[idx_a];
- Float3 angvel_a = MAKE_FLOAT3(angvel4_a.x, angvel4_a.y, angvel4_a…
+ Float3 angvel_a = MAKE_FLOAT3(
+ angvel4_a.x, angvel4_a.y, angvel4_a.z);
// Loop over all possible contacts, and remove contacts
// that no longer are valid (delta_n > 0.0)
t@@ -608,8 +609,9 @@ __global__ void interact(
w_n = MAKE_FLOAT3(w_0_nx.x, w_0_nx.y, w_0_nx.z);
if (delta_w < 0.0f) {
w_0_force = contactLinear_wall(&F, &T, &es_dot, &ev_dot, &p, idx_a,
- radius_a, dev_vel_sorted, dev_angve…
- w_0_mvfd.y);
+ radius_a, dev_vel_sorted,
+ dev_angvel_sorted, w_n,
+ delta_w, w_0_mvfd.y);
}
// Lower wall (force on wall not stored)
t@@ -617,8 +619,9 @@ __global__ void interact(
w_n = MAKE_FLOAT3(0.0f, 0.0f, 1.0f);
if (delta_w < 0.0f) {
(void)contactLinear_wall(&F, &T, &es_dot, &ev_dot, &p, idx_a,
- radius_a, dev_vel_sorted, dev_angvel_sort…
- w_n, delta_w, 0.0f);
+ radius_a, dev_vel_sorted,
+ dev_angvel_sorted, w_n, delta_w,
+ 0.0);
}
t@@ -629,7 +632,9 @@ __global__ void interact(
w_n = MAKE_FLOAT3(w_1_nx.x, w_1_nx.y, w_1_nx.z);
if (delta_w < 0.0f) {
w_1_force = contactLinear_wall(&F, &T, &es_dot, &ev_dot, &p,
- idx_a, radius_a, dev_vel_sorted…
+ idx_a, radius_a,
+ dev_vel_sorted,
+ dev_angvel_sorted, w_n,
delta_w, w_1_mvfd.y);
}
t@@ -638,7 +643,9 @@ __global__ void interact(
w_n = MAKE_FLOAT3(w_2_nx.x, w_2_nx.y, w_2_nx.z);
if (delta_w < 0.0f) {
w_2_force = contactLinear_wall(&F, &T, &es_dot, &ev_dot, &p,
- idx_a, radius_a, dev_vel_sorted…
+ idx_a, radius_a,
+ dev_vel_sorted,
+ dev_angvel_sorted, w_n,
delta_w, w_2_mvfd.y);
}
t@@ -647,7 +654,9 @@ __global__ void interact(
w_n = MAKE_FLOAT3(w_3_nx.x, w_3_nx.y, w_3_nx.z);
if (delta_w < 0.0f) {
w_3_force = contactLinear_wall(&F, &T, &es_dot, &ev_dot, &p,
- idx_a, radius_a, dev_vel_sorted…
+ idx_a, radius_a,
+ dev_vel_sorted,
+ dev_angvel_sorted, w_n,
delta_w, w_3_mvfd.y);
}
t@@ -656,7 +665,9 @@ __global__ void interact(
w_n = MAKE_FLOAT3(w_4_nx.x, w_4_nx.y, w_4_nx.z);
if (delta_w < 0.0f) {
w_4_force = contactLinear_wall(&F, &T, &es_dot, &ev_dot, &p,
- idx_a, radius_a, dev_vel_sorted…
+ idx_a, radius_a,
+ dev_vel_sorted,
+ dev_angvel_sorted, w_n,
delta_w, w_4_mvfd.y);
}
t@@ -667,7 +678,9 @@ __global__ void interact(
w_n = MAKE_FLOAT3(w_3_nx.x, w_3_nx.y, w_3_nx.z);
if (delta_w < 0.0f) {
w_3_force = contactLinear_wall(&F, &T, &es_dot, &ev_dot, &p,
- idx_a, radius_a, dev_vel_sorted…
+ idx_a, radius_a,
+ dev_vel_sorted,
+ dev_angvel_sorted, w_n,
delta_w, w_3_mvfd.y);
}
t@@ -676,7 +689,9 @@ __global__ void interact(
w_n = MAKE_FLOAT3(w_4_nx.x, w_4_nx.y, w_4_nx.z);
if (delta_w < 0.0f) {
w_4_force = contactLinear_wall(&F, &T, &es_dot, &ev_dot, &p,
- idx_a, radius_a, dev_vel_sorted…
+ idx_a, radius_a,
+ dev_vel_sorted,
+ dev_angvel_sorted, w_n,
delta_w, w_4_mvfd.y);
}
}
diff --git a/tests/fluid_particle_interaction.py b/tests/fluid_particle_interac…
t@@ -50,4 +50,4 @@ test(sim.vel[0,0] > 0.0, 'Particle 0 velocity:')
test(sim.vel[1,0] > 0.0, 'Particle 1 velocity:')
test(sim.vel[2,0] > 0.0, 'Particle 2 velocity:')
-#sim.cleanup()
+sim.cleanup()
You are viewing proxied material from mx1.adamsgaard.dk. The copyright of proxied material belongs to its original authors. Any comments or complaints in relation to proxied material should be directed to the original authors of the content concerned. Please see the disclaimer for more details.