112  CHARACTER(LEN=3) :: donst
   113  INTEGER :: idim, jdim, lsm, lsoil, lugb, iy, im, id, ih, ialb
   114  INTEGER :: isot, ivegsrc, lensfc, zsea1_mm, zsea2_mm, ierr
   115  INTEGER :: nprocs, myrank, num_threads, 
num_parthds, max_tasks
   116  REAL    :: fh, deltsfc, zsea1, zsea2
   117  LOGICAL :: use_ufo, do_nsst, do_lndinc, do_sfccycle
   119  NAMELIST/namcyc/ idim,jdim,lsm,lsoil,lugb,iy,im,id,ih,fh,&
   120                   deltsfc,ialb,use_ufo,donst,             &
   121                   do_sfccycle,isot,ivegsrc,zsea1_mm,      &
   122                   zsea2_mm, max_tasks, do_lndinc
   124  DATA idim,jdim,lsm,lsoil/96,96,1,4/
   125  DATA iy,im,id,ih,fh/1997,8,2,0,0./
   126  DATA lugb/51/, deltsfc/0.0/, ialb/1/, max_tasks/99999/
   127  DATA isot/1/, ivegsrc/2/, zsea1_mm/0/, zsea2_mm/0/
   130  CALL mpi_comm_size(mpi_comm_world, nprocs, ierr)
   131  CALL mpi_comm_rank(mpi_comm_world, myrank, ierr)
   133  if (myrank==0) 
call w3tagb(
'GLOBAL_CYCLE',2018,0179,0055,
'NP20')
   138  print*,
"STARTING CYCLE PROGRAM ON RANK ", myrank
   139  print*,
"RUNNING WITH ", nprocs, 
"TASKS"   140  print*,
"AND WITH ", num_threads, 
" THREADS."   148  print*,
"READ NAMCYC NAMELIST."   150  CALL baopenr(36, 
"fort.36", ierr)
   154  IF (max_tasks < 99999 .AND. myrank > (max_tasks - 1)) 
THEN   155    print*,
"USER SPECIFIED MAX NUMBER OF TASKS: ", max_tasks
   156    print*,
"WILL NOT RUN CYCLE PROGRAM ON RANK: ", myrank
   162  zsea1 = float(zsea1_mm) / 1000.0  
   163  zsea2 = float(zsea2_mm) / 1000.0
   165  IF (donst == 
"YES") 
THEN   172  IF (myrank==0) print*,
"LUGB,IDIM,JDIM,LSM,ISOT,IVEGSRC,LSOIL,DELTSFC,IY,IM,ID,IH,FH: ", &
   173               lugb,idim,jdim,lsm,isot,ivegsrc,lsoil,deltsfc,iy,im,id,ih,fh
   175  CALL sfcdrv(lugb,idim,jdim,lsm,lensfc,lsoil,deltsfc,  &
   176              iy,im,id,ih,fh,ialb,                  &
   177              use_ufo,do_nsst,do_sfccycle,do_lndinc, &
   178              zsea1,zsea2,isot,ivegsrc,myrank)
   181  print*,
'CYCLE PROGRAM COMPLETED NORMALLY ON RANK: ', myrank
   185  CALL mpi_barrier(mpi_comm_world, ierr)
   187  if (myrank==0) 
call w3tage(
'GLOBAL_CYCLE')
   189  CALL mpi_finalize(ierr)
   304  SUBROUTINE sfcdrv(LUGB, IDIM,JDIM,LSM,LENSFC,LSOIL,DELTSFC,  &
   305                    IY,IM,ID,IH,FH,IALB,                  &
   306                    USE_UFO,DO_NSST,DO_SFCCYCLE,DO_LNDINC,&
   307                    ZSEA1,ZSEA2,ISOT,IVEGSRC,MYRANK)
   312  USE land_increments, 
ONLY: add_increment_soil,     &
   313                             add_increment_snow,     &
   314                             calculate_landinc_mask, &
   315                             apply_land_da_adjustments_stc, &
   316                             apply_land_da_adjustments_snd
   320  INTEGER, 
INTENT(IN) :: IDIM, JDIM, LSM,LENSFC, LSOIL, IALB
   321  INTEGER, 
INTENT(IN) :: LUGB, IY, IM, ID, IH
   322  INTEGER, 
INTENT(IN) :: ISOT, IVEGSRC, MYRANK
   324  LOGICAL, 
INTENT(IN) :: USE_UFO, DO_NSST,DO_SFCCYCLE
   325  LOGICAL, 
INTENT(IN) :: DO_LNDINC
   327  REAL, 
INTENT(IN)    :: FH, DELTSFC, ZSEA1, ZSEA2
   329  INTEGER, 
PARAMETER  :: NLUNIT=35
   330  INTEGER, 
PARAMETER  :: SZ_NML=1
   332  CHARACTER(LEN=5)    :: TILE_NUM
   333  CHARACTER(LEN=500)  :: NST_FILE
   334  CHARACTER(LEN=500)  :: LND_SOI_FILE
   335  CHARACTER(LEN=4)    :: INPUT_NML_FILE(SZ_NML)
   338  INTEGER             :: I_INDEX(LENSFC), J_INDEX(LENSFC)
   339  INTEGER             :: IDUM(IDIM,JDIM)
   340  integer             :: num_parthds, num_threads
   342  real(kind=kind_io8) :: min_ice(lensfc)
   344  REAL                :: SLMASK(LENSFC), OROG(LENSFC)
   345  REAL                :: SIHFCS(LENSFC), SICFCS(LENSFC)
   346  REAL                :: SITFCS(LENSFC), TSFFCS(LENSFC)
   347  REAL                :: SWEFCS(LENSFC), ZORFCS(LENSFC)
   348  REAL                :: ALBFCS(LENSFC,4), TG3FCS(LENSFC)
   349  REAL                :: CNPFCS(LENSFC), SMCFCS(LENSFC,LSOIL)
   350  REAL                :: STCFCS(LENSFC,LSOIL), SLIFCS(LENSFC)
   351  REAL                :: AISFCS(LENSFC), F10M(LENSFC)
   352  REAL                :: VEGFCS(LENSFC), VETFCS(LENSFC)
   353  REAL                :: SOTFCS(LENSFC), ALFFCS(LENSFC,2)
   354  REAL                :: CVFCS(LENSFC), CVTFCS(LENSFC)
   355  REAL                :: CVBFCS(LENSFC), TPRCP(LENSFC)
   356  REAL                :: SRFLAG(LENSFC), SNDFCS(LENSFC)
   357  REAL                :: SLCFCS(LENSFC,LSOIL), VMXFCS(LENSFC)
   358  REAL                :: VMNFCS(LENSFC), T2M(LENSFC)
   359  REAL                :: Q2M(LENSFC), SLPFCS(LENSFC)
   360  REAL                :: ABSFCS(LENSFC), OROG_UF(LENSFC)
   361  REAL                :: USTAR(LENSFC)
   362  REAL                :: FMM(LENSFC), FHH(LENSFC)
   363  REAL                :: RLA(LENSFC), RLO(LENSFC)
   364  REAL(KIND=4)        :: ZSOIL(LSOIL)
   365  REAL                :: SIG1T(LENSFC)
   368  REAL, 
