FV3 Bundle
fv_pressure.F90
Go to the documentation of this file.
1 ! (C) Copyright 2018 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
6 !> Variable transforms on pressure variables for fv3
7 !> Daniel Holdaway, NASA/JCSDA
8 
10 
11 use mapl_constantsmod, only: kappa => mapl_kappa
12 
13 implicit none
14 public
15 
16 contains
17 
18 !----------------------------------------------------------------------------
19 ! Compute all pressures needed as input to fv3 ------------------------------
20 !----------------------------------------------------------------------------
21 
22 subroutine compute_fv3_pressures( is, ie, js, je, isd, ied, jsd, jed, npz, &
23  kappa, ptop, delp, pe, pk, pkz, peln )
24 
25  implicit none
26  integer, intent(in) :: is, ie, js, je, isd, ied, jsd, jed, npz
27  real, intent(in) :: kappa, ptop
28  real, intent(in) :: delp(isd:ied, jsd:jed, npz)
29  real, intent(out) :: pe(is-1:ie+1, npz+1, js-1:je+1)
30  real, intent(out) :: pk(is:ie, js:je, npz+1)
31  real, intent(out) :: peln(is:ie, npz+1, js:je)
32  real, intent(out) :: pkz(is:ie, js:je, npz)
33  integer :: i, j, k
34 
35  pe(:, :, :) = 0.0
36  pe(:, 1, :) = ptop
37  do k=2,npz+1
38  do j=js,je
39  do i=is,ie
40  pe(i, k, j) = pe(i, k-1, j) + delp(i, j, k-1)
41  end do
42  end do
43  end do
44 
45  do k=1,npz+1
46  do j=js,je
47  do i=is,ie
48  peln(i, k, j) = log(pe(i, k, j))
49  end do
50  end do
51  end do
52 
53  do k=1,npz+1
54  do j=js,je
55  do i=is,ie
56  pk(i, j, k) = exp(kappa*peln(i, k, j))
57  end do
58  end do
59  end do
60 
61  do k=1,npz
62  do j=js,je
63  do i=is,ie
64  pkz(i, j, k) = (pk(i, j, k+1)-pk(i, j, k))/(kappa*(peln(i, k+1, j)-peln(i, k, j)))
65  end do
66  end do
67  end do
68 
69 end subroutine compute_fv3_pressures
70 
71 subroutine compute_fv3_pressures_tlm( is, ie, js, je, isd, ied, jsd, jed, npz, &
72  kappa, ptop, delp, delp_tl, &
73  pe, pe_tl, pk, pk_tl, pkz, pkz_tl, peln, peln_tl )
74 
75  implicit none
76  integer, intent(in) :: is, ie, js, je, isd, ied, jsd, jed, npz
77  real, intent(in) :: kappa, ptop
78  real, intent(in) :: delp(isd:ied, jsd:jed, npz)
79  real, intent(in) :: delp_tl(isd:ied, jsd:jed, npz)
80  real, intent(out) :: pe(is-1:ie+1, npz+1, js-1:je+1)
81  real, intent(out) :: pe_tl(is-1:ie+1, npz+1, js-1:je+1)
82  real, intent(out) :: pk(is:ie, js:je, npz+1)
83  real, intent(out) :: pk_tl(is:ie, js:je, npz+1)
84  real, intent(out) :: peln(is:ie, npz+1, js:je)
85  real, intent(out) :: peln_tl(is:ie, npz+1, js:je)
86  real, intent(out) :: pkz(is:ie, js:je, npz)
87  real, intent(out) :: pkz_tl(is:ie, js:je, npz)
88  integer :: i, j, k
89 
90  pe(:, :, :) = 0.0
91  pe(:, 1, :) = ptop
92  pe_tl = 0.0
93  do k=2,npz+1
94  do j=js,je
95  do i=is,ie
96  pe_tl(i, k, j) = pe_tl(i, k-1, j) + delp_tl(i, j, k-1)
97  pe(i, k, j) = pe(i, k-1, j) + delp(i, j, k-1)
98  end do
99  end do
100  end do
101 
102  peln_tl = 0.0
103  do k=1,npz+1
104  do j=js,je
105  do i=is,ie
106  peln_tl(i, k, j) = pe_tl(i, k, j)/pe(i, k, j)
107  peln(i, k, j) = log(pe(i, k, j))
108  end do
109  end do
110  end do
111 
112  pk_tl = 0.0
113  do k=1,npz+1
114  do j=js,je
115  do i=is,ie
116  pk_tl(i, j, k) = kappa*peln_tl(i, k, j)*exp(kappa*peln(i, k, j))
117  pk(i, j, k) = exp(kappa*peln(i, k, j))
118  end do
119  end do
120  end do
121 
122  pkz_tl = 0.0
123  do k=1,npz
124  do j=js,je
125  do i=is,ie
126  pkz_tl(i, j, k) = ((pk_tl(i, j, k+1)-pk_tl(i, j, k))*kappa*(peln(i, k+1, j)-peln(i, k, j)) &
127  -(pk(i, j, k+1)-pk(i, j, k))*kappa*(peln_tl(i, k+1, j)-peln_tl(i, k, j))) &
128  /(kappa*(peln(i, k+1, j)-peln(i, k, j)))**2
129  pkz(i, j, k) = (pk(i, j, k+1)-pk(i, j, k))/(kappa*(peln(i, k+1, j)-peln(i, k, j)))
130  end do
131  end do
132  end do
133 
134 end subroutine compute_fv3_pressures_tlm
135 
136 subroutine compute_fv3_pressures_bwd( is, ie, js, je, isd, ied, jsd, jed, npz, &
137  kappa, ptop, delp, delp_ad, &
138  pe, pe_ad, pk, pk_ad, pkz, pkz_ad, peln, peln_ad)
140  implicit none
141  integer, intent(in) :: is, ie, js, je, isd, ied, jsd, jed, npz
142  real, intent(in) :: kappa, ptop
143  real, intent(in) :: delp(isd:ied, jsd:jed, npz)
144  real, intent(inout) :: delp_ad(isd:ied, jsd:jed, npz)
145  real, intent(in) :: pe(is-1:ie+1, npz+1, js-1:je+1)
146  real, intent(inout) :: pe_ad(is-1:ie+1, npz+1, js-1:je+1)
147  real, intent(in) :: pk(is:ie, js:je, npz+1)
148  real, intent(inout) :: pk_ad(is:ie, js:je, npz+1)
149  real, intent(in) :: peln(is:ie, npz+1, js:je)
150  real, intent(inout) :: peln_ad(is:ie, npz+1, js:je)
151  real, intent(in) :: pkz(is:ie, js:je, npz)
152  real, intent(inout) :: pkz_ad(is:ie, js:je, npz)
153  integer :: i, j, k
154  real :: temp_tj1, temp_ad1, temp_ad2
155 
156  do k=npz,1,-1
157  do j=je,js,-1
158  do i=ie,is,-1
159  temp_tj1 = kappa*(peln(i, k+1, j)-peln(i, k, j))
160  temp_ad1 = pkz_ad(i, j, k)/temp_tj1
161  temp_ad2 = -((pk(i, j, k+1)-pk(i, j, k))*kappa*temp_ad1/temp_tj1)
162  pk_ad(i, j, k+1) = pk_ad(i, j, k+1) + temp_ad1
163  pk_ad(i, j, k) = pk_ad(i, j, k) - temp_ad1
164  peln_ad(i, k+1, j) = peln_ad(i, k+1, j) + temp_ad2
165  peln_ad(i, k, j) = peln_ad(i, k, j) - temp_ad2
166  pkz_ad(i, j, k) = 0.0
167  end do
168  end do
169  end do
170 
171  do k=npz+1,1,-1
172  do j=je,js,-1
173  do i=ie,is,-1
174  peln_ad(i, k, j) = peln_ad(i, k, j) + exp(kappa*peln(i, k, j))*kappa*pk_ad(i, j, k)
175  pk_ad(i, j, k) = 0.0
176  end do
177  end do
178  end do
179 
180  do k=npz+1,1,-1
181  do j=je,js,-1
182  do i=ie,is,-1
183  pe_ad(i, k, j) = pe_ad(i, k, j) + peln_ad(i, k, j)/pe(i, k, j)
184  peln_ad(i, k, j) = 0.0
185  end do
186  end do
187  end do
188 
189  do k=npz+1,2,-1
190  do j=je,js,-1
191  do i=ie,is,-1
192  pe_ad(i, k-1, j) = pe_ad(i, k-1, j) + pe_ad(i, k, j)
193  delp_ad(i, j, k-1) = delp_ad(i, j, k-1) + pe_ad(i, k, j)
194  pe_ad(i, k, j) = 0.0
195  end do
196  end do
197  end do
198 
199 end subroutine compute_fv3_pressures_bwd
200 
201 !----------------------------------------------------------------------------
202 
203 end module fv_pressure_mod
subroutine compute_fv3_pressures(is, ie, js, je, isd, ied, jsd, jed, npz, kappa, ptop, delp, pe, pk, pkz, peln)
Definition: fv_pressure.F90:24
real, parameter, public mapl_kappa
Variable transforms on pressure variables for fv3 Daniel Holdaway, NASA/JCSDA.
Definition: fv_pressure.F90:9
subroutine compute_fv3_pressures_bwd(is, ie, js, je, isd, ied, jsd, jed, npz, kappa, ptop, delp, delp_ad, pe, pe_ad, pk, pk_ad, pkz, pkz_ad, peln, peln_ad)
subroutine compute_fv3_pressures_tlm(is, ie, js, je, isd, ied, jsd, jed, npz, kappa, ptop, delp, delp_tl, pe, pe_tl, pk, pk_tl, pkz, pkz_tl, peln, peln_tl)
Definition: fv_pressure.F90:74