tset missing ghost nodes for F_pf and div_phi_vi_v - 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 2e8ed2ac27a09fd1d8cd08617169d289e1808081
 (DIR) parent 03a1e5cfa91cabac2788716303dc592da5d261f1
 (HTM) Author: Anders Damsgaard <anders.damsgaard@geo.au.dk>
       Date:   Mon,  9 Jun 2014 12:57:58 +0200
       
       set missing ghost nodes for F_pf and div_phi_vi_v
       
       Diffstat:
         M src/debug.h                         |       4 ++--
         M src/device.cu                       |      12 ++++++++++++
         M src/navierstokes.cuh                |       6 +++---
       
       3 files changed, 17 insertions(+), 5 deletions(-)
       ---
 (DIR) diff --git a/src/debug.h b/src/debug.h
       t@@ -40,10 +40,10 @@ const int conv_log_interval = 4;
        #define CHECK_NS_FINITE
        
        // Enable reporting of velocity prediction components to stdout
       -//#define REPORT_V_P_COMPONENTS
       +#define REPORT_V_P_COMPONENTS
        
        // Enable reporting of forcing function terms to stdout
       -//#define REPORT_FORCING_TERMS
       +#define REPORT_FORCING_TERMS
        
        // Choose solver model (see Zhou et al. 2010 "Discrete particle simulation of
        // particle-fluid flow: model formulations and their applicability", table. 1.
 (DIR) diff --git a/src/device.cu b/src/device.cu
       t@@ -970,6 +970,11 @@ __host__ void DEM::startTime()
                        cudaThreadSynchronize();
                        checkForCudaErrorsIter("Post applyInteractionForceToFluid",
                                iter);
       +
       +                setNSghostNodes<Float3><<<dimGridFluid, dimBlockFluid>>>(
       +                        dev_ns_F_pf, ns.bc_bot, ns.bc_top);
       +                cudaThreadSynchronize();
       +                checkForCudaErrorsIter("Post setNSghostNodes(F_pf)", iter);
                    }
        #endif
        
       t@@ -1075,6 +1080,13 @@ __host__ void DEM::startTime()
                                    &t_findNSdivphiviv);
                        checkForCudaErrorsIter("Post findNSdivphiviv", iter);
        
       +                // Set cell-center ghost nodes
       +                setNSghostNodes<Float3><<<dimGridFluid, dimBlockFluid>>>(
       +                    dev_ns_div_phi_vi_v, ns.bc_bot, ns.bc_top);
       +                cudaThreadSynchronize();
       +                checkForCudaErrorsIter("Post setNSghostNodes(div_phi_vi_v)",
       +                                       iter);
       +
                        // Predict the fluid velocities on the base of the old pressure
                        // field and ignoring the incompressibility constraint
                        if (PROFILING == 1)
 (DIR) diff --git a/src/navierstokes.cuh b/src/navierstokes.cuh
       t@@ -2139,8 +2139,8 @@ __global__ void findPredNSvelocities(
            const unsigned int cellidx = idx(x,y,z);
        
            // Check that we are not outside the fluid grid
       -    if (x <= nx && y <= ny && z <= nz) {
       -    //if (x < nx && y < ny && z < nz) {
       +    //if (x <= nx && y <= ny && z <= nz) {
       +    if (x < nx && y < ny && z < nz) {
        
                // Values that are needed for calculating the predicted velocity
                __syncthreads();
       t@@ -3005,7 +3005,7 @@ __global__ void findInteractionForce(
        
        #ifdef CHECK_NS_FINITE
                //*
       -        printf("\n%d [%d,%d,%d]\n"
       +        printf("\nfindInteractionForce %d [%d,%d,%d]\n"
                       "\tV_p = %f Re=%f Cd=%f chi=%f\n"
                       "\tf_d = %+e %+e %+e\n"
                       "\tf_p = %+e %+e %+e\n"