ALLOCATABLE   :: STC_BCK(:,:), SMC_BCK(:,:), SLC_BCK(:,:)
   369  REAL, 
ALLOCATABLE   :: SLIFCS_FG(:)
   370  INTEGER, 
ALLOCATABLE :: LANDINC_MASK_FG(:), LANDINC_MASK(:)
   371  REAL, 
ALLOCATABLE   :: SND_BCK(:), SND_INC(:), SWE_BCK(:)
   372  REAL(KIND=KIND_IO8), 
ALLOCATABLE :: SLMASKL(:), SLMASKW(:)
   374  TYPE(NSST_DATA)     :: NSST
   375  real, 
dimension(idim,jdim) :: tf_clm,tf_trd,sal_clm
   376  real, 
dimension(lensfc)    :: tf_clm_tile,tf_trd_tile,sal_clm_tile
   377  INTEGER             :: veg_type_landice
   379  LOGICAL :: FILE_EXISTS, DO_SOI_INC, DO_SNO_INC
   385  NAMELIST/namsfcd/ nst_file, lnd_soi_file, do_sno_inc
   387  DATA nst_file/
'NULL'/
   388  DATA lnd_soi_file/
'NULL'/
   396  input_nml_file = 
"NULL"   398  CALL baopenr(37, 
"fort.37", ierr)
   399  READ (37, nml=namsfcd)
   402  print*,
'IN ROUTINE SFCDRV,IDIM=',idim,
'JDIM=',jdim,
'FH=',fh
   408  CALL read_lat_lon_orog(rla,rlo,orog,orog_uf,tile_num,idim,jdim,lensfc)
   414  i_index = reshape(idum, (/lensfc/))
   420  j_index = reshape(idum, (/lensfc/) )
   424    print*,
"WILL PROCESS NSST RECORDS."   425    ALLOCATE(nsst%C_0(lensfc))
   426    ALLOCATE(nsst%C_D(lensfc))
   427    ALLOCATE(nsst%D_CONV(lensfc))
   428    ALLOCATE(nsst%DT_COOL(lensfc))
   429    ALLOCATE(nsst%IFD(lensfc))
   430    ALLOCATE(nsst%QRAIN(lensfc))
   431    ALLOCATE(nsst%TREF(lensfc))
   432    ALLOCATE(nsst%TFINC(lensfc))
   433    ALLOCATE(nsst%W_0(lensfc))
   434    ALLOCATE(nsst%W_D(lensfc))
   435    ALLOCATE(nsst%XS(lensfc))
   436    ALLOCATE(nsst%XT(lensfc))
   437    ALLOCATE(nsst%XTTS(lensfc))
   438    ALLOCATE(nsst%XU(lensfc))
   439    ALLOCATE(nsst%XV(lensfc))
   440    ALLOCATE(nsst%XZ(lensfc))
   441    ALLOCATE(nsst%XZTS(lensfc))
   442    ALLOCATE(nsst%Z_C(lensfc))
   443    ALLOCATE(nsst%ZM(lensfc))
   444    ALLOCATE(slifcs_fg(lensfc))
   449    IF  (trim(lnd_soi_file) .NE. 
"NULL") 
THEN   452        print*,
" APPLYING SOIL INCREMENTS FROM THE GSI"   453        ALLOCATE(stc_bck(lensfc, lsoil), smc_bck(lensfc, lsoil), slc_bck(lensfc,lsoil))
   454        ALLOCATE(landinc_mask_fg(lensfc))
   462        print*,
" APPLYING SNOW INCREMENTS FROM JEDI"   463        ALLOCATE(snd_bck(lensfc), snd_inc(lensfc), swe_bck(lensfc))
   466    ALLOCATE(landinc_mask(lensfc))
   467    if (ivegsrc == 2) 
then      478  CALL read_data(lsoil,lensfc,do_nsst,.false.,tsffcs=tsffcs,smcfcs=smcfcs,   &
   479                 swefcs=swefcs,stcfcs=stcfcs,tg3fcs=tg3fcs,zorfcs=zorfcs,  &
   480                 cvfcs=cvfcs,  cvbfcs=cvbfcs,cvtfcs=cvtfcs,albfcs=albfcs,  &
   481                 vegfcs=vegfcs,slifcs=slifcs,cnpfcs=cnpfcs,f10m=f10m    ,  &
   482                 vetfcs=vetfcs,sotfcs=sotfcs,alffcs=alffcs,ustar=ustar  ,  &
   483                 fmm=fmm      ,fhh=fhh      ,sihfcs=sihfcs,sicfcs=sicfcs,  &
   484                 sitfcs=sitfcs,tprcp=tprcp  ,srflag=srflag,sndfcs=sndfcs,  &
   485                 vmnfcs=vmnfcs,vmxfcs=vmxfcs,slcfcs=slcfcs,slpfcs=slpfcs,  &
   486                 absfcs=absfcs,t2m=t2m      ,q2m=q2m      ,slmask=slmask,  &
   487                 zsoil=zsoil,   nsst=nsst)
   491    print*,
'USE UNFILTERED OROGRAPHY.'   498    IF(nint(slifcs(i)).EQ.2) aisfcs(i) = 1.
   502    IF (.NOT. do_sfccycle ) 
THEN   504      print*,
"FIRST GUESS MASK ADJUSTED BY IFD RECORD"   506      WHERE(nint(nsst%IFD) == 3) slifcs_fg = 2.0
   509      print*,
"SAVE FIRST GUESS MASK"   516     CALL calculate_landinc_mask(slcfcs(:,1),swefcs, vetfcs,  &
   517                     lensfc,veg_type_landice,  landinc_mask)
   527  IF (do_sfccycle) 
THEN   528    ALLOCATE(slmaskl(lensfc), slmaskw(lensfc))
   531      IF(nint(slmask(i)) == 1) 
THEN   532        slmaskl(i) = 1.0_kind_io8
   533        slmaskw(i) = 1.0_kind_io8
   535        slmaskl(i) = 0.0_kind_io8
   536        slmaskw(i) = 0.0_kind_io8
   538      if(nint(slmask(i)) == 0) 
then   539        min_ice(i) = 0.15_kind_io8
   541        min_ice(i) = 0.0_kind_io8
   544    num_threads = num_parthds()
   546    print*,
