95 CHARACTER(LEN=3) :: donst
96 INTEGER :: idim, jdim, lsoil, lugb, iy, im, id, ih, ialb
97 INTEGER :: isot, ivegsrc, lensfc, zsea1_mm, zsea2_mm, ierr
98 INTEGER :: nprocs, myrank, num_threads,
num_parthds, max_tasks
99 REAL :: fh, deltsfc, zsea1, zsea2
100 LOGICAL :: use_ufo, do_nsst, adjt_nst_only
102 namelist/namcyc/ idim,jdim,lsoil,lugb,iy,im,id,ih,fh, &
103 deltsfc,ialb,use_ufo,donst, &
104 adjt_nst_only,isot,ivegsrc,zsea1_mm, &
107 DATA idim,jdim,lsoil/96,96,4/
108 DATA iy,im,id,ih,fh/1997,8,2,0,0./
109 DATA lugb/51/, deltsfc/0.0/, ialb/1/, max_tasks/99999/
110 DATA isot/1/, ivegsrc/2/, zsea1_mm/0/, zsea2_mm/0/
113 CALL mpi_comm_size(mpi_comm_world, nprocs, ierr)
114 CALL mpi_comm_rank(mpi_comm_world, myrank, ierr)
116 if (myrank==0) call w3tagb(
'GLOBAL_CYCLE',2018,0179,0055,
'NP20')
121 print*,
"STARTING CYCLE PROGRAM ON RANK ", myrank
122 print*,
"RUNNING WITH ", nprocs,
"TASKS"
123 print*,
"AND WITH ", num_threads,
" THREADS."
127 adjt_nst_only = .false.
130 print*,
"READ NAMCYC NAMELIST."
132 CALL baopenr(36,
"fort.36", ierr)
134 IF (myrank==0)
WRITE(6,namcyc)
136 IF (max_tasks < 99999 .AND. myrank > (max_tasks - 1))
THEN
137 print*,
"USER SPECIFIED MAX NUMBER OF TASKS: ", max_tasks
138 print*,
"WILL NOT RUN CYCLE PROGRAM ON RANK: ", myrank
144 zsea1 = float(zsea1_mm) / 1000.0
145 zsea2 = float(zsea2_mm) / 1000.0
147 IF (donst ==
"YES")
THEN
154 IF (myrank==0) print*,
"LUGB,IDIM,JDIM,LSOIL,DELTSFC,IY,IM,ID,IH,FH: ", &
155 lugb,idim,jdim,lsoil,deltsfc,iy,im,id,ih,fh
157 CALL
sfcdrv(lugb,idim,jdim,lensfc,lsoil,deltsfc, &
158 iy,im,id,ih,fh,ialb, &
159 use_ufo,do_nsst,adjt_nst_only, &
160 zsea1,zsea2,isot,ivegsrc,myrank)
163 print*,
'CYCLE PROGRAM COMPLETED NORMALLY ON RANK: ', myrank
167 CALL mpi_barrier(mpi_comm_world, ierr)
169 if (myrank==0) call w3tage(
'GLOBAL_CYCLE')
171 CALL mpi_finalize(ierr)
283 SUBROUTINE sfcdrv(LUGB, IDIM,JDIM,LENSFC,LSOIL,DELTSFC, &
284 iy,im,id,ih,fh,ialb, &
285 use_ufo,do_nsst,adjt_nst_only, &
286 zsea1,zsea2,isot,ivegsrc,myrank)
292 INTEGER,
INTENT(IN) :: idim, jdim, lensfc, lsoil, ialb
293 INTEGER,
INTENT(IN) :: lugb, iy, im, id, ih
294 INTEGER,
INTENT(IN) :: isot, ivegsrc, myrank
296 LOGICAL,
INTENT(IN) :: use_ufo, do_nsst, adjt_nst_only
298 REAL,
INTENT(IN) :: fh, deltsfc, zsea1, zsea2
300 INTEGER,
PARAMETER :: nlunit=35
301 INTEGER,
PARAMETER :: sz_nml=1
303 CHARACTER(LEN=5) :: tile_num
304 CHARACTER(LEN=500) :: gsi_file
305 CHARACTER(LEN=4) :: input_nml_file(sz_nml)
308 INTEGER :: i_index(lensfc), j_index(lensfc)
309 INTEGER :: idum(idim,jdim)
311 REAL :: slmask(lensfc), orog(lensfc)
312 REAL :: sihfcs(lensfc), sicfcs(lensfc)
313 REAL :: sitfcs(lensfc), tsffcs(lensfc)
314 REAL :: snofcs(lensfc), zorfcs(lensfc)
315 REAL :: albfcs(lensfc,4), tg3fcs(lensfc)
316 REAL :: cnpfcs(lensfc), smcfcs(lensfc,lsoil)
317 REAL :: stcfcs(lensfc,lsoil), slifcs(lensfc)
318 REAL :: aisfcs(lensfc), f10m(lensfc)
319 REAL :: vegfcs(lensfc), vetfcs(lensfc)
320 REAL :: sotfcs(lensfc), alffcs(lensfc,2)
321 REAL :: cvfcs(lensfc), cvtfcs(lensfc)
322 REAL :: cvbfcs(lensfc), tprcp(lensfc)
323 REAL :: srflag(lensfc), swdfcs(lensfc)
324 REAL :: slcfcs(lensfc,lsoil), vmxfcs(lensfc)
325 REAL :: vmnfcs(lensfc), t2m(lensfc)
326 REAL :: q2m(lensfc), slpfcs(lensfc)
327 REAL :: absfcs(lensfc), orog_uf(lensfc)
328 REAL :: ustar(lensfc)
329 REAL :: fmm(lensfc), fhh(lensfc)
330 REAL :: rla(lensfc), rlo(lensfc)
331 REAL(KIND=4) :: zsoil(lsoil)
332 REAL :: sig1t(lensfc)
335 REAL,
ALLOCATABLE :: slifcs_fg(:)
338 real,
dimension(idim,jdim) :: tf_clm,tf_trd,sal_clm
339 real,
dimension(lensfc) :: tf_clm_tile,tf_trd_tile,sal_clm_tile
346 DATA gsi_file/
'NULL'/
348 namelist/namsfcd/ gsi_file
352 input_nml_file =
"NULL"
354 CALL baopenr(37,
"fort.37", ierr)
355 READ (37, nml=namsfcd)
359 print*,
'IN ROUTINE SFCDRV,IDIM=',idim,
'JDIM=',jdim,
'FH=',fh
371 i_index = reshape(idum, (/lensfc/))
377 j_index = reshape(idum, (/lensfc/) )
381 print*,
"WILL PROCESS NSST RECORDS."
382 ALLOCATE(nsst%C_0(lensfc))
383 ALLOCATE(nsst%C_D(lensfc))
384 ALLOCATE(nsst%D_CONV(lensfc))
385 ALLOCATE(nsst%DT_COOL(lensfc))
386 ALLOCATE(nsst%IFD(lensfc))
387 ALLOCATE(nsst%QRAIN(lensfc))
388 ALLOCATE(nsst%TREF(lensfc))
389 ALLOCATE(nsst%TFINC(lensfc))
390 ALLOCATE(nsst%W_0(lensfc))
391 ALLOCATE(nsst%W_D(lensfc))
392 ALLOCATE(nsst%XS(lensfc))
393 ALLOCATE(nsst%XT(lensfc))
394 ALLOCATE(nsst%XTTS(lensfc))
395 ALLOCATE(nsst%XU(lensfc))
396 ALLOCATE(nsst%XV(lensfc))
397 ALLOCATE(nsst%XZ(lensfc))
398 ALLOCATE(nsst%XZTS(lensfc))
399 ALLOCATE(nsst%Z_C(lensfc))
400 ALLOCATE(nsst%ZM(lensfc))
401 ALLOCATE(slifcs_fg(lensfc))
408 CALL
read_data(tsffcs,smcfcs,snofcs,stcfcs,tg3fcs,zorfcs, &
409 cvfcs,cvbfcs,cvtfcs,albfcs,slifcs, &
410 vegfcs,cnpfcs,f10m,vetfcs,sotfcs, &
411 alffcs,ustar,fmm,fhh,sihfcs,sicfcs, &
412 sitfcs,tprcp,srflag,swdfcs,vmnfcs, &
413 vmxfcs,slcfcs,slpfcs,absfcs,t2m,q2m, &
414 slmask,zsoil,lsoil,lensfc,do_nsst,nsst)
418 print*,
'USE UNFILTERED OROGRAPHY.'
425 IF(nint(slifcs(i)).EQ.2) aisfcs(i) = 1.
429 IF (adjt_nst_only)
THEN
431 print*,
"FIRST GUESS MASK ADJUSTED BY IFD RECORD"
433 WHERE(nint(nsst%IFD) == 3) slifcs_fg = 2.0
436 print*,
"SAVE FIRST GUESS MASK"
445 IF (.NOT. adjt_nst_only)
THEN
447 print*,
"CALL SFCCYCLE TO UPDATE SURFACE FIELDS."
448 CALL
sfccycle(lugb,lensfc,lsoil,sig1t,deltsfc, &
449 iy,im,id,ih,fh,rla,rlo, &
450 slmask,orog, orog_uf, use_ufo, do_nsst, &
451 sihfcs,sicfcs,sitfcs,swdfcs,slcfcs, &
452 vmnfcs,vmxfcs,slpfcs,absfcs, &
453 tsffcs,snofcs,zorfcs,albfcs,tg3fcs, &
454 cnpfcs,smcfcs,stcfcs,slifcs,aisfcs, &
455 vegfcs,vetfcs,sotfcs,alffcs, &
456 cvfcs,cvbfcs,cvtfcs,myrank,nlunit, &
457 sz_nml, input_nml_file, &
458 ialb,isot,ivegsrc,tile_num,i_index,j_index)
468 IF (gsi_file ==
"NULL")
THEN
470 print*,
"NO GSI FILE. ADJUST IFD FOR FORMER ICE POINTS."
472 IF (nint(slifcs_fg(i)) == 2 .AND. nint(slifcs(i)) == 0)
THEN
479 print*,
"ADJUST TREF FROM GSI INCREMENT"
483 call
get_tf_clm(rla,rlo,jdim,idim,iy,im,id,ih,tf_clm,tf_trd)
484 tf_clm_tile(:) = reshape(tf_clm, (/lensfc/) )
485 tf_trd_tile(:) = reshape(tf_trd, (/lensfc/) )
489 call
get_sal_clm(rla,rlo,jdim,idim,iy,im,id,ih,sal_clm)
490 sal_clm_tile(:) = reshape(sal_clm, (/lensfc/) )
498 CALL
adjust_nsst(rla,rlo,slifcs,slifcs_fg,tsffcs,sitfcs,sicfcs,stcfcs, &
499 nsst,lensfc,lsoil,idim,jdim,zsea1,zsea2,im,id,deltsfc, &
500 tf_clm_tile,tf_trd_tile,sal_clm_tile)
508 CALL
write_data(slifcs,tsffcs,snofcs,tg3fcs,zorfcs, &
509 albfcs,alffcs,vegfcs,cnpfcs,f10m, &
510 t2m,q2m,vetfcs,sotfcs,ustar,fmm,fhh, &
511 sicfcs,sihfcs,sitfcs, &
512 tprcp,srflag,swdfcs, &
513 vmnfcs,vmxfcs,slpfcs,absfcs, &
514 slcfcs,smcfcs,stcfcs, &
515 idim,jdim,lensfc,lsoil,do_nsst,nsst)
520 DEALLOCATE(nsst%D_CONV)
521 DEALLOCATE(nsst%DT_COOL)
523 DEALLOCATE(nsst%QRAIN)
524 DEALLOCATE(nsst%TREF)
525 DEALLOCATE(nsst%TFINC)
530 DEALLOCATE(nsst%XTTS)
534 DEALLOCATE(nsst%XZTS)
537 DEALLOCATE(slifcs_fg)
575 SUBROUTINE adjust_nsst(RLA,RLO,SLMSK_TILE,SLMSK_FG_TILE,SKINT_TILE,&
576 sicet_tile,sice_tile,soilt_tile,nsst,lensfc,lsoil, &
577 idim,jdim,zsea1,zsea2,mon,day,deltsfc, &
578 tf_clm_tile,tf_trd_tile,sal_clm_tile)
582 slmsk_gaus, dtref_gaus, &
589 INTEGER,
INTENT(IN) :: lensfc, lsoil, idim, jdim, mon, day
591 REAL,
INTENT(IN) :: slmsk_tile(lensfc), slmsk_fg_tile(lensfc)
592 real,
intent(in) :: tf_clm_tile(lensfc),tf_trd_tile(lensfc),sal_clm_tile(lensfc)
593 REAL,
INTENT(IN) :: zsea1, zsea2, deltsfc
594 REAL,
INTENT(INOUT) :: rla(lensfc), rlo(lensfc), skint_tile(lensfc)
595 REAL,
INTENT(INOUT) :: sicet_tile(lensfc),sice_tile(lensfc),soilt_tile(lensfc,lsoil)
599 REAL,
PARAMETER :: tmax=313.0,tzero=273.16
601 INTEGER :: iopt, nret, kgds_gaus(200)
602 INTEGER :: igaus, jgaus, ij, ii, jj, iii, jjj, krad
603 INTEGER :: istart, iend, jstart, jend
604 INTEGER :: mask_tile, mask_fg_tile
605 INTEGER :: itile, jtile
606 INTEGER :: max_search, j, ierr
607 INTEGER :: igausp1, jgausp1
608 integer :: nintp,nsearched,nice,nland
609 integer :: nfill,nfill_tice,nfill_clm
610 integer :: nset_thaw,nset_thaw_s,nset_thaw_i,nset_thaw_c
612 INTEGER,
ALLOCATABLE :: id1(:,:), id2(:,:), jdc(:,:)
617 REAL :: tref_save,wsum,tf_ice,tf_thaw
618 REAL :: fill, dtzm, gaus_res_km, dtref
619 REAL,
ALLOCATABLE :: xpts(:), ypts(:), lats(:), lons(:)
620 REAL,
ALLOCATABLE :: dum2d(:,:), lats_rad(:), lons_rad(:)
621 REAL,
ALLOCATABLE :: agrid(:,:,:), s2c(:,:,:)
625 kgds_gaus(2) = idim_gaus
626 kgds_gaus(3) = jdim_gaus
630 kgds_gaus(7) = -90000
631 kgds_gaus(8) = nint(-360000./float(idim_gaus))
632 kgds_gaus(9) = nint((360.0 / float(idim_gaus))*1000.0)
634 kgds_gaus(10) = jdim_gaus/2
639 print*,
'ADJUST NSST USING GSI INCREMENTS ON GAUSSIAN GRID'
647 ALLOCATE(xpts(idim_gaus*jdim_gaus))
648 ALLOCATE(ypts(idim_gaus*jdim_gaus))
649 ALLOCATE(lats(idim_gaus*jdim_gaus))
650 ALLOCATE(lons(idim_gaus*jdim_gaus))
656 CALL gdswzd(kgds_gaus,iopt,(idim_gaus*jdim_gaus),fill,xpts,ypts,lons,lats,nret)
658 IF (nret /= (idim_gaus*jdim_gaus))
THEN
659 print*,
'FATAL ERROR: PROBLEM IN GDSWZD. STOP.'
660 CALL mpi_abort(mpi_comm_world, 12, ierr)
663 DEALLOCATE (xpts, ypts)
665 ALLOCATE(dum2d(idim_gaus,jdim_gaus))
666 dum2d = reshape(lats, (/idim_gaus,jdim_gaus/) )
669 ALLOCATE(lats_rad(jdim_gaus))
671 lats_rad(j) = dum2d(1,jdim_gaus-j+1) * 3.1415926 / 180.0
674 dum2d = reshape(lons, (/idim_gaus,jdim_gaus/) )
676 ALLOCATE(lons_rad(idim_gaus))
677 lons_rad = dum2d(:,1) * 3.1415926 / 180.0
681 ALLOCATE(agrid(idim,jdim,2))
682 agrid(:,:,1) = reshape(rlo, (/idim,jdim/) )
683 agrid(:,:,2) = reshape(rla, (/idim,jdim/) )
684 agrid = agrid * 3.1415926 / 180.0
686 ALLOCATE(id1(idim,jdim))
687 ALLOCATE(id2(idim,jdim))
688 ALLOCATE(jdc(idim,jdim))
689 ALLOCATE(s2c(idim,jdim,4))
697 CALL
remap_coef( 1, idim, 1, jdim, idim_gaus, jdim_gaus, &
698 lons_rad, lats_rad, id1, id2, jdc, s2c, agrid )
700 DEALLOCATE(lons_rad, lats_rad, agrid)
707 gaus_res_km = 360.0 / idim_gaus * 111.0
708 max_search = ceiling(500.0/gaus_res_km)
711 print*,
'MAXIMUM SEARCH IS ',max_search,
' GAUSSIAN POINTS.'
734 ij_loop :
DO ij = 1, lensfc
736 mask_tile = nint(slmsk_tile(ij))
737 mask_fg_tile = nint(slmsk_fg_tile(ij))
742 tf_ice =
tfreez(sal_clm_tile(ij)) + tzero
747 IF (mask_tile == 1)
THEN
755 if (mask_tile == 2)
then
757 skint_tile(ij)=(1.0-sice_tile(ij))*nsst%tref(ij)+sice_tile(ij)*sicet_tile(ij)
765 jtile = (ij-1) / idim + 1
767 IF (itile==0) itile = idim
775 IF (mask_fg_tile == 2 .AND. mask_tile == 0)
THEN
779 call
tf_thaw_set(nsst%tref,nint(slmsk_fg_tile),itile,jtile,tf_ice,tf_clm_tile(ij),tf_thaw,idim,jdim, &
780 nset_thaw_s,nset_thaw_i,nset_thaw_c)
782 nset_thaw = nset_thaw + 1
798 igaus = id1(itile,jtile)
799 jgaus = jdc(itile,jtile)
800 igausp1 = id2(itile,jtile)
801 jgausp1 = jdc(itile,jtile)+1
803 IF (slmsk_gaus(igaus,jgaus) == 0 .OR. &
804 slmsk_gaus(igausp1,jgaus) == 0 .OR. &
805 slmsk_gaus(igausp1,jgausp1) == 0 .OR. &
806 slmsk_gaus(igaus,jgausp1) == 0)
THEN
811 IF (slmsk_gaus(igaus,jgaus) == 0)
THEN
812 dtref = dtref + (s2c(itile,jtile,1) * dtref_gaus(igaus,jgaus))
813 wsum = wsum + s2c(itile,jtile,1)
816 IF (slmsk_gaus(igausp1,jgaus) == 0)
THEN
817 dtref = dtref + (s2c(itile,jtile,2) * dtref_gaus(igausp1,jgaus))
818 wsum = wsum + s2c(itile,jtile,2)
821 IF (slmsk_gaus(igausp1,jgausp1) == 0)
THEN
822 dtref = dtref + (s2c(itile,jtile,3) * dtref_gaus(igausp1,jgausp1))
823 wsum = wsum + s2c(itile,jtile,3)
826 IF (slmsk_gaus(igaus,jgausp1) == 0)
THEN
827 dtref = dtref + (s2c(itile,jtile,4) * dtref_gaus(igaus,jgausp1))
828 wsum = wsum + s2c(itile,jtile,4)
834 tref_save = nsst%TREF(ij)
835 nsst%TREF(ij) = nsst%TREF(ij) + dtref
836 nsst%TREF(ij) = max(nsst%TREF(ij), tf_ice)
837 nsst%TREF(ij) = min(nsst%TREF(ij), tmax)
838 nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
840 CALL
dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij), &
841 nsst%Z_C(ij),zsea1,zsea2,dtzm)
843 skint_tile(ij) = nsst%TREF(ij) + dtzm
844 skint_tile(ij) = max(skint_tile(ij), tf_ice)
845 skint_tile(ij) = min(skint_tile(ij), tmax)
847 sicet_tile(ij) = skint_tile(ij)
848 soilt_tile(ij,:) = skint_tile(ij)
859 DO krad = 1, max_search
861 istart = igaus - krad
863 jstart = jgaus - krad
869 IF((jj == jstart) .OR. (jj == jend) .OR. &
870 (ii == istart) .OR. (ii == iend))
THEN
872 IF ((jj >= 1) .AND. (jj <= jdim_gaus))
THEN
877 ELSE IF (ii >= (idim_gaus+1))
THEN
889 IF (krad <= 2 .AND. slmsk_gaus(iii,jjj) == 2) is_ice = .true.
891 IF (slmsk_gaus(iii,jjj) == 0)
THEN
895 nsearched = nsearched + 1
897 tref_save = nsst%TREF(ij)
898 nsst%TREF(ij ) = nsst%TREF(ij) + dtref_gaus(iii,jjj)
899 nsst%TREF(ij) = max(nsst%TREF(ij), tf_ice)
900 nsst%TREF(ij) = min(nsst%TREF(ij), tmax)
901 nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
903 CALL
dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij), &
904 nsst%Z_C(ij),zsea1,zsea2,dtzm)
906 skint_tile(ij) = nsst%TREF(ij) + dtzm
907 skint_tile(ij) = max(skint_tile(ij), tf_ice)
908 skint_tile(ij) = min(skint_tile(ij), tmax)
910 sicet_tile(ij) = skint_tile(ij)
911 soilt_tile(ij,:) = skint_tile(ij)
934 nsst%TREF(ij) = tf_ice
936 nfill_tice = nfill_tice + 1
938 tref_save = nsst%TREF(ij)
939 nsst%TREF(ij) = nsst%TREF(ij) + tf_trd_tile(ij)
940 nsst%TREF(ij) = max(nsst%TREF(ij), tf_ice)
941 nsst%TREF(ij) = min(nsst%TREF(ij), tmax)
942 nsst%TFINC(ij) = nsst%TREF(ij) - tref_save
944 nfill_clm = nfill_clm + 1
947 CALL
dtzm_point(nsst%XT(ij),nsst%XZ(ij),nsst%DT_COOL(ij), &
948 nsst%Z_C(ij),zsea1,zsea2,dtzm)
950 skint_tile(ij) = nsst%TREF(ij) + dtzm
951 skint_tile(ij) = max(skint_tile(ij), tf_ice)
952 skint_tile(ij) = min(skint_tile(ij), tmax)
954 sicet_tile(ij) = skint_tile(ij)
955 soilt_tile(ij,:) = skint_tile(ij)
961 write(*,
'(a)')
'statistics of grids number processed for tile : '
962 write(*,
'(a,I8)')
' nintp = ',nintp
963 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
964 write(*,
'(a,I8)')
' nsearched = ',nsearched
965 write(*,
'(a,3I6)')
' nfill,nfill_tice,nfill_clm = ',nfill,nfill_tice,nfill_clm
966 write(*,
'(a,I8)')
' nice = ',nice
967 write(*,
'(a,I8)')
' nland = ',nland
969 DEALLOCATE(id1, id2, jdc, s2c)
986 INTEGER,
INTENT(IN) :: mon, day
988 REAL,
INTENT(IN) :: latitude, deltsfc
989 REAL,
INTENT(OUT) :: dtref
991 INTEGER :: num_days(12), mon2, mon1
993 REAL,
TARGET :: sst_80_90(12)
994 REAL,
TARGET :: sst_70_80(12)
995 REAL,
TARGET :: sst_60_70(12)
996 REAL,
TARGET :: sst_50_60(12)
997 REAL,
TARGET :: sst_40_50(12)
998 REAL,
TARGET :: sst_30_40(12)
999 REAL,
TARGET :: sst_20_30(12)
1000 REAL,
TARGET :: sst_10_20(12)
1001 REAL,
TARGET :: sst_00_10(12)
1002 REAL,
TARGET :: sst_m10_00(12)
1003 REAL,
TARGET :: sst_m20_m10(12)
1004 REAL,
TARGET :: sst_m30_m20(12)
1005 REAL,
TARGET :: sst_m40_m30(12)
1006 REAL,
TARGET :: sst_m50_m40(12)
1007 REAL,
TARGET :: sst_m60_m50(12)
1008 REAL,
TARGET :: sst_m70_m60(12)
1009 REAL,
TARGET :: sst_m80_m70(12)
1010 REAL,
TARGET :: sst_m90_m80(12)
1012 REAL,
POINTER :: sst(:)
1014 DATA num_days /31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31/
1016 DATA sst_80_90 /271.466, 271.458, 271.448, 271.445, 271.519, 271.636, &
1017 272.023, 272.066, 272.001, 271.698, 271.510, 271.472/
1019 DATA sst_70_80 /272.149, 272.103, 272.095, 272.126, 272.360, 272.988, &
1020 274.061, 274.868, 274.415, 273.201, 272.468, 272.268/
1022 DATA sst_60_70 /274.240, 274.019, 273.988, 274.185, 275.104, 276.875, &
1023 279.005, 280.172, 279.396, 277.586, 275.818, 274.803/
1025 DATA sst_50_60 /277.277, 276.935, 277.021, 277.531, 279.100, 281.357, &
1026 283.735, 285.171, 284.399, 282.328, 279.918, 278.199/
1028 DATA sst_40_50 /281.321, 280.721, 280.850, 281.820, 283.958, 286.588, &
1029 289.195, 290.873, 290.014, 287.652, 284.898, 282.735/
1031 DATA sst_30_40 /289.189, 288.519, 288.687, 289.648, 291.547, 293.904, &
1032 296.110, 297.319, 296.816, 295.225, 292.908, 290.743/
1034 DATA sst_20_30 /294.807, 294.348, 294.710, 295.714, 297.224, 298.703, &
1035 299.682, 300.127, 300.099, 299.455, 297.953, 296.177/
1037 DATA sst_10_20 /298.878, 298.720, 299.033, 299.707, 300.431, 300.709, &
1038 300.814, 300.976, 301.174, 301.145, 300.587, 299.694/
1040 DATA sst_00_10 /300.415, 300.548, 300.939, 301.365, 301.505, 301.141, &
1041 300.779, 300.660, 300.818, 300.994, 300.941, 300.675/
1043 DATA sst_m10_00 /300.226, 300.558, 300.914, 301.047, 300.645, 299.870, &
1044 299.114, 298.751, 298.875, 299.294, 299.721, 299.989/
1046 DATA sst_m20_m10 /299.547, 299.985, 300.056, 299.676, 298.841, 297.788, &
1047 296.893, 296.491, 296.687, 297.355, 298.220, 298.964/
1049 DATA sst_m30_m20 /297.524, 298.073, 297.897, 297.088, 295.846, 294.520, &
1050 293.525, 293.087, 293.217, 293.951, 295.047, 296.363/
1052 DATA sst_m40_m30 /293.054, 293.765, 293.468, 292.447, 291.128, 289.781, &
1053 288.773, 288.239, 288.203, 288.794, 289.947, 291.553/
1055 DATA sst_m50_m40 /285.052, 285.599, 285.426, 284.681, 283.761, 282.826, &
1056 282.138, 281.730, 281.659, 281.965, 282.768, 283.961/
1058 DATA sst_m60_m50 /277.818, 278.174, 277.991, 277.455, 276.824, 276.229, &
1059 275.817, 275.585, 275.560, 275.687, 276.142, 276.968/
1061 DATA sst_m70_m60 /273.436, 273.793, 273.451, 272.813, 272.349, 272.048, &
1062 271.901, 271.838, 271.845, 271.889, 272.080, 272.607/
1064 DATA sst_m80_m70 /271.579, 271.578, 271.471, 271.407, 271.392, 271.391, &
1065 271.390, 271.391, 271.394, 271.401, 271.422, 271.486/
1067 DATA sst_m90_m80 /271.350, 271.350, 271.350, 271.350, 271.350, 271.350, &
1068 271.350, 271.350, 271.350, 271.350, 271.350, 271.350/
1071 IF (latitude > 80.0)
THEN
1073 ELSEIF (latitude > 70.0)
THEN
1075 ELSEIF (latitude > 60.0)
THEN
1077 ELSEIF (latitude > 50.0)
THEN
1079 ELSEIF (latitude > 40.0)
THEN
1081 ELSEIF (latitude > 30.0)
THEN
1083 ELSEIF (latitude > 20.0)
THEN
1085 ELSEIF (latitude > 10.0)
THEN
1087 ELSEIF (latitude > 0.0)
THEN
1089 ELSEIF (latitude > -10.0)
THEN
1091 ELSEIF (latitude > -20.0)
THEN
1093 ELSEIF (latitude > -30.0)
THEN
1095 ELSEIF (latitude > -40.0)
THEN
1097 ELSEIF (latitude > -50.0)
THEN
1099 ELSEIF (latitude > -60.0)
THEN
1101 ELSEIF (latitude > -70.0)
THEN
1103 ELSEIF (latitude > -80.0)
THEN
1111 IF(mon2 == 13) mon2 = 1
1113 dtref = (sst(mon2) - sst(mon1)) / num_days(mon1)
1116 IF (mon1 == 0) mon1=12
1118 dtref = (sst(mon2) - sst(mon1)) / num_days(mon1)
1121 dtref = dtref * (deltsfc / 24.0)
1139 real,
intent(in) :: xt,xz,dt_cool,zc,z1,z2
1140 real,
intent(out) :: dtzm
1142 real,
parameter :: zero = 0.0
1143 real,
parameter :: one = 1.0
1144 real,
parameter :: half = 0.5
1145 real :: dt_warm,dtw,dtc
1150 if ( xt > zero )
then
1151 dt_warm = (xt+xt)/xz
1154 dtw = dt_warm*(one-(z1+z2)/(xz+xz))
1155 elseif ( z1 < xz .and. z2 >= xz )
then
1156 dtw = half*(one-z1/xz)*dt_warm*(xz-z1)/(z2-z1)
1158 elseif ( z1 == z2 )
then
1160 dtw = dt_warm*(one-z1/xz)
1168 if ( zc > zero )
then
1171 dtc = dt_cool*(one-(z1+z2)/(zc+zc))
1172 elseif ( z1 < zc .and. z2 >= zc )
then
1173 dtc = half*(one-z1/zc)*dt_cool*(zc-z1)/(z2-z1)
1175 elseif ( z1 == z2 )
then
1177 dtc = dt_cool*(one-z1/zc)
1215 im, jm, lon, lat, id1, id2, jdc, s2c, agrid )
1218 integer,
intent(in):: is, ie, js, je
1219 integer,
intent(in):: im, jm
1220 real,
intent(in):: lon(im), lat(jm)
1221 real,
intent(out):: s2c(is:ie,js:je,4)
1222 integer,
intent(out),
dimension(is:ie,js:je):: id1, id2, jdc
1223 real,
intent(in):: agrid(is:ie,js:je,2)
1228 real,
parameter :: pi = 3.1415926
1229 integer i,j, i1, i2, jc, i0, j0
1231 rdlon(i) = 1. / (lon(i+1) - lon(i))
1233 rdlon(im) = 1. / (lon(1) + 2.*pi - lon(im))
1236 rdlat(j) = 1. / (lat(j+1) - lat(j))
1244 if ( agrid(i,j,1)>lon(im) )
then
1246 a1 = (agrid(i,j,1)-lon(im)) * rdlon(im)
1247 elseif ( agrid(i,j,1)<lon(1) )
then
1249 a1 = (agrid(i,j,1)+2.*pi-lon(im)) * rdlon(im)
1252 if ( agrid(i,j,1)>=lon(i0) .and. agrid(i,j,1)<=lon(i0+1) )
then
1254 a1 = (agrid(i,j,1)-lon(i1)) * rdlon(i0)
1261 if ( agrid(i,j,2)<lat(1) )
then
1264 elseif ( agrid(i,j,2)>lat(jm) )
then
1269 if ( agrid(i,j,2)>=lat(j0) .and. agrid(i,j,2)<=lat(j0+1) )
then
1271 b1 = (agrid(i,j,2)-lat(jc)) * rdlat(jc)
1278 if ( a1<0.0 .or. a1>1.0 .or. b1<0.0 .or. b1>1.0 )
then
1279 write(*,*)
'gid=', i,j,a1, b1
1282 s2c(i,j,1) = (1.-a1) * (1.-b1)
1283 s2c(i,j,2) = a1 * (1.-b1)
1284 s2c(i,j,3) = a1 * b1
1285 s2c(i,j,4) = (1.-a1) * b1
1313 subroutine tf_thaw_set(tf_ij,mask_ij,itile,jtile,tice,tclm,tf_thaw,nx,ny, &
1314 nset_thaw_s,nset_thaw_i,nset_thaw_c)
1316 real,
dimension(nx*ny),
intent(in) :: tf_ij
1317 integer,
dimension(nx*ny),
intent(in) :: mask_ij
1318 real,
intent(in) :: tice,tclm
1319 integer,
intent(in) :: itile,jtile,nx,ny
1320 real,
intent(out) :: tf_thaw
1321 integer,
intent(inout) :: nset_thaw_s,nset_thaw_i,nset_thaw_c
1323 real,
parameter :: bmiss = -999.0
1324 real,
dimension(nx,ny) :: tf
1325 integer,
dimension(nx,ny) :: mask
1326 integer :: krad,max_search,istart,iend,jstart,jend
1327 integer :: ii,jj,iii,jjj
1332 mask(:,:) = reshape(mask_ij,(/nx,ny/) )
1333 tf(:,:) = reshape(tf_ij,(/nx,ny/) )
1337 do krad = 1, max_search
1339 istart = itile - krad
1341 jstart = jtile - krad
1344 do jj = jstart, jend
1345 do ii = istart, iend
1348 if ((jj == jstart) .or. (jj == jend) .or. &
1349 (ii == istart) .or. (ii == iend))
then
1351 if ((jj >= 1) .and. (jj <= ny))
then
1355 else if (ii >= (nx+1))
then
1366 if (krad <= 2 .and. mask(iii,jjj) == 2) is_ice = .true.
1368 if (mask(iii,jjj) == 0)
then
1369 tf_thaw = tf(iii,jjj)
1370 nset_thaw_s = nset_thaw_s + 1
1371 write(*,
'(a,I4,2F9.3)')
'nset_thaw_s,tf(iii,jjj),tclm : ',nset_thaw_s,tf(iii,jjj),tclm
1383 if ( tf_thaw == bmiss )
then
1386 nset_thaw_i = nset_thaw_i + 1
1387 write(*,
'(a,I4,F9.3)')
'nset_thaw_i,tf_ice : ',nset_thaw_i,tice
1389 tf_thaw = 0.8*tice+0.2*tclm
1390 nset_thaw_c = nset_thaw_c + 1
1391 write(*,
'(a,I4,2F9.3)')
'nset_thaw_c,tf_ice,tclm : ',nset_thaw_c,tice,tclm
1408 integer,
intent(in) :: ij
1410 real,
intent(in) :: tf_thaw
1416 nsst%d_conv(ij) = 0.0
1417 nsst%dt_cool(ij) = 0.0
1419 nsst%qrain(ij) = 0.0
1420 nsst%tref(ij) = tf_thaw
1450 subroutine get_tf_clm(xlats_ij,xlons_ij,ny,nx,iy,im,id,ih,tf_clm,tf_trd)
1455 real,
dimension(nx*ny),
intent(in) :: xlats_ij
1456 real,
dimension(nx*ny),
intent(in) :: xlons_ij
1457 real,
dimension(nx,ny),
intent(out) :: tf_clm
1458 real,
dimension(nx,ny),
intent(out) :: tf_trd
1459 integer,
intent(in) :: iy,im,id,ih,nx,ny
1461 real,
allocatable,
dimension(:,:) :: tf_clm0
1462 real,
allocatable,
dimension(:,:) :: tf_trd0
1463 real,
allocatable,
dimension(:) :: cxlats
1464 real,
allocatable,
dimension(:) :: cxlons
1466 real,
dimension(nx*ny) :: tf_clm_ij
1467 real,
dimension(nx*ny) :: tf_trd_ij
1469 integer :: nxc,nyc,mon1,mon2,i,j
1470 character (len=6),
parameter :: fin_tf_clm=
'sstclm'
1479 allocate( tf_clm0(nxc,nyc),tf_trd0(nxc,nyc),cxlats(nyc),cxlons(nxc) )
1483 call
get_tf_clm_ta(tf_clm0,tf_trd0,cxlats,cxlons,nyc,nxc,mon1,mon2,wei1,wei2)
1487 if ( nx == nxc .and. ny == nyc )
then
1488 tf_clm(:,:) = tf_clm0(:,:)
1489 tf_trd(:,:) = tf_trd0(:,:)
1493 call
intp_tile(tf_clm0, cxlats, cxlons, nyc, nxc, &
1494 tf_clm_ij,xlats_ij,xlons_ij,ny, nx)
1495 call
intp_tile(tf_trd0, cxlats, cxlons, nyc, nxc, &
1496 tf_trd_ij,xlats_ij,xlons_ij,ny, nx)
1499 tf_clm(:,:) = reshape(tf_clm_ij, (/nx,ny/) )
1500 tf_trd(:,:) = reshape(tf_trd_ij, (/nx,ny/) )
1519 subroutine get_tf_clm_ta(tf_clm_ta,tf_clm_trend,xlats,xlons,nlat,nlon,mon1,mon2,wei1,wei2)
1524 integer,
intent(in) :: nlat,nlon,mon1,mon2
1525 real,
intent(in) :: wei1,wei2
1527 real,
dimension(nlon,nlat),
intent(out) :: tf_clm_ta,tf_clm_trend
1528 real,
dimension(nlat),
intent(out) :: xlats
1529 real,
dimension(nlon),
intent(out) :: xlons
1532 character (len=6),
parameter :: fin_tf_clm=
'sstclm'
1535 real,
dimension(nlon,nlat) :: tf_clm1,tf_clm2
1545 tf_clm_ta(:,:) = wei1*tf_clm1(:,:)+wei2*tf_clm2(:,:)
1549 tf_clm_trend(:,:) = (tf_clm2(:,:)-tf_clm1(:,:))/120.0
1551 write(*,
'(a,2f9.3)')
'tf_clm_ta, min, max : ',minval(tf_clm_ta),maxval(tf_clm_ta)
1552 write(*,
'(a,2f9.3)')
'tf_clm_trend, min, max : ',minval(tf_clm_trend),maxval(tf_clm_trend)
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) :: sal_clm
1574 integer,
intent(in) :: iy,im,id,ih,nx,ny
1576 real,
allocatable,
dimension(:,:) :: sal_clm0
1577 real,
allocatable,
dimension(:) :: cxlats
1578 real,
allocatable,
dimension(:) :: cxlons
1580 real,
dimension(nx*ny) :: sal_clm_ij
1582 integer :: nxc,nyc,mon1,mon2,i,j
1583 character (len=6),
parameter :: fin_sal_clm=
'salclm'
1592 allocate( sal_clm0(nxc,nyc),cxlats(nyc),cxlons(nxc) )
1596 call
get_sal_clm_ta(sal_clm0,cxlats,cxlons,nyc,nxc,mon1,mon2,wei1,wei2)
1600 if ( nx == nxc .and. ny == nyc )
then
1601 sal_clm(:,:) = sal_clm0(:,:)
1605 call
intp_tile(sal_clm0, cxlats, cxlons, nyc, nxc, &
1606 sal_clm_ij,xlats_ij,xlons_ij,ny, nx)
1610 sal_clm(:,:) = reshape(sal_clm_ij, (/nx,ny/) )
1633 integer,
intent(in) :: nlat,nlon,mon1,mon2
1634 real,
intent(in) :: wei1,wei2
1636 real,
dimension(nlon,nlat),
intent(out) :: sal_clm_ta
1637 real,
dimension(nlat),
intent(out) :: xlats
1638 real,
dimension(nlon),
intent(out) :: xlons
1641 character (len=6),
parameter :: fin_sal_clm=
'salclm'
1644 real,
dimension(nlon,nlat) :: sal_clm1,sal_clm2
1654 sal_clm_ta(:,:) = wei1*sal_clm1(:,:)+wei2*sal_clm2(:,:)
1655 write(*,
'(a,2f9.3)')
'sal_clm_ta, min, max : ',minval(sal_clm_ta),maxval(sal_clm_ta)
1672 subroutine intp_tile(tf_lalo,dlats_lalo,dlons_lalo,jdim_lalo,idim_lalo, &
1673 tf_tile,xlats_tile,xlons_tile,jdim_tile,idim_tile)
1678 real,
dimension(idim_lalo,jdim_lalo),
intent(in) :: tf_lalo
1679 real,
dimension(jdim_lalo),
intent(in) :: dlats_lalo
1680 real,
dimension(idim_lalo),
intent(in) :: dlons_lalo
1681 real,
dimension(jdim_tile*idim_tile),
intent(in) :: xlats_tile
1682 real,
dimension(jdim_tile*idim_tile),
intent(in) :: xlons_tile
1683 integer,
intent(in) :: jdim_lalo,idim_lalo,jdim_tile,idim_tile
1684 real,
dimension(jdim_tile*idim_tile),
intent(out) :: tf_tile
1687 real,
parameter :: deg2rad=3.1415926/180.0
1688 real,
dimension(jdim_lalo) :: xlats_lalo
1689 real,
dimension(idim_lalo) :: xlons_lalo
1690 real :: tf,wsum,res_km
1691 integer :: itile,jtile
1692 integer :: ii,jj,ij,iii,jjj
1693 integer :: ilalo,jlalo,ilalop1,jlalop1
1694 integer :: istart,iend,jstart,jend,krad
1696 integer,
allocatable,
dimension(:,:) :: id1,id2,jdc
1697 real,
allocatable,
dimension(:,:,:) :: agrid,s2c
1700 print*,
'interpolate from lat/lon grids to any one grid with known lat/lon'
1702 xlats_lalo = dlats_lalo*deg2rad
1703 xlons_lalo = dlons_lalo*deg2rad
1705 allocate(agrid(idim_tile,jdim_tile,2))
1706 agrid(:,:,1) = reshape(xlons_tile, (/idim_tile,jdim_tile/) )
1707 agrid(:,:,2) = reshape(xlats_tile, (/idim_tile,jdim_tile/) )
1708 agrid = agrid*deg2rad
1710 allocate(id1(idim_tile,jdim_tile))
1711 allocate(id2(idim_tile,jdim_tile))
1712 allocate(jdc(idim_tile,jdim_tile))
1713 allocate(s2c(idim_tile,jdim_tile,4))
1721 call
remap_coef( 1, idim_tile, 1, jdim_tile, idim_lalo, jdim_lalo, &
1722 xlons_lalo, xlats_lalo, id1, id2, jdc, s2c, agrid )
1724 do ij = 1, jdim_tile*idim_tile
1726 jtile = (ij-1)/idim_tile + 1
1727 itile = mod(ij,idim_tile)
1728 if (itile==0) itile = idim_tile
1730 ilalo = id1(itile,jtile)
1731 ilalop1 = id2(itile,jtile)
1732 jlalo = jdc(itile,jtile)
1733 jlalop1 = jdc(itile,jtile) + 1
1735 wsum = s2c(itile,jtile,1) + s2c(itile,jtile,2) + &
1736 s2c(itile,jtile,3) + s2c(itile,jtile,4)
1738 tf_tile(ij) = ( s2c(itile,jtile,1)*tf_lalo(ilalo,jlalo) + &
1739 s2c(itile,jtile,2)*tf_lalo(ilalop1,jlalo) + &
1740 s2c(itile,jtile,3)*tf_lalo(ilalop1,jlalop1) + &
1741 s2c(itile,jtile,4)*tf_lalo(ilalo,jlalop1) )/wsum
1744 deallocate(id1, id2, jdc, s2c)
1764 integer,
intent(in) :: iy,im,id,ih
1766 integer,
intent(out) :: mon1,mon2
1767 real,
intent(out) :: wei1,wei2
1771 integer :: mon,monend,monm,monp,jdow,jdoy,jday
1776 real,
dimension(13) :: dayhf
1777 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/
1790 call w3doxdat(jda,jdow,jdoy,jday)
1791 rjday=jdoy+jda(5)/24.
1792 if(rjday.lt.dayhf(1)) rjday=rjday+365.
1798 if( rjday >= dayhf(monm) .and. rjday < dayhf(monp) )
then
1805 print *,
'FATAL ERROR in get_tim_wei, wrong rjday',rjday
1809 wei1 = (dayhf(mon2)-rjday)/(dayhf(mon2)-dayhf(mon1))
1810 wei2 = (rjday-dayhf(mon1))/(dayhf(mon2)-dayhf(mon1))
1812 if( mon2 == 13 ) mon2=1
1814 write(*,
'(a,2i4,3f9.3)')
'mon1,mon2,rjday,wei1,wei2=',mon1,mon2,rjday,wei1,wei2
1833 parameter(a1 = -0.0575)
1834 parameter(a2 = 1.710523e-3)
1835 parameter(a3 = -2.154996e-4)
1837 IF (salinity .LT. 0.)
THEN
1842 tfreez = sal*(a1+a2*sqrt(sal)+a3*sal)
integer function num_parthds()
Return the number of omp threads.
subroutine adjust_nsst(RLA, RLO, SLMSK_TILE, SLMSK_FG_TILE, SKINT_TILE, SICET_TILE, sice_tile, SOILT_TILE, NSST, LENSFC, LSOIL, IDIM, JDIM, ZSEA1, ZSEA2, MON, DAY, DELTSFC, tf_clm_tile, tf_trd_tile, sal_clm_tile)
Read in gsi file with the updated reference temperature increments (on the gaussian grid)...
subroutine, public read_gsi_data(GSI_FILE)
Read file from the GSI containing the foundation temperature increments and mask. ...
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, public read_lat_lon_orog(RLA, RLO, OROG, OROG_UF, TILE_NUM, IDIM, JDIM, IJDIM)
Read latitude and longitude for the cubed-sphere tile from the 'grid' file.
subroutine, public read_data(TSFFCS, SMCFCS, SNOFCS, STCFCS, TG3FCS, ZORFCS, CVFCS, CVBFCS, CVTFCS, ALBFCS, SLIFCS, VEGFCS, CNPFCS, F10M, VETFCS, SOTFCS, ALFFCS, USTAR, FMM, FHH, SIHFCS, SICFCS, SITFCS, TPRCP, SRFLAG, SWDFCS, VMNFCS, VMXFCS, SLCFCS, SLPFCS, ABSFCS, T2M, Q2M, SLMASK, ZSOIL, LSOIL, LENSFC, DO_NSST, NSST)
Read the first guess surface records and nsst records (if selected) for a single cubed-sphere tile...
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 dtzm_point(XT, XZ, DT_COOL, ZC, Z1, Z2, DTZM)
Compute the vertical mean of the NSST t-profile.
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, public get_dim_nc(filename, nlat, nlon)
Get the i/j dimensions of the data from a NetCDF file.
subroutine, public read_salclm_gfs_nc(filename, sal, xlats, xlons, nlat, nlon, itime)
Read the woa05 salinity monthly climatology file.
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).
subroutine remap_coef(is, ie, js, je, im, jm, lon, lat, id1, id2, jdc, s2c, agrid)
Generate the weights and index of the grids used in the bilinear interpolation.
subroutine climo_trend(LATITUDE, MON, DAY, DELTSFC, DTREF)
If the tile point is an isolated water point that has no corresponding gsi water point, then tref is updated using the rtg sst climo trend.
subroutine, public get_tf_clm_dim(file_sst, mlat_sst, mlon_sst)
Get the i/j dimensions of RTG SST climatology file.
subroutine, public write_data(slifcs, tsffcs, snofcs, tg3fcs, zorfcs, albfcs, alffcs, vegfcs, cnpfcs, f10m, t2m, q2m, vetfcs, sotfcs, ustar, fmm, fhh, sicfcs, sihfcs, sitfcs, tprcp, srflag, swdfcs, vmnfcs, vmxfcs, slpfcs, absfcs, slcfcs, smcfcs, stcfcs, idim, jdim, lensfc, lsoil, do_nsst, nsst)
Write out all surface records - and nsst records if selected - on a single cubed-sphere tile to a mod...
program sfc_drv
Stand alone surface/NSST cycle driver for the cubed-sphere grid.
subroutine sfccycle(lugb, len, lsoil, sig1t, deltsfc, iy, im, id, ih, fh, rla, rlo, slmask, orog, orog_uf, use_ufo, nst_anl, sihfcs, sicfcs, sitfcs, swdfcs, slcfcs, vmnfcs, vmxfcs, slpfcs, absfcs, tsffcs, snofcs, zorfcs, albfcs, tg3fcs, cnpfcs, smcfcs, stcfcs, slifcs, aisfcs, vegfcs, vetfcs, sotfcs, alffcs, cvfcs, cvbfcs, cvtfcs, me, nlunit, sz_nml, input_nml_file, ialb, isot, ivegsrc, tile_num_ch, i_index, j_index)
Surface cycling driver routine.
subroutine sfcdrv(LUGB, IDIM, JDIM, LENSFC, LSOIL, DELTSFC, IY, IM, ID, IH, FH, IALB, USE_UFO, DO_NSST, ADJT_NST_ONLY, ZSEA1, ZSEA2, ISOT, IVEGSRC, MYRANK)
Driver routine for updating the surface fields.
real function tfreez(salinity)
Compute the freezing point of water as a function of salinity.
This module contains routines that read and write data.
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.
subroutine, public read_tf_clim_grb(file_sst, sst, rlats_sst, rlons_sst, mlat_sst, mlon_sst, mon)
Read a GRIB1 sst climatological analysis file.