ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ERF_TerrainMetrics.H
Go to the documentation of this file.
1 #ifndef ERF_TERRAIN_METRIC_H_
2 #define ERF_TERRAIN_METRIC_H_
3 
4 #include <AMReX.H>
5 #include <AMReX_Geometry.H>
6 #include <AMReX_MultiFab.H>
7 #include <ERF_IndexDefines.H>
8 
9 /**
10  * Routine to define default z_phys_nd and z_phys_cc
11  */
12 void
13 init_default_zphys (int lev, const amrex::Geometry& geom,
14  amrex::MultiFab& z_phys_nd, amrex::MultiFab& z_phys_cc);
15 
16 /**
17  * Utility routines for constructing terrain metric terms
18  */
19 
20 // Declare functions for ERF.cpp
21 void init_zlevels (amrex::Vector<amrex::Vector<amrex::Real>> & zlevels_stag,
22  amrex::Vector<amrex::Vector<amrex::Real>> & stretched_dz_h,
23  amrex::Vector<amrex::Gpu::DeviceVector<amrex::Real>> & stretched_dz_d,
24  amrex::Vector<amrex::Geometry> const& geom,
25  amrex::Vector<amrex::IntVect> const& ref_ratio,
26  const amrex::Real grid_stretching_ratio,
27  const amrex::Real zsurf,
28  const amrex::Real dz0);
29 
30 void make_terrain_fitted_coords (int lev, const amrex::Geometry& geom,
31  amrex::MultiFab& z_phys_nd,
32  amrex::Vector<amrex::Real> const& z_levels_h,
33  amrex::GpuArray<ERF_BC, AMREX_SPACEDIM*2>& phys_bc_type);
34 
35 void init_which_terrain_grid (int lev, const amrex::Geometry& geom,
36  amrex::MultiFab& z_phys_nd,
37  amrex::Vector<amrex::Real> const& z_levels_h);
38 
39 
40 //*****************************************************************************************
41 // Compute terrain metric terms at cell-center
42 //*****************************************************************************************
43 // Metric is at cell center
44 AMREX_FORCE_INLINE
45 AMREX_GPU_DEVICE
46 amrex::Real
47 Compute_h_zeta_AtCellCenter (const int &i, const int &j, const int &k,
48  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
49  const amrex::Array4<const amrex::Real>& z_nd)
50 {
51  amrex::Real dzInv = cellSizeInv[2];
52  amrex::Real met_h_zeta = 0.25 * dzInv *
53  ( z_nd(i+1,j,k+1) + z_nd(i+1,j+1,k+1) + z_nd(i,j,k+1) + z_nd(i,j+1,k+1)
54  -z_nd(i+1,j,k ) - z_nd(i+1,j+1,k ) - z_nd(i,j,k ) - z_nd(i,j+1,k ) );
55  return met_h_zeta;
56 }
57 
58 // Metric is at cell center
59 AMREX_GPU_DEVICE
60 AMREX_FORCE_INLINE
61 amrex::Real
62 Compute_h_xi_AtCellCenter (const int &i, const int &j, const int &k,
63  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
64  const amrex::Array4<const amrex::Real>& z_nd)
65 {
66  amrex::Real dxInv = cellSizeInv[0];
67  amrex::Real met_h_xi = 0.25 * dxInv *
68  ( z_nd(i+1,j,k) + z_nd(i+1,j+1,k) + z_nd(i+1,j,k+1) + z_nd(i+1,j+1,k+1)
69  -z_nd(i ,j,k) - z_nd(i ,j+1,k) - z_nd(i ,j,k+1) - z_nd(i ,j+1,k+1) );
70  return met_h_xi;
71 }
72 
73 // Metric is at cell center
74 AMREX_GPU_DEVICE
75 AMREX_FORCE_INLINE
76 amrex::Real
77 Compute_h_eta_AtCellCenter (const int &i, const int &j, const int &k,
78  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
79  const amrex::Array4<const amrex::Real>& z_nd)
80 {
81  amrex::Real dyInv = cellSizeInv[1];
82  amrex::Real met_h_eta = 0.25 * dyInv *
83  ( z_nd(i,j+1,k) + z_nd(i+1,j+1,k) + z_nd(i,j+1,k+1) + z_nd(i+1,j+1,k+1)
84  -z_nd(i,j ,k) - z_nd(i+1,j ,k) - z_nd(i,j ,k+1) - z_nd(i+1,j ,k+1) );
85  return met_h_eta;
86 }
87 
88 
89 //*****************************************************************************************
90 // Compute terrain metric terms at face-centers
91 //*****************************************************************************************
92 // Metric coincides with U location
93 AMREX_GPU_DEVICE
94 AMREX_FORCE_INLINE
95 amrex::Real
96 Compute_h_zeta_AtIface (const int &i, const int &j, const int &k,
97  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
98  const amrex::Array4<const amrex::Real>& z_nd)
99 {
100  amrex::Real dzInv = cellSizeInv[2];
101  amrex::Real met_h_zeta = 0.5 * dzInv * ( z_nd(i,j,k+1) + z_nd(i,j+1,k+1)
102  - z_nd(i,j,k ) - z_nd(i,j+1,k ) );
103  return met_h_zeta;
104 }
105 
106 // Metric coincides with U location
107 AMREX_GPU_DEVICE
108 AMREX_FORCE_INLINE
109 amrex::Real
110 Compute_h_xi_AtIface (const int &i, const int &j, const int &k,
111  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
112  const amrex::Array4<const amrex::Real>& z_nd)
113 {
114  amrex::Real dxInv = cellSizeInv[0];
115  amrex::Real met_h_xi = 0.125 * dxInv *
116  ( z_nd(i+1,j,k) + z_nd(i+1,j,k+1) + z_nd(i+1,j+1,k) + z_nd(i+1,j+1,k+1)
117  -z_nd(i-1,j,k) - z_nd(i-1,j,k+1) - z_nd(i-1,j+1,k) - z_nd(i-1,j+1,k+1) );
118  return met_h_xi;
119 }
120 
121 // Metric coincides with U location
122 AMREX_GPU_DEVICE
123 AMREX_FORCE_INLINE
124 amrex::Real
125 Compute_h_eta_AtIface (const int &i, const int &j, const int &k,
126  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
127  const amrex::Array4<const amrex::Real>& z_nd)
128 {
129  amrex::Real dyInv = cellSizeInv[1];
130  amrex::Real met_h_eta = 0.5 * dyInv * ( z_nd(i,j+1,k ) + z_nd(i,j+1,k+1)
131  - z_nd(i,j ,k ) - z_nd(i,j ,k+1) );
132  return met_h_eta;
133 }
134 
135 // Metric coincides with V location
136 AMREX_GPU_DEVICE
137 AMREX_FORCE_INLINE
138 amrex::Real
139 Compute_h_zeta_AtJface (const int &i, const int &j, const int &k,
140  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
141  const amrex::Array4<const amrex::Real>& z_nd)
142 {
143  amrex::Real dzInv = cellSizeInv[2];
144  amrex::Real met_h_zeta = 0.5 * dzInv * ( z_nd(i,j,k+1) + z_nd(i+1,j,k+1)
145  - z_nd(i,j,k ) - z_nd(i+1,j,k ) );
146  return met_h_zeta;
147 }
148 
149 // Metric coincides with V location
150 AMREX_GPU_DEVICE
151 AMREX_FORCE_INLINE
152 amrex::Real
153 Compute_h_xi_AtJface (const int &i, const int &j, const int &k,
154  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
155  const amrex::Array4<const amrex::Real>& z_nd)
156 {
157  amrex::Real dxInv = cellSizeInv[0];
158  amrex::Real met_h_xi = 0.5 * dxInv * ( z_nd(i+1,j,k) + z_nd(i+1,j,k+1)
159  -z_nd(i ,j,k) - z_nd(i ,j,k+1) );
160  return met_h_xi;
161 }
162 
163 // Metric coincides with V location
164 AMREX_GPU_DEVICE
165 AMREX_FORCE_INLINE
166 amrex::Real
167 Compute_h_eta_AtJface (const int &i, const int &j, const int &k,
168  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
169  const amrex::Array4<const amrex::Real>& z_nd)
170 {
171  amrex::Real dyInv = cellSizeInv[1];
172  amrex::Real met_h_eta = 0.125 * dyInv *
173  ( z_nd(i,j+1,k) + z_nd(i,j+1,k+1) + z_nd(i+1,j+1,k) + z_nd(i+1,j+1,k+1)
174  -z_nd(i,j-1,k) - z_nd(i,j-1,k+1) - z_nd(i+1,j-1,k) - z_nd(i+1,j-1,k+1) );
175  return met_h_eta;
176 }
177 
178 // Metric coincides with K location
179 AMREX_GPU_DEVICE
180 AMREX_FORCE_INLINE
181 amrex::Real
182 Compute_h_zeta_AtKface (const int &i, const int &j, const int &k,
183  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
184  const amrex::Array4<const amrex::Real>& z_nd)
185 {
186  amrex::Real dzInv = cellSizeInv[2];
187  amrex::Real met_h_zeta = 0.125 * dzInv *
188  ( z_nd(i,j,k+1) + z_nd(i+1,j,k+1) + z_nd(i,j+1,k+1) + z_nd(i+1,j+1,k+1)
189  -z_nd(i,j,k-1) - z_nd(i+1,j,k-1) - z_nd(i,j+1,k-1) - z_nd(i+1,j+1,k-1) );
190  return met_h_zeta;
191 }
192 
193 // Metric coincides with K location
194 AMREX_GPU_DEVICE
195 AMREX_FORCE_INLINE
196 amrex::Real
197 Compute_h_xi_AtKface (const int &i, const int &j, const int &k,
198  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
199  const amrex::Array4<const amrex::Real>& z_nd)
200 {
201  amrex::Real dxInv = cellSizeInv[0];
202  amrex::Real met_h_xi = 0.5 * dxInv * ( z_nd(i+1,j,k) + z_nd(i+1,j+1,k)
203  -z_nd(i ,j,k) - z_nd(i ,j+1,k) );
204  return met_h_xi;
205 }
206 
207 // Metric coincides with K location
208 AMREX_GPU_DEVICE
209 AMREX_FORCE_INLINE
210 amrex::Real
211 Compute_h_eta_AtKface (const int &i, const int &j, const int &k,
212  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
213  const amrex::Array4<const amrex::Real>& z_nd)
214 {
215  amrex::Real dyInv = cellSizeInv[1];
216  amrex::Real met_h_eta = 0.5 * dyInv * ( z_nd(i,j+1,k) + z_nd(i+1,j+1,k)
217  -z_nd(i,j ,k) - z_nd(i+1,j ,k) );
218  return met_h_eta;
219 }
220 
221 
222 //*****************************************************************************************
223 // Compute terrain metric terms at edge-centers
224 //*****************************************************************************************
225 // -- EdgeCenterK --
226 
227 // Metric is at edge and center Z (red pentagon)
228 AMREX_GPU_DEVICE
229 AMREX_FORCE_INLINE
230 amrex::Real
231 Compute_h_zeta_AtEdgeCenterK (const int &i, const int &j, const int &k,
232  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
233  const amrex::Array4<const amrex::Real>& z_nd)
234 {
235  amrex::Real dzInv = cellSizeInv[2];
236  amrex::Real met_h_zeta = dzInv * (z_nd(i,j,k+1) - z_nd(i,j,k));
237  return met_h_zeta;
238 }
239 
240 // Metric is at edge and center Z (red pentagon)
241 AMREX_GPU_DEVICE
242 AMREX_FORCE_INLINE
243 amrex::Real
244 Compute_h_xi_AtEdgeCenterK (const int &i, const int &j, const int &k,
245  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
246  const amrex::Array4<const amrex::Real>& z_nd)
247 {
248  amrex::Real dxInv = cellSizeInv[0];
249  amrex::Real met_h_xi = 0.25 * dxInv *
250  ( z_nd(i+1,j,k) + z_nd(i+1,j,k+1)
251  -z_nd(i-1,j,k) - z_nd(i-1,j,k+1) );
252  return met_h_xi;
253 }
254 
255 // Metric is at edge and center Z (red pentagon)
256 AMREX_GPU_DEVICE
257 AMREX_FORCE_INLINE
258 amrex::Real
259 Compute_h_eta_AtEdgeCenterK (const int &i, const int &j, const int &k,
260  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
261  const amrex::Array4<const amrex::Real>& z_nd)
262 {
263  amrex::Real dyInv = cellSizeInv[1];
264  amrex::Real met_h_eta = 0.25 * dyInv *
265  ( z_nd(i,j+1,k) + z_nd(i,j+1,k+1)
266  -z_nd(i,j-1,k) - z_nd(i,j-1,k+1) );
267  return met_h_eta;
268 }
269 
270 // -- EdgeCenterJ --
271 
272 // Metric is at edge and center Y (magenta cross)
273 AMREX_GPU_DEVICE
274 AMREX_FORCE_INLINE
275 amrex::Real
276 Compute_h_zeta_AtEdgeCenterJ (const int &i, const int &j, const int &k,
277  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
278  const amrex::Array4<const amrex::Real>& z_nd)
279 {
280  amrex::Real dzInv = cellSizeInv[2];
281  amrex::Real met_h_zeta = 0.25 * dzInv * ( z_nd(i,j,k+1) + z_nd(i,j+1,k+1)
282  -z_nd(i,j,k-1) - z_nd(i,j+1,k-1) );
283  return met_h_zeta;
284 }
285 
286 // Metric is at edge and center Y (magenta cross)
287 AMREX_GPU_DEVICE
288 AMREX_FORCE_INLINE
289 amrex::Real
290 Compute_h_xi_AtEdgeCenterJ (const int &i, const int &j, const int &k,
291  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
292  const amrex::Array4<const amrex::Real>& z_nd)
293 {
294  amrex::Real dxInv = cellSizeInv[0];
295  amrex::Real met_h_xi = 0.25 * dxInv *
296  ( z_nd(i+1,j+1,k) + z_nd(i+1,j ,k)
297  -z_nd(i-1,j+1,k) - z_nd(i-1,j ,k) );
298  return met_h_xi;
299 }
300 
301 // Metric is at edge and center Y (magenta cross)
302 AMREX_GPU_DEVICE
303 AMREX_FORCE_INLINE
304 amrex::Real
305 Compute_h_eta_AtEdgeCenterJ (const int &i, const int &j, const int &k,
306  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
307  const amrex::Array4<const amrex::Real>& z_nd)
308 {
309  amrex::Real dyInv = cellSizeInv[1];
310  amrex::Real met_h_eta = dyInv * ( z_nd(i,j+1,k) - z_nd(i,j,k) );
311  return met_h_eta;
312 }
313 
314 // -- EdgeCenterI --
315 
316 // Metric is at edge and center X (purple hexagon)
317 AMREX_GPU_DEVICE
318 AMREX_FORCE_INLINE
319 amrex::Real
320 Compute_h_zeta_AtEdgeCenterI (const int &i, const int &j, const int &k,
321  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
322  const amrex::Array4<const amrex::Real>& z_nd)
323 {
324  amrex::Real dzInv = cellSizeInv[2];
325  amrex::Real met_h_zeta = 0.25 * dzInv * ( z_nd(i,j,k+1) + z_nd(i+1,j,k+1)
326  -z_nd(i,j,k-1) - z_nd(i+1,j,k-1) );
327  return met_h_zeta;
328 }
329 
330 // Metric is at edge and center X (purple hexagon)
331 AMREX_GPU_DEVICE
332 AMREX_FORCE_INLINE
333 amrex::Real
334 Compute_h_xi_AtEdgeCenterI (const int &i, const int &j, const int &k,
335  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
336  const amrex::Array4<const amrex::Real>& z_nd)
337 {
338  amrex::Real dxInv = cellSizeInv[0];
339  amrex::Real met_h_xi = dxInv * ( z_nd(i+1,j,k) - z_nd(i,j,k) );
340  return met_h_xi;
341 }
342 
343 // Metric is at edge and center X (purple hexagon)
344 AMREX_GPU_DEVICE
345 AMREX_FORCE_INLINE
346 amrex::Real
347 Compute_h_eta_AtEdgeCenterI (const int &i, const int &j, const int &k,
348  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& cellSizeInv,
349  const amrex::Array4<const amrex::Real>& z_nd)
350 {
351  amrex::Real dyInv = cellSizeInv[1];
352  amrex::Real met_h_eta = 0.25 * dyInv *
353  ( z_nd(i+1,j+1,k) + z_nd(i,j+1,k)
354  -z_nd(i+1,j-1,k) - z_nd(i,j-1,k) );
355  return met_h_eta;
356 }
357 
358 // Relative height above terrain surface at cell center from z_nd (nodal absolute height)
359 AMREX_GPU_DEVICE
360 AMREX_FORCE_INLINE
361 amrex::Real
362 Compute_Z_AtCellCenter (const int &i, const int &j, const int &k,
363  const amrex::Array4<const amrex::Real>& z_nd)
364 {
365  const amrex::Real z_cc = 0.125*( z_nd(i ,j ,k ) + z_nd(i ,j ,k+1) +
366  + z_nd(i+1,j ,k ) + z_nd(i ,j ,k+1)
367  + z_nd(i ,j+1,k ) + z_nd(i ,j+1,k+1)
368  + z_nd(i+1,j+1,k ) + z_nd(i ,j+1,k+1));
369 
370  return z_cc;
371 }
372 
373 // Relative height above terrain surface at cell center from z_nd (nodal absolute height)
374 AMREX_GPU_DEVICE
375 AMREX_FORCE_INLINE
376 amrex::Real
377 Compute_Z_AtWFace (const int &i, const int &j, const int &k,
378  const amrex::Array4<const amrex::Real>& z_nd)
379 {
380  const amrex::Real z_wf = 0.25*( z_nd(i ,j ,k ) + z_nd(i+1,j ,k )
381  + z_nd(i ,j+1,k ) + z_nd(i+1,j+1,k ) );
382 
383  return z_wf;
384 }
385 
386 // Relative height above terrain surface at cell center from z_nd (nodal absolute height)
387 AMREX_GPU_DEVICE
388 AMREX_FORCE_INLINE
389 amrex::Real
390 Compute_Zrel_AtCellCenter (const int &i, const int &j, const int &k,
391  const amrex::Array4<const amrex::Real>& z_nd)
392 {
393  const amrex::Real z_cc = 0.125*( z_nd(i ,j ,k ) + z_nd(i ,j ,k+1) +
394  + z_nd(i+1,j ,k ) + z_nd(i ,j ,k+1)
395  + z_nd(i ,j+1,k ) + z_nd(i ,j+1,k+1)
396  + z_nd(i+1,j+1,k ) + z_nd(i ,j+1,k+1));
397 
398  // Note: we assume the z_nd array spans from the bottom to top of the domain
399  // i.e. no domain decomposition across processors in vertical direction
400  const amrex::Real z0_cc = 0.25*( z_nd(i ,j ,0) + z_nd(i ,j+1,0)
401  + z_nd(i+1,j ,0) + z_nd(i+1,j+1,0));
402 
403  return (z_cc - z0_cc);
404 }
405 
406 
407 //*****************************************************************************************
408 // Contravariant velocity computation
409 //*****************************************************************************************
410 
411 // Define omega given u,v and w (4th order)
412 AMREX_GPU_DEVICE
413 AMREX_FORCE_INLINE
414 amrex::Real
415 OmegaFromW (int& i, int& j, int& k, amrex::Real w,
416  const amrex::Array4<const amrex::Real>& u_arr,
417  const amrex::Array4<const amrex::Real>& v_arr,
418  const amrex::Array4<const amrex::Real>& mf_u,
419  const amrex::Array4<const amrex::Real>& mf_v,
420  const amrex::Array4<const amrex::Real>& z_nd,
421  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& dxInv)
422 {
423  // This is dh/dxi at hi and lo edges
424  amrex::Real z_x_p2 = Compute_Z_AtWFace(i+2,j,k,z_nd);
425  amrex::Real z_x_p1 = Compute_Z_AtWFace(i+1,j,k,z_nd);
426  amrex::Real z_x_m1 = Compute_Z_AtWFace(i-1,j,k,z_nd);
427  amrex::Real z_x_m2 = Compute_Z_AtWFace(i-2,j,k,z_nd);
428  amrex::Real met_xi = (amrex::Real(1.0/12.0)*(z_x_m2 - z_x_p2)
429  + amrex::Real(8.0/12.0)*(z_x_p1 - z_x_m1)) * dxInv[0];
430 
431  // This is dh/deta at hi and lo edges
432  amrex::Real z_y_p2 = Compute_Z_AtWFace(i,j+2,k,z_nd);
433  amrex::Real z_y_p1 = Compute_Z_AtWFace(i,j+1,k,z_nd);
434  amrex::Real z_y_m1 = Compute_Z_AtWFace(i,j-1,k,z_nd);
435  amrex::Real z_y_m2 = Compute_Z_AtWFace(i,j-2,k,z_nd);
436  amrex::Real met_eta = (amrex::Real(1.0/12.0)*(z_y_m2 - z_y_p2)
437  + amrex::Real(8.0/12.0)*(z_y_p1 - z_y_m1)) * dxInv[1];
438 
439  // Use extrapolation instead of interpolation if at the bottom boundary
440  amrex::Real u_hi = (k == 0) ? 1.5 * u_arr(i+1,j ,k ) - 0.5 * u_arr(i+1,j ,k+1) :
441  0.5 * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
442  amrex::Real u_lo = (k == 0) ? 1.5 * u_arr(i ,j ,k ) - 0.5 * u_arr(i ,j ,k+1) :
443  0.5 * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
444  amrex::Real mf_u_hi = mf_u(i+1,j,0);
445  amrex::Real mf_u_lo = mf_u(i ,j,0);
446 
447  amrex::Real v_hi = (k == 0) ? 1.5 * v_arr(i ,j+1,k ) - 0.5 * v_arr(i ,j+1,k+1) :
448  0.5 * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
449  amrex::Real v_lo = (k == 0) ? 1.5 * v_arr(i ,j ,k ) - 0.5 * v_arr(i ,j ,k+1) :
450  0.5 * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
451  amrex::Real mf_v_hi = mf_v(i,j+1,0);
452  amrex::Real mf_v_lo = mf_v(i,j ,0);
453 
454  amrex::Real u_met = met_xi * 0.5 * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
455  amrex::Real v_met = met_eta * 0.5 * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
456 
457  amrex::Real omega = w - u_met - v_met;
458  return omega;
459 }
460 
461 // Define w given u and v arrays and scalar omega (4th order)
462 AMREX_GPU_DEVICE
463 AMREX_FORCE_INLINE
464 amrex::Real
465 WFromOmega (int& i, int& j, int& k, amrex::Real omega,
466  const amrex::Array4<const amrex::Real>& u_arr,
467  const amrex::Array4<const amrex::Real>& v_arr,
468  const amrex::Array4<const amrex::Real>& mf_u,
469  const amrex::Array4<const amrex::Real>& mf_v,
470  const amrex::Array4<const amrex::Real>& z_nd,
471  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& dxInv)
472 {
473  // This is dh/dxi at hi and lo edges
474  amrex::Real z_x_p2 = Compute_Z_AtWFace(i+2,j,k,z_nd);
475  amrex::Real z_x_p1 = Compute_Z_AtWFace(i+1,j,k,z_nd);
476  amrex::Real z_x_m1 = Compute_Z_AtWFace(i-1,j,k,z_nd);
477  amrex::Real z_x_m2 = Compute_Z_AtWFace(i-2,j,k,z_nd);
478  amrex::Real met_xi = (amrex::Real(1.0/12.0)*(z_x_m2 - z_x_p2)
479  + amrex::Real(8.0/12.0)*(z_x_p1 - z_x_m1)) * dxInv[0];
480 
481  // This is dh/deta at hi and lo edges
482  amrex::Real z_y_p2 = Compute_Z_AtWFace(i,j+2,k,z_nd);
483  amrex::Real z_y_p1 = Compute_Z_AtWFace(i,j+1,k,z_nd);
484  amrex::Real z_y_m1 = Compute_Z_AtWFace(i,j-1,k,z_nd);
485  amrex::Real z_y_m2 = Compute_Z_AtWFace(i,j-2,k,z_nd);
486  amrex::Real met_eta = (amrex::Real(1.0/12.0)*(z_y_m2 - z_y_p2)
487  + amrex::Real(8.0/12.0)*(z_y_p1 - z_y_m1)) * dxInv[1];
488 
489  // Use extrapolation instead of interpolation if at the bottom boundary
490  amrex::Real u_hi = (k == 0) ? 1.5 * u_arr(i+1,j ,k ) - 0.5 * u_arr(i+1,j ,k+1) :
491  0.5 * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
492  amrex::Real u_lo = (k == 0) ? 1.5 * u_arr(i ,j ,k ) - 0.5 * u_arr(i ,j ,k+1) :
493  0.5 * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
494  amrex::Real mf_u_hi = mf_u(i+1,j,0);
495  amrex::Real mf_u_lo = mf_u(i ,j,0);
496 
497  amrex::Real v_hi = (k == 0) ? 1.5 * v_arr(i ,j+1,k ) - 0.5 * v_arr(i ,j+1,k+1) :
498  0.5 * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
499  amrex::Real v_lo = (k == 0) ? 1.5 * v_arr(i ,j ,k ) - 0.5 * v_arr(i ,j ,k+1) :
500  0.5 * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
501  amrex::Real mf_v_hi = mf_v(i,j+1,0);
502  amrex::Real mf_v_lo = mf_v(i,j ,0);
503 
504  amrex::Real u_met = met_xi * 0.5 * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
505  amrex::Real v_met = met_eta * 0.5 * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
506 
507  amrex::Real w = omega + u_met + v_met;
508  return w;
509 }
510 
511 //
512 // NOTE: 2nd order kept here for debugging purposes
513 //
514 
515 # if 0
516 // Define omega given u,v and w (2nd order)
517 AMREX_GPU_DEVICE
518 AMREX_FORCE_INLINE
519 amrex::Real
520 OmegaFromW (int& i, int& j, int& k, amrex::Real w,
521  const amrex::Array4<const amrex::Real>& u_arr,
522  const amrex::Array4<const amrex::Real>& v_arr,
523  const amrex::Array4<const amrex::Real>& mf_u,
524  const amrex::Array4<const amrex::Real>& mf_v,
525  const amrex::Array4<const amrex::Real>& z_nd,
526  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& dxInv)
527 {
528  // This is dh/dxi at hi and lo edges
529  amrex::Real z_x_h = Compute_Z_AtWFace(i+1,j,k,z_nd);
530  amrex::Real z_x_l = Compute_Z_AtWFace(i-1,j,k,z_nd);
531  amrex::Real met_xi = 0.5 * (z_x_h - z_x_l) * dxInv[0];
532 
533  // This is dh/deta at hi and lo edges
534  amrex::Real z_y_h = Compute_Z_AtWFace(i,j+1,k,z_nd);
535  amrex::Real z_y_l = Compute_Z_AtWFace(i,j-1,k,z_nd);
536  amrex::Real met_eta = 0.5 * (z_y_h - z_y_l) * dxInv[1];
537 
538  // Use extrapolation instead of interpolation if at the bottom boundary
539  amrex::Real u_hi = (k == 0) ? 1.5 * u_arr(i+1,j ,k ) - 0.5 * u_arr(i+1,j ,k+1) :
540  0.5 * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
541  amrex::Real u_lo = (k == 0) ? 1.5 * u_arr(i ,j ,k ) - 0.5 * u_arr(i ,j ,k+1) :
542  0.5 * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
543  amrex::Real mf_u_hi = mf_u(i+1,j,0);
544  amrex::Real mf_u_lo = mf_u(i ,j,0);
545 
546  amrex::Real v_hi = (k == 0) ? 1.5 * v_arr(i ,j+1,k ) - 0.5 * v_arr(i ,j+1,k+1) :
547  0.5 * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
548  amrex::Real v_lo = (k == 0) ? 1.5 * v_arr(i ,j ,k ) - 0.5 * v_arr(i ,j ,k+1) :
549  0.5 * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
550  amrex::Real mf_v_hi = mf_v(i,j+1,0);
551  amrex::Real mf_v_lo = mf_v(i,j ,0);
552 
553  amrex::Real u_met = met_xi * 0.5 * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
554  amrex::Real v_met = met_eta * 0.5 * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
555 
556  amrex::Real omega = w - u_met - v_met;
557  return omega;
558 }
559 
560 // Define w given u and v arrays and scalar omega (2nd order)
561 AMREX_GPU_DEVICE
562 AMREX_FORCE_INLINE
563 amrex::Real
564 WFromOmega (int& i, int& j, int& k, amrex::Real omega,
565  const amrex::Array4<const amrex::Real>& u_arr,
566  const amrex::Array4<const amrex::Real>& v_arr,
567  const amrex::Array4<const amrex::Real>& mf_u,
568  const amrex::Array4<const amrex::Real>& mf_v,
569  const amrex::Array4<const amrex::Real>& z_nd,
570  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& dxInv)
571 {
572  // This is dh/dxi at hi and lo edges
573  amrex::Real z_x_h = Compute_Z_AtWFace(i+1,j,k,z_nd);
574  amrex::Real z_x_l = Compute_Z_AtWFace(i-1,j,k,z_nd);
575  amrex::Real met_xi = 0.5 * (z_x_h - z_x_l) * dxInv[0];
576 
577  // This is dh/deta at hi and lo edges
578  amrex::Real z_y_h = Compute_Z_AtWFace(i,j+1,k,z_nd);
579  amrex::Real z_y_l = Compute_Z_AtWFace(i,j-1,k,z_nd);
580  amrex::Real met_eta = 0.5 * (z_y_h - z_y_l) * dxInv[1];
581 
582  // Use extrapolation instead of interpolation if at the bottom boundary
583  amrex::Real u_hi = (k == 0) ? 1.5 * u_arr(i+1,j ,k ) - 0.5 * u_arr(i+1,j ,k+1) :
584  0.5 * ( u_arr(i+1,j ,k-1) + u_arr(i+1,j ,k ) );
585  amrex::Real u_lo = (k == 0) ? 1.5 * u_arr(i ,j ,k ) - 0.5 * u_arr(i ,j ,k+1) :
586  0.5 * ( u_arr(i ,j ,k-1) + u_arr(i ,j ,k ) );
587  amrex::Real mf_u_hi = mf_u(i+1,j,0);
588  amrex::Real mf_u_lo = mf_u(i ,j,0);
589 
590  amrex::Real v_hi = (k == 0) ? 1.5 * v_arr(i ,j+1,k ) - 0.5 * v_arr(i ,j+1,k+1) :
591  0.5 * ( v_arr(i ,j+1,k-1) + v_arr(i ,j+1,k ) );
592  amrex::Real v_lo = (k == 0) ? 1.5 * v_arr(i ,j ,k ) - 0.5 * v_arr(i ,j ,k+1) :
593  0.5 * ( v_arr(i ,j ,k-1) + v_arr(i ,j ,k ) );
594  amrex::Real mf_v_hi = mf_v(i,j+1,0);
595  amrex::Real mf_v_lo = mf_v(i,j ,0);
596 
597  amrex::Real u_met = met_xi * 0.5 * ( u_hi*mf_u_hi + u_lo*mf_u_lo );
598  amrex::Real v_met = met_eta * 0.5 * ( v_hi*mf_v_hi + v_lo*mf_v_lo );
599 
600  amrex::Real w = omega + u_met + v_met;
601  return w;
602 }
603 #endif
604 
605 
606 //*****************************************************************************************
607 // Rotate scalar flux vector and stress tensor for MOST
608 //*****************************************************************************************
609 AMREX_GPU_DEVICE
610 AMREX_FORCE_INLINE
611 void
612 rotate_scalar_flux (const int& i,
613  const int& j,
614  const int& klo,
615  const amrex::Real& flux,
616  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& dxInv,
617  const amrex::Array4<const amrex::Real>& zphys_arr,
618  const amrex::Array4<amrex::Real>& phi1_arr,
619  const amrex::Array4<amrex::Real>& phi2_arr,
620  const amrex::Array4<amrex::Real>& phi3_arr)
621 {
622  amrex::Real h_xi = Compute_h_xi_AtCellCenter(i, j, klo, dxInv, zphys_arr);
623  amrex::Real h_eta = Compute_h_eta_AtCellCenter(i, j, klo, dxInv, zphys_arr);
624  amrex::Real InvNorm = 1.0 / std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
625  phi1_arr(i,j,klo) = -h_xi * flux * InvNorm;
626  phi2_arr(i,j,klo) = -h_eta * flux * InvNorm;
627  phi3_arr(i,j,klo) = flux * InvNorm;
628 }
629 
630 AMREX_GPU_DEVICE
631 AMREX_FORCE_INLINE
632 void
633 rotate_stress_tensor (const int& i,
634  const int& j,
635  const int& klo,
636  const amrex::Real& flux,
637  const amrex::GpuArray<amrex::Real, AMREX_SPACEDIM>& dxInv,
638  const amrex::Array4<const amrex::Real>& zphys_arr,
639  const amrex::Array4<const amrex::Real>& u_arr,
640  const amrex::Array4<const amrex::Real>& v_arr,
641  const amrex::Array4<const amrex::Real>& w_arr,
642  const amrex::Array4<amrex::Real>& tau11_arr,
643  const amrex::Array4<amrex::Real>& tau22_arr,
644  const amrex::Array4<amrex::Real>& tau33_arr,
645  const amrex::Array4<amrex::Real>& tau12_arr,
646  const amrex::Array4<amrex::Real>& tau21_arr,
647  const amrex::Array4<amrex::Real>& tau13_arr,
648  const amrex::Array4<amrex::Real>& tau31_arr,
649  const amrex::Array4<amrex::Real>& tau23_arr,
650  const amrex::Array4<amrex::Real>& tau32_arr)
651 {
652  // Unit-normal vector
653  amrex::Array1D<amrex::Real,0,2> n_hat;
654 
655  // Unit-tangent vectors
656  amrex::Array1D<amrex::Real,0,2> t_hat_1;
657  amrex::Array1D<amrex::Real,0,2> t_hat_2;
658  amrex::Array1D<amrex::Real,0,2> u_t_hat;
659 
660  // Final basis vector for right-hand coord sys
661  amrex::Array1D<amrex::Real,0,2> a_hat;
662 
663  // Rotation matrix
664  amrex::Array2D<amrex::Real,0,2,0,2> R_mat;
665 
666  // Metric data
667  amrex::Real h_xi = Compute_h_xi_AtIface(i, j, klo, dxInv, zphys_arr);
668  amrex::Real h_eta = Compute_h_eta_AtIface(i, j, klo, dxInv, zphys_arr);
669 
670  // Populate the normal vector
671  amrex::Real Inormn = 1./std::sqrt(1.0 + h_xi*h_xi + h_eta*h_eta);
672  n_hat(0) = -Inormn*h_xi; n_hat(1) = -Inormn*h_eta; n_hat(2) = Inormn;
673 
674  // Populate the tangent vectors
675  amrex::Real Inorm1 = 1./std::sqrt(1.0 + h_xi*h_xi);
676  amrex::Real Inorm2 = 1./std::sqrt(1.0 + h_eta*h_eta);
677  t_hat_1(0) = Inorm1; t_hat_2(1) = Inorm2;
678  t_hat_1(2) = Inorm1*h_xi; t_hat_2(2) = Inorm2*h_eta;
679 
680  // Populate the u_t vector
681  amrex::Real Norm_u_t = 0.0;
682  amrex::Real mag1 = (u_arr(i,j,klo) + h_xi *w_arr(i,j,klo))*Inorm1;
683  amrex::Real mag2 = (v_arr(i,j,klo) + h_eta*w_arr(i,j,klo))*Inorm2;
684  for (int icol(0); icol<3; ++icol) {
685  u_t_hat(icol) = mag1*t_hat_1(icol) + mag2*t_hat_2(icol);
686  Norm_u_t += u_t_hat(icol)*u_t_hat(icol);
687  }
688  for (int icol(0); icol<3; ++icol) {
689  u_t_hat(icol) /= std::sqrt(Norm_u_t);
690  }
691 
692  // Populate the a_hat vector
693  a_hat(0) = n_hat(1)*u_t_hat(2) - n_hat(2)*u_t_hat(1);
694  a_hat(1) = -(n_hat(0)*u_t_hat(2) - n_hat(2)*u_t_hat(0));
695  a_hat(2) = n_hat(1)*u_t_hat(1) - n_hat(1)*u_t_hat(0);
696 
697  // Copy column vectors into R_mat
698  int jrow;
699  jrow = 0;
700  for (int icol(0); icol<3; ++icol) {
701  R_mat(icol,jrow) = u_t_hat(icol);
702  }
703  jrow = 1;
704  for (int icol(0); icol<3; ++icol) {
705  R_mat(icol,jrow) = a_hat(icol);
706  }
707  jrow = 2;
708  for (int icol(0); icol<3; ++icol) {
709  R_mat(icol,jrow) = n_hat(icol);
710  }
711 
712  // Assign the stresses
713  tau11_arr(i,j,klo) = 2.0*R_mat(0,0)*R_mat(2,0)*flux;
714  tau22_arr(i,j,klo) = 2.0*R_mat(0,1)*R_mat(2,1)*flux;
715  tau33_arr(i,j,klo) = 2.0*R_mat(0,2)*R_mat(2,2)*flux;
716 
717  tau12_arr(i,j,klo) = (R_mat(0,0)*R_mat(2,1) + R_mat(2,0)*R_mat(0,1))*flux;
718  tau21_arr(i,j,klo) = tau12_arr(i,j,klo);
719 
720  tau13_arr(i,j,klo) = (R_mat(0,0)*R_mat(2,2) + R_mat(0,2)*R_mat(2,0))*flux;
721  tau31_arr(i,j,klo) = tau13_arr(i,j,klo);
722 
723  tau23_arr(i,j,klo) = (R_mat(0,1)*R_mat(2,2) + R_mat(0,2)*R_mat(2,1))*flux;
724  tau32_arr(i,j,klo) = tau23_arr(i,j,klo);
725 }
726 #endif
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real OmegaFromW(int &i, int &j, int &k, amrex::Real w, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:415
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterJ(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:305
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterK(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:244
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterJ(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:290
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:110
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterJ(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:276
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterK(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:231
void init_zlevels(amrex::Vector< amrex::Vector< amrex::Real >> &zlevels_stag, amrex::Vector< amrex::Vector< amrex::Real >> &stretched_dz_h, amrex::Vector< amrex::Gpu::DeviceVector< amrex::Real >> &stretched_dz_d, amrex::Vector< amrex::Geometry > const &geom, amrex::Vector< amrex::IntVect > const &ref_ratio, const amrex::Real grid_stretching_ratio, const amrex::Real zsurf, const amrex::Real dz0)
AMREX_FORCE_INLINE AMREX_GPU_DEVICE amrex::Real Compute_h_zeta_AtCellCenter(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:47
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtEdgeCenterI(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:320
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Zrel_AtCellCenter(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:390
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtCellCenter(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:77
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtKface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:182
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtCellCenter(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:362
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void rotate_scalar_flux(const int &i, const int &j, const int &klo, const amrex::Real &flux, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const amrex::Array4< const amrex::Real > &zphys_arr, const amrex::Array4< amrex::Real > &phi1_arr, const amrex::Array4< amrex::Real > &phi2_arr, const amrex::Array4< amrex::Real > &phi3_arr)
Definition: ERF_TerrainMetrics.H:612
void init_default_zphys(int lev, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::MultiFab &z_phys_cc)
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:96
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void rotate_stress_tensor(const int &i, const int &j, const int &klo, const amrex::Real &flux, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv, const amrex::Array4< const amrex::Real > &zphys_arr, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &w_arr, const amrex::Array4< amrex::Real > &tau11_arr, const amrex::Array4< amrex::Real > &tau22_arr, const amrex::Array4< amrex::Real > &tau33_arr, const amrex::Array4< amrex::Real > &tau12_arr, const amrex::Array4< amrex::Real > &tau21_arr, const amrex::Array4< amrex::Real > &tau13_arr, const amrex::Array4< amrex::Real > &tau31_arr, const amrex::Array4< amrex::Real > &tau23_arr, const amrex::Array4< amrex::Real > &tau32_arr)
Definition: ERF_TerrainMetrics.H:633
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:153
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_Z_AtWFace(const int &i, const int &j, const int &k, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:377
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtKface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:197
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterK(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:259
void init_which_terrain_grid(int lev, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::Vector< amrex::Real > const &z_levels_h)
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtIface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:125
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtEdgeCenterI(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:334
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_zeta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:139
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_xi_AtCellCenter(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:62
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtJface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:167
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtEdgeCenterI(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:347
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real Compute_h_eta_AtKface(const int &i, const int &j, const int &k, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &cellSizeInv, const amrex::Array4< const amrex::Real > &z_nd)
Definition: ERF_TerrainMetrics.H:211
void make_terrain_fitted_coords(int lev, const amrex::Geometry &geom, amrex::MultiFab &z_phys_nd, amrex::Vector< amrex::Real > const &z_levels_h, amrex::GpuArray< ERF_BC, AMREX_SPACEDIM *2 > &phys_bc_type)
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real WFromOmega(int &i, int &j, int &k, amrex::Real omega, const amrex::Array4< const amrex::Real > &u_arr, const amrex::Array4< const amrex::Real > &v_arr, const amrex::Array4< const amrex::Real > &mf_u, const amrex::Array4< const amrex::Real > &mf_v, const amrex::Array4< const amrex::Real > &z_nd, const amrex::GpuArray< amrex::Real, AMREX_SPACEDIM > &dxInv)
Definition: ERF_TerrainMetrics.H:465
@ omega
Definition: ERF_Morrison.H:53