"CALL SFCCYCLE TO UPDATE SURFACE FIELDS."   547    CALL sfccycle(lugb,lensfc,lsoil,sig1t,deltsfc,          &
   548                iy,im,id,ih,fh,rla,rlo,                   &
   549                slmaskl,slmaskw, orog, orog_uf, use_ufo, do_nsst,   &
   550                sihfcs,sicfcs,sitfcs,sndfcs,slcfcs,       &
   551                vmnfcs,vmxfcs,slpfcs,absfcs,              &
   552                tsffcs,swefcs,zorfcs,albfcs,tg3fcs,       &
   553                cnpfcs,smcfcs,stcfcs,slifcs,aisfcs,       &
   554                vegfcs,vetfcs,sotfcs,alffcs,              &
   555                cvfcs,cvbfcs,cvtfcs,myrank,num_threads, nlunit,        &
   556                sz_nml, input_nml_file,                   &
   558                ialb,isot,ivegsrc,tile_num,i_index,j_index)
   559    DEALLOCATE(slmaskl, slmaskw)
   569    IF (nst_file == 
"NULL") 
THEN   571      print*,
"NO GSI FILE.  ADJUST IFD FOR FORMER ICE POINTS."   573        IF (nint(slifcs_fg(i)) == 2 .AND. nint(slifcs(i)) == 0) 
THEN   580      print*,
"ADJUST TREF FROM GSI INCREMENT"   584      call get_tf_clm(rla,rlo,jdim,idim,iy,im,id,ih,tf_clm,tf_trd)
   585      tf_clm_tile(:) = reshape(tf_clm, (/lensfc/) )
   586      tf_trd_tile(:) = reshape(tf_trd, (/lensfc/) )
   590      call get_sal_clm(rla,rlo,jdim,idim,iy,im,id,ih,sal_clm)
   591      sal_clm_tile(:) = reshape(sal_clm, (/lensfc/) )
   595      CALL read_gsi_data(nst_file, 
'NST')
   599      CALL adjust_nsst(rla,rlo,slifcs,slifcs_fg,tsffcs,sitfcs,sicfcs,stcfcs, &
   600                     nsst,lensfc,lsoil,idim,jdim,zsea1,zsea2,im,id,deltsfc,  &
   601                     tf_clm_tile,tf_trd_tile,sal_clm_tile)
   623        CALL read_data(lsoil,lensfc,.false.,.true.,sndfcs=snd_inc)
   633        CALL add_increment_snow(snd_inc,landinc_mask,lensfc,sndfcs)
   639        CALL apply_land_da_adjustments_snd(lsm, lensfc, landinc_mask, swe_bck, snd_bck, &
   651         landinc_mask_fg = landinc_mask
   653         IF (do_sfccycle .OR. do_sno_inc)  
THEN   654             CALL calculate_landinc_mask(slcfcs(:,1),swefcs, vetfcs, lensfc, &
   655                                                         veg_type_landice, landinc_mask )
   662         INQUIRE(file=trim(lnd_soi_file), exist=file_exists)
   663         IF (.not. file_exists) 
then   664             print *, 
'FATAL ERROR: land increment update requested, but file does not exist: ', &
   666             call mpi_abort(mpi_comm_world, 10, ierr)
   669         CALL read_gsi_data(lnd_soi_file, 
'LND', lsoil=lsoil)
   684         CALL add_increment_soil(rla,rlo,stcfcs,landinc_mask,landinc_mask_fg,  &
   685             lensfc,lsoil,idim,jdim, myrank)
   691         CALL apply_land_da_adjustments_stc(lsm, isot, ivegsrc,lensfc, lsoil, &
   692             sotfcs, landinc_mask_fg, stc_bck, stcfcs, smcfcs, slcfcs)
   701    IF(
ALLOCATED(landinc_mask_fg)) 
DEALLOCATE(landinc_mask_fg)
   702    IF(
ALLOCATED(landinc_mask)) 
DEALLOCATE(landinc_mask)
   703    IF(
ALLOCATED(stc_bck)) 
DEALLOCATE(stc_bck)
   704    IF(
ALLOCATED(smc_bck)) 
DEALLOCATE(smc_bck)
   705    IF(
ALLOCATED(slc_bck)) 
DEALLOCATE(slc_bck)
   706    IF(
ALLOCATED(snd_bck)) 
DEALLOCATE(snd_bck)
   707    IF(
ALLOCATED(swe_bck)) 
DEALLOCATE(swe_bck)
   708    IF(
ALLOCATED(snd_inc)) 
DEALLOCATE(snd_inc)
   715  CALL write_data(slifcs,tsffcs,swefcs,tg3fcs,zorfcs,         &
   716                  albfcs,alffcs,vegfcs,cnpfcs,f10m,           &
   717                  t2m,q2m,vetfcs,sotfcs,ustar,fmm,fhh,        &
   718                  sicfcs,sihfcs,sitfcs,                       &
   719                  tprcp,srflag,sndfcs,                        &
   720                  vmnfcs,vmxfcs,slpfcs,absfcs,                &
   721                  slcfcs,smcfcs,stcfcs,                       &
   722                  idim,jdim,lensfc,lsoil,do_nsst,nsst)
   727    DEALLOCATE(nsst%D_CONV)
   728    DEALLOCATE(nsst%DT_COOL)
   730    DEALLOCATE(nsst%QRAIN)
   731    DEALLOCATE(nsst%TREF)
   732    DEALLOCATE(nsst%TFINC)
   737    DEALLOCATE(nsst%XTTS)
   741    DEALLOCATE(nsst%XZTS)
   744    DEALLOCATE(slifcs_fg)
   749  END SUBROUTINE sfcdrv
   782  SUBROUTINE adjust_nsst(RLA,RLO,SLMSK_TILE,SLMSK_FG_TILE,SKINT_TILE,&
   783                         SICET_TILE,sice_tile,SOILT_TILE,NSST,LENSFC,LSOIL,    &
   784                         IDIM,JDIM,ZSEA1,ZSEA2,MON,DAY,DELTSFC, &
   785                         tf_clm_tile,tf_trd_tile,sal_clm_tile)
   789  USE read_write_data, 
ONLY : idim_gaus, jdim_gaus, &
   790                              slmsk_gaus, dtref_gaus, &
   797  INTEGER, 
INTENT(IN)      :: LENSFC, LSOIL, IDIM, JDIM, MON, DAY
   799  REAL, 
INTENT(IN)         :: SLMSK_TILE(LENSFC), SLMSK_FG_TILE(LENSFC)
   800  real, 
intent(in)         :: tf_clm_tile(lensfc),tf_trd_tile(lensfc),sal_clm_tile(lensfc)
   801  REAL, 
INTENT(IN)         :: ZSEA1, ZSEA2, DELTSFC
   802  REAL, 
INTENT(INOUT)      :: RLA(LENSFC), RLO(LENSFC), SKINT_TILE(LENSFC)
   803  REAL, 
INTENT(INOUT)      :: SICET_TILE(LENSFC),sice_tile(lensfc),SOILT_TILE(LENSFC,LSOIL)
   805  TYPE(NSST_DATA)          :: NSST
   807  REAL, 
PARAMETER          :: TMAX=313.0,tzero=273.16
   809  INTEGER                  :: IOPT, NRET, KGDS_GAUS(200)
   810  INTEGER                  :: IGAUS, JGAUS, IJ, II, JJ, III, JJJ, KRAD
   811  INTEGER                  :: ISTART, IEND, JSTART, JEND
   812  INTEGER                  :: MASK_TILE, MASK_FG_TILE
   813  INTEGER                  :: ITILE, JTILE
   814  INTEGER                  :: MAX_SEARCH, J, IERR
   815  INTEGER                  :: IGAUSP1, JGAUSP1
   816  integer                  :: nintp,nsearched,nice,nland
   817  integer                  :: nfill,nfill_tice,nfill_clm
   818  integer                  :: nset_thaw,nset_thaw_s,nset_thaw_i,nset_thaw_c
   820  INTEGER, 
ALLOCATABLE     :: ID1(:,:), ID2(:,:), JDC(:,:)
   825  REAL                     :: TREF_SAVE,WSUM,tf_ice,tf_thaw
   826  REAL                     :: FILL, DTZM, GAUS_RES_KM, DTREF
   827  REAL, 
ALLOCATABLE        :: XPTS(:), YPTS(:), LATS(:), LONS(:)
   828  REAL, 
ALLOCATABLE        :: DUM2D(:,:), LATS_RAD(:), LONS_RAD(:)
   829  REAL, 
ALLOCATABLE        :: AGRID(:,:,:), S2C(:,:,:)
   833  kgds_gaus(2)  = idim_gaus  
   834  kgds_gaus(3)  = jdim_gaus
   838  kgds_gaus(7)  = -90000     
   839  kgds_gaus(8)  = nint(-360000./float(idim_gaus))  
   840  kgds_gaus(9)  = nint((360.0 / float(idim_gaus))*1000.0)
   842  kgds_gaus(10) = jdim_gaus/2     
   847  print*,
'ADJUST NSST USING GSI INCREMENTS ON GAUSSIAN GRID'   855  ALLOCATE(xpts(idim_gaus*jdim_gaus))
   856  ALLOCATE(ypts(idim_gaus*jdim_gaus))
   857  ALLOCATE(lats(idim_gaus*jdim_gaus))
   858  ALLOCATE(lons(idim_gaus*jdim_gaus))
   864  CALL gdswzd(kgds_gaus,iopt,(idim_gaus*jdim_gaus),fill,xpts,ypts,lons,lats,nret)
   866  IF (nret /= (idim_gaus*jdim_gaus)) 
THEN   867    print*,
'FATAL ERROR: PROBLEM IN GDSWZD. STOP.'   868    CALL mpi_abort(mpi_comm_world, 12, ierr)
   871  DEALLOCATE (xpts, ypts)
   873  ALLOCATE(dum2d(idim_gaus,jdim_gaus))
   874  dum2d = reshape(lats, (/idim_gaus,jdim_gaus/) )
   877  ALLOCATE(lats_rad(jdim_gaus))
   879    lats_rad(j) = dum2d(1,jdim_gaus-j+1) * 3.1415926 / 180.0
   882  dum2d = reshape(lons, (/idim_gaus,jdim_gaus/) )
   884  ALLOCATE(lons_rad(idim_gaus))
   885  lons_rad = dum2d(:,1) * 3.1415926 / 180.0
   889  ALLOCATE(agrid(idim,jdim,2))
   890  agrid(:,:,1) = reshape(rlo, (/idim,jdim/) )
   891  agrid(:,:,2) = reshape(rla, (/idim,jdim/) )
   892  agrid        = agrid * 3.1415926 / 180.0
   894  ALLOCATE(id1(idim,jdim))
   895  ALLOCATE(id2(idim,jdim))
   896  ALLOCATE(jdc(idim,jdim))
   897  ALLOCATE(s2c(idim,jdim,4))
   905  CALL remap_coef( 1, idim, 1, jdim, idim_gaus, jdim_gaus, &
   906                   lons_rad, lats_rad, id1, id2, jdc, s2c, agrid )
   908  DEALLOCATE(lons_rad, lats_rad, agrid)
   915  gaus_res_km = 360.0 / idim_gaus * 111.0
   916  max_search  = ceiling(500.0/gaus_res_km)
   919  print*,
'MAXIMUM SEARCH IS ',max_search, 
' GAUSSIAN POINTS.'   942  ij_loop : 
DO ij = 1, lensfc
   944    mask_tile    = nint(slmsk_tile(ij))
   945    mask_fg_tile = nint(slmsk_fg_tile(ij))
   950    tf_ice = tfreez(sal_clm_tile(ij)) + tzero
   955    IF (mask_tile == 1) 
THEN   963    if (mask_tile == 2) 
then   965      skint_tile(ij)=(1.0-sice_tile(ij))*nsst%tref(ij)+sice_tile(ij)*sicet_tile(ij)
   973    jtile = (ij-1) / idim + 1
   975    IF (itile==0) itile = idim
   983    IF (mask_fg_tile == 2 .AND. mask_tile == 0) 
THEN   987      call tf_thaw_set(nsst%tref,nint(slmsk_fg_tile),itile,jtile,tf_ice,tf_clm_tile(ij),tf_thaw,idim,jdim, &
   988                       nset_thaw_s,nset_thaw_i,nset_thaw_c)
   990      nset_thaw = nset_thaw + 1
  1006    igaus   = id1(itile,jtile)
  1007    jgaus   = jdc(itile,jtile)
  1008    igausp1 = id2(itile,jtile)
  1009    jgausp1 = jdc(itile,jtile)+1
  1011    IF (slmsk_gaus(igaus,jgaus)     == 0 .OR. &
  1012        slmsk_gaus(igausp1,jgaus)   == 0 .OR. &
  1013        slmsk_gaus(igausp1,jgausp1) == 0 .OR. &
  1014        slmsk_gaus(igaus,jgausp1)   == 0) 
THEN  1019      IF (slmsk_gaus(igaus,jgaus) == 0) 
THEN  1020        dtref = dtref + (s2c(itile,jtile,1) * dtref_gaus(igaus,jgaus))
  1021        wsum  = wsum + s2c(itile,jtile,1)
  1024      IF (slmsk_gaus(igausp1,jgaus) == 0) 
THEN  1025        dtref = dtref + (s2c(itile,jtile,2) * dtref_gaus(igausp1,jgaus))
  1026        wsum  = wsum + s2c(itile,jtile,2)
  1029      IF (slmsk_gaus(igausp1,jgausp1) == 0) 
THEN  1030        dtref = dtref + (s2c(itile,jtile,3) * dtref_gaus(igausp1,jgausp1))
  1031        wsum  = wsum + s2c(itile,jtile,3)
  1034      IF (slmsk_gaus(igaus,jgausp1) == 0) 
THEN  1035        dtref = dtref + (s2c(itile,jtile,4) * dtref_gaus(igaus,jgausp1))
  1036        wsum  = wsum + s2c(itile,jtile,4)
  1040      dtref = dtref / wsum
  1042      tref_save      = nsst%TREF(ij)
  1043      nsst%TREF(ij)  = nsst%TREF(ij) + dtref
  1044      nsst%TREF(ij)  = max(nsst%TREF(ij), tf_ice)
  1045      nsst%TREF(ij)  = min(nsst%TREF(ij), tmax)
  1046      nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
  1048      CALL dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij),  &
  1049                      nsst%Z_C(ij),zsea1,zsea2,dtzm)
  1051      skint_tile(ij) = nsst%TREF(ij) + dtzm
  1052      skint_tile(ij) = max(skint_tile(ij), tf_ice)
  1053      skint_tile(ij) = min(skint_tile(ij), tmax)
  1055      sicet_tile(ij)   = skint_tile(ij)
  1056      soilt_tile(ij,:) = skint_tile(ij)
  1067      DO krad = 1, max_search
  1069        istart = igaus - krad
  1071        jstart = jgaus - krad
  1074        DO jj = jstart, jend
  1075        DO ii = istart, iend
  1077          IF((jj == jstart) .OR. (jj == jend) .OR.   &
  1078             (ii == istart) .OR. (ii == iend))  
