OOPS
qg_geom_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 
9 module qg_geom_mod
10 
11 use atlas_module, only: atlas_field, atlas_fieldset, atlas_real, atlas_functionspace_pointcloud
12 use fckit_configuration_module, only: fckit_configuration
13 use fckit_log_module,only: fckit_log
14 use kinds
15 use iso_c_binding
18 
19 implicit none
20 
21 private
22 public :: qg_geom
23 public :: qg_geom_registry
25 ! ------------------------------------------------------------------------------
26 type :: qg_geom
27  integer :: nx !< Number of points in the zonal direction
28  integer :: ny !< Number of points in the meridional direction
29  integer :: nz !< Number of vertical levels
30  real(kind_real) :: deltax !< Zonal cell size
31  real(kind_real) :: deltay !< Meridional cell size
32  real(kind_real),allocatable :: x(:) !< Zonal coordinate
33  real(kind_real),allocatable :: y(:) !< Meridional coordinate
34  real(kind_real),allocatable :: z(:) !< Altitude
35  real(kind_real),allocatable :: lat(:,:) !< Latitude
36  real(kind_real),allocatable :: lon(:,:) !< Longitude
37  real(kind_real),allocatable :: area(:,:) !< Area
38  real(kind_real),allocatable :: f(:,:) !< Coefficients of PV operator
39  real(kind_real),allocatable :: f_p(:,:) !< Coefficients of PV operator, right eigenvectors
40  real(kind_real),allocatable :: f_pinv(:,:) !< Coefficients of PV operator, right eigenvectors inverse
41  real(kind_real),allocatable :: f_d(:) !< Coefficients of PV operator, eigenvalues
42  real(kind_real),allocatable :: bet(:) !< Beta coefficient
43  real(kind_real),allocatable :: heat(:,:) !< Heating term
44  type(atlas_functionspace_pointcloud) :: afunctionspace !< ATLAS function space
45 end type qg_geom
46 
47 #define LISTED_TYPE qg_geom
48 
49 !> Linked list interface - defines registry_t type
50 #include "oops/util/linkedList_i.f"
51 
52 !> Global registry
53 type(registry_t) :: qg_geom_registry
54 ! ------------------------------------------------------------------------------
55 contains
56 ! ------------------------------------------------------------------------------
57 ! Public
58 ! ------------------------------------------------------------------------------
59 !> Linked list implementation
60 #include "oops/util/linkedList_c.f"
61 ! ------------------------------------------------------------------------------
62 !> Setup geometry
63 subroutine qg_geom_setup(self,f_conf)
64 
65 ! Passed variables
66 type(qg_geom),intent(inout) :: self !< Geometry
67 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
68 
69 ! Local variables
70 integer :: ix,iy,iz,ix_c,iy_c,lwork,info
71 integer,allocatable :: ipiv(:),ipivsave(:)
72 real(kind_real) :: mapfac,distx,disty,f
73 real(kind_real),allocatable :: real_array(:),depths(:),wi(:),vl(:,:),work(:)
74 real(kind_real),allocatable :: fsave(:,:),vrlu(:,:),vrlusave(:,:)
75 character(len=1024) :: record
76 logical :: htype
77 character(len=:),allocatable :: str
78 
79 ! Get horizontal resolution data
80 call f_conf%get_or_die("nx",self%nx)
81 call f_conf%get_or_die("ny",self%ny)
82 self%nz = f_conf%get_size("depths")
83 
84 ! Allocation
85 allocate(depths(self%nz))
86 allocate(self%x(self%nx))
87 allocate(self%y(self%ny))
88 allocate(self%z(self%nz))
89 allocate(self%lon(self%nx,self%ny))
90 allocate(self%lat(self%nx,self%ny))
91 allocate(self%area(self%nx,self%ny))
92 allocate(self%f(self%nz,self%nz))
93 allocate(self%f_p(self%nz,self%nz))
94 allocate(self%f_pinv(self%nz,self%nz))
95 allocate(self%f_d(self%nz))
96 allocate(self%bet(self%ny))
97 allocate(self%heat(self%nx,self%ny))
98 allocate(wi(self%nz))
99 allocate(vl(self%nz,self%nz))
100 allocate(fsave(self%nz,self%nz))
101 allocate(vrlu(self%nz,self%nz))
102 allocate(vrlusave(self%nz,self%nz))
103 allocate(ipiv(self%nz))
104 allocate(ipivsave(self%nz))
105 
106 ! Get depths
107 call f_conf%get_or_die("depths",real_array)
108 depths = real_array
109 
110 ! Define dx/dy
111 self%deltax = domain_zonal/real(self%nx,kind_real)
112 self%deltay = domain_meridional/real(self%ny+1,kind_real)
113 
114 ! Print sizes and dx/dy
115 write(record,'(a,i3,a,i3,a,i3)') 'qg_geom_create: nx/ny/nz = ',self%nx,' /',self%ny,' /',self%nz
116 call fckit_log%info(record)
117 write(record,'(a,f7.2,a,f7.2,a)') ' deltax/deltay = ',self%deltax*1.0e-3,' km / ',self%deltay*1.0e-3,' km'
118 call fckit_log%info(record)
119 
120 ! Define x/y
121 do ix=1,self%nx
122  self%x(ix) = (real(ix,kind_real)-0.5)*self%deltax
123 enddo
124 do iy=1,self%ny
125  self%y(iy) = real(iy,kind_real)*self%deltay
126 enddo
127 
128 ! Define lon/lat/area
129 do iy=1,self%ny
130  do ix=1,self%nx
131  call xy_to_lonlat(self%x(ix),self%y(iy),self%lon(ix,iy),self%lat(ix,iy),mapfac)
132  self%area(ix,iy) = self%deltax*self%deltay/mapfac
133  end do
134 end do
135 
136 ! Set heights
137 self%z(1) = 0.5*depths(1)
138 do iz=2,self%nz
139  self%z(iz) = sum(depths(1:iz-1))+0.5*depths(iz)
140 end do
141 
142 ! Coefficients of PV operator
143 self%f = 0.0
144 do iz=1,self%nz
145  f = f0**2/(g*dlogtheta*depths(iz))
146  if (iz>1) then
147  self%f(iz,iz-1) = f
148  self%f(iz,iz) = self%f(iz,iz)-f
149  end if
150  if (iz<self%nz) then
151  self%f(iz,iz+1) = f
152  self%f(iz,iz) = self%f(iz,iz)-f
153  end if
154 enddo
155 
156 ! Compute eigendecomposition of ff
157 fsave = self%f
158 allocate(work(1))
159 call dgeev('V','V',self%nz,fsave,self%nz,self%f_d,wi,vl,self%nz,self%f_p,self%nz,work,-1,info)
160 if (info/=0) call abor1_ftn('error in dgeev, first pass')
161 lwork = int(work(1))
162 deallocate(work)
163 allocate(work(lwork))
164 fsave = self%f
165 call dgeev('V','V',self%nz,fsave,self%nz,self%f_d,wi,vl,self%nz,self%f_p,self%nz,work,lwork,info)
166 if (info/=0) call abor1_ftn('error in dgeev, second pass')
167 deallocate(work)
168 
169 ! Compute inverse of right eigenvectors of ff
170 vrlu = self%f_p
171 call dgetrf(self%nz,self%nz,vrlu,self%nz,ipiv,info)
172 if (info/=0) call abor1_ftn('error in dgetrf')
173 allocate(work(1))
174 vrlusave = vrlu
175 ipivsave = ipiv
176 call dgetri(self%nz,vrlusave,self%nz,ipivsave,work,-1,info)
177 if (info/=0) call abor1_ftn('error in dgetri, first pass')
178 lwork = int(work(1))
179 deallocate(work)
180 allocate(work(lwork))
181 self%f_pinv = vrlu
182 ipivsave = ipiv
183 call dgetri(self%nz,self%f_pinv,self%nz,ipivsave,work,lwork,info)
184 if (info/=0) call abor1_ftn('error in dgetri, second pass')
185 deallocate(work)
186 
187 ! Beta coefficient
188 do iy=1,self%ny
189  self%bet(iy) = real(iy-(self%ny+1)/2,kind_real)*self%deltay*bet0
190 enddo
191 
192 ! Set heating term
193 call f_conf%get_or_die("heating",htype)
194 if (.not. htype) then
195  ! No heating term
196  call fckit_log%info('qg_geom_setup: heating off')
197  self%heat = 0.0
198 else
199  call fckit_log%info('qg_geom_setup: Gaussian heating on')
200  ! Gaussian source
201  ix_c = self%nx/4
202  iy_c = 3*self%ny/4
203  do iy=1,self%ny
204  do ix=1,self%nx
205  distx = abs(self%x(ix)-self%x(ix_c))
206  if (distx>0.5*domain_zonal) distx = domain_zonal-distx
207  disty = abs(self%y(iy)-self%y(iy_c))
208  self%heat(ix,iy) = heating_amplitude*exp(-(distx**2+disty**2)/heating_scale**2)
209  enddo
210  enddo
211 endif
212 
213 end subroutine qg_geom_setup
214 ! ------------------------------------------------------------------------------
215 !> Set ATLAS lon/lat field
216 subroutine qg_geom_set_atlas_lonlat(self,afieldset)
217 
218 ! Passed variables
219 type(qg_geom),intent(inout) :: self !< Geometry
220 type(atlas_fieldset),intent(inout) :: afieldset !< ATLAS fieldset
221 
222 ! Local variables
223 integer :: ix,iy,inode
224 real(kind_real), pointer :: real_ptr(:,:)
225 type(atlas_field) :: afield
226 
227 ! Create lon/lat field
228 afield = atlas_field(name="lonlat", kind=atlas_real(kind_real), shape=(/2,self%nx*self%ny/))
229 call afield%data(real_ptr)
230 inode = 0
231 do iy=1,self%ny
232  do ix=1,self%nx
233  inode = inode+1
234  real_ptr(1,inode) = self%lon(ix,iy)
235  real_ptr(2,inode) = self%lat(ix,iy)
236  end do
237 end do
238 call afieldset%add(afield)
239 call afield%final()
240 
241 end subroutine qg_geom_set_atlas_lonlat
242 ! ------------------------------------------------------------------------------
243 !> Fill ATLAS fieldset
244 subroutine qg_geom_fill_atlas_fieldset(self,afieldset)
245 
246 ! Passed variables
247 type(qg_geom),intent(inout) :: self !< Geometry
248 type(atlas_fieldset),intent(inout) :: afieldset !< ATLAS fieldset
249 
250 ! Local variables
251 integer :: ix,iy,iz,inode
252 real(kind_real),pointer :: real_ptr_1(:),real_ptr_2(:,:)
253 type(atlas_field) :: afield
254 
255 ! Add area
256 afield = self%afunctionspace%create_field(name='area',kind=atlas_real(kind_real),levels=0)
257 call afield%data(real_ptr_1)
258 inode = 0
259 do iy=1,self%ny
260  do ix=1,self%nx
261  inode = inode+1
262  real_ptr_1(inode) = self%area(ix,iy)
263  enddo
264 enddo
265 call afieldset%add(afield)
266 call afield%final()
267 
268 ! Add vertical unit
269 afield = self%afunctionspace%create_field(name='vunit',kind=atlas_real(kind_real),levels=self%nz)
270 call afield%data(real_ptr_2)
271 do iz=1,self%nz
272  real_ptr_2(iz,1:self%nx*self%ny) = self%z(iz)
273 end do
274 call afieldset%add(afield)
275 call afield%final()
276 
277 end subroutine qg_geom_fill_atlas_fieldset
278 ! ------------------------------------------------------------------------------
279 !> Clone geometry
280 subroutine qg_geom_clone(self,other)
281 
282 ! Passed variables
283 type(qg_geom),intent(inout) :: self !< Geometry
284 type(qg_geom),intent(in) :: other !< Other geometry
285 
286 ! Copy dimensions
287 self%nx = other%nx
288 self%ny = other%ny
289 self%nz = other%nz
290 
291 ! Allocation
292 allocate(self%x(self%nx))
293 allocate(self%y(self%ny))
294 allocate(self%z(self%nz))
295 allocate(self%lon(self%nx,self%ny))
296 allocate(self%lat(self%nx,self%ny))
297 allocate(self%area(self%nx,self%ny))
298 allocate(self%f(self%nz,self%nz))
299 allocate(self%f_p(self%nz,self%nz))
300 allocate(self%f_pinv(self%nz,self%nz))
301 allocate(self%f_d(self%nz))
302 allocate(self%bet(self%ny))
303 allocate(self%heat(self%nx,self%ny))
304 
305 ! Copy data
306 self%deltax = other%deltax
307 self%deltay = other%deltay
308 self%x = other%x
309 self%y = other%y
310 self%z = other%z
311 self%lon = other%lon
312 self%lat = other%lat
313 self%area = other%area
314 self%f = other%f
315 self%f_p = other%f_p
316 self%f_pinv = other%f_pinv
317 self%f_d = other%f_d
318 self%bet = other%bet
319 self%heat = other%heat
320 self%afunctionspace = atlas_functionspace_pointcloud(other%afunctionspace%c_ptr())
321 
322 end subroutine qg_geom_clone
323 ! ------------------------------------------------------------------------------
324 !> Delete geometry
325 subroutine qg_geom_delete(self)
326 
327 ! Passed variables
328 type(qg_geom),intent(inout) :: self !< Geometry
329 
330 ! Release memory
331 deallocate(self%x)
332 deallocate(self%y)
333 deallocate(self%z)
334 deallocate(self%lon)
335 deallocate(self%lat)
336 deallocate(self%area)
337 deallocate(self%f)
338 deallocate(self%f_p)
339 deallocate(self%f_pinv)
340 deallocate(self%f_d)
341 deallocate(self%bet)
342 deallocate(self%heat)
343 call self%afunctionspace%final()
344 
345 end subroutine qg_geom_delete
346 ! ------------------------------------------------------------------------------
347 !> Get geometry info
348 subroutine qg_geom_info(self,nx,ny,nz,deltax,deltay)
349 
350 ! Passed variables
351 type(qg_geom),intent(in) :: self !< Geometry
352 integer,intent(out) :: nx !< Number of points in the zonal direction
353 integer,intent(out) :: ny !< Number of points in the meridional direction
354 integer,intent(out) :: nz !< Number of vertical levels
355 real(kind_real),intent(out) :: deltax !< Zonal cell size
356 real(kind_real),intent(out) :: deltay !< Meridional cell size
357 
358 ! Copy data
359 nx = self%nx
360 ny = self%ny
361 nz = self%nz
362 deltax = self%deltax
363 deltay = self%deltay
364 
365 end subroutine qg_geom_info
366 ! ------------------------------------------------------------------------------
367 end module qg_geom_mod
qg_geom_mod
Definition: qg_geom_mod.F90:9
qg_projection_mod
Definition: qg_projection_mod.F90:9
qg_constants_mod::domain_zonal
real(kind_real), parameter, public domain_zonal
Model domain in zonal direction (m)
Definition: qg_constants_mod.F90:32
qg_geom_mod::qg_geom_info
subroutine, public qg_geom_info(self, nx, ny, nz, deltax, deltay)
Get geometry info.
Definition: qg_geom_mod.F90:349
qg_geom_mod::qg_geom_clone
subroutine, public qg_geom_clone(self, other)
Clone geometry.
Definition: qg_geom_mod.F90:281
qg_constants_mod::bet0
real(kind_real), parameter, public bet0
Meridional gradient of f (s^{-1} m^{-1})
Definition: qg_constants_mod.F90:36
qg_constants_mod::g
real(kind_real), parameter, public g
Gravity (m^2 s^{-2})
Definition: qg_constants_mod.F90:29
qg_constants_mod::heating_scale
real(kind_real), parameter, public heating_scale
Heating term scale (m)
Definition: qg_constants_mod.F90:41
qg_geom_mod::qg_geom_set_atlas_lonlat
subroutine, public qg_geom_set_atlas_lonlat(self, afieldset)
Set ATLAS lon/lat field.
Definition: qg_geom_mod.F90:217
qg_geom_mod::qg_geom
Definition: qg_geom_mod.F90:26
qg_geom_mod::qg_geom_registry
type(registry_t), public qg_geom_registry
Linked list interface - defines registry_t type.
Definition: qg_geom_mod.F90:53
qg_geom_mod::qg_geom_setup
subroutine, public qg_geom_setup(self, f_conf)
Linked list implementation.
Definition: qg_geom_mod.F90:64
qg_constants_mod
Definition: qg_constants_mod.F90:9
qg_geom_mod::qg_geom_fill_atlas_fieldset
subroutine, public qg_geom_fill_atlas_fieldset(self, afieldset)
Fill ATLAS fieldset.
Definition: qg_geom_mod.F90:245
qg_projection_mod::xy_to_lonlat
subroutine, public xy_to_lonlat(x, y, lon, lat, mapfac)
Convert x/y to lon/lat.
Definition: qg_projection_mod.F90:23
qg_constants_mod::dlogtheta
real(kind_real), parameter, public dlogtheta
Difference in log(pot. T) across boundary.
Definition: qg_constants_mod.F90:37
qg_geom_mod::qg_geom_delete
subroutine, public qg_geom_delete(self)
Delete geometry.
Definition: qg_geom_mod.F90:326
qg_constants_mod::f0
real(kind_real), parameter, public f0
Coriolis parameter at the center of the domain.
Definition: qg_constants_mod.F90:35
qg_constants_mod::domain_meridional
real(kind_real), parameter, public domain_meridional
Model domain in meridional direction (m)
Definition: qg_constants_mod.F90:33
qg_constants_mod::heating_amplitude
real(kind_real), parameter, public heating_amplitude
Heating term amplitude.
Definition: qg_constants_mod.F90:40