Line data Source code
1 : module dadadj_cam
2 :
3 : ! CAM interfaces for the dry adiabatic adjustment parameterization
4 :
5 : use shr_kind_mod, only: r8=>shr_kind_r8, cs=>shr_kind_cs
6 : use ppgrid, only: pcols, pver, pverp
7 : use constituents, only: pcnst
8 : use air_composition, only: cappav, cpairv
9 : use physconst, only: pi
10 : use physics_types, only: physics_state, physics_ptend, physics_ptend_init
11 : use phys_control, only: use_simple_phys
12 : use cam_abortutils, only: endrun
13 : use cam_logfile, only: iulog
14 : use error_messages, only: handle_errmsg
15 :
16 : use spmd_utils, only: masterproc, masterprocid, mpicom, mpi_integer
17 : use namelist_utils, only: find_group_name
18 : use units, only: getunit, freeunit
19 :
20 : use dadadj, only: dadadj_initial, dadadj_calc
21 :
22 : implicit none
23 : private
24 : save
25 :
26 : public :: &
27 : dadadj_readnl, &
28 : dadadj_init, &
29 : dadadj_tend
30 :
31 : ! Namelist variables
32 : integer :: dadadj_nlvdry = 3 ! number of layers from top of model to apply the adjustment
33 : integer :: dadadj_niter = 15 ! number of iterations for convergence
34 :
35 : !===============================================================================
36 : contains
37 : !===============================================================================
38 :
39 1536 : subroutine dadadj_readnl(filein)
40 :
41 : character(len=cs), intent(in) :: filein ! Input namelist filename
42 :
43 : namelist /dadadj_nl/ dadadj_nlvdry, dadadj_niter
44 :
45 : integer :: unitn, ierr
46 : character(len=*), parameter :: sub='dadadj_readnl'
47 : !------------------------------------------------------------------
48 :
49 : ! Read namelist
50 1536 : if (masterproc) then
51 2 : unitn = getunit()
52 2 : open(unitn, file=trim(filein), status='old')
53 2 : call find_group_name(unitn, 'dadadj_nl', status=ierr)
54 2 : if (ierr == 0) then
55 0 : read(unitn, dadadj_nl, iostat=ierr)
56 0 : if (ierr /= 0) then
57 0 : call endrun( sub//':: ERROR reading namelist')
58 : end if
59 : end if
60 2 : close(unitn)
61 2 : call freeunit(unitn)
62 : end if
63 :
64 : #ifdef SPMD
65 : ! Broadcast namelist variables
66 1536 : call mpibcast(dadadj_nlvdry, 1, mpi_integer, masterprocid, mpicom)
67 1536 : call mpibcast(dadadj_niter, 1, mpi_integer, masterprocid, mpicom)
68 : #endif
69 :
70 1536 : call dadadj_initial(dadadj_nlvdry, dadadj_niter)
71 :
72 1536 : if (masterproc .and. .not. use_simple_phys) then
73 2 : write(iulog,*)'Dry adiabatic adjustment applied to top N layers; N=', &
74 4 : dadadj_nlvdry
75 2 : write(iulog,*)'Dry adiabatic adjustment number of iterations for convergence =', &
76 4 : dadadj_niter
77 : end if
78 :
79 1536 : end subroutine dadadj_readnl
80 :
81 :
82 : !===============================================================================
83 :
84 1536 : subroutine dadadj_init()
85 : use cam_history, only: addfld
86 :
87 3072 : call addfld('DADADJ_PD', (/ 'lev' /), 'A', 'probability', 'dry adiabatic adjustment probability')
88 :
89 1536 : end subroutine dadadj_init
90 :
91 :
92 : !===============================================================================
93 :
94 2730672 : subroutine dadadj_tend(dt, state, ptend)
95 1536 : use cam_history, only: outfld
96 :
97 : real(r8), intent(in) :: dt ! Time step [s]
98 : type(physics_state), intent(in) :: state ! Physics state variables
99 : type(physics_ptend), intent(out) :: ptend ! parameterization tendencies
100 :
101 : logical :: lq(pcnst)
102 : real(r8) :: dadpdf(pcols, pver)
103 : integer :: ncol, lchnk, icol_err
104 : character(len=128) :: errstring ! Error string
105 :
106 65016 : ncol = state%ncol
107 65016 : lchnk = state%lchnk
108 65016 : lq(:) = .FALSE.
109 65016 : lq(1) = .TRUE.
110 65016 : call physics_ptend_init(ptend, state%psetcols, 'dadadj', ls=.true., lq=lq)
111 :
112 : ! use the ptend components for temporary storate and copy state info for input to
113 : ! dadadj_calc which directly updates the temperature and moisture input arrays.
114 :
115 101027304 : ptend%s(:ncol,:pver) = state%t(:ncol,:pver)
116 101027304 : ptend%q(:ncol,:pver,1) = state%q(:ncol,:pver,1)
117 :
118 : call dadadj_calc( &
119 : ncol, state%pmid, state%pint, state%pdel, cappav(:,:,lchnk), ptend%s, &
120 65016 : ptend%q(:,:,1), dadpdf, icol_err)
121 :
122 23735376 : call outfld('DADADJ_PD', dadpdf(:ncol,:), ncol, lchnk)
123 :
124 65016 : if (icol_err > 0) then
125 : ! error exit
126 : write(errstring, *) &
127 0 : 'dadadj_calc: No convergence in column at lat,lon:', &
128 0 : state%lat(icol_err)*180._r8/pi, state%lon(icol_err)*180._r8/pi
129 0 : call handle_errmsg(errstring, subname="dadadj_tend")
130 : end if
131 :
132 101027304 : ptend%s(:ncol,:) = (ptend%s(:ncol,:) - state%t(:ncol,:) )/dt * cpairv(:ncol,:,lchnk)
133 101027304 : ptend%q(:ncol,:,1) = (ptend%q(:ncol,:,1) - state%q(:ncol,:,1))/dt
134 :
135 65016 : end subroutine dadadj_tend
136 :
137 : !===============================================================================
138 : end module dadadj_cam
|