FMS 2025.01-dev
Flexible Modeling System
Loading...
Searching...
No Matches
horiz_interp_bicubic.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_bicubic_mod horiz_interp_bicubic_mod
20!> @ingroup horiz_interp
21!> @brief Delivers methods for bicubic interpolation from a coarse regular grid
22!! on a fine regular grid
23!!
24!> This module delivers methods for bicubic interpolation from a
25!! coarse regular grid on a fine regular grid.
26!! Subroutines
27!!
28!! - @ref bcuint
29!! - @ref bcucof
30!!
31!! are methods taken from
32!!
33!! W. H. Press, S. A. Teukolski, W. T. Vetterling and B. P. Flannery,
34!! Numerical Recipies in FORTRAN, The Art of Scientific Computing.
35!! Cambridge University Press, 1992
36!!
37!! written by
38!! martin.schmidt@io-warnemuende.de (2004)
39!! revised by
40!! martin.schmidt@io-warnemuende.de (2004)
41!!
42!! Version 1.0.0.2005-07-06
43!! The module is thought to interact with MOM-4.
44!! Alle benotigten Felder werden extern von MOM verwaltet, da sie
45!! nicht fur alle interpolierten Daten die gleiche Dimension haben mussen.
46module horiz_interp_bicubic_mod
47
48 use mpp_mod, only: mpp_error, fatal, stdout, mpp_pe, mpp_root_pe
49 use fms_mod, only: write_version_number
50 use horiz_interp_type_mod, only: horiz_interp_type, bicubic
51 use constants_mod, only: pi
52 use platform_mod, only: r4_kind, r8_kind
53
54
55 implicit none
56
57 private
58
61
62 !> Creates a new @ref horiz_interp_type for bicubic interpolation.
63 !! Allocates space and initializes a derived-type variable
64 !! that contains pre-computed interpolation indices and weights.
65 !> @ingroup horiz_interp_bicubic_mod
67 module procedure horiz_interp_bicubic_new_1d_r8
68 module procedure horiz_interp_bicubic_new_1d_s_r8
69 module procedure horiz_interp_bicubic_new_1d_r4
70 module procedure horiz_interp_bicubic_new_1d_s_r4
71 end interface
72
73 !> @brief Perform bicubic horizontal interpolation
75 module procedure horiz_interp_bicubic_r4
76 module procedure horiz_interp_bicubic_r8
77 end interface
78
79!> @addtogroup horiz_interp_bicubic_mod
80!> @{
81
82! Include variable "version" to be written to log file.
83#include<file_version.h>
84 logical :: module_is_initialized = .false.
85 integer :: verbose_bicubic = 0
86
87! Grid variables
88! xc, yc : co-ordinates of the coarse grid
89! xf, yf : co-ordinates of the fine grid
90! fc : variable to be interpolated at the coarse grid
91! dfc_x : x-derivative of fc at the coarse grid
92! dfc_y : y-derivative of fc at the coarse grid
93! dfc_xy : x-y-derivative of fc at the coarse grid
94! ff : variable to be interpolated at the fine grid
95! dff_x : x-derivative of fc at the fine grid
96! dff_y : y-derivative of fc at the fine grid
97! dff_xy : x-y-derivative of fc at the fine grid
98
99
100 real(r8_kind) :: tpi
101
102 !! Private interfaces for mixed precision helper routines
103
104 interface fill_xy
105 module procedure fill_xy_r4
106 module procedure fill_xy_r8
107 end interface
108
109 interface bcuint
110 module procedure bcuint_r4
111 module procedure bcuint_r8
112 end interface
113
114 interface bcucof
115 module procedure bcucof_r4
116 module procedure bcucof_r8
117 end interface
118
119 !> find the lower neighbour of xf in field xc, return is the index
120 interface indl
121 module procedure indl_r4
122 module procedure indl_r8
123 end interface
124
125 !> find the upper neighbour of xf in field xc, return is the index
126 interface indu
127 module procedure indu_r4
128 module procedure indu_r8
129 end interface
130
131 contains
132
133 !> @brief Initializes module and writes version number to logfile.out
135
136 if(module_is_initialized) return
137 call write_version_number("HORIZ_INTERP_BICUBIC_MOD", version)
138 module_is_initialized = .true.
139 tpi = real(2.0_r8_kind*pi, r8_kind)
140
141 end subroutine horiz_interp_bicubic_init
142
143 !> Free memory from a horiz_interp_type used for bicubic interpolation
144 !! (allocated via @ref horiz_bicubic_new)
145 subroutine horiz_interp_bicubic_del( Interp )
146 type(horiz_interp_type), intent(inout) :: interp
147
148 if(interp%horizInterpReals8_type%is_allocated) then
149 if(allocated(interp%horizInterpReals8_type%rat_x)) deallocate ( interp%horizInterpReals8_type%rat_x )
150 if(allocated(interp%horizInterpReals8_type%rat_y)) deallocate ( interp%horizInterpReals8_type%rat_y )
151 if(allocated(interp%horizInterpReals8_type%lon_in)) deallocate ( interp%horizInterpReals8_type%lon_in )
152 if(allocated(interp%horizInterpReals8_type%lat_in)) deallocate ( interp%horizInterpReals8_type%lat_in )
153 if(allocated(interp%horizInterpReals8_type%wti)) deallocate ( interp%horizInterpReals8_type%wti )
154 else if(interp%horizInterpReals4_type%is_allocated) then
155 if(allocated(interp%horizInterpReals4_type%rat_x)) deallocate ( interp%horizInterpReals4_type%rat_x )
156 if(allocated(interp%horizInterpReals4_type%rat_y)) deallocate ( interp%horizInterpReals4_type%rat_y )
157 if(allocated(interp%horizInterpReals4_type%lon_in)) deallocate ( interp%horizInterpReals4_type%lon_in )
158 if(allocated(interp%horizInterpReals4_type%lat_in)) deallocate ( interp%horizInterpReals4_type%lat_in )
159 if(allocated(interp%horizInterpReals4_type%wti)) deallocate ( interp%horizInterpReals4_type%wti )
160 endif
161 if( allocated(interp%i_lon) ) deallocate( interp%i_lon )
162 if( allocated(interp%j_lat) ) deallocate( interp%j_lat )
163
164 interp%horizInterpReals8_type%is_allocated = .false.
165 interp%horizInterpReals4_type%is_allocated = .false.
166
167 end subroutine horiz_interp_bicubic_del
168
169#include "horiz_interp_bicubic_r4.fh"
170#include "horiz_interp_bicubic_r8.fh"
171
172end module horiz_interp_bicubic_mod
173!> @}
174! close documentation
subroutine, public horiz_interp_bicubic_init
Initializes module and writes version number to logfile.out.
subroutine, public horiz_interp_bicubic_del(interp)
Free memory from a horiz_interp_type used for bicubic interpolation (allocated via horiz_bicubic_new)
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...
Error handler.
Definition mpp.F90:382
Perform bicubic horizontal interpolation.