tdouble precision comparison - sphere - GPU-based 3D discrete element method algorithm with optional fluid coupling
 (HTM) git clone git://src.adamsgaard.dk/sphere
 (DIR) Log
 (DIR) Files
 (DIR) Refs
 (DIR) LICENSE
       ---
 (DIR) commit 6de1181f3e695a80fe56a139e3384cba82034404
 (DIR) parent a3f8391f6d7702842e46e8ead63ad008d89106f3
 (HTM) Author: Anders Damsgaard <anders.damsgaard@geo.au.dk>
       Date:   Mon, 30 Jun 2014 10:20:05 +0200
       
       double precision comparison
       
       Diffstat:
         M src/contactsearch.cuh               |      24 ++++++++++++------------
       
       1 file changed, 12 insertions(+), 12 deletions(-)
       ---
 (DIR) diff --git a/src/contactsearch.cuh b/src/contactsearch.cuh
       t@@ -381,19 +381,19 @@ __global__ void interact(
                Float4* dev_angvel_sorted,           // in
                Float4* dev_vel,                     // in
                Float4* dev_angvel,                  // in
       -        Float4* dev_force,         // out
       -        Float4* dev_torque,        // out
       -        Float* dev_es_dot,         // out
       -        Float* dev_ev_dot,         // out
       -        Float* dev_es,             // out
       -        Float* dev_ev,             // out
       -        Float* dev_p,              // out
       +        Float4* dev_force,          // out
       +        Float4* dev_torque,         // out
       +        Float* dev_es_dot,          // out
       +        Float* dev_ev_dot,          // out
       +        Float* dev_es,              // out
       +        Float* dev_ev,              // out
       +        Float* dev_p,               // out
                Float4* dev_walls_nx,                // in
                Float4* dev_walls_mvfd,              // in
       -        Float* dev_walls_force_pp, // out
       -        unsigned int* dev_contacts,          // in
       +        Float* dev_walls_force_pp,  // out
       +        unsigned int* dev_contacts, // out
                Float4* dev_distmod,                 // in
       -        Float4* dev_delta_t)                 // in
       +        Float4* dev_delta_t)        // out
        {
            // Thread index equals index of particle A
            unsigned int idx_a = blockIdx.x * blockDim.x + threadIdx.x;
       t@@ -491,7 +491,7 @@ __global__ void interact(
                        if (idx_b_orig != (unsigned int)devC_np) {
        
                            // Read inter-particle distance correction vector
       -                    distmod    = dev_distmod[mempos];
       +                    distmod = dev_distmod[mempos];
        
                            // Read particle b position and radius
                            x_b = dev_x[idx_b_orig];
       t@@ -506,7 +506,7 @@ __global__ void interact(
                            delta_n = x_ab_length - (radius_a + radius_b);
        
                            // Process collision if the particles are overlapping
       -                    if (delta_n < 0.0f) {
       +                    if (delta_n < 0.0) {
                                if (devC_params.contactmodel == 2) {
                                    contactLinear(&F, &T, &es_dot, &ev_dot, &p, 
                                            idx_a_orig,