OOPS
qg_gom_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_gom_mod
10 
11 use atlas_module, only: atlas_field
12 use fckit_configuration_module, only: fckit_configuration
13 use fckit_log_module, only: fckit_log
14 use iso_c_binding
15 use kinds
16 use netcdf
19 use qg_geom_mod
20 use qg_locs_mod
22 use qg_tools_mod
23 use random_mod
24 
25 implicit none
26 private
27 public :: qg_gom
28 public :: qg_gom_registry
32 ! ------------------------------------------------------------------------------
33 type :: qg_gom
34  integer :: nobs !< Number of observations
35  real(kind_real), allocatable :: x(:) !< Streamfunction observations values
36  real(kind_real), allocatable :: q(:) !< Potential vorticity observations values
37  real(kind_real), allocatable :: u(:) !< Zonal wind observations values
38  real(kind_real), allocatable :: v(:) !< Meridian wind observations values
39  logical :: lalloc = .false. !< Allocation flag
40  type(oops_variables) :: vars !< Variables
41 end type qg_gom
42 
43 #define LISTED_TYPE qg_gom
44 
45 !> Linked list interface - defines registry_t type
46 #include "oops/util/linkedList_i.f"
47 
48 !> Global registry
49 type(registry_t) :: qg_gom_registry
50 ! ------------------------------------------------------------------------------
51 contains
52 ! ------------------------------------------------------------------------------
53 ! Public
54 ! ------------------------------------------------------------------------------
55 !> Linked list implementation
56 #include "oops/util/linkedList_c.f"
57 ! ------------------------------------------------------------------------------
58 !> Setup GOM
59 subroutine qg_gom_setup(self,nobs)
60 
61 implicit none
62 
63 ! Passed variables
64 type(qg_gom),intent(inout) :: self !< GOM
65 integer, intent(in) :: nobs !< Number of observations
66 
67 ! Set attributes
68 self%nobs = nobs
69 
70 ! Allocation
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))
75 self%lalloc = .true.
76 
77 end subroutine qg_gom_setup
78 ! ------------------------------------------------------------------------------
79 !> Delete GOM
80 subroutine qg_gom_delete(self)
81 
82 implicit none
83 
84 ! Passed variables
85 type(qg_gom),intent(inout) :: self !< GOM
86 
87 ! Release memory
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)
92 self%lalloc = .false.
93 
94 end subroutine qg_gom_delete
95 ! ------------------------------------------------------------------------------
96 !> Copy GOM
97 subroutine qg_gom_copy(self,other)
98 
99 implicit none
100 
101 ! Passed variables
102 type(qg_gom),intent(inout) :: self !< GOM
103 type(qg_gom),intent(in) :: other !< Other GOM
104 
105 ! Copy attributes
106 self%nobs = other%nobs
107 
108 ! Allocation
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))
114  self%lalloc = .true.
115 endif
116 
117 ! Copy
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
122 
123 end subroutine qg_gom_copy
124 ! ------------------------------------------------------------------------------
125 !> Set GOM to zero
126 subroutine qg_gom_zero(self)
127 
128 implicit none
129 
130 ! Passed variables
131 type(qg_gom),intent(inout) :: self !< GOM
132 
133 ! Set to zero
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
138 
139 end subroutine qg_gom_zero
140 ! ------------------------------------------------------------------------------
141 !> Get GOM absolute value
142 subroutine qg_gom_abs(self)
143 
144 implicit none
145 
146 ! Passed variables
147 type(qg_gom),intent(inout) :: self !< GOM
148 
149 ! Get absolute value
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)
154 
155 end subroutine qg_gom_abs
156 ! ------------------------------------------------------------------------------
157 !> Generate random GOM values
158 subroutine qg_gom_random(self)
159 
160 implicit none
161 
162 ! Passed variables
163 type(qg_gom),intent(inout) :: self !< GOM
164 
165 ! Local variables
166 integer :: nv
167 real(kind_real),allocatable :: values(:,:)
168 
169 ! TODO(Benjamin): change that in a following PR
170 nv = 0
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))
176 
177 ! Generate random GOM values
178 call normal_distribution(values,0.0_kind_real,1.0_kind_real)
179 
180 ! Split random values
181 nv = 0
182 if (self%vars%has('x')) then
183  nv = nv+1
184  self%x = values(nv,:)
185 endif
186 if (self%vars%has('q')) then
187  nv = nv+1
188  self%q = values(nv,:)
189 endif
190 if (self%vars%has('u')) then
191  nv = nv+1
192  self%u = values(nv,:)
193 endif
194 if (self%vars%has('v')) then
195  nv = nv+1
196  self%v = values(nv,:)
197 endif
198 
199 end subroutine qg_gom_random
200 ! ------------------------------------------------------------------------------
201 !> Multiply GOM with a scalar
202 subroutine qg_gom_mult(self,zz)
203 
204 implicit none
205 
206 ! Passed variables
207 type(qg_gom),intent(inout) :: self !< GOM
208 real(kind_real),intent(in) :: zz !< Multiplier
209 
210 ! Multiply GOM with a scalar
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
215 
216 end subroutine qg_gom_mult
217 ! ------------------------------------------------------------------------------
218 !> Add GOM
219 subroutine qg_gom_add(self,other)
220 
221 implicit none
222 
223 ! Passed variables
224 type(qg_gom),intent(inout) :: self !< GOM
225 type(qg_gom),intent(in) :: other !< Other GOM
226 
227 ! Add GOM
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
232 
233 end subroutine qg_gom_add
234 ! ------------------------------------------------------------------------------
235 !> Subtract GOM
236 subroutine qg_gom_diff(self,other)
237 
238 implicit none
239 
240 ! Passed variables
241 type(qg_gom),intent(inout) :: self !< GOM
242 type(qg_gom),intent(in) :: other !< Other GOM
243 
244 ! Subtract GOM
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
249 
250 end subroutine qg_gom_diff
251 ! ------------------------------------------------------------------------------
252 !> Schur product for GOM
253 subroutine qg_gom_schurmult(self,other)
254 
255 implicit none
256 
257 ! Passed variables
258 type(qg_gom),intent(inout) :: self !< GOM
259 type(qg_gom),intent(in) :: other !< Other GOM
260 
261 ! Multiply GOM
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
266 
267 end subroutine qg_gom_schurmult
268 ! ------------------------------------------------------------------------------
269 !> Schur division for GOM
270 subroutine qg_gom_divide(self,other)
271 
272 implicit none
273 
274 ! Passed variables
275 type(qg_gom),intent(inout) :: self !< GOM
276 type(qg_gom),intent(in) :: other !< Other GOM
277 
278 ! Local variables
279 real(kind_real) :: tol
280 integer :: jloc
281 
282 ! Set tolerance
283 tol = epsilon(tol)
284 
285 ! Conditional division
286 do jloc=1,self%nobs
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)
290  else
291  self%x(jloc) = 0.0
292  endif
293  endif
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)
297  else
298  self%q(jloc) = 0.0
299  endif
300  endif
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)
304  else
305  self%u(jloc) = 0.0
306  endif
307  endif
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)
311  else
312  self%v(jloc) = 0.0
313  endif
314  endif
315 enddo
316 
317 end subroutine qg_gom_divide
318 ! ------------------------------------------------------------------------------
319 !> Compute GOM RMS
320 subroutine qg_gom_rms(self,rms)
321 
322 implicit none
323 
324 ! Passed variables
325 type(qg_gom),intent(inout) :: self !< GOM
326 real(kind_real),intent(inout) :: rms !< RMS
327 
328 ! Local variables
329 integer :: nv
330 
331 ! Initialization
332 rms = 0.0
333 nv = 0
334 
335 ! Loop over values
336 if (self%vars%has('x')) then
337  rms = rms+sum(self%x**2)
338  nv = nv+1
339 endif
340 if (self%vars%has('q')) then
341  rms = rms+sum(self%q**2)
342  nv = nv+1
343 endif
344 if (self%vars%has('u')) then
345  rms = rms+sum(self%u**2)
346  nv = nv+1
347 endif
348 if (self%vars%has('v')) then
349  rms = rms+sum(self%v**2)
350  nv = nv+1
351 endif
352 
353 ! Normalize and take square-root
354 rms = sqrt(rms/real(self%nobs*nv,kind_real))
355 
356 end subroutine qg_gom_rms
357 ! ------------------------------------------------------------------------------
358 !> GOM dot product
359 subroutine qg_gom_dotprod(gom1,gom2,prod)
360 
361 implicit none
362 
363 ! Passed variables
364 type(qg_gom),intent(in) :: gom1 !< GOM 1
365 type(qg_gom),intent(in) :: gom2 !< GOM 2
366 real(kind_real),intent(inout) :: prod !< Dot product
367 
368 ! Local variables
369 integer :: jo,jv
370 
371 ! Check
372 if (gom1%nobs/=gom2%nobs) call abor1_ftn('qg_gom_dotprod: inconsistent GOM sizes')
373 
374 ! Initialization
375 prod = 0.0
376 
377 ! Dot product
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)
382 
383 end subroutine qg_gom_dotprod
384 ! ------------------------------------------------------------------------------
385 !> Compute GOM stats
386 subroutine qg_gom_stats(self,kobs,pmin,pmax,prms)
387 
388 implicit none
389 
390 ! Passed variables
391 type(qg_gom),intent(inout) :: self !< GOM
392 integer,intent(inout) :: kobs !< Number of observations
393 real(kind_real),intent(inout) :: pmin !< Minimum value
394 real(kind_real),intent(inout) :: pmax !< Maximum value
395 real(kind_real),intent(inout) :: prms !< RMS
396 
397 ! Local variables
398 integer :: nv
399 
400 ! Compute GOM stats
401 kobs = self%nobs
402 if (self%nobs>0) then
403  pmin = huge(1.0)
404  pmax = -huge(1.0)
405  prms = 0.0
406  nv = 0
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)
411  nv = nv+1
412  endif
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)
417  nv = nv+1
418  endif
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)
423  nv = nv+1
424  endif
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)
429  nv = nv+1
430  endif
431  prms = sqrt(prms/real(self%nobs*nv,kind_real))
432 else
433  pmin = 0.0
434  pmax = 0.0
435  prms = 0.0
436 end if
437 
438 end subroutine qg_gom_stats
439 ! ------------------------------------------------------------------------------
440 !> Find and locate GOM max. value
441 subroutine qg_gom_maxloc(self,mxval,mxloc,mxvar)
442 
443 implicit none
444 
445 ! Passed variables
446 type(qg_gom),intent(inout) :: self !< GOM
447 real(kind_real),intent(inout) :: mxval !< Maximum value
448 integer,intent(inout) :: mxloc !< Location of maximum value
449 type(oops_variables),intent(inout) :: mxvar !< Variable of maximum value
450 
451 ! Local variables
452 integer :: mxloc_arr(1),mxval_tmp
453 character(len=1) :: var
454 
455 ! Initialization
456 mxval = -huge(1.0)
457 
458 ! Find GOM max. value
459 if (self%vars%has('x')) then
460  mxval_tmp = maxval(self%x)
461  if (mxval_tmp>mxval) then
462  mxval = mxval
463  mxloc_arr = maxloc(self%x)
464  var = 'x'
465  endif
466 endif
467 if (self%vars%has('q')) then
468  mxval_tmp = maxval(self%q)
469  if (mxval_tmp>mxval) then
470  mxval = mxval
471  mxloc_arr = maxloc(self%q)
472  var = 'q'
473  endif
474 endif
475 if (self%vars%has('u')) then
476  mxval_tmp = maxval(self%u)
477  if (mxval_tmp>mxval) then
478  mxval = mxval
479  mxloc_arr = maxloc(self%u)
480  var = 'u'
481  endif
482 endif
483 if (self%vars%has('v')) then
484  mxval_tmp = maxval(self%v)
485  if (mxval_tmp>mxval) then
486  mxval = mxval
487  mxloc_arr = maxloc(self%v)
488  var = 'v'
489  endif
490 endif
491 
492 ! Locate GOM max. value
493 mxloc = mxloc_arr(1)
494 
495 ! Set GOM max. variable
496 call mxvar%push_back(var)
497 
498 end subroutine qg_gom_maxloc
499 ! ------------------------------------------------------------------------------
500 !> Read GOM from file
501 subroutine qg_gom_read_file(self,f_conf)
502 
503 implicit none
504 
505 ! Passed variables
506 type(qg_gom),intent(inout) :: self !< GOM
507 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
508 
509 ! Local variables
510 integer :: ncid,nobs_id,nobs,x_id,q_id,u_id,v_id
511 character(len=1024) :: filename
512 character(len=:),allocatable :: str
513 
514 ! Get filename
515 call f_conf%get_or_die("filename",str)
516 filename = str
517 call fckit_log%info('qg_gom_read_file: reading '//trim(filename))
518 
519 ! Open NetCDF file
520 call ncerr(nf90_open(trim(filename)//'.nc',nf90_nowrite,ncid))
521 
522 ! Get dimension id
523 call ncerr(nf90_inq_dimid(ncid,'nobs',nobs_id))
524 
525 ! Get dimension
526 call ncerr(nf90_inquire_dimension(ncid,nobs_id,len=nobs))
527 
528 ! GOM setup
529 call qg_gom_setup(self,nobs)
530 
531 ! Get variables ids
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))
536 
537 ! Get variables
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))
542 
543 ! Close NetCDF file
544 call ncerr(nf90_close(ncid))
545 
546 end subroutine qg_gom_read_file
547 ! ------------------------------------------------------------------------------
548 !> Write GOM to file
549 subroutine qg_gom_write_file(self,f_conf)
550 
551 implicit none
552 
553 ! Passed variables
554 type(qg_gom),intent(inout) :: self !< GOM
555 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
556 
557 ! Local variables
558 integer :: ncid,nobs_id,x_id,q_id,u_id,v_id
559 character(len=1024) :: filename
560 character(len=:),allocatable :: str
561 
562 ! Check allocation
563 if (.not.self%lalloc) call abor1_ftn('qg_gom_write_file: gom not allocated')
564 
565 ! Set filename
566 call f_conf%get_or_die("filename",str)
567 filename = str
568 call fckit_log%info('qg_gom_write_file: writing '//trim(filename))
569 
570 ! Create NetCDF file
571 call ncerr(nf90_create(trim(filename)//'.nc',or(nf90_clobber,nf90_64bit_offset),ncid))
572 
573 ! Define dimensions
574 call ncerr(nf90_def_dim(ncid,'nobs',self%nobs,nobs_id))
575 
576 ! Define variables
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))
581 
582 ! End definitions
583 call ncerr(nf90_enddef(ncid))
584 
585 ! Put variables
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))
590 
591 ! Close NetCDF file
592 call ncerr(nf90_close(ncid))
593 
594 end subroutine qg_gom_write_file
595 ! ------------------------------------------------------------------------------
596 !> GOM analytic initialization
597 subroutine qg_gom_analytic_init(self,locs,f_conf)
598 
599 implicit none
600 
601 ! Passed variables
602 type(qg_gom),intent(inout) :: self !< GOM
603 type(qg_locs),intent(inout) :: locs !< Locations
604 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
605 
606 ! Local variables
607 integer :: iloc
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
613 
614 ! get locations
615 lonlat_field = locs%lonlat()
616 call lonlat_field%data(lonlat)
617 
618 z_field = locs%altitude()
619 call z_field%data(z)
620 
621 ! Check allocation
622 if (.not.self%lalloc) call abor1_ftn('qg_gom_analytic init: gom not allocated')
623 
624 ! Get analytic configuration
625 call f_conf%get_or_die("analytic_init",str)
626 ic = 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')
631  ! Go to cartesian coordinates
632  call lonlat_to_xy(lonlat(1,iloc),lonlat(2,iloc),x,y)
633 
634  ! Compute values for baroclinic instability
635  if (self%vars%has('x')) call baroclinic_instability(x,y,z(iloc),'x',self%x(iloc))
636  if (self%vars%has('q')) call baroclinic_instability(x,y,z(iloc),'q',self%q(iloc))
637  if (self%vars%has('u')) call baroclinic_instability(x,y,z(iloc),'u',self%u(iloc))
638  if (self%vars%has('v')) call baroclinic_instability(x,y,z(iloc),'v',self%v(iloc))
639  case ('large-vortices')
640  ! Go to cartesian coordinates
641  call lonlat_to_xy(lonlat(1,iloc),lonlat(2,iloc),x,y)
642 
643  ! Compute values for 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))
648  case default
649  call abor1_ftn('qg_gom_analytic_init: unknown initialization')
650  endselect
651 enddo
652 
653 call lonlat_field%final()
654 call z_field%final()
655 
656 end subroutine qg_gom_analytic_init
657 ! ------------------------------------------------------------------------------
658 end module qg_gom_mod
Fortran interface to Variables.
subroutine, public qg_gom_rms(self, rms)
Compute GOM RMS.
Definition: qg_gom_mod.F90:321
subroutine, public qg_gom_write_file(self, f_conf)
Write GOM to file.
Definition: qg_gom_mod.F90:550
subroutine, public qg_gom_add(self, other)
Add GOM.
Definition: qg_gom_mod.F90:220
type(registry_t), public qg_gom_registry
Linked list interface - defines registry_t type.
Definition: qg_gom_mod.F90:49
subroutine, public qg_gom_schurmult(self, other)
Schur product for GOM.
Definition: qg_gom_mod.F90:254
subroutine, public qg_gom_dotprod(gom1, gom2, prod)
GOM dot product.
Definition: qg_gom_mod.F90:360
subroutine, public qg_gom_random(self)
Generate random GOM values.
Definition: qg_gom_mod.F90:159
subroutine, public qg_gom_setup(self, nobs)
Linked list implementation.
Definition: qg_gom_mod.F90:60
subroutine, public qg_gom_maxloc(self, mxval, mxloc, mxvar)
Find and locate GOM max. value.
Definition: qg_gom_mod.F90:442
subroutine, public qg_gom_diff(self, other)
Subtract GOM.
Definition: qg_gom_mod.F90:237
subroutine, public qg_gom_abs(self)
Get GOM absolute value.
Definition: qg_gom_mod.F90:143
subroutine, public qg_gom_mult(self, zz)
Multiply GOM with a scalar.
Definition: qg_gom_mod.F90:203
subroutine, public qg_gom_read_file(self, f_conf)
Read GOM from file.
Definition: qg_gom_mod.F90:502
subroutine, public qg_gom_delete(self)
Delete GOM.
Definition: qg_gom_mod.F90:81
subroutine, public qg_gom_analytic_init(self, locs, f_conf)
GOM analytic initialization.
Definition: qg_gom_mod.F90:598
subroutine, public qg_gom_divide(self, other)
Schur division for GOM.
Definition: qg_gom_mod.F90:271
subroutine, public qg_gom_copy(self, other)
Copy GOM.
Definition: qg_gom_mod.F90:98
subroutine, public qg_gom_zero(self)
Set GOM to zero.
Definition: qg_gom_mod.F90:127
subroutine, public qg_gom_stats(self, kobs, pmin, pmax, prms)
Compute GOM stats.
Definition: qg_gom_mod.F90:387
subroutine, public lonlat_to_xy(lon, lat, x, y)
Convert lon/lat to x/y.
subroutine, public large_vortices(x, y, z, var, res)
Generate values for large vortices.
subroutine, public baroclinic_instability(x, y, z, var, res)
Generate values for baroclinic instability.
subroutine, public ncerr(info)
Check NetCDF status.