SCALE-RM
scale_interpolation_nest.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 !-------------------------------------------------------------------------------
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  use scale_prof
22  use scale_debug
25  use scale_index
26  use scale_tracer
27  use scale_const, only: &
28  r_in_m => const_radius
29  !-----------------------------------------------------------------------------
30  implicit none
31  private
32  !-----------------------------------------------------------------------------
33  !
34  !++ Public procedure
35  !
36  public :: intrpnest_setup
40 
41  !-----------------------------------------------------------------------------
42  !
43  !++ Public parameters & variables
44  !
45  !-----------------------------------------------------------------------------
46  !
47  !++ Private procedure
48  !
49  private :: intrpnest_search_nearest_block
50  private :: intrpnest_search_horiz_1points
51  private :: intrpnest_search_horiz_3points
52  private :: intrpnest_search_horiz_4points
53  private :: intrpnest_search_horiz_8points
54  private :: intrpnest_search_horiz_12points
55  private :: intrpnest_search_vert_offline
56  private :: intrpnest_search_vert_online
57  private :: intrpnest_interp_2d_1points
58  private :: intrpnest_interp_3d_1points
59  private :: intrpnest_interp_2d_3points
60  private :: intrpnest_interp_3d_3points
61  private :: intrpnest_interp_2d_4points
62  private :: intrpnest_interp_3d_4points
63  private :: intrpnest_interp_2d_8points
64  private :: intrpnest_interp_3d_8points
65  private :: intrpnest_interp_2d_12points
66  private :: intrpnest_interp_3d_12points
67  private :: intrpnest_haversine
68 
69  !-----------------------------------------------------------------------------
70  abstract interface
71  subroutine intrpnest_intfc_search_h( &
72  hfact, & ! (out)
73  igrd, & ! (out)
74  jgrd, & ! (out)
75  mylat, & ! (in)
76  mylon, & ! (in)
77  inlat, & ! (in)
78  inlon, & ! (in)
79  is, & ! (in)
80  ie, & ! (in)
81  js, & ! (in)
82  je ) ! (in)
83  use scale_precision
84  implicit none
85 
86  real(RP), intent(out) :: hfact(:) ! horizontal interp factor
87  integer, intent(out) :: igrd (:) ! grid points of interp target
88  integer, intent(out) :: jgrd (:) ! grid points of interp target
89  real(RP), intent(in) :: mylat ! latitude data of mine
90  real(RP), intent(in) :: mylon ! longitude data of mine
91  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
92  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
93  integer, intent(in) :: is ! start index for x-direction
94  integer, intent(in) :: ie ! end index for x-direction
95  integer, intent(in) :: js ! start index for y-direction
96  integer, intent(in) :: je ! end index for y-direction
97  end subroutine intrpnest_intfc_search_h
98  end interface
99  procedure(intrpnest_intfc_search_h), pointer :: intrpnest_search_horiz => null()
100  private :: intrpnest_search_horiz
101 
102  !-----------------------------------------------------------------------------
103  abstract interface
104  subroutine intrpnest_intfc_search_v( &
105  vfact, & ! (out)
106  kgrd, & ! (out)
107  ncopy, & ! (out)
108  igrd, & ! (in)
109  jgrd, & ! (in)
110  myhgt, & ! (in)
111  inhgt, & ! (in)
112  iloc, & ! (in)
113  jloc, & ! (in)
114  ks, & ! (in)
115  ke, & ! (in)
116  inka, & ! (in)
117  lndgrd ) ! (in)
118  use scale_precision
119  implicit none
120 
121  real(RP), intent(inout) :: vfact(:,:,:,:,:)! vertical interp factor
122  integer, intent(inout) :: kgrd (:,:,:,:,:)! grid points of interp target
123  integer, intent(out) :: ncopy(:) ! number of daughter's layers below the parent's lowest layer
124  integer, intent(in) :: igrd(:) ! grid points of interp target
125  integer, intent(in) :: jgrd(:) ! grid points of interp target
126  real(RP), intent(in) :: myhgt(:) ! height data of mine
127  real(RP), intent(in) :: inhgt(:,:,:) ! height data of you (input)
128  integer, intent(in) :: iloc ! locator index for x-direction
129  integer, intent(in) :: jloc ! locator index for y-direction
130  integer, intent(in) :: ks ! start index for z-direction
131  integer, intent(in) :: ke ! end index for z-direction
132  integer, intent(in) :: inka ! grid number of you (input)
133  logical, intent(in) :: lndgrd ! flag of land grid
134  end subroutine intrpnest_intfc_search_v
135  end interface
136  procedure(intrpnest_intfc_search_v), pointer :: intrpnest_search_vert => null()
137  private :: intrpnest_search_vert
138 
139  !-----------------------------------------------------------------------------
140  abstract interface
141  subroutine intrpnest_intfc_interp_2d( &
142  intp, & ! (out)
143  ref, & ! (in)
144  hfact, & ! (in)
145  igrd, & ! (in)
146  jgrd, & ! (in)
147  ia, & ! (in)
148  ja ) ! (in)
149  use scale_precision
150  implicit none
151 
152  real(RP), intent(out) :: intp(:,:) ! interpolated data
153  real(RP), intent(in) :: ref (:,:) ! reference data
154  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
155  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
156  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
157  integer, intent(in) :: ia ! grid number of mine
158  integer, intent(in) :: ja ! grid number of mine
159  end subroutine intrpnest_intfc_interp_2d
160  end interface
161  procedure(intrpnest_intfc_interp_2d), pointer :: intrpnest_interp_2d => null()
162  public :: intrpnest_interp_2d
163 
164  !-----------------------------------------------------------------------------
165  abstract interface
166  subroutine intrpnest_intfc_interp_3d( &
167  intp, & ! (out)
168  ref, & ! (in)
169  hfact, & ! (in)
170  vfact, & ! (in)
171  kgrd, & ! (in)
172  igrd, & ! (in)
173  jgrd, & ! (in)
174  ia, & ! (in)
175  ja, & ! (in)
176  ks, & ! (in)
177  ke, & ! (in)
178  logwegt ) ! (in)
179  use scale_precision
180  implicit none
181 
182  real(RP), intent(out) :: intp(:,:,:) ! interpolated data
183  real(RP), intent(in) :: ref (:,:,:) ! reference data
184  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
185  real(RP), intent(in) :: vfact(:,:,:,:,:) ! vertical interp factor
186  integer, intent(in) :: kgrd (:,:,:,:,:) ! grid points of interp target
187  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
188  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
189  integer, intent(in) :: ia ! grid number of mine
190  integer, intent(in) :: ja ! grid number of mine
191  integer, intent(in) :: ks ! start grid number of mine
192  integer, intent(in) :: ke ! end grid number of mine
193  logical, intent(in), optional :: logwegt
194  end subroutine intrpnest_intfc_interp_3d
195  end interface
196  procedure(intrpnest_intfc_interp_3d), pointer :: intrpnest_interp_3d => null()
197  public :: intrpnest_interp_3d
198 
199  !-----------------------------------------------------------------------------
200  !
201  !++ Private parameters & variables
202  !
203  real(RP), private, parameter :: large_number_1 = 9.999e+15_rp
204  real(RP), private, parameter :: large_number_2 = 9.888e+15_rp
205  real(RP), private, parameter :: large_number_3 = 9.777e+15_rp
206  real(RP), private, parameter :: large_number_4 = 9.666e+15_rp
207  real(RP), private, parameter :: large_number_5 = 9.555e+15_rp
208  real(RP), private, parameter :: large_number_6 = 9.444e+15_rp
209  real(RP), private, parameter :: large_number_7 = 9.333e+15_rp
210  real(RP), private, parameter :: large_number_8 = 9.222e+15_rp
211  real(RP), private, parameter :: large_number_9 = 9.111e+15_rp
212  real(RP), private, parameter :: large_number_10 = 9.000e+15_rp
213  real(RP), private, parameter :: large_number_11 = 8.999e+15_rp
214  real(RP), private, parameter :: large_number_12 = 8.888e+15_rp
215 
216  integer, private :: divnum
217  integer, private :: itp_nh
218 
219  !-----------------------------------------------------------------------------
220 contains
221  !-----------------------------------------------------------------------------
223  subroutine intrpnest_setup ( &
224  interp_search_divnum, &
225  NEST_INTERP_LEVEL, &
226  OFFLINE )
227  use scale_process, only: &
229  implicit none
230 
231  integer, intent(in) :: interp_search_divnum
232  integer, intent(in) :: NEST_INTERP_LEVEL
233  logical, intent(in) :: OFFLINE
234 
235  character(7) :: select_type
236  !---------------------------------------------------------------------------
237 
238  if( io_l ) write(io_fid_log,*)
239  if( io_l ) write(io_fid_log,*) '+++ Module[NEST]/Categ[GRID INTERP]'
240 
241  divnum = interp_search_divnum
242 
243  select case ( nest_interp_level )
244  case ( 1 )
245  intrpnest_search_horiz => intrpnest_search_horiz_1points
246  intrpnest_interp_2d => intrpnest_interp_2d_1points
247  intrpnest_interp_3d => intrpnest_interp_3d_1points
248  itp_nh = 1
249 
250  case ( 3 )
251  intrpnest_search_horiz => intrpnest_search_horiz_3points
252  intrpnest_interp_2d => intrpnest_interp_2d_3points
253  intrpnest_interp_3d => intrpnest_interp_3d_3points
254  itp_nh = 3
255 
256  case ( 4 )
257  intrpnest_search_horiz => intrpnest_search_horiz_4points
258  intrpnest_interp_2d => intrpnest_interp_2d_4points
259  intrpnest_interp_3d => intrpnest_interp_3d_4points
260  itp_nh = 4
261 
262  case ( 8 )
263  intrpnest_search_horiz => intrpnest_search_horiz_8points
264  intrpnest_interp_2d => intrpnest_interp_2d_8points
265  intrpnest_interp_3d => intrpnest_interp_3d_8points
266  itp_nh = 8
267 
268  case ( 12 )
269  intrpnest_search_horiz => intrpnest_search_horiz_12points
270  intrpnest_interp_2d => intrpnest_interp_2d_12points
271  intrpnest_interp_3d => intrpnest_interp_3d_12points
272  itp_nh = 12
273 
274  case default
275  write(*,*) 'xxx invarid NEST_INTERP_LEVEL (', nest_interp_level, &
276  ') [setup: nest/interp]'
277  call prc_mpistop
278  end select
279 
280  if ( offline ) then
281  select_type = "offline"
282  intrpnest_search_vert => intrpnest_search_vert_offline
283  else
284  select_type = "online"
285  intrpnest_search_vert => intrpnest_search_vert_online
286  endif
287 
288  if( io_l ) write(io_fid_log,*) '+++ horizontal interpolation with ', &
289  nest_interp_level, " points."
290  if( io_l ) write(io_fid_log,*) '+++ vertical interpolation for ', &
291  trim(select_type)
292 
293  return
294  end subroutine intrpnest_setup
295 
296 
297  !-----------------------------------------------------------------------------
298  ! make interpolation factor using Lat-Lon
299  subroutine intrpnest_interp_fact_latlon( &
300  hfact, & ! (out)
301  igrd, & ! (out)
302  jgrd, & ! (out)
303  mylat, & ! (in)
304  mylon, & ! (in)
305  myia, & ! (in)
306  myja, & ! (in)
307  inlat, & ! (in)
308  inlon, & ! (in)
309  inia, & ! (in)
310  inja ) ! (in)
311  use scale_process, only: &
313  implicit none
314 
315  real(RP), intent(out) :: hfact(:,:,:) ! horizontal interp factor
316  integer, intent(out) :: igrd (:,:,:) ! grid points of interp target
317  integer, intent(out) :: jgrd (:,:,:) ! grid points of interp target
318 
319  real(RP), intent(in) :: mylat(:,:) ! latitude data of mine
320  real(RP), intent(in) :: mylon(:,:) ! longitude data of mine
321  integer, intent(in) :: myIA ! grid number of mine
322  integer, intent(in) :: myJA ! grid number of mine
323 
324  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
325  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
326  integer, intent(in) :: inIA ! grid number of you (input)
327  integer, intent(in) :: inJA ! grid number of you (input)
328 
329  integer :: i, j
330  integer :: is, ie
331  integer :: js, je
332  !---------------------------------------------------------------------------
333 
334  hfact(:,:,:) = 0.0_rp
335 
336  do j = 1, myja
337  do i = 1, myia
338  ! nearest block search
339  call intrpnest_search_nearest_block( is, ie, js, je, &
340  mylat(i,j), mylon(i,j), &
341  inlat(:,:), inlon(:,:), &
342  inia, inja )
343 
344  ! main search
345  call intrpnest_search_horiz( hfact(i,j,:), &
346  igrd(i,j,:), &
347  jgrd(i,j,:), &
348  mylat(i,j), &
349  mylon(i,j), &
350  inlat, &
351  inlon, &
352  is, ie, &
353  js, je )
354  enddo
355  enddo
356 
357  return
358  end subroutine intrpnest_interp_fact_latlon
359 
360 
361  !-----------------------------------------------------------------------------
362  ! make interpolation factor using Lat-Lon and Z-Height information
363  subroutine intrpnest_interp_fact_llz( &
364  hfact, & ! (out)
365  vfact, & ! (out)
366  kgrd, & ! (out)
367  igrd, & ! (out)
368  jgrd, & ! (out)
369  ncopy, & ! (out)
370  myhgt, & ! (in)
371  mylat, & ! (in)
372  mylon, & ! (in)
373  myks, & ! (in)
374  myke, & ! (in)
375  myia, & ! (in)
376  myja, & ! (in)
377  inhgt, & ! (in)
378  inlat, & ! (in)
379  inlon, & ! (in)
380  inka, & ! (in)
381  inia, & ! (in)
382  inja, & ! (in)
383  landgrid ) ! (in)
384  use scale_process, only: &
386  implicit none
387 
388  real(RP), intent(out) :: hfact(:,:,:) ! horizontal interp factor
389  real(RP), intent(out) :: vfact(:,:,:,:,:) ! vertical interp factor
390  integer, intent(out) :: kgrd (:,:,:,:,:) ! grid points of interp target
391  integer, intent(out) :: igrd (:,:,:) ! grid points of interp target
392  integer, intent(out) :: jgrd (:,:,:) ! grid points of interp target
393  integer, intent(out) :: ncopy(:,:,:) ! number of daughter's layers below parent lowest layer
394 
395  real(RP), intent(in) :: myhgt(:,:,:) ! height data of mine
396  real(RP), intent(in) :: mylat(:,:) ! latitude data of mine
397  real(RP), intent(in) :: mylon(:,:) ! longitude data of mine
398  integer, intent(in) :: myKS ! start grid number of mine
399  integer, intent(in) :: myKE ! end grid number of mine
400  integer, intent(in) :: myIA ! grid number of mine
401  integer, intent(in) :: myJA ! grid number of mine
402 
403  real(RP), intent(in) :: inhgt(:,:,:) ! height data of you (input)
404  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
405  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
406  integer, intent(in) :: inKA ! grid number of you (input)
407  integer, intent(in) :: inIA ! grid number of you (input)
408  integer, intent(in) :: inJA ! grid number of you (input)
409 
410  logical, intent(in), optional :: landgrid
411 
412  integer :: i, j
413  integer :: is, ie
414  integer :: js, je
415  logical :: lndgrd
416  !---------------------------------------------------------------------------
417 
418  lndgrd = .false.
419  if ( present(landgrid) ) then
420  if ( landgrid ) then
421  lndgrd = .true.
422  endif
423  endif
424 
425  hfact(:,:,:) = 0.0_rp
426  vfact(:,:,:,:,:) = 0.0_rp
427  ncopy(:,:,:) = 0
428 
429  do j = 1, myja
430  do i = 1, myia
431 
432  ! nearest block search
433  call intrpnest_search_nearest_block( is, ie, js, je, &
434  mylat(i,j), mylon(i,j), &
435  inlat(:,:), inlon(:,:), &
436  inia, inja )
437 
438  ! main search
439  call intrpnest_search_horiz( hfact(i,j,:), &
440  igrd(i,j,:), &
441  jgrd(i,j,:), &
442  mylat(i,j), &
443  mylon(i,j), &
444  inlat, &
445  inlon, &
446  is, ie, &
447  js, je )
448 
449 
450  call intrpnest_search_vert( vfact, &
451  kgrd, &
452  ncopy(i,j,:), &
453  igrd(i,j,:), &
454  jgrd(i,j,:), &
455  myhgt(:,i,j), &
456  inhgt, &
457  i, j, &
458  myks, myke, &
459  inka, &
460  lndgrd )
461 
462  enddo
463  enddo
464 
465  return
466  end subroutine intrpnest_interp_fact_llz
467 
468 
469  !-----------------------------------------------------------------------------
470  ! search of nearest region for speed up of interpolation
471  subroutine intrpnest_search_nearest_block( &
472  is, & ! (out)
473  ie, & ! (out)
474  js, & ! (out)
475  je, & ! (out)
476  mylat, & ! (in)
477  mylon, & ! (in)
478  inlat, & ! (in)
479  inlon, & ! (in)
480  inia, & ! (in)
481  inja ) ! (in)
482  implicit none
483 
484  integer, intent(out) :: is ! start index for x-direction
485  integer, intent(out) :: ie ! end index for x-direction
486  integer, intent(out) :: js ! start index for y-direction
487  integer, intent(out) :: je ! end index for y-direction
488 
489  real(RP), intent(in) :: mylat ! latitude data of mine
490  real(RP), intent(in) :: mylon ! longitude data of mine
491  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
492  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
493  integer, intent(in) :: inIA ! grid number of you (input)
494  integer, intent(in) :: inJA ! grid number of you (input)
495 
496  real(RP) :: distance, dist
497  integer :: ii, jj
498  integer :: iinc, jinc
499  integer :: blk_i, blk_j
500  !---------------------------------------------------------------------------
501 
502  iinc = max( (inia + 1) / divnum, 1 )
503  jinc = max( (inja + 1) / divnum, 1 )
504  dist = large_number_1
505 
506  jj = 1 + (jinc/2)
507  do while (jj <= inja)
508  ii = 1 + (iinc/2)
509  do while (ii <= inia)
510  distance = intrpnest_haversine( mylat, mylon, &
511  inlat(ii,jj), inlon(ii,jj) )
512  if( distance < dist )then
513  dist = distance
514  blk_i = ii
515  blk_j = jj
516  endif
517  ii = ii + iinc
518  enddo
519  jj = jj + jinc
520  enddo
521 
522  ! +- 3 is buffer for 12 points
523  is = blk_i - (iinc/2) - 3
524  if( is < 1 ) is = 1
525  ie = blk_i + (iinc/2) + 3
526  if( ie > inia ) ie = inia
527  js = blk_j - (jinc/2) - 3
528  if( js < 1 ) js = 1
529  je = blk_j + (jinc/2) + 3
530  if( je > inja ) je = inja
531 
532  return
533  end subroutine intrpnest_search_nearest_block
534 
535 
536  !-----------------------------------------------------------------------------
537  ! horizontal search of interpolation points for one-points (nearest-neighbor)
538  subroutine intrpnest_search_horiz_1points( &
539  hfact, & ! (out)
540  igrd, & ! (out)
541  jgrd, & ! (out)
542  mylat, & ! (in)
543  mylon, & ! (in)
544  inlat, & ! (in)
545  inlon, & ! (in)
546  is, & ! (in)
547  ie, & ! (in)
548  js, & ! (in)
549  je ) ! (in)
550  implicit none
551 
552  real(RP), intent(out) :: hfact(:) ! horizontal interp factor
553  integer, intent(out) :: igrd (:) ! grid points of interp target
554  integer, intent(out) :: jgrd (:) ! grid points of interp target
555 
556  real(RP), intent(in) :: mylat ! latitude data of mine
557  real(RP), intent(in) :: mylon ! longitude data of mine
558  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
559  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
560 
561  integer, intent(in) :: is ! start index for x-direction
562  integer, intent(in) :: ie ! end index for x-direction
563  integer, intent(in) :: js ! start index for y-direction
564  integer, intent(in) :: je ! end index for y-direction
565 
566  real(RP) :: distance
567  real(RP) :: dist
568  integer :: ii, jj
569  !---------------------------------------------------------------------------
570 
571  dist = large_number_1
572  igrd(:) = -1
573  jgrd(:) = -1
574 
575  do jj = js, je
576  do ii = is, ie
577  distance = intrpnest_haversine( mylat,mylon,inlat(ii,jj),inlon(ii,jj) )
578  if ( distance <= dist ) then
579  dist = distance; igrd(1) = ii; jgrd(1) = jj
580  endif
581  enddo
582  enddo
583 
584  hfact(1) = 1.0_rp
585 
586  return
587  end subroutine intrpnest_search_horiz_1points
588 
589 
590  !-----------------------------------------------------------------------------
591  ! horizontal search of interpolation points for three-points
592  subroutine intrpnest_search_horiz_3points( &
593  hfact, & ! (out)
594  igrd, & ! (out)
595  jgrd, & ! (out)
596  mylat, & ! (in)
597  mylon, & ! (in)
598  inlat, & ! (in)
599  inlon, & ! (in)
600  is, & ! (in)
601  ie, & ! (in)
602  js, & ! (in)
603  je ) ! (in)
604  use scale_const, only: &
605  eps => const_eps
606  implicit none
607 
608  real(RP), intent(out) :: hfact(:) ! horizontal interp factor
609  integer, intent(out) :: igrd (:) ! grid points of interp target
610  integer, intent(out) :: jgrd (:) ! grid points of interp target
611 
612  real(RP), intent(in) :: mylat ! latitude data of mine
613  real(RP), intent(in) :: mylon ! longitude data of mine
614  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
615  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
616 
617  integer, intent(in) :: is ! start index for x-direction
618  integer, intent(in) :: ie ! end index for x-direction
619  integer, intent(in) :: js ! start index for y-direction
620  integer, intent(in) :: je ! end index for y-direction
621 
622  real(RP) :: distance
623  real(RP) :: denom
624  real(RP) :: dist(3)
625  integer :: ii, jj
626  !---------------------------------------------------------------------------
627 
628  dist(1) = large_number_3
629  dist(2) = large_number_2
630  dist(3) = large_number_1
631  igrd(:) = -1
632  jgrd(:) = -1
633 
634  do jj = js, je
635  do ii = is, ie
636  distance = intrpnest_haversine( mylat,mylon,inlat(ii,jj),inlon(ii,jj) )
637  if ( distance <= dist(1) ) then
638  dist(3) = dist(2); igrd(3) = igrd(2); jgrd(3) = jgrd(2)
639  dist(2) = dist(1); igrd(2) = igrd(1); jgrd(2) = jgrd(1)
640  dist(1) = distance; igrd(1) = ii; jgrd(1) = jj
641  elseif ( dist(1) < distance .AND. distance <= dist(2) ) then
642  dist(3) = dist(2); igrd(3) = igrd(2); jgrd(3) = jgrd(2)
643  dist(2) = distance; igrd(2) = ii; jgrd(2) = jj
644  elseif ( dist(2) < distance .AND. distance <= dist(3) ) then
645  dist(3) = distance; igrd(3) = ii; jgrd(3) = jj
646  endif
647  enddo
648  enddo
649 
650  if ( abs(dist(1)) < eps ) then
651  hfact(1) = 1.0_rp
652  hfact(2) = 0.0_rp
653  hfact(3) = 0.0_rp
654  else
655  denom = 1.0_rp / ( (1.0_rp/dist(1)) + (1.0_rp/dist(2)) + (1.0_rp/dist(3)) )
656  hfact(1) = ( 1.0_rp/dist(1) ) * denom
657  hfact(2) = ( 1.0_rp/dist(2) ) * denom
658  hfact(3) = ( 1.0_rp/dist(3) ) * denom
659  endif
660 
661  return
662  end subroutine intrpnest_search_horiz_3points
663 
664 
665  !-----------------------------------------------------------------------------
666  ! horizontal search of interpolation points for four-points
667  subroutine intrpnest_search_horiz_4points( &
668  hfact, & ! (out)
669  igrd, & ! (out)
670  jgrd, & ! (out)
671  mylat, & ! (in)
672  mylon, & ! (in)
673  inlat, & ! (in)
674  inlon, & ! (in)
675  is, & ! (in)
676  ie, & ! (in)
677  js, & ! (in)
678  je ) ! (in)
679  use scale_const, only: &
680  eps => const_eps
681  implicit none
682 
683  real(RP), intent(out) :: hfact(:) ! horizontal interp factor
684  integer, intent(out) :: igrd (:) ! grid points of interp target
685  integer, intent(out) :: jgrd (:) ! grid points of interp target
686 
687  real(RP), intent(in) :: mylat ! latitude data of mine
688  real(RP), intent(in) :: mylon ! longitude data of mine
689  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
690  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
691 
692  integer, intent(in) :: is ! start index for x-direction
693  integer, intent(in) :: ie ! end index for x-direction
694  integer, intent(in) :: js ! start index for y-direction
695  integer, intent(in) :: je ! end index for y-direction
696 
697  real(RP) :: distance
698  real(RP) :: denom
699  real(RP) :: dist(4)
700  integer :: ii, jj
701  !---------------------------------------------------------------------------
702 
703  dist(1) = large_number_4
704  dist(2) = large_number_3
705  dist(3) = large_number_2
706  dist(4) = large_number_1
707  igrd(:) = -1
708  jgrd(:) = -1
709 
710  do jj = js, je
711  do ii = is, ie
712  distance = intrpnest_haversine( mylat,mylon,inlat(ii,jj),inlon(ii,jj) )
713  if ( distance <= dist(1) ) then
714  dist(4) = dist(3); igrd(4) = igrd(3); jgrd(4) = jgrd(3)
715  dist(3) = dist(2); igrd(3) = igrd(2); jgrd(3) = jgrd(2)
716  dist(2) = dist(1); igrd(2) = igrd(1); jgrd(2) = jgrd(1)
717  dist(1) = distance; igrd(1) = ii; jgrd(1) = jj
718  elseif ( dist(1) < distance .AND. distance <= dist(2) ) then
719  dist(4) = dist(3); igrd(4) = igrd(3); jgrd(4) = jgrd(3)
720  dist(3) = dist(2); igrd(3) = igrd(2); jgrd(3) = jgrd(2)
721  dist(2) = distance; igrd(2) = ii; jgrd(2) = jj
722  elseif ( dist(2) < distance .AND. distance <= dist(3) ) then
723  dist(4) = dist(3); igrd(4) = igrd(3); jgrd(4) = jgrd(3)
724  dist(3) = distance; igrd(3) = ii; jgrd(3) = jj
725  elseif ( dist(3) < distance .AND. distance <= dist(4) ) then
726  dist(4) = distance; igrd(4) = ii; jgrd(4) = jj
727  endif
728  enddo
729  enddo
730 
731  if ( abs(dist(1)) < eps ) then
732  hfact(1) = 1.0_rp
733  hfact(2) = 0.0_rp
734  hfact(3) = 0.0_rp
735  hfact(4) = 0.0_rp
736  else
737  denom = 1.0_rp / ( (1.0_rp/dist(1)) + (1.0_rp/dist(2)) &
738  + (1.0_rp/dist(3)) + (1.0_rp/dist(4)) )
739  hfact(1) = ( 1.0_rp/dist(1) ) * denom
740  hfact(2) = ( 1.0_rp/dist(2) ) * denom
741  hfact(3) = ( 1.0_rp/dist(3) ) * denom
742  hfact(4) = ( 1.0_rp/dist(4) ) * denom
743  endif
744 
745  return
746  end subroutine intrpnest_search_horiz_4points
747 
748 
749  !-----------------------------------------------------------------------------
750  ! horizontal search of interpolation points for eight-points
751  subroutine intrpnest_search_horiz_8points( &
752  hfact, & ! (out)
753  igrd, & ! (out)
754  jgrd, & ! (out)
755  mylat, & ! (in)
756  mylon, & ! (in)
757  inlat, & ! (in)
758  inlon, & ! (in)
759  is, & ! (in)
760  ie, & ! (in)
761  js, & ! (in)
762  je ) ! (in)
763  use scale_const, only: &
764  eps => const_eps
765  implicit none
766 
767  real(RP), intent(out) :: hfact(:) ! horizontal interp factor
768  integer, intent(out) :: igrd (:) ! grid points of interp target
769  integer, intent(out) :: jgrd (:) ! grid points of interp target
770 
771  real(RP), intent(in) :: mylat ! latitude data of mine
772  real(RP), intent(in) :: mylon ! longitude data of mine
773  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
774  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
775 
776  integer, intent(in) :: is ! start index for x-direction
777  integer, intent(in) :: ie ! end index for x-direction
778  integer, intent(in) :: js ! start index for y-direction
779  integer, intent(in) :: je ! end index for y-direction
780 
781  real(RP) :: distance
782  real(RP) :: denom
783  real(RP) :: dist(8)
784  integer :: ii, jj
785  !---------------------------------------------------------------------------
786 
787  dist(1) = large_number_8
788  dist(2) = large_number_7
789  dist(3) = large_number_6
790  dist(4) = large_number_5
791  dist(5) = large_number_4
792  dist(6) = large_number_3
793  dist(7) = large_number_2
794  dist(8) = large_number_1
795  igrd(:) = -1
796  jgrd(:) = -1
797 
798  do jj = js, je
799  do ii = is, ie
800  distance = intrpnest_haversine( mylat,mylon,inlat(ii,jj),inlon(ii,jj) )
801  if ( distance <= dist(1) ) then
802  dist(8) = dist(7); igrd(8) = igrd(7); jgrd(8) = jgrd(7)
803  dist(7) = dist(6); igrd(7) = igrd(6); jgrd(7) = jgrd(6)
804  dist(6) = dist(5); igrd(6) = igrd(5); jgrd(6) = jgrd(5)
805  dist(5) = dist(4); igrd(5) = igrd(4); jgrd(5) = jgrd(4)
806  dist(4) = dist(3); igrd(4) = igrd(3); jgrd(4) = jgrd(3)
807  dist(3) = dist(2); igrd(3) = igrd(2); jgrd(3) = jgrd(2)
808  dist(2) = dist(1); igrd(2) = igrd(1); jgrd(2) = jgrd(1)
809  dist(1) = distance; igrd(1) = ii; jgrd(1) = jj
810  elseif ( dist(1) < distance .AND. distance <= dist(2) ) then
811  dist(8) = dist(7); igrd(8) = igrd(7); jgrd(8) = jgrd(7)
812  dist(7) = dist(6); igrd(7) = igrd(6); jgrd(7) = jgrd(6)
813  dist(6) = dist(5); igrd(6) = igrd(5); jgrd(6) = jgrd(5)
814  dist(5) = dist(4); igrd(5) = igrd(4); jgrd(5) = jgrd(4)
815  dist(4) = dist(3); igrd(4) = igrd(3); jgrd(4) = jgrd(3)
816  dist(3) = dist(2); igrd(3) = igrd(2); jgrd(3) = jgrd(2)
817  dist(2) = distance; igrd(2) = ii; jgrd(2) = jj
818  elseif ( dist(2) < distance .AND. distance <= dist(3) ) then
819  dist(8) = dist(7); igrd(8) = igrd(7); jgrd(8) = jgrd(7)
820  dist(7) = dist(6); igrd(7) = igrd(6); jgrd(7) = jgrd(6)
821  dist(6) = dist(5); igrd(6) = igrd(5); jgrd(6) = jgrd(5)
822  dist(5) = dist(4); igrd(5) = igrd(4); jgrd(5) = jgrd(4)
823  dist(4) = dist(3); igrd(4) = igrd(3); jgrd(4) = jgrd(3)
824  dist(3) = distance; igrd(3) = ii; jgrd(3) = jj
825  elseif ( dist(3) < distance .AND. distance <= dist(4) ) then
826  dist(8) = dist(7); igrd(8) = igrd(7); jgrd(8) = jgrd(7)
827  dist(7) = dist(6); igrd(7) = igrd(6); jgrd(7) = jgrd(6)
828  dist(6) = dist(5); igrd(6) = igrd(5); jgrd(6) = jgrd(5)
829  dist(5) = dist(4); igrd(5) = igrd(4); jgrd(5) = jgrd(4)
830  dist(4) = distance; igrd(4) = ii; jgrd(4) = jj
831  elseif ( dist(4) < distance .AND. distance <= dist(5) ) then
832  dist(8) = dist(7); igrd(8) = igrd(7); jgrd(8) = jgrd(7)
833  dist(7) = dist(6); igrd(7) = igrd(6); jgrd(7) = jgrd(6)
834  dist(6) = dist(5); igrd(6) = igrd(5); jgrd(6) = jgrd(5)
835  dist(5) = distance; igrd(5) = ii; jgrd(5) = jj
836  elseif ( dist(5) < distance .AND. distance <= dist(6) ) then
837  dist(8) = dist(7); igrd(8) = igrd(7); jgrd(8) = jgrd(7)
838  dist(7) = dist(6); igrd(7) = igrd(6); jgrd(7) = jgrd(6)
839  dist(6) = distance; igrd(6) = ii; jgrd(6) = jj
840  elseif ( dist(6) < distance .AND. distance <= dist(7) ) then
841  dist(8) = dist(7); igrd(8) = igrd(7); jgrd(8) = jgrd(7)
842  dist(7) = distance; igrd(7) = ii; jgrd(7) = jj
843  elseif ( dist(7) < distance .AND. distance <= dist(8) ) then
844  dist(8) = distance; igrd(8) = ii; jgrd(8) = jj
845  endif
846  enddo
847  enddo
848 
849  if ( abs(dist(1)) < eps ) then
850  hfact(:) = 0.0_rp
851  hfact(1) = 1.0_rp
852  else
853  denom = 1.0_rp / ( (1.0_rp/dist(1)) + (1.0_rp/dist(2)) &
854  + (1.0_rp/dist(3)) + (1.0_rp/dist(4)) &
855  + (1.0_rp/dist(5)) + (1.0_rp/dist(6)) &
856  + (1.0_rp/dist(7)) + (1.0_rp/dist(8)) )
857  hfact(1) = ( 1.0_rp/dist(1) ) * denom
858  hfact(2) = ( 1.0_rp/dist(2) ) * denom
859  hfact(3) = ( 1.0_rp/dist(3) ) * denom
860  hfact(4) = ( 1.0_rp/dist(4) ) * denom
861  hfact(5) = ( 1.0_rp/dist(5) ) * denom
862  hfact(6) = ( 1.0_rp/dist(6) ) * denom
863  hfact(7) = ( 1.0_rp/dist(7) ) * denom
864  hfact(8) = ( 1.0_rp/dist(8) ) * denom
865  endif
866 
867  return
868  end subroutine intrpnest_search_horiz_8points
869 
870 
871  !-----------------------------------------------------------------------------
872  ! horizontal search of interpolation points for twelve-points
873  subroutine intrpnest_search_horiz_12points( &
874  hfact, & ! (out)
875  igrd, & ! (out)
876  jgrd, & ! (out)
877  mylat, & ! (in)
878  mylon, & ! (in)
879  inlat, & ! (in)
880  inlon, & ! (in)
881  is, & ! (in)
882  ie, & ! (in)
883  js, & ! (in)
884  je ) ! (in)
885  use scale_const, only: &
886  eps => const_eps
887  implicit none
888 
889  real(RP), intent(out) :: hfact(:) ! horizontal interp factor
890  integer, intent(out) :: igrd (:) ! grid points of interp target
891  integer, intent(out) :: jgrd (:) ! grid points of interp target
892 
893  real(RP), intent(in) :: mylat ! latitude data of mine
894  real(RP), intent(in) :: mylon ! longitude data of mine
895  real(RP), intent(in) :: inlat(:,:) ! latitude data of you (input)
896  real(RP), intent(in) :: inlon(:,:) ! longitude data of you (input)
897 
898  integer, intent(in) :: is ! start index for x-direction
899  integer, intent(in) :: ie ! end index for x-direction
900  integer, intent(in) :: js ! start index for y-direction
901  integer, intent(in) :: je ! end index for y-direction
902 
903  real(RP) :: distance
904  real(RP) :: denom
905  real(RP) :: dist(12)
906  integer :: ii, jj
907  !---------------------------------------------------------------------------
908 
909  dist(1 ) = large_number_12
910  dist(2 ) = large_number_11
911  dist(3 ) = large_number_10
912  dist(4 ) = large_number_9
913  dist(5 ) = large_number_8
914  dist(6 ) = large_number_7
915  dist(7 ) = large_number_6
916  dist(8 ) = large_number_5
917  dist(9 ) = large_number_4
918  dist(10) = large_number_3
919  dist(11) = large_number_2
920  dist(12) = large_number_1
921  igrd(:) = -1
922  jgrd(:) = -1
923 
924  do jj = js, je
925  do ii = is, ie
926  distance = intrpnest_haversine( mylat,mylon,inlat(ii,jj),inlon(ii,jj) )
927  if ( distance <= dist(1) ) then
928  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
929  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
930  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
931  dist(9 ) = dist(8 ); igrd(9 ) = igrd(8 ); jgrd(9 ) = jgrd(8 )
932  dist(8 ) = dist(7 ); igrd(8 ) = igrd(7 ); jgrd(8 ) = jgrd(7 )
933  dist(7 ) = dist(6 ); igrd(7 ) = igrd(6 ); jgrd(7 ) = jgrd(6 )
934  dist(6 ) = dist(5 ); igrd(6 ) = igrd(5 ); jgrd(6 ) = jgrd(5 )
935  dist(5 ) = dist(4 ); igrd(5 ) = igrd(4 ); jgrd(5 ) = jgrd(4 )
936  dist(4 ) = dist(3 ); igrd(4 ) = igrd(3 ); jgrd(4 ) = jgrd(3 )
937  dist(3 ) = dist(2 ); igrd(3 ) = igrd(2 ); jgrd(3 ) = jgrd(2 )
938  dist(2 ) = dist(1 ); igrd(2 ) = igrd(1 ); jgrd(2 ) = jgrd(1 )
939  dist(1 ) = distance; igrd(1 ) = ii; jgrd(1 ) = jj
940  elseif ( dist(1) < distance .AND. distance <= dist(2) ) then
941  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
942  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
943  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
944  dist(9 ) = dist(8 ); igrd(9 ) = igrd(8 ); jgrd(9 ) = jgrd(8 )
945  dist(8 ) = dist(7 ); igrd(8 ) = igrd(7 ); jgrd(8 ) = jgrd(7 )
946  dist(7 ) = dist(6 ); igrd(7 ) = igrd(6 ); jgrd(7 ) = jgrd(6 )
947  dist(6 ) = dist(5 ); igrd(6 ) = igrd(5 ); jgrd(6 ) = jgrd(5 )
948  dist(5 ) = dist(4 ); igrd(5 ) = igrd(4 ); jgrd(5 ) = jgrd(4 )
949  dist(4 ) = dist(3 ); igrd(4 ) = igrd(3 ); jgrd(4 ) = jgrd(3 )
950  dist(3 ) = dist(2 ); igrd(3 ) = igrd(2 ); jgrd(3 ) = jgrd(2 )
951  dist(2 ) = distance; igrd(2 ) = ii; jgrd(2 ) = jj
952  elseif ( dist(2) < distance .AND. distance <= dist(3) ) then
953  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
954  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
955  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
956  dist(9 ) = dist(8 ); igrd(9 ) = igrd(8 ); jgrd(9 ) = jgrd(8 )
957  dist(8 ) = dist(7 ); igrd(8 ) = igrd(7 ); jgrd(8 ) = jgrd(7 )
958  dist(7 ) = dist(6 ); igrd(7 ) = igrd(6 ); jgrd(7 ) = jgrd(6 )
959  dist(6 ) = dist(5 ); igrd(6 ) = igrd(5 ); jgrd(6 ) = jgrd(5 )
960  dist(5 ) = dist(4 ); igrd(5 ) = igrd(4 ); jgrd(5 ) = jgrd(4 )
961  dist(4 ) = dist(3 ); igrd(4 ) = igrd(3 ); jgrd(4 ) = jgrd(3 )
962  dist(3 ) = distance; igrd(3 ) = ii; jgrd(3 ) = jj
963  elseif ( dist(3) < distance .AND. distance <= dist(4) ) then
964  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
965  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
966  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
967  dist(9 ) = dist(8 ); igrd(9 ) = igrd(8 ); jgrd(9 ) = jgrd(8 )
968  dist(8 ) = dist(7 ); igrd(8 ) = igrd(7 ); jgrd(8 ) = jgrd(7 )
969  dist(7 ) = dist(6 ); igrd(7 ) = igrd(6 ); jgrd(7 ) = jgrd(6 )
970  dist(6 ) = dist(5 ); igrd(6 ) = igrd(5 ); jgrd(6 ) = jgrd(5 )
971  dist(5 ) = dist(4 ); igrd(5 ) = igrd(4 ); jgrd(5 ) = jgrd(4 )
972  dist(4 ) = distance; igrd(4 ) = ii; jgrd(4 ) = jj
973  elseif ( dist(4) < distance .AND. distance <= dist(5) ) then
974  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
975  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
976  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
977  dist(9 ) = dist(8 ); igrd(9 ) = igrd(8 ); jgrd(9 ) = jgrd(8 )
978  dist(8 ) = dist(7 ); igrd(8 ) = igrd(7 ); jgrd(8 ) = jgrd(7 )
979  dist(7 ) = dist(6 ); igrd(7 ) = igrd(6 ); jgrd(7 ) = jgrd(6 )
980  dist(6 ) = dist(5 ); igrd(6 ) = igrd(5 ); jgrd(6 ) = jgrd(5 )
981  dist(5 ) = distance; igrd(5 ) = ii; jgrd(5 ) = jj
982  elseif ( dist(5) < distance .AND. distance <= dist(6) ) then
983  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
984  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
985  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
986  dist(9 ) = dist(8 ); igrd(9 ) = igrd(8 ); jgrd(9 ) = jgrd(8 )
987  dist(8 ) = dist(7 ); igrd(8 ) = igrd(7 ); jgrd(8 ) = jgrd(7 )
988  dist(7 ) = dist(6 ); igrd(7 ) = igrd(6 ); jgrd(7 ) = jgrd(6 )
989  dist(6 ) = distance; igrd(6 ) = ii; jgrd(6 ) = jj
990  elseif ( dist(6) < distance .AND. distance <= dist(7) ) then
991  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
992  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
993  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
994  dist(9 ) = dist(8 ); igrd(9 ) = igrd(8 ); jgrd(9 ) = jgrd(8 )
995  dist(8 ) = dist(7 ); igrd(8 ) = igrd(7 ); jgrd(8 ) = jgrd(7 )
996  dist(7 ) = distance; igrd(7 ) = ii; jgrd(7 ) = jj
997  elseif ( dist(7) < distance .AND. distance <= dist(8) ) then
998  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
999  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
1000  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
1001  dist(9 ) = dist(8 ); igrd(9 ) = igrd(8 ); jgrd(9 ) = jgrd(8 )
1002  dist(8 ) = distance; igrd(8 ) = ii; jgrd(8 ) = jj
1003  elseif ( dist(8) < distance .AND. distance <= dist(9) ) then
1004  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
1005  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
1006  dist(10) = dist(9 ); igrd(10) = igrd(9 ); jgrd(10) = jgrd(9 )
1007  dist(9 ) = distance; igrd(9 ) = ii; jgrd(9 ) = jj
1008  elseif ( dist(9) < distance .AND. distance <= dist(10) ) then
1009  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
1010  dist(11) = dist(10); igrd(11) = igrd(10); jgrd(11) = jgrd(10)
1011  dist(10) = distance; igrd(10) = ii; jgrd(10) = jj
1012  elseif ( dist(10) < distance .AND. distance <= dist(11) ) then
1013  dist(12) = dist(11); igrd(12) = igrd(11); jgrd(12) = jgrd(11)
1014  dist(11) = distance; igrd(11) = ii; jgrd(11) = jj
1015  elseif ( dist(11) < distance .AND. distance <= dist(12) ) then
1016  dist(12) = distance; igrd(12) = ii; jgrd(12) = jj
1017  endif
1018  enddo
1019  enddo
1020 
1021  if ( abs(dist(1)) < eps ) then
1022  hfact(:) = 0.0_rp
1023  hfact(1) = 1.0_rp
1024  else
1025  denom = 1.0_rp / ( (1.0_rp/dist(1 )) + (1.0_rp/dist(2 )) &
1026  + (1.0_rp/dist(3 )) + (1.0_rp/dist(4 )) &
1027  + (1.0_rp/dist(5 )) + (1.0_rp/dist(6 )) &
1028  + (1.0_rp/dist(7 )) + (1.0_rp/dist(8 )) &
1029  + (1.0_rp/dist(9 )) + (1.0_rp/dist(10)) &
1030  + (1.0_rp/dist(11)) + (1.0_rp/dist(12)) )
1031  hfact(1 ) = ( 1.0_rp/dist(1 ) ) * denom
1032  hfact(2 ) = ( 1.0_rp/dist(2 ) ) * denom
1033  hfact(3 ) = ( 1.0_rp/dist(3 ) ) * denom
1034  hfact(4 ) = ( 1.0_rp/dist(4 ) ) * denom
1035  hfact(5 ) = ( 1.0_rp/dist(5 ) ) * denom
1036  hfact(6 ) = ( 1.0_rp/dist(6 ) ) * denom
1037  hfact(7 ) = ( 1.0_rp/dist(7 ) ) * denom
1038  hfact(8 ) = ( 1.0_rp/dist(8 ) ) * denom
1039  hfact(9 ) = ( 1.0_rp/dist(9 ) ) * denom
1040  hfact(10) = ( 1.0_rp/dist(10) ) * denom
1041  hfact(11) = ( 1.0_rp/dist(11) ) * denom
1042  hfact(12) = ( 1.0_rp/dist(12) ) * denom
1043  endif
1044 
1045  return
1046  end subroutine intrpnest_search_horiz_12points
1047 
1048 
1049  !-----------------------------------------------------------------------------
1050  ! vertical search of interpolation points for two-points (online)
1051  subroutine intrpnest_search_vert_online( &
1052  vfact, & ! (inout)
1053  kgrd, & ! (inout)
1054  ncopy, & ! (out)
1055  igrd, & ! (in)
1056  jgrd, & ! (in)
1057  myhgt, & ! (in)
1058  inhgt, & ! (in)
1059  iloc, & ! (in)
1060  jloc, & ! (in)
1061  ks, & ! (in)
1062  ke, & ! (in)
1063  inka, & ! (in)
1064  lndgrd ) ! (in)
1065  use scale_process, only: &
1066  prc_mpistop
1067  implicit none
1068 
1069  real(RP), intent(inout) :: vfact(:,:,:,:,:) ! vertical interp factor
1070  integer, intent(inout) :: kgrd (:,:,:,:,:) ! grid points of interp target
1071  integer, intent(out) :: ncopy(:) ! number of daughter's layers below inKS
1072 
1073  integer, intent(in) :: igrd(:) ! grid points of interp target
1074  integer, intent(in) :: jgrd(:) ! grid points of interp target
1075  real(RP), intent(in) :: myhgt(:) ! height data of mine
1076  real(RP), intent(in) :: inhgt(:,:,:) ! height data of you (input)
1077  integer, intent(in) :: iloc ! locator index for x-direction
1078  integer, intent(in) :: jloc ! locator index for y-direction
1079  integer, intent(in) :: ks ! start index for z-direction
1080  integer, intent(in) :: ke ! end index for z-direction
1081  integer, intent(in) :: inKA ! grid number of you (input)
1082  logical, intent(in) :: lndgrd ! flag of land grid
1083 
1084  real(RP) :: dist(2)
1085  logical :: dflag ! flag: data found or not
1086  integer :: ii, jj, idx
1087  integer :: k, kk
1088  integer :: inKS, inKE
1089  logical :: copy
1090  !---------------------------------------------------------------------------
1091 
1092  if ( lndgrd ) then
1093  write(*,*) 'xxx internal error [interporation: nest/interp]'
1094  write(*,*) ' land grid is not araviable in online'
1095  call prc_mpistop
1096  endif
1097 
1098  inks = 1 + khalo
1099  inke = inka - khalo
1100 
1101  do idx = 1, itp_nh
1102  ii = igrd(idx)
1103  jj = jgrd(idx)
1104 
1105  do k = ks, ke
1106  dist(1) = large_number_2
1107  dist(2) = large_number_1
1108  copy = .false.
1109 
1110  kgrd(k,iloc,jloc,idx,:) = -1
1111 
1112  dflag = .false.
1113 
1114  if( myhgt(k) < inhgt(inks,ii,jj) ) then
1115  copy = .true.
1116  ncopy(idx) = ncopy(idx) + 1
1117 
1118  kgrd(k,iloc,jloc,idx,:) = inks
1119 
1120  vfact(k,iloc,jloc,idx,1) = 1.0_rp
1121  vfact(k,iloc,jloc,idx,2) = 0.0_rp
1122 
1123  dflag = .true.
1124  else
1125  do kk = inks, inke
1126  dist(1) = myhgt(k) - inhgt(kk ,ii,jj)
1127  dist(2) = myhgt(k) - inhgt(kk+1,ii,jj)
1128 
1129  if( dist(1) >= 0.0_rp .AND. dist(2) < 0.0_rp ) then
1130  kgrd(k,iloc,jloc,idx,1) = kk
1131  kgrd(k,iloc,jloc,idx,2) = kk+1
1132 
1133  vfact(k,iloc,jloc,idx,1) = abs(dist(2)) / ( abs(dist(1)) + abs(dist(2)) )
1134  vfact(k,iloc,jloc,idx,2) = abs(dist(1)) / ( abs(dist(1)) + abs(dist(2)) )
1135 
1136  dflag = .true.
1137 
1138  exit ! loop end
1139  endif
1140  enddo
1141  endif
1142 
1143  if( .NOT. dflag ) then
1144  write(*,*) 'xxx internal error [INTRPNEST_search_vert_online]'
1145  write(*,*) 'xxx data for interpolation was not found.'
1146  write(*,*) 'xxx iloc=',iloc,' jloc=',jloc,' k=',k,' idx=',idx
1147  call prc_mpistop
1148  endif
1149  enddo
1150 
1151  enddo
1152 
1153  return
1154  end subroutine intrpnest_search_vert_online
1155 
1156  !-----------------------------------------------------------------------------
1157  ! vertical search of interpolation points for two-points (offline)
1158  subroutine intrpnest_search_vert_offline( &
1159  vfact, & ! (inout)
1160  kgrd, & ! (inout)
1161  ncopy, & ! (out)
1162  igrd, & ! (in)
1163  jgrd, & ! (in)
1164  myhgt, & ! (in)
1165  inhgt, & ! (in)
1166  iloc, & ! (in)
1167  jloc, & ! (in)
1168  ks, & ! (in)
1169  ke, & ! (in)
1170  inka, & ! (in)
1171  lndgrd ) ! (in)
1172  use scale_const, only: &
1173  eps => const_eps
1174  use scale_process, only: &
1175  prc_mpistop
1176  implicit none
1177 
1178  real(RP), intent(inout) :: vfact(:,:,:,:,:) ! vertical interp factor
1179  integer, intent(inout) :: kgrd (:,:,:,:,:) ! grid points of interp target
1180  integer, intent(out) :: ncopy(:) ! number of daughter's layers below inKS
1181 
1182  integer, intent(in) :: igrd(:) ! grid points of interp target
1183  integer, intent(in) :: jgrd(:) ! grid points of interp target
1184  real(RP), intent(in) :: myhgt(:) ! height data of mine
1185  real(RP), intent(in) :: inhgt(:,:,:) ! height data of you (input)
1186  integer, intent(in) :: iloc ! locator index for x-direction
1187  integer, intent(in) :: jloc ! locator index for y-direction
1188  integer, intent(in) :: ks ! start index for z-direction
1189  integer, intent(in) :: ke ! end index for z-direction
1190  integer, intent(in) :: inKA ! grid number of you (input)
1191  logical, intent(in) :: lndgrd ! flag of land grid
1192 
1193  real(RP) :: denom
1194  real(RP) :: dist(2)
1195  logical :: dflag ! flag: data found or not
1196  integer :: ii, jj
1197  integer :: k, kk, kks, kke
1198  integer :: idx
1199  !---------------------------------------------------------------------------
1200 
1201  if ( lndgrd ) then
1202  kks = 1; kke = lkmax
1203  else
1204  kks = ks; kke = ke
1205  endif
1206 
1207  ncopy = 0 ! dummy
1208  do idx = 1, itp_nh
1209  ii = igrd(idx)
1210  jj = jgrd(idx)
1211 
1212  do k = kks, kke
1213  dist(1) = large_number_2
1214  dist(2) = large_number_1
1215  kgrd(k,iloc,jloc,idx,:) = -1
1216  dflag = .false.
1217 
1218  if( myhgt(k) < inhgt(1,ii,jj) )then ! copy
1219  kgrd(k,iloc,jloc,idx,:) = 1
1220  vfact(k,iloc,jloc,idx,1) = 1.0_rp
1221  vfact(k,iloc,jloc,idx,2) = 0.0_rp
1222  dflag = .true.
1223  else if( abs(inhgt(inka,ii,jj)-myhgt(k))<eps )then
1224  kgrd(k,iloc,jloc,idx,:) = inka
1225  vfact(k,iloc,jloc,idx,1) = 1.0_rp
1226  vfact(k,iloc,jloc,idx,2) = 0.0_rp
1227  dflag = .true.
1228  else if( inhgt(inka,ii,jj) < myhgt(k) )then
1229  if ( lndgrd ) then ! copy
1230  kgrd(k,iloc,jloc,idx,:) = inka
1231  vfact(k,iloc,jloc,idx,1) = 1.0_rp
1232  vfact(k,iloc,jloc,idx,2) = 0.0_rp
1233  dflag = .true.
1234  else
1235  write(*,*) 'xxx internal error [INTRPNEST_search_vert_offline]'
1236  write(*,*) 'xxx data level is beyond parent data'
1237  write(*,*) 'in',ii,jj,inka,inhgt(inka,ii,jj),'my',iloc,jloc,k,myhgt(k)
1238  call prc_mpistop
1239  end if
1240  else
1241 
1242  do kk = 1, inka-1
1243  if( (inhgt(kk,ii,jj)<=myhgt(k)).AND.(myhgt(k)<inhgt(kk+1,ii,jj)) )then
1244  kgrd(k,iloc,jloc,idx,1) = kk
1245  kgrd(k,iloc,jloc,idx,2) = kk+1
1246  dist(1) = abs( myhgt(k) - inhgt(kk,ii,jj) )
1247  dist(2) = abs( myhgt(k) - inhgt(kk+1,ii,jj) )
1248  dflag = .true.
1249  if ( abs(dist(1))<eps )then
1250  vfact(k,iloc,jloc,idx,1) = 1.0_rp
1251  vfact(k,iloc,jloc,idx,2) = 0.0_rp
1252  else
1253  denom = 1.0_rp / ( (1.0_rp/dist(1)) + (1.0_rp/dist(2)) )
1254  vfact(k,iloc,jloc,idx,1) = ( 1.0_rp/dist(1) ) * denom
1255  vfact(k,iloc,jloc,idx,2) = ( 1.0_rp/dist(2) ) * denom
1256  endif
1257  exit
1258  endif
1259  enddo
1260  endif
1261 
1262  if( .NOT. dflag )then
1263  write(*,*) 'xxx internal error [INTRPNEST_search_vert_offline]'
1264  write(*,*) 'xxx data for interpolation was not found.'
1265  write(*,*) 'xxx iloc=',iloc,' jloc=',jloc,' k=',k,' idx=',idx
1266  call prc_mpistop
1267  endif
1268  enddo
1269 
1270  enddo
1271 
1272  return
1273  end subroutine intrpnest_search_vert_offline
1274 
1275  !-----------------------------------------------------------------------------
1276  ! interpolation using one-points for 2D data (nearest-neighbor)
1277  subroutine intrpnest_interp_2d_1points( &
1278  intp, & ! (out)
1279  ref, & ! (in)
1280  hfact, & ! (in)
1281  igrd, & ! (in)
1282  jgrd, & ! (in)
1283  ia, & ! (in)
1284  ja ) ! (in)
1285  implicit none
1286 
1287  real(RP), intent(out) :: intp(:,:) ! interpolated data
1288 
1289  real(RP), intent(in) :: ref (:,:) ! reference data
1290  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1291  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1292  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1293  integer, intent(in) :: ia ! grid number of mine
1294  integer, intent(in) :: ja ! grid number of mine
1295 
1296  integer :: i, j
1297  !---------------------------------------------------------------------------
1298 
1299 !OCL PREFETCH
1300  do j = 1, ja
1301  do i = 1, ia
1302  intp(i,j) = ref(igrd(i,j,1),jgrd(i,j,1)) * hfact(i,j,1)
1303  end do
1304  end do
1305 
1306  return
1307  end subroutine intrpnest_interp_2d_1points
1308 
1309 
1310  !-----------------------------------------------------------------------------
1311  ! interpolation using one-points for 3D data (nearest-neighbor)
1312  subroutine intrpnest_interp_3d_1points( &
1313  intp, & ! (out)
1314  ref, & ! (in)
1315  hfact, & ! (in)
1316  vfact, & ! (in)
1317  kgrd, & ! (in)
1318  igrd, & ! (in)
1319  jgrd, & ! (in)
1320  ia, & ! (in)
1321  ja, & ! (in)
1322  ks, & ! (in)
1323  ke, & ! (in)
1324  logwegt ) ! (in)
1325  implicit none
1326 
1327  real(RP), intent(out) :: intp(:,:,:) ! interpolated data
1328 
1329  real(RP), intent(in) :: ref (:,:,:) ! reference data
1330  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1331  real(RP), intent(in) :: vfact(:,:,:,:,:) ! vertical interp factor
1332  integer, intent(in) :: kgrd (:,:,:,:,:) ! grid points of interp target
1333  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1334  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1335  integer, intent(in) :: ia ! grid number of mine
1336  integer, intent(in) :: ja ! grid number of mine
1337  integer, intent(in) :: ks ! start grid number of mine
1338  integer, intent(in) :: ke ! end grid number of mine
1339 
1340  logical, intent(in), optional :: logwegt
1341 
1342  integer :: i, j, k
1343  logical :: logarithmic
1344  !---------------------------------------------------------------------------
1345 
1346  logarithmic = .false.
1347  if ( present(logwegt) ) then
1348  if ( logwegt ) then
1349  logarithmic = .true.
1350  endif
1351  endif
1352 
1353  ! linear interpolation
1354 !OCL PREFETCH
1355  do j = 1, ja
1356  do i = 1, ia
1357  do k = ks, ke
1358  intp(k,i,j) = ref(kgrd(k,i,j,1,1),igrd(i,j,1),jgrd(i,j,1)) &
1359  * hfact(i,j,1) * vfact(k,i,j,1,1) &
1360  + ref(kgrd(k,i,j,1,2),igrd(i,j,1),jgrd(i,j,1)) &
1361  * hfact(i,j,1) * vfact(k,i,j,1,2)
1362  end do
1363  end do
1364  end do
1365 
1366  ! logarithmic weighting (for pres, dens)
1367  if ( logarithmic ) then
1368  do j = 1, ja
1369  do i = 1, ia
1370  do k = ks, ke
1371  intp(k,i,j) = exp( intp(k,i,j) )
1372  end do
1373  end do
1374  end do
1375  endif
1376 
1377  return
1378  end subroutine intrpnest_interp_3d_1points
1379 
1380 
1381  !-----------------------------------------------------------------------------
1382  ! interpolation using three-points for 2D data
1383  subroutine intrpnest_interp_2d_3points( &
1384  intp, & ! (out)
1385  ref, & ! (in)
1386  hfact, & ! (in)
1387  igrd, & ! (in)
1388  jgrd, & ! (in)
1389  ia, & ! (in)
1390  ja ) ! (in)
1391  implicit none
1392 
1393  real(RP), intent(out) :: intp(:,:) ! interpolated data
1394 
1395  real(RP), intent(in) :: ref (:,:) ! reference data
1396  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1397  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1398  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1399  integer, intent(in) :: ia ! grid number of mine
1400  integer, intent(in) :: ja ! grid number of mine
1401 
1402  integer :: i, j
1403  !---------------------------------------------------------------------------
1404 
1405 !OCL PREFETCH
1406  do j = 1, ja
1407  do i = 1, ia
1408  intp(i,j) = ref(igrd(i,j,1),jgrd(i,j,1)) * hfact(i,j,1) &
1409  + ref(igrd(i,j,2),jgrd(i,j,2)) * hfact(i,j,2) &
1410  + ref(igrd(i,j,3),jgrd(i,j,3)) * hfact(i,j,3)
1411  end do
1412  end do
1413 
1414  return
1415  end subroutine intrpnest_interp_2d_3points
1416 
1417 
1418  !-----------------------------------------------------------------------------
1419  ! interpolation using three-points for 3D data
1420  subroutine intrpnest_interp_3d_3points( &
1421  intp, & ! (out)
1422  ref, & ! (in)
1423  hfact, & ! (in)
1424  vfact, & ! (in)
1425  kgrd, & ! (in)
1426  igrd, & ! (in)
1427  jgrd, & ! (in)
1428  ia, & ! (in)
1429  ja, & ! (in)
1430  ks, & ! (in)
1431  ke, & ! (in)
1432  logwegt ) ! (in)
1433  implicit none
1434 
1435  real(RP), intent(out) :: intp(:,:,:) ! interpolated data
1436 
1437  real(RP), intent(in) :: ref (:,:,:) ! reference data
1438  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1439  real(RP), intent(in) :: vfact(:,:,:,:,:) ! vertical interp factor
1440  integer, intent(in) :: kgrd (:,:,:,:,:) ! grid points of interp target
1441  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1442  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1443  integer, intent(in) :: ia ! grid number of mine
1444  integer, intent(in) :: ja ! grid number of mine
1445  integer, intent(in) :: ks ! start grid number of mine
1446  integer, intent(in) :: ke ! end grid number of mine
1447 
1448  logical, intent(in), optional :: logwegt
1449 
1450  integer :: i, j, k
1451  logical :: logarithmic
1452  !---------------------------------------------------------------------------
1453 
1454  logarithmic = .false.
1455  if ( present(logwegt) ) then
1456  if ( logwegt ) then
1457  logarithmic = .true.
1458  endif
1459  endif
1460 
1461  ! linear interpolation
1462 !OCL PREFETCH
1463  do j = 1, ja
1464  do i = 1, ia
1465  do k = ks, ke
1466  intp(k,i,j) = ref(kgrd(k,i,j,1,1),igrd(i,j,1),jgrd(i,j,1)) &
1467  * hfact(i,j,1) * vfact(k,i,j,1,1) &
1468  + ref(kgrd(k,i,j,2,1),igrd(i,j,2),jgrd(i,j,2)) &
1469  * hfact(i,j,2) * vfact(k,i,j,2,1) &
1470  + ref(kgrd(k,i,j,3,1),igrd(i,j,3),jgrd(i,j,3)) &
1471  * hfact(i,j,3) * vfact(k,i,j,3,1) &
1472  + ref(kgrd(k,i,j,1,2),igrd(i,j,1),jgrd(i,j,1)) &
1473  * hfact(i,j,1) * vfact(k,i,j,1,2) &
1474  + ref(kgrd(k,i,j,2,2),igrd(i,j,2),jgrd(i,j,2)) &
1475  * hfact(i,j,2) * vfact(k,i,j,2,2) &
1476  + ref(kgrd(k,i,j,3,2),igrd(i,j,3),jgrd(i,j,3)) &
1477  * hfact(i,j,3) * vfact(k,i,j,3,2)
1478  end do
1479  end do
1480  end do
1481 
1482  ! logarithmic weighting (for pres, dens)
1483  if ( logarithmic ) then
1484  do j = 1, ja
1485  do i = 1, ia
1486  do k = ks, ke
1487  intp(k,i,j) = exp( intp(k,i,j) )
1488  end do
1489  end do
1490  end do
1491  endif
1492 
1493  return
1494  end subroutine intrpnest_interp_3d_3points
1495 
1496 
1497  !-----------------------------------------------------------------------------
1498  ! interpolation using four-points for 2D data
1499  subroutine intrpnest_interp_2d_4points( &
1500  intp, & ! (out)
1501  ref, & ! (in)
1502  hfact, & ! (in)
1503  igrd, & ! (in)
1504  jgrd, & ! (in)
1505  ia, & ! (in)
1506  ja ) ! (in)
1507  implicit none
1508 
1509  real(RP), intent(out) :: intp(:,:) ! interpolated data
1510 
1511  real(RP), intent(in) :: ref (:,:) ! reference data
1512  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1513  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1514  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1515  integer, intent(in) :: ia ! grid number of mine
1516  integer, intent(in) :: ja ! grid number of mine
1517 
1518  integer :: i, j
1519  !---------------------------------------------------------------------------
1520 
1521 !OCL PREFETCH
1522  do j = 1, ja
1523  do i = 1, ia
1524  intp(i,j) = ref(igrd(i,j,1),jgrd(i,j,1)) * hfact(i,j,1) &
1525  + ref(igrd(i,j,2),jgrd(i,j,2)) * hfact(i,j,2) &
1526  + ref(igrd(i,j,3),jgrd(i,j,3)) * hfact(i,j,3) &
1527  + ref(igrd(i,j,4),jgrd(i,j,4)) * hfact(i,j,4)
1528  end do
1529  end do
1530 
1531  return
1532  end subroutine intrpnest_interp_2d_4points
1533 
1534 
1535  !-----------------------------------------------------------------------------
1536  ! interpolation using four-points for 3D data
1537  subroutine intrpnest_interp_3d_4points( &
1538  intp, & ! (out)
1539  ref, & ! (in)
1540  hfact, & ! (in)
1541  vfact, & ! (in)
1542  kgrd, & ! (in)
1543  igrd, & ! (in)
1544  jgrd, & ! (in)
1545  ia, & ! (in)
1546  ja, & ! (in)
1547  ks, & ! (in)
1548  ke, & ! (in)
1549  logwegt ) ! (in)
1550  implicit none
1551 
1552  real(RP), intent(out) :: intp(:,:,:) ! interpolated data
1553 
1554  real(RP), intent(in) :: ref (:,:,:) ! reference data
1555  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1556  real(RP), intent(in) :: vfact(:,:,:,:,:) ! vertical interp factor
1557  integer, intent(in) :: kgrd (:,:,:,:,:) ! grid points of interp target
1558  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1559  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1560  integer, intent(in) :: ia ! grid number of mine
1561  integer, intent(in) :: ja ! grid number of mine
1562  integer, intent(in) :: ks ! start grid number of mine
1563  integer, intent(in) :: ke ! end grid number of mine
1564 
1565  logical, intent(in), optional :: logwegt
1566 
1567  integer :: i, j, k
1568  logical :: logarithmic
1569  !---------------------------------------------------------------------------
1570 
1571  logarithmic = .false.
1572  if ( present(logwegt) ) then
1573  if ( logwegt ) then
1574  logarithmic = .true.
1575  endif
1576  endif
1577 
1578  ! linear interpolation
1579 !OCL PREFETCH
1580  do j = 1, ja
1581  do i = 1, ia
1582  do k = ks, ke
1583  intp(k,i,j) = ref(kgrd(k,i,j,1,1),igrd(i,j,1),jgrd(i,j,1)) &
1584  * hfact(i,j,1) * vfact(k,i,j,1,1) &
1585  + ref(kgrd(k,i,j,2,1),igrd(i,j,2),jgrd(i,j,2)) &
1586  * hfact(i,j,2) * vfact(k,i,j,2,1) &
1587  + ref(kgrd(k,i,j,3,1),igrd(i,j,3),jgrd(i,j,3)) &
1588  * hfact(i,j,3) * vfact(k,i,j,3,1) &
1589  + ref(kgrd(k,i,j,4,1),igrd(i,j,4),jgrd(i,j,4)) &
1590  * hfact(i,j,4) * vfact(k,i,j,4,1) &
1591  + ref(kgrd(k,i,j,1,2),igrd(i,j,1),jgrd(i,j,1)) &
1592  * hfact(i,j,1) * vfact(k,i,j,1,2) &
1593  + ref(kgrd(k,i,j,2,2),igrd(i,j,2),jgrd(i,j,2)) &
1594  * hfact(i,j,2) * vfact(k,i,j,2,2) &
1595  + ref(kgrd(k,i,j,3,2),igrd(i,j,3),jgrd(i,j,3)) &
1596  * hfact(i,j,3) * vfact(k,i,j,3,2) &
1597  + ref(kgrd(k,i,j,4,2),igrd(i,j,4),jgrd(i,j,4)) &
1598  * hfact(i,j,4) * vfact(k,i,j,4,2)
1599  end do
1600  end do
1601  end do
1602 
1603  ! logarithmic weighting (for pres, dens)
1604  if ( logarithmic ) then
1605  do j = 1, ja
1606  do i = 1, ia
1607  do k = ks, ke
1608  intp(k,i,j) = exp( intp(k,i,j) )
1609  end do
1610  end do
1611  end do
1612  endif
1613 
1614  return
1615  end subroutine intrpnest_interp_3d_4points
1616 
1617 
1618  !-----------------------------------------------------------------------------
1619  ! interpolation using eight-points for 2D data
1620  subroutine intrpnest_interp_2d_8points( &
1621  intp, & ! (out)
1622  ref, & ! (in)
1623  hfact, & ! (in)
1624  igrd, & ! (in)
1625  jgrd, & ! (in)
1626  ia, & ! (in)
1627  ja ) ! (in)
1628  implicit none
1629 
1630  real(RP), intent(out) :: intp(:,:) ! interpolated data
1631 
1632  real(RP), intent(in) :: ref (:,:) ! reference data
1633  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1634  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1635  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1636  integer, intent(in) :: ia ! grid number of mine
1637  integer, intent(in) :: ja ! grid number of mine
1638 
1639  integer :: i, j
1640  !---------------------------------------------------------------------------
1641 
1642 !OCL PREFETCH
1643  do j = 1, ja
1644  do i = 1, ia
1645  intp(i,j) = ref(igrd(i,j,1),jgrd(i,j,1)) * hfact(i,j,1) &
1646  + ref(igrd(i,j,2),jgrd(i,j,2)) * hfact(i,j,2) &
1647  + ref(igrd(i,j,3),jgrd(i,j,3)) * hfact(i,j,3) &
1648  + ref(igrd(i,j,4),jgrd(i,j,4)) * hfact(i,j,4) &
1649  + ref(igrd(i,j,5),jgrd(i,j,5)) * hfact(i,j,5) &
1650  + ref(igrd(i,j,6),jgrd(i,j,6)) * hfact(i,j,6) &
1651  + ref(igrd(i,j,7),jgrd(i,j,7)) * hfact(i,j,7) &
1652  + ref(igrd(i,j,8),jgrd(i,j,8)) * hfact(i,j,8)
1653  end do
1654  end do
1655 
1656  return
1657  end subroutine intrpnest_interp_2d_8points
1658 
1659 
1660  !-----------------------------------------------------------------------------
1661  ! interpolation using eight-points for 3D data
1662  subroutine intrpnest_interp_3d_8points( &
1663  intp, & ! (out)
1664  ref, & ! (in)
1665  hfact, & ! (in)
1666  vfact, & ! (in)
1667  kgrd, & ! (in)
1668  igrd, & ! (in)
1669  jgrd, & ! (in)
1670  ia, & ! (in)
1671  ja, & ! (in)
1672  ks, & ! (in)
1673  ke, & ! (in)
1674  logwegt ) ! (in)
1675  implicit none
1676 
1677  real(RP), intent(out) :: intp(:,:,:) ! interpolated data
1678 
1679  real(RP), intent(in) :: ref (:,:,:) ! reference data
1680  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1681  real(RP), intent(in) :: vfact(:,:,:,:,:) ! vertical interp factor
1682  integer, intent(in) :: kgrd (:,:,:,:,:) ! grid points of interp target
1683  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1684  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1685  integer, intent(in) :: ia ! grid number of mine
1686  integer, intent(in) :: ja ! grid number of mine
1687  integer, intent(in) :: ks ! start grid number of mine
1688  integer, intent(in) :: ke ! end grid number of mine
1689 
1690  logical, intent(in), optional :: logwegt
1691 
1692  integer :: i, j, k
1693  logical :: logarithmic
1694  !---------------------------------------------------------------------------
1695 
1696  logarithmic = .false.
1697  if ( present(logwegt) ) then
1698  if ( logwegt ) then
1699  logarithmic = .true.
1700  endif
1701  endif
1702 
1703  ! linear interpolation
1704 !OCL PREFETCH
1705  do j = 1, ja
1706  do i = 1, ia
1707  do k = ks, ke
1708  intp(k,i,j) = ref(kgrd(k,i,j,1,1),igrd(i,j,1),jgrd(i,j,1)) &
1709  * hfact(i,j,1) * vfact(k,i,j,1,1) &
1710  + ref(kgrd(k,i,j,2,1),igrd(i,j,2),jgrd(i,j,2)) &
1711  * hfact(i,j,2) * vfact(k,i,j,2,1) &
1712  + ref(kgrd(k,i,j,3,1),igrd(i,j,3),jgrd(i,j,3)) &
1713  * hfact(i,j,3) * vfact(k,i,j,3,1) &
1714  + ref(kgrd(k,i,j,4,1),igrd(i,j,4),jgrd(i,j,4)) &
1715  * hfact(i,j,4) * vfact(k,i,j,4,1) &
1716  + ref(kgrd(k,i,j,5,1),igrd(i,j,5),jgrd(i,j,5)) &
1717  * hfact(i,j,5) * vfact(k,i,j,5,1) &
1718  + ref(kgrd(k,i,j,6,1),igrd(i,j,6),jgrd(i,j,6)) &
1719  * hfact(i,j,6) * vfact(k,i,j,6,1) &
1720  + ref(kgrd(k,i,j,7,1),igrd(i,j,7),jgrd(i,j,7)) &
1721  * hfact(i,j,7) * vfact(k,i,j,7,1) &
1722  + ref(kgrd(k,i,j,8,1),igrd(i,j,8),jgrd(i,j,8)) &
1723  * hfact(i,j,8) * vfact(k,i,j,8,1) &
1724  + ref(kgrd(k,i,j,1,2),igrd(i,j,1),jgrd(i,j,1)) &
1725  * hfact(i,j,1) * vfact(k,i,j,1,2) &
1726  + ref(kgrd(k,i,j,2,2),igrd(i,j,2),jgrd(i,j,2)) &
1727  * hfact(i,j,2) * vfact(k,i,j,2,2) &
1728  + ref(kgrd(k,i,j,3,2),igrd(i,j,3),jgrd(i,j,3)) &
1729  * hfact(i,j,3) * vfact(k,i,j,3,2) &
1730  + ref(kgrd(k,i,j,4,2),igrd(i,j,3),jgrd(i,j,3)) &
1731  * hfact(i,j,4) * vfact(k,i,j,4,2) &
1732  + ref(kgrd(k,i,j,5,2),igrd(i,j,3),jgrd(i,j,3)) &
1733  * hfact(i,j,5) * vfact(k,i,j,5,2) &
1734  + ref(kgrd(k,i,j,6,2),igrd(i,j,3),jgrd(i,j,3)) &
1735  * hfact(i,j,6) * vfact(k,i,j,6,2) &
1736  + ref(kgrd(k,i,j,7,2),igrd(i,j,3),jgrd(i,j,3)) &
1737  * hfact(i,j,7) * vfact(k,i,j,7,2) &
1738  + ref(kgrd(k,i,j,8,2),igrd(i,j,8),jgrd(i,j,8)) &
1739  * hfact(i,j,8) * vfact(k,i,j,8,2)
1740  end do
1741  end do
1742  end do
1743 
1744  ! logarithmic weighting (for pres, dens)
1745  if ( logarithmic ) then
1746  do j = 1, ja
1747  do i = 1, ia
1748  do k = ks, ke
1749  intp(k,i,j) = exp( intp(k,i,j) )
1750  end do
1751  end do
1752  end do
1753  endif
1754 
1755  return
1756  end subroutine intrpnest_interp_3d_8points
1757 
1758 
1759  !-----------------------------------------------------------------------------
1760  ! interpolation using twelve-points for 2D data
1761  subroutine intrpnest_interp_2d_12points( &
1762  intp, & ! (out)
1763  ref, & ! (in)
1764  hfact, & ! (in)
1765  igrd, & ! (in)
1766  jgrd, & ! (in)
1767  ia, & ! (in)
1768  ja ) ! (in)
1769  implicit none
1770 
1771  real(RP), intent(out) :: intp(:,:) ! interpolated data
1772 
1773  real(RP), intent(in) :: ref (:,:) ! reference data
1774  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1775  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1776  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1777  integer, intent(in) :: ia ! grid number of mine
1778  integer, intent(in) :: ja ! grid number of mine
1779 
1780  integer :: i, j
1781  !---------------------------------------------------------------------------
1782 
1783 !OCL PREFETCH
1784  do j = 1, ja
1785  do i = 1, ia
1786  intp(i,j) = ref(igrd(i,j,1), jgrd(i,j,1)) * hfact(i,j,1) &
1787  + ref(igrd(i,j,2), jgrd(i,j,2)) * hfact(i,j,2) &
1788  + ref(igrd(i,j,3), jgrd(i,j,3)) * hfact(i,j,3) &
1789  + ref(igrd(i,j,4), jgrd(i,j,4)) * hfact(i,j,4) &
1790  + ref(igrd(i,j,5), jgrd(i,j,5)) * hfact(i,j,5) &
1791  + ref(igrd(i,j,6), jgrd(i,j,6)) * hfact(i,j,6) &
1792  + ref(igrd(i,j,7), jgrd(i,j,7)) * hfact(i,j,7) &
1793  + ref(igrd(i,j,8), jgrd(i,j,8)) * hfact(i,j,8) &
1794  + ref(igrd(i,j,9), jgrd(i,j,9)) * hfact(i,j,9) &
1795  + ref(igrd(i,j,10),jgrd(i,j,10)) * hfact(i,j,10) &
1796  + ref(igrd(i,j,11),jgrd(i,j,11)) * hfact(i,j,11) &
1797  + ref(igrd(i,j,12),jgrd(i,j,12)) * hfact(i,j,12)
1798  end do
1799  end do
1800 
1801  return
1802  end subroutine intrpnest_interp_2d_12points
1803 
1804 
1805  !-----------------------------------------------------------------------------
1806  ! interpolation using twelve-points for 3D data
1807  subroutine intrpnest_interp_3d_12points( &
1808  intp, & ! (out)
1809  ref, & ! (in)
1810  hfact, & ! (in)
1811  vfact, & ! (in)
1812  kgrd, & ! (in)
1813  igrd, & ! (in)
1814  jgrd, & ! (in)
1815  ia, & ! (in)
1816  ja, & ! (in)
1817  ks, & ! (in)
1818  ke, & ! (in)
1819  logwegt ) ! (in)
1820  implicit none
1821 
1822  real(RP), intent(out) :: intp(:,:,:) ! interpolated data
1823 
1824  real(RP), intent(in) :: ref (:,:,:) ! reference data
1825  real(RP), intent(in) :: hfact(:,:,:) ! horizontal interp factor
1826  real(RP), intent(in) :: vfact(:,:,:,:,:) ! vertical interp factor
1827  integer, intent(in) :: kgrd (:,:,:,:,:) ! grid points of interp target
1828  integer, intent(in) :: igrd (:,:,:) ! grid points of interp target
1829  integer, intent(in) :: jgrd (:,:,:) ! grid points of interp target
1830  integer, intent(in) :: ia ! grid number of mine
1831  integer, intent(in) :: ja ! grid number of mine
1832  integer, intent(in) :: ks ! start grid number of mine
1833  integer, intent(in) :: ke ! end grid number of mine
1834 
1835  logical, intent(in), optional :: logwegt
1836 
1837  integer :: i, j, k
1838  logical :: logarithmic
1839  !---------------------------------------------------------------------------
1840 
1841  logarithmic = .false.
1842  if ( present(logwegt) ) then
1843  if ( logwegt ) then
1844  logarithmic = .true.
1845  endif
1846  endif
1847 
1848  ! linear interpolation
1849 !OCL PREFETCH
1850  do j = 1, ja
1851  do i = 1, ia
1852  do k = ks, ke
1853  intp(k,i,j) = ref(kgrd(k,i,j,1, 1),igrd(i,j,1 ),jgrd(i,j,1 )) &
1854  * hfact(i,j,1 ) * vfact(k,i,j,1, 1) &
1855  + ref(kgrd(k,i,j,2, 1),igrd(i,j,2 ),jgrd(i,j,2 )) &
1856  * hfact(i,j,2 ) * vfact(k,i,j,2, 1) &
1857  + ref(kgrd(k,i,j,3, 1),igrd(i,j,3 ),jgrd(i,j,3 )) &
1858  * hfact(i,j,3 ) * vfact(k,i,j,3, 1) &
1859  + ref(kgrd(k,i,j,4, 1),igrd(i,j,4 ),jgrd(i,j,4 )) &
1860  * hfact(i,j,4 ) * vfact(k,i,j,4, 1) &
1861  + ref(kgrd(k,i,j,5, 1),igrd(i,j,5 ),jgrd(i,j,5 )) &
1862  * hfact(i,j,5 ) * vfact(k,i,j,5, 1) &
1863  + ref(kgrd(k,i,j,6, 1),igrd(i,j,6 ),jgrd(i,j,6 )) &
1864  * hfact(i,j,6 ) * vfact(k,i,j,6, 1) &
1865  + ref(kgrd(k,i,j,7, 1),igrd(i,j,7 ),jgrd(i,j,7 )) &
1866  * hfact(i,j,7 ) * vfact(k,i,j,7, 1) &
1867  + ref(kgrd(k,i,j,8, 1),igrd(i,j,8 ),jgrd(i,j,8 )) &
1868  * hfact(i,j,8 ) * vfact(k,i,j,8, 1) &
1869  + ref(kgrd(k,i,j,9, 1),igrd(i,j,9 ),jgrd(i,j,9 )) &
1870  * hfact(i,j,9 ) * vfact(k,i,j,9, 1) &
1871  + ref(kgrd(k,i,j,10,1),igrd(i,j,10),jgrd(i,j,10)) &
1872  * hfact(i,j,10) * vfact(k,i,j,10,1) &
1873  + ref(kgrd(k,i,j,11,1),igrd(i,j,11),jgrd(i,j,11)) &
1874  * hfact(i,j,11) * vfact(k,i,j,11,1) &
1875  + ref(kgrd(k,i,j,12,1),igrd(i,j,12),jgrd(i,j,12)) &
1876  * hfact(i,j,12) * vfact(k,i,j,12,1) &
1877  + ref(kgrd(k,i,j,1, 2),igrd(i,j,1 ),jgrd(i,j,1 )) &
1878  * hfact(i,j,1 ) * vfact(k,i,j,1, 2) &
1879  + ref(kgrd(k,i,j,2, 2),igrd(i,j,2 ),jgrd(i,j,2 )) &
1880  * hfact(i,j,2 ) * vfact(k,i,j,2, 2) &
1881  + ref(kgrd(k,i,j,3, 2),igrd(i,j,3 ),jgrd(i,j,3 )) &
1882  * hfact(i,j,3 ) * vfact(k,i,j,3, 2) &
1883  + ref(kgrd(k,i,j,4, 2),igrd(i,j,4 ),jgrd(i,j,4 )) &
1884  * hfact(i,j,4 ) * vfact(k,i,j,4, 2) &
1885  + ref(kgrd(k,i,j,5, 2),igrd(i,j,5 ),jgrd(i,j,5 )) &
1886  * hfact(i,j,5 ) * vfact(k,i,j,5, 2) &
1887  + ref(kgrd(k,i,j,6, 2),igrd(i,j,6 ),jgrd(i,j,6 )) &
1888  * hfact(i,j,6 ) * vfact(k,i,j,6, 2) &
1889  + ref(kgrd(k,i,j,7, 2),igrd(i,j,7 ),jgrd(i,j,7 )) &
1890  * hfact(i,j,7 ) * vfact(k,i,j,7, 2) &
1891  + ref(kgrd(k,i,j,8, 2),igrd(i,j,8 ),jgrd(i,j,8 )) &
1892  * hfact(i,j,8 ) * vfact(k,i,j,8, 2) &
1893  + ref(kgrd(k,i,j,9, 2),igrd(i,j,9 ),jgrd(i,j,9 )) &
1894  * hfact(i,j,9 ) * vfact(k,i,j,9, 2) &
1895  + ref(kgrd(k,i,j,10,2),igrd(i,j,10),jgrd(i,j,10)) &
1896  * hfact(i,j,10) * vfact(k,i,j,10,2) &
1897  + ref(kgrd(k,i,j,11,2),igrd(i,j,11),jgrd(i,j,11)) &
1898  * hfact(i,j,11) * vfact(k,i,j,11,2) &
1899  + ref(kgrd(k,i,j,12,2),igrd(i,j,12),jgrd(i,j,12)) &
1900  * hfact(i,j,12) * vfact(k,i,j,12,2)
1901  end do
1902  end do
1903  end do
1904 
1905  ! logarithmic weighting (for pres, dens)
1906  if ( logarithmic ) then
1907  do j = 1, ja
1908  do i = 1, ia
1909  do k = ks, ke
1910  intp(k,i,j) = exp( intp(k,i,j) )
1911  end do
1912  end do
1913  end do
1914  endif
1915 
1916  return
1917  end subroutine intrpnest_interp_3d_12points
1918 
1919 
1920  !-----------------------------------------------------------------------------
1921  subroutine intrpnest_domain_compatibility( &
1922  lon_org, & ! (in)
1923  lat_org, & ! (in)
1924  lev_org, & ! (in)
1925  lon_loc, & ! (in)
1926  lat_loc, & ! (in)
1927  lev_loc, & ! (in)
1928  skip_x, & ! (in)
1929  skip_y, & ! (in)
1930  skip_z ) ! (in)
1931  use scale_process, only: &
1932  prc_mpistop
1933  use scale_const, only: &
1934  d2r => const_d2r
1935  implicit none
1936  real(RP), intent(in) :: lon_org(:,:)
1937  real(RP), intent(in) :: lat_org(:,:)
1938  real(RP), intent(in) :: lev_org(:,:,:)
1939  real(RP), intent(in) :: lon_loc(:,:)
1940  real(RP), intent(in) :: lat_loc(:,:)
1941  real(RP), intent(in) :: lev_loc(:,:,:)
1942  logical, intent(in), optional :: skip_x
1943  logical, intent(in), optional :: skip_y
1944  logical, intent(in), optional :: skip_z
1945 
1946  real(RP) :: max_ref, min_ref
1947  real(RP) :: max_loc, min_loc
1948 
1949  logical :: do_xdirec
1950  logical :: do_ydirec
1951  logical :: do_zdirec
1952 
1953  intrinsic size
1954  !---------------------------------------------------------------------------
1955 
1956  do_xdirec = .true.
1957  if ( present(skip_x) ) then
1958  if ( skip_x ) then
1959  do_xdirec = .false.
1960  endif
1961  endif
1962 
1963  do_ydirec = .true.
1964  if ( present(skip_y) ) then
1965  if ( skip_y ) then
1966  do_ydirec = .false.
1967  endif
1968  endif
1969 
1970  do_zdirec = .true.
1971  if ( present(skip_z) ) then
1972  if ( skip_z ) then
1973  do_zdirec = .false.
1974  endif
1975  endif
1976 
1977  if ( do_xdirec ) then
1978  max_ref = maxval( lon_org(:,:) / d2r )
1979  min_ref = minval( lon_org(:,:) / d2r )
1980  max_loc = maxval( lon_loc(:,:) / d2r )
1981  min_loc = minval( lon_loc(:,:) / d2r )
1982 
1983  if ( (min_ref+360.0_rp-max_ref) < 360.0_rp / size(lon_org,1) * 2.0_rp) then
1984  ! cyclic OK
1985  else if ( max_ref < max_loc .OR. min_ref > min_loc ) then
1986  write(*,*) 'xxx ERROR: REQUESTED DOMAIN IS TOO MUCH BROAD'
1987  write(*,*) 'xxx -- LONGITUDINAL direction over the limit'
1988  write(*,*) 'xxx -- reference max: ', max_ref
1989  write(*,*) 'xxx -- reference min: ', min_ref
1990  write(*,*) 'xxx -- local max: ', max_loc
1991  write(*,*) 'xxx -- local min: ', min_loc
1992  call prc_mpistop
1993  endif
1994  endif
1995 
1996  if ( do_ydirec ) then
1997  max_ref = maxval( lat_org(:,:) / d2r )
1998  min_ref = minval( lat_org(:,:) / d2r )
1999  max_loc = maxval( lat_loc(:,:) / d2r )
2000  min_loc = minval( lat_loc(:,:) / d2r )
2001 
2002  if ( max_ref < max_loc .OR. min_ref > min_loc ) then
2003  write(*,*) 'xxx ERROR: REQUESTED DOMAIN IS TOO MUCH BROAD'
2004  write(*,*) 'xxx -- LATITUDINAL direction over the limit'
2005  write(*,*) 'xxx -- reference max: ', max_ref
2006  write(*,*) 'xxx -- reference min: ', min_ref
2007  write(*,*) 'xxx -- local max: ', max_loc
2008  write(*,*) 'xxx -- local min: ', min_loc
2009  call prc_mpistop
2010  endif
2011  endif
2012 
2013  if ( do_zdirec ) then
2014  max_ref = maxval( lev_org(:,:,:) )
2015  !max_loc = maxval( lev_loc(KS-1:KE+1,:,:) ) ! HALO + 1
2016  max_loc = maxval( lev_loc(:,:,:) ) ! HALO + 1
2017  !min_ref = minval( lev_org(:,:,:) )
2018  !min_loc = minval( lev_loc(3:KA,:,:) ) ! HALO + 1
2019 
2020  if ( max_ref < max_loc ) then
2021  !if ( max_ref < max_loc .OR. min_ref > min_loc ) then
2022  write(*,*) 'xxx ERROR: REQUESTED DOMAIN IS TOO MUCH BROAD'
2023  write(*,*) 'xxx -- VERTICAL direction over the limit'
2024  write(*,*) 'xxx -- reference max: ', max_ref
2025  !write(*,*) 'xxx -- reference min: ', min_ref
2026  write(*,*) 'xxx -- local max: ', max_loc
2027  !write(*,*) 'xxx -- local min: ', min_loc
2028  call prc_mpistop
2029  endif
2030  endif
2031 
2032  return
2033  end subroutine intrpnest_domain_compatibility
2034 
2035  !-----------------------------------------------------------------------------
2036  ! Haversine Formula (from R.W. Sinnott, "Virtues of the Haversine",
2037  ! Sky and Telescope, vol. 68, no. 2, 1984, p. 159):
2038  function intrpnest_haversine( &
2039  la0, &
2040  lo0, &
2041  la, &
2042  lo ) &
2043  result( d )
2044  implicit none
2045  real(RP), intent(in) :: la0, lo0, la, lo ! la,la0: Lat, lo,lo0: Lon; [rad]
2046  real(RP) :: d, dlon, dlat, work1, work2
2047  !---------------------------------------------------------------------------
2048 
2049  ! output unit : [m]
2050  dlon = lo0 - lo
2051  dlat = la0 - la
2052  work1 = (sin(dlat/2.0_rp))**2.0_rp + &
2053  cos(la0) * cos(la) * (sin(dlon/2.0_rp))**2.0_rp
2054  work2 = 2.0_rp * asin(min( 1.0_rp, sqrt(work1) ))
2055  d = r_in_m * work2
2056 
2057  end function intrpnest_haversine
2058 
2059 end module scale_interpolation_nest
integer, public is
start point of inner domain: x, local
module DEBUG
Definition: scale_debug.F90:13
subroutine, public intrpnest_domain_compatibility(lon_org, lat_org, lev_org, lon_loc, lat_loc, lev_loc, skip_x, skip_y, skip_z)
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), public const_radius
radius of the planet [m]
Definition: scale_const.F90:46
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
procedure(intrpnest_intfc_interp_2d), pointer, public intrpnest_interp_2d
subroutine, public intrpnest_interp_fact_latlon(hfact, igrd, jgrd, mylat, mylon, myIA, myJA, inlat, inlon, inIA, inJA)
subroutine, public intrpnest_interp_fact_llz(hfact, vfact, kgrd, igrd, jgrd, ncopy, myhgt, mylat, mylon, myKS, myKE, myIA, myJA, inhgt, inlat, inlon, inKA, inIA, inJA, landgrid)
module grid index
module TRACER
module Index
Definition: scale_index.F90:14
integer, public ia
of x whole cells (local, with HALO)
integer, public js
start point of inner domain: y, local
module PROCESS
procedure(intrpnest_intfc_interp_3d), pointer, public intrpnest_interp_3d
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
integer, parameter, public khalo
of halo cells: z
module INTERPOLATION (nesting system)
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
real(rp), public const_eps
small number
Definition: scale_const.F90:36
module PRECISION
module land grid index
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine, public intrpnest_setup(interp_search_divnum, NEST_INTERP_LEVEL, OFFLINE)
Setup.
integer, public ja
of y whole cells (local, with HALO)