| tcorrect velocity gradient at bottom neumann BC in staggered grid, show v_p for… | |
| git clone git://src.adamsgaard.dk/sphere | |
| Log | |
| Files | |
| Refs | |
| LICENSE | |
| --- | |
| commit ff67b9fd391caf9576d5c437bb92a14d914d63bd | |
| parent bbc594bbe6b8c9dcbfa4e229fa28a6a1f6248959 | |
| Author: Anders Damsgaard <[email protected]> | |
| Date: Wed, 22 Oct 2014 10:35:56 +0200 | |
| correct velocity gradient at bottom neumann BC in staggered grid, show v_p for … | |
| Diffstat: | |
| M src/debug.h | 6 +++--- | |
| M src/navierstokes.cuh | 36 ++++++++++++++++-------------… | |
| 2 files changed, 21 insertions(+), 21 deletions(-) | |
| --- | |
| diff --git a/src/debug.h b/src/debug.h | |
| t@@ -21,8 +21,8 @@ const unsigned int nijacnorm = 10; | |
| const int write_res_log = 0; | |
| // Report epsilon values during Jacobi iterations to stdout | |
| -#define REPORT_EPSILON | |
| -#define REPORT_MORE_EPSILON | |
| +//#define REPORT_EPSILON | |
| +//#define REPORT_MORE_EPSILON | |
| // Report the number of iterations it took before convergence to logfile | |
| // 'output/<sid>-conv.dat' | |
| t@@ -47,7 +47,7 @@ const int conv_log_interval = 1; | |
| #define REPORT_V_C_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. | |
| diff --git a/src/navierstokes.cuh b/src/navierstokes.cuh | |
| t@@ -692,7 +692,8 @@ __global__ void setNSghostNodes( | |
| // z | |
| if (z == 0 && bc_bot == 0) | |
| dev_scalarfield[idx(x,y,-1)] = val; // Dirichlet | |
| - if (z == 1 && bc_bot == 1) | |
| + //if (z == 1 && bc_bot == 1) | |
| + if (z == 0 && bc_bot == 1) | |
| dev_scalarfield[idx(x,y,-1)] = val; // Neumann | |
| if (z == 0 && bc_bot == 2) | |
| dev_scalarfield[idx(x,y,nz)] = val; // Periodic -z | |
| t@@ -2450,23 +2451,22 @@ __global__ void findPredNSvelocities( | |
| #ifdef REPORT_V_P_COMPONENTS | |
| // Report velocity components to stdout for debugging | |
| - if (z==0) | |
| - printf("\n[%d,%d,%d]" | |
| - "\tv_p = %+e %+e %+e\n" | |
| - "\tpres = %+e %+e %+e\n" | |
| - "\tinteract = %+e %+e %+e\n" | |
| - "\tdiff = %+e %+e %+e\n" | |
| - "\tgrav = %+e %+e %+e\n" | |
| - "\tporos = %+e %+e %+e\n" | |
| - "\tadv = %+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, | |
| - interaction_term.x, interaction_term.y, interaction_term.z, | |
| - 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); | |
| + printf("\n[%d,%d,%d]" | |
| + "\tv_p = %+e %+e %+e\n" | |
| + "\tpres = %+e %+e %+e\n" | |
| + "\tinteract = %+e %+e %+e\n" | |
| + "\tdiff = %+e %+e %+e\n" | |
| + "\tgrav = %+e %+e %+e\n" | |
| + "\tporos = %+e %+e %+e\n" | |
| + "\tadv = %+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, | |
| + interaction_term.x, interaction_term.y, interaction_term.z, | |
| + 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); | |
| #endif | |
| // Save the predicted velocity |