MPAS-JEDI
mpasjedi_unstructured_interp_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2020 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
7 
8 use fckit_log_module, only: fckit_log
9 use fckit_mpi_module, only: fckit_mpi_comm, &
10  & fckit_mpi_sum, &
11  & fckit_mpi_max, &
12  & fckit_mpi_min
13 ! oops
14 use unstructured_interpolation_mod, only: unstrc_interp
15 
16 !mpas-jedi
18 use mpas_geom_mod
19 
20 private
21 
22 public :: unsinterp_integer_apply, &
24 
25 contains
26 
27 ! ------------------------------------------------------------------------------
28 
29 ! This subroutine populates field_out(n) with the field_in variable value that
30 ! has the greatest sum of weights contributed from all the neighbors of the
31 ! output location. (Only works with a discrete-valued field_in)
32 ! Note: assumes wtype='barycent'.
33 subroutine unsinterp_integer_apply(unsinterp, field_in, field_out)
34 
35  type(unstrc_interp), intent(in) :: unsinterp
36  real(kind=kind_real), intent(in) :: field_in(:) !Integer field in
37  real(kind=kind_real), intent(inout) :: field_out(:) !Integer field out
38 
39  integer :: maxtypel, mintypel, maxtype, mintype, ngrid_out
40  integer :: i, j, k, n, index
41  real(kind=kind_real), allocatable :: interp_w(:,:)
42  real(kind=kind_real), allocatable :: field_out_tmp(:)
43  real(kind=kind_real), allocatable :: field_neighbors(:,:)
44  real(kind=kind_real), allocatable :: field_types(:)
45 
46  ! Inteprolation of integer fields
47 
48  ! Size of output
49  ngrid_out = size(field_out)
50 
51  ! Get nearest neighbors
52  allocate(field_neighbors(unsinterp%nn,ngrid_out))
53  allocate(field_out_tmp(ngrid_out))
54  call unsinterp%apply(field_in, field_out_tmp, field_neighbors)
55 
56  ! Global min and max integers in field
57  maxtypel = int(maxval(field_in))
58  mintypel = int(minval(field_in))
59  call unsinterp%comm%allreduce(maxtypel,maxtype,fckit_mpi_max())
60  call unsinterp%comm%allreduce(mintypel,mintype,fckit_mpi_min())
61 
62  ! Put weights into field type array and pick max for interpolated value
63  allocate(field_types(mintype:maxtype))
64 
65  field_out = 0.0_kind_real
66  do i = 1,ngrid_out
67  field_types = 0.0
68  do n = 1, unsinterp%nn
69  index = int(field_neighbors(n,i))
70  field_types(index) = field_types(index) + unsinterp%interp_w(n,i)
71  enddo
72  field_out(i) = real(maxloc(field_types,1)+(mintype-1),kind_real)
73  enddo
74 
75 end subroutine unsinterp_integer_apply
76 
77 ! ------------------------------------------------------------------------------
78 
79 ! This subroutine populates field_out(n) with the field_in variable value contained
80 ! in the nearest neighbor (i.e. neighbor with greatest weight) of the output location.
81 ! (Would work with either discrete- or continuous-valued field_in)
82 ! Note: assumes wtype='barycent'.
83 subroutine unsinterp_nearest_apply(unsinterp, field_in, field_out)
84 
85  type(unstrc_interp), intent(in) :: unsinterp
86  real(kind=kind_real), intent(in) :: field_in(:) !Integer field in
87  real(kind=kind_real), intent(inout) :: field_out(:) !Integer field out
88 
89  integer :: n, ngrid_out
90  real(kind=kind_real), allocatable :: field_out_tmp(:)
91  real(kind=kind_real), allocatable :: field_neighbors(:,:)
92 
93  ! Inteprolation using the nearest neighbor
94 
95  ! Size of output
96  ngrid_out = size(field_out)
97 
98  ! Get nearest neighbors
99  allocate(field_neighbors(unsinterp%nn,ngrid_out))
100  allocate(field_out_tmp(ngrid_out))
101  call unsinterp%apply(field_in, field_out_tmp, field_neighbors)
102 
103  ! Find nearest neighbor
104  do n = 1, ngrid_out
105  field_out(n) = field_neighbors(maxloc(unsinterp%interp_w(:,n),1),n)
106  enddo
107 
108 end subroutine unsinterp_nearest_apply
109 
110 ! ------------------------------------------------------------------------------
subroutine, public unsinterp_integer_apply(unsinterp, field_in, field_out)
subroutine, public unsinterp_nearest_apply(unsinterp, field_in, field_out)