115 module amip_interp_mod
119 use time_manager_mod,
only:
time_type,
operator(+),
operator(>), &
129 horiz_interp_type,
assignment(=)
133 mpp_pe, lowercase, mpp_root_pe, &
136 use constants_mod,
only: tfreeze, pi
137 use platform_mod,
only: r4_kind, r8_kind, i2_kind, fms_file_len
138 use mpp_mod,
only: input_nml_file
139 use fms2_io_mod,
only: fmsnetcdffile_t, fms2_io_file_exists=>file_exists,
open_file,
close_file, &
140 get_dimension_size, fms2_io_read_data=>
read_data
141 use netcdf,
only: nf90_max_name
154 integer :: i_sst = 1200
155 integer :: j_sst = 600
156 real(r8_kind),
parameter:: big_number = 1.e30_r8_kind
157 logical :: forecast_mode = .false.
167 #include<file_version.h>
170 real(r8_kind),
allocatable,
dimension(:,:) :: tempamip
181 integer :: year, month, day
186 interface assignment(=)
192 interface operator (==)
198 interface operator (/=)
204 interface operator (>)
210 module procedure get_amip_sst_r4, get_amip_sst_r8
215 module procedure get_amip_ice_r4, get_amip_ice_r8
275 module procedure amip_interp_new_1d_r4, amip_interp_new_1d_r8
276 module procedure amip_interp_new_2d_r4, amip_interp_new_2d_r8
286 type (horiz_interp_type) :: Hintrp, Hintrp2
287 real(r4_kind),
dimension(:,:),
allocatable :: data1_r4, data2_r4
288 real(r8_kind),
dimension(:,:),
allocatable :: data1_r8, data2_r8
289 type (date_type) :: Date1, Date2
290 logical :: use_climo, use_annual
291 logical :: I_am_initialized=.false.
299 integer :: mobs, nobs
300 real(r8_kind),
allocatable,
dimension(:) :: lon_bnd, lat_bnd
305 character(len=FMS_FILE_LEN) :: file_name_sst, file_name_ice
306 type(FmsNetcdfFile_t),
target :: fileobj_sst, fileobj_ice
308 type (date_type) :: Curr_date =
date_type( -99, -99, -99 )
309 type (date_type) :: Date_end =
date_type( -99, -99, -99 )
311 real(r8_kind) :: tice_crit_k
312 integer(i2_kind) :: ice_crit
314 logical :: module_is_initialized = .false.
330 real(r8_kind) ::
teq = 305._r8_kind
331 real(r8_kind) ::
tdif = 50._r8_kind
332 real(r8_kind) ::
tann = 20._r8_kind
333 real(r8_kind) ::
tlag = 0.875_r8_kind
341 logical :: do_sst_pert = .false.
359 verbose, i_sst, j_sst, forecast_mode, &
368 integer :: iunit,io,ierr
376 read (input_nml_file, amip_interp_nml, iostat=io)
384 write (iunit,nml=amip_interp_nml)
390 'MPP_IO is no longer supported. Please remove use_mpp_io from amip_interp_nml',&
398 if ( tice_crit_k < 200._r8_kind )
then
399 tice_crit_k = tice_crit_k + tfreeze
401 ice_crit = nint((tice_crit_k-tfreeze)*100._r8_kind, i2_kind)
408 if (lowercase(trim(
data_set)) ==
'amip1')
then
409 file_name_sst =
'INPUT/' //
'amip1_sst.data'
410 file_name_ice =
'INPUT/' //
'amip1_sst.data'
411 mobs = 180; nobs = 91
412 call set_sst_grid_edges_amip1
414 call error_mesg (
'amip_interp_init',
'using AMIP 1 sst', note)
416 else if (lowercase(trim(
data_set)) ==
'amip2')
then
417 file_name_sst =
'INPUT/' //
'amip2_sst.data'
418 file_name_ice =
'INPUT/' //
'amip2_ice.data'
419 mobs = 360; nobs = 180
420 call set_sst_grid_edges_oi
422 tice_crit_k = 271.38_r8_kind
424 call error_mesg (
'amip_interp_init',
'using AMIP 2 sst', note)
426 else if (lowercase(trim(
data_set)) ==
'hurrell')
then
427 file_name_sst =
'INPUT/' //
'hurrell_sst.data'
428 file_name_ice =
'INPUT/' //
'hurrell_ice.data'
429 mobs = 360; nobs = 180
430 call set_sst_grid_edges_oi
432 tice_crit_k = 271.38_r8_kind
434 call error_mesg (
'amip_interp_init',
'using HURRELL sst', note)
437 else if (lowercase(trim(
data_set)) ==
'daily')
then
438 file_name_sst =
'INPUT/' //
'hurrell_sst.data'
439 file_name_ice =
'INPUT/' //
'hurrell_ice.data'
440 mobs = 360; nobs = 180
441 call set_sst_grid_edges_oi
443 call error_mesg (
'amip_interp_init',
'using AVHRR daily sst', note)
446 else if (lowercase(trim(
data_set)) ==
'reynolds_eof')
then
447 file_name_sst =
'INPUT/' //
'reynolds_sst.data'
448 file_name_ice =
'INPUT/' //
'reynolds_sst.data'
449 mobs = 180; nobs = 90
450 call set_sst_grid_edges_oi
453 'using NCEP Reynolds Historical Reconstructed SST', note)
455 else if (lowercase(trim(
data_set)) ==
'reynolds_oi')
then
456 file_name_sst =
'INPUT/' //
'reyoi_sst.data'
457 file_name_ice =
'INPUT/' //
'reyoi_sst.data'
460 mobs = i_sst; nobs = j_sst
462 mobs = 360; nobs = 180
465 call set_sst_grid_edges_oi
467 call error_mesg (
'amip_interp_init',
'using Reynolds OI SST', &
471 call error_mesg (
'amip_interp_init',
'the value of the &
472 &namelist parameter DATA_SET being used is not allowed', fatal)
476 print *,
'ice_crit,tice_crit_k=',ice_crit,tice_crit_k
479 file_name_sst = trim(file_name_sst)//
'.nc'
480 file_name_ice = trim(file_name_ice)//
'.nc'
482 if (.not. fms2_io_file_exists(trim(file_name_sst)) )
then
484 'file '//trim(file_name_sst)//
' does not exist', fatal)
486 if (.not. fms2_io_file_exists(trim(file_name_ice)) )
then
488 'file '//trim(file_name_ice)//
' does not exist', fatal)
491 if (.not.
open_file(fileobj_sst, trim(file_name_sst),
'read')) &
492 call error_mesg (
'amip_interp_init',
'Error in opening file '//trim(file_name_sst), fatal)
493 if (.not.
open_file(fileobj_ice, trim(file_name_ice),
'read')) &
494 call error_mesg (
'amip_interp_init',
'Error in opening file '//trim(file_name_ice), fatal)
495 module_is_initialized = .true.
498 subroutine set_sst_grid_edges_amip1
500 real(r8_kind) :: hpie, dlon, dlat, wb, sb
502 allocate(lon_bnd(mobs+1))
503 allocate(lat_bnd(nobs+1))
507 hpie = pi / 2._r8_kind
509 dlon = 4._r8_kind*hpie/real(mobs, r8_kind)
510 wb = -0.5_r8_kind*dlon
513 lon_bnd(i) = wb + dlon*real(i-1, r8_kind)
515 lon_bnd(mobs+1) = lon_bnd(1) + 4._r8_kind*hpie
517 dlat = 2._r8_kind*hpie/real(nobs-1, r8_kind)
518 sb = -hpie + 0.5_r8_kind*dlat
521 lat_bnd(nobs+1) = hpie
523 lat_bnd(j) = sb + dlat * real(j-2, r8_kind)
525 end subroutine set_sst_grid_edges_amip1
527 subroutine set_sst_grid_edges_oi
529 real(r8_kind) :: hpie, dlon, dlat, wb, sb
532 if(
allocated(lon_bnd))
deallocate(lon_bnd)
533 if(
allocated(lat_bnd))
deallocate(lat_bnd)
536 allocate(lon_bnd(mobs+1))
537 allocate(lat_bnd(nobs+1))
541 hpie = pi / 2._r8_kind
542 dlon = 4._r8_kind*hpie/real(mobs, r8_kind)
547 lon_bnd(i) = wb + dlon * real(i-1, r8_kind)
549 lon_bnd(mobs+1) = lon_bnd(1) + 4._r8_kind*hpie
551 dlat = 2._r8_kind*hpie/real(nobs, r8_kind)
555 lat_bnd(nobs+1) = hpie
557 lat_bnd(j) = sb + dlat * real(j-1, r8_kind)
559 end subroutine set_sst_grid_edges_oi
568 if(
allocated(interp%data1_r4))
deallocate(interp%data1_r4)
569 if(
allocated(interp%data1_r8))
deallocate(interp%data1_r8)
570 if(
allocated(interp%data2_r4))
deallocate(interp%data2_r4)
571 if(
allocated(interp%data2_r8))
deallocate(interp%data2_r8)
573 if(
allocated(lon_bnd))
deallocate(lon_bnd)
574 if(
allocated(lat_bnd))
deallocate(lat_bnd)
578 interp%I_am_initialized = .false.
586 integer,
intent(out) :: nlon
588 integer,
intent(out) :: nlat
593 nlon = mobs; nlat = nobs
598 type (
date_type),
intent(in) :: left, right
601 if (left % year == right % year .and. &
602 left % month == right % month .and. &
603 left % day == right % day )
then
612 type (
date_type),
intent(in) :: left, right
615 if (left % year == right % year .and. &
616 left % month == right % month .and. &
617 left % day == right % day )
then
625 function date_gt (Left, Right)
result (answer)
626 type (
date_type),
intent(in) :: left, right
630 dif(1) = left%year - right%year
631 dif(2) = left%month - right%month
632 dif(3) = left%day - right%day
635 if (dif(i) == 0) cycle
648 if(.not.amip_interp_in%I_am_initialized)
then
649 call mpp_error(fatal,
'amip_interp_type_eq: amip_interp_type variable on right hand side is unassigned')
652 amip_interp_out%Hintrp = amip_interp_in%Hintrp
653 amip_interp_out%Hintrp2 = amip_interp_in%Hintrp2
654 amip_interp_out%data1_r4 = amip_interp_in%data1_r4
655 amip_interp_out%data1_r8 = amip_interp_in%data1_r8
656 amip_interp_out%data2_r4 = amip_interp_in%data2_r4
657 amip_interp_out%data2_r8 = amip_interp_in%data2_r8
658 amip_interp_out%Date1 = amip_interp_in%Date1
659 amip_interp_out%Date2 = amip_interp_in%Date2
660 amip_interp_out%Date1 = amip_interp_in%Date1
661 amip_interp_out%Date2 = amip_interp_in%Date2
662 amip_interp_out%use_climo = amip_interp_in%use_climo
663 amip_interp_out%use_annual = amip_interp_in%use_annual
664 amip_interp_out%I_am_initialized = .true.
667 #include "amip_interp_r4.fh"
668 #include "amip_interp_r8.fh"
670 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...