THEN  1080            IF ((jj >= 1) .AND. (jj <= jdim_gaus)) 
THEN  1084                iii = idim_gaus + ii
  1085              ELSE IF (ii >= (idim_gaus+1)) 
THEN  1086                iii = ii - idim_gaus
  1097              IF (krad <= 2 .AND. slmsk_gaus(iii,jjj) == 2) is_ice = .true.
  1099              IF (slmsk_gaus(iii,jjj) == 0) 
THEN  1103                nsearched = nsearched + 1
  1105                tref_save      = nsst%TREF(ij)
  1106                nsst%TREF(ij ) = nsst%TREF(ij) + dtref_gaus(iii,jjj)
  1107                nsst%TREF(ij)  = max(nsst%TREF(ij), tf_ice)
  1108                nsst%TREF(ij)  = min(nsst%TREF(ij), tmax)
  1109                nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
  1111                CALL dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij),  &
  1112                      nsst%Z_C(ij),zsea1,zsea2,dtzm)
  1114                skint_tile(ij) = nsst%TREF(ij) + dtzm
  1115                skint_tile(ij) = max(skint_tile(ij), tf_ice)
  1116                skint_tile(ij) = min(skint_tile(ij), tmax)
  1118                sicet_tile(ij)   = skint_tile(ij)
  1119                soilt_tile(ij,:) = skint_tile(ij)
  1142        nsst%TREF(ij) = tf_ice
  1144        nfill_tice = nfill_tice + 1
  1146        tref_save      = nsst%TREF(ij)
  1147        nsst%TREF(ij)  = nsst%TREF(ij) + tf_trd_tile(ij)
  1148        nsst%TREF(ij)  = max(nsst%TREF(ij), tf_ice)
  1149        nsst%TREF(ij)  = min(nsst%TREF(ij), tmax)
  1150        nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
  1152        nfill_clm = nfill_clm + 1
  1155      CALL dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij),  &
  1156                      nsst%Z_C(ij),zsea1,zsea2,dtzm)
  1158      skint_tile(ij) = nsst%TREF(ij) + dtzm
  1159      skint_tile(ij) = max(skint_tile(ij), tf_ice)
  1160      skint_tile(ij) = min(skint_tile(ij), tmax)
  1162      sicet_tile(ij)   = skint_tile(ij)
  1163      soilt_tile(ij,:) = skint_tile(ij)
  1169  write(*,
'(a)') 
'statistics of grids number processed for tile : '  1170  write(*,
'(a,I8)') 
' nintp = ',nintp
  1171  write(*,
'(a,4I8)') 
'nset_thaw,nset_thaw_s,nset_thaw_i,nset_thaw_c =',nset_thaw,nset_thaw_s,nset_thaw_i,nset_thaw_c
  1172  write(*,
'(a,I8)') 
' nsearched = ',nsearched
  1173  write(*,
'(a,3I6)') 
' nfill,nfill_tice,nfill_clm = ',nfill,nfill_tice,nfill_clm
  1174  write(*,
'(a,I8)') 
' nice = ',nice
  1175  write(*,
'(a,I8)') 
' nland = ',nland
  1177  DEALLOCATE(id1, id2, jdc, s2c)
  1179  END SUBROUTINE adjust_nsst
  1191  SUBROUTINE climo_trend(LATITUDE, MON, DAY, DELTSFC, DTREF)
  1194  INTEGER, 
