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 de5a49343d0628cabc393a22841501d2666073af
parent 03437c82b542d8cf41ae59bad59dc5e405ceee52
Author: Anders Damsgaard Christensen <[email protected]>
Date: Tue, 20 Sep 2016 09:57:20 -0700
Merge branch 'master' of github.com:anders-dc/sphere
Diffstat:
M src/darcy.cuh | 36 +++++++++++++++++++++++++++++…
M src/device.cu | 11 +++++++++++
2 files changed, 46 insertions(+), 1 deletion(-)
---
diff --git a/src/darcy.cuh b/src/darcy.cuh
t@@ -761,7 +761,7 @@ __global__ void copyDarcyPorositiesToEdges(
// 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…
+// bottom from the grid interior at the frictionless Y boundaries (grid.period…
// == 2).
__global__ void copyDarcyPorositiesToBottom(
Float* __restrict__ dev_darcy_phi, // in + out
t@@ -794,6 +794,40 @@ __global__ void copyDarcyPorositiesToBottom(
}
+// Copy the porosity, porosity change, div_v_p and vp_avg values to the grid t…
+// from the grid interior if the grid is adaptive (grid.adaptive == 1).
+__global__ void copyDarcyPorositiesToTop(
+ 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.adaptive == 1 &&
+ x < nx && y < ny && z == nz - 1) {
+
+ const unsigned int readidx = d_idx(x, y, nz - 2);
+ 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@@ -2037,6 +2037,17 @@ __host__ void DEM::startTime()
cudaThreadSynchronize();
}
+ // copy porosities to the frictionless lower Z boundary
+ if (grid.adaptive == 1) {
+ copyDarcyPorositiesToTop<<<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) {
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.