Line data Source code
1 : 2 : module iop_forcing 3 : 4 : ! ----------------------------------------------- ! 5 : ! User-defined manipulation of forcings for SCAM. ! 6 : ! ! 7 : ! Author : Sungsu Park and John Truesdale ! 8 : ! ! 9 : ! ----------------------------------------------- ! 10 : 11 : use shr_kind_mod, only: r8 => shr_kind_r8 12 : implicit none 13 : private 14 : 15 : public scam_use_iop_srf 16 : public scam_set_iop_Tg 17 : 18 : contains 19 : 20 0 : subroutine scam_use_iop_srf( cam_in ) 21 : ! ------------------------------------------------------- ! 22 : ! Use the SCAM-IOP specified surface LHFLX/SHFLX/ustar/Tg ! 23 : ! instead of using internally-computed values. ! 24 : ! ! 25 : ! Author : John Truesdale and Sungsu Park ! 26 : ! ------------------------------------------------------- ! 27 : use ppgrid, only: begchunk, endchunk 28 : use camsrfexch, only: cam_in_t 29 : use physconst, only: stebol, latvap 30 : use scamMod 31 : use cam_abortutils, only: endrun 32 : 33 : implicit none 34 : save 35 : 36 : type(cam_in_t), intent(INOUT) :: cam_in(begchunk:endchunk) 37 : integer :: c ! Chunk index 38 : integer :: ncol ! Number of columns 39 : 40 0 : if( scm_iop_lhflxshflxTg .and. scm_iop_Tg ) then 41 0 : call endrun( 'scam_use_iop_srf : scm_iop_lhflxshflxTg and scm_iop_Tg must not be specified at the same time.') 42 : end if 43 : 44 0 : if( scm_iop_lhflxshflxTg ) then 45 0 : do c = begchunk, endchunk 46 0 : ncol = cam_in(c)%ncol 47 0 : if( have_lhflx ) then 48 0 : cam_in(c)%lhf(1) = lhflxobs(1) 49 0 : cam_in(c)%cflx(1,1) = lhflxobs(1) / latvap 50 : endif 51 0 : if( have_shflx ) cam_in(c)%shf(1) = shflxobs(1) 52 0 : if( have_Tg ) then 53 0 : cam_in(c)%ts(1) = tground(1) 54 0 : cam_in(c)%lwup(1) = stebol * tground(1)**4 55 : endif 56 : end do 57 : endif 58 : 59 0 : if( scm_iop_Tg .or. scm_crm_mode) then 60 0 : do c = begchunk, endchunk 61 0 : ncol = cam_in(c)%ncol 62 0 : if( have_Tg ) then 63 0 : cam_in(c)%ts(1) = tground(1) 64 0 : cam_in(c)%lwup(1) = stebol * tground(1)**4 65 : endif 66 : end do 67 : endif 68 : 69 0 : end subroutine scam_use_iop_srf 70 : 71 0 : subroutine scam_set_iop_Tg( cam_out ) 72 : ! ----------------------------- ! 73 : ! USE the SCAM-IOP specified Tg ! 74 : ! ----------------------------- ! 75 0 : use ppgrid, only: begchunk, endchunk 76 : use camsrfexch, only: cam_out_t 77 : use scamMod 78 : use cam_abortutils, only: endrun 79 : 80 : implicit none 81 : save 82 : 83 : type(cam_out_t), intent(INOUT) :: cam_out(begchunk:endchunk) 84 : integer :: c ! Chunk index 85 : 86 0 : if( scm_iop_lhflxshflxTg .and. scm_iop_Tg ) then 87 0 : call endrun( 'scam_use_iop_srf : scm_iop_lhflxshflxTg and scm_iop_Tg must not be specified at the same time.') 88 : end if 89 0 : if( scm_iop_Tg .or. scm_crm_mode ) then 90 0 : do c = begchunk, endchunk 91 0 : cam_out(c)%tbot(1) = tground(1) 92 : end do 93 : endif 94 : 95 0 : end subroutine scam_set_iop_Tg 96 : 97 : end module iop_forcing