FMS 2025.01-dev
Flexible Modeling System
|
Delivers methods for bicubic interpolation from a coarse regular grid on a fine regular grid. More...
Data Types | |
interface | bcucof |
interface | bcuint |
interface | fill_xy |
interface | horiz_interp_bicubic_new |
Creates a new horiz_interp_type for bicubic interpolation. Allocates space and initializes a derived-type variable that contains pre-computed interpolation indices and weights. More... | |
interface | indl |
find the lower neighbour of xf in field xc, return is the index More... | |
interface | indu |
find the upper neighbour of xf in field xc, return is the index More... | |
Functions/Subroutines | |
subroutine | bcucof_ (y, y1, y2, y12, d1, d2, c) |
bcucof_r4 | |
bcucof_r8 | |
subroutine | bcuint_ (y, y1, y2, y12, x1l, x1u, x2l, x2u, t, u, ansy, ansy1, ansy2) |
bcuint_r4 | |
bcuint_r8 | |
subroutine | fill_xy_ (fi, ics, ice, jcs, jce, mask, maxpass) |
fill_xy_r4 | |
fill_xy_r8 | |
subroutine | horiz_interp_bicubic_ (interp, data_in, data_out, verbose, mask_in, mask_out, missing_value, missing_permit) |
Perform bicubic horizontal interpolation. | |
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. | |
subroutine | horiz_interp_bicubic_new_1d_ (interp, lon_in, lat_in, lon_out, lat_out, verbose, src_modulo) |
Creates a new horiz_interp_type. | |
subroutine | horiz_interp_bicubic_new_1d_s_ (interp, lon_in, lat_in, lon_out, lat_out, verbose, src_modulo) |
Creates a new horiz_interp_type. | |
integer function | indl_ (xc, xf) |
find the lower neighbour of xf in field xc, return is the index | |
indl_r4 | |
indl_r8 | |
integer function | indu_ (xc, xf) |
find the upper neighbour of xf in field xc, return is the index | |
indu_r4 | |
indu_r8 | |
Variables | |
logical | module_is_initialized = .FALSE. |
real(r8_kind) | tpi |
integer | verbose_bicubic = 0 |
Delivers methods for bicubic interpolation from a coarse regular grid on a fine regular grid.
This module delivers methods for bicubic interpolation from a coarse regular grid on a fine regular grid. Subroutines
are methods taken from
W. H. Press, S. A. Teukolski, W. T. Vetterling and B. P. Flannery, Numerical Recipies in FORTRAN, The Art of Scientific Computing. Cambridge University Press, 1992
written by marti.nosp@m.n.sc.nosp@m.hmidt.nosp@m.@io-.nosp@m.warne.nosp@m.muen.nosp@m.de.de (2004) revised by marti.nosp@m.n.sc.nosp@m.hmidt.nosp@m.@io-.nosp@m.warne.nosp@m.muen.nosp@m.de.de (2004)
Version 1.0.0.2005-07-06 The module is thought to interact with MOM-4. Alle benotigten Felder werden extern von MOM verwaltet, da sie nicht fur alle interpolierten Daten die gleiche Dimension haben mussen.
interface horiz_interp_bicubic_mod::bcucof |
Definition at line 114 of file horiz_interp_bicubic.F90.
Public Member Functions | |
bcucof_r4 | |
bcucof_r8 | |
interface horiz_interp_bicubic_mod::bcuint |
Definition at line 109 of file horiz_interp_bicubic.F90.
Public Member Functions | |
bcuint_r4 | |
bcuint_r8 | |
interface horiz_interp_bicubic_mod::fill_xy |
Definition at line 104 of file horiz_interp_bicubic.F90.
Public Member Functions | |
fill_xy_r4 | |
fill_xy_r8 | |
interface horiz_interp_bicubic_mod::horiz_interp_bicubic_new |
Creates a new horiz_interp_type for bicubic interpolation. Allocates space and initializes a derived-type variable that contains pre-computed interpolation indices and weights.
Definition at line 66 of file horiz_interp_bicubic.F90.
Public Member Functions | |
horiz_interp_bicubic_new_1d_r4 | |
horiz_interp_bicubic_new_1d_r8 | |
horiz_interp_bicubic_new_1d_s_r4 | |
horiz_interp_bicubic_new_1d_s_r8 | |
interface horiz_interp_bicubic_mod::indl |
find the lower neighbour of xf in field xc, return is the index
Definition at line 120 of file horiz_interp_bicubic.F90.
Public Member Functions | |
indl_r4 | |
indl_r8 | |
interface horiz_interp_bicubic_mod::indu |
find the upper neighbour of xf in field xc, return is the index
Definition at line 126 of file horiz_interp_bicubic.F90.
Public Member Functions | |
indu_r4 | |
indu_r8 | |
subroutine bcucof_ | ( | real(fms_hi_kind_), dimension(4) | y, |
real(fms_hi_kind_), dimension(4) | y1, | ||
real(fms_hi_kind_), dimension(4) | y2, | ||
real(fms_hi_kind_), dimension(4) | y12, | ||
real(fms_hi_kind_) | d1, | ||
real(fms_hi_kind_) | d2, | ||
real(fms_hi_kind_), dimension(4,4) | c | ||
) |
Definition at line 460 of file horiz_interp_bicubic.inc.
subroutine bcuint_ | ( | real(fms_hi_kind_), dimension(4) | y, |
real(fms_hi_kind_), dimension(4) | y1, | ||
real(fms_hi_kind_), dimension(4) | y2, | ||
real(fms_hi_kind_), dimension(4) | y12, | ||
real(fms_hi_kind_) | x1l, | ||
real(fms_hi_kind_) | x1u, | ||
real(fms_hi_kind_) | x2l, | ||
real(fms_hi_kind_) | x2u, | ||
real(fms_hi_kind_) | t, | ||
real(fms_hi_kind_) | u, | ||
real(fms_hi_kind_) | ansy, | ||
real(fms_hi_kind_) | ansy1, | ||
real(fms_hi_kind_) | ansy2 | ||
) |
Definition at line 438 of file horiz_interp_bicubic.inc.
subroutine fill_xy_ | ( | real(fms_hi_kind_), dimension(ics:ice,jcs:jce), intent(inout) | fi, |
integer, intent(in) | ics, | ||
integer, intent(in) | ice, | ||
integer, intent(in) | jcs, | ||
integer, intent(in) | jce, | ||
real(fms_hi_kind_), dimension(ics:ice,jcs:jce), intent(in), optional | mask, | ||
integer, intent(in) | maxpass | ||
) |
Definition at line 548 of file horiz_interp_bicubic.inc.
subroutine horiz_interp_bicubic_ | ( | type (horiz_interp_type), intent(in) | interp, |
real(fms_hi_kind_), dimension(:,:), intent(in) | data_in, | ||
real(fms_hi_kind_), dimension(:,:), intent(out) | data_out, | ||
integer, intent(in), optional | verbose, | ||
real(fms_hi_kind_), dimension(:,:), intent(in), optional | mask_in, | ||
real(fms_hi_kind_), dimension(:,:), intent(out), optional | mask_out, | ||
real(fms_hi_kind_), intent(in), optional | missing_value, | ||
integer, intent(in), optional | missing_permit | ||
) |
Perform bicubic horizontal interpolation.
Definition at line 354 of file horiz_interp_bicubic.inc.
subroutine, public horiz_interp_bicubic_del | ( | type(horiz_interp_type), intent(inout) | interp | ) |
Free memory from a horiz_interp_type used for bicubic interpolation (allocated via horiz_bicubic_new)
Definition at line 145 of file horiz_interp_bicubic.F90.
subroutine, public horiz_interp_bicubic_init |
Initializes module and writes version number to logfile.out.
Definition at line 134 of file horiz_interp_bicubic.F90.
subroutine horiz_interp_bicubic_new_1d_ | ( | type(horiz_interp_type), intent(inout) | interp, |
real(fms_hi_kind_), dimension(:), intent(in) | lon_in, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lat_in, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lon_out, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lat_out, | ||
integer, intent(in), optional | verbose, | ||
logical, intent(in), optional | src_modulo | ||
) |
Creates a new horiz_interp_type.
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indices and weights.
Definition at line 201 of file horiz_interp_bicubic.inc.
subroutine horiz_interp_bicubic_new_1d_s_ | ( | type(horiz_interp_type), intent(inout) | interp, |
real(fms_hi_kind_), dimension(:), intent(in) | lon_in, | ||
real(fms_hi_kind_), dimension(:), intent(in) | lat_in, | ||
real(fms_hi_kind_), dimension(:,:), intent(in) | lon_out, | ||
real(fms_hi_kind_), dimension(:,:), intent(in) | lat_out, | ||
integer, intent(in), optional | verbose, | ||
logical, intent(in), optional | src_modulo | ||
) |
Creates a new horiz_interp_type.
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indices and weights.
[in,out] | interp | A derived-type variable containing indices and weights used for subsequent interpolations. To reinitialize this variable for a different grid-to-grid interpolation you must first use the HORIZ_INTERP_BICUBIC_NEW__del interface. |
[in] | lon_in | Longitude (radians) for source data grid |
[in] | lat_in | Latitude (radians) for source data grid |
[in] | lon_out | Longitude (radians) for output data grid |
[in] | lat_out | Latitude (radians) for output data grid |
[in] | verbose | flag for print output amount |
[in] | src_modulo | indicates if the boundary condition along zonal boundary is cyclic or not. Zonal boundary condition is cyclic when true |
Definition at line 26 of file horiz_interp_bicubic.inc.
integer function indl_ | ( | real(fms_hi_kind_), dimension(1:), intent(in) | xc, |
real(fms_hi_kind_), intent(in) | xf | ||
) |
find the lower neighbour of xf in field xc, return is the index
Definition at line 516 of file horiz_interp_bicubic.inc.
integer function indu_ | ( | real(fms_hi_kind_), dimension(1:), intent(in) | xc, |
real(fms_hi_kind_), intent(in) | xf | ||
) |
find the upper neighbour of xf in field xc, return is the index
Definition at line 533 of file horiz_interp_bicubic.inc.
|
private |
Definition at line 84 of file horiz_interp_bicubic.F90.
|
private |
Definition at line 100 of file horiz_interp_bicubic.F90.
|
private |
Definition at line 85 of file horiz_interp_bicubic.F90.