ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Interpolation_Bilinear.H
Go to the documentation of this file.
1 #ifndef ERF_INTERPOLATE_BILINEAR_H_
2 #define ERF_INTERPOLATE_BILINEAR_H_
3 
4 #include "AMReX_GpuPrint.H"
5 #include "ERF_DataStruct.H"
6 
7 /**
8  * Bilinear interpolation routines. There are different versions for
9  * different scenarios that take in a different argument list
10  */
11 
12 AMREX_FORCE_INLINE
13 AMREX_GPU_HOST_DEVICE
14 int get_single_index(int i, int j, int k,
15  int nx, int ny)
16 {
17 
18  int si = k*nx*ny + j*nx + i;
19  return si;
20 }
21 
22 AMREX_FORCE_INLINE
23 AMREX_GPU_HOST_DEVICE
24 int clamp_interp_index (int idx, int npts)
25 {
26  if (npts <= 1) {
27  return 0;
28  }
29 
30  return amrex::max(0, amrex::min(idx, npts - 2));
31 }
32 
33 AMREX_FORCE_INLINE
34 AMREX_GPU_HOST_DEVICE
36 {
37  return (xhi > xlo) ? (x - xlo) / (xhi - xlo) : amrex::Real(0);
38 }
39 
40 AMREX_FORCE_INLINE
41 AMREX_GPU_HOST_DEVICE
42 void bilinear_interpolation(const amrex::Real* xvec, const amrex::Real* yvec, const amrex::Real* zvec,
43  const amrex::Real dxvec, const amrex::Real dyvec,
44  const int nx, const int ny, const int nz,
45  const amrex::Real x, const amrex::Real y, const amrex::Real z,
46  const amrex::Real* varvec,
47  amrex::Real& tmp_var)
48 {
49 
50  if(z < amrex::Real(0)){
51  tmp_var = amrex::Real(0);
52  return;
53  }
54 
55  amrex::Real x_interp = x;
56  amrex::Real y_interp = y;
57  amrex::Real z_interp = z;
58 
59  if (x_interp < xvec[0]) x_interp = xvec[0];
60  if (x_interp > xvec[nx-1]) x_interp = xvec[nx-1];
61 
62  if (y_interp < yvec[0]) y_interp = yvec[0];
63  if (y_interp > yvec[ny-1]) y_interp = yvec[ny-1];
64 
65  if (z_interp < zvec[0]) z_interp = zvec[0];
66  if (z_interp > zvec[nz-1]) z_interp = zvec[nz-1];
67 
68  int iloc = clamp_interp_index(static_cast<int>(std::floor((x_interp - xvec[0]) / dxvec)), nx);
69  int jloc = clamp_interp_index(static_cast<int>(std::floor((y_interp - yvec[0]) / dyvec)), ny);
70 
71  int kloc = 0;
72  if (nz > 1) {
73  kloc = nz - 2;
74  for (int k = 1; k < nz; ++k) {
75  if (zvec[k] >= z_interp) {
76  kloc = k - 1;
77  break;
78  }
79  }
80  }
81 
82  int ihi = (nx > 1) ? iloc + 1 : iloc;
83  int jhi = (ny > 1) ? jloc + 1 : jloc;
84  int khi = (nz > 1) ? kloc + 1 : kloc;
85 
86  amrex::Real xlo = xvec[iloc];
87  amrex::Real ylo = yvec[jloc];
88  amrex::Real zlo = zvec[kloc];
89 
90  amrex::Real w_x = get_interp_weight(x_interp, xlo, xvec[ihi]);
91  amrex::Real w_y = get_interp_weight(y_interp, ylo, yvec[jhi]);
92  amrex::Real w_z = get_interp_weight(z_interp, zlo, zvec[khi]);
93 
94  int ind0 = get_single_index(iloc,jloc,kloc,nx,ny);
95  int ind1 = get_single_index(ihi ,jloc,kloc,nx,ny);
96  int ind2 = get_single_index(ihi ,jhi ,kloc,nx,ny);
97  int ind3 = get_single_index(iloc,jhi ,kloc,nx,ny);
98  int ind4 = get_single_index(iloc,jloc,khi ,nx,ny);
99  int ind5 = get_single_index(ihi ,jloc,khi ,nx,ny);
100  int ind6 = get_single_index(ihi ,jhi ,khi ,nx,ny);
101  int ind7 = get_single_index(iloc,jhi ,khi ,nx,ny);
102 
103  tmp_var = (1-w_x)*(1-w_y)*(1-w_z)*varvec[ind0] + w_x*(1-w_y)*(1-w_z)*varvec[ind1] +
104  w_x*w_y*(1-w_z)*varvec[ind2] + (1-w_x)*w_y*(1-w_z)*varvec[ind3] +
105  (1-w_x)*(1-w_y)*w_z*varvec[ind4] + w_x*(1-w_y)*w_z*varvec[ind5] +
106  w_x*w_y*w_z*varvec[ind6] + (1-w_x)*w_y*w_z*varvec[ind7];
107 
108  //std::cout << "Variable value is " << tmp_var << "\n";
109 }
110 
111 AMREX_GPU_DEVICE
112 inline
113 amrex::Real interpolate_from_coarse(const amrex::Array4<const amrex::Real>& crse,
114  int n,
116  const amrex::Real* prob_lo_crse,
117  const amrex::Real* dx_crse)
118 {
119  // Convert to coarse grid index space
120  amrex::Real ic_real = (x - prob_lo_crse[0]) / dx_crse[0];
121  amrex::Real jc_real = (y - prob_lo_crse[1]) / dx_crse[1];
122  amrex::Real kc_real = (z - prob_lo_crse[2]) / dx_crse[2];
123 
124  int ic = static_cast<int>(amrex::Math::floor(ic_real));
125  int jc = static_cast<int>(amrex::Math::floor(jc_real));
126  int kc = static_cast<int>(amrex::Math::floor(kc_real));
127 
128  amrex::Real dx = ic_real - ic;
129  amrex::Real dy = jc_real - jc;
130  amrex::Real dz = kc_real - kc;
131 
132  // Trilinear interpolation
133  amrex::Real v000 = crse(ic, jc, kc, n);
134  amrex::Real v100 = crse(ic + 1, jc, kc, n);
135  amrex::Real v010 = crse(ic, jc + 1, kc, n);
136  amrex::Real v110 = crse(ic + 1, jc + 1, kc, n);
137  amrex::Real v001 = crse(ic, jc, kc + 1, n);
138  amrex::Real v101 = crse(ic + 1, jc, kc + 1, n);
139  amrex::Real v011 = crse(ic, jc + 1, kc + 1, n);
140  amrex::Real v111 = crse(ic + 1, jc + 1, kc + 1, n);
141 
142  amrex::Real v = v000 * (1 - dx) * (1 - dy) * (1 - dz) +
143  v100 * dx * (1 - dy) * (1 - dz) +
144  v010 * (1 - dx) * dy * (1 - dz) +
145  v110 * dx * dy * (1 - dz) +
146  v001 * (1 - dx) * (1 - dy) * dz +
147  v101 * dx * (1 - dy) * dz +
148  v011 * (1 - dx) * dy * dz +
149  v111 * dx * dy * dz;
150 
151  return v;
152 }
153 
154 AMREX_FORCE_INLINE
155 AMREX_GPU_HOST_DEVICE
156 void bilinear_interpolation_2d(const amrex::Real* xvec, const amrex::Real* yvec,
157  const amrex::Real dxvec, const amrex::Real dyvec,
158  const int nx, const int ny,
160  const amrex::Real* varvec,
161  amrex::Real& tmp_var)
162 {
163  if(x<xvec[0]) {
164  x = xvec[0];
165  }
166  if(x>xvec[nx-1]) {
167  x = xvec[nx-1];
168  }
169 
170  if(y<yvec[0]) {
171  y = yvec[0];
172  }
173  if(y>yvec[ny-1]) {
174  y = yvec[ny-1];
175  }
176 
177  int iloc = clamp_interp_index(static_cast<int>(std::floor((x-xvec[0])/dxvec)), nx);
178  int jloc = clamp_interp_index(static_cast<int>(std::floor((y-yvec[0])/dyvec)), ny);
179  int ihi = (nx > 1) ? iloc + 1 : iloc;
180  int jhi = (ny > 1) ? jloc + 1 : jloc;
181 
182  amrex::Real xlo = xvec[iloc];
183  amrex::Real ylo = yvec[jloc];
184 
185  amrex::Real w_x = get_interp_weight(x, xlo, xvec[ihi]);
186  amrex::Real w_y = get_interp_weight(y, ylo, yvec[jhi]);
187 
188  int ind0 = get_single_index(iloc,jloc,0,nx,ny);
189  int ind1 = get_single_index(ihi ,jloc,0,nx,ny);
190  int ind2 = get_single_index(ihi ,jhi ,0,nx,ny);
191  int ind3 = get_single_index(iloc,jhi ,0,nx,ny);
192 
193  tmp_var = (1-w_x)*(1-w_y)*varvec[ind0] + w_x*(1-w_y)*varvec[ind1] +
194  w_x*w_y*varvec[ind2] + (1-w_x)*w_y*varvec[ind3];
195 
196  //std::cout << "Variable value is " << tmp_var << "\n";
197 }
198 
199 
200 #endif
const Real dy
Definition: ERF_InitCustomPert_ABL.H:24
const Real dx
Definition: ERF_InitCustomPert_ABL.H:23
const int khi
Definition: ERF_InitCustomPert_Bubble.H:21
const Real dz
Definition: ERF_InitCustomPert_Bubble.H:25
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE void bilinear_interpolation_2d(const amrex::Real *xvec, const amrex::Real *yvec, const amrex::Real dxvec, const amrex::Real dyvec, const int nx, const int ny, amrex::Real x, amrex::Real y, const amrex::Real *varvec, amrex::Real &tmp_var)
Definition: ERF_Interpolation_Bilinear.H:156
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE int get_single_index(int i, int j, int k, int nx, int ny)
Definition: ERF_Interpolation_Bilinear.H:14
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE int clamp_interp_index(int idx, int npts)
Definition: ERF_Interpolation_Bilinear.H:24
AMREX_GPU_DEVICE amrex::Real interpolate_from_coarse(const amrex::Array4< const amrex::Real > &crse, int n, amrex::Real x, amrex::Real y, amrex::Real z, const amrex::Real *prob_lo_crse, const amrex::Real *dx_crse)
Definition: ERF_Interpolation_Bilinear.H:113
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE void bilinear_interpolation(const amrex::Real *xvec, const amrex::Real *yvec, const amrex::Real *zvec, const amrex::Real dxvec, const amrex::Real dyvec, const int nx, const int ny, const int nz, const amrex::Real x, const amrex::Real y, const amrex::Real z, const amrex::Real *varvec, amrex::Real &tmp_var)
Definition: ERF_Interpolation_Bilinear.H:42
AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE amrex::Real get_interp_weight(amrex::Real x, amrex::Real xlo, amrex::Real xhi)
Definition: ERF_Interpolation_Bilinear.H:35
amrex::Real Real
Definition: ERF_ShocInterface.H:19