116 module amip_interp_mod
120 use time_manager_mod,
only:
time_type,
operator(+),
operator(>), &
130 horiz_interp_type,
assignment(=)
134 mpp_pe, lowercase, mpp_root_pe, &
137 use constants_mod,
only: tfreeze, pi
138 use platform_mod,
only: r4_kind, r8_kind, i2_kind, fms_file_len
139 use mpp_mod,
only: input_nml_file
140 use fms2_io_mod,
only: fmsnetcdffile_t, fms2_io_file_exists=>file_exists,
open_file,
close_file, &
141 get_dimension_size, fms2_io_read_data=>
read_data
142 use netcdf,
only: nf90_max_name
155 integer :: i_sst = 1200
156 integer :: j_sst = 600
157 real(r8_kind),
parameter:: big_number = 1.e30_r8_kind
158 logical :: forecast_mode = .false.
159 real(r8_kind),
allocatable,
dimension(:,:) :: sst_ncep, sst_anom
161 public i_sst, j_sst, sst_ncep, sst_anom, forecast_mode,
use_ncep_sst
169 #include<file_version.h>
172 real(r8_kind),
allocatable,
dimension(:,:) :: tempamip
183 integer :: year, month, day
188 interface assignment(=)
194 interface operator (==)
200 interface operator (/=)
206 interface operator (>)
212 module procedure get_amip_sst_r4, get_amip_sst_r8
217 module procedure get_amip_ice_r4, get_amip_ice_r8
277 module procedure amip_interp_new_1d_r4, amip_interp_new_1d_r8
278 module procedure amip_interp_new_2d_r4, amip_interp_new_2d_r8
288 type (horiz_interp_type) :: Hintrp, Hintrp2
289 real(r4_kind),
dimension(:,:),
allocatable :: data1_r4, data2_r4
290 real(r8_kind),
dimension(:,:),
allocatable :: data1_r8, data2_r8
291 type (date_type) :: Date1, Date2
292 logical :: use_climo, use_annual
293 logical :: I_am_initialized=.false.
301 integer :: mobs, nobs
302 real(r8_kind),
allocatable,
dimension(:) :: lon_bnd, lat_bnd
307 character(len=FMS_FILE_LEN) :: file_name_sst, file_name_ice
308 type(FmsNetcdfFile_t),
target :: fileobj_sst, fileobj_ice
310 type (date_type) :: Curr_date =
date_type( -99, -99, -99 )
311 type (date_type) :: Date_end =
date_type( -99, -99, -99 )
313 real(r8_kind) :: tice_crit_k
314 integer(i2_kind) :: ice_crit
316 logical :: module_is_initialized = .false.
332 real(r8_kind) ::
teq = 305._r8_kind
333 real(r8_kind) ::
tdif = 50._r8_kind
334 real(r8_kind) ::
tann = 20._r8_kind
335 real(r8_kind) ::
tlag = 0.875_r8_kind
343 logical :: do_sst_pert = .false.
361 verbose, i_sst, j_sst, forecast_mode, &
370 integer :: iunit,io,ierr
378 read (input_nml_file, amip_interp_nml, iostat=io)
386 write (iunit,nml=amip_interp_nml)
392 'MPP_IO is no longer supported. Please remove use_mpp_io from amip_interp_nml',&
400 if ( tice_crit_k < 200._r8_kind )
then
401 tice_crit_k = tice_crit_k + tfreeze
403 ice_crit = nint((tice_crit_k-tfreeze)*100._r8_kind, i2_kind)
410 if (lowercase(trim(
data_set)) ==
'amip1')
then
411 file_name_sst =
'INPUT/' //
'amip1_sst.data'
412 file_name_ice =
'INPUT/' //
'amip1_sst.data'
413 mobs = 180; nobs = 91
414 call set_sst_grid_edges_amip1
416 call error_mesg (
'amip_interp_init',
'using AMIP 1 sst', note)
418 else if (lowercase(trim(
data_set)) ==
'amip2')
then
419 file_name_sst =
'INPUT/' //
'amip2_sst.data'
420 file_name_ice =
'INPUT/' //
'amip2_ice.data'
421 mobs = 360; nobs = 180
422 call set_sst_grid_edges_oi
424 tice_crit_k = 271.38_r8_kind
426 call error_mesg (
'amip_interp_init',
'using AMIP 2 sst', note)
428 else if (lowercase(trim(
data_set)) ==
'hurrell')
then
429 file_name_sst =
'INPUT/' //
'hurrell_sst.data'
430 file_name_ice =
'INPUT/' //
'hurrell_ice.data'
431 mobs = 360; nobs = 180
432 call set_sst_grid_edges_oi
434 tice_crit_k = 271.38_r8_kind
436 call error_mesg (
'amip_interp_init',
'using HURRELL sst', note)
439 else if (lowercase(trim(
data_set)) ==
'daily')
then
440 file_name_sst =
'INPUT/' //
'hurrell_sst.data'
441 file_name_ice =
'INPUT/' //
'hurrell_ice.data'
442 mobs = 360; nobs = 180
443 call set_sst_grid_edges_oi
445 call error_mesg (
'amip_interp_init',
'using AVHRR daily sst', note)
448 else if (lowercase(trim(
data_set)) ==
'reynolds_eof')
then
449 file_name_sst =
'INPUT/' //
'reynolds_sst.data'
450 file_name_ice =
'INPUT/' //
'reynolds_sst.data'
451 mobs = 180; nobs = 90
452 call set_sst_grid_edges_oi
455 'using NCEP Reynolds Historical Reconstructed SST', note)
457 else if (lowercase(trim(
data_set)) ==
'reynolds_oi')
then
458 file_name_sst =
'INPUT/' //
'reyoi_sst.data'
459 file_name_ice =
'INPUT/' //
'reyoi_sst.data'
462 mobs = i_sst; nobs = j_sst
463 if (.not.
allocated (sst_ncep))
then
464 allocate (sst_ncep(i_sst,j_sst))
465 sst_ncep(:,:) = big_number
467 if (.not.
allocated (sst_anom))
then
468 allocate (sst_anom(i_sst,j_sst))
469 sst_anom(:,:) = big_number
472 mobs = 360; nobs = 180
475 call set_sst_grid_edges_oi
477 call error_mesg (
'amip_interp_init',
'using Reynolds OI SST', &
481 call error_mesg (
'amip_interp_init',
'the value of the &
482 &namelist parameter DATA_SET being used is not allowed', fatal)
486 print *,
'ice_crit,tice_crit_k=',ice_crit,tice_crit_k
489 file_name_sst = trim(file_name_sst)//
'.nc'
490 file_name_ice = trim(file_name_ice)//
'.nc'
492 if (.not. fms2_io_file_exists(trim(file_name_sst)) )
then
494 'file '//trim(file_name_sst)//
' does not exist', fatal)
496 if (.not. fms2_io_file_exists(trim(file_name_ice)) )
then
498 'file '//trim(file_name_ice)//
' does not exist', fatal)
501 if (.not.
open_file(fileobj_sst, trim(file_name_sst),
'read')) &
502 call error_mesg (
'amip_interp_init',
'Error in opening file '//trim(file_name_sst), fatal)
503 if (.not.
open_file(fileobj_ice, trim(file_name_ice),
'read')) &
504 call error_mesg (
'amip_interp_init',
'Error in opening file '//trim(file_name_ice), fatal)
505 module_is_initialized = .true.
508 subroutine set_sst_grid_edges_amip1
510 real(r8_kind) :: hpie, dlon, dlat, wb, sb
512 allocate(lon_bnd(mobs+1))
513 allocate(lat_bnd(nobs+1))
517 hpie = pi / 2._r8_kind
519 dlon = 4._r8_kind*hpie/real(mobs, r8_kind)
520 wb = -0.5_r8_kind*dlon
523 lon_bnd(i) = wb + dlon*real(i-1, r8_kind)
525 lon_bnd(mobs+1) = lon_bnd(1) + 4._r8_kind*hpie
527 dlat = 2._r8_kind*hpie/real(nobs-1, r8_kind)
528 sb = -hpie + 0.5_r8_kind*dlat
531 lat_bnd(nobs+1) = hpie
533 lat_bnd(j) = sb + dlat * real(j-2, r8_kind)
535 end subroutine set_sst_grid_edges_amip1
537 subroutine set_sst_grid_edges_oi
539 real(r8_kind) :: hpie, dlon, dlat, wb, sb
542 if(
allocated(lon_bnd))
deallocate(lon_bnd)
543 if(
allocated(lat_bnd))
deallocate(lat_bnd)
546 allocate(lon_bnd(mobs+1))
547 allocate(lat_bnd(nobs+1))
551 hpie = pi / 2._r8_kind
552 dlon = 4._r8_kind*hpie/real(mobs, r8_kind)
557 lon_bnd(i) = wb + dlon * real(i-1, r8_kind)
559 lon_bnd(mobs+1) = lon_bnd(1) + 4._r8_kind*hpie
561 dlat = 2._r8_kind*hpie/real(nobs, r8_kind)
565 lat_bnd(nobs+1) = hpie
567 lat_bnd(j) = sb + dlat * real(j-1, r8_kind)
569 end subroutine set_sst_grid_edges_oi
578 if(
allocated(interp%data1_r4))
deallocate(interp%data1_r4)
579 if(
allocated(interp%data1_r8))
deallocate(interp%data1_r8)
580 if(
allocated(interp%data2_r4))
deallocate(interp%data2_r4)
581 if(
allocated(interp%data2_r8))
deallocate(interp%data2_r8)
583 if(
allocated(lon_bnd))
deallocate(lon_bnd)
584 if(
allocated(lat_bnd))
deallocate(lat_bnd)
588 interp%I_am_initialized = .false.
596 integer,
intent(out) :: nlon
598 integer,
intent(out) :: nlat
603 nlon = mobs; nlat = nobs
608 type (
date_type),
intent(in) :: left, right
611 if (left % year == right % year .and. &
612 left % month == right % month .and. &
613 left % day == right % day )
then
622 type (
date_type),
intent(in) :: left, right
625 if (left % year == right % year .and. &
626 left % month == right % month .and. &
627 left % day == right % day )
then
635 function date_gt (Left, Right)
result (answer)
636 type (
date_type),
intent(in) :: left, right
640 dif(1) = left%year - right%year
641 dif(2) = left%month - right%month
642 dif(3) = left%day - right%day
645 if (dif(i) == 0) cycle
658 if(.not.amip_interp_in%I_am_initialized)
then
659 call mpp_error(fatal,
'amip_interp_type_eq: amip_interp_type variable on right hand side is unassigned')
662 amip_interp_out%Hintrp = amip_interp_in%Hintrp
663 amip_interp_out%Hintrp2 = amip_interp_in%Hintrp2
664 amip_interp_out%data1_r4 = amip_interp_in%data1_r4
665 amip_interp_out%data1_r8 = amip_interp_in%data1_r8
666 amip_interp_out%data2_r4 = amip_interp_in%data2_r4
667 amip_interp_out%data2_r8 = amip_interp_in%data2_r8
668 amip_interp_out%Date1 = amip_interp_in%Date1
669 amip_interp_out%Date2 = amip_interp_in%Date2
670 amip_interp_out%Date1 = amip_interp_in%Date1
671 amip_interp_out%Date2 = amip_interp_in%Date2
672 amip_interp_out%use_climo = amip_interp_in%use_climo
673 amip_interp_out%use_annual = amip_interp_in%use_annual
674 amip_interp_out%I_am_initialized = .true.
677 #include "amip_interp_r4.fh"
678 #include "amip_interp_r8.fh"
680 end module amip_interp_mod
real(r8_kind) tice_crit
in degC or degK
character(len=24) data_set
use 'amip1', 'amip2', 'reynolds_eof' 'reynolds_oi', 'hurrell', or 'daily', when "use_daily=....
real(r8_kind) teq
parameters for prescribed zonal sst option
logical function date_not_equals(Left, Right)
logical interp_oi_sst
changed to false for regular runs
subroutine, public amip_interp_init
initialize amip_interp_mod for use
logical use_mpp_io
Set to .true. to use mpp_io, otherwise fms2io is used.
real(r8_kind) tann
parameters for prescribed zonal sst option
logical no_anom_sst
SJL: During nudging: use_ncep_sst = .T.; no_anom_sst = .T. during forecast: use_ncep_sst = ....
logical use_ncep_ice
For seasonal forecast: use_ncep_ice = .F.
logical function date_equals(Left, Right)
real(r8_kind) sst_pert
global temperature perturbation used for sensitivity experiments
integer verbose
0 <= verbose <= 3
character(len=16) date_out_of_range
use 'fail', 'initclimo', or 'climo'
subroutine, public amip_interp_del(Interp)
Frees data associated with a amip_interp_type variable. Should be used for any variables initialized ...
logical, public use_ncep_sst
SJL: During nudging: use_ncep_sst = .T.; no_anom_sst = .T. during forecast: use_ncep_sst = ....
real(r8_kind) tlag
parameters for prescribed zonal sst option
integer, dimension(3) amip_date
amip date for repeating single day (rsd) option
logical use_zonal
parameters for prescribed zonal sst option
logical use_daily
if '.true.', give 'data_set = 'daily''
logical function date_gt(Left, Right)
subroutine get_sst_grid_size(nlon, nlat)
Returns the size (i.e., number of longitude and latitude points) of the observed data grid.
subroutine amip_interp_type_eq(amip_interp_out, amip_interp_in)
character(len=6) sst_pert_type
use 'random' or 'fixed'
real(r8_kind) tdif
parameters for prescribed zonal sst option
Initializes data needed for the horizontal interpolation between the sst data and model grid.
Contains information needed by the interpolation module (exchange_mod) and buffers data (r4_kind flav...
Private data type for representing a calendar date.
Close a netcdf or domain file opened with open_file or open_virtual_file.
Opens a given netcdf or domain file.
Read data from a defined field in a file.
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 error_mesg(routine, message, level)
Print notes, warnings and error messages; terminates program for warning and error messages....
Added for mixed precision support. Updates force time_manager math to be done with kind=8 reals _wrap...
subroutine, public horiz_interp_del(Interp)
Deallocates memory used by "horiz_interp_type" variables. Must be called before reinitializing with h...
subroutine, public horiz_interp_init
Initialize module and writes version number to logfile.out.
Subroutine for performing the horizontal interpolation between two grids.
integer function stdlog()
This function returns the current standard fortran unit numbers for log messages. Log messages,...
integer function mpp_pe()
Returns processor ID.
real(r8_kind) function, public fraction_of_year(Time)
Wrapper function to return the fractional time into the current year Always returns an r8_kind,...
Returns a weight and dates or indices for interpolating between two dates. The interface fraction_of_...
subroutine, public get_date(time, year, month, day, hour, minute, second, tick, err_msg)
Gets the date for different calendar types. Given a time_interval, returns the corresponding date und...
Given an input date in year, month, days, etc., creates a time_type that represents this time interva...
Given some number of seconds and days, returns the corresponding time_type.
Type to represent amounts of time. Implemented as seconds and days to allow for larger intervals.
AMIP interpolation for ice.
Retrieve sea surface temperature data and interpolated grid.
Allocates space and initializes a derived-type variable that contains pre-computed interpolation indi...