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
30 real(kind_real) :: deltax
31 real(kind_real) :: deltay
32 real(kind_real),
allocatable :: x(:)
33 real(kind_real),
allocatable :: y(:)
34 real(kind_real),
allocatable :: z(:)
35 real(kind_real),
allocatable :: lat(:,:)
36 real(kind_real),
allocatable :: lon(:,:)
37 real(kind_real),
allocatable :: area(:,:)
38 real(kind_real),
allocatable :: f(:,:)
39 real(kind_real),
allocatable :: f_p(:,:)
40 real(kind_real),
allocatable :: f_pinv(:,:)
41 real(kind_real),
allocatable :: f_d(:)
42 real(kind_real),
allocatable :: bet(:)
43 real(kind_real),
allocatable :: heat(:,:)
44 type(atlas_functionspace_pointcloud) :: afunctionspace
47 #define LISTED_TYPE qg_geom
50 #include "oops/util/linkedList_i.f"
60 #include "oops/util/linkedList_c.f"
66 type(
qg_geom),
intent(inout) :: self
67 type(fckit_configuration),
intent(in) :: f_conf
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
78 character(len=:),
allocatable :: str
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")
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))
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))
108 call f_conf%get_or_die(
"depths",real_array)
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)
123 self%x(ix) = (real(ix,kind_real)-0.5)*self%deltax
126 self%y(iy) = real(iy,kind_real)*self%deltay
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
138 self%z(1) = 0.5*depths(1)
140 self%z(iz) = sum(depths(1:iz-1))+0.5*depths(iz)
149 self%f(iz,iz) = self%f(iz,iz)-f
153 self%f(iz,iz) = self%f(iz,iz)-f
158 norm=maxval(abs(self%f))
159 fsave = self%f / norm
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')
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')
173 call dgetrf(self%nz,self%nz,vrlu,self%nz,ipiv,info)
174 if (info/=0)
call abor1_ftn(
'error in dgetrf')
178 call dgetri(self%nz,vrlusave,self%nz,ipivsave,work,-1,info)
179 if (info/=0)
call abor1_ftn(
'error in dgetri, first pass')
182 allocate(work(lwork))
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')
190 self%f_d = self%f_d * norm
191 self%f_p = self%f_p * norm
192 self%f_pinv = self%f_pinv / norm
196 self%bet(iy) = real(iy-(self%ny+1)/2,kind_real)*self%deltay*
bet0
200 call f_conf%get_or_die(
"heating",htype)
201 if (.not. htype)
then
203 call fckit_log%info(
'qg_geom_setup: heating off')
206 call fckit_log%info(
'qg_geom_setup: Gaussian heating on')
212 distx = abs(self%x(ix)-self%x(ix_c))
214 disty = abs(self%y(iy)-self%y(iy_c))
226 type(
qg_geom),
intent(inout) :: self
227 type(atlas_fieldset),
intent(inout) :: afieldset
230 integer :: ix,iy,inode
231 real(kind_real),
pointer :: real_ptr(:,:)
232 type(atlas_field) :: afield
235 afield = atlas_field(name=
"lonlat", kind=atlas_real(kind_real), shape=(/2,self%nx*self%ny/))
236 call afield%data(real_ptr)
241 real_ptr(1,inode) = self%lon(ix,iy)
242 real_ptr(2,inode) = self%lat(ix,iy)
245 call afieldset%add(afield)
254 type(
qg_geom),
intent(inout) :: self
255 type(atlas_fieldset),
intent(inout) :: afieldset
258 integer :: ix,iy,iz,inode
259 real(kind_real),
pointer :: real_ptr_1(:),real_ptr_2(:,:)
260 type(atlas_field) :: afield
263 afield = self%afunctionspace%create_field(name=
'area',kind=atlas_real(kind_real),levels=0)
264 call afield%data(real_ptr_1)
269 real_ptr_1(inode) = self%area(ix,iy)
272 call afieldset%add(afield)
276 afield = self%afunctionspace%create_field(name=
'vunit',kind=atlas_real(kind_real),levels=self%nz)
277 call afield%data(real_ptr_2)
279 real_ptr_2(iz,1:self%nx*self%ny) = self%z(iz)
281 call afieldset%add(afield)
290 type(
qg_geom),
intent(inout) :: self
291 type(
qg_geom),
intent(in) :: other
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))
313 self%deltax = other%deltax
314 self%deltay = other%deltay
320 self%area = other%area
323 self%f_pinv = other%f_pinv
326 self%heat = other%heat
327 self%afunctionspace = atlas_functionspace_pointcloud(other%afunctionspace%c_ptr())
335 type(
qg_geom),
intent(inout) :: self
343 deallocate(self%area)
346 deallocate(self%f_pinv)
349 deallocate(self%heat)
350 call self%afunctionspace%final()
358 type(
qg_geom),
intent(in) :: self
359 integer,
intent(out) :: nx
360 integer,
intent(out) :: ny
361 integer,
intent(out) :: nz
362 real(kind_real),
intent(out) :: deltax
363 real(kind_real),
intent(out) :: deltay
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.
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.
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.