URI:
       tstgarted implementing slab on a slope, no slip boundaries missing - sphere - GPU-based 3D discrete element method algorithm with optional fluid coupling
  HTML git clone git://src.adamsgaard.dk/sphere
   DIR Log
   DIR Files
   DIR Refs
   DIR LICENSE
       ---
   DIR commit e6024e9e3bc8ef33c96b2d38c7e23a62dc190e4d
   DIR parent 983b9baa877f7dac3ef84c42ab16f60feb5e5de4
  HTML Author: Anders Damsgaard <anders.damsgaard@geo.au.dk>
       Date:   Wed, 18 Jun 2014 11:43:59 +0200
       
       stgarted implementing slab on a slope, no slip boundaries missing
       
       Diffstat:
         M src/debug.h                         |       2 +-
         M src/navierstokes.cuh                |       4 +---
         A tests/cfd_inclined.py               |      43 ++++++++++++++++++++++++++++++
         M tests/cfd_simple.py                 |      13 ++++++++++---
         M tests/cfd_tests_neumann.py          |       2 +-
       
       5 files changed, 56 insertions(+), 8 deletions(-)
       ---
   DIR diff --git a/src/debug.h b/src/debug.h
       t@@ -51,7 +51,7 @@ const int conv_log_interval = 4;
        // simulation of particulate systems: Theoretical developments".
        // SET_2 corresponds approximately to Model A in Zhu et al. 2007.
        // Choose exactly one.
       -//#define SET_1
       +//#define SET_1 // set 1 not functional!
        #define SET_2
        
        #endif
   DIR diff --git a/src/navierstokes.cuh b/src/navierstokes.cuh
       t@@ -2270,7 +2270,6 @@ __global__ void findPredNSvelocities(
                       "\tgrav     = %+e %+e %+e\n"
                       "\tporos    = %+e %+e %+e\n"
                       "\tadv      = %+e %+e %+e\n"
       -               "\tdiv_tau  = %+e %+e %+e\n",
                       x, y, z,
                       v_p.x, v_p.y, v_p.z,
                       pressure_term.x, pressure_term.y, pressure_term.z, 
       t@@ -2278,8 +2277,7 @@ __global__ void findPredNSvelocities(
                       diffusion_term.x, diffusion_term.y, diffusion_term.z, 
                       gravity_term.x, gravity_term.y, gravity_term.z, 
                       porosity_term.x, porosity_term.y, porosity_term.z, 
       -               advection_term.x, advection_term.y, advection_term.z,
       -               div_tau.x, div_tau.y, div_tau.z);
       +               advection_term.x, advection_term.y, advection_term.z);
        #endif
        
                // Enforce Neumann BC if specified
   DIR diff --git a/tests/cfd_inclined.py b/tests/cfd_inclined.py
       t@@ -0,0 +1,43 @@
       +#!/usr/bin/env python
       +import sphere
       +from pytestutils import *
       +
       +orig = sphere.sim('cfd_incl', fluid=True)
       +orig.cleanup()
       +#orig.defineWorldBoundaries([0.3, 0.3, 0.3], dx = 0.1)
       +orig.defineWorldBoundaries([0.3, 0.3, 0.3], dx = 0.06)
       +orig.initFluid(mu=8.9e-4) # inviscid "fluids" (mu=0) won't work!
       +#orig.initTemporal(total = 0.5, file_dt = 0.05, dt = 1.0e-4)
       +orig.initTemporal(total = 1.0e-0, file_dt = 1.0e-1, dt = 1.0e-3)
       +orig.bc_bot[0] = 1 # No-flow BC at bottom (Neumann)
       +orig.bc_top[0] = 1 # No-flow BC at top (Neumann)
       +
       +angle = 10.0 # slab inclination in degrees
       +g_magnitude = 10.0
       +orig.g[0] =  numpy.sin(numpy.radians(angle))*g_magnitude
       +orig.g[2] = -numpy.cos(numpy.radians(angle))*g_magnitude
       +
       +tau_d = orig.g * orig.rho_f * orig.L[3] # analytical driving stress
       +v_sur = tau_d * orig.L[3] / orig.mu     # analytical surface velocity
       +
       +# increase the max iterations for first step
       +orig.setMaxIterations(1e5)
       +
       +# Homogeneous pressure, no gravity
       +orig.run(verbose=False)
       +orig.writeVTKall()
       +
       +py = sphere.sim(sid=orig.sid, fluid=True)
       +py.readlast(verbose = False)
       +ones = numpy.ones((orig.num))
       +zeros = numpy.zeros((orig.num[0], orig.num[1], orig.num[2], 3))
       +#compareNumpyArraysClose(ones, py.p_f, "Conservation of pressure:",
       +        #tolerance = 1.0e-5)
       +#compareNumpyArraysClose(zeros, py.v_f, "Flow field:              ",
       +        #tolerance = 1.0e-5)
       +#ideal_grad_p_z = numpy.linspace(
       +        #orig.p_f[0,0,0] + orig.L[2]*orig.rho_f*numpy.abs(orig.g[2]),
       +        #orig.p_f[0,0,-1], orig.num[2])
       +#compareNumpyArraysClose(ideal_grad_p_z, py.p_f[0,0,:],
       +        #"Pressure gradient:\t", tolerance=1.0e2)
       +#orig.cleanup()
   DIR diff --git a/tests/cfd_simple.py b/tests/cfd_simple.py
       t@@ -4,12 +4,14 @@ from pytestutils import *
        
        orig = sphere.sim('cfd_simple', fluid=True)
        orig.cleanup()
       -orig.defineWorldBoundaries([0.3, 0.3, 0.3], dx = 0.1)
       +#orig.defineWorldBoundaries([0.3, 0.3, 0.3], dx = 0.1)
       +orig.defineWorldBoundaries([0.3, 0.3, 0.3], dx = 0.06)
       +#orig.initFluid(mu=0.0)
        orig.initFluid(mu=8.9e-4)
        #orig.initTemporal(total = 0.5, file_dt = 0.05, dt = 1.0e-4)
       -orig.initTemporal(total = 1.0e-3, file_dt = 1.0e-3, dt = 1.0e-3)
       +orig.initTemporal(total = 1.0e-0, file_dt = 1.0e-1, dt = 1.0e-3)
        orig.bc_bot[0] = 1      # No-flow BC at bottom (Neumann)
       -#orig.g[2] = -10.0
       +orig.g[2] = -10.0
        
        # Homogeneous pressure, no gravity
        orig.run(verbose=False)
       t@@ -23,4 +25,9 @@ zeros = numpy.zeros((orig.num[0], orig.num[1], orig.num[2], 3))
                #tolerance = 1.0e-5)
        #compareNumpyArraysClose(zeros, py.v_f, "Flow field:              ",
                #tolerance = 1.0e-5)
       +ideal_grad_p_z = numpy.linspace(
       +        orig.p_f[0,0,0] + orig.L[2]*orig.rho_f*numpy.abs(orig.g[2]),
       +        orig.p_f[0,0,-1], orig.num[2])
       +compareNumpyArraysClose(ideal_grad_p_z, py.p_f[0,0,:],
       +        "Pressure gradient:\t", tolerance=1.0e2)
        #orig.cleanup()
   DIR diff --git a/tests/cfd_tests_neumann.py b/tests/cfd_tests_neumann.py
       t@@ -26,7 +26,7 @@ py.readlast(verbose = False)
        ones = numpy.ones((orig.num))
        py.readlast(verbose = False)
        compareNumpyArraysClose(ones, py.p_f, "Conservation of pressure:",
       -        tolerance = 1.0e-3)
       +        tolerance = 1.0e-5)
        
        # Fluid flow along z should be very small
        if ((numpy.abs(py.v_f[:,:,:,:]) < 1.0e-6).all()):