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, cm=>shr_kind_cm
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_init, dadadj_run
21 :
22 : implicit none
23 : private
24 : save
25 :
26 : public :: &
27 : dadadj_readnl, &
28 : dadadj_cam_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 : integer :: errflg ! CCPP physics scheme error flag
47 : character(len=512) :: errmsg ! CCPP physics scheme error message
48 : character(len=*), parameter :: sub='dadadj_readnl'
49 : !------------------------------------------------------------------
50 :
51 : ! Read namelist
52 1536 : if (masterproc) then
53 2 : unitn = getunit()
54 2 : open(unitn, file=trim(filein), status='old')
55 2 : call find_group_name(unitn, 'dadadj_nl', status=ierr)
56 2 : if (ierr == 0) then
57 0 : read(unitn, dadadj_nl, iostat=ierr)
58 0 : if (ierr /= 0) then
59 0 : call endrun( sub//':: ERROR reading namelist')
60 : end if
61 : end if
62 2 : close(unitn)
63 2 : call freeunit(unitn)
64 : end if
65 :
66 : #ifdef SPMD
67 : ! Broadcast namelist variables
68 1536 : call mpibcast(dadadj_nlvdry, 1, mpi_integer, masterprocid, mpicom)
69 1536 : call mpibcast(dadadj_niter, 1, mpi_integer, masterprocid, mpicom)
70 : #endif
71 :
72 1536 : call dadadj_init(dadadj_nlvdry, dadadj_niter, pver, errmsg, errflg)
73 1536 : if (errflg /=0) then
74 0 : call endrun('dadadj_readnl: Error returned from dadadj_init: '//trim(errmsg))
75 : end if
76 :
77 1536 : if (masterproc .and. .not. use_simple_phys) then
78 2 : write(iulog,*)'Dry adiabatic adjustment applied to top N layers; N=', &
79 4 : dadadj_nlvdry
80 2 : write(iulog,*)'Dry adiabatic adjustment number of iterations for convergence =', &
81 4 : dadadj_niter
82 : end if
83 :
84 1536 : end subroutine dadadj_readnl
85 :
86 :
87 : !===============================================================================
88 :
89 1536 : subroutine dadadj_cam_init()
90 : use cam_history, only: addfld
91 :
92 3072 : call addfld('DADADJ_PD', (/ 'lev' /), 'A', 'probability', 'dry adiabatic adjustment probability')
93 :
94 1536 : end subroutine dadadj_cam_init
95 :
96 :
97 : !===============================================================================
98 :
99 5981472 : subroutine dadadj_tend(dt, state, ptend)
100 1536 : use cam_history, only: outfld
101 :
102 : real(r8), intent(in) :: dt ! Time step [s]
103 : type(physics_state), intent(in) :: state ! Physics state variables
104 : type(physics_ptend), intent(out) :: ptend ! parameterization tendencies
105 :
106 : character(len=512) :: errstring ! Error string
107 : character(len=512) :: errmsg ! CCPP physics scheme error message
108 : character(len=64) :: scheme_name! CCPP physics scheme name (not used in CAM)
109 : integer :: icol_err
110 : integer :: lchnk
111 : integer :: ncol
112 : integer :: errflg ! CCPP physics scheme error flag
113 : logical :: lq(pcnst)
114 : real(r8) :: dadpdf(pcols, pver)
115 :
116 : !------------------------------------------------------------------
117 1495368 : ncol = state%ncol
118 1495368 : lchnk = state%lchnk
119 1495368 : lq(:) = .FALSE.
120 1495368 : lq(1) = .TRUE.
121 1495368 : call physics_ptend_init(ptend, state%psetcols, 'dadadj', ls=.true., lq=lq)
122 :
123 : !REMOVECAM - no longer need these when CAM is retired and pcols no longer exists
124 1495368 : dadpdf = 0._r8
125 662448024 : ptend%s = 0._r8
126 1988839440 : ptend%q = 0._r8
127 : !REMOVECAM_END
128 :
129 : ! dadadj_run returns t tend, we are passing the ptend%s array to receive the t tendency and will convert it to s
130 : ! before it is returned to CAM..
131 : call dadadj_run( &
132 0 : ncol, pver, dt, state%pmid(:ncol,:), state%pint(:ncol,:), state%pdel(:ncol,:), &
133 0 : state%t(:ncol,:), state%q(:ncol,:,1), cappav(:ncol,:,lchnk), cpairv(:ncol,:,lchnk), ptend%s(:ncol,:), &
134 1495368 : ptend%q(:ncol,:,1), dadpdf(:ncol,:), scheme_name, errmsg, errflg)
135 :
136 : ! error exit
137 1495368 : if (errflg /= 0) then
138 : ! If this is a Convergence error then output lat lon of problem column using column index (errflg)
139 0 : if(index('Convergence', errmsg) /= 0)then
140 0 : write(errstring, *) trim(adjustl(errmsg)),' lat:',state%lat(errflg)*180._r8/pi,' lon:', &
141 0 : state%lon(errflg)*180._r8/pi
142 : else
143 0 : errstring=trim(errmsg)
144 : end if
145 0 : call endrun('Error dadadj_tend:'//trim(errstring))
146 : end if
147 :
148 153698328 : call outfld('DADADJ_PD', dadpdf(:ncol,:), ncol, lchnk)
149 :
150 1495368 : end subroutine dadadj_tend
151 :
152 : !===============================================================================
153 : end module dadadj_cam
|