ERF
Energy Research and Forecasting: An Atmospheric Modeling Code
ERF_Interpolation_Bilinear.H File Reference
#include "ERF_DataStruct.H"
Include dependency graph for ERF_Interpolation_Bilinear.H:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE int get_single_index (int i, int j, int k, int nx, int ny)
 
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)
 
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)
 

Function Documentation

◆ bilinear_interpolation()

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 
)
29 {
30 
31  if(z < 0.0){
32  tmp_var = 0.0;
33  return;
34  }
35 
36  int iloc=-1, jloc=-1, kloc=-1;
37  for(int k=0;k<nz;k++){
38  if(zvec[k] > z){
39  kloc = k-1;
40  break;
41  }
42  else if (zvec[k] == z) {
43  kloc = k;
44  }
45  }
46 
47  iloc = static_cast<int>(std::floor((x - xvec[0]) / dxvec));
48  jloc = static_cast<int>(std::floor((y - yvec[0]) / dyvec));
49 
50  if(iloc > nx-1 or iloc < 0 or
51  jloc > ny-1 or iloc < 0 or
52  kloc > nz-1 or kloc < 0){
53  printf("There is a problem inside here %d %d %d %d %d %d\n", iloc, jloc, kloc, nx-1, ny-1, nz-1);
54  for(int k=0;k<nz;k++){
55  printf("Value of z is %d %0.15g %0.15g\n", k, zvec[k], z);
56  }
57  }
58 
59  amrex::Real xlo = xvec[0] + iloc*dxvec;
60  amrex::Real ylo = yvec[0] + jloc*dyvec;
61  amrex::Real zlo = zvec[kloc];
62 
63  amrex::Real w_x = (x - xlo)/dxvec;
64  amrex::Real w_y = (y - ylo)/dyvec;
65  amrex::Real w_z = (z - zlo)/(zvec[kloc+1] - zvec[kloc]);
66 
67  int ind0 = get_single_index(iloc,jloc,kloc,nx,ny);
68  int ind1 = get_single_index(iloc+1,jloc,kloc,nx,ny);
69  int ind2 = get_single_index(iloc+1,jloc+1,kloc,nx,ny);
70  int ind3 = get_single_index(iloc,jloc+1,kloc,nx,ny);
71  int ind4 = get_single_index(iloc,jloc,kloc+1,nx,ny);
72  int ind5 = get_single_index(iloc+1,jloc,kloc+1,nx,ny);
73  int ind6 = get_single_index(iloc+1,jloc+1,kloc+1,nx,ny);
74  int ind7 = get_single_index(iloc,jloc+1,kloc+1,nx,ny);
75 
76  tmp_var = (1-w_x)*(1-w_y)*(1-w_z)*varvec[ind0] + w_x*(1-w_y)*(1-w_z)*varvec[ind1] +
77  w_x*w_y*(1-w_z)*varvec[ind2] + (1-w_x)*w_y*(1-w_z)*varvec[ind3] +
78  (1-w_x)*(1-w_y)*w_z*varvec[ind4] + w_x*(1-w_y)*w_z*varvec[ind5] +
79  w_x*w_y*w_z*varvec[ind6] + (1-w_x)*w_y*w_z*varvec[ind7];
80 
81  //std::cout << "Variable value is " << tmp_var << "\n";
82 }
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:13
amrex::Real Real
Definition: ERF_ShocInterface.H:19

Referenced by ERF::FillForecastStateMultiFabs().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ get_single_index()

AMREX_FORCE_INLINE AMREX_GPU_HOST_DEVICE int get_single_index ( int  i,
int  j,
int  k,
int  nx,
int  ny 
)

Bilinear interpolation routines. There are different versions for different scenarios that take in a different argument list

15 {
16 
17  int si = k*nx*ny + j*nx + i;
18  return si;
19 }

Referenced by bilinear_interpolation().

Here is the caller graph for this function:

◆ interpolate_from_coarse()

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 
)
inline
91 {
92  // Convert to coarse grid index space
93  amrex::Real ic_real = (x - prob_lo_crse[0]) / dx_crse[0];
94  amrex::Real jc_real = (y - prob_lo_crse[1]) / dx_crse[1];
95  amrex::Real kc_real = (z - prob_lo_crse[2]) / dx_crse[2];
96 
97  int ic = static_cast<int>(amrex::Math::floor(ic_real));
98  int jc = static_cast<int>(amrex::Math::floor(jc_real));
99  int kc = static_cast<int>(amrex::Math::floor(kc_real));
100 
101  amrex::Real dx = ic_real - ic;
102  amrex::Real dy = jc_real - jc;
103  amrex::Real dz = kc_real - kc;
104 
105  // Trilinear interpolation
106  amrex::Real v000 = crse(ic, jc, kc, n);
107  amrex::Real v100 = crse(ic + 1, jc, kc, n);
108  amrex::Real v010 = crse(ic, jc + 1, kc, n);
109  amrex::Real v110 = crse(ic + 1, jc + 1, kc, n);
110  amrex::Real v001 = crse(ic, jc, kc + 1, n);
111  amrex::Real v101 = crse(ic + 1, jc, kc + 1, n);
112  amrex::Real v011 = crse(ic, jc + 1, kc + 1, n);
113  amrex::Real v111 = crse(ic + 1, jc + 1, kc + 1, n);
114 
115  amrex::Real v = v000 * (1 - dx) * (1 - dy) * (1 - dz) +
116  v100 * dx * (1 - dy) * (1 - dz) +
117  v010 * (1 - dx) * dy * (1 - dz) +
118  v110 * dx * dy * (1 - dz) +
119  v001 * (1 - dx) * (1 - dy) * dz +
120  v101 * dx * (1 - dy) * dz +
121  v011 * (1 - dx) * dy * dz +
122  v111 * dx * dy * dz;
123 
124  return v;
125 }