INTENT(IN)    :: MON, DAY
  1196  REAL, 
INTENT(IN)       :: LATITUDE, DELTSFC
  1197  REAL, 
INTENT(OUT)      :: DTREF
  1199  INTEGER                :: NUM_DAYS(12), MON2, MON1
  1201  REAL, 
TARGET           :: SST_80_90(12)
  1202  REAL, 
TARGET           :: SST_70_80(12)
  1203  REAL, 
TARGET           :: SST_60_70(12)
  1204  REAL, 
TARGET           :: SST_50_60(12)
  1205  REAL, 
TARGET           :: SST_40_50(12)
  1206  REAL, 
TARGET           :: SST_30_40(12)
  1207  REAL, 
TARGET           :: SST_20_30(12)
  1208  REAL, 
TARGET           :: SST_10_20(12)
  1209  REAL, 
TARGET           :: SST_00_10(12)
  1210  REAL, 
TARGET           :: SST_M10_00(12)
  1211  REAL, 
TARGET           :: SST_M20_M10(12)
  1212  REAL, 
TARGET           :: SST_M30_M20(12)
  1213  REAL, 
TARGET           :: SST_M40_M30(12)
  1214  REAL, 
TARGET           :: SST_M50_M40(12)
  1215  REAL, 
TARGET           :: SST_M60_M50(12)
  1216  REAL, 
TARGET           :: SST_M70_M60(12)
  1217  REAL, 
TARGET           :: SST_M80_M70(12)
  1218  REAL, 
TARGET           :: SST_M90_M80(12)
  1220  REAL, 
POINTER          :: SST(:)
  1222  DATA num_days  /31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31/
  1224  DATA sst_80_90 /271.466, 271.458, 271.448, 271.445, 271.519, 271.636, & 
  1225                  272.023, 272.066, 272.001, 271.698, 271.510, 271.472/
  1227  DATA sst_70_80 /272.149, 272.103, 272.095, 272.126, 272.360, 272.988, &
  1228                  274.061, 274.868, 274.415, 273.201, 272.468, 272.268/
  1230  DATA sst_60_70 /274.240, 274.019, 273.988, 274.185, 275.104, 276.875, &
  1231                  279.005, 280.172, 279.396, 277.586, 275.818, 274.803/
  1233  DATA sst_50_60 /277.277, 276.935, 277.021, 277.531, 279.100, 281.357, &
  1234                  283.735, 285.171, 284.399, 282.328, 279.918, 278.199/
  1236  DATA sst_40_50 /281.321, 280.721, 280.850, 281.820, 283.958, 286.588, &
  1237                  289.195, 290.873, 290.014, 287.652, 284.898, 282.735/
  1239  DATA sst_30_40 /289.189, 288.519, 288.687, 289.648, 291.547, 293.904, &
  1240                  296.110, 297.319, 296.816, 295.225, 292.908, 290.743/
  1242  DATA sst_20_30 /294.807, 294.348, 294.710, 295.714, 297.224, 298.703, &
  1243                  299.682, 300.127, 300.099, 299.455, 297.953, 296.177/
  1245  DATA sst_10_20 /298.878, 298.720, 299.033, 299.707, 300.431, 300.709, &
  1246                  300.814, 300.976, 301.174, 301.145, 300.587, 299.694/
  1248  DATA sst_00_10 /300.415, 300.548, 300.939, 301.365, 301.505, 301.141, &
  1249                  300.779, 300.660, 300.818, 300.994, 300.941, 300.675/
  1251  DATA sst_m10_00 /300.226, 300.558, 300.914, 301.047, 300.645, 299.870, &
  1252                   299.114, 298.751, 298.875, 299.294, 299.721, 299.989/
  1254  DATA sst_m20_m10 /299.547, 299.985, 300.056, 299.676, 298.841, 297.788, &
  1255                    296.893, 296.491, 296.687, 297.355, 298.220, 298.964/
  1257  DATA sst_m30_m20 /297.524, 298.073, 297.897, 297.088, 295.846, 294.520, &
  1258                    293.525, 293.087, 293.217, 293.951, 295.047, 296.363/
  1260  DATA sst_m40_m30 /293.054, 293.765, 293.468, 292.447, 291.128, 289.781, &
  1261                    288.773, 288.239, 288.203, 288.794, 289.947, 291.553/
  1263  DATA sst_m50_m40 /285.052, 285.599, 285.426, 284.681, 283.761, 282.826, &
  1264                    282.138, 281.730, 281.659, 281.965, 282.768, 283.961/
  1266  DATA sst_m60_m50 /277.818, 278.174, 277.991, 277.455, 276.824, 276.229, &
  1267                    275.817, 275.585, 275.560, 275.687, 276.142, 276.968/
  1269  DATA sst_m70_m60 /273.436, 273.793, 273.451, 272.813, 272.349, 272.048, &
  1270                    271.901, 271.838, 271.845, 271.889, 272.080, 272.607/
  1272  DATA sst_m80_m70 /271.579, 271.578, 271.471, 271.407, 271.392, 271.391, &
  1273                    271.390, 271.391, 271.394, 271.401, 271.422, 271.486/
  1275  DATA sst_m90_m80 /271.350, 271.350, 271.350, 271.350, 271.350, 271.350, &
  1276                    271.350, 271.350, 271.350, 271.350, 271.350, 271.350/
  1279  IF (latitude > 80.0) 
