global_cycle  1.10.0
cycle.f90
Go to the documentation of this file.
1 
4 
106  PROGRAM sfc_drv
107 
108  use mpi
109 
110  IMPLICIT NONE
111 !
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
118 !
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
123 !
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/
128 !
129  CALL mpi_init(ierr)
130  CALL mpi_comm_size(mpi_comm_world, nprocs, ierr)
131  CALL mpi_comm_rank(mpi_comm_world, myrank, ierr)
132 
133  if (myrank==0) call w3tagb('GLOBAL_CYCLE',2018,0179,0055,'NP20')
134 
135  num_threads = num_parthds()
136 
137  print*
138  print*,"STARTING CYCLE PROGRAM ON RANK ", myrank
139  print*,"RUNNING WITH ", nprocs, "TASKS"
140  print*,"AND WITH ", num_threads, " THREADS."
141 
142  use_ufo = .false.
143  donst = "NO"
144  do_lndinc = .false.
145  do_sfccycle = .true.
146 
147  print*
148  print*,"READ NAMCYC NAMELIST."
149 
150  CALL baopenr(36, "fort.36", ierr)
151  READ(36, nml=namcyc)
152  !IF (MYRANK==0) WRITE(6,NAMCYC)
153 
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
157  GOTO 333
158  ENDIF
159 
160  lensfc = idim*jdim ! TOTAL NUMBER OF POINTS FOR THE CUBED-SPHERE TILE
161 
162  zsea1 = float(zsea1_mm) / 1000.0 ! CONVERT FROM MM TO METERS
163  zsea2 = float(zsea2_mm) / 1000.0
164 
165  IF (donst == "YES") THEN
166  do_nsst=.true.
167  ELSE
168  do_nsst=.false.
169  ENDIF
170 
171  print*
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
174 
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)
179 
180  print*
181  print*,'CYCLE PROGRAM COMPLETED NORMALLY ON RANK: ', myrank
182 
183  333 CONTINUE
184 
185  CALL mpi_barrier(mpi_comm_world, ierr)
186 
187  if (myrank==0) call w3tage('GLOBAL_CYCLE')
188 
189  CALL mpi_finalize(ierr)
190 
191  stop
192 
193  END PROGRAM sfc_drv
194 
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)
308 !
309  USE read_write_data
310  use machine
311  USE mpi
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
317 
318  IMPLICIT NONE
319 
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
323 
324  LOGICAL, INTENT(IN) :: USE_UFO, DO_NSST,DO_SFCCYCLE
325  LOGICAL, INTENT(IN) :: DO_LNDINC
326 
327  REAL, INTENT(IN) :: FH, DELTSFC, ZSEA1, ZSEA2
328 
329  INTEGER, PARAMETER :: NLUNIT=35
330  INTEGER, PARAMETER :: SZ_NML=1
331 
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)
336 
337  INTEGER :: I, IERR
338  INTEGER :: I_INDEX(LENSFC), J_INDEX(LENSFC)
339  INTEGER :: IDUM(IDIM,JDIM)
340  integer :: num_parthds, num_threads
341 
342  LOGICAL :: IS_NOAHMP=.false.
343 
344  real(kind=kind_io8) :: min_ice(lensfc)
345 
346  REAL :: SLMASK(LENSFC), OROG(LENSFC)
347  REAL :: SIHFCS(LENSFC), SICFCS(LENSFC)
348  REAL :: SITFCS(LENSFC), TSFFCS(LENSFC)
349  REAL :: SWEFCS(LENSFC), ZORFCS(LENSFC)
350  REAL :: ALBFCS(LENSFC,4), TG3FCS(LENSFC)
351  REAL :: CNPFCS(LENSFC), SMCFCS(LENSFC,LSOIL)
352  REAL :: STCFCS(LENSFC,LSOIL), SLIFCS(LENSFC)
353  REAL :: AISFCS(LENSFC), F10M(LENSFC)
354  REAL :: VEGFCS(LENSFC), VETFCS(LENSFC)
355  REAL :: SOTFCS(LENSFC), ALFFCS(LENSFC,2)
356  REAL :: CVFCS(LENSFC), CVTFCS(LENSFC)
357  REAL :: CVBFCS(LENSFC), TPRCP(LENSFC)
358  REAL :: SRFLAG(LENSFC), SNDFCS(LENSFC)
359  REAL :: SLCFCS(LENSFC,LSOIL), VMXFCS(LENSFC)
360  REAL :: VMNFCS(LENSFC), T2M(LENSFC)
361  REAL :: Q2M(LENSFC), SLPFCS(LENSFC)
362  REAL :: ABSFCS(LENSFC), OROG_UF(LENSFC)
363  REAL :: USTAR(LENSFC)
364  REAL :: FMM(LENSFC), FHH(LENSFC)
365  REAL :: RLA(LENSFC), RLO(LENSFC)
366  REAL(KIND=4) :: ZSOIL(LSOIL)
367  REAL :: SIG1T(LENSFC)
370  REAL, ALLOCATABLE :: STC_BCK(:,:), SMC_BCK(:,:), SLC_BCK(:,:)
371  REAL, ALLOCATABLE :: SLIFCS_FG(:)
372  INTEGER, ALLOCATABLE :: LANDINC_MASK_FG(:), LANDINC_MASK(:)
373  REAL, ALLOCATABLE :: SND_BCK(:), SND_INC(:), SWE_BCK(:)
374  REAL(KIND=KIND_IO8), ALLOCATABLE :: SLMASKL(:), SLMASKW(:)
375 
376  TYPE(NSST_DATA) :: NSST
377  real, dimension(idim,jdim) :: tf_clm,tf_trd,sal_clm
378  real, dimension(lensfc) :: tf_clm_tile,tf_trd_tile,sal_clm_tile
379  INTEGER :: veg_type_landice
380 
381  LOGICAL :: FILE_EXISTS, DO_SOI_INC, DO_SNO_INC
382 !--------------------------------------------------------------------------------
383 ! NST_FILE is the path/name of the gaussian GSI file which contains NSST
384 ! increments.
385 !--------------------------------------------------------------------------------
386 
387  NAMELIST/namsfcd/ nst_file, lnd_soi_file, do_sno_inc
388 
389  DATA nst_file/'NULL'/
390  DATA lnd_soi_file/'NULL'/
391 
392  do_sno_inc = .false.
393  do_soi_inc = .false.
394 
395 
396  sig1t = 0.0 ! Not a dead start!
397 
398  input_nml_file = "NULL"
399 
400  CALL baopenr(37, "fort.37", ierr)
401  READ (37, nml=namsfcd)
402 
403  print*
404  print*,'IN ROUTINE SFCDRV,IDIM=',idim,'JDIM=',jdim,'FH=',fh
405 
406 !--------------------------------------------------------------------------------
407 ! READ THE OROGRAPHY AND GRID POINT LAT/LONS FOR THE CUBED-SPHERE TILE.
408 !--------------------------------------------------------------------------------
409 
410  CALL read_lat_lon_orog(rla,rlo,orog,orog_uf,tile_num,idim,jdim,lensfc)
411 
412  DO i = 1, idim
413  idum(i,:) = i
414  ENDDO
415 
416  i_index = reshape(idum, (/lensfc/))
417 
418  DO i = 1, jdim
419  idum(:,i) = i
420  ENDDO
421 
422  j_index = reshape(idum, (/lensfc/) )
423 
424  IF (do_nsst) THEN
425  print*
426  print*,"WILL PROCESS NSST RECORDS."
427  ALLOCATE(nsst%C_0(lensfc))
428  ALLOCATE(nsst%C_D(lensfc))
429  ALLOCATE(nsst%D_CONV(lensfc))
430  ALLOCATE(nsst%DT_COOL(lensfc))
431  ALLOCATE(nsst%IFD(lensfc))
432  ALLOCATE(nsst%QRAIN(lensfc))
433  ALLOCATE(nsst%TREF(lensfc))
434  ALLOCATE(nsst%TFINC(lensfc))
435  ALLOCATE(nsst%W_0(lensfc))
436  ALLOCATE(nsst%W_D(lensfc))
437  ALLOCATE(nsst%XS(lensfc))
438  ALLOCATE(nsst%XT(lensfc))
439  ALLOCATE(nsst%XTTS(lensfc))
440  ALLOCATE(nsst%XU(lensfc))
441  ALLOCATE(nsst%XV(lensfc))
442  ALLOCATE(nsst%XZ(lensfc))
443  ALLOCATE(nsst%XZTS(lensfc))
444  ALLOCATE(nsst%Z_C(lensfc))
445  ALLOCATE(nsst%ZM(lensfc))
446  ALLOCATE(slifcs_fg(lensfc))
447  ENDIF
448 
449 IF (do_lndinc) THEN
450  ! identify variables to be updates, and allocate arrays.
451  IF (trim(lnd_soi_file) .NE. "NULL") THEN
452  do_soi_inc = .true.
453  print*
454  print*," APPLYING SOIL INCREMENTS FROM THE GSI"
455  ALLOCATE(stc_bck(lensfc, lsoil), smc_bck(lensfc, lsoil), slc_bck(lensfc,lsoil))
456  ALLOCATE(landinc_mask_fg(lensfc))
457  ENDIF
458  ! FOR NOW, CODE SO CAN DO BOTH, BUT MIGHT NEED TO THINK ABOUT THIS SOME MORE.
459  IF (do_sno_inc) THEN
460  ! ideally, would check here that sfcsub snow DA update is not also requested
461  ! but latter is controlled by fnsol, which is read in within that routine.
462  ! should be done at script level.
463  print*
464  print*," APPLYING SNOW INCREMENTS FROM JEDI"
465  ALLOCATE(snd_bck(lensfc), snd_inc(lensfc), swe_bck(lensfc))
466  ENDIF
467  ! set-up land mask info
468  ALLOCATE(landinc_mask(lensfc))
469  if (ivegsrc == 2) then ! sib
470  veg_type_landice=13
471  else
472  veg_type_landice=15
473  endif
474 ENDIF
475 
476 !--------------------------------------------------------------------------------
477 ! READ THE INPUT SURFACE DATA ON THE CUBED-SPHERE TILE.
478 !--------------------------------------------------------------------------------
479 
480  CALL read_data(lsoil,lensfc,do_nsst,.false.,is_noahmp=is_noahmp, &
481  tsffcs=tsffcs,smcfcs=smcfcs, &
482  swefcs=swefcs,stcfcs=stcfcs,tg3fcs=tg3fcs,zorfcs=zorfcs, &
483  cvfcs=cvfcs, cvbfcs=cvbfcs,cvtfcs=cvtfcs,albfcs=albfcs, &
484  vegfcs=vegfcs,slifcs=slifcs,cnpfcs=cnpfcs,f10m=f10m , &
485  vetfcs=vetfcs,sotfcs=sotfcs,alffcs=alffcs,ustar=ustar , &
486  fmm=fmm ,fhh=fhh ,sihfcs=sihfcs,sicfcs=sicfcs, &
487  sitfcs=sitfcs,tprcp=tprcp ,srflag=srflag,sndfcs=sndfcs, &
488  vmnfcs=vmnfcs,vmxfcs=vmxfcs,slcfcs=slcfcs,slpfcs=slpfcs, &
489  absfcs=absfcs,t2m=t2m ,q2m=q2m ,slmask=slmask, &
490  zsoil=zsoil, nsst=nsst)
491 
492  IF (use_ufo) THEN
493  print*
494  print*,'USE UNFILTERED OROGRAPHY.'
495  ELSE
496  orog_uf = 0.0
497  ENDIF
498 
499  DO i=1,lensfc
500  aisfcs(i) = 0.
501  IF(nint(slifcs(i)).EQ.2) aisfcs(i) = 1.
502  ENDDO
503 
504  IF (do_nsst) THEN
505  IF (.NOT. do_sfccycle ) THEN
506  print*
507  print*,"FIRST GUESS MASK ADJUSTED BY IFD RECORD"
508  slifcs_fg = slifcs
509  WHERE(nint(nsst%IFD) == 3) slifcs_fg = 2.0
510  ELSE
511  print*
512  print*,"SAVE FIRST GUESS MASK"
513  slifcs_fg = slifcs
514  ENDIF
515  ENDIF
516 
517  ! CALCULATE MASK FOR LAND INCREMENTS
518  IF (do_lndinc) &
519  CALL calculate_landinc_mask(slcfcs(:,1),swefcs, vetfcs, &
520  lensfc,veg_type_landice, landinc_mask)
521 
522 !--------------------------------------------------------------------------------
523 ! UPDATE SURFACE FIELDS.
524 !
525 ! FIRST, SET WATER AND LAND MASKS - SLMASKW/SLMASKL. FOR UNCOUPLED
526 ! (NON-FRACTIONAL) MODE, THESE ARE IDENTICAL TO THE CURRENT
527 ! MASK - '0' WATER; '1' LAND.
528 !--------------------------------------------------------------------------------
529 
530  IF (do_sfccycle) THEN
531  ALLOCATE(slmaskl(lensfc), slmaskw(lensfc))
532 ! for running uncoupled (non-fractional grid)
533  DO i=1,lensfc
534  IF(nint(slmask(i)) == 1) THEN
535  slmaskl(i) = 1.0_kind_io8
536  slmaskw(i) = 1.0_kind_io8
537  ELSE
538  slmaskl(i) = 0.0_kind_io8
539  slmaskw(i) = 0.0_kind_io8
540  ENDIF
541  if(nint(slmask(i)) == 0) then
542  min_ice(i) = 0.15_kind_io8
543  else
544  min_ice(i) = 0.0_kind_io8
545  endif
546  ENDDO
547  num_threads = num_parthds()
548  print*
549  print*,"CALL SFCCYCLE TO UPDATE SURFACE FIELDS."
550  CALL sfccycle(lugb,lensfc,lsoil,sig1t,deltsfc, &
551  iy,im,id,ih,fh,rla,rlo, &
552  slmaskl,slmaskw, orog, orog_uf, use_ufo, do_nsst, &
553  sihfcs,sicfcs,sitfcs,sndfcs,slcfcs, &
554  vmnfcs,vmxfcs,slpfcs,absfcs, &
555  tsffcs,swefcs,zorfcs,albfcs,tg3fcs, &
556  cnpfcs,smcfcs,stcfcs,slifcs,aisfcs, &
557  vegfcs,vetfcs,sotfcs,alffcs, &
558  cvfcs,cvbfcs,cvtfcs,myrank,num_threads, nlunit, &
559  sz_nml, input_nml_file, &
560  min_ice, &
561  ialb,isot,ivegsrc,tile_num,i_index,j_index)
562  DEALLOCATE(slmaskl, slmaskw)
563  ENDIF
564 
565 !--------------------------------------------------------------------------------
566 ! IF RUNNING WITH NSST, READ IN GSI FILE WITH THE UPDATED INCREMENTS (ON THE
567 ! GAUSSIAN GRID), INTERPOLATE INCREMENTS TO THE CUBED-SPHERE TILE, AND PERFORM
568 ! REQUIRED ADJUSTMENTS AND QC.
569 !--------------------------------------------------------------------------------
570 
571  IF (do_nsst) THEN
572  IF (nst_file == "NULL") THEN
573  print*
574  print*,"NO GSI FILE. ADJUST IFD FOR FORMER ICE POINTS."
575  DO i = 1, lensfc
576  IF (nint(slifcs_fg(i)) == 2 .AND. nint(slifcs(i)) == 0) THEN
577  nsst%IFD(i) = 3.0
578  ENDIF
579  ENDDO
580  nsst%TFINC = 0.0
581  ELSE
582  print*
583  print*,"ADJUST TREF FROM GSI INCREMENT"
584 !
585 ! Get tf climatology at the time
586 !
587  call get_tf_clm(rla,rlo,jdim,idim,iy,im,id,ih,tf_clm,tf_trd)
588  tf_clm_tile(:) = reshape(tf_clm, (/lensfc/) )
589  tf_trd_tile(:) = reshape(tf_trd, (/lensfc/) )
590 !
591 ! Get salinity climatology at the time
592 !
593  call get_sal_clm(rla,rlo,jdim,idim,iy,im,id,ih,sal_clm)
594  sal_clm_tile(:) = reshape(sal_clm, (/lensfc/) )
595 !
596 ! read tf analysis increment generated by GSI
597 !
598  CALL read_gsi_data(nst_file, 'NST')
599 !
600 ! update foundation & surface temperature for NSST
601 !
602  CALL adjust_nsst(rla,rlo,slifcs,slifcs_fg,tsffcs,sitfcs,sicfcs,stcfcs, &
603  nsst,lensfc,lsoil,idim,jdim,zsea1,zsea2,im,id,deltsfc, &
604  tf_clm_tile,tf_trd_tile,sal_clm_tile)
605  ENDIF
606  ENDIF
607 
608 !--------------------------------------------------------------------------------
609 ! READ IN AND APPLY LAND INCREMENTS FROM THE GSI
610 !--------------------------------------------------------------------------------
611 
612  IF (do_lndinc) THEN
613 
614  ! SNOW INCREMENTS
615  ! do snow first, as temperature updates will use snow analaysis
616  IF (do_sno_inc) THEN
617  ! updates are made to snow depth only over land (and not-land ice).
618  ! SWE is then updated from the snow depth analysis, using the model
619  ! forecast density
620 
621  !--------------------------------------------------------------------------------
622  ! read increments in
623  !--------------------------------------------------------------------------------
624 
625  ! Only coded for DA on native model grid (would always be the case for cycling DA)
626  CALL read_data(lsoil,lensfc,.false.,.true.,sndfcs=snd_inc)
627 
628  !--------------------------------------------------------------------------------
629  ! add increments to state vars
630  !--------------------------------------------------------------------------------
631 
632  ! store background states
633  snd_bck = sndfcs
634  swe_bck = swefcs
635 
636  CALL add_increment_snow(snd_inc,landinc_mask,lensfc,sndfcs)
637 
638  !--------------------------------------------------------------------------------
639  ! make any necessary adjustments to dependent variables
640  !--------------------------------------------------------------------------------
641 
642  CALL apply_land_da_adjustments_snd(lsm, lensfc, landinc_mask, swe_bck, snd_bck, &
643  sndfcs, swefcs)
644 
645  ENDIF
646 
647  ! SOIL INCREMENTS
648  IF (do_soi_inc) THEN
649 
650  !--------------------------------------------------------------------------------
651  ! re-calculate soilsnow mask if snow has been updated.
652  !--------------------------------------------------------------------------------
653 
654  landinc_mask_fg = landinc_mask
655 
656  IF (do_sfccycle .OR. do_sno_inc) THEN
657  CALL calculate_landinc_mask(slcfcs(:,1),swefcs, vetfcs, lensfc, &
658  veg_type_landice, landinc_mask )
659  ENDIF
660 
661  !--------------------------------------------------------------------------------
662  ! read increments in
663  !--------------------------------------------------------------------------------
664 
665  INQUIRE(file=trim(lnd_soi_file), exist=file_exists)
666  IF (.not. file_exists) then
667  print *, 'FATAL ERROR: land increment update requested, but file does not exist: ', &
668  trim(lnd_soi_file)
669  call mpi_abort(mpi_comm_world, 10, ierr)
670  ENDIF
671 
672  CALL read_gsi_data(lnd_soi_file, 'LND', lsoil=lsoil)
673 
674  !--------------------------------------------------------------------------------
675  ! add increments to state vars
676  !--------------------------------------------------------------------------------
677  ! when applying increments, will often need to adjust other land states in response
678  ! to the changes made. Need to store bacground, apply the increments, then make
679  ! secondart adjustments. When updating more than one state, be careful about the
680  ! order if increments and secondary adjustments.
681 
682  ! store background states
683  stc_bck = stcfcs
684  smc_bck = smcfcs ! not used yet.
685  slc_bck = slcfcs ! not used yet.
686 
687  CALL add_increment_soil(rla,rlo,stcfcs,landinc_mask,landinc_mask_fg, &
688  lensfc,lsoil,idim,jdim, myrank)
689 
690  !--------------------------------------------------------------------------------
691  ! make any necessary adjustments to dependent variables
692  !--------------------------------------------------------------------------------
693 
694  CALL apply_land_da_adjustments_stc(lsm, isot, ivegsrc,lensfc, lsoil, &
695  sotfcs, landinc_mask_fg, stc_bck, stcfcs, smcfcs, slcfcs)
696 
697  ENDIF ! soil increments
698 
699 !--------------------------------------------------------------------------------
700 ! clean up
701 !--------------------------------------------------------------------------------
702 
703  ! to do - save and write out STC_INC? (soil temperature increments)
704  IF(ALLOCATED(landinc_mask_fg)) DEALLOCATE(landinc_mask_fg)
705  IF(ALLOCATED(landinc_mask)) DEALLOCATE(landinc_mask)
706  IF(ALLOCATED(stc_bck)) DEALLOCATE(stc_bck)
707  IF(ALLOCATED(smc_bck)) DEALLOCATE(smc_bck)
708  IF(ALLOCATED(slc_bck)) DEALLOCATE(slc_bck)
709  IF(ALLOCATED(snd_bck)) DEALLOCATE(snd_bck)
710  IF(ALLOCATED(swe_bck)) DEALLOCATE(swe_bck)
711  IF(ALLOCATED(snd_inc)) DEALLOCATE(snd_inc)
712 
713  ENDIF
714 !--------------------------------------------------------------------------------
715 ! WRITE OUT UPDATED SURFACE DATA ON THE CUBED-SPHERE TILE.
716 !--------------------------------------------------------------------------------
717 
718  IF (is_noahmp) THEN
719 
720  CALL write_data(lensfc,idim,jdim,lsoil,do_nsst,nsst,vegfcs=vegfcs)
721 
722  ELSE
723 
724  CALL write_data(lensfc,idim,jdim,lsoil, &
725  do_nsst,nsst,slifcs=slifcs,tsffcs=tsffcs,vegfcs=vegfcs, &
726  swefcs=swefcs,tg3fcs=tg3fcs,zorfcs=zorfcs, &
727  albfcs=albfcs,alffcs=alffcs,cnpfcs=cnpfcs, &
728  f10m=f10m,t2m=t2m,q2m=q2m,vetfcs=vetfcs, &
729  sotfcs=sotfcs,ustar=ustar,fmm=fmm,fhh=fhh, &
730  sicfcs=sicfcs,sihfcs=sihfcs,sitfcs=sitfcs,tprcp=tprcp, &
731  srflag=srflag,swdfcs=sndfcs,vmnfcs=vmnfcs, &
732  vmxfcs=vmxfcs,slpfcs=slpfcs,absfcs=absfcs, &
733  slcfcs=slcfcs,smcfcs=smcfcs,stcfcs=stcfcs)
734 
735  ENDIF
736 
737  IF (do_nsst) THEN
738  DEALLOCATE(nsst%C_0)
739  DEALLOCATE(nsst%C_D)
740  DEALLOCATE(nsst%D_CONV)
741  DEALLOCATE(nsst%DT_COOL)
742  DEALLOCATE(nsst%IFD)
743  DEALLOCATE(nsst%QRAIN)
744  DEALLOCATE(nsst%TREF)
745  DEALLOCATE(nsst%TFINC)
746  DEALLOCATE(nsst%W_0)
747  DEALLOCATE(nsst%W_D)
748  DEALLOCATE(nsst%XS)
749  DEALLOCATE(nsst%XT)
750  DEALLOCATE(nsst%XTTS)
751  DEALLOCATE(nsst%XU)
752  DEALLOCATE(nsst%XV)
753  DEALLOCATE(nsst%XZ)
754  DEALLOCATE(nsst%XZTS)
755  DEALLOCATE(nsst%Z_C)
756  DEALLOCATE(nsst%ZM)
757  DEALLOCATE(slifcs_fg)
758  ENDIF
759 
760  RETURN
761 
762  END SUBROUTINE sfcdrv
763 
795  SUBROUTINE adjust_nsst(RLA,RLO,SLMSK_TILE,SLMSK_FG_TILE,SKINT_TILE,&
796  SICET_TILE,sice_tile,SOILT_TILE,NSST,LENSFC,LSOIL, &
797  IDIM,JDIM,ZSEA1,ZSEA2,MON,DAY,DELTSFC, &
798  tf_clm_tile,tf_trd_tile,sal_clm_tile)
800  USE utils
801  USE gdswzd_mod
802  USE read_write_data, ONLY : idim_gaus, jdim_gaus, &
803  slmsk_gaus, dtref_gaus, &
804  nsst_data
805 
806  USE mpi
807 
808  IMPLICIT NONE
809 
810  INTEGER, INTENT(IN) :: LENSFC, LSOIL, IDIM, JDIM, MON, DAY
811 
812  REAL, INTENT(IN) :: SLMSK_TILE(LENSFC), SLMSK_FG_TILE(LENSFC)
813  real, intent(in) :: tf_clm_tile(lensfc),tf_trd_tile(lensfc),sal_clm_tile(lensfc)
814  REAL, INTENT(IN) :: ZSEA1, ZSEA2, DELTSFC
815  REAL, INTENT(INOUT) :: RLA(LENSFC), RLO(LENSFC), SKINT_TILE(LENSFC)
816  REAL, INTENT(INOUT) :: SICET_TILE(LENSFC),sice_tile(lensfc),SOILT_TILE(LENSFC,LSOIL)
817 
818  TYPE(NSST_DATA) :: NSST
819 
820  REAL, PARAMETER :: TMAX=313.0,tzero=273.16
821 
822  INTEGER :: IOPT, NRET, KGDS_GAUS(200)
823  INTEGER :: IGAUS, JGAUS, IJ, II, JJ, III, JJJ, KRAD
824  INTEGER :: ISTART, IEND, JSTART, JEND
825  INTEGER :: MASK_TILE, MASK_FG_TILE
826  INTEGER :: ITILE, JTILE
827  INTEGER :: MAX_SEARCH, J, IERR
828  INTEGER :: IGAUSP1, JGAUSP1
829  integer :: nintp,nsearched,nice,nland
830  integer :: nfill,nfill_tice,nfill_clm
831  integer :: nset_thaw,nset_thaw_s,nset_thaw_i,nset_thaw_c
832 
833  INTEGER, ALLOCATABLE :: ID1(:,:), ID2(:,:), JDC(:,:)
834 
835  LOGICAL :: IS_ICE
836 
837  real :: tfreez
838  REAL :: TREF_SAVE,WSUM,tf_ice,tf_thaw
839  REAL :: FILL, DTZM, GAUS_RES_KM, DTREF
840  REAL, ALLOCATABLE :: XPTS(:), YPTS(:), LATS(:), LONS(:)
841  REAL, ALLOCATABLE :: DUM2D(:,:), LATS_RAD(:), LONS_RAD(:)
842  REAL, ALLOCATABLE :: AGRID(:,:,:), S2C(:,:,:)
843 
844  kgds_gaus = 0
845  kgds_gaus(1) = 4 ! OCT 6 - TYPE OF GRID (GAUSSIAN)
846  kgds_gaus(2) = idim_gaus ! OCT 7-8 - # PTS ON LATITUDE CIRCLE
847  kgds_gaus(3) = jdim_gaus
848  kgds_gaus(4) = 90000 ! OCT 11-13 - LAT OF ORIGIN
849  kgds_gaus(5) = 0 ! OCT 14-16 - LON OF ORIGIN
850  kgds_gaus(6) = 128 ! OCT 17 - RESOLUTION FLAG
851  kgds_gaus(7) = -90000 ! OCT 18-20 - LAT OF EXTREME POINT
852  kgds_gaus(8) = nint(-360000./float(idim_gaus)) ! OCT 21-23 - LON OF EXTREME POINT
853  kgds_gaus(9) = nint((360.0 / float(idim_gaus))*1000.0)
854  ! OCT 24-25 - LONGITUDE DIRECTION INCR.
855  kgds_gaus(10) = jdim_gaus/2 ! OCT 26-27 - NUMBER OF CIRCLES POLE TO EQUATOR
856  kgds_gaus(12) = 255 ! OCT 29 - RESERVED
857  kgds_gaus(20) = 255 ! OCT 5 - NOT USED, SET TO 255
858 
859  print*
860  print*,'ADJUST NSST USING GSI INCREMENTS ON GAUSSIAN GRID'
861 
862 !----------------------------------------------------------------------
863 ! CALL GDSWZD TO COMPUTE THE LAT/LON OF EACH GSI GAUSSIAN GRID POINT.
864 !----------------------------------------------------------------------
865 
866  iopt = 0
867  fill = -9999.
868  ALLOCATE(xpts(idim_gaus*jdim_gaus))
869  ALLOCATE(ypts(idim_gaus*jdim_gaus))
870  ALLOCATE(lats(idim_gaus*jdim_gaus))
871  ALLOCATE(lons(idim_gaus*jdim_gaus))
872  xpts = fill
873  ypts = fill
874  lats = fill
875  lons = fill
876 
877  CALL gdswzd(kgds_gaus,iopt,(idim_gaus*jdim_gaus),fill,xpts,ypts,lons,lats,nret)
878 
879  IF (nret /= (idim_gaus*jdim_gaus)) THEN
880  print*,'FATAL ERROR: PROBLEM IN GDSWZD. STOP.'
881  CALL mpi_abort(mpi_comm_world, 12, ierr)
882  ENDIF
883 
884  DEALLOCATE (xpts, ypts)
885 
886  ALLOCATE(dum2d(idim_gaus,jdim_gaus))
887  dum2d = reshape(lats, (/idim_gaus,jdim_gaus/) )
888  DEALLOCATE(lats)
889 
890  ALLOCATE(lats_rad(jdim_gaus))
891  DO j = 1, jdim_gaus
892  lats_rad(j) = dum2d(1,jdim_gaus-j+1) * 3.1415926 / 180.0
893  ENDDO
894 
895  dum2d = reshape(lons, (/idim_gaus,jdim_gaus/) )
896  DEALLOCATE(lons)
897  ALLOCATE(lons_rad(idim_gaus))
898  lons_rad = dum2d(:,1) * 3.1415926 / 180.0
899 
900  DEALLOCATE(dum2d)
901 
902  ALLOCATE(agrid(idim,jdim,2))
903  agrid(:,:,1) = reshape(rlo, (/idim,jdim/) )
904  agrid(:,:,2) = reshape(rla, (/idim,jdim/) )
905  agrid = agrid * 3.1415926 / 180.0
906 
907  ALLOCATE(id1(idim,jdim))
908  ALLOCATE(id2(idim,jdim))
909  ALLOCATE(jdc(idim,jdim))
910  ALLOCATE(s2c(idim,jdim,4))
911 
912 !----------------------------------------------------------------------
913 ! COMPUTE BILINEAR WEIGHTS FOR EACH MODEL POINT FROM THE NEAREST
914 ! FOUR GSI/GAUSSIAN POINTS. DOES NOT ACCOUNT FOR MASK. THAT
915 ! HAPPENS LATER.
916 !----------------------------------------------------------------------
917 
918  CALL remap_coef( 1, idim, 1, jdim, idim_gaus, jdim_gaus, &
919  lons_rad, lats_rad, id1, id2, jdc, s2c, agrid )
920 
921  DEALLOCATE(lons_rad, lats_rad, agrid)
922 
923 !----------------------------------------------------------------------
924 ! THE MAXIMUM DISTANCE TO SEARCH IS 500 KM. HOW MANY GAUSSIAN
925 ! GRID LENGTHS IS THAT?
926 !----------------------------------------------------------------------
927 
928  gaus_res_km = 360.0 / idim_gaus * 111.0
929  max_search = ceiling(500.0/gaus_res_km)
930 
931  print*
932  print*,'MAXIMUM SEARCH IS ',max_search, ' GAUSSIAN POINTS.'
933  print*
934 
935 !
936 ! Initialize variables for counts statitics to be zeros
937 !
938  nintp = 0
939  nset_thaw = 0
940  nset_thaw_s = 0
941  nset_thaw_i = 0
942  nset_thaw_c = 0
943  nsearched = 0
944  nfill = 0
945  nfill_tice = 0
946  nfill_clm = 0
947  nice = 0
948  nland = 0
949 !----------------------------------------------------------------------
950 ! TREF INCREMENT WILL BE OUTPUT. INITIALIZE TO ZERO.
951 !----------------------------------------------------------------------
952 
953  nsst%TFINC = 0.0
954 
955  ij_loop : DO ij = 1, lensfc
956 
957  mask_tile = nint(slmsk_tile(ij))
958  mask_fg_tile = nint(slmsk_fg_tile(ij))
959 
960 !
961 ! when sea ice exists, get salinity dependent water temperature
962 !
963  tf_ice = tfreez(sal_clm_tile(ij)) + tzero
964 !----------------------------------------------------------------------
965 ! SKIP LAND POINTS. NSST NOT APPLIED AT LAND.
966 !----------------------------------------------------------------------
967 
968  IF (mask_tile == 1) THEN
969  nland = nland + 1
970  cycle ij_loop
971  ENDIF
972 
973 !
974 ! these are ice points. set tref to tf_ice and update tmpsfc.
975 !
976  if (mask_tile == 2) then
977  nsst%tref(ij)=tf_ice ! water part tmp set
978  skint_tile(ij)=(1.0-sice_tile(ij))*nsst%tref(ij)+sice_tile(ij)*sicet_tile(ij)
979  nice = nice + 1
980  cycle ij_loop
981  endif
982 
983 !
984 ! Get i,j index on array of (idim,jdim) from known ij
985 !
986  jtile = (ij-1) / idim + 1
987  itile = mod(ij,idim)
988  IF (itile==0) itile = idim
989 
990 !----------------------------------------------------------------------
991 ! IF THE MODEL POINT WAS ICE COVERED, BUT IS NOW OPEN WATER, SET
992 ! TREF TO searched adjascent open water onea, if failed the search, set to
993 ! weighted average of tf_ice and tf_clm. For NSST vars, set xz TO '30' AND ALL OTHER FIELDS TO ZERO.
994 !----------------------------------------------------------------------
995 
996  IF (mask_fg_tile == 2 .AND. mask_tile == 0) THEN
997 !
998 ! set background for the thaw (just melted water) situation
999 !
1000  call tf_thaw_set(nsst%tref,nint(slmsk_fg_tile),itile,jtile,tf_ice,tf_clm_tile(ij),tf_thaw,idim,jdim, &
1001  nset_thaw_s,nset_thaw_i,nset_thaw_c)
1002  call nsst_water_reset(nsst,ij,tf_thaw)
1003  nset_thaw = nset_thaw + 1
1004  ENDIF
1005 
1006 !----------------------------------------------------------------------
1007 ! THESE ARE POINTS THAT ARE OPEN WATER AND WERE OPEN WATER PRIOR
1008 ! TO ANY ICE UPDATE BY SFCCYCLE. UPDATE TREF AND SKIN TEMP.
1009 ! AT OPEN WATER POINTS, THE SEA ICE TEMPERATURE (SICET_TILE) AND
1010 ! SOIL COLUMN TEMPERATURE (SOILT_TILE) ARE SET TO THE SKIN TEMP.
1011 ! IT IS SIMPLY A FILLER VALUE. THESE FIELDS ARE NOT USED AT
1012 ! OPEN WATER POINTS.
1013 !----------------------------------------------------------------------
1014 !----------------------------------------------------------------------
1015 ! SEE IF ANY OF THE NEAREST GSI POINTS MASK AREA OPEN WATER.
1016 ! IF SO, APPLY NSST INCREMENT USING BILINEAR INTERPOLATION.
1017 !----------------------------------------------------------------------
1018 
1019  igaus = id1(itile,jtile)
1020  jgaus = jdc(itile,jtile)
1021  igausp1 = id2(itile,jtile)
1022  jgausp1 = jdc(itile,jtile)+1
1023 
1024  IF (slmsk_gaus(igaus,jgaus) == 0 .OR. &
1025  slmsk_gaus(igausp1,jgaus) == 0 .OR. &
1026  slmsk_gaus(igausp1,jgausp1) == 0 .OR. &
1027  slmsk_gaus(igaus,jgausp1) == 0) THEN
1028 
1029  dtref = 0.0
1030  wsum = 0.0
1031 
1032  IF (slmsk_gaus(igaus,jgaus) == 0) THEN
1033  dtref = dtref + (s2c(itile,jtile,1) * dtref_gaus(igaus,jgaus))
1034  wsum = wsum + s2c(itile,jtile,1)
1035  ENDIF
1036 
1037  IF (slmsk_gaus(igausp1,jgaus) == 0) THEN
1038  dtref = dtref + (s2c(itile,jtile,2) * dtref_gaus(igausp1,jgaus))
1039  wsum = wsum + s2c(itile,jtile,2)
1040  ENDIF
1041 
1042  IF (slmsk_gaus(igausp1,jgausp1) == 0) THEN
1043  dtref = dtref + (s2c(itile,jtile,3) * dtref_gaus(igausp1,jgausp1))
1044  wsum = wsum + s2c(itile,jtile,3)
1045  ENDIF
1046 
1047  IF (slmsk_gaus(igaus,jgausp1) == 0) THEN
1048  dtref = dtref + (s2c(itile,jtile,4) * dtref_gaus(igaus,jgausp1))
1049  wsum = wsum + s2c(itile,jtile,4)
1050  ENDIF
1051 
1052  nintp = nintp + 1
1053  dtref = dtref / wsum
1054 
1055  tref_save = nsst%TREF(ij)
1056  nsst%TREF(ij) = nsst%TREF(ij) + dtref
1057  nsst%TREF(ij) = max(nsst%TREF(ij), tf_ice)
1058  nsst%TREF(ij) = min(nsst%TREF(ij), tmax)
1059  nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
1060 
1061  CALL dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij), &
1062  nsst%Z_C(ij),zsea1,zsea2,dtzm)
1063 
1064  skint_tile(ij) = nsst%TREF(ij) + dtzm
1065  skint_tile(ij) = max(skint_tile(ij), tf_ice)
1066  skint_tile(ij) = min(skint_tile(ij), tmax)
1067 
1068  sicet_tile(ij) = skint_tile(ij)
1069  soilt_tile(ij,:) = skint_tile(ij)
1070 
1071 !----------------------------------------------------------------------
1072 ! NO NEARBY GSI/GAUSSIAN OPEN WATER POINTS. PERFORM A SPIRAL SEARCH TO
1073 ! FIND NEAREST NON-LAND POINT ON GSI/GAUSSIAN GRID.
1074 !----------------------------------------------------------------------
1075 
1076  ELSE
1077 
1078  is_ice = .false.
1079 
1080  DO krad = 1, max_search
1081 
1082  istart = igaus - krad
1083  iend = igaus + krad
1084  jstart = jgaus - krad
1085  jend = jgaus + krad
1086 
1087  DO jj = jstart, jend
1088  DO ii = istart, iend
1089 
1090  IF((jj == jstart) .OR. (jj == jend) .OR. &
1091  (ii == istart) .OR. (ii == iend)) THEN
1092 
1093  IF ((jj >= 1) .AND. (jj <= jdim_gaus)) THEN
1094 
1095  jjj = jj
1096  IF (ii <= 0) THEN
1097  iii = idim_gaus + ii
1098  ELSE IF (ii >= (idim_gaus+1)) THEN
1099  iii = ii - idim_gaus
1100  ELSE
1101  iii = ii
1102  END IF
1103 
1104 !----------------------------------------------------------------------
1105 ! SEE IF NEARBY POINTS ARE SEA ICE. IF THEY ARE, AND THE SEARCH FOR
1106 ! A GAUSSIAN GRID OPEN WATER POINT FAILS, THEN TREF WILL BE SET TO
1107 ! FREEZING BELOW.
1108 !----------------------------------------------------------------------
1109 
1110  IF (krad <= 2 .AND. slmsk_gaus(iii,jjj) == 2) is_ice = .true.
1111 
1112  IF (slmsk_gaus(iii,jjj) == 0) THEN
1113 
1114 ! PRINT*,'MISMATCH AT TILE POINT ',ITILE,JTILE
1115 ! PRINT*,'UPDATE TREF USING GSI INCREMENT AT ',III,JJJ,DTREF_GAUS(III,JJJ)
1116  nsearched = nsearched + 1
1117 
1118  tref_save = nsst%TREF(ij)
1119  nsst%TREF(ij ) = nsst%TREF(ij) + dtref_gaus(iii,jjj)
1120  nsst%TREF(ij) = max(nsst%TREF(ij), tf_ice)
1121  nsst%TREF(ij) = min(nsst%TREF(ij), tmax)
1122  nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
1123 
1124  CALL dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij), &
1125  nsst%Z_C(ij),zsea1,zsea2,dtzm)
1126 
1127  skint_tile(ij) = nsst%TREF(ij) + dtzm
1128  skint_tile(ij) = max(skint_tile(ij), tf_ice)
1129  skint_tile(ij) = min(skint_tile(ij), tmax)
1130 
1131  sicet_tile(ij) = skint_tile(ij)
1132  soilt_tile(ij,:) = skint_tile(ij)
1133  cycle ij_loop
1134 
1135  ENDIF ! GSI/Gaussian mask is open water
1136 
1137  ENDIF
1138 
1139  ENDIF
1140 
1141  ENDDO
1142  ENDDO
1143 
1144  ENDDO ! KRAD LOOP
1145 
1146 !----------------------------------------------------------------------
1147 ! THE SEARCH FAILED. IF THERE IS NEARBY ICE, SET TREF TO FREEZING.
1148 ! ELSE UPDATE TREF BASED ON THE ANNUAL SST CYCLE.
1149 !----------------------------------------------------------------------
1150 
1151 ! PRINT*,'WARNING !!!!!! SEARCH FAILED AT TILE POINT ',ITILE,JTILE
1152 
1153  nfill = nfill + 1
1154  IF (is_ice) THEN
1155  nsst%TREF(ij) = tf_ice
1156 ! PRINT*,"NEARBY ICE. SET TREF TO FREEZING"
1157  nfill_tice = nfill_tice + 1
1158  ELSE
1159  tref_save = nsst%TREF(ij)
1160  nsst%TREF(ij) = nsst%TREF(ij) + tf_trd_tile(ij)
1161  nsst%TREF(ij) = max(nsst%TREF(ij), tf_ice)
1162  nsst%TREF(ij) = min(nsst%TREF(ij), tmax)
1163  nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
1164 ! PRINT*,'UPDATE TREF FROM SST CLIMO ',DTREF
1165  nfill_clm = nfill_clm + 1
1166  ENDIF
1167 
1168  CALL dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij), &
1169  nsst%Z_C(ij),zsea1,zsea2,dtzm)
1170 
1171  skint_tile(ij) = nsst%TREF(ij) + dtzm
1172  skint_tile(ij) = max(skint_tile(ij), tf_ice)
1173  skint_tile(ij) = min(skint_tile(ij), tmax)
1174 
1175  sicet_tile(ij) = skint_tile(ij)
1176  soilt_tile(ij,:) = skint_tile(ij)
1177 
1178  ENDIF ! NEARBY GAUSSIAN POINTS ARE OPEN WATER?
1179 
1180  ENDDO ij_loop
1181 
1182  write(*,'(a)') 'statistics of grids number processed for tile : '
1183  write(*,'(a,I8)') ' nintp = ',nintp
1184  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
1185  write(*,'(a,I8)') ' nsearched = ',nsearched
1186  write(*,'(a,3I6)') ' nfill,nfill_tice,nfill_clm = ',nfill,nfill_tice,nfill_clm
1187  write(*,'(a,I8)') ' nice = ',nice
1188  write(*,'(a,I8)') ' nland = ',nland
1189 
1190  DEALLOCATE(id1, id2, jdc, s2c)
1191 
1192  END SUBROUTINE adjust_nsst
1193 
1204  SUBROUTINE climo_trend(LATITUDE, MON, DAY, DELTSFC, DTREF)
1205  IMPLICIT NONE
1206 
1207  INTEGER, INTENT(IN) :: MON, DAY
1208 
1209  REAL, INTENT(IN) :: LATITUDE, DELTSFC
1210  REAL, INTENT(OUT) :: DTREF
1211 
1212  INTEGER :: NUM_DAYS(12), MON2, MON1
1213 
1214  REAL, TARGET :: SST_80_90(12)
1215  REAL, TARGET :: SST_70_80(12)
1216  REAL, TARGET :: SST_60_70(12)
1217  REAL, TARGET :: SST_50_60(12)
1218  REAL, TARGET :: SST_40_50(12)
1219  REAL, TARGET :: SST_30_40(12)
1220  REAL, TARGET :: SST_20_30(12)
1221  REAL, TARGET :: SST_10_20(12)
1222  REAL, TARGET :: SST_00_10(12)
1223  REAL, TARGET :: SST_M10_00(12)
1224  REAL, TARGET :: SST_M20_M10(12)
1225  REAL, TARGET :: SST_M30_M20(12)
1226  REAL, TARGET :: SST_M40_M30(12)
1227  REAL, TARGET :: SST_M50_M40(12)
1228  REAL, TARGET :: SST_M60_M50(12)
1229  REAL, TARGET :: SST_M70_M60(12)
1230  REAL, TARGET :: SST_M80_M70(12)
1231  REAL, TARGET :: SST_M90_M80(12)
1232 
1233  REAL, POINTER :: SST(:)
1234 
1235  DATA num_days /31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31/
1236 
1237  DATA sst_80_90 /271.466, 271.458, 271.448, 271.445, 271.519, 271.636, &
1238  272.023, 272.066, 272.001, 271.698, 271.510, 271.472/
1239 
1240  DATA sst_70_80 /272.149, 272.103, 272.095, 272.126, 272.360, 272.988, &
1241  274.061, 274.868, 274.415, 273.201, 272.468, 272.268/
1242 
1243  DATA sst_60_70 /274.240, 274.019, 273.988, 274.185, 275.104, 276.875, &
1244  279.005, 280.172, 279.396, 277.586, 275.818, 274.803/
1245 
1246  DATA sst_50_60 /277.277, 276.935, 277.021, 277.531, 279.100, 281.357, &
1247  283.735, 285.171, 284.399, 282.328, 279.918, 278.199/
1248 
1249  DATA sst_40_50 /281.321, 280.721, 280.850, 281.820, 283.958, 286.588, &
1250  289.195, 290.873, 290.014, 287.652, 284.898, 282.735/
1251 
1252  DATA sst_30_40 /289.189, 288.519, 288.687, 289.648, 291.547, 293.904, &
1253  296.110, 297.319, 296.816, 295.225, 292.908, 290.743/
1254 
1255  DATA sst_20_30 /294.807, 294.348, 294.710, 295.714, 297.224, 298.703, &
1256  299.682, 300.127, 300.099, 299.455, 297.953, 296.177/
1257 
1258  DATA sst_10_20 /298.878, 298.720, 299.033, 299.707, 300.431, 300.709, &
1259  300.814, 300.976, 301.174, 301.145, 300.587, 299.694/
1260 
1261  DATA sst_00_10 /300.415, 300.548, 300.939, 301.365, 301.505, 301.141, &
1262  300.779, 300.660, 300.818, 300.994, 300.941, 300.675/
1263 
1264  DATA sst_m10_00 /300.226, 300.558, 300.914, 301.047, 300.645, 299.870, &
1265  299.114, 298.751, 298.875, 299.294, 299.721, 299.989/
1266 
1267  DATA sst_m20_m10 /299.547, 299.985, 300.056, 299.676, 298.841, 297.788, &
1268  296.893, 296.491, 296.687, 297.355, 298.220, 298.964/
1269 
1270  DATA sst_m30_m20 /297.524, 298.073, 297.897, 297.088, 295.846, 294.520, &
1271  293.525, 293.087, 293.217, 293.951, 295.047, 296.363/
1272 
1273  DATA sst_m40_m30 /293.054, 293.765, 293.468, 292.447, 291.128, 289.781, &
1274  288.773, 288.239, 288.203, 288.794, 289.947, 291.553/
1275 
1276  DATA sst_m50_m40 /285.052, 285.599, 285.426, 284.681, 283.761, 282.826, &
1277  282.138, 281.730, 281.659, 281.965, 282.768, 283.961/
1278 
1279  DATA sst_m60_m50 /277.818, 278.174, 277.991, 277.455, 276.824, 276.229, &
1280  275.817, 275.585, 275.560, 275.687, 276.142, 276.968/
1281 
1282  DATA sst_m70_m60 /273.436, 273.793, 273.451, 272.813, 272.349, 272.048, &
1283  271.901, 271.838, 271.845, 271.889, 272.080, 272.607/
1284 
1285  DATA sst_m80_m70 /271.579, 271.578, 271.471, 271.407, 271.392, 271.391, &
1286  271.390, 271.391, 271.394, 271.401, 271.422, 271.486/
1287 
1288  DATA sst_m90_m80 /271.350, 271.350, 271.350, 271.350, 271.350, 271.350, &
1289  271.350, 271.350, 271.350, 271.350, 271.350, 271.350/
1290 
1291  NULLIFY(sst)
1292  IF (latitude > 80.0) THEN
1293  sst => sst_80_90
1294  ELSEIF (latitude > 70.0) THEN
1295  sst => sst_70_80
1296  ELSEIF (latitude > 60.0) THEN
1297  sst => sst_60_70
1298  ELSEIF (latitude > 50.0) THEN
1299  sst => sst_50_60
1300  ELSEIF (latitude > 40.0) THEN
1301  sst => sst_40_50
1302  ELSEIF (latitude > 30.0) THEN
1303  sst => sst_30_40
1304  ELSEIF (latitude > 20.0) THEN
1305  sst => sst_20_30
1306  ELSEIF (latitude > 10.0) THEN
1307  sst => sst_10_20
1308  ELSEIF (latitude > 0.0) THEN
1309  sst => sst_00_10
1310  ELSEIF (latitude > -10.0) THEN
1311  sst => sst_m10_00
1312  ELSEIF (latitude > -20.0) THEN
1313  sst => sst_m20_m10
1314  ELSEIF (latitude > -30.0) THEN
1315  sst => sst_m30_m20
1316  ELSEIF (latitude > -40.0) THEN
1317  sst => sst_m40_m30
1318  ELSEIF (latitude > -50.0) THEN
1319  sst => sst_m50_m40
1320  ELSEIF (latitude > -60.0) THEN
1321  sst => sst_m60_m50
1322  ELSEIF (latitude > -70.0) THEN
1323  sst => sst_m70_m60
1324  ELSEIF (latitude > -80.0) THEN
1325  sst => sst_m80_m70
1326  ELSE
1327  sst => sst_m90_m80
1328  END IF
1329 
1330  IF (day >= 15) THEN
1331  mon2 = mon+1
1332  IF(mon2 == 13) mon2 = 1
1333  mon1 = mon
1334  dtref = (sst(mon2) - sst(mon1)) / num_days(mon1)
1335  ELSE
1336  mon1 = mon - 1
1337  IF (mon1 == 0) mon1=12
1338  mon2 = mon
1339  dtref = (sst(mon2) - sst(mon1)) / num_days(mon1)
1340  ENDIF
1341 
1342  dtref = dtref * (deltsfc / 24.0)
1343 
1344  END SUBROUTINE climo_trend
1345 
1357  SUBROUTINE dtzm_point(XT,XZ,DT_COOL,ZC,Z1,Z2,DTZM)
1358  implicit none
1359 
1360  real, intent(in) :: xt,xz,dt_cool,zc,z1,z2
1361  real, intent(out) :: dtzm
1362 
1363  real, parameter :: zero = 0.0
1364  real, parameter :: one = 1.0
1365  real, parameter :: half = 0.5
1366  real :: dt_warm,dtw,dtc
1367 !
1368 ! get the mean warming in the range of z=z1 to z=z2
1369 !
1370  dtw = zero
1371  if ( xt > zero ) then
1372  dt_warm = (xt+xt)/xz ! Tw(0)
1373  if ( z1 < z2) then
1374  if ( z2 < xz ) then
1375  dtw = dt_warm*(one-(z1+z2)/(xz+xz))
1376  elseif ( z1 < xz .and. z2 >= xz ) then
1377  dtw = half*(one-z1/xz)*dt_warm*(xz-z1)/(z2-z1)
1378  endif
1379  elseif ( z1 == z2 ) then
1380  if ( z1 < xz ) then
1381  dtw = dt_warm*(one-z1/xz)
1382  endif
1383  endif
1384  endif
1385 !
1386 ! get the mean cooling in the range of z=z1 to z=z2
1387 !
1388  dtc = zero
1389  if ( zc > zero ) then
1390  if ( z1 < z2) then
1391  if ( z2 < zc ) then
1392  dtc = dt_cool*(one-(z1+z2)/(zc+zc))
1393  elseif ( z1 < zc .and. z2 >= zc ) then
1394  dtc = half*(one-z1/zc)*dt_cool*(zc-z1)/(z2-z1)
1395  endif
1396  elseif ( z1 == z2 ) then
1397  if ( z1 < zc ) then
1398  dtc = dt_cool*(one-z1/zc)
1399  endif
1400  endif
1401  endif
1402 
1403 !
1404 ! get the mean T departure from Tf in the range of z=z1 to z=z2
1405 !
1406  dtzm = dtw - dtc
1407 
1408  END SUBROUTINE dtzm_point
1409 
1429  subroutine tf_thaw_set(tf_ij,mask_ij,itile,jtile,tice,tclm,tf_thaw,nx,ny, &
1430  nset_thaw_s,nset_thaw_i,nset_thaw_c)
1432  real, dimension(nx*ny), intent(in) :: tf_ij
1433  integer, dimension(nx*ny), intent(in) :: mask_ij
1434  real, intent(in) :: tice,tclm
1435  integer, intent(in) :: itile,jtile,nx,ny
1436  real, intent(out) :: tf_thaw
1437  integer, intent(inout) :: nset_thaw_s,nset_thaw_i,nset_thaw_c
1438 ! Local
1439  real, parameter :: bmiss = -999.0
1440  real, dimension(nx,ny) :: tf
1441  integer, dimension(nx,ny) :: mask
1442  integer :: krad,max_search,istart,iend,jstart,jend
1443  integer :: ii,jj,iii,jjj
1444  logical :: is_ice
1445 
1446  max_search = 2
1447 
1448  mask(:,:) = reshape(mask_ij,(/nx,ny/) )
1449  tf(:,:) = reshape(tf_ij,(/nx,ny/) )
1450 
1451  tf_thaw = bmiss
1452 
1453  do krad = 1, max_search
1454 
1455  istart = itile - krad
1456  iend = itile + krad
1457  jstart = jtile - krad
1458  jend = jtile + krad
1459 
1460  do jj = jstart, jend
1461  do ii = istart, iend
1462 
1463 
1464  if ((jj == jstart) .or. (jj == jend) .or. &
1465  (ii == istart) .or. (ii == iend)) then
1466 
1467  if ((jj >= 1) .and. (jj <= ny)) then
1468  jjj = jj
1469  if (ii <= 0) then
1470  iii = nx + ii
1471  else if (ii >= (nx+1)) then
1472  iii = ii - nx
1473  else
1474  iii = ii
1475  endif
1476 
1477 !----------------------------------------------------------------------
1478 ! SEE IF NEARBY POINTS ARE SEA ICE. IF THEY ARE, AND THE SEARCH FOR
1479 ! A GAUSSIAN GRID OPEN WATER POINT FAILS, THEN TREF WILL BE SET TO
1480 ! FREEZING BELOW.
1481 !----------------------------------------------------------------------
1482  if (krad <= 2 .and. mask(iii,jjj) == 2) is_ice = .true.
1483 
1484  if (mask(iii,jjj) == 0) then
1485  tf_thaw = tf(iii,jjj)
1486  nset_thaw_s = nset_thaw_s + 1
1487  write(*,'(a,I4,2F9.3)') 'nset_thaw_s,tf(iii,jjj),tclm : ',nset_thaw_s,tf(iii,jjj),tclm
1488  go to 100
1489  endif ! tile mask is open water
1490 
1491  endif
1492  endif
1493  enddo
1494  enddo
1495  enddo ! krad loop
1496 
1497  100 continue
1498 
1499  if ( tf_thaw == bmiss ) then
1500  if (is_ice) then
1501  tf_thaw = tice
1502  nset_thaw_i = nset_thaw_i + 1
1503  write(*,'(a,I4,F9.3)') 'nset_thaw_i,tf_ice : ',nset_thaw_i,tice
1504  else
1505  tf_thaw = 0.8*tice+0.2*tclm
1506  nset_thaw_c = nset_thaw_c + 1
1507  write(*,'(a,I4,2F9.3)') 'nset_thaw_c,tf_ice,tclm : ',nset_thaw_c,tice,tclm
1508  endif
1509  endif
1510 
1511  end subroutine tf_thaw_set
1512 
1520  subroutine nsst_water_reset(nsst,ij,tf_thaw)
1521  use read_write_data, only : nsst_data
1522  implicit none
1523 
1524  integer, intent(in) :: ij
1525 
1526  real, intent(in) :: tf_thaw
1527 
1528  type(nsst_data), intent(inout) :: nsst
1529 
1530  nsst%c_0(ij) = 0.0
1531  nsst%c_d(ij) = 0.0
1532  nsst%d_conv(ij) = 0.0
1533  nsst%dt_cool(ij) = 0.0
1534  nsst%ifd(ij) = 0.0
1535  nsst%qrain(ij) = 0.0
1536  nsst%tref(ij) = tf_thaw
1537  nsst%w_0(ij) = 0.0
1538  nsst%w_d(ij) = 0.0
1539  nsst%xs(ij) = 0.0
1540  nsst%xt(ij) = 0.0
1541  nsst%xtts(ij) = 0.0
1542  nsst%xu(ij) = 0.0
1543  nsst%xv(ij) = 0.0
1544  nsst%xz(ij) = 30.0
1545  nsst%xzts(ij) = 0.0
1546  nsst%z_c(ij) = 0.0
1547  nsst%zm(ij) = 0.0
1548 
1549  end subroutine nsst_water_reset
1550 
1566 subroutine get_tf_clm(xlats_ij,xlons_ij,ny,nx,iy,im,id,ih,tf_clm,tf_trd)
1567  use read_write_data, only : get_tf_clm_dim
1568 
1569  implicit none
1570 
1571  real, dimension(nx*ny), intent(in) :: xlats_ij
1572  real, dimension(nx*ny), intent(in) :: xlons_ij
1573  real, dimension(nx,ny), intent(out) :: tf_clm
1574  real, dimension(nx,ny), intent(out) :: tf_trd
1575  integer, intent(in) :: iy,im,id,ih,nx,ny
1576 ! local declare
1577  real, allocatable, dimension(:,:) :: tf_clm0 ! sst climatology at the valid time (nxc,nyc)
1578  real, allocatable, dimension(:,:) :: tf_trd0 ! 6-hourly sst climatology tendency valid at atime (nxc,nyc)
1579  real, allocatable, dimension(:) :: cxlats ! latitudes of sst climatology
1580  real, allocatable, dimension(:) :: cxlons ! longitudes of sst climatology
1581 
1582  real, dimension(nx*ny) :: tf_clm_ij ! sst climatology at target grids (nx*ny)
1583  real, dimension(nx*ny) :: tf_trd_ij ! 6-hourly sst climatology tendency
1584  real :: wei1,wei2
1585  integer :: nxc,nyc,mon1,mon2,i,j
1586  character (len=6), parameter :: fin_tf_clm='sstclm' ! sst climatology file name
1587 !
1588 ! get which two months used and their weights from atime
1589 !
1590  call get_tim_wei(iy,im,id,ih,mon1,mon2,wei1,wei2)
1591 !
1592 ! get the dimensions of the sst climatology & allocate the related arrays
1593 !
1594  call get_tf_clm_dim(fin_tf_clm,nyc,nxc)
1595  allocate( tf_clm0(nxc,nyc),tf_trd0(nxc,nyc),cxlats(nyc),cxlons(nxc) )
1596 !
1597 ! get tf_clm at the analysis time from monthly climatology & cxlats, cxlons
1598 !
1599  call get_tf_clm_ta(tf_clm0,tf_trd0,cxlats,cxlons,nyc,nxc,mon1,mon2,wei1,wei2)
1600 !
1601 ! get tf_clm (nx by ny lat/lon) valid at atime
1602 !
1603  if ( nx == nxc .and. ny == nyc ) then
1604  tf_clm(:,:) = tf_clm0(:,:)
1605  tf_trd(:,:) = tf_trd0(:,:)
1606 ! write(*,'(a,2f9.3)') 'same dimensions, tf_clm, min : ',minval(tf_clm),maxval(tf_clm)
1607  else
1608 ! write(*,'(a,4i8)') 'different dimensions,nx,ny,nxc,nyc : ',nx,ny,nxc,nyc
1609  call intp_tile(tf_clm0, cxlats, cxlons, nyc, nxc, &
1610  tf_clm_ij,xlats_ij,xlons_ij,ny, nx)
1611  call intp_tile(tf_trd0, cxlats, cxlons, nyc, nxc, &
1612  tf_trd_ij,xlats_ij,xlons_ij,ny, nx)
1613 ! write(*,'(a,2f9.3)') 'tf_clm0, min, max : ',minval(tf_clm0),maxval(tf_clm0)
1614 
1615  tf_clm(:,:) = reshape(tf_clm_ij, (/nx,ny/) )
1616  tf_trd(:,:) = reshape(tf_trd_ij, (/nx,ny/) )
1617  endif
1618 
1619 end subroutine get_tf_clm
1620 
1635 subroutine get_tf_clm_ta(tf_clm_ta,tf_clm_trend,xlats,xlons,nlat,nlon,mon1,mon2,wei1,wei2)
1636  use read_write_data, only : read_tf_clim_grb
1637  implicit none
1638 
1639 ! input
1640  integer, intent(in) :: nlat,nlon,mon1,mon2
1641  real, intent(in) :: wei1,wei2
1642 ! output
1643  real, dimension(nlon,nlat), intent(out) :: tf_clm_ta,tf_clm_trend
1644  real, dimension(nlat), intent(out) :: xlats
1645  real, dimension(nlon), intent(out) :: xlons
1646 
1647 !input/output data file names
1648  character (len=6), parameter :: fin_tf_clm='sstclm'
1649 
1650 ! local declare
1651  real, dimension(nlon,nlat) :: tf_clm1,tf_clm2
1652 
1653 !
1654 ! read in rtg sst climatology without bitmap (surface mask) for mon1 and mon2
1655 !
1656  call read_tf_clim_grb(trim(fin_tf_clm),tf_clm1,xlats,xlons,nlat,nlon,mon1)
1657  call read_tf_clim_grb(trim(fin_tf_clm),tf_clm2,xlats,xlons,nlat,nlon,mon2)
1658 !
1659 ! tf_clim at the analysis time
1660 !
1661  tf_clm_ta(:,:) = wei1*tf_clm1(:,:)+wei2*tf_clm2(:,:)
1662 !
1663 ! tf_clim trend at the analysis time (month to month)
1664 !
1665  tf_clm_trend(:,:) = (tf_clm2(:,:)-tf_clm1(:,:))/120.0
1666 
1667  write(*,'(a,2f9.3)') 'tf_clm_ta, min, max : ',minval(tf_clm_ta),maxval(tf_clm_ta)
1668  write(*,'(a,2f9.3)') 'tf_clm_trend, min, max : ',minval(tf_clm_trend),maxval(tf_clm_trend)
1669  end subroutine get_tf_clm_ta
1670 
1683 subroutine get_sal_clm(xlats_ij,xlons_ij,ny,nx,iy,im,id,ih,sal_clm)
1684  use read_write_data, only : get_dim_nc
1685  implicit none
1686 
1687  real, dimension(nx*ny), intent(in) :: xlats_ij !
1688  real, dimension(nx*ny), intent(in) :: xlons_ij !
1689  real, dimension(nx,ny), intent(out) :: sal_clm !
1690  integer, intent(in) :: iy,im,id,ih,nx,ny
1691 ! local declare
1692  real, allocatable, dimension(:,:) :: sal_clm0 ! salinity climatology at the valid time
1693  real, allocatable, dimension(:) :: cxlats ! latitudes of sst climatology
1694  real, allocatable, dimension(:) :: cxlons ! longitudes of sst climatology
1695 
1696  real, dimension(nx*ny) :: sal_clm_ij ! salinity climatology at target grids (nx*ny)
1697  real :: wei1,wei2
1698  integer :: nxc,nyc,mon1,mon2,i,j
1699  character (len=6), parameter :: fin_sal_clm='salclm' ! salinity climatology file name
1700 !
1701 ! get which two months used and their weights from atime
1702 !
1703  call get_tim_wei(iy,im,id,ih,mon1,mon2,wei1,wei2)
1704 !
1705 ! get the dimensions of the sst climatology & allocate the related arrays
1706 !
1707  call get_dim_nc(fin_sal_clm,nyc,nxc)
1708  allocate( sal_clm0(nxc,nyc),cxlats(nyc),cxlons(nxc) )
1709 !
1710 ! get sal_clm at the analysis time from monthly climatology & cxlats, cxlons
1711 !
1712  call get_sal_clm_ta(sal_clm0,cxlats,cxlons,nyc,nxc,mon1,mon2,wei1,wei2)
1713 !
1714 ! get sal_clm (nx by ny lat/lon) valid at atime
1715 !
1716  if ( nx == nxc .and. ny == nyc ) then
1717  sal_clm(:,:) = sal_clm0(:,:)
1718 ! write(*,'(a,2f9.3)') 'same dimensions, sal_clm, min : ',minval(sal_clm),maxval(sal_clm)
1719  else
1720 ! write(*,'(a,4i8)') 'different dimensions,nx,ny,nxc,nyc : ',nx,ny,nxc,nyc
1721  call intp_tile(sal_clm0, cxlats, cxlons, nyc, nxc, &
1722  sal_clm_ij,xlats_ij,xlons_ij,ny, nx)
1723 ! write(*,'(a,2f9.3)') 'sal_clm0, min, max : ',minval(sal_clm0),maxval(sal_clm0)
1724 ! write(*,'(a,2f9.3)') 'done with intp_tile for sal_clm, min, max : ',minval(sal_clm_ij),maxval(sal_clm_ij)
1725 
1726  sal_clm(:,:) = reshape(sal_clm_ij, (/nx,ny/) )
1727  endif
1728 
1729 end subroutine get_sal_clm
1730 
1743 subroutine get_sal_clm_ta(sal_clm_ta,xlats,xlons,nlat,nlon,mon1,mon2,wei1,wei2)
1745  use read_write_data, only : read_salclm_gfs_nc
1746  implicit none
1747 
1748 ! input
1749  integer, intent(in) :: nlat,nlon,mon1,mon2
1750  real, intent(in) :: wei1,wei2
1751 ! output
1752  real, dimension(nlon,nlat), intent(out) :: sal_clm_ta
1753  real, dimension(nlat), intent(out) :: xlats
1754  real, dimension(nlon), intent(out) :: xlons
1755 
1756 !input/output data file names
1757  character (len=6), parameter :: fin_sal_clm='salclm'
1758 
1759 ! local declare
1760  real, dimension(nlon,nlat) :: sal_clm1,sal_clm2
1761 
1762 !
1763 ! read in rtg sst climatology without bitmap (surface mask) for mon1 and mon2
1764 !
1765  call read_salclm_gfs_nc(trim(fin_sal_clm),sal_clm1,xlats,xlons,nlat,nlon,mon1)
1766  call read_salclm_gfs_nc(trim(fin_sal_clm),sal_clm2,xlats,xlons,nlat,nlon,mon2)
1767 !
1768 ! sal_clim at the analysis time
1769 !
1770  sal_clm_ta(:,:) = wei1*sal_clm1(:,:)+wei2*sal_clm2(:,:)
1771  write(*,'(a,2f9.3)') 'sal_clm_ta, min, max : ',minval(sal_clm_ta),maxval(sal_clm_ta)
1772  end subroutine get_sal_clm_ta
1773 
1788 subroutine intp_tile(tf_lalo,dlats_lalo,dlons_lalo,jdim_lalo,idim_lalo, &
1789  tf_tile,xlats_tile,xlons_tile,jdim_tile,idim_tile)
1791  use utils
1792 
1793  implicit none
1794 
1795 ! input/output
1796  real, dimension(idim_lalo,jdim_lalo), intent(in) :: tf_lalo
1797  real, dimension(jdim_lalo), intent(in) :: dlats_lalo
1798  real, dimension(idim_lalo), intent(in) :: dlons_lalo
1799  real, dimension(jdim_tile*idim_tile), intent(in) :: xlats_tile
1800  real, dimension(jdim_tile*idim_tile), intent(in) :: xlons_tile
1801  integer, intent(in) :: jdim_lalo,idim_lalo,jdim_tile,idim_tile
1802  real, dimension(jdim_tile*idim_tile), intent(out) :: tf_tile
1803 
1804 ! local
1805  real, parameter :: deg2rad=3.1415926/180.0
1806  real, dimension(jdim_lalo) :: xlats_lalo
1807  real, dimension(idim_lalo) :: xlons_lalo
1808  real :: tf,wsum,res_km
1809  integer :: itile,jtile
1810  integer :: ii,jj,ij,iii,jjj
1811  integer :: ilalo,jlalo,ilalop1,jlalop1
1812  integer :: istart,iend,jstart,jend,krad
1813 
1814  integer, allocatable, dimension(:,:) :: id1,id2,jdc
1815  real, allocatable, dimension(:,:,:) :: agrid,s2c
1816 
1817  print*
1818  print*,'interpolate from lat/lon grids to any one grid with known lat/lon'
1819 
1820  xlats_lalo = dlats_lalo*deg2rad
1821  xlons_lalo = dlons_lalo*deg2rad
1822 
1823  allocate(agrid(idim_tile,jdim_tile,2))
1824  agrid(:,:,1) = reshape(xlons_tile, (/idim_tile,jdim_tile/) )
1825  agrid(:,:,2) = reshape(xlats_tile, (/idim_tile,jdim_tile/) )
1826  agrid = agrid*deg2rad
1827 
1828  allocate(id1(idim_tile,jdim_tile))
1829  allocate(id2(idim_tile,jdim_tile))
1830  allocate(jdc(idim_tile,jdim_tile))
1831  allocate(s2c(idim_tile,jdim_tile,4))
1832 
1833 !----------------------------------------------------------------------
1834 ! compute bilinear weights for each model point from the nearest
1835 ! four lalo points. does not account for mask. that
1836 ! happens later.
1837 !----------------------------------------------------------------------
1838 
1839  call remap_coef( 1, idim_tile, 1, jdim_tile, idim_lalo, jdim_lalo, &
1840  xlons_lalo, xlats_lalo, id1, id2, jdc, s2c, agrid )
1841 
1842  do ij = 1, jdim_tile*idim_tile
1843 
1844  jtile = (ij-1)/idim_tile + 1
1845  itile = mod(ij,idim_tile)
1846  if (itile==0) itile = idim_tile
1847 
1848  ilalo = id1(itile,jtile)
1849  ilalop1 = id2(itile,jtile)
1850  jlalo = jdc(itile,jtile)
1851  jlalop1 = jdc(itile,jtile) + 1
1852 
1853  wsum = s2c(itile,jtile,1) + s2c(itile,jtile,2) + &
1854  s2c(itile,jtile,3) + s2c(itile,jtile,4)
1855 
1856  tf_tile(ij) = ( s2c(itile,jtile,1)*tf_lalo(ilalo,jlalo) + &
1857  s2c(itile,jtile,2)*tf_lalo(ilalop1,jlalo) + &
1858  s2c(itile,jtile,3)*tf_lalo(ilalop1,jlalop1) + &
1859  s2c(itile,jtile,4)*tf_lalo(ilalo,jlalop1) )/wsum
1860  enddo
1861 
1862  deallocate(id1, id2, jdc, s2c)
1863 
1864 end subroutine intp_tile
1865 
1878 subroutine get_tim_wei(iy,im,id,ih,mon1,mon2,wei1,wei2)
1879  implicit none
1880 
1881 ! input
1882  integer, intent(in) :: iy,im,id,ih
1883 ! output
1884  integer, intent(out) :: mon1,mon2
1885  real, intent(out) :: wei1,wei2
1886 
1887 ! local declare
1888  real :: rjday
1889  integer :: mon,monend,monm,monp,jdow,jdoy,jday
1890  integer :: jda(8)
1891 !
1892 !dayhf : julian day of the middle of each month
1893 !
1894  real, dimension(13) :: dayhf
1895  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/
1896 
1897 ! 15, 44, 73.5, 104, 134.5, 165, 195.5, 226.5, 257, 287.5, 318.5, 349 ! from
1898 ! woa05
1899 
1900  jda=0
1901  jda(1)=iy
1902  jda(2)=im
1903  jda(3)=id
1904  jda(5)=ih
1905  jdow = 0
1906  jdoy = 0
1907  jday = 0
1908  call w3doxdat(jda,jdow,jdoy,jday)
1909  rjday=jdoy+jda(5)/24.
1910  if(rjday.lt.dayhf(1)) rjday=rjday+365.
1911 
1912  monend = 12
1913  do mon = 1, monend
1914  monm = mon
1915  monp = mon + 1
1916  if( rjday >= dayhf(monm) .and. rjday < dayhf(monp) ) then
1917  mon1 = monm
1918  mon2 = monp
1919  go to 10
1920  endif
1921  enddo
1922 
1923  print *,'FATAL ERROR in get_tim_wei, wrong rjday',rjday
1924  call abort
1925  10 continue
1926 
1927  wei1 = (dayhf(mon2)-rjday)/(dayhf(mon2)-dayhf(mon1))
1928  wei2 = (rjday-dayhf(mon1))/(dayhf(mon2)-dayhf(mon1))
1929 
1930  if( mon2 == 13 ) mon2=1
1931 
1932  write(*,'(a,2i4,3f9.3)') 'mon1,mon2,rjday,wei1,wei2=',mon1,mon2,rjday,wei1,wei2
1933 
1934  end subroutine get_tim_wei
1935 
1945  real function tfreez(salinity)
1947  implicit none
1948 
1949  REAL salinity,sal
1950  REAL a1, a2, a3
1951  parameter(a1 = -0.0575)
1952  parameter(a2 = 1.710523e-3)
1953  parameter(a3 = -2.154996e-4)
1954 
1955  IF (salinity .LT. 0.) THEN
1956  sal = 0.
1957  ELSE
1958  sal = salinity
1959  ENDIF
1960  tfreez = sal*(a1+a2*sqrt(sal)+a3*sal)
1961 
1962  return
1963  end
integer function num_parthds()
Return the number of omp threads.
Definition: num_parthds.f90:10
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.
Definition: cycle.f90:1567
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.
Definition: cycle.f90:1636
subroutine get_sal_clm_ta(sal_clm_ta, xlats, xlons, nlat, nlon, mon1, mon2, wei1, wei2)
Get climatological salinity at the analysis time.
Definition: cycle.f90:1744
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.
Definition: cycle.f90:1431
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).
Definition: cycle.f90:1790
Module containing utility routines.
Definition: utils.F90:7
real function tfreez(salinity)
Compute the freezing point of water as a function of salinity.
Definition: cycle.f90:1946
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...
Definition: cycle.f90:1879
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.
Definition: cycle.f90:1684
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.
Definition: cycle.f90:1521