tBugfix-related changes - 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 c5711cbfa0d4ccc6b229dfbc50404a22966fe418
 (DIR) parent 6d5161287b1acd8c09dddf393ecc32e30cdfa2b3
 (HTM) Author: Anders Damsgaard <adc@geo.au.dk>
       Date:   Tue, 27 Nov 2012 13:45:08 +0100
       
       Bugfix-related changes
       
       Diffstat:
         M src/device.cu                       |      14 ++++++++++++++
         M src/integration.cuh                 |       4 ++--
         M src/sphere.cpp                      |       6 ++++--
       
       3 files changed, 20 insertions(+), 4 deletions(-)
       ---
 (DIR) diff --git a/src/device.cu b/src/device.cu
       t@@ -400,6 +400,19 @@ __host__ void DEM::transferToGlobalDeviceMemory()
          cudaMemcpy( dev_walls_mvfd,  walls.mvfd,
              sizeof(Float4)*walls.nw, cudaMemcpyHostToDevice);
        
       +  std::cout << "walls.nw = " << walls.nw << '\n';
       +  std::cout << "walls.nx[0] = "
       +      << walls.nx[0].x << '\t'
       +      << walls.nx[0].y << '\t'
       +      << walls.nx[0].w << '\t'
       +      << walls.nx[0].z << '\n';
       +
       +  std::cout << "walls.mvfd[0] = "
       +      << walls.mvfd[0].x << '\t'
       +      << walls.mvfd[0].y << '\t'
       +      << walls.mvfd[0].w << '\t'
       +      << walls.mvfd[0].z << '\n';
       +
          checkForCudaErrors("End of transferToGlobalDeviceMemory");
          if (verbose == 1)
            std::cout << "Done\n";
       t@@ -822,6 +835,7 @@ __host__ void DEM::startTime()
        
              filetimeclock = 0.0;
            }
       +    //break; // Stop after the first iteration
          }
        
          // Stop clock and display calculation time spent
 (DIR) diff --git a/src/integration.cuh b/src/integration.cuh
       t@@ -242,7 +242,7 @@ __global__ void integrateWalls(Float4* dev_walls_nx,
              return;
        
            // Find the final sum of forces on wall
       -    w_mvfd.z = 0.0f;
       +    w_mvfd.z = 0.0;
            for (int i=0; i<blocksPerGrid; ++i) {
              w_mvfd.z += dev_walls_force_partial[i];
            }
       t@@ -259,7 +259,7 @@ __global__ void integrateWalls(Float4* dev_walls_nx,
        
            // If Wall BC is controlled by velocity, it should not change
            if (wmode == 2) { 
       -      acc = 0.0f;
       +      acc = 0.0;
            }
        
            // Update position. Second-order scheme based on Taylor expansion 
 (DIR) diff --git a/src/sphere.cpp b/src/sphere.cpp
       t@@ -224,9 +224,11 @@ void DEM::reportValues()
            if (walls.wmode[0] == 0)
              cout << "Fixed\n";
            else if (walls.wmode[0] == 1)
       -      cout << "Deviatoric stress\n";
       +      cout << "Deviatoric stress, "
       +        << walls.mvfd[0].w << " Pa\n";
            else if (walls.wmode[0] == 2)
       -      cout << "Velocity\n";
       +      cout << "Velocity, "
       +        << walls.mvfd[0].y << " m/s\n";
            else {
              cerr << "Top BC not recognized!\n";
              exit(1);