THEN  1281  ELSEIF (latitude > 70.0) 
THEN  1283  ELSEIF (latitude > 60.0) 
THEN  1285  ELSEIF (latitude > 50.0) 
THEN  1287  ELSEIF (latitude > 40.0) 
THEN  1289  ELSEIF (latitude > 30.0) 
THEN  1291  ELSEIF (latitude > 20.0) 
THEN  1293  ELSEIF (latitude > 10.0) 
THEN  1295  ELSEIF (latitude > 0.0) 
THEN  1297  ELSEIF (latitude > -10.0) 
THEN  1299  ELSEIF (latitude > -20.0) 
THEN  1301  ELSEIF (latitude > -30.0) 
THEN  1303  ELSEIF (latitude > -40.0) 
THEN  1305  ELSEIF (latitude > -50.0) 
THEN  1307  ELSEIF (latitude > -60.0) 
THEN  1309  ELSEIF (latitude > -70.0) 
THEN  1311  ELSEIF (latitude > -80.0) 
THEN  1319    IF(mon2 == 13) mon2 = 1
  1321    dtref = (sst(mon2) - sst(mon1)) / num_days(mon1)
  1324    IF (mon1 == 0) mon1=12
  1326    dtref = (sst(mon2) - sst(mon1)) / num_days(mon1)
  1329  dtref = dtref * (deltsfc / 24.0)
  1331  END SUBROUTINE climo_trend
  1344  SUBROUTINE dtzm_point(XT,XZ,DT_COOL,ZC,Z1,Z2,DTZM)
  1347   real, 
intent(in)  :: xt,xz,dt_cool,zc,z1,z2
  1348   real, 
intent(out) :: dtzm
  1350   real, 
parameter   :: zero = 0.0
  1351   real, 
parameter   :: one  = 1.0
  1352   real, 
parameter   :: half = 0.5
  1353   real              :: dt_warm,dtw,dtc
  1358   if ( xt > zero ) 
then  1359     dt_warm = (xt+xt)/xz      
  1362         dtw = dt_warm*(one-(z1+z2)/(xz+xz))
  1363       elseif ( z1 < xz .and. z2 >= xz ) 
then  1364         dtw = half*(one-z1/xz)*dt_warm*(xz-z1)/(z2-z1)
  1366     elseif ( z1 == z2 ) 
then  1368         dtw = dt_warm*(one-z1/xz)
  1376   if ( zc > zero ) 
then  1379         dtc = dt_cool*(one-(z1+z2)/(zc+zc))
  1380       elseif ( z1 < zc .and. z2 >= zc ) 
then  1381         dtc = half*(one-z1/zc)*dt_cool*(zc-z1)/(z2-z1)
  1383     elseif ( z1 == z2 ) 
then  1385         dtc = dt_cool*(one-z1/zc)
  1395  END SUBROUTINE dtzm_point
  1416  subroutine tf_thaw_set(tf_ij,mask_ij,itile,jtile,tice,tclm,tf_thaw,nx,ny, &
  1417                         nset_thaw_s,nset_thaw_i,nset_thaw_c)
  1419  real,    
dimension(nx*ny), 
intent(in)    :: tf_ij
  1420  integer, 
dimension(nx*ny), 
intent(in)    :: mask_ij
  1421  real,                      
intent(in)    :: tice,tclm
  1422  integer,                   
intent(in)    :: itile,jtile,nx,ny
  1423  real,                      
intent(out)   :: tf_thaw
  1424  integer,                   
intent(inout) :: nset_thaw_s,nset_thaw_i,nset_thaw_c
  1426  real, 
parameter :: bmiss = -999.0
  1427  real,    
dimension(nx,ny) :: tf
  1428  integer, 
dimension(nx,ny) :: mask
  1429  integer :: krad,max_search,istart,iend,jstart,jend
  1430  integer :: ii,jj,iii,jjj
  1435  mask(:,:) = reshape(mask_ij,(/nx,ny/) )
  1436  tf(:,:)   = reshape(tf_ij,(/nx,ny/) )
  1440  do krad = 1, max_search
  1442     istart = itile - krad
  1444     jstart = jtile - krad
  1447     do jj = jstart, jend
  1448        do ii = istart, iend
  1451           if ((jj == jstart) .or. (jj == jend) .or.   &
  1452               (ii == istart) .or. (ii == iend))  
then  1454              if ((jj >= 1) .and. (jj <= ny)) 
then  1458                 else if (ii >= (nx+1)) 
then  1469                 if (krad <= 2 .and. mask(iii,jjj) == 2) is_ice = .true.
  1471                 if (mask(iii,jjj) == 0) 
then  1472                    tf_thaw = tf(iii,jjj)
  1473                    nset_thaw_s = nset_thaw_s + 1
  1474                    write(*,
'(a,I4,2F9.3)') 
'nset_thaw_s,tf(iii,jjj),tclm : ',nset_thaw_s,tf(iii,jjj),tclm
  1486  if ( tf_thaw == bmiss ) 
then  1489        nset_thaw_i = nset_thaw_i + 1
  1490        write(*,
'(a,I4,F9.3)') 
'nset_thaw_i,tf_ice : ',nset_thaw_i,tice
  1492        tf_thaw = 0.8*tice+0.2*tclm
  1493        nset_thaw_c = nset_thaw_c + 1
  1494        write(*,
'(a,I4,2F9.3)') 
'nset_thaw_c,tf_ice,tclm : ',nset_thaw_c,tice,tclm
  1508  use read_write_data, 
only : nsst_data
  1511  integer, 
intent(in)  :: ij
  1513  real, 
intent(in)     :: tf_thaw
  1515  type(nsst_data), 
intent(inout)      :: nsst
  1519  nsst%d_conv(ij)  = 0.0
  1520  nsst%dt_cool(ij) = 0.0
  1522  nsst%qrain(ij)   = 0.0
  1523  nsst%tref(ij)    = tf_thaw
  1553 subroutine get_tf_clm(xlats_ij,xlons_ij,ny,nx,iy,im,id,ih,tf_clm,tf_trd)
  1554  use read_write_data, 
