FMS  2025.04
Flexible Modeling System
horiz_interp_bicubic.F90
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @defgroup horiz_interp_bicubic_mod horiz_interp_bicubic_mod
19 !> @ingroup horiz_interp
20 !> @brief Delivers methods for bicubic interpolation from a coarse regular grid
21 !! on a fine regular grid
22 !!
23 !> This module delivers methods for bicubic interpolation from a
24 !! coarse regular grid on a fine regular grid.
25 !! Subroutines
26 !!
27 !! - @ref bcuint
28 !! - @ref bcucof
29 !!
30 !! are methods taken from
31 !!
32 !! W. H. Press, S. A. Teukolski, W. T. Vetterling and B. P. Flannery,
33 !! Numerical Recipies in FORTRAN, The Art of Scientific Computing.
34 !! Cambridge University Press, 1992
35 !!
36 !! written by
37 !! martin.schmidt@io-warnemuende.de (2004)
38 !! revised by
39 !! martin.schmidt@io-warnemuende.de (2004)
40 !!
41 !! Version 1.0.0.2005-07-06
42 !! The module is thought to interact with MOM-4.
43 !! Alle benotigten Felder werden extern von MOM verwaltet, da sie
44 !! nicht fur alle interpolierten Daten die gleiche Dimension haben mussen.
45 module horiz_interp_bicubic_mod
46 
47  use mpp_mod, only: mpp_error, fatal, stdout, mpp_pe, mpp_root_pe
48  use fms_mod, only: write_version_number
49  use horiz_interp_type_mod, only: horiz_interp_type, bicubic
50  use constants_mod, only: pi
51  use platform_mod, only: r4_kind, r8_kind
52 
53 
54  implicit none
55 
56  private
57 
60 
61  !> Creates a new @ref horiz_interp_type for bicubic interpolation.
62  !! Allocates space and initializes a derived-type variable
63  !! that contains pre-computed interpolation indices and weights.
64  !> @ingroup horiz_interp_bicubic_mod
66  module procedure horiz_interp_bicubic_new_1d_r8
67  module procedure horiz_interp_bicubic_new_1d_s_r8
68  module procedure horiz_interp_bicubic_new_1d_r4
69  module procedure horiz_interp_bicubic_new_1d_s_r4
70  end interface
71 
72  !> @brief Perform bicubic horizontal interpolation
74  module procedure horiz_interp_bicubic_r4
75  module procedure horiz_interp_bicubic_r8
76  end interface
77 
78 !> @addtogroup horiz_interp_bicubic_mod
79 !> @{
80 
81 ! Include variable "version" to be written to log file.
82 #include<file_version.h>
83  logical :: module_is_initialized = .false.
84  integer :: verbose_bicubic = 0
85 
86 ! Grid variables
87 ! xc, yc : co-ordinates of the coarse grid
88 ! xf, yf : co-ordinates of the fine grid
89 ! fc : variable to be interpolated at the coarse grid
90 ! dfc_x : x-derivative of fc at the coarse grid
91 ! dfc_y : y-derivative of fc at the coarse grid
92 ! dfc_xy : x-y-derivative of fc at the coarse grid
93 ! ff : variable to be interpolated at the fine grid
94 ! dff_x : x-derivative of fc at the fine grid
95 ! dff_y : y-derivative of fc at the fine grid
96 ! dff_xy : x-y-derivative of fc at the fine grid
97 
98 
99  real(r8_kind) :: tpi
100 
101  !! Private interfaces for mixed precision helper routines
102 
103  interface fill_xy
104  module procedure fill_xy_r4
105  module procedure fill_xy_r8
106  end interface
107 
108  interface bcuint
109  module procedure bcuint_r4
110  module procedure bcuint_r8
111  end interface
112 
113  interface bcucof
114  module procedure bcucof_r4
115  module procedure bcucof_r8
116  end interface
117 
118  !> find the lower neighbour of xf in field xc, return is the index
119  interface indl
120  module procedure indl_r4
121  module procedure indl_r8
122  end interface
123 
124  !> find the upper neighbour of xf in field xc, return is the index
125  interface indu
126  module procedure indu_r4
127  module procedure indu_r8
128  end interface
129 
130  contains
131 
132  !> @brief Initializes module and writes version number to logfile.out
134 
135  if(module_is_initialized) return
136  call write_version_number("HORIZ_INTERP_BICUBIC_MOD", version)
137  module_is_initialized = .true.
138  tpi = real(2.0_r8_kind*pi, r8_kind)
139 
140  end subroutine horiz_interp_bicubic_init
141 
142  !> Free memory from a horiz_interp_type used for bicubic interpolation
143  !! (allocated via @ref horiz_bicubic_new)
144  subroutine horiz_interp_bicubic_del( Interp )
145  type(horiz_interp_type), intent(inout) :: interp
146 
147  if(interp%horizInterpReals8_type%is_allocated) then
148  if(allocated(interp%horizInterpReals8_type%rat_x)) deallocate ( interp%horizInterpReals8_type%rat_x )
149  if(allocated(interp%horizInterpReals8_type%rat_y)) deallocate ( interp%horizInterpReals8_type%rat_y )
150  if(allocated(interp%horizInterpReals8_type%lon_in)) deallocate ( interp%horizInterpReals8_type%lon_in )
151  if(allocated(interp%horizInterpReals8_type%lat_in)) deallocate ( interp%horizInterpReals8_type%lat_in )
152  if(allocated(interp%horizInterpReals8_type%wti)) deallocate ( interp%horizInterpReals8_type%wti )
153  else if(interp%horizInterpReals4_type%is_allocated) then
154  if(allocated(interp%horizInterpReals4_type%rat_x)) deallocate ( interp%horizInterpReals4_type%rat_x )
155  if(allocated(interp%horizInterpReals4_type%rat_y)) deallocate ( interp%horizInterpReals4_type%rat_y )
156  if(allocated(interp%horizInterpReals4_type%lon_in)) deallocate ( interp%horizInterpReals4_type%lon_in )
157  if(allocated(interp%horizInterpReals4_type%lat_in)) deallocate ( interp%horizInterpReals4_type%lat_in )
158  if(allocated(interp%horizInterpReals4_type%wti)) deallocate ( interp%horizInterpReals4_type%wti )
159  endif
160  if( allocated(interp%i_lon) ) deallocate( interp%i_lon )
161  if( allocated(interp%j_lat) ) deallocate( interp%j_lat )
162 
163  interp%horizInterpReals8_type%is_allocated = .false.
164  interp%horizInterpReals4_type%is_allocated = .false.
165 
166  end subroutine horiz_interp_bicubic_del
167 
168 #include "horiz_interp_bicubic_r4.fh"
169 #include "horiz_interp_bicubic_r8.fh"
170 
171 end module horiz_interp_bicubic_mod
172 !> @}
173 ! close documentation
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:701
subroutine, public horiz_interp_bicubic_del(Interp)
Free memory from a horiz_interp_type used for bicubic interpolation (allocated via horiz_bicubic_new)
subroutine, public horiz_interp_bicubic_init
Initializes module and writes version number to logfile.out.
Creates a new horiz_interp_type for bicubic interpolation. Allocates space and initializes a derived-...
find the lower neighbour of xf in field xc, return is the index
find the upper neighbour of xf in field xc, return is the index
Holds data pointers and metadata for horizontal interpolations, passed between the horiz_interp modul...
integer function stdout()
This function returns the current standard fortran unit numbers for output.
Definition: mpp_util.inc:42
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406
Error handler.
Definition: mpp.F90:381
Perform bicubic horizontal interpolation.