ewCudaKernels.cu 5.84 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24
/*
 * EasyWave - A realtime tsunami simulation program with GPU support.
 * Copyright (C) 2014  Andrey Babeyko, Johannes Spazier
 * GFZ German Research Centre for Geosciences (http://www.gfz-potsdam.de)
 *
 * Parts of this program (especially the GPU extension) were developed
 * within the context of the following publicly funded project:
 * - TRIDEC, EU 7th Framework Programme, Grant Agreement 258723
 *   (http://www.tridec-online.eu)
 *
 * This program is free software: you can redistribute it and/or modify
 * it under the terms of the GNU Affero General Public License as
 * published by the Free Software Foundation, either version 3 of the
 * License, or (at your option) any later version.
 *
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU Affero General Public License for more details.
 *
 * You should have received a copy of the GNU Affero General Public License
 * along with this program.  If not, see <http://www.gnu.org/licenses/>.
 */

25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93
#include "ewGpuNode.cuh"
#include "ewCudaKernels.cuh"

__global__ void runWaveUpdateKernel( KernelData data ) {

  Params& dp = data.params;

  int i = blockIdx.y * blockDim.y + threadIdx.y + dp.iMin;
  int j = blockIdx.x * blockDim.x + threadIdx.x + dp.jMin;
  int ij = data.idx(i,j);
  float absH;

  /* maybe unnecessary if controlled from outside */
  if( i <= dp.iMax && j <= dp.jMax && data.d[ij] != 0  ) {

	  float hh = data.h[ij] - data.cR1[ij] * ( data.fM[ij] - data.fM[data.le(ij)] + data.fN[ij] * data.cR6[j] - data.fN[data.dn(ij)]*data.cR6[j-1] );

	  absH = fabs(hh);

	  if( absH < dp.sshZeroThreshold ) {
		  hh = 0.f;
	  } else if( hh > data.hMax[ij] ) {
		  data.hMax[ij] = hh;
		  //hMax[ij] = fmaxf(hMax[ij],h[ij]);
	  }

	  if( dp.sshArrivalThreshold && data.tArr[ij] < 0 && absH > dp.sshArrivalThreshold )
		  data.tArr[ij] = dp.mTime;

	  data.h[ij] = hh;

  }

}

__global__ void runFluxUpdateKernel( KernelData data ) {

	Params& dp = data.params;

	int i = blockIdx.y * blockDim.y + threadIdx.y + dp.iMin;
	int j = blockIdx.x * blockDim.x + threadIdx.x + dp.jMin;
	int ij = data.idx(i,j);

	if( i <= dp.iMax && j <= dp.jMax && data.d[ij] != 0  ) {

	  float hh = data.h[ij];

	  if( data.d[data.ri(ij)] != 0 )
		  data.fM[ij] = data.fM[ij] - data.cR2[ij]*(data.h[data.ri(ij)] - hh);

	  if( data.d[data.up(ij)] != 0 )
		  data.fN[ij] = data.fN[ij] - data.cR4[ij]*(data.h[data.up(ij)] - hh);

	}

}

