11 use fckit_log_module,
only: fckit_log
31 type(
qg_geom),
intent(in) :: geom
32 real(kind_real),
intent(in) :: lon
33 real(kind_real),
intent(in) :: lat
34 real(kind_real),
intent(in) :: z
35 real(kind_real),
intent(in) :: field(geom%nx,geom%ny,geom%nz)
36 real(kind_real),
intent(out) :: val
39 integer :: jxm1,jxo,jxp1,jxp2
40 integer :: jym1,jyo,jyp1,jyp2
42 real(kind_real) :: x,y
43 real(kind_real) :: ax,ay,az
44 real(kind_real) :: oo,op1,p1o,p1p1,o,p1
56 ay = ay+real(jyo-1,kind_real)
60 if (jyp1>geom%ny)
then
61 ay = ay+real(jyp1-geom%ny,kind_real)
67 oo = (1.0-ax)*field(jxo,jyo,jzo)+ax*field(jxp1,jyo,jzo)
68 op1 = (1.0-ax)*field(jxo,jyo,jzp1)+ax*field(jxp1,jyo,jzp1)
69 p1o = (1.0-ax)*field(jxo,jyp1,jzo)+ax*field(jxp1,jyp1,jzo)
70 p1p1 = (1.0-ax)*field(jxo,jyp1,jzp1)+ax*field(jxp1,jyp1,jzp1)
73 o = (1.0-ay)*oo+ay*p1o
74 p1 = (1.0-ay)*op1+ay*p1p1
77 val = (1.0-az)*o+az*p1
85 type(
qg_geom),
intent(in) :: geom
86 real(kind_real),
intent(in) :: lon
87 real(kind_real),
intent(in) :: lat
88 real(kind_real),
intent(in) :: z
89 real(kind_real),
intent(in) :: val
90 real(kind_real),
intent(inout) :: field(geom%nx,geom%ny,geom%nz)
93 integer :: jxm1,jxo,jxp1,jxp2
94 integer :: jym1,jyo,jyp1,jyp2
96 real(kind_real) :: x,y
97 real(kind_real) :: ax,ay,az
98 real(kind_real) :: oo,op1,p1o,p1p1,o,p1
110 ay = ay+real(jyo-1,kind_real)
114 if (jyp1>geom%ny)
then
115 ay = ay+real(jyp1-geom%ny,kind_real)
131 field(jxo,jyo,jzo) = field(jxo,jyo,jzo)+(1.0-ax)*oo
132 field(jxp1,jyo,jzo) = field(jxp1,jyo,jzo)+ax*oo
133 field(jxo,jyo,jzp1) = field(jxo,jyo,jzp1)+(1.0-ax)*op1
134 field(jxp1,jyo,jzp1) = field(jxp1,jyo,jzp1)+ax*op1
135 field(jxo,jyp1,jzo) = field(jxo,jyp1,jzo)+(1.0-ax)*p1o
136 field(jxp1,jyp1,jzo) = field(jxp1,jyp1,jzo)+ax*p1o
137 field(jxo,jyp1,jzp1) = field(jxo,jyp1,jzp1)+(1.0-ax)*p1p1
138 field(jxp1,jyp1,jzp1) = field(jxp1,jyp1,jzp1)+ax*p1p1
146 type(
qg_geom),
intent(in) :: geom
147 real(kind_real),
intent(in) :: x
148 real(kind_real),
intent(in) :: y
149 real(kind_real),
intent(in) :: gfld2dext(geom%nx,-1:geom%ny+2)
150 real(kind_real),
intent(out) :: val
153 integer :: jxm1,jxo,jxp1,jxp2
154 integer :: jym1,jyo,jyp1,jyp2
155 real(kind_real) :: ax,ay,m1,o,p1,p2
162 call cubic(ax,gfld2dext(jxm1,jym1),gfld2dext(jxo,jym1),gfld2dext(jxp1,jym1),gfld2dext(jxp2,jym1),m1)
163 call cubic(ax,gfld2dext(jxm1,jyo),gfld2dext(jxo,jyo),gfld2dext(jxp1,jyo),gfld2dext(jxp2,jyo),o)
164 call cubic(ax,gfld2dext(jxm1,jyp1),gfld2dext(jxo,jyp1),gfld2dext(jxp1,jyp1),gfld2dext(jxp2,jyp1),p1)
165 call cubic(ax,gfld2dext(jxm1,jyp2),gfld2dext(jxo,jyp2),gfld2dext(jxp1,jyp2),gfld2dext(jxp2,jyp2),p2)
168 call cubic(ay,m1,o,p1,p2,val)
176 type(
qg_geom),
intent(in) :: geom
177 real(kind_real),
intent(in) :: x
178 real(kind_real),
intent(in) :: y
179 real(kind_real),
intent(in) :: gfld2dext_traj(geom%nx,-1:geom%ny+2)
180 real(kind_real),
intent(in) :: dx
181 real(kind_real),
intent(in) :: dy
182 real(kind_real),
intent(in) :: gfld2dext(geom%nx,-1:geom%ny+2)
183 real(kind_real),
intent(out) :: val
186 integer :: jxm1,jxo,jxp1,jxp2
187 integer :: jym1,jyo,jyp1,jyp2
188 real(kind_real) :: ax_traj,ay_traj,m1_traj,o_traj,p1_traj,p2_traj,ax,ay,m1,o,p1,p2
195 call cubic(ax_traj,gfld2dext_traj(jxm1,jym1),gfld2dext_traj(jxo,jym1),gfld2dext_traj(jxp1,jym1),gfld2dext_traj(jxp2,jym1),m1_traj)
196 call cubic(ax_traj,gfld2dext_traj(jxm1,jyo),gfld2dext_traj(jxo,jyo),gfld2dext_traj(jxp1,jyo),gfld2dext_traj(jxp2,jyo),o_traj)
197 call cubic(ax_traj,gfld2dext_traj(jxm1,jyp1),gfld2dext_traj(jxo,jyp1),gfld2dext_traj(jxp1,jyp1),gfld2dext_traj(jxp2,jyp1),p1_traj)
198 call cubic(ax_traj,gfld2dext_traj(jxm1,jyp2),gfld2dext_traj(jxo,jyp2),gfld2dext_traj(jxp1,jyp2),gfld2dext_traj(jxp2,jyp2),p2_traj)
205 call cubic_tl(ax_traj,gfld2dext_traj(jxm1,jym1),gfld2dext_traj(jxo,jym1),gfld2dext_traj(jxp1,jym1),gfld2dext_traj(jxp2,jym1), &
206 & ax,gfld2dext(jxm1,jym1),gfld2dext(jxo,jym1),gfld2dext(jxp1,jym1),gfld2dext(jxp2,jym1),m1)
207 call cubic_tl(ax_traj,gfld2dext_traj(jxm1,jyo),gfld2dext_traj(jxo,jyo),gfld2dext_traj(jxp1,jyo),gfld2dext_traj(jxp2,jyo), &
208 & ax,gfld2dext(jxm1,jyo),gfld2dext(jxo,jyo),gfld2dext(jxp1,jyo),gfld2dext(jxp2,jyo),o)
209 call cubic_tl(ax_traj,gfld2dext_traj(jxm1,jyp1),gfld2dext_traj(jxo,jyp1),gfld2dext_traj(jxp1,jyp1),gfld2dext_traj(jxp2,jyp1), &
210 & ax,gfld2dext(jxm1,jyp1),gfld2dext(jxo,jyp1),gfld2dext(jxp1,jyp1),gfld2dext(jxp2,jyp1),p1)
211 call cubic_tl(ax_traj,gfld2dext_traj(jxm1,jyp2),gfld2dext_traj(jxo,jyp2),gfld2dext_traj(jxp1,jyp2),gfld2dext_traj(jxp2,jyp2), &
212 & ax,gfld2dext(jxm1,jyp2),gfld2dext(jxo,jyp2),gfld2dext(jxp1,jyp2),gfld2dext(jxp2,jyp2),p2)
215 call cubic_tl(ay_traj,m1_traj,o_traj,p1_traj,p2_traj,ay,m1,o,p1,p2,val)
223 type(
qg_geom),
intent(in) :: geom
224 real(kind_real),
intent(in) :: x
225 real(kind_real),
intent(in) :: y
226 real(kind_real),
intent(in) :: gfld2dext_traj(geom%nx,-1:geom%ny+2)
227 real(kind_real),
intent(in) :: val
228 real(kind_real),
intent(inout) :: dx
229 real(kind_real),
intent(inout) :: dy
230 real(kind_real),
intent(inout) :: gfld2dext(geom%nx,-1:geom%ny+2)
233 integer :: jxm1,jxo,jxp1,jxp2
234 integer :: jym1,jyo,jyp1,jyp2
235 real(kind_real) :: ax_traj,ay_traj,m1_traj,o_traj,p1_traj,p2_traj,ax,ay,m1,o,p1,p2
242 call cubic(ax_traj,gfld2dext_traj(jxm1,jym1),gfld2dext_traj(jxo,jym1),gfld2dext_traj(jxp1,jym1),gfld2dext_traj(jxp2,jym1),m1_traj)
243 call cubic(ax_traj,gfld2dext_traj(jxm1,jyo),gfld2dext_traj(jxo,jyo),gfld2dext_traj(jxp1,jyo),gfld2dext_traj(jxp2,jyo),o_traj)
244 call cubic(ax_traj,gfld2dext_traj(jxm1,jyp1),gfld2dext_traj(jxo,jyp1),gfld2dext_traj(jxp1,jyp1),gfld2dext_traj(jxp2,jyp1),p1_traj)
245 call cubic(ax_traj,gfld2dext_traj(jxm1,jyp2),gfld2dext_traj(jxo,jyp2),gfld2dext_traj(jxp1,jyp2),gfld2dext_traj(jxp2,jyp2),p2_traj)
256 call cubic_ad(ay_traj,m1_traj,o_traj,p1_traj,p2_traj,val,ay,m1,o,p1,p2)
259 call cubic_ad(ax_traj,gfld2dext_traj(jxm1,jym1),gfld2dext_traj(jxo,jym1),gfld2dext_traj(jxp1,jym1),gfld2dext_traj(jxp2,jym1), &
260 & m1,ax,gfld2dext(jxm1,jym1),gfld2dext(jxo,jym1),gfld2dext(jxp1,jym1),gfld2dext(jxp2,jym1))
261 call cubic_ad(ax_traj,gfld2dext_traj(jxm1,jyo),gfld2dext_traj(jxo,jyo),gfld2dext_traj(jxp1,jyo),gfld2dext_traj(jxp2,jyo), &
262 & o,ax,gfld2dext(jxm1,jyo),gfld2dext(jxo,jyo),gfld2dext(jxp1,jyo),gfld2dext(jxp2,jyo))
263 call cubic_ad(ax_traj,gfld2dext_traj(jxm1,jyp1),gfld2dext_traj(jxo,jyp1),gfld2dext_traj(jxp1,jyp1),gfld2dext_traj(jxp2,jyp1), &
264 & p1,ax,gfld2dext(jxm1,jyp1),gfld2dext(jxo,jyp1),gfld2dext(jxp1,jyp1),gfld2dext(jxp2,jyp1))
265 call cubic_ad(ax_traj,gfld2dext_traj(jxm1,jyp2),gfld2dext_traj(jxo,jyp2),gfld2dext_traj(jxp1,jyp2),gfld2dext_traj(jxp2,jyp2), &
266 & p2,ax,gfld2dext(jxm1,jyp2),gfld2dext(jxo,jyp2),gfld2dext(jxp1,jyp2),gfld2dext(jxp2,jyp2))
269 dx = dx+ax/geom%deltax
270 dy = dy+ay/geom%deltay
278 type(
qg_geom),
intent(in) :: geom
279 real(kind_real),
intent(in) :: x
280 integer,
intent(out) :: jxm1
281 integer,
intent(out) :: jxo
282 integer,
intent(out) :: jxp1
283 integer,
intent(out) :: jxp2
284 real(kind_real),
intent(out) :: ax
288 real(kind_real) :: xx
300 if ((0.0<=xx).and.(xx<geom%x(1)))
then
303 elseif ((geom%x(geom%nx)<=xx).and.(xx<=
domain_zonal))
then
305 ax = (xx-geom%x(jxo))/geom%deltax
309 if ((geom%x(ix)<=xx).and.(xx<geom%x(ix+1)))
then
314 if (jxo==-1)
call abor1_ftn(
'find_x_indices: cannot find jxo')
315 ax = (xx-geom%x(jxo))/geom%deltax
320 if (jxm1<1) jxm1 = jxm1+geom%nx
322 if (jxp1>geom%nx) jxp1 = jxp1-geom%nx
324 if (jxp2>geom%nx) jxp2 = jxp2-geom%nx
332 type(
qg_geom),
intent(in) :: geom
333 real(kind_real),
intent(in) :: y
334 integer,
intent(out) :: jym1
335 integer,
intent(out) :: jyo
336 integer,
intent(out) :: jyp1
337 integer,
intent(out) :: jyp2
338 real(kind_real),
intent(out) :: ay
342 real(kind_real) :: yy
348 if ((0.0<=yy).and.(yy<geom%y(1)))
then
353 ay = (yy-geom%y(jyo))/geom%deltay
357 if ((geom%y(iy)<=yy).and.(yy<geom%y(iy+1)))
then
362 if (jyo==-1)
call abor1_ftn(
'find_y_indices: cannot find jyo')
363 ay = (yy-geom%y(jyo))/geom%deltay
377 type(
qg_geom),
intent(in) :: geom
378 real(kind_real),
intent(in) :: z
379 integer,
intent(out) :: jzo
380 integer,
intent(out) :: jzp1
381 real(kind_real),
intent(out) :: az
385 real(kind_real) :: zz
391 if ((0.0<=zz).and.(zz<geom%z(1)))
then
393 elseif ((geom%z(geom%nz)<=zz).and.(zz<=
domain_depth))
then
398 if ((geom%z(iz)<=zz).and.(zz<geom%z(iz+1)))
then
403 if (jzo==-1)
call abor1_ftn(
'find_z_indices: cannot find jzo')
405 az = (zz-geom%z(jzo))/(geom%z(jzo+1)-geom%z(jzo))
416 real(kind_real),
intent(in) :: a
417 real(kind_real),
intent(in) :: m1
418 real(kind_real),
intent(in) :: o
419 real(kind_real),
intent(in) :: p1
420 real(kind_real),
intent(in) :: p2
421 real(kind_real),
intent(out) :: res
423 res = -a*(a-1.0)*(a-2.0)/6.0*m1 &
424 & +(a+1.0)*(a-1.0)*(a-2.0)/2.0*o &
425 & -(a+1.0)*a*(a-2.0)/2.0*p1 &
426 & +(a+1.0)*a*(a-1.0)/6.0*p2
431 subroutine cubic_tl(a_traj,m1_traj,o_traj,p1_traj,p2_traj,a,m1,o,p1,p2,res)
434 real(kind_real),
intent(in) :: a_traj
435 real(kind_real),
intent(in) :: m1_traj
436 real(kind_real),
intent(in) :: o_traj
437 real(kind_real),
intent(in) :: p1_traj
438 real(kind_real),
intent(in) :: p2_traj
439 real(kind_real),
intent(in) :: a
440 real(kind_real),
intent(in) :: m1
441 real(kind_real),
intent(in) :: o
442 real(kind_real),
intent(in) :: p1
443 real(kind_real),
intent(in) :: p2
444 real(kind_real),
intent(out) :: res
446 res = -a*(a_traj-1.0)*(a_traj-2.0)/6.0*m1_traj &
447 & -a_traj*a*(a_traj-2.0)/6.0*m1_traj &
448 & -a_traj*(a_traj-1.0)*a/6.0*m1_traj &
449 & -a_traj*(a_traj-1.0)*(a_traj-2.0)/6.0*m1 &
450 & +a*(a_traj-1.0)*(a_traj-2.0)/2.0*o_traj &
451 & +(a_traj+1.0)*a*(a_traj-2.0)/2.0*o_traj &
452 & +(a_traj+1.0)*(a_traj-1.0)*a/2.0*o_traj &
453 & +(a_traj+1.0)*(a_traj-1.0)*(a_traj-2.0)/2.0*o &
454 & -a*a_traj*(a_traj-2.0)/2.0*p1_traj &
455 & -(a_traj+1.0)*a*(a_traj-2.0)/2.0*p1_traj &
456 & -(a_traj+1.0)*a_traj*a/2.0*p1_traj &
457 & -(a_traj+1.0)*a_traj*(a_traj-2.0)/2.0*p1 &
458 & +a*a_traj*(a_traj-1.0)/6.0*p2_traj &
459 & +(a_traj+1.0)*a*(a_traj-1.0)/6.0*p2_traj &
460 & +(a_traj+1.0)*a_traj*a/6.0*p2_traj &
461 & +(a_traj+1.0)*a_traj*(a_traj-1.0)/6.0*p2
466 subroutine cubic_ad(a_traj,m1_traj,o_traj,p1_traj,p2_traj,res,a,m1,o,p1,p2)
469 real(kind_real),
intent(in) :: a_traj
470 real(kind_real),
intent(in) :: m1_traj
471 real(kind_real),
intent(in) :: o_traj
472 real(kind_real),
intent(in) :: p1_traj
473 real(kind_real),
intent(in) :: p2_traj
474 real(kind_real),
intent(in) :: res
475 real(kind_real),
intent(inout) :: a
476 real(kind_real),
intent(inout) :: m1
477 real(kind_real),
intent(inout) :: o
478 real(kind_real),
intent(inout) :: p1
479 real(kind_real),
intent(inout) :: p2
481 a = a+(-(a_traj-1.0)*(a_traj-2.0)/6.0*m1_traj &
482 & -a_traj*(a_traj-2.0)/6.0*m1_traj &
483 & -a_traj*(a_traj-1.0)/6.0*m1_traj &
484 & +(a_traj-1.0)*(a_traj-2.0)/2.0*o_traj &
485 & +(a_traj+1.0)*(a_traj-2.0)/2.0*o_traj &
486 & +(a_traj+1.0)*(a_traj-1.0)/2.0*o_traj &
487 & -a_traj*(a_traj-2.0)/2.0*p1_traj &
488 & -(a_traj+1.0)*(a_traj-2.0)/2.0*p1_traj &
489 & -(a_traj+1.0)*a_traj/2.0*p1_traj &
490 & +a_traj*(a_traj-1.0)/6.0*p2_traj &
491 & +(a_traj+1.0)*(a_traj-1.0)/6.0*p2_traj &
492 & +(a_traj+1.0)*a_traj/6.0*p2_traj) &
494 m1 = m1-a_traj*(a_traj-1.0)*(a_traj-2.0)/6.0*res
495 o = o+(a_traj+1.0)*(a_traj-1.0)*(a_traj-2.0)/2.0*res
496 p1 = p1-(a_traj+1.0)*a_traj*(a_traj-2.0)/2.0*res
497 p2 = p2+(a_traj+1.0)*a_traj*(a_traj-1.0)/6.0*res
real(kind_real), parameter, public domain_depth
Model depth (m)
real(kind_real), parameter, public domain_meridional
Model domain in meridional direction (m)
real(kind_real), parameter, public domain_zonal
Model domain in zonal direction (m)
subroutine, public find_x_indices(geom, x, jxm1, jxo, jxp1, jxp2, ax)
Find x indices.
subroutine, public qg_interp_bicubic_ad(geom, x, y, gfld2dext_traj, val, dx, dy, gfld2dext)
Bicubic horizontal interpolation - adjoint.
subroutine find_z_indices(geom, z, jzo, jzp1, az)
Find z indices.
subroutine cubic(a, m1, o, p1, p2, res)
Cubic interpolation.
subroutine, public qg_interp_bicubic_tl(geom, x, y, gfld2dext_traj, dx, dy, gfld2dext, val)
Bicubic horizontal interpolation - tangent linear.
subroutine cubic_ad(a_traj, m1_traj, o_traj, p1_traj, p2_traj, res, a, m1, o, p1, p2)
Cubic interpolation - adjoint.
subroutine, public qg_interp_trilinear_ad(geom, lon, lat, z, val, field)
Interpolation - adjoint.
subroutine, public find_y_indices(geom, y, jym1, jyo, jyp1, jyp2, ay)
Find y indices.
subroutine, public qg_interp_trilinear(geom, lon, lat, z, field, val)
Trilinear interpolation.
subroutine, public qg_interp_bicubic(geom, x, y, gfld2dext, val)
Bicubic horizontal interpolation.
subroutine cubic_tl(a_traj, m1_traj, o_traj, p1_traj, p2_traj, a, m1, o, p1, p2, res)
Cubic interpolation - tangent linear.
subroutine, public lonlat_to_xy(lon, lat, x, y)
Convert lon/lat to x/y.