Newer
Older
#include "sfx_esm.h"
#include "sfx_model_compute_subfunc.cuh"
#include "sfx_surface.cuh"
#include "sfx_memory_processing.cuh"
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
namespace sfx_kernel
{
template<typename T>
__global__ void compute_flux(sfxDataVecTypeC sfx,
meteoDataVecTypeC meteo,
const sfx_esm_param model_param,
const sfx_surface_param surface_param,
const sfx_esm_numericsTypeC numerics,
const sfx_phys_constants phys_constants,
const int grid_size);
}
template<typename T>
__global__ void sfx_kernel::compute_flux(sfxDataVecTypeC sfx,
meteoDataVecTypeC meteo,
const sfx_esm_param model_param,
const sfx_surface_param surface_param,
const sfx_esm_numericsTypeC numerics,
const sfx_phys_constants phys_constants,
const int grid_size)
{
const int index = blockIdx.x * blockDim.x + threadIdx.x;
T h, U, dT, Tsemi, dQ, z0_m;
T Re, z0_t, B, h0_m, h0_t, u_dyn0, zeta, Rib, zeta_conv_lim, Rib_conv_lim, f_m_conv_lim, f_h_conv_lim, psi_m, psi_h, phi_m, phi_h, Km, Pr_t_inv, Cm, Ct;
int surface_type;
T fval;
if(index < grid_size)
{
U = meteo.U[index];
Tsemi = meteo.Tsemi[index];
dT = meteo.dT[index];
dQ = meteo.dQ[index];
h = meteo.h[index];
z0_m = meteo.z0_m[index];
surface_type = z0_m < 0.0 ? surface_param.surface_ocean : surface_param.surface_land;
if (surface_type == surface_param.surface_ocean)
{
get_charnock_roughness(z0_m, u_dyn0, U, h, surface_param, numerics.maxiters_charnock);
h0_m = h / z0_m;
}
if (surface_type == surface_param.surface_land)
{
h0_m = h / z0_m;
u_dyn0 = U * model_param.kappa / log(h0_m);
}
Re = u_dyn0 * z0_m / phys_constants.nu_air;
get_thermal_roughness(z0_t, B, z0_m, Re, surface_param, surface_type);
h0_t = h / z0_t;
Rib = (phys_constants.g / Tsemi) * h * (dT + 0.61e0 * Tsemi * dQ) / (U*U);
get_convection_lim(zeta_conv_lim, Rib_conv_lim, f_m_conv_lim, f_h_conv_lim,
h0_m, h0_t, B,
model_param);
if (Rib > 0.0)
{
Rib = sfx_math::min(Rib, model_param.Rib_max);
get_psi_stable(psi_m, psi_h, zeta, Rib, h0_m, h0_t, B, model_param);
fval = model_param.beta_m * zeta;
phi_m = 1.0 + fval;
phi_h = 1.0/model_param.Pr_t_0_inv + fval;
}
else if (Rib < Rib_conv_lim)
{
get_psi_convection(psi_m, psi_h, zeta, Rib, h0_m, h0_t, B, zeta_conv_lim, f_m_conv_lim, f_h_conv_lim, model_param, numerics.maxiters_convection);
fval = pow(zeta_conv_lim / zeta, 1.0/3.0);
phi_m = fval / f_m_conv_lim;
phi_h = fval / (model_param.Pr_t_0_inv * f_h_conv_lim);
}
else if (Rib > -0.001)
{
get_psi_neutral(psi_m, psi_h, zeta, h0_m, h0_t, B, model_param);
phi_m = 1.0;
phi_h = 1.0 / model_param.Pr_t_0_inv;
}
else
{
get_psi_semi_convection(psi_m, psi_h, zeta, Rib, h0_m, h0_t, B, model_param, numerics.maxiters_convection);
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
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
phi_m = pow(1.0 - model_param.alpha_m * zeta, -0.25);
phi_h = 1.0 / (model_param.Pr_t_0_inv * sqrt(1.0 - model_param.alpha_h_fix * zeta));
}
Cm = model_param.kappa / psi_m;
Ct = model_param.kappa / psi_h;
Km = model_param.kappa * Cm * U * h / phi_m;
Pr_t_inv = phi_m / phi_h;
sfx.zeta[index] = zeta;
sfx.Rib[index] = Rib;
sfx.Re[index] = Re;
sfx.B[index] = B;
sfx.z0_m[index] = z0_m;
sfx.z0_t[index] = z0_t;
sfx.Rib_conv_lim[index] = Rib_conv_lim;
sfx.Cm[index] = Cm;
sfx.Ct[index] = Ct;
sfx.Km[index] = Km;
sfx.Pr_t_inv[index] = Pr_t_inv;
}
}
template<typename T, MemType memIn, MemType memOut >
void FluxEsm<T, memIn, memOut, MemType::GPU>::compute_flux()
{
const int BlockCount = int(ceil(float(grid_size) / 1024.0));
dim3 cuBlock = dim3(1024, 1, 1);
dim3 cuGrid = dim3(BlockCount, 1, 1);
sfx_kernel::compute_flux<T><<<cuGrid, cuBlock>>>(sfx, meteo, model_param,
surface_param, numerics, phys_constants, grid_size);
if(MemType::GPU != memOut)
{
const size_t new_size = grid_size * sizeof(T);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->zeta, (void*&)sfx.zeta, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->Rib, (void*&)sfx.Rib, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->Re, (void*&)sfx.Re, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->B, (void*&)sfx.B, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->z0_m, (void*&)sfx.z0_m, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->z0_t, (void*&)sfx.z0_t, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->Rib_conv_lim, (void*&)sfx.Rib_conv_lim, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->Cm, (void*&)sfx.Cm, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->Ct, (void*&)sfx.Ct, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->Km, (void*&)sfx.Km, new_size);
memproc::memcopy<memOut, MemType::GPU>((void*&)res_sfx->Pr_t_inv, (void*&)sfx.Pr_t_inv, new_size);
}
}
template class FluxEsmBase<float, MemType::GPU, MemType::GPU, MemType::GPU>;
template class FluxEsmBase<float, MemType::GPU, MemType::GPU, MemType::CPU>;
template class FluxEsmBase<float, MemType::GPU, MemType::CPU, MemType::GPU>;
template class FluxEsmBase<float, MemType::CPU, MemType::GPU, MemType::GPU>;
template class FluxEsmBase<float, MemType::CPU, MemType::CPU, MemType::GPU>;
template class FluxEsmBase<float, MemType::CPU, MemType::GPU, MemType::CPU>;
template class FluxEsmBase<float, MemType::GPU, MemType::CPU, MemType::CPU>;
template class FluxEsm<float, MemType::GPU, MemType::GPU, MemType::GPU>;
template class FluxEsm<float, MemType::GPU, MemType::GPU, MemType::CPU>;
template class FluxEsm<float, MemType::GPU, MemType::CPU, MemType::GPU>;
template class FluxEsm<float, MemType::CPU, MemType::GPU, MemType::GPU>;
template class FluxEsm<float, MemType::CPU, MemType::CPU, MemType::GPU>;
template class FluxEsm<float, MemType::CPU, MemType::GPU, MemType::CPU>;
template class FluxEsm<float, MemType::GPU, MemType::CPU, MemType::CPU>;