tenabled debugging, small error in printf - 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 d843fd95007da8e532a444c0e3aca8df9010fd6f
(DIR) parent 2184bb21b72c8c0b61d8fee27347fdc2ba40a48e
(HTM) Author: Anders Damsgaard <anders.damsgaard@geo.au.dk>
Date: Fri, 23 May 2014 11:11:05 +0200
enabled debugging, small error in printf
Diffstat:
M src/debug.h | 12 ++++++------
M src/navierstokes.cuh | 12 ++++++------
2 files changed, 12 insertions(+), 12 deletions(-)
---
(DIR) 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@@ -40,10 +40,10 @@ const int conv_log_interval = 10;
#define CHECK_NS_FINITE
// Enable reporting of forcing function terms to stdout
-//#define REPORT_FORCING_TERMS
+#define REPORT_FORCING_TERMS
// Enable reporting of velocity prediction components to stdout
-//#define REPORT_V_P_COMPONENTS
+#define REPORT_V_P_COMPONENTS
// Choose solver model (see Zhou et al. 2010 "Discrete particle simulation of
// particle-fluid flow: model formulations and their applicability", table. 1.
t@@ -51,7 +51,7 @@ const int conv_log_interval = 10;
// 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_2
+//#define SET_1
+#define SET_2
#endif
(DIR) diff --git a/src/navierstokes.cuh b/src/navierstokes.cuh
t@@ -2108,19 +2108,19 @@ __global__ void findPredNSvelocities(
const Float3 diffusion_term = dt/(rho*phi)*div_tau;
const Float3 gravity_term = MAKE_FLOAT3(
devC_params.g[0], devC_params.g[1], devC_params.g[2])*dt;
- const Float3 porosity_term = -v*dphi/phi;
- const Float3 advection_term = -div_phi_vi_v*dt/phi;
+ const Float3 porosity_term = -1.0*v*dphi/phi;
+ const Float3 advection_term = -1.0*div_phi_vi_v*dt/phi;
#endif
#ifdef SET_2
const Float3 interaction_term = -dt/(rho*phi)*f_i;
const Float3 diffusion_term = dt/rho*div_tau;
const Float3 gravity_term = MAKE_FLOAT3(
devC_params.g[0], devC_params.g[1], devC_params.g[2])*dt;
- const Float3 porosity_term = -v*dphi/phi;
- const Float3 advection_term = -div_phi_vi_v*dt/phi;
+ const Float3 porosity_term = -1.0*v*dphi/phi;
+ const Float3 advection_term = -1.0*div_phi_vi_v*dt/phi;
#endif
- Float v_p = v
+ Float3 v_p = v
+ pressure_term
+ interaction_term
+ diffusion_term
t@@ -2137,7 +2137,7 @@ __global__ void findPredNSvelocities(
"diff = %e\t%e\t%e\t"
"grav = %e\t%e\t%e\t"
"poros = %e\t%e\t%e\t"
- "adv = %e\t%e\t%e\n"
+ "adv = %e\t%e\t%e\n",
x, y, z,
v_p.x, v_p.y, v_p.z,
pressure_term.x, pressure_term.y, pressure_term.z,