SCALE-RM
scale_atmos_sub_profile.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
13 !-------------------------------------------------------------------------------
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use scale_precision
20  use scale_stdio
21  use scale_prof
23 
24  use scale_const, only: &
25  grav => const_grav, &
26  cpdry => const_cpdry, &
27  rdry => const_rdry, &
28  p00 => const_pre00
29  !-----------------------------------------------------------------------------
30  implicit none
31  private
32  !-----------------------------------------------------------------------------
33  !
34  !++ Public procedure
35  !
36  public :: atmos_profile_isa
37 
38  interface atmos_profile_isa
39  module procedure atmos_profile_isa_1d
40  module procedure atmos_profile_isa_3d
41  end interface atmos_profile_isa
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Public parameters & variables
46  !
47  !-----------------------------------------------------------------------------
48  !
49  !++ Private procedure
50  !
51  !-----------------------------------------------------------------------------
52  !
53  !++ Private parameters & variables
54  !
55  ! ISA profile
56  integer, private, parameter :: nref = 8
57  real(RP), private, parameter :: z_isa(nref) = (/ 0.0_rp, &
58  11000.0_rp, &
59  20000.0_rp, &
60  32000.0_rp, &
61  47000.0_rp, &
62  51000.0_rp, &
63  71000.0_rp, &
64  84852.0_rp /)
65  real(RP), private, parameter :: gamma(nref) = (/ -6.5e-3_rp, &
66  0.0_rp, &
67  1.0e-3_rp, &
68  2.8e-3_rp, &
69  0.0e-3_rp, &
70  -2.8e-3_rp, &
71  -2.0e-3_rp, &
72  0.0_rp /)
73 
74  !-----------------------------------------------------------------------------
75 contains
76  !-----------------------------------------------------------------------------
78  subroutine atmos_profile_isa_1d( &
79  KA, KS, KE, &
80  temp_sfc, &
81  pres_sfc, &
82  z, &
83  pott )
84  implicit none
85 
86  integer, intent(in) :: ka, ks, ke
87  real(RP), intent(in) :: temp_sfc
88  real(RP), intent(in) :: pres_sfc
89  real(RP), intent(in) :: z (ka)
90  real(RP), intent(out) :: pott(ka)
91 
92  real(RP) :: temp_isa(nref)
93  real(RP) :: pres_isa(nref)
94  real(RP) :: temp(ka)
95  real(RP) :: pres(ka)
96 
97  real(RP) :: gmr ! grav / Rdry
98  real(RP) :: rovcp ! CPdry / Rdry
99  integer :: k, n
100  !---------------------------------------------------------------------------
101 
102  gmr = grav / rdry
103  rovcp = rdry / cpdry
104 
105  !--- ISA profile
106  temp_isa(1) = temp_sfc
107  pres_isa(1) = pres_sfc
108 
109  do n = 2, nref
110  temp_isa(n) = temp_isa(n-1) + gamma(n-1) * ( z_isa(n)-z_isa(n-1) )
111 
112  if ( gamma(n-1) == 0.0_rp ) then
113  pres_isa(n) = pres_isa(n-1) * exp( -gmr / temp_isa(n) * ( z_isa(n)-z_isa(n-1) ) )
114  else
115  pres_isa(n) = pres_isa(n-1) * ( temp_isa(n)/temp_isa(n-1) ) ** ( -gmr/gamma(n-1) )
116  endif
117  enddo
118 
119  if( io_l ) write(io_fid_log,*)
120  if( io_l ) write(io_fid_log,*) '###### ICAO International Standard Atmosphere ######'
121  if( io_l ) write(io_fid_log,*) ' height: lapse rate: pressure: temperature'
122  do n = 1, nref
123  if( io_l ) write(io_fid_log,'(4F13.5)') z_isa(n), gamma(n), pres_isa(n), temp_isa(n)
124  enddo
125  if( io_l ) write(io_fid_log,*) '####################################################'
126 
127  !--- make reference state
128  do k = ks, ke
129  if ( z(k) <= z_isa(1) ) then
130 
131  temp(k) = temp_isa(1) + gamma(1) * ( z(k)-z_isa(1) )
132  pres(k) = pres_isa(1) * ( temp(k)/temp_isa(1) ) ** ( -gmr/gamma(1) )
133 
134  elseif( z(k) > z_isa(nref) ) then
135 
136  temp(k) = temp_isa(nref)
137  pres(k) = pres_isa(nref) * exp( -gmr/temp_isa(nref) * ( z(k)-z_isa(nref) ) )
138 
139  else
140  do n = 2, nref
141  if ( z(k) > z_isa(n-1) .AND. z(k) <= z_isa(n) ) then
142 
143  temp(k) = temp_isa(n-1) + gamma(n-1) * ( z(k)-z_isa(n-1) )
144  if ( gamma(n-1) == 0.0_rp ) then
145  pres(k) = pres_isa(n-1) * exp( -gmr/temp_isa(n-1) * ( z(k)-z_isa(n-1) ) )
146  else
147  pres(k) = pres_isa(n-1) * ( temp(k)/temp_isa(n-1) ) ** ( -gmr/gamma(n-1) )
148  endif
149 
150  endif
151  enddo
152  endif
153 
154  pott(k) = temp(k) * ( p00/pres(k) )**rovcp
155  enddo
156 
157  return
158  end subroutine atmos_profile_isa_1d
159 
160  !-----------------------------------------------------------------------------
162  subroutine atmos_profile_isa_3d( &
163  KA, KS, KE, &
164  IA, IS, IE, &
165  JA, JS, JE, &
166  temp_sfc, &
167  pres_sfc, &
168  z, &
169  pott )
170  implicit none
171 
172  integer, intent(in) :: ka, ks, ke
173  integer, intent(in) :: ia, is, ie
174  integer, intent(in) :: ja, js, je
175  real(RP), intent(in) :: temp_sfc(ia,ja)
176  real(RP), intent(in) :: pres_sfc(ia,ja)
177  real(RP), intent(in) :: z (ka,ia,ja)
178  real(RP), intent(out) :: pott (ka,ia,ja)
179 
180  real(RP) :: temp_isa(nref,ia,ja)
181  real(RP) :: pres_isa(nref,ia,ja)
182  real(RP) :: temp(ka)
183  real(RP) :: pres(ka)
184 
185  real(RP) :: gmr ! grav / Rdry
186  real(RP) :: rovcp ! CPdry / Rdry
187  integer :: k, i, j, n
188  !---------------------------------------------------------------------------
189 
190  if( io_l ) write(io_fid_log,*)
191  if( io_l ) write(io_fid_log,*) '++++++ Module[PROFILE] / Categ[ATMOS SHARE] / Origin[SCALElib]'
192 
193  gmr = grav / rdry
194  rovcp = rdry / cpdry
195 
196  !--- ISA profile
197  do j = js, je
198  do i = is, ie
199  temp_isa(1,i,j) = temp_sfc(i,j)
200  pres_isa(1,i,j) = pres_sfc(i,j)
201  enddo
202  enddo
203 
204  do j = js, je
205  do i = is, ie
206  do n = 2, nref
207  temp_isa(n,i,j) = temp_isa(n-1,i,j) + gamma(n-1) * ( z_isa(n)-z_isa(n-1) )
208 
209  if ( gamma(n-1) == 0.0_rp ) then
210  pres_isa(n,i,j) = pres_isa(n-1,i,j) * exp( -gmr / temp_isa(n,i,j) * ( z_isa(n)-z_isa(n-1) ) )
211  else
212  pres_isa(n,i,j) = pres_isa(n-1,i,j) * ( temp_isa(n,i,j)/temp_isa(n-1,i,j) ) ** ( -gmr/gamma(n-1) )
213  endif
214  enddo
215  enddo
216  enddo
217 
218  if( io_l ) write(io_fid_log,*)
219  if( io_l ) write(io_fid_log,*) '###### ICAO International Standard Atmosphere ######'
220  if( io_l ) write(io_fid_log,*) ' height: lapse rate: pressure: temperature'
221  do n = 1, nref
222  if( io_l ) write(io_fid_log,'(4F13.5)') z_isa(n), gamma(n), pres_isa(n,is,js), temp_isa(n,is,js)
223  enddo
224  if( io_l ) write(io_fid_log,*) '####################################################'
225 
226  !--- make reference state
227  do j = js, je
228  do i = is, ie
229  do k = ks, ke
230  if ( z(k,i,j) <= z_isa(1) ) then
231 
232  temp(k) = temp_isa(1,i,j) + gamma(1) * ( z(k,i,j)-z_isa(1) )
233  pres(k) = pres_isa(1,i,j) * ( temp(k)/temp_isa(1,i,j) ) ** ( -gmr/gamma(1) )
234 
235  elseif ( z(k,i,j) > z_isa(nref) ) then
236 
237  temp(k) = temp_isa(nref,i,j)
238  pres(k) = pres_isa(nref,i,j) * exp( -gmr/temp_isa(nref,i,j) * ( z(k,i,j)-z_isa(nref) ) )
239 
240  else
241  do n = 2, nref
242  if ( z(k,i,j) > z_isa(n-1) .AND. z(k,i,j) <= z_isa(n) ) then
243 
244  temp(k) = temp_isa(n-1,i,j) + gamma(n-1) * ( z(k,i,j)-z_isa(n-1) )
245  if ( gamma(n-1) == 0.0_rp ) then
246  pres(k) = pres_isa(n-1,i,j) * exp( -gmr/temp_isa(n-1,i,j) * ( z(k,i,j)-z_isa(n-1) ) )
247  else
248  pres(k) = pres_isa(n-1,i,j) * ( temp(k)/temp_isa(n-1,i,j) ) ** ( -gmr/gamma(n-1) )
249  endif
250 
251  endif
252  enddo
253  endif
254 
255  pott(k,i,j) = temp(k) * ( p00/pres(k) )**rovcp
256  enddo
257  enddo
258  enddo
259 
260  return
261  end subroutine atmos_profile_isa_3d
262 
263 end module scale_atmos_profile
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
real(rp), public const_cpdry
specific heat (dry air,constant pressure) [J/kg/K]
Definition: scale_const.F90:58
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), public const_rdry
specific gas constant (dry air) [J/kg/K]
Definition: scale_const.F90:57
module grid index
module ATMOSPHERE / Typical vertical profile
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
real(rp), public const_pre00
pressure reference [Pa]
Definition: scale_const.F90:93
real(rp), public const_grav
standard acceleration of gravity [m/s2]
Definition: scale_const.F90:48
integer, public js
start point of inner domain: y, local
module CONSTANT
Definition: scale_const.F90:14
integer, public ks
start point of inner domain: z, local
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
module PRECISION
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public ja
of y whole cells (local, with HALO)