-
Notifications
You must be signed in to change notification settings - Fork 10
/
fv_regridding_utils.F90
341 lines (297 loc) · 11.2 KB
/
fv_regridding_utils.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
332
333
334
335
336
337
338
339
340
#include "MAPL_Generic.h"
module fv_regridding_utils
use ESMF
use fv_arrays_mod, only: fv_atmos_type, fv_grid_type, fv_grid_bounds_type, FVPRC, REAL4, REAL8
use fv_diagnostics_mod,only: prt_maxmin
use fv_mp_mod, only: is_master, ng
use fv_mapz_mod, only: mappm
use mpp_mod, only: mpp_error, FATAL, NOTE, mpp_broadcast,mpp_npes
use MAPL
implicit none
private
public remap_scalar
public fv_rst
public copy_fv_rst
real(FVPRC), parameter :: PI = MAPL_PI_R8
real(FVPRC), parameter :: OMEGA = MAPL_OMEGA
real(FVPRC), parameter :: GRAV = MAPL_GRAV
real(FVPRC), parameter :: KAPPA = MAPL_KAPPA
real(FVPRC), parameter :: RDGAS = MAPL_RGAS
real(FVPRC), parameter :: RVGAS = MAPL_RVAP
real(FVPRC), parameter :: CP_AIR = MAPL_CP
real(FVPRC), parameter:: zvir = rvgas/rdgas - 1.
type fv_var
character(len=128) :: name
integer :: nlev
integer :: n_ungrid
integer :: rank
real(FVPRC), allocatable :: ptr2d(:,:)
real(FVPRC), allocatable :: ptr3d(:,:,:)
real(FVPRC), allocatable :: ptr4d(:,:,:,:)
contains
procedure :: alloc_var
procedure :: dealloc_var
end type fv_var
type fv_rst
integer :: ungrid_size
logical :: has_edge
logical :: has_center
character(len=1024) :: file_name
logical :: have_descriptor
type(fv_var), allocatable :: vars(:)
end type fv_rst
contains
subroutine dealloc_var(this)
class(fv_var), intent(inout) :: this
if (this%rank==2) then
deallocate(this%ptr2d)
else if (this%rank==3) then
deallocate(this%ptr3d)
else if (this%rank==4) then
deallocate(this%ptr4d)
end if
end subroutine dealloc_var
subroutine alloc_var(this,is,ie,js,je,km,rc)
class(fv_var), intent(inout) :: this
integer, intent(in) :: is,ie,js,je
integer, intent(in), optional :: km
integer, intent(out), optional :: rc
integer :: status
integer :: km_use
if (this%rank==2) then
allocate(this%ptr2d(is:ie,js:je),source=0.0_FVPRC)
else if (this%rank==3) then
if (this%n_ungrid > 0) then
allocate(this%ptr3d(is:ie,js:je,this%n_ungrid),source=0.0_FVPRC)
else if (this%n_ungrid == 0) then
if (present(km)) then
km_use = km
else
km_use = this%nlev
end if
allocate(this%ptr3d(is:ie,js:je,km_use),source=0.0_FVPRC)
end if
else if (this%rank == 4) then
if (present(km)) then
km_use = km
else
km_use = this%nlev
end if
allocate(this%ptr4d(is:ie,js:je,km_use,this%n_ungrid),source=0.0_FVPRC)
end if
_RETURN(_SUCCESS)
end subroutine alloc_var
subroutine copy_fv_rst(in_rst,out_rst)
type(fv_rst), pointer, intent(inout) :: in_rst(:)
type(fv_rst), pointer, intent(inout) :: out_rst(:)
integer :: ifile,ivar
allocate(out_rst(size(in_rst)) )
do ifile=1,size(in_rst)
allocate( out_rst(ifile)%vars(size(in_rst(ifile)%vars) ) )
out_rst(ifile)%file_name=in_rst(ifile)%file_name
out_rst(ifile)%have_descriptor=in_rst(ifile)%have_descriptor
out_rst(ifile)%ungrid_size=in_rst(ifile)%ungrid_size
out_rst(ifile)%has_center=in_rst(ifile)%has_center
out_rst(ifile)%has_edge=in_rst(ifile)%has_edge
do ivar=1,size(in_rst(ifile)%vars)
out_rst(ifile)%vars(ivar)%name=in_rst(ifile)%vars(ivar)%name
out_rst(ifile)%vars(ivar)%nlev=in_rst(ifile)%vars(ivar)%nlev
out_rst(ifile)%vars(ivar)%n_ungrid=in_rst(ifile)%vars(ivar)%n_ungrid
out_rst(ifile)%vars(ivar)%rank=in_rst(ifile)%vars(ivar)%rank
enddo
enddo
end subroutine copy_fv_rst
subroutine remap_scalar(im, jm, km, npz, nq, ncnst, ak0, bk0, psc, gzc, ta, qa, Atm, in_fv_rst,out_fv_rst)
type(fv_atmos_type), intent(inout) :: Atm
integer, intent(in):: im, jm, km, npz, nq, ncnst
real(FVPRC), intent(in):: ak0(km+1), bk0(km+1)
real(FVPRC), intent(in), dimension(Atm%bd%is:Atm%bd%ie,Atm%bd%js:Atm%bd%je):: psc, gzc
real(FVPRC), intent(in), dimension(Atm%bd%is:Atm%bd%ie,Atm%bd%js:Atm%bd%je,km):: ta
real(FVPRC), intent(in), dimension(Atm%bd%is:Atm%bd%ie,Atm%bd%js:Atm%bd%je,km,ncnst):: qa
type(fv_rst), pointer, intent(inout) :: in_fv_rst(:)
type(fv_rst), pointer, intent(inout) :: out_fv_rst(:)
! local:
real(FVPRC), dimension(Atm%bd%is:Atm%bd%ie,km):: tp
real(FVPRC), dimension(Atm%bd%is:Atm%bd%ie,km+1):: pe0, pn0
real(FVPRC), dimension(Atm%bd%is:Atm%bd%ie,npz):: qn1
real(FVPRC), dimension(Atm%bd%is:Atm%bd%ie,npz+1):: pe1, pn1
real(FVPRC), dimension(Atm%bd%is:Atm%bd%ie,npz+1):: q_edge_old,q_edge_new
real(FVPRC) pt0(km), gz(km+1), pk0(km+1)
real(FVPRC) qp( Atm%bd%is:Atm%bd%ie,km,ncnst)
real(FVPRC) qp1(Atm%bd%is:Atm%bd%ie,km)
real(FVPRC) pst, p1, p2, alpha, rdg
integer i,j,k, iq
integer sphum,ifile,ivar,n_ungrid
integer :: is, ie, js, je
integer :: isd, ied, jsd, jed
logical :: doVert
is = Atm%bd%is
ie = Atm%bd%ie
js = Atm%bd%js
je = Atm%bd%je
isd = Atm%bd%isd
ied = Atm%bd%ied
jsd = Atm%bd%jsd
jed = Atm%bd%jed
sphum = 1
if ( sphum/=1 ) then
call mpp_error(FATAL,'SPHUM must be 1st tracer')
endif
do j=js,je
do i=is,ie
do iq=1,ncnst
do k=1,km
qp(i,k,iq) = qa(i,j,k,iq)
enddo
enddo
do k=1,km
tp(i,k) = ta(i,j,k)*(1.+zvir*qp(i,k,sphum))
enddo
! Tracers:
do k=1,km+1
pe0(i,k) = ak0(k) + bk0(k)*psc(i,j)
pn0(i,k) = log(pe0(i,k))
pk0(k) = pe0(i,k)**kappa
enddo
! * Adjust interpolated ps to model terrain
gz(km+1) = gzc(i,j)
do k=km,1,-1
gz(k) = gz(k+1) + rdgas*tp(i,k)*(pn0(i,k+1)-pn0(i,k))
enddo
! Only lowest layer potential temp is needed
pt0(km) = tp(i,km)/(pk0(km+1)-pk0(km))*(kappa*(pn0(i,km+1)-pn0(i,km)))
if( Atm%phis(i,j)>gzc(i,j) ) then
do k=km,1,-1
if( Atm%phis(i,j) < gz(k) .and. &
Atm%phis(i,j) >= gz(k+1) ) then
pst = pk0(k) + (pk0(k+1)-pk0(k))*(gz(k)-Atm%phis(i,j))/(gz(k)-gz(k+1))
go to 123
endif
enddo
else
! Extrapolation into the ground
pst = pk0(km+1) + (gzc(i,j)-Atm%phis(i,j))/(cp_air*pt0(km))
endif
123 Atm%ps(i,j) = pst**(1./kappa)
enddo !i-loop
do i=is,ie
pe1(i,1) = Atm%ak(1)
pn1(i,1) = log(pe1(i,1))
enddo
do k=2,npz+1
do i=is,ie
pe1(i,k) = Atm%ak(k) + Atm%bk(k)*Atm%ps(i,j)
pn1(i,k) = log(pe1(i,k))
enddo
enddo
! * Compute delp
do k=1,npz
do i=is,ie
Atm%delp(i,j,k) = pe1(i,k+1) - pe1(i,k)
enddo
enddo
!---------------
! map tracers
!----------------
do iq=1,ncnst
call mappm(km, pe0, qp(is,1,iq), npz, pe1, qn1, is,ie, 0, 11, Atm%ptop)
do k=1,npz
do i=is,ie
Atm%q(i,j,k,iq) = qn1(i,k)
enddo
enddo
enddo
!---------------
! map extra 3d variables
!----------------
do ifile=1,size(out_fv_rst)
if (out_fv_rst(ifile)%have_descriptor) then
do ivar=1,size(out_fv_rst(ifile)%vars)
if (out_fv_rst(ifile)%vars(ivar)%rank==2) then
out_fv_rst(ifile)%vars(ivar)%ptr2d(is:ie,j)=in_fv_rst(ifile)%vars(ivar)%ptr2d(is:ie,j)
else if (out_fv_rst(ifile)%vars(ivar)%rank==3) then
if (out_fv_rst(ifile)%vars(ivar)%nLev==npz) then
do k=1,in_fv_rst(ifile)%vars(ivar)%nLev
qp1(is:ie,k)=in_fv_rst(ifile)%vars(ivar)%ptr3d(is:ie,j,k)
enddo
call mappm(km, pe0, qp1, npz, pe1, qn1, is,ie, 0, 11, Atm%ptop)
do k=1,npz
do i=is,ie
out_fv_rst(ifile)%vars(ivar)%ptr3d(i,j,k) = qn1(i,k)
enddo
enddo
else if (out_fv_rst(ifile)%vars(ivar)%nLev==npz+1) then
do k=1,in_fv_rst(ifile)%vars(ivar)%nLev
q_edge_old(is:ie,k)=in_fv_rst(ifile)%vars(ivar)%ptr3d(is:ie,j,k)
enddo
call remap_edge(q_edge_old,q_edge_new,is,ie,km,npz,Atm%ak,Atm%bk)
do k=1,npz+1
do i=is,ie
out_fv_rst(ifile)%vars(ivar)%ptr3d(i,j,k) = q_edge_new(i,k)
enddo
enddo
end if
else if (out_fv_rst(ifile)%vars(ivar)%rank==4) then
do n_ungrid=1,out_fv_rst(ifile)%vars(ivar)%n_ungrid
if (out_fv_rst(ifile)%vars(ivar)%nLev==npz) then
do k=1,in_fv_rst(ifile)%vars(ivar)%nLev
qp1(is:ie,k)=in_fv_rst(ifile)%vars(ivar)%ptr4d(is:ie,j,k,n_ungrid)
enddo
call mappm(km, pe0, qp1, npz, pe1, qn1, is,ie, 0, 11, Atm%ptop)
do k=1,npz
do i=is,ie
out_fv_rst(ifile)%vars(ivar)%ptr4d(i,j,k,n_ungrid) = qn1(i,k)
enddo
enddo
else if (out_fv_rst(ifile)%vars(ivar)%nLev==npz+1) then
do k=1,in_fv_rst(ifile)%vars(ivar)%nLev
q_edge_old(is:ie,k)=in_fv_rst(ifile)%vars(ivar)%ptr4d(is:ie,j,k,n_ungrid)
enddo
call remap_edge(q_edge_old,q_edge_new,is,ie,km,npz,Atm%ak,Atm%bk)
do k=1,npz+1
do i=is,ie
out_fv_rst(ifile)%vars(ivar)%ptr4d(i,j,k,n_ungrid) = q_edge_new(i,k)
enddo
enddo
end if
end do
end if
enddo
else
do ivar=1,size(out_fv_rst(ifile)%vars)
out_fv_rst(ifile)%vars(ivar)%ptr3d(is:ie,j,:)=in_fv_rst(ifile)%vars(ivar)%ptr3d(is:ie,j,:)
end do
end if
enddo
!-------------------------------------------------------------
! map virtual temperature using geopotential conserving scheme.
!-------------------------------------------------------------
call mappm(km, pn0, tp, npz, pn1, qn1, is,ie, 1, 9, Atm%ptop)
do k=1,npz
do i=is,ie
Atm%pt(i,j,k) = qn1(i,k)/(1.+zvir*Atm%q(i,j,k,sphum))
enddo
enddo
enddo
call prt_maxmin('PS_model', Atm%ps, is, ie, js, je, ng, 1, 0.01_FVPRC)
if (is_master()) write(*,*) 'done remap_scalar'
end subroutine remap_scalar
subroutine remap_edge(q1,q2,is,ie,km,kn,ak,bk)
! q1,km - old levels
! q2,kn - new levels
integer, intent(in) :: is,ie,km,kn
real(FVPRC),intent(in) :: ak(kn+1), bk(kn+1)
real(FVPRC),intent(in) :: q1(is:ie,km+1)
real(FVPRC),intent(out) :: q2(is:ie,kn+1)
integer i,k
do i=is,ie
do k=1,kn+1
if (bk(k)==0.0) then
q2(i,k)=0.0
else
q2(i,k)=bk(k)*q1(i,km+1)
end if
enddo
enddo
end subroutine remap_edge
end module fv_regridding_utils