only : get_tf_clm_dim
  1558  real,    
dimension(nx*ny), 
intent(in)  :: xlats_ij 
  1559  real,    
dimension(nx*ny), 
intent(in)  :: xlons_ij 
  1560  real,    
dimension(nx,ny), 
intent(out) :: tf_clm   
  1561  real,    
dimension(nx,ny), 
intent(out) :: tf_trd   
  1562  integer, 
intent(in) :: iy,im,id,ih,nx,ny
  1564  real,    
allocatable, 
dimension(:,:)   :: tf_clm0    
  1565  real,    
allocatable, 
dimension(:,:)   :: tf_trd0    
  1566  real,    
allocatable, 
dimension(:)     :: cxlats     
  1567  real,    
allocatable, 
dimension(:)     :: cxlons     
  1569  real,    
dimension(nx*ny)  :: tf_clm_ij  
  1570  real,    
dimension(nx*ny)  :: tf_trd_ij  
  1572  integer :: nxc,nyc,mon1,mon2,i,j
  1573  character (len=6), 
parameter :: fin_tf_clm=
'sstclm'   1581  call get_tf_clm_dim(fin_tf_clm,nyc,nxc)
  1582  allocate( tf_clm0(nxc,nyc),tf_trd0(nxc,nyc),cxlats(nyc),cxlons(nxc) )
  1586  call get_tf_clm_ta(tf_clm0,tf_trd0,cxlats,cxlons,nyc,nxc,mon1,mon2,wei1,wei2)
  1590  if ( nx == nxc .and. ny == nyc ) 
then  1591     tf_clm(:,:) = tf_clm0(:,:)
  1592     tf_trd(:,:) = tf_trd0(:,:)
  1596     call intp_tile(tf_clm0,  cxlats,  cxlons,  nyc, nxc, &
  1597                    tf_clm_ij,xlats_ij,xlons_ij,ny,  nx)
  1598     call intp_tile(tf_trd0,  cxlats,  cxlons,  nyc, nxc, &
  1599                    tf_trd_ij,xlats_ij,xlons_ij,ny,  nx)
  1602     tf_clm(:,:) = reshape(tf_clm_ij, (/nx,ny/) )
  1603     tf_trd(:,:) = reshape(tf_trd_ij, (/nx,ny/) )
  1622 subroutine get_tf_clm_ta(tf_clm_ta,tf_clm_trend,xlats,xlons,nlat,nlon,mon1,mon2,wei1,wei2)
  1623  use read_write_data, 
only : read_tf_clim_grb
  1627  integer, 
intent(in) :: nlat,nlon,mon1,mon2
  1628  real,    
intent(in) :: wei1,wei2
  1630  real, 
dimension(nlon,nlat), 
intent(out)   :: tf_clm_ta,tf_clm_trend
  1631  real, 
dimension(nlat),      
intent(out)   :: xlats
  1632  real, 
dimension(nlon),      
intent(out)   :: xlons
  1635  character (len=6),  
parameter :: fin_tf_clm=
'sstclm'  1638  real, 
dimension(nlon,nlat) :: tf_clm1,tf_clm2
  1643   call read_tf_clim_grb(trim(fin_tf_clm),tf_clm1,xlats,xlons,nlat,nlon,mon1)
  1644   call read_tf_clim_grb(trim(fin_tf_clm),tf_clm2,xlats,xlons,nlat,nlon,mon2)
  1648    tf_clm_ta(:,:) = wei1*tf_clm1(:,:)+wei2*tf_clm2(:,:)
  1652    tf_clm_trend(:,:) = (tf_clm2(:,:)-tf_clm1(:,:))/120.0
  1654    write(*,
'(a,2f9.3)') 
'tf_clm_ta, min, max : ',minval(tf_clm_ta),maxval(tf_clm_ta)
  1655    write(*,
'(a,2f9.3)') 
'tf_clm_trend, min, max : ',minval(tf_clm_trend),maxval(tf_clm_trend)
  1670 subroutine get_sal_clm(xlats_ij,xlons_ij,ny,nx,iy,im,id,ih,sal_clm)
  1671  use read_write_data, 
only : get_dim_nc
  1674  real,    
dimension(nx*ny), 
intent(in)  :: xlats_ij   
  1675  real,    
dimension(nx*ny), 
intent(in)  :: xlons_ij   
  1676  real,    
dimension(nx,ny), 
intent(out) :: sal_clm    
  1677  integer, 
intent(in) :: iy,im,id,ih,nx,ny
  1679  real,    
allocatable, 
dimension(:,:)   :: sal_clm0   
  1680  real,    
allocatable, 
dimension(:)     :: cxlats     
  1681  real,    
allocatable, 
dimension(:)     :: cxlons     
  1683  real,    
dimension(nx*ny)  :: sal_clm_ij  
  1685  integer :: nxc,nyc,mon1,mon2,i,j
  1686  character (len=6), 
parameter :: fin_sal_clm=
'salclm'   1694  call get_dim_nc(fin_sal_clm,nyc,nxc)
  1695  allocate( sal_clm0(nxc,nyc),cxlats(nyc),cxlons(nxc) )
  1699  call get_sal_clm_ta(sal_clm0,cxlats,cxlons,nyc,nxc,mon1,mon2,wei1,wei2)
  1703  if ( nx == nxc .and. ny == nyc ) 
then  1704     sal_clm(:,:) = sal_clm0(:,:)
  1708     call intp_tile(sal_clm0,  cxlats,  cxlons,  nyc, nxc, &
  1709                    sal_clm_ij,xlats_ij,xlons_ij,ny,  nx)
  1713     sal_clm(:,:) = reshape(sal_clm_ij, (/nx,ny/) )
  1730 subroutine get_sal_clm_ta(sal_clm_ta,xlats,xlons,nlat,nlon,mon1,mon2,wei1,wei2)
  1732  use read_write_data, 
only : read_salclm_gfs_nc
  1736  integer, 
intent(in) :: nlat,nlon,mon1,mon2
  1737  real,    
intent(in) :: wei1,wei2
  1739  real, 
dimension(nlon,nlat), 
intent(out)   :: sal_clm_ta
  1740  real, 
dimension(nlat),      
intent(out)   :: xlats
  1741  real, 
dimension(nlon),      
intent(out)   :: xlons
  1744  character (len=6),  
parameter :: fin_sal_clm=
'salclm'  1747  real, 
dimension(nlon,nlat) :: sal_clm1,sal_clm2
  1752   call read_salclm_gfs_nc(trim(fin_sal_clm),sal_clm1,xlats,xlons,nlat,nlon,mon1)
  1753   call read_salclm_gfs_nc(trim(fin_sal_clm),sal_clm2,xlats,xlons,nlat,nlon,mon2)
  1757    sal_clm_ta(:,:) = wei1*sal_clm1(:,:)+wei2*sal_clm2(:,:)
  1758    write(*,
'(a,2f9.3)') 
'sal_clm_ta, min, max : ',minval(sal_clm_ta),maxval(sal_clm_ta)
  1775 subroutine intp_tile(tf_lalo,dlats_lalo,dlons_lalo,jdim_lalo,idim_lalo, &
  1776                      tf_tile,xlats_tile,xlons_tile,jdim_tile,idim_tile)
  1783  real,    
