Line data Source code
1 : module beljaars_drag_cam
2 :
3 : use shr_kind_mod, only: r8 => shr_kind_r8
4 : use spmd_utils, only: masterproc
5 : use cam_abortutils, only: endrun
6 : use shr_log_mod, only: errMsg => shr_log_errMsg
7 : use cam_logfile, only: iulog
8 : use ppgrid, only: pcols, pver
9 :
10 : implicit none
11 : private
12 :
13 : public :: beljaars_drag_readnl
14 : public :: beljaars_drag_register
15 : public :: beljaars_drag_init
16 : public :: beljaars_drag_tend
17 :
18 : ! Is this module on at all?
19 : logical, public, protected :: do_beljaars = .false.
20 :
21 : ! Tuning parameters for TMS.
22 : real(r8) :: blj_orocnst
23 : real(r8) :: blj_z0fac
24 :
25 : ! pbuf field indices
26 : integer :: &
27 : sgh30_idx = -1, &
28 : dragblj_idx = -1, &
29 : taubljx_idx = -1, &
30 : taubljy_idx = -1
31 :
32 : contains
33 :
34 1490712 : subroutine beljaars_drag_readnl(nlfile)
35 : use namelist_utils, only: find_group_name
36 : use units, only: getunit, freeunit
37 : use spmd_utils, only: masterprocid, mpi_logical, mpi_real8, mpicom
38 :
39 : ! filepath for file containing namelist input
40 : character(len=*), intent(in) :: nlfile
41 :
42 : ! file unit and error code
43 : integer :: unitn, ierr
44 :
45 : character(len=*), parameter :: subname = "beljaars_drag_readnl"
46 :
47 : namelist /blj_nl/ do_beljaars
48 :
49 1536 : ierr = 0
50 :
51 1536 : if (masterproc) then
52 2 : unitn = getunit()
53 2 : open( unitn, file=trim(nlfile), status='old' )
54 2 : call find_group_name(unitn, 'blj_nl', status=ierr)
55 2 : if (ierr == 0) then
56 2 : read(unitn, blj_nl, iostat=ierr)
57 2 : if (ierr /= 0) then
58 0 : call endrun(subname // ':: ERROR reading namelist')
59 : end if
60 : end if
61 2 : close(unitn)
62 2 : call freeunit(unitn)
63 : end if
64 :
65 1536 : call mpi_bcast(do_beljaars, 1, mpi_logical, masterprocid, mpicom, ierr)
66 1536 : if (ierr /= 0) call endrun(errMsg(__FILE__, __LINE__)//" mpi_bcast error")
67 :
68 1536 : end subroutine beljaars_drag_readnl
69 :
70 1536 : subroutine beljaars_drag_register()
71 : use physics_buffer, only: pbuf_add_field, dtype_r8
72 :
73 1536 : call pbuf_add_field("dragblj", "physpkg", dtype_r8, (/pcols,pver/), dragblj_idx)
74 1536 : call pbuf_add_field("taubljx", "physpkg", dtype_r8, (/pcols/), taubljx_idx)
75 1536 : call pbuf_add_field("taubljy", "physpkg", dtype_r8, (/pcols/), taubljy_idx)
76 :
77 1536 : end subroutine beljaars_drag_register
78 :
79 1536 : subroutine beljaars_drag_init()
80 :
81 1536 : use cam_history, only: addfld, add_default, horiz_only
82 : use error_messages, only: handle_errmsg
83 : use phys_control, only: phys_getopts
84 : use physconst, only: karman, gravit, rair
85 : use physics_buffer, only: pbuf_get_index
86 : use beljaars_drag, only: init_blj
87 :
88 : logical :: history_amwg
89 :
90 : character(len=128) :: errstring
91 :
92 1536 : if (.not. do_beljaars) return
93 :
94 1536 : call phys_getopts(history_amwg_out=history_amwg)
95 :
96 1536 : call init_blj( r8, gravit, rair, errstring )
97 1536 : call handle_errmsg(errstring, subname="init_blj")
98 :
99 3072 : call addfld('DRAGBLJ', (/ 'lev' /) , 'A', '1/s', 'Drag profile from Beljaars SGO ')
100 1536 : call addfld('TAUBLJX', horiz_only, 'A', 'N/m2', 'Zonal integrated drag from Beljaars SGO')
101 1536 : call addfld('TAUBLJY', horiz_only, 'A', 'N/m2', 'Meridional integrated drag from Beljaars SGO')
102 1536 : if (history_amwg) then
103 1536 : call add_default( 'TAUBLJX ', 1, ' ' )
104 1536 : call add_default( 'TAUBLJY ', 1, ' ' )
105 : end if
106 :
107 1536 : if (masterproc) then
108 2 : write(iulog,*)'Using Beljaars SGO drag module'
109 : end if
110 :
111 1536 : sgh30_idx = pbuf_get_index("SGH30")
112 :
113 1536 : end subroutine beljaars_drag_init
114 :
115 1489176 : subroutine beljaars_drag_tend(state, pbuf, cam_in)
116 1536 : use physics_buffer, only: physics_buffer_desc, pbuf_get_field
117 : use physics_types, only: physics_state
118 : use camsrfexch, only: cam_in_t
119 : use cam_history, only: outfld
120 : use beljaars_drag, only: compute_blj
121 :
122 : type(physics_state), intent(in) :: state
123 : type(physics_buffer_desc), pointer, intent(in) :: pbuf(:)
124 : type(cam_in_t), intent(in) :: cam_in
125 :
126 1489176 : real(r8), pointer :: sgh30(:)
127 1489176 : real(r8), pointer :: dragblj(:,:)
128 1489176 : real(r8), pointer :: taubljx(:), taubljy(:)
129 :
130 1489176 : call pbuf_get_field(pbuf, dragblj_idx, dragblj)
131 1489176 : call pbuf_get_field(pbuf, taubljx_idx, taubljx)
132 1489176 : call pbuf_get_field(pbuf, taubljy_idx, taubljy)
133 :
134 1489176 : if (.not. do_beljaars) then
135 0 : dragblj = 0._r8
136 0 : taubljx = 0._r8
137 0 : taubljy = 0._r8
138 : return
139 : end if
140 :
141 1489176 : call pbuf_get_field(pbuf, sgh30_idx, sgh30)
142 :
143 : call compute_blj( pcols , pver , state%ncol , &
144 : state%u , state%v , state%t , state%pmid , &
145 : state%pdel , state%zm , sgh30 , dragblj , &
146 1489176 : taubljx , taubljy , cam_in%landfrac )
147 :
148 1489176 : call outfld("TAUBLJX", taubljx, pcols, state%lchnk)
149 1489176 : call outfld("TAUBLJY", taubljy, pcols, state%lchnk)
150 1489176 : call outfld("DRAGBLJ", dragblj, pcols, state%lchnk)
151 :
152 2978352 : end subroutine beljaars_drag_tend
153 :
154 : end module beljaars_drag_cam
|