OOPS
qg_tools_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
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 ! In applying this licence, ECMWF does not waive the privileges and immunities
6 ! granted to it by virtue of its status as an intergovernmental organisation nor
7 ! does it submit to any jurisdiction.
8 
10 
11 use fckit_configuration_module, only: fckit_configuration
12 use datetime_mod
13 use duration_mod
14 use iso_c_binding
15 use kinds
16 use netcdf
18 
19 implicit none
20 
21 private
23 ! ------------------------------------------------------------------------------
24 real(kind_real),parameter :: ubot = -2.0_kind_real !< Zonal wind at the surface (m/s)
25 real(kind_real),parameter :: utop = 58.0_kind_real !< Zonal wind at the top (m/s)
26 ! ------------------------------------------------------------------------------
27 contains
28 ! ------------------------------------------------------------------------------
29 !> Generate filename
30 function genfilename(f_conf,length,vdate)
31 use string_utils
32 
33 implicit none
34 
35 ! Passed variables
36 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
37 integer,intent(in) :: length !< Length
38 type(datetime),intent(in) :: vdate !< Date and time
39 
40 ! Result
41 character(len=2*length) :: genfilename
42 
43 ! Local variables
44 integer :: lenfn
45 character(len=length) :: fdbdir,expver,typ,validitydate,referencedate,sstep,mmb,iter
46 character(len=2*length) :: prefix
47 character(len=:),allocatable :: str
48 
49 type(datetime) :: rdate
50 type(duration) :: step
51 
52 ! Get configuration parameters
53 call f_conf%get_or_die("datadir",str)
54 call swap_name_member(f_conf, str)
55 fdbdir = str
56 call f_conf%get_or_die("exp",str)
57 call swap_name_member(f_conf, str)
58 expver = str
59 call f_conf%get_or_die("type",str)
60 typ = str
61 
62 ! Ensemble case
63 if (typ=='ens') then
64  call f_conf%get_or_die("member",str)
65  mmb = str
66  lenfn = len_trim(fdbdir) + 1 + len_trim(expver) + 1 + len_trim(typ) + 1 + len_trim(mmb)
67  prefix = trim(fdbdir) // '/' // trim(expver) // '.' // trim(typ) // '.' // trim(mmb)
68 else
69  lenfn = len_trim(fdbdir) + 1 + len_trim(expver) + 1 + len_trim(typ)
70  prefix = trim(fdbdir) // '/' // trim(expver) // '.' // trim(typ)
71 endif
72 
73 ! Forecast / ensemble cases
74 if ((typ=='fc').or.(typ=='ens')) then
75  call f_conf%get_or_die("date",str)
76  referencedate = str
77  call datetime_to_string(vdate,validitydate)
78  call datetime_create(trim(referencedate),rdate)
79  call datetime_diff(vdate,rdate,step)
80  call duration_to_string(step,sstep)
81  lenfn = lenfn+1+len_trim(referencedate)+1+len_trim(sstep)
82  genfilename = trim(prefix)//'.'//trim(referencedate)//'.'// trim(sstep)//'.nc'
83 endif
84 
85 ! Analysis or increment case
86 if ((typ=='an').or.(typ=='in')) then
87  call datetime_to_string(vdate,validitydate)
88  lenfn = lenfn+1+len_trim(validitydate)
89  genfilename = trim(prefix)//'.'//trim(validitydate)//'.nc'
90 endif
91 
92 if (typ=='krylov') then
93  call f_conf%get_or_die("iteration",str)
94  iter = str
95  call datetime_to_string(vdate,validitydate)
96  lenfn = lenfn+1+len_trim(iter)+1+len_trim(validitydate)
97  genfilename = trim(prefix)//'.'// trim(iter)//'.'//trim(validitydate)//'.nc'
98 endif
99 
100 ! Check filename length
101 if (lenfn>length) call abor1_ftn('genfilename: filename too long')
102 
103 end function genfilename
104 ! ------------------------------------------------------------------------------
105 !> Check NetCDF status
106 subroutine ncerr(info)
107 
108 implicit none
109 
110 ! Passed variables
111 integer,intent(in) :: info !< Info index
112 
113 ! Check status
114 if (info/=nf90_noerr) call abor1_ftn(trim(nf90_strerror(info)))
115 
116 end subroutine ncerr
117 ! ------------------------------------------------------------------------------
118 !> Generate values for baroclinic instability
119 subroutine baroclinic_instability(x,y,z,var,res)
120 
121 implicit none
122 
123 ! Passed variables
124 real(kind_real),intent(in) :: x !< X value
125 real(kind_real),intent(in) :: y !< Y value
126 real(kind_real),intent(in) :: z !< Z value
127 character(len=1),intent(in) :: var !< Variable
128 real(kind_real),intent(out) :: res !< Results
129 
130 ! Local variable
131 real(kind_real) :: u
132 
133 ! Define zonal wind
134 u = ubot+(utop-ubot)*z/domain_depth
135 
136 select case (var)
137 case ('x')
138  ! Streamfunction
139  res = -u*y
140 case ('q')
141  call abor1_ftn('baroclinic_instability: cannot define q')
142 case ('u')
143  ! Zonal wind
144  res = u
145 case ('v')
146  ! Meridional wind
147  res = 0.0
148 case default
149  call abor1_ftn('baroclinic_instability: wrong variable')
150 end select
151 
152 end subroutine baroclinic_instability
153 ! ------------------------------------------------------------------------------
154 !> Generate values for large vortices
155 subroutine large_vortices(x,y,z,var,res)
156 
157 implicit none
158 
159 ! Passed variables
160 real(kind_real),intent(in) :: x !< X value
161 real(kind_real),intent(in) :: y !< Y value
162 real(kind_real),intent(in) :: z !< Z value
163 character(len=1),intent(in) :: var !< Variable
164 real(kind_real),intent(out) :: res !< Results
165 
166 ! Local variable
167 real(kind_real) :: ff
168 
169 ! Define wind speed
170 ff = ubot+(utop-ubot)*z/domain_depth
171 
172 select case (var)
173 case ('x')
174  ! Streamfunction
175  res = ff*domain_meridional*cos(2.0*pi*x/domain_zonal)*sin(pi*y/domain_meridional)
176 case ('q')
177  call abor1_ftn('large_vortices: cannot define q')
178 case ('u')
179  ! Zonal wind
180  res = -ff*pi*cos(2.0*pi*x/domain_zonal)*cos(pi*y/domain_meridional)
181 case ('v')
182  ! Meridional wind
183  res = -2.0*ff*pi*domain_meridional/domain_zonal*sin(2.0*pi*x/domain_zonal)*sin(pi*y/domain_meridional)
184 case default
185  call abor1_ftn('large_vortices: wrong variable')
186 end select
187 
188 end subroutine large_vortices
189 ! ------------------------------------------------------------------------------
190 end module qg_tools_mod
real(kind_real), parameter, public pi
Pi.
real(kind_real), parameter, public domain_depth
Model depth (m)
real(kind_real), parameter, public domain_meridional
Model domain in meridional direction (m)
real(kind_real), parameter, public domain_zonal
Model domain in zonal direction (m)
subroutine, public large_vortices(x, y, z, var, res)
Generate values for large vortices.
real(kind_real), parameter ubot
Zonal wind at the surface (m/s)
subroutine, public baroclinic_instability(x, y, z, var, res)
Generate values for baroclinic instability.
subroutine, public ncerr(info)
Check NetCDF status.
character(len=2 *length) function, public genfilename(f_conf, length, vdate)
Generate filename.
real(kind_real), parameter utop
Zonal wind at the top (m/s)