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