-
Notifications
You must be signed in to change notification settings - Fork 1
/
local_physics_types.F90
331 lines (255 loc) · 12.5 KB
/
local_physics_types.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
!-------------------------------------------------------------------------------
!physics data types module
!-------------------------------------------------------------------------------
module physics_types
use shr_kind_mod, only: r8 => shr_kind_r8
implicit none
private ! Make default type private to the module
logical, parameter :: adjust_te = .FALSE.
integer, parameter :: pcnst_max = 2000
! Public types:
public physics_state
public physics_tend
public physics_ptend
public physics_ptend_init
!-------------------------------------------------------------------------------
type physics_state
integer :: &
lchnk, &! chunk index
ngrdcol, &! -- Grid -- number of active columns (on the grid)
psetcols=0, &! -- -- max number of columns set - if subcols = pcols*psubcols, else = pcols
ncol=0 ! -- -- sum of nsubcol for all ngrdcols - number of active columns
real(r8), dimension(:), allocatable :: &
lat, &! latitude (radians)
lon, &! longitude (radians)
ps, &! surface pressure
psdry, &! dry surface pressure
phis, &! surface geopotential
ulat, &! unique latitudes (radians)
ulon ! unique longitudes (radians)
real(r8), dimension(:,:),allocatable :: &
t, &! temperature (K)
u, &! zonal wind (m/s)
v, &! meridional wind (m/s)
s, &! dry static energy
omega, &! vertical pressure velocity (Pa/s)
pmid, &! midpoint pressure (Pa)
pmiddry, &! midpoint pressure dry (Pa)
pdel, &! layer thickness (Pa)
pdeldry, &! layer thickness dry (Pa)
rpdel, &! reciprocal of layer thickness (Pa)
rpdeldry,&! recipricol layer thickness dry (Pa)
lnpmid, &! ln(pmid)
lnpmiddry,&! log midpoint pressure dry (Pa)
exner, &! inverse exner function w.r.t. surface pressure (ps/p)^(R/cp)
zm ! geopotential height above surface at midpoints (m)
real(r8), dimension(:,:,:),allocatable :: &
q ! constituent mixing ratio (kg/kg moist or dry air depending on type)
real(r8), dimension(:,:),allocatable :: &
pint, &! interface pressure (Pa)
pintdry, &! interface pressure dry (Pa)
lnpint, &! ln(pint)
lnpintdry,&! log interface pressure dry (Pa)
zi ! geopotential height above surface at interfaces (m)
real(r8), dimension(:),allocatable :: &
te_ini, &! vertically integrated total (kinetic + static) energy of initial state
te_cur, &! vertically integrated total (kinetic + static) energy of current state
tw_ini, &! vertically integrated total water of initial state
tw_cur ! vertically integrated total water of new state
integer :: count ! count of values with significant energy or water imbalances
integer, dimension(:),allocatable :: &
latmapback, &! map from column to unique lat for that column
lonmapback, &! map from column to unique lon for that column
cid ! unique column id
integer :: ulatcnt, &! number of unique lats in chunk
uloncnt ! number of unique lons in chunk
end type physics_state
!-------------------------------------------------------------------------------
type physics_tend
integer :: psetcols=0 ! max number of columns set- if subcols = pcols*psubcols, else = pcols
real(r8), dimension(:,:),allocatable :: dtdt, dudt, dvdt
real(r8), dimension(:), allocatable :: flx_net
real(r8), dimension(:), allocatable :: &
te_tnd, &! cumulative boundary flux of total energy
tw_tnd ! cumulative boundary flux of total water
end type physics_tend
!-------------------------------------------------------------------------------
! This is for tendencies returned from individual parameterizations
type physics_ptend
integer :: psetcols=0 ! max number of columns set- if subcols = pcols*psubcols, else = pcols
character*24 :: name ! name of parameterization which produced tendencies.
logical :: &
ls = .false., &! true if dsdt is returned
lu = .false., &! true if dudt is returned
lv = .false. ! true if dvdt is returned
logical, dimension(pcnst_max) :: lq = .false. ! true if dqdt() is returned
!logical,dimension(1) :: lq = .false. ! true if dqdt() is returned
integer :: &
top_level, &! top level index for which nonzero tendencies have been set
bot_level ! bottom level index for which nonzero tendencies have been set
real(r8), dimension(:,:),allocatable :: &
s, &! heating rate (J/kg/s)
u, &! u momentum tendency (m/s/s)
v ! v momentum tendency (m/s/s)
real(r8), dimension(:,:,:),allocatable :: &
q ! consituent tendencies (kg/kg/s)
! boundary fluxes
real(r8), dimension(:),allocatable ::&
hflux_srf, &! net heat flux at surface (W/m2)
hflux_top, &! net heat flux at top of model (W/m2)
taux_srf, &! net zonal stress at surface (Pa)
taux_top, &! net zonal stress at top of model (Pa)
tauy_srf, &! net meridional stress at surface (Pa)
tauy_top ! net meridional stress at top of model (Pa)
real(r8), dimension(:,:),allocatable ::&
cflx_srf, &! constituent flux at surface (kg/m2/s)
cflx_top ! constituent flux top of model (kg/m2/s)
end type physics_ptend
contains
!===============================================================================
subroutine physics_ptend_init(ptend, psetcols, pver, pcnst, name, ls, lu, lv, lq)
!-----------------------------------------------------------------------
! Allocate the fields in the structure which are specified.
! Initialize the parameterization tendency structure to "empty"
!-----------------------------------------------------------------------
!------------------------------Arguments--------------------------------
type(physics_ptend), intent(out) :: ptend ! Parameterization tendencies
integer, intent(in) :: psetcols ! maximum number of columns
integer, intent(in) :: pver
integer, intent(in) :: pcnst
character(len=*) :: name ! optional name of parameterization which produced tendencies.
logical, optional :: ls ! if true, then fields to support dsdt are allocated
logical, optional :: lu ! if true, then fields to support dudt are allocated
logical, optional :: lv ! if true, then fields to support dvdt are allocated
logical, dimension(pcnst),optional :: lq ! if true, then fields to support dqdt are allocated
!-----------------------------------------------------------------------
if (allocated(ptend%s)) then
call endrun(' physics_ptend_init: ptend should not be allocated before calling this routine')
end if
ptend%name = name
ptend%psetcols = psetcols
! If no fields being stored, initialize all values to appropriate nulls and return
if (.not. present(ls) .and. .not. present(lu) .and. .not. present(lv) .and. .not. present(lq) ) then
ptend%ls = .false.
ptend%lu = .false.
ptend%lv = .false.
ptend%lq(:) = .false.
ptend%top_level = 1
ptend%bot_level = pver
return
end if
if (present(ls)) then
ptend%ls = ls
else
ptend%ls = .false.
end if
if (present(lu)) then
ptend%lu = lu
else
ptend%lu = .false.
end if
if (present(lv)) then
ptend%lv = lv
else
ptend%lv = .false.
end if
if (present(lq)) then
ptend%lq(1:pcnst) = lq(1:pcnst)
else
ptend%lq(:) = .false.
end if
call physics_ptend_alloc(ptend, psetcols , pver, pcnst )
call physics_ptend_reset(ptend, pver )
return
end subroutine physics_ptend_init
!===============================================================================
!===============================================================================
subroutine physics_ptend_reset(ptend , pver )
!-----------------------------------------------------------------------
! Reset the parameterization tendency structure to "empty"
!-----------------------------------------------------------------------
!------------------------------Arguments--------------------------------
type(physics_ptend), intent(inout) :: ptend ! Parameterization tendencies
integer, intent(in) :: pver
!-----------------------------------------------------------------------
integer :: m ! Index for constiuent
!-----------------------------------------------------------------------
if(ptend%ls) then
ptend%s = 0._r8
ptend%hflux_srf = 0._r8
ptend%hflux_top = 0._r8
endif
if(ptend%lu) then
ptend%u = 0._r8
ptend%taux_srf = 0._r8
ptend%taux_top = 0._r8
endif
if(ptend%lv) then
ptend%v = 0._r8
ptend%tauy_srf = 0._r8
ptend%tauy_top = 0._r8
endif
if(any (ptend%lq(:))) then
ptend%q = 0._r8
ptend%cflx_srf = 0._r8
ptend%cflx_top = 0._r8
end if
ptend%top_level = 1
ptend%bot_level = pver
return
end subroutine physics_ptend_reset
!===============================================================================
!===============================================================================
subroutine physics_ptend_alloc(ptend,psetcols,pver,pcnst)
! allocate the individual ptend components
type(physics_ptend), intent(inout) :: ptend
integer, intent(in) :: psetcols, pver, pcnst
integer :: ierr = 0
ptend%psetcols = psetcols
if (ptend%ls) then
allocate(ptend%s(psetcols,pver), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%s')
allocate(ptend%hflux_srf(psetcols), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%hflux_srf')
allocate(ptend%hflux_top(psetcols), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%hflux_top')
end if
if (ptend%lu) then
allocate(ptend%u(psetcols,pver), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%u')
allocate(ptend%taux_srf(psetcols), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%taux_srf')
allocate(ptend%taux_top(psetcols), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%taux_top')
end if
if (ptend%lv) then
allocate(ptend%v(psetcols,pver), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%v')
allocate(ptend%tauy_srf(psetcols), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%tauy_srf')
allocate(ptend%tauy_top(psetcols), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%tauy_top')
end if
if (any(ptend%lq)) then
allocate(ptend%q(psetcols,pver,pcnst), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%q')
allocate(ptend%cflx_srf(psetcols,pcnst), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%cflx_srf')
allocate(ptend%cflx_top(psetcols,pcnst), stat=ierr)
if ( ierr /= 0 ) call endrun('physics_ptend_alloc error: allocation error for ptend%cflx_top')
end if
end subroutine physics_ptend_alloc
!===============================================================================
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
subroutine endrun(msg)
integer :: iulog
character(len=*), intent(in), optional :: msg ! string to be printed
iulog=6
if (present (msg)) then
write(iulog,*)'ENDRUN:', msg
else
write(iulog,*)'ENDRUN: called without a message string'
end if
stop
end subroutine endrun
end module physics_types