trename wall shear stress from tau to tau_x - sphere - GPU-based 3D discrete el… | |
git clone git://src.adamsgaard.dk/sphere | |
Log | |
Files | |
Refs | |
LICENSE | |
--- | |
commit cbdda5df5cf0c5d7ce15213aa4591a78737b6b4e | |
parent 3c5a8efe14ab698d5695f0ec74d9f0e5074cfdd6 | |
Author: Anders Damsgaard <[email protected]> | |
Date: Wed, 14 Jan 2015 12:22:46 +0100 | |
rename wall shear stress from tau to tau_x | |
Diffstat: | |
M python/sphere.py | 12 +++++------- | |
M src/datatypes.h | 2 +- | |
M src/device.cu | 7 +++++++ | |
M src/file_io.cpp | 6 +++--- | |
M src/sphere.h | 1 + | |
5 files changed, 17 insertions(+), 11 deletions(-) | |
--- | |
diff --git a/python/sphere.py b/python/sphere.py | |
t@@ -251,7 +251,7 @@ class sim: | |
self.w_sigma0_f = numpy.zeros(1, dtype=numpy.float64) | |
# Wall shear stress, enforced when wmode == 3 | |
- self.w_tau = numpy.zeros(1, dtype=numpy.float64) | |
+ self.w_tau_x = numpy.zeros(1, dtype=numpy.float64) | |
## Bond parameters | |
# Radius multiplier to the parallel-bond radii | |
t@@ -564,7 +564,7 @@ class sim: | |
elif (self.w_sigma0_f != other.w_sigma0_f): | |
print(52) | |
return 52 | |
- elif (self.w_tau != other.w_tau): | |
+ elif (self.w_tau_x != other.w_tau_x): | |
print(52.5) | |
return 52.5 | |
elif (self.gamma_wn != other.gamma_wn): | |
t@@ -1041,7 +1041,7 @@ class sim: | |
self.w_sigma0_A = numpy.fromfile(fh, dtype=numpy.float64, coun… | |
self.w_sigma0_f = numpy.fromfile(fh, dtype=numpy.float64, coun… | |
if self.version >= 2.1: | |
- self.w_tau = numpy.fromfile(fh, dtype=numpy.float64, count=1) | |
+ self.w_tau_x = numpy.fromfile(fh, dtype=numpy.float64, count=1) | |
if bonds: | |
# Inter-particle bonds | |
t@@ -1330,7 +1330,7 @@ class sim: | |
fh.write(self.w_sigma0[i].astype(numpy.float64)) | |
fh.write(self.w_sigma0_A.astype(numpy.float64)) | |
fh.write(self.w_sigma0_f.astype(numpy.float64)) | |
- fh.write(self.w_tau.astype(numpy.float64)) | |
+ fh.write(self.w_tau_x.astype(numpy.float64)) | |
fh.write(self.lambda_bar.astype(numpy.float64)) | |
fh.write(self.nb0.astype(numpy.uint32)) | |
t@@ -5349,15 +5349,13 @@ class sim: | |
def shearStress(self): | |
''' | |
Calculates the sum of shear stress values measured on any moving | |
- particles. | |
+ particles with a finite and fixed velocity. | |
:returns: The shear stress in Pa | |
:return type: numpy.array | |
''' | |
fixvel = numpy.nonzero(self.fixvel > 0.0) | |
- shearvel = self.shearVel() | |
- | |
force = numpy.zeros(3) | |
# Summation of shear stress contributions | |
diff --git a/src/datatypes.h b/src/datatypes.h | |
t@@ -104,7 +104,7 @@ struct Walls { | |
int wmode[MAXWALLS]; // Wall modes | |
Float4* nx; // Wall normal and position | |
Float4* mvfd; // Wall mass, velocity, force and normal stress | |
- Float* tau; // Wall shear stress | |
+ Float* tau_x; // Wall shear stress | |
}; | |
// Structures containing fluid parameters | |
diff --git a/src/device.cu b/src/device.cu | |
t@@ -16,6 +16,7 @@ | |
#include "utility.h" | |
#include "constants.cuh" | |
#include "debug.h" | |
+#include "version.h" | |
#include "sorting.cuh" | |
#include "contactmodels.cuh" | |
t@@ -380,6 +381,7 @@ __host__ void DEM::allocateGlobalDeviceMemory(void) | |
cudaMalloc((void**)&dev_walls_wmode, sizeof(int)*walls.nw); | |
cudaMalloc((void**)&dev_walls_nx, sizeof(Float4)*walls.nw); | |
cudaMalloc((void**)&dev_walls_mvfd, sizeof(Float4)*walls.nw); | |
+ cudaMalloc((void**)&dev_walls_tau_x, sizeof(Float)*walls.nw); | |
cudaMalloc((void**)&dev_walls_force_pp, sizeof(Float)*walls.nw*np); | |
cudaMalloc((void**)&dev_walls_acc, sizeof(Float)*walls.nw); | |
// dev_walls_force_partial allocated later | |
t@@ -545,6 +547,7 @@ __host__ void DEM::freeGlobalDeviceMemory() | |
// Wall arrays | |
cudaFree(dev_walls_nx); | |
cudaFree(dev_walls_mvfd); | |
+ cudaFree(dev_walls_tau_x); | |
cudaFree(dev_walls_force_partial); | |
cudaFree(dev_walls_force_pp); | |
cudaFree(dev_walls_acc); | |
t@@ -631,6 +634,8 @@ __host__ void DEM::transferToGlobalDeviceMemory(int status… | |
sizeof(Float4)*walls.nw, cudaMemcpyHostToDevice); | |
cudaMemcpy( dev_walls_mvfd, walls.mvfd, | |
sizeof(Float4)*walls.nw, cudaMemcpyHostToDevice); | |
+ cudaMemcpy( dev_walls_tau_x, walls.tau_x, | |
+ sizeof(Float)*walls.nw, cudaMemcpyHostToDevice); | |
// Fluid arrays | |
if (fluid == 1) { | |
t@@ -711,6 +716,8 @@ __host__ void DEM::transferFromGlobalDeviceMemory() | |
sizeof(Float4)*walls.nw, cudaMemcpyDeviceToHost); | |
cudaMemcpy( walls.mvfd, dev_walls_mvfd, | |
sizeof(Float4)*walls.nw, cudaMemcpyDeviceToHost); | |
+ cudaMemcpy( walls.tau_x, dev_walls_tau_x, | |
+ sizeof(Float)*walls.nw, cudaMemcpyDeviceToHost); | |
// Fluid arrays | |
if (fluid == 1 && cfd_solver == 0) { | |
diff --git a/src/file_io.cpp b/src/file_io.cpp | |
t@@ -202,7 +202,7 @@ void DEM::readbin(const char *target) | |
// Wall mass (x), velocity (y), force (z), and deviatoric stress (w) | |
walls.nx = new Float4[walls.nw]; | |
walls.mvfd = new Float4[walls.nw]; | |
- walls.tau = new Float[walls.nw]; | |
+ walls.tau_x = new Float[walls.nw]; | |
for (i = 0; i<walls.nw; ++i) | |
ifs.read(as_bytes(walls.wmode[i]), sizeof(walls.wmode[i])); | |
t@@ -220,7 +220,7 @@ void DEM::readbin(const char *target) | |
} | |
ifs.read(as_bytes(params.sigma0_A), sizeof(params.sigma0_A)); | |
ifs.read(as_bytes(params.sigma0_f), sizeof(params.sigma0_f)); | |
- ifs.read(as_bytes(walls.tau), sizeof(walls.tau)); | |
+ ifs.read(as_bytes(walls.tau_x), sizeof(walls.tau_x)); | |
// Read bond parameters | |
ifs.read(as_bytes(params.lambda_bar), sizeof(params.lambda_bar)); | |
t@@ -522,7 +522,7 @@ void DEM::writebin(const char *target) | |
} | |
ofs.write(as_bytes(params.sigma0_A), sizeof(params.sigma0_A)); | |
ofs.write(as_bytes(params.sigma0_f), sizeof(params.sigma0_f)); | |
- ofs.write(as_bytes(walls.tau), sizeof(walls.tau)); | |
+ ofs.write(as_bytes(walls.tau_x), sizeof(walls.tau_x)); | |
// Write bond parameters | |
ofs.write(as_bytes(params.lambda_bar), sizeof(params.lambda_bar)); | |
diff --git a/src/sphere.h b/src/sphere.h | |
t@@ -97,6 +97,7 @@ class DEM { | |
int *dev_walls_wmode; | |
Float4 *dev_walls_nx; // normal, pos. | |
Float4 *dev_walls_mvfd; // mass, velo., force, dev. stress | |
+ Float *dev_walls_tau_x; // wall shear stress | |
Float *dev_walls_force_partial; // Pre-sum per wall | |
Float *dev_walls_force_pp; // Force per particle per wall | |
Float *dev_walls_acc; // Wall acceleration |