FMS  2024.03
Flexible Modeling System
horiz_interp_conserve.F90
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 !> @defgroup horiz_interp_conserve_mod horiz_interp_conserve_mod
20 !> @ingroup horiz_interp
21 !> @brief Performs spatial interpolation between grids using conservative interpolation
22 !!
23 !> @author Bruce Wyman, Zhi Liang
24 !!
25 !> This module can conservatively interpolate data from any logically rectangular grid
26 !! to any rectangular grid. The interpolation scheme is area-averaging
27 !! conservative scheme. There is an optional mask field for missing input data in both
28 !! horiz_interp__conserveinit and horiz_interp_conserve. For efficiency purpose, mask should only be
29 !! kept in horiz_interp_init (will remove the mask in horiz_interp in the future).
30 !! There are 1-D and 2-D version of horiz_interp_conserve_init for 1-D and 2-D grid.
31 !! There is a optional argument mask in horiz_interp_conserve_init_2d and no mask should
32 !! to passed into horiz_interp_conserv. optional argument mask will not be passed into
33 !! horiz_interp_conserve_init_1d and optional argument mask may be passed into
34 !! horiz_interp_conserve (For the purpose of reproduce Memphis??? results).
35 !! An optional output mask field may be used in conjunction with the input mask to show
36 !! where output data exists.
37 
38 module horiz_interp_conserve_mod
39 
40  use platform_mod, only: r4_kind, r8_kind
41  use mpp_mod, only: mpp_send, mpp_recv, mpp_pe, mpp_root_pe, mpp_npes
42  use mpp_mod, only: mpp_error, fatal, mpp_sync_self
43  use mpp_mod, only: comm_tag_1, comm_tag_2
44  use fms_mod, only: write_version_number
45  use grid2_mod, only: get_great_circle_algorithm
46  use constants_mod, only: pi
47  use horiz_interp_type_mod, only: horiz_interp_type, conserve
48 
49 
50  implicit none
51  private
52 
53  ! public interface
54 
55 
56  !> @brief Allocates space and initializes a derived-type variable
57  !! that contains pre-computed interpolation indices and weights.
58  !!
59  !> Allocates space and initializes a derived-type variable
60  !! that contains pre-computed interpolation indices and weights
61  !! for improved performance of multiple interpolations between
62  !! the same grids.
63  !! @param lon_in
64  !! Longitude (in radians) for source data grid.
65  !!
66  !! @param lat_in
67  !! Latitude (in radians) for source data grid.
68  !!
69  !! @param lon_out
70  !! Longitude (in radians) for destination data grid.
71  !!
72  !! @param lat_out
73  !! Latitude (in radians) for destination data grid.
74  !!
75  !! @param verbose
76  !! flag for the amount of print output.
77  !!
78  !! @param mask_in
79  !! Input mask. must be the size (size(lon_in)-1, size(lon. The real value of
80  !! mask_in must be in the range (0.,1.). Set mask_in=0.0 for data points
81  !! that should not be used or have missing data.
82  !!
83  !! @param mask_out
84  !! Output mask that specifies whether data was computed.
85  !!
86  !! @param Interp
87  !! A derived-type variable containing indices and weights used for subsequent
88  !! interpolations. To reinitialize this variable for a different grid-to-grid
89  !! interpolation you must first use the "horiz_interp_del" interface.
90  !!
91  !> @ingroup horiz_interp_conserve_mod
93  module procedure horiz_interp_conserve_new_1dx1d_r4
94  module procedure horiz_interp_conserve_new_1dx2d_r4
95  module procedure horiz_interp_conserve_new_2dx1d_r4
96  module procedure horiz_interp_conserve_new_2dx2d_r4
97  module procedure horiz_interp_conserve_new_1dx1d_r8
98  module procedure horiz_interp_conserve_new_1dx2d_r8
99  module procedure horiz_interp_conserve_new_2dx1d_r8
100  module procedure horiz_interp_conserve_new_2dx2d_r8
101  end interface
102 
104  module procedure horiz_interp_conserve_r4
105  module procedure horiz_interp_conserve_r8
106  end interface
107 
108 !> private helper routines
109  interface data_sum
110  module procedure data_sum_r4
111  module procedure data_sum_r8
112  end interface
113 
114  interface stats
115  module procedure stats_r4
116  module procedure stats_r8
117  end interface
118 
120  module procedure horiz_interp_conserve_version1_r8
121  module procedure horiz_interp_conserve_version1_r4
122  end interface
123 
125  module procedure horiz_interp_conserve_version2_r8
126  module procedure horiz_interp_conserve_version2_r4
127  end interface
128 
129 
130 
131  !> @addtogroup horiz_interp_conserve_mod
132  !> @{
135 
136  integer :: pe, root_pe
137  !-----------------------------------------------------------------------
138  ! Include variable "version" to be written to log file.
139 #include<file_version.h>
140  logical :: module_is_initialized = .false.
141 
142  logical :: great_circle_algorithm = .false.
143 
144 contains
145 
146  !> Writes version number to logfile.
148 
149  if(module_is_initialized) return
150  call write_version_number("HORIZ_INTERP_CONSERVE_MOD", version)
151 
152  great_circle_algorithm = get_great_circle_algorithm()
153 
154  module_is_initialized = .true.
155 
156  end subroutine horiz_interp_conserve_init
157 
158  !> Deallocates memory used by "HI_KIND_TYPE" variables.
159  !! Must be called before reinitializing with horiz_interp_new.
160  subroutine horiz_interp_conserve_del ( Interp )
161 
162  type (horiz_interp_type), intent(inout) :: interp !< A derived-type variable returned by
163  !! previous call to horiz_interp_new. The input variable must have
164  !! allocated arrays. The returned variable will contain deallocated arrays.
165 
166  select case(interp%version)
167  case (1)
168  if( interp%horizInterpReals8_type%is_allocated) then
169  if(allocated(interp%horizInterpReals8_type%area_src)) deallocate(interp%horizInterpReals8_type%area_src)
170  if(allocated(interp%horizInterpReals8_type%area_dst)) deallocate(interp%horizInterpReals8_type%area_dst)
171  if(allocated(interp%horizInterpReals8_type%facj)) deallocate(interp%horizInterpReals8_type%facj)
172  if(allocated(interp%jlat)) deallocate(interp%jlat)
173  if(allocated(interp%horizInterpReals8_type%faci)) deallocate(interp%horizInterpReals8_type%faci)
174  if(allocated(interp%ilon)) deallocate(interp%ilon)
175  else if( interp%horizInterpReals4_type%is_allocated) then
176  if(allocated(interp%horizInterpReals4_type%area_src)) deallocate(interp%horizInterpReals4_type%area_src)
177  if(allocated(interp%horizInterpReals4_type%area_dst)) deallocate(interp%horizInterpReals4_type%area_dst)
178  if(allocated(interp%horizInterpReals4_type%facj)) deallocate(interp%horizInterpReals4_type%facj)
179  if(allocated(interp%jlat)) deallocate(interp%jlat)
180  if(allocated(interp%horizInterpReals4_type%faci)) deallocate(interp%horizInterpReals4_type%faci)
181  if(allocated(interp%ilon)) deallocate(interp%ilon)
182  endif
183  case (2)
184  if( interp%horizInterpReals8_type%is_allocated) then
185  if(allocated(interp%i_src)) deallocate(interp%i_src)
186  if(allocated(interp%j_src)) deallocate(interp%j_src)
187  if(allocated(interp%i_dst)) deallocate(interp%i_dst)
188  if(allocated(interp%j_dst)) deallocate(interp%j_dst)
189  if(allocated(interp%horizInterpReals8_type%area_frac_dst)) &
190  deallocate(interp%horizInterpReals8_type%area_frac_dst)
191  else if( interp%horizInterpReals4_type%is_allocated ) then
192  if(allocated(interp%i_src)) deallocate(interp%i_src)
193  if(allocated(interp%j_src)) deallocate(interp%j_src)
194  if(allocated(interp%i_dst)) deallocate(interp%i_dst)
195  if(allocated(interp%j_dst)) deallocate(interp%j_dst)
196  if(allocated(interp%horizInterpReals4_type%area_frac_dst)) &
197  deallocate(interp%horizInterpReals4_type%area_frac_dst)
198  endif
199  end select
200  interp%horizInterpReals4_type%is_allocated = .false.
201  interp%horizInterpReals8_type%is_allocated = .false.
202 
203  end subroutine horiz_interp_conserve_del
204 
205 #include "horiz_interp_conserve_r4.fh"
206 #include "horiz_interp_conserve_r8.fh"
207 
208 end module horiz_interp_conserve_mod
209 !> @}
210 ! close documentation grouping
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
Definition: fms.F90:758
logical function, public get_great_circle_algorithm()
Determine if we are using the great circle algorithm.
Definition: grid2.F90:171
subroutine, public horiz_interp_conserve_init
Writes version number to logfile.
subroutine, public horiz_interp_conserve_del(Interp)
Deallocates memory used by "HI_KIND_TYPE" variables. Must be called before reinitializing with horiz_...
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indi...
Holds data pointers and metadata for horizontal interpolations, passed between the horiz_interp modul...
subroutine mpp_sync_self(pelist, check, request, msg_size, msg_type)
This is to check if current PE's outstanding puts are complete but we can't use shmem_fence because w...
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:421
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407
Error handler.
Definition: mpp.F90:382
Recieve data from another PE.
Definition: mpp.F90:937
Send data to a receiving PE.
Definition: mpp.F90:1004