Introduction
Introduction Statistics Contact Development Disclaimer Help
tMerge branch 'master' of github.com:anders-dc/sphere - sphere - GPU-based 3D d…
git clone git://src.adamsgaard.dk/sphere
Log
Files
Refs
LICENSE
---
commit 7651aa467789fd8a7e2b916d002ca92555d332e5
parent d1f1f06606be74e7eaed6f43e868b4515df28c43
Author: Anders Damsgaard Christensen <[email protected]>
Date: Fri, 12 Aug 2016 10:47:16 -0700
Merge branch 'master' of github.com:anders-dc/sphere
Diffstat:
M src/darcy.cuh | 79 +++++++++++++++++++++++++++++…
M src/device.cu | 34 ++++++++++++++++++++++++++---…
2 files changed, 108 insertions(+), 5 deletions(-)
---
diff --git a/src/darcy.cuh b/src/darcy.cuh
t@@ -716,6 +716,85 @@ __global__ void findDarcyPorositiesLinear(
}
+// Copy the porosity, porosity change, div_v_p and vp_avg values to the grid
+// edges from the grid interior at the frictionless Y boundaries (grid.periodi…
+// == 2).
+__global__ void copyDarcyPorositiesToEdges(
+ Float* __restrict__ dev_darcy_phi, // in + out
+ Float* __restrict__ dev_darcy_dphi, // in + out
+ Float* __restrict__ dev_darcy_div_v_p, // in + out
+ Float3* __restrict__ dev_darcy_vp_avg) // in + out
+{
+ // 3D thread index
+ const unsigned int x = blockDim.x * blockIdx.x + threadIdx.x;
+ const unsigned int y = blockDim.y * blockIdx.y + threadIdx.y;
+ const unsigned int z = blockDim.z * blockIdx.z + threadIdx.z;
+
+ // Grid dimensions
+ const unsigned int nx = devC_grid.num[0];
+ const unsigned int ny = devC_grid.num[1];
+ const unsigned int nz = devC_grid.num[2];
+
+ // check that we are not outside the fluid grid
+ if (devC_grid.periodic == 2 &&
+ x < nx && (y == 0 || y == ny - 1) && z < nz) {
+
+ // Read porosities from this cell
+ int y_read;
+
+ // read values from inside cells
+ if (y == 0)
+ y_read = 1;
+ if (y == ny - 1)
+ y_read = ny - 2;
+
+ const unsigned int readidx = d_idx(x, y_read, z);
+ const unsigned int writeidx = d_idx(x, y, z);
+
+ __syncthreads();
+ dev_darcy_phi[writeidx] = dev_darcy_phi[readidx];
+ dev_darcy_dphi[writeidx] = dev_darcy_dphi[readidx];
+ dev_darcy_div_v_p[writeidx] = dev_darcy_div_v_p[readidx];
+ dev_darcy_vp_avg[writeidx] = dev_darcy_vp_avg[readidx];
+ }
+}
+
+
+// Copy the porosity, porosity change, div_v_p and vp_avg values to the grid
+// edges from the grid interior at the frictionless Y boundaries (grid.periodi…
+// == 2).
+__global__ void copyDarcyPorositiesToBottom(
+ Float* __restrict__ dev_darcy_phi, // in + out
+ Float* __restrict__ dev_darcy_dphi, // in + out
+ Float* __restrict__ dev_darcy_div_v_p, // in + out
+ Float3* __restrict__ dev_darcy_vp_avg) // in + out
+{
+ // 3D thread index
+ const unsigned int x = blockDim.x * blockIdx.x + threadIdx.x;
+ const unsigned int y = blockDim.y * blockIdx.y + threadIdx.y;
+ const unsigned int z = blockDim.z * blockIdx.z + threadIdx.z;
+
+ // Grid dimensions
+ const unsigned int nx = devC_grid.num[0];
+ const unsigned int ny = devC_grid.num[1];
+ const unsigned int nz = devC_grid.num[2];
+
+ // check that we are not outside the fluid grid
+ if (devC_grid.periodic == 2 &&
+ x < nx && y < ny && z == 0) {
+
+ const unsigned int readidx = d_idx(x, y, 1);
+ const unsigned int writeidx = d_idx(x, y, z);
+
+ __syncthreads();
+ dev_darcy_phi[writeidx] = dev_darcy_phi[readidx];
+ dev_darcy_dphi[writeidx] = dev_darcy_dphi[readidx];
+ dev_darcy_div_v_p[writeidx] = dev_darcy_div_v_p[readidx];
+ dev_darcy_vp_avg[writeidx] = dev_darcy_vp_avg[readidx];
+ }
+}
+
+
// Find the porosity in each cell on the base of a sphere, centered at the cell
// center.
__global__ void findDarcyPorosities(
diff --git a/src/device.cu b/src/device.cu
t@@ -37,11 +37,13 @@ int cudaCoresPerSM(int major, int minor)
return 32;
else if (major == 2 && minor == 1)
return 48;
- else if (major == 3 && minor == 0)
+ else if (major == 3)
return 192;
- else if (major == 3 && minor == 5)
- return 192;
- else if (major == 5 && minor == 0)
+ else if (major == 5)
+ return 128;
+ else if (major == 6 && minor == 0)
+ return 64;
+ else if (major == 6 && minor == 1)
return 128;
else
printf("Error in cudaCoresPerSM",
t@@ -1963,6 +1965,28 @@ __host__ void DEM::startTime()
&t_findDarcyPorosities);
checkForCudaErrorsIter("Post findDarcyPorosities", iter);
+ // copy porosities to the frictionless Y boundaries
+ if (grid.periodic == 2) {
+ copyDarcyPorositiesToEdges<<<dimGridFluid,
+ dimBlockFluid>>>(
+ dev_darcy_phi,
+ dev_darcy_dphi,
+ dev_darcy_div_v_p,
+ dev_darcy_vp_avg);
+ cudaThreadSynchronize();
+ }
+
+ // copy porosities to the frictionless lower Z boundary
+ if (grid.periodic == 2) {
+ copyDarcyPorositiesToBottom<<<dimGridFluid,
+ dimBlockFluid>>>(
+ dev_darcy_phi,
+ dev_darcy_dphi,
+ dev_darcy_div_v_p,
+ dev_darcy_vp_avg);
+ cudaThreadSynchronize();
+ }
+
// Modulate the pressures at the upper boundary cells
if ((darcy.p_mod_A > 1.0e-5 || darcy.p_mod_A < -1.0e-5) &&
darcy.p_mod_f > 1.0e-7) {
t@@ -2377,7 +2401,7 @@ __host__ void DEM::startTime()
velocity_state == 1) {
change_velocity_state = 1.0;
velocity_state = 2;
- } else if (time.current >= 10.0 && velocity_state == 2) {
+ } else if (time.current >= v2_end && velocity_state == 2) {
change_velocity_state = -1.0;
velocity_state = 1;
}
You are viewing proxied material from mx1.adamsgaard.dk. The copyright of proxied material belongs to its original authors. Any comments or complaints in relation to proxied material should be directed to the original authors of the content concerned. Please see the disclaimer for more details.