dimension(idim_lalo,jdim_lalo), 
intent(in)  :: tf_lalo
  1784  real,    
dimension(jdim_lalo),           
intent(in)  :: dlats_lalo
  1785  real,    
dimension(idim_lalo),           
intent(in)  :: dlons_lalo
  1786  real,    
dimension(jdim_tile*idim_tile), 
intent(in)  :: xlats_tile
  1787  real,    
dimension(jdim_tile*idim_tile), 
intent(in)  :: xlons_tile
  1788  integer,                                 
intent(in)  :: jdim_lalo,idim_lalo,jdim_tile,idim_tile
  1789  real,    
dimension(jdim_tile*idim_tile), 
intent(out) :: tf_tile
  1792  real, 
parameter :: deg2rad=3.1415926/180.0
  1793  real,    
dimension(jdim_lalo) :: xlats_lalo
  1794  real,    
dimension(idim_lalo) :: xlons_lalo
  1795  real    :: tf,wsum,res_km
  1796  integer :: itile,jtile
  1797  integer :: ii,jj,ij,iii,jjj
  1798  integer :: ilalo,jlalo,ilalop1,jlalop1
  1799  integer :: istart,iend,jstart,jend,krad
  1801  integer, 
allocatable, 
dimension(:,:)   :: id1,id2,jdc
  1802  real,    
allocatable, 
dimension(:,:,:) :: agrid,s2c
  1805  print*,
'interpolate from lat/lon grids to any one grid with known lat/lon'  1807  xlats_lalo = dlats_lalo*deg2rad
  1808  xlons_lalo = dlons_lalo*deg2rad
  1810  allocate(agrid(idim_tile,jdim_tile,2))
  1811  agrid(:,:,1) = reshape(xlons_tile, (/idim_tile,jdim_tile/) )
  1812  agrid(:,:,2) = reshape(xlats_tile, (/idim_tile,jdim_tile/) )
  1813  agrid        = agrid*deg2rad
  1815  allocate(id1(idim_tile,jdim_tile))
  1816  allocate(id2(idim_tile,jdim_tile))
  1817  allocate(jdc(idim_tile,jdim_tile))
  1818  allocate(s2c(idim_tile,jdim_tile,4))
  1826  call remap_coef( 1, idim_tile, 1, jdim_tile, idim_lalo, jdim_lalo, &
  1827                   xlons_lalo, xlats_lalo, id1, id2, jdc, s2c, agrid )
  1829  do ij = 1, jdim_tile*idim_tile
  1831     jtile = (ij-1)/idim_tile + 1
  1832     itile = mod(ij,idim_tile)
  1833     if (itile==0) itile = idim_tile
  1835     ilalo   = id1(itile,jtile)
  1836     ilalop1 = id2(itile,jtile)
  1837     jlalo   = jdc(itile,jtile)
  1838     jlalop1 = jdc(itile,jtile) + 1
  1840     wsum = s2c(itile,jtile,1) + s2c(itile,jtile,2) + &
  1841            s2c(itile,jtile,3) + s2c(itile,jtile,4)
  1843     tf_tile(ij)  = ( s2c(itile,jtile,1)*tf_lalo(ilalo,jlalo)     + &
  1844                      s2c(itile,jtile,2)*tf_lalo(ilalop1,jlalo)   + &
  1845                      s2c(itile,jtile,3)*tf_lalo(ilalop1,jlalop1) + &
  1846                      s2c(itile,jtile,4)*tf_lalo(ilalo,jlalop1) )/wsum
  1849  deallocate(id1, id2, jdc, s2c)
  1865 subroutine get_tim_wei(iy,im,id,ih,mon1,mon2,wei1,wei2)
  1869  integer, 
intent(in) :: iy,im,id,ih
  1871  integer, 
intent(out) :: mon1,mon2
  1872  real,    
intent(out) :: wei1,wei2
  1876  integer :: mon,monend,monm,monp,jdow,jdoy,jday
  1881  real, 
dimension(13) ::  dayhf
  1882  data dayhf/15.5,45.0,74.5,105.0,135.5,166.0,196.5,227.5,258.0,288.5,319.0,349.5,380.5/
  1895  call w3doxdat(jda,jdow,jdoy,jday)
  1896  rjday=jdoy+jda(5)/24.
  1897  if(rjday.lt.dayhf(1)) rjday=rjday+365.
  1903     if( rjday >= dayhf(monm) .and. rjday < dayhf(monp) ) 
then  1910  print *,
'FATAL ERROR in get_tim_wei, wrong rjday',rjday
  1914  wei1 = (dayhf(mon2)-rjday)/(dayhf(mon2)-dayhf(mon1))
  1915  wei2 = (rjday-dayhf(mon1))/(dayhf(mon2)-dayhf(mon1))
  1917  if( mon2 == 13 ) mon2=1
  1919  write(*,
'(a,2i4,3f9.3)') 
'mon1,mon2,rjday,wei1,wei2=',mon1,mon2,rjday,wei1,wei2
  1932  real function tfreez(salinity)
  1938  parameter(a1 = -0.0575)
  1939  parameter(a2 =  1.710523e-3)
  1940  parameter(a3 = -2.154996e-4)
  1942  IF (salinity .LT. 0.) 
THEN  1947  tfreez = sal*(a1+a2*sqrt(sal)+a3*sal)
 integer function num_parthds()
Return the number of omp threads. 
 
subroutine get_tf_clm(xlats_ij, xlons_ij, ny, nx, iy, im, id, ih, tf_clm, tf_trd)
Get the sst climatology at the valid time and on the target grid. 
 
subroutine get_tf_clm_ta(tf_clm_ta, tf_clm_trend, xlats, xlons, nlat, nlon, mon1, mon2, wei1, wei2)
Get the reference temperature/sst climatology and its trend at analysis time. 
 
subroutine get_sal_clm_ta(sal_clm_ta, xlats, xlons, nlat, nlon, mon1, mon2, wei1, wei2)
Get climatological salinity at the analysis time. 
 
subroutine tf_thaw_set(tf_ij, mask_ij, itile, jtile, tice, tclm, tf_thaw, nx, ny, nset_thaw_s, nset_thaw_i, nset_thaw_c)
Set the background reference temperature (tf) for points where the ice has just melted. 
 
subroutine intp_tile(tf_lalo, dlats_lalo, dlons_lalo, jdim_lalo, idim_lalo, tf_tile, xlats_tile, xlons_tile, jdim_tile, idim_tile)
Interpolate lon/lat grid data to the fv3 native grid (tf_lalo => tf_tile). 
 
Module containing utility routines. 
 
real function tfreez(salinity)
Compute the freezing point of water as a function of salinity. 
 
subroutine get_tim_wei(iy, im, id, ih, mon1, mon2, wei1, wei2)
For a given date, determine the bounding months and the linear time interpolation weights...
 
subroutine get_sal_clm(xlats_ij, xlons_ij, ny, nx, iy, im, id, ih, sal_clm)
Get salinity climatology at the valid time on the target grid. 
 
subroutine nsst_water_reset(nsst, ij, tf_thaw)
If the first guess was sea ice, but the analysis is open water, reset all nsst variables.