tenforce dirichlet pressure BC when bc_top == 3 - sphere - GPU-based 3D discret… | |
git clone git://src.adamsgaard.dk/sphere | |
Log | |
Files | |
Refs | |
LICENSE | |
--- | |
commit 4b3e865a3c097b32ea168764b5953c9b93df2cee | |
parent 456f57069bb6aa097d784fa4bff2075a5baa14c7 | |
Author: Anders Damsgaard <[email protected]> | |
Date: Mon, 19 Jan 2015 09:31:01 +0100 | |
enforce dirichlet pressure BC when bc_top == 3 | |
Diffstat: | |
M python/halfshear-darcy-stress-star… | 4 ++-- | |
M src/darcy.cuh | 6 +++--- | |
M src/device.cu | 2 +- | |
3 files changed, 6 insertions(+), 6 deletions(-) | |
--- | |
diff --git a/python/halfshear-darcy-stress-starter.py b/python/halfshear-darcy-… | |
t@@ -71,8 +71,8 @@ sim.setDampingTangential(0.0) | |
#sim.deleteAllParticles() | |
#sim.fixvel[:] = -1.0 | |
-#sim.initTemporal(total = 20.0, file_dt = 0.01, epsilon=0.07) | |
-sim.initTemporal(total = 20.0, file_dt = 0.00001, epsilon=0.07) | |
+sim.initTemporal(total = 20.0, file_dt = 0.01, epsilon=0.07) | |
+#sim.initTemporal(total = 20.0, file_dt = 0.00001, epsilon=0.07) | |
#sim.time_dt[0] *= 1.0e-2 | |
#sim.initTemporal(total = 1.0e-4, file_dt = 1.0e-5, epsilon=0.07) | |
diff --git a/src/darcy.cuh b/src/darcy.cuh | |
t@@ -240,7 +240,7 @@ __global__ void setDarcyGhostNodes( | |
if (z == 0 && bc_bot == 2) | |
dev_scalarfield[idx(x,y,nz)] = val; // Periodic -z | |
- if (z == nz-1 && bc_top == 0) | |
+ if (z == nz-1 && (bc_top == 0 || bc_top == 3)) | |
dev_scalarfield[idx(x,y,nz)] = val; // Dirichlet | |
if (z == nz-2 && bc_top == 1) | |
dev_scalarfield[idx(x,y,nz)] = val; // Neumann | |
t@@ -938,7 +938,7 @@ __global__ void firstDarcySolution( | |
if (z == 0 && bc_bot == 1) | |
p_zn = p; | |
//if ((z == nz-1 && bc_top == 1) || z >= wall0_iz) | |
- if (z == nz-1 && bc_top == 1) | |
+ if (z == nz-1 && (bc_top == 1 || bc_top == 3)) | |
p_zp = p; | |
// upwind coefficients for grad(p) determined from values of k | |
t@@ -1089,7 +1089,7 @@ __global__ void updateDarcySolution( | |
if (z == 0 && bc_bot == 1) | |
p_zn = p; | |
//if ((z == nz-1 && bc_top == 1) || z >= wall0_iz) | |
- if (z == nz-1 && bc_top == 1) | |
+ if (z == nz-1 && (bc_top == 1 || bc_top == 3)) | |
p_zp = p; | |
// upwind coefficients for grad(p) determined from values of k | |
diff --git a/src/device.cu b/src/device.cu | |
t@@ -2037,7 +2037,7 @@ __host__ void DEM::startTime() | |
checkForCudaErrorsIter("Post updateDarcySolution", | |
iter); | |
- if (darcy.bc_top == 1) { | |
+ if (darcy.bc_top == 1 || darcy.bc_top == 3) { | |
if (PROFILING == 1) | |
startTimer(&kernel_tic); | |
setDarcyTopWallFixedFlow |