11 use atlas_module,
only: atlas_field
12 use fckit_configuration_module,
only: fckit_configuration
13 use fckit_log_module,
only: fckit_log
35 real(kind_real),
allocatable :: x(:)
36 real(kind_real),
allocatable :: q(:)
37 real(kind_real),
allocatable :: u(:)
38 real(kind_real),
allocatable :: v(:)
39 logical :: lalloc = .false.
43 #define LISTED_TYPE qg_gom
46 #include "oops/util/linkedList_i.f"
56 #include "oops/util/linkedList_c.f"
64 type(
qg_gom),
intent(inout) :: self
65 integer,
intent(in) :: nobs
71 if (self%vars%has(
'x'))
allocate(self%x(self%nobs))
72 if (self%vars%has(
'q'))
allocate(self%q(self%nobs))
73 if (self%vars%has(
'u'))
allocate(self%u(self%nobs))
74 if (self%vars%has(
'v'))
allocate(self%v(self%nobs))
85 type(
qg_gom),
intent(inout) :: self
88 if (
allocated(self%x))
deallocate(self%x)
89 if (
allocated(self%q))
deallocate(self%q)
90 if (
allocated(self%u))
deallocate(self%u)
91 if (
allocated(self%v))
deallocate(self%v)
102 type(
qg_gom),
intent(inout) :: self
103 type(
qg_gom),
intent(in) :: other
106 self%nobs = other%nobs
109 if (.not.self%lalloc)
then
110 if (self%vars%has(
'x'))
allocate(self%x(self%nobs))
111 if (self%vars%has(
'q'))
allocate(self%q(self%nobs))
112 if (self%vars%has(
'u'))
allocate(self%u(self%nobs))
113 if (self%vars%has(
'v'))
allocate(self%v(self%nobs))
118 if (self%vars%has(
'x')) self%x = other%x
119 if (self%vars%has(
'q')) self%q = other%q
120 if (self%vars%has(
'u')) self%u = other%u
121 if (self%vars%has(
'v')) self%v = other%v
131 type(
qg_gom),
intent(inout) :: self
134 if (self%vars%has(
'x')) self%x = 0.0
135 if (self%vars%has(
'q')) self%q = 0.0
136 if (self%vars%has(
'u')) self%u = 0.0
137 if (self%vars%has(
'v')) self%v = 0.0
147 type(
qg_gom),
intent(inout) :: self
150 if (self%vars%has(
'x')) self%x = abs(self%x)
151 if (self%vars%has(
'q')) self%q = abs(self%q)
152 if (self%vars%has(
'u')) self%u = abs(self%u)
153 if (self%vars%has(
'v')) self%v = abs(self%v)
163 type(
qg_gom),
intent(inout) :: self
167 real(kind_real),
allocatable :: values(:,:)
171 if (self%vars%has(
'x')) nv = nv+1
172 if (self%vars%has(
'q')) nv = nv+1
173 if (self%vars%has(
'u')) nv = nv+1
174 if (self%vars%has(
'v')) nv = nv+1
175 allocate(values(nv,self%nobs))
178 call normal_distribution(values,0.0_kind_real,1.0_kind_real)
182 if (self%vars%has(
'x'))
then
184 self%x = values(nv,:)
186 if (self%vars%has(
'q'))
then
188 self%q = values(nv,:)
190 if (self%vars%has(
'u'))
then
192 self%u = values(nv,:)
194 if (self%vars%has(
'v'))
then
196 self%v = values(nv,:)
207 type(
qg_gom),
intent(inout) :: self
208 real(kind_real),
intent(in) :: zz
211 if (self%vars%has(
'x')) self%x = zz*self%x
212 if (self%vars%has(
'q')) self%q = zz*self%q
213 if (self%vars%has(
'u')) self%u = zz*self%u
214 if (self%vars%has(
'v')) self%v = zz*self%v
224 type(
qg_gom),
intent(inout) :: self
225 type(
qg_gom),
intent(in) :: other
228 if (self%vars%has(
'x')) self%x = self%x+other%x
229 if (self%vars%has(
'q')) self%q = self%q+other%q
230 if (self%vars%has(
'u')) self%u = self%u+other%u
231 if (self%vars%has(
'v')) self%v = self%v+other%v
241 type(
qg_gom),
intent(inout) :: self
242 type(
qg_gom),
intent(in) :: other
245 if (self%vars%has(
'x')) self%x = self%x-other%x
246 if (self%vars%has(
'q')) self%q = self%q-other%q
247 if (self%vars%has(
'u')) self%u = self%u-other%u
248 if (self%vars%has(
'v')) self%v = self%v-other%v
258 type(
qg_gom),
intent(inout) :: self
259 type(
qg_gom),
intent(in) :: other
262 if (self%vars%has(
'x')) self%x = self%x*other%x
263 if (self%vars%has(
'q')) self%q = self%q*other%q
264 if (self%vars%has(
'u')) self%u = self%u*other%u
265 if (self%vars%has(
'v')) self%v = self%v*other%v
275 type(
qg_gom),
intent(inout) :: self
276 type(
qg_gom),
intent(in) :: other
279 real(kind_real) :: tol
287 if (self%vars%has(
'x'))
then
288 if (abs(other%x(jloc))>tol)
then
289 self%x(jloc) = self%x(jloc)/other%x(jloc)
294 if (self%vars%has(
'q'))
then
295 if (abs(other%q(jloc))>tol)
then
296 self%q(jloc) = self%q(jloc)/other%q(jloc)
301 if (self%vars%has(
'u'))
then
302 if (abs(other%u(jloc))>tol)
then
303 self%u(jloc) = self%u(jloc)/other%u(jloc)
308 if (self%vars%has(
'v'))
then
309 if (abs(other%v(jloc))>tol)
then
310 self%v(jloc) = self%v(jloc)/other%v(jloc)
325 type(
qg_gom),
intent(inout) :: self
326 real(kind_real),
intent(inout) :: rms
336 if (self%vars%has(
'x'))
then
337 rms = rms+sum(self%x**2)
340 if (self%vars%has(
'q'))
then
341 rms = rms+sum(self%q**2)
344 if (self%vars%has(
'u'))
then
345 rms = rms+sum(self%u**2)
348 if (self%vars%has(
'v'))
then
349 rms = rms+sum(self%v**2)
354 rms = sqrt(rms/real(self%nobs*nv,kind_real))
364 type(
qg_gom),
intent(in) :: gom1
365 type(
qg_gom),
intent(in) :: gom2
366 real(kind_real),
intent(inout) :: prod
372 if (gom1%nobs/=gom2%nobs)
call abor1_ftn(
'qg_gom_dotprod: inconsistent GOM sizes')
378 if (gom1%vars%has(
'x').and.gom2%vars%has(
'x')) prod = prod+sum(gom1%x*gom2%x)
379 if (gom1%vars%has(
'q').and.gom2%vars%has(
'q')) prod = prod+sum(gom1%q*gom2%q)
380 if (gom1%vars%has(
'u').and.gom2%vars%has(
'u')) prod = prod+sum(gom1%u*gom2%u)
381 if (gom1%vars%has(
'v').and.gom2%vars%has(
'v')) prod = prod+sum(gom1%v*gom2%v)
391 type(
qg_gom),
intent(inout) :: self
392 integer,
intent(inout) :: kobs
393 real(kind_real),
intent(inout) :: pmin
394 real(kind_real),
intent(inout) :: pmax
395 real(kind_real),
intent(inout) :: prms
402 if (self%nobs>0)
then
407 if (self%vars%has(
'x'))
then
408 pmin = min(pmin,minval(self%x))
409 pmax = max(pmax,maxval(self%x))
410 prms = prms+sum(self%x**2)
413 if (self%vars%has(
'q'))
then
414 pmin = min(pmin,minval(self%q))
415 pmax = max(pmax,maxval(self%q))
416 prms = prms+sum(self%q**2)
419 if (self%vars%has(
'u'))
then
420 pmin = min(pmin,minval(self%u))
421 pmax = max(pmax,maxval(self%u))
422 prms = prms+sum(self%u**2)
425 if (self%vars%has(
'v'))
then
426 pmin = min(pmin,minval(self%v))
427 pmax = max(pmax,maxval(self%v))
428 prms = prms+sum(self%v**2)
431 prms = sqrt(prms/real(self%nobs*nv,kind_real))
446 type(
qg_gom),
intent(inout) :: self
447 real(kind_real),
intent(inout) :: mxval
448 integer,
intent(inout) :: mxloc
452 integer :: mxloc_arr(1),mxval_tmp
453 character(len=1) :: var
459 if (self%vars%has(
'x'))
then
460 mxval_tmp = maxval(self%x)
461 if (mxval_tmp>mxval)
then
463 mxloc_arr = maxloc(self%x)
467 if (self%vars%has(
'q'))
then
468 mxval_tmp = maxval(self%q)
469 if (mxval_tmp>mxval)
then
471 mxloc_arr = maxloc(self%q)
475 if (self%vars%has(
'u'))
then
476 mxval_tmp = maxval(self%u)
477 if (mxval_tmp>mxval)
then
479 mxloc_arr = maxloc(self%u)
483 if (self%vars%has(
'v'))
then
484 mxval_tmp = maxval(self%v)
485 if (mxval_tmp>mxval)
then
487 mxloc_arr = maxloc(self%v)
496 call mxvar%push_back(var)
506 type(
qg_gom),
intent(inout) :: self
507 type(fckit_configuration),
intent(in) :: f_conf
510 integer :: ncid,nobs_id,nobs,x_id,q_id,u_id,v_id
511 character(len=1024) :: filename
512 character(len=:),
allocatable :: str
515 call f_conf%get_or_die(
"filename",str)
517 call fckit_log%info(
'qg_gom_read_file: reading '//trim(filename))
520 call ncerr(nf90_open(trim(filename)//
'.nc',nf90_nowrite,ncid))
523 call ncerr(nf90_inq_dimid(ncid,
'nobs',nobs_id))
526 call ncerr(nf90_inquire_dimension(ncid,nobs_id,len=nobs))
532 if (self%vars%has(
'x'))
call ncerr(nf90_inq_varid(ncid,
'x',x_id))
533 if (self%vars%has(
'q'))
call ncerr(nf90_inq_varid(ncid,
'q',q_id))
534 if (self%vars%has(
'u'))
call ncerr(nf90_inq_varid(ncid,
'u',u_id))
535 if (self%vars%has(
'v'))
call ncerr(nf90_inq_varid(ncid,
'v',v_id))
538 if (self%vars%has(
'x'))
call ncerr(nf90_get_var(ncid,x_id,self%x))
539 if (self%vars%has(
'q'))
call ncerr(nf90_get_var(ncid,q_id,self%q))
540 if (self%vars%has(
'u'))
call ncerr(nf90_get_var(ncid,u_id,self%u))
541 if (self%vars%has(
'v'))
call ncerr(nf90_get_var(ncid,v_id,self%v))
544 call ncerr(nf90_close(ncid))
554 type(
qg_gom),
intent(inout) :: self
555 type(fckit_configuration),
intent(in) :: f_conf
558 integer :: ncid,nobs_id,x_id,q_id,u_id,v_id
559 character(len=1024) :: filename
560 character(len=:),
allocatable :: str
563 if (.not.self%lalloc)
call abor1_ftn(
'qg_gom_write_file: gom not allocated')
566 call f_conf%get_or_die(
"filename",str)
568 call fckit_log%info(
'qg_gom_write_file: writing '//trim(filename))
571 call ncerr(nf90_create(trim(filename)//
'.nc',or(nf90_clobber,nf90_64bit_offset),ncid))
574 call ncerr(nf90_def_dim(ncid,
'nobs',self%nobs,nobs_id))
577 if (self%vars%has(
'x'))
call ncerr(nf90_def_var(ncid,
'x',nf90_double,(/nobs_id/),x_id))
578 if (self%vars%has(
'q'))
call ncerr(nf90_def_var(ncid,
'q',nf90_double,(/nobs_id/),q_id))
579 if (self%vars%has(
'u'))
call ncerr(nf90_def_var(ncid,
'u',nf90_double,(/nobs_id/),u_id))
580 if (self%vars%has(
'v'))
call ncerr(nf90_def_var(ncid,
'v',nf90_double,(/nobs_id/),v_id))
583 call ncerr(nf90_enddef(ncid))
586 if (self%vars%has(
'x'))
call ncerr(nf90_put_var(ncid,x_id,self%x))
587 if (self%vars%has(
'q'))
call ncerr(nf90_put_var(ncid,q_id,self%q))
588 if (self%vars%has(
'u'))
call ncerr(nf90_put_var(ncid,u_id,self%u))
589 if (self%vars%has(
'v'))
call ncerr(nf90_put_var(ncid,v_id,self%v))
592 call ncerr(nf90_close(ncid))
602 type(
qg_gom),
intent(inout) :: self
603 type(
qg_locs),
intent(inout) :: locs
604 type(fckit_configuration),
intent(in) :: f_conf
608 real(kind_real) :: x,y
609 character(len=30) :: ic
610 character(len=:),
allocatable :: str
611 real(kind_real),
pointer :: lonlat(:,:), z(:)
612 type(atlas_field) :: lonlat_field, z_field
615 lonlat_field = locs%lonlat()
616 call lonlat_field%data(lonlat)
618 z_field = locs%altitude()
622 if (.not.self%lalloc)
call abor1_ftn(
'qg_gom_analytic init: gom not allocated')
625 call f_conf%get_or_die(
"analytic_init",str)
627 call fckit_log%info(
'qg_gom_analytic_init: ic = '//trim(ic))
628 do iloc=1,locs%nlocs()
629 select case (trim(ic))
630 case (
'baroclinic-instability')
639 case (
'large-vortices')
644 if (self%vars%has(
'x'))
call large_vortices(x,y,z(iloc),
'x',self%x(iloc))
645 if (self%vars%has(
'q'))
call large_vortices(x,y,z(iloc),
'q',self%q(iloc))
646 if (self%vars%has(
'u'))
call large_vortices(x,y,z(iloc),
'u',self%u(iloc))
647 if (self%vars%has(
'v'))
call large_vortices(x,y,z(iloc),
'v',self%v(iloc))
649 call abor1_ftn(
'qg_gom_analytic_init: unknown initialization')
653 call lonlat_field%final()
Fortran interface to Variables.
subroutine, public qg_gom_rms(self, rms)
Compute GOM RMS.
subroutine, public qg_gom_write_file(self, f_conf)
Write GOM to file.
subroutine, public qg_gom_add(self, other)
Add GOM.
type(registry_t), public qg_gom_registry
Linked list interface - defines registry_t type.
subroutine, public qg_gom_schurmult(self, other)
Schur product for GOM.
subroutine, public qg_gom_dotprod(gom1, gom2, prod)
GOM dot product.
subroutine, public qg_gom_random(self)
Generate random GOM values.
subroutine, public qg_gom_setup(self, nobs)
Linked list implementation.
subroutine, public qg_gom_maxloc(self, mxval, mxloc, mxvar)
Find and locate GOM max. value.
subroutine, public qg_gom_diff(self, other)
Subtract GOM.
subroutine, public qg_gom_abs(self)
Get GOM absolute value.
subroutine, public qg_gom_mult(self, zz)
Multiply GOM with a scalar.
subroutine, public qg_gom_read_file(self, f_conf)
Read GOM from file.
subroutine, public qg_gom_delete(self)
Delete GOM.
subroutine, public qg_gom_analytic_init(self, locs, f_conf)
GOM analytic initialization.
subroutine, public qg_gom_divide(self, other)
Schur division for GOM.
subroutine, public qg_gom_copy(self, other)
Copy GOM.
subroutine, public qg_gom_zero(self)
Set GOM to zero.
subroutine, public qg_gom_stats(self, kobs, pmin, pmax, prms)
Compute GOM stats.
subroutine, public lonlat_to_xy(lon, lat, x, y)
Convert lon/lat to x/y.