37 module horiz_interp_mod
48 use mpp_mod,
only: input_nml_file, warning,
mpp_pe, mpp_root_pe
49 use constants_mod,
only: pi
51 use horiz_interp_type_mod,
only: conserve, bilinear, spherical, bicubic
61 use platform_mod,
only: r4_kind, r8_kind
127 module procedure horiz_interp_new_1d_r4
128 module procedure horiz_interp_new_1d_r8
130 module procedure horiz_interp_new_1d_src_r4
131 module procedure horiz_interp_new_1d_src_r8
133 module procedure horiz_interp_new_2d_r4
134 module procedure horiz_interp_new_2d_r8
136 module procedure horiz_interp_new_1d_dst_r4
137 module procedure horiz_interp_new_1d_dst_r8
143 module procedure horiz_interp_read_weights_r4
144 module procedure horiz_interp_read_weights_r8
204 module procedure horiz_interp_base_2d_r4
205 module procedure horiz_interp_base_2d_r8
206 module procedure horiz_interp_base_3d_r4
207 module procedure horiz_interp_base_3d_r8
208 module procedure horiz_interp_solo_1d_r4
209 module procedure horiz_interp_solo_1d_r8
210 module procedure horiz_interp_solo_1d_src_r4
211 module procedure horiz_interp_solo_1d_src_r8
212 module procedure horiz_interp_solo_2d_r4
213 module procedure horiz_interp_solo_2d_r8
214 module procedure horiz_interp_solo_1d_dst_r4
215 module procedure horiz_interp_solo_1d_dst_r8
216 module procedure horiz_interp_solo_old_r4
217 module procedure horiz_interp_solo_old_r8
222 module procedure is_lat_lon_r4
223 module procedure is_lat_lon_r8
227 module procedure horiz_interp_solo_1d_r4
228 module procedure horiz_interp_solo_1d_r8
245 #include<file_version.h>
246 logical :: module_is_initialized = .false.
255 integer :: iunit, ierr, io
257 if(module_is_initialized)
return
258 call write_version_number(
"HORIZ_INTERP_MOD", version)
260 read (input_nml_file, horiz_interp_nml, iostat=io)
262 if (
mpp_pe() == mpp_root_pe() )
then
264 write (iunit, nml=horiz_interp_nml)
268 call mpp_error(fatal,
"horiz_interp_mod: You have overridden the default value of " // &
269 "reproduce_siena and set it to .true. in horiz_interp_nml. This was a temporary workaround to " // &
270 "allow for consistency in continuing experiments and is no longer supported. " // &
271 "Please remove this namelist.")
274 call horiz_interp_conserve_init
275 call horiz_interp_bilinear_init
276 call horiz_interp_bicubic_init
277 call horiz_interp_spherical_init
279 module_is_initialized = .true.
296 select case(interp % interp_method)
300 call horiz_interp_bilinear_del(interp )
307 interp%I_am_initialized = .false.
319 #include "horiz_interp_r4.fh"
320 #include "horiz_interp_r8.fh"
322 end module horiz_interp_mod
integer function, public check_nml_error(IOSTAT, NML_NAME)
Checks the iostat argument that is returned after reading a namelist and determines if the error code...
subroutine, public write_version_number(version, tag, unit)
Prints to the log file (or a specified unit) the version id string and tag name.
logical function, public fms_error_handler(routine, message, err_msg)
Facilitates the control of fatal error conditions.
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-...
subroutine, public horiz_interp_bilinear_del(Interp)
Deallocates memory used by "horiz_interp_type" variables.
subroutine, public horiz_interp_bilinear_init
Initialize this module and writes version number to logfile.
Creates a horiz_interp_type for bilinear interpolation.
Subroutines for reading in weight files and using that to fill in the horiz_interp type instead calcu...
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...
subroutine, public horiz_interp_del(Interp)
Deallocates memory used by "horiz_interp_type" variables. Must be called before reinitializing with h...
logical reproduce_siena
Set reproduce_siena = .true. to reproduce siena results. Set reproduce_siena = .false....
subroutine, public horiz_interp_init
Initialize module and writes version number to logfile.out.
subroutine, public horiz_interp_end
Dummy routine.
Subroutine for performing the horizontal interpolation between two grids.
subroutine, public horiz_interp_spherical_init
Initializes module and writes version number to logfile.out.
subroutine, public horiz_interp_spherical_del(Interp)
Deallocates memory used by "HI_KIND_TYPE" variables. Must be called before reinitializing with horiz_...
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.
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
integer function mpp_pe()
Returns processor ID.
Reduction operations. Find the min of scalar a from the PEs in pelist result is also automatically br...
Perform bicubic horizontal interpolation.
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indi...
Subroutines for reading in weight files and using that to fill in the horiz_interp type instead calcu...