__global__ void runWaveBoundaryKernel( KernelData data ) {

	KernelData& dt = data;
	Params& dp = data.params;

	int id = blockIdx.x * blockDim.x + threadIdx.x + 2;
	int ij;

	if( dt.devID == 0 ) {

		if( id <= dp.nJ-1 ) {
		  ij = dt.idx(1,id);
94
		  dt.h[ij] = sqrtf( powf(dt.fM[ij],2.0f) + 0.25f*powf((dt.fN[ij] + dt.fN[dt.dn(ij)]),2.0f) )*dt.cB2[id-1];
95 96 97 98 99
		  if( dt.fM[ij] > 0 ) dt.h[ij] = -dt.h[ij];
		}

		if( id == 2 ) {
		  ij = dt.idx(1,1);
100
		  dt.h[ij] = sqrtf( powf(dt.fM[ij],2.0f) + powf(dt.fN[ij],2.0f) )*dt.cB1[0];
101 102 103
		  if( dt.fN[ij] > 0 ) dt.h[ij] = -dt.h[ij];

		  ij = dt.idx(1,dp.nJ);
104
		  dt.h[ij] = sqrtf( powf(dt.fM[ij],2.0f) + powf(dt.fN[dt.dn(ij)],2.0f) )*dt.cB3[0];
105 106 107 108 109 110 111
		  if( dt.fN[dt.dn(ij)] < 0 ) dt.h[ij] = -dt.h[ij];
		}

	} else if( dt.devID == dt.devNum - 1 ) {

		if( id <= dp.nJ-1 ) {
		  ij = dt.idx(dp.nI,id);
112
		  dt.h[ij] = sqrtf( powf(dt.fM[dt.le(ij)],2.0f) + 0.25f*powf((dt.fN[ij] + dt.fN[dt.dn(ij)]),2.0f) )*dt.cB4[id-1];
113 114 115 116 117
		  if( dt.fM[dt.le(ij)] < 0 ) dt.h[ij] = -dt.h[ij];
		}

		if( id == 2 ) {
		  ij = dt.idx(dp.nI,1);
118
		  dt.h[ij] = sqrtf( powf(dt.fM[dt.le(ij)],2.0f) + powf(dt.fN[ij],2.0f) )*dt.cB1[dp.nI-1];
119 120 121
		  if( dt.fN[ij] > 0 ) dt.h[ij] = -dt.h[ij];

		  ij = dt.idx(dp.nI,dp.nJ);
122
		  dt.h[ij] = sqrtf( powf(dt.fM[dt.le(ij)],2.0f) + powf(dt.fN[dt.dn(ij)],2.0f) )*dt.cB3[dp.nI-1];
123 124 125 126 127 128 129
		  if( dt.fN[dt.dn(ij)] < 0 ) dt.h[ij] = -dt.h[ij];
		}

	}

	if( id <= dp.nI - 1 ) {
	  ij = dt.idx(id,1);
130
	  dt.h[ij] = sqrtf( powf(dt.fN[ij],2.0f) + 0.25f*powf((dt.fM[ij] + dt.fM[dt.le(ij)]),2.0f) )*dt.cB1[id-1];
131 132 133 134 135
	  if( dt.fN[ij] > 0 ) dt.h[ij] = -dt.h[ij];
	}

	if( id <= dp.nI - 1 ) {
	  ij = dt.idx(id,dp.nJ);
136
	  dt.h[ij] = sqrtf( powf(dt.fN[dt.dn(ij)],2.0f) + 0.25f*powf((dt.fM[ij] + dt.fM[dt.dn(ij)]),2.0f) )*dt.cB3[id-1];
137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212
	  if( dt.fN[dt.dn(ij)] < 0 ) dt.h[ij] = -dt.h[ij];
	}

}

__global__ void runFluxBoundaryKernel( KernelData data ) {

	KernelData& dt = data;
	Params& dp = data.params;

	int id = blockIdx.x * blockDim.x + threadIdx.x + 1;
	int ij;

	if( dt.devID == 0 ) {

		if( id <= dp.nJ ) {
		  ij = dt.idx(1,id);
		  dt.fM[ij] = dt.fM[ij] - dt.cR2[ij]*(dt.h[dt.ri(ij)] - dt.h[ij]);
		}

		if( id <= dp.nJ-1 ) {
		  ij = dt.idx(1,id);
		  dt.fN[ij] = dt.fN[ij] - dt.cR4[ij]*(dt.h[dt.up(ij)] - dt.h[ij]);
		}

	} else if( dt.devID == dt.devNum - 1 ) {

		if( id <= dp.nJ-1 ) {
		  ij = dt.idx(dp.nI,id);
		  dt.fN[ij] = dt.fN[ij] - dt.cR4[ij]*(dt.h[dt.up(ij)] - dt.h[ij]);
		}

	}

	if( id <= dp.nI - 1 ) {
	  ij = dt.idx(id,1);
	  dt.fM[ij] = dt.fM[ij] - dt.cR2[ij]*(dt.h[dt.ri(ij)] - dt.h[ij]);
	}

	if( id <= dp.nI - 1 ) {
	  ij = dt.idx(id,dp.nJ);
	  dt.fM[ij] = dt.fM[ij] - dt.cR2[ij]*(dt.h[dt.ri(ij)] - dt.h[ij]);
	}

	if( id <= dp.nI ) {
	  ij = dt.idx(id,1);
	  dt.fN[ij] = dt.fN[ij] - dt.cR4[ij]*(dt.h[dt.up(ij)] - dt.h[ij]);
	}

}

__global__ void runGridExtendKernel( KernelData data ) {

	Params& dp = data.params;

	int id = blockIdx.x * blockDim.x + threadIdx.x;

	if( id >= dp.jMin && id <= dp.jMax ) {

		if( fabsf(data.h[data.idx(dp.iMin+2,id)]) > dp.sshClipThreshold )
			atomicAdd( &(data.extend->x), 1 );

		if( fabsf(data.h[data.idx(dp.iMax-2,id)]) > dp.sshClipThreshold )
			atomicAdd( &(data.extend->y), 1 );
	}

	if( id >= dp.iMin && id <= dp.iMax ) {

	  if( fabsf(data.h[data.idx(id,dp.jMin+2)]) > dp.sshClipThreshold )
		  atomicAdd( &(data.extend->z), 1 );

	  if( fabsf(data.h[data.idx(id,dp.jMax-2)]) > dp.sshClipThreshold )
		  atomicAdd( &(data.extend->w), 1 );
	}

}