OOPS
qg_fields_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 
10 
11 use atlas_module, only: atlas_field, atlas_fieldset, atlas_real
12 use fckit_configuration_module, only: fckit_configuration
13 use datetime_mod
14 use duration_mod
15 use fckit_log_module,only: fckit_log
16 use iso_c_binding
17 use kinds
18 use missing_values_mod
19 use netcdf
20 !$ use omp_lib
27 use qg_geom_mod
29 use qg_interp_mod
30 use qg_locs_mod
31 use qg_tools_mod
32 use random_mod
33 
34 implicit none
35 
36 private
37 public :: qg_fields
38 public :: qg_fields_registry
47 ! ------------------------------------------------------------------------------
48 integer,parameter :: rseed = 7 !< Random seed (for reproducibility)
49 
50 type :: qg_fields
51  type(qg_geom),pointer :: geom !< Geometry
52  logical :: lbc !< Boundaries are present
53  real(kind_real),allocatable :: x(:,:,:) !< Streamfunction
54  real(kind_real),allocatable :: q(:,:,:) !< Potential vorticity
55  real(kind_real),allocatable :: u(:,:,:) !< U wind
56  real(kind_real),allocatable :: v(:,:,:) !< V wind
57  real(kind_real),allocatable :: x_north(:) !< Streamfunction on northern wall
58  real(kind_real),allocatable :: x_south(:) !< Streamfunction on southern wall
59  real(kind_real),allocatable :: q_north(:,:) !< q on northern wall
60  real(kind_real),allocatable :: q_south(:,:) !< q on southern wall
61 end type qg_fields
62 
63 #define LISTED_TYPE qg_fields
64 
65 !> Linked list interface - defines registry_t type
66 #include "oops/util/linkedList_i.f"
67 
68 !> Global registry
69 type(registry_t) :: qg_fields_registry
70 ! ------------------------------------------------------------------------------
71 contains
72 ! ------------------------------------------------------------------------------
73 ! Public
74 ! ------------------------------------------------------------------------------
75 !> Linked list implementation
76 #include "oops/util/linkedList_c.f"
77 ! ------------------------------------------------------------------------------
78 !> Create fields from geometry and variables
79 subroutine qg_fields_create(self,geom,vars,lbc)
80 
81 implicit none
82 
83 ! Passed variables
84 type(qg_fields),intent(inout) :: self !< Fields
85 type(qg_geom),target,intent(in) :: geom !< Geometry
86 type(oops_variables),intent(in) :: vars !< List of variables
87 logical,intent(in) :: lbc !< Boundaries flag
88 
89 ! Local variables
90 character(len=1024) :: record
91 
92 ! Associate geometry
93 self%geom => geom
94 
95 ! Set boundaries
96 self%lbc = lbc
97 
98 ! Allocate 3d fields
99 if (vars%has('x')) allocate(self%x(self%geom%nx,self%geom%ny,self%geom%nz))
100 if (vars%has('q')) allocate(self%q(self%geom%nx,self%geom%ny,self%geom%nz))
101 if (vars%has('u')) allocate(self%u(self%geom%nx,self%geom%ny,self%geom%nz))
102 if (vars%has('v')) allocate(self%v(self%geom%nx,self%geom%ny,self%geom%nz))
103 
104 ! Allocate boundaries
105 if (self%lbc) then
106  ! Allocation
107  allocate(self%x_north(self%geom%nz))
108  allocate(self%x_south(self%geom%nz))
109  allocate(self%q_north(self%geom%nx,self%geom%nz))
110  allocate(self%q_south(self%geom%nx,self%geom%nz))
111 endif
112 
113 ! Initialize
114 call qg_fields_zero(self)
115 
116 end subroutine qg_fields_create
117 ! ------------------------------------------------------------------------------
118 !> Create fields from another one
119 subroutine qg_fields_create_from_other(self,other,geom)
120 
121 implicit none
122 
123 ! Passed variables
124 type(qg_fields),intent(inout) :: self !< Fields
125 type(qg_fields),intent(in) :: other !< Other fields
126 type(qg_geom),target,intent(in) :: geom !< Geometry
127 
128 ! Associate geometry
129 self%geom => geom
130 
131 ! Copy attributes
132 self%lbc = other%lbc
133 
134 ! Allocate 3d fields
135 if (allocated(other%x)) allocate(self%x(self%geom%nx,self%geom%ny,self%geom%nz))
136 if (allocated(other%q)) allocate(self%q(self%geom%nx,self%geom%ny,self%geom%nz))
137 if (allocated(other%u)) allocate(self%u(self%geom%nx,self%geom%ny,self%geom%nz))
138 if (allocated(other%v)) allocate(self%v(self%geom%nx,self%geom%ny,self%geom%nz))
139 
140 ! Allocate boundaries
141 if (self%lbc) then
142  ! Allocation
143  allocate(self%x_north(self%geom%nz))
144  allocate(self%x_south(self%geom%nz))
145  allocate(self%q_north(self%geom%nx,self%geom%nz))
146  allocate(self%q_south(self%geom%nx,self%geom%nz))
147 endif
148 
149 ! Initialize
150 call qg_fields_zero(self)
151 
152 end subroutine qg_fields_create_from_other
153 ! ------------------------------------------------------------------------------
154 !> Delete fields
155 subroutine qg_fields_delete(self)
156 
157 implicit none
158 
159 ! Passed variables
160 type(qg_fields),intent(inout) :: self !< Fields
161 
162 ! Release memory
163 if (allocated(self%x)) deallocate(self%x)
164 if (allocated(self%q)) deallocate(self%q)
165 if (allocated(self%u)) deallocate(self%u)
166 if (allocated(self%v)) deallocate(self%v)
167 if (allocated(self%x_north)) deallocate(self%x_north)
168 if (allocated(self%x_south)) deallocate(self%x_south)
169 if (allocated(self%q_north)) deallocate(self%q_north)
170 if (allocated(self%q_south)) deallocate(self%q_south)
171 
172 end subroutine qg_fields_delete
173 ! ------------------------------------------------------------------------------
174 !> Set fields to zero
175 subroutine qg_fields_zero(self)
176 
177 implicit none
178 
179 ! Passed variables
180 type(qg_fields),intent(inout) :: self
181 
182 ! Check field
183 call qg_fields_check(self)
184 
185 ! Set fields to zero
186 if (allocated(self%x)) self%x = 0.0_kind_real
187 if (allocated(self%q)) self%q = 0.0_kind_real
188 if (allocated(self%u)) self%u = 0.0_kind_real
189 if (allocated(self%v)) self%v = 0.0_kind_real
190 if (self%lbc) then
191  self%x_north = 0.0_kind_real
192  self%x_south = 0.0_kind_real
193  self%q_north = 0.0_kind_real
194  self%q_south = 0.0_kind_real
195 endif
196 
197 end subroutine qg_fields_zero
198 ! ------------------------------------------------------------------------------
199 !> Set fields to ones
200 subroutine qg_fields_ones(self)
201 
202 implicit none
203 
204 ! Passed variables
205 type(qg_fields),intent(inout) :: self
206 
207 ! Check field
208 call qg_fields_check(self)
209 
210 ! Set fields to ones
211 if (allocated(self%x)) self%x = 1.0
212 if (allocated(self%q)) self%q = 1.0
213 if (allocated(self%u)) self%u = 1.0
214 if (allocated(self%v)) self%v = 1.0
215 if (self%lbc) then
216  self%x_north = 1.0
217  self%x_south = 1.0
218  self%q_north = 1.0
219  self%q_south = 1.0
220 endif
221 
222 end subroutine qg_fields_ones
223 ! ------------------------------------------------------------------------------
224 !> Set fields to Diracs
225 subroutine qg_fields_dirac(self,f_conf)
226 
227 implicit none
228 
229 ! Passed variables
230 type(qg_fields),intent(inout) :: self !< Fields
231 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
232 
233 ! Local variables
234 integer :: ndir,idir
235 integer,allocatable :: ixdir(:),iydir(:),izdir(:)
236 character(len=1) :: var
237 character(len=:),allocatable :: str
238 
239 ! Check field
240 call qg_fields_check(self)
241 
242 ! Get Diracs size
243 ndir = f_conf%get_size('ixdir')
244 if ((f_conf%get_size('iydir')/=ndir).or.(f_conf%get_size('izdir')/=ndir)) &
245  & call abor1_ftn('qg_fields_dirac: inconsistent sizes for ixdir, iydir and izdir')
246 
247 ! Allocation
248 allocate(ixdir(ndir))
249 allocate(iydir(ndir))
250 allocate(izdir(ndir))
251 
252 ! Get Diracs positions
253 call f_conf%get_or_die("ixdir",ixdir)
254 call f_conf%get_or_die("iydir",iydir)
255 call f_conf%get_or_die("izdir",izdir)
256 call f_conf%get_or_die("var",str)
257 var = str
258 
259 ! Check Diracs positions
260 if (any(ixdir<1).or.any(ixdir>self%geom%nx)) call abor1_ftn('qg_fields_dirac: invalid ixdir')
261 if (any(iydir<1).or.any(iydir>self%geom%ny)) call abor1_ftn('qg_fields_dirac: invalid iydir')
262 if (any(izdir<1).or.any(izdir>self%geom%nz)) call abor1_ftn('qg_fields_dirac: invalid izdir')
263 
264 ! Setup Diracs
265 call qg_fields_zero(self)
266 do idir=1,ndir
267  select case (var)
268  case ('x')
269  if (.not.allocated(self%x)) call abor1_ftn('qg_fields_dirac: x should be allocated')
270  self%x(ixdir(idir),iydir(idir),izdir(idir)) = 1.0
271  case ('q')
272  if (.not.allocated(self%q)) call abor1_ftn('qg_fields_dirac: q should be allocated')
273  self%q(ixdir(idir),iydir(idir),izdir(idir)) = 1.0
274  case default
275  call abor1_ftn('qg_fields_dirac: wrong variable')
276  endselect
277 enddo
278 
279 ! Complete other fields
280 call qg_fields_complete(self,var)
281 
282 end subroutine qg_fields_dirac
283 ! ------------------------------------------------------------------------------
284 !> Generate random fields
285 subroutine qg_fields_random(self,var)
286 
287 implicit none
288 
289 ! Passed variables
290 type(qg_fields),intent(inout) :: self !< Fields
291 character(len=1),intent(in) :: var !< Variable to randomize ('x' or 'q')
292 
293 ! Check field
294 call qg_fields_check(self)
295 
296 ! Set at random value
297 select case (var)
298 case ('x')
299  if (.not.allocated(self%x)) call abor1_ftn('qg_fields_random: x should be allocated')
300  call normal_distribution(self%x,0.0_kind_real,1.0_kind_real,rseed)
301 case ('q')
302  if (.not.allocated(self%q)) call abor1_ftn('qg_fields_random: q should be allocated')
303  call normal_distribution(self%q,0.0_kind_real,1.0_kind_real,rseed)
304 case default
305  call abor1_ftn('qg_fields_random: wrong variable')
306 endselect
307 
308 ! Complete other fields
309 call qg_fields_complete(self,var)
310 
311 end subroutine qg_fields_random
312 ! ------------------------------------------------------------------------------
313 !> Copy fields
314 subroutine qg_fields_copy(self,other)
315 
316 implicit none
317 
318 ! Passed variables
319 type(qg_fields),intent(inout) :: self !< Fields
320 type(qg_fields),intent(in) :: other !< Other fields
321 
322 ! Check resolution
323 call qg_fields_check_resolution(self,other)
324 
325 ! Copy 3D field
326 if (allocated(self%x).and.allocated(other%x)) self%x = other%x
327 if (allocated(self%q).and.allocated(other%q)) self%q = other%q
328 if (allocated(self%u).and.allocated(other%u)) self%u = other%u
329 if (allocated(self%v).and.allocated(other%v)) self%v = other%v
330 
331 ! Copy LBC
332 call qg_fields_copy_lbc(self,other)
333 
334 end subroutine qg_fields_copy
335 ! ------------------------------------------------------------------------------
336 !> Copy fields LBC
337 subroutine qg_fields_copy_lbc(self,other)
338 
339 implicit none
340 
341 ! Passed variables
342 type(qg_fields),intent(inout) :: self !< Fields
343 type(qg_fields),intent(in) :: other !< Other fields
344 
345 ! Check resolution
346 call qg_fields_check_resolution(self,other)
347 
348 if (self%lbc) then
349  if (other%lbc) then
350  self%x_north = other%x_north
351  self%x_south = other%x_south
352  self%q_north = other%q_north
353  self%q_south = other%q_south
354  else
355  self%x_north = 0.0_kind_real
356  self%x_south = 0.0_kind_real
357  self%q_north = 0.0_kind_real
358  self%q_south = 0.0_kind_real
359  endif
360 endif
361 
362 end subroutine qg_fields_copy_lbc
363 ! ------------------------------------------------------------------------------
364 !> Add fields
365 subroutine qg_fields_self_add(self,rhs)
366 
367 implicit none
368 
369 ! Passed variables
370 type(qg_fields),intent(inout) :: self !< Fields
371 type(qg_fields),intent(in) :: rhs !< Right-hand side
372 
373 ! Check resolution
374 call qg_fields_check_resolution(self,rhs)
375 
376 ! Add field
377 if (allocated(self%x).and.allocated(rhs%x)) self%x = self%x + rhs%x
378 if (allocated(self%q).and.allocated(rhs%q)) self%q = self%q + rhs%q
379 if (allocated(self%u).and.allocated(rhs%u)) self%u = self%u + rhs%u
380 if (allocated(self%v).and.allocated(rhs%v)) self%v = self%v + rhs%v
381 if (self%lbc.and.rhs%lbc) then
382  self%x_north = self%x_north+rhs%x_north
383  self%x_south = self%x_south+rhs%x_south
384  self%q_north = self%q_north+rhs%q_north
385  self%q_south = self%q_south+rhs%q_south
386 endif
387 
388 end subroutine qg_fields_self_add
389 ! ------------------------------------------------------------------------------
390 !> Subtract fields
391 subroutine qg_fields_self_sub(self,rhs)
392 
393 implicit none
394 
395 ! Passed variables
396 type(qg_fields),intent(inout) :: self !< Fields
397 type(qg_fields),intent(in) :: rhs !< Right-hand side
398 
399 ! Check resolution
400 call qg_fields_check_resolution(self,rhs)
401 
402 ! Subtract field
403 if (allocated(self%x).and.allocated(rhs%x)) self%x = self%x - rhs%x
404 if (allocated(self%q).and.allocated(rhs%q)) self%q = self%q - rhs%q
405 if (allocated(self%u).and.allocated(rhs%u)) self%u = self%u - rhs%u
406 if (allocated(self%v).and.allocated(rhs%v)) self%v = self%v - rhs%v
407 if (self%lbc.and.rhs%lbc) then
408  self%x_north = self%x_north-rhs%x_north
409  self%x_south = self%x_south-rhs%x_south
410  self%q_north = self%q_north-rhs%q_north
411  self%q_south = self%q_south-rhs%q_south
412 endif
413 
414 end subroutine qg_fields_self_sub
415 ! ------------------------------------------------------------------------------
416 !> Multiply fields by a scalar
417 subroutine qg_fields_self_mul(self,zz)
418 
419 implicit none
420 
421 ! Passed variables
422 type(qg_fields),intent(inout) :: self !< Fields
423 real(kind_real),intent(in) :: zz !< Multiplier
424 
425 ! Check field
426 call qg_fields_check(self)
427 
428 ! Multiply with a scalar
429 if (allocated(self%x)) self%x = zz * self%x
430 if (allocated(self%q)) self%q = zz * self%q
431 if (allocated(self%u)) self%u = zz * self%u
432 if (allocated(self%v)) self%v = zz * self%v
433 if (self%lbc) then
434  self%x_north = zz*self%x_north
435  self%x_south = zz*self%x_south
436  self%q_north = zz*self%q_north
437  self%q_south = zz*self%q_south
438 endif
439 
440 end subroutine qg_fields_self_mul
441 ! ------------------------------------------------------------------------------
442 !> Apply axpy operator to fields
443 subroutine qg_fields_axpy(self,zz,rhs)
444 
445 implicit none
446 
447 ! Passed variables
448 type(qg_fields),intent(inout) :: self !< Fields
449 real(kind_real),intent(in) :: zz !< Multiplier
450 type(qg_fields),intent(in) :: rhs !< Right-hand side
451 
452 ! Check resolution
453 call qg_fields_check_resolution(self,rhs)
454 
455 ! Apply apxy
456 if (allocated(self%x).and.allocated(rhs%x)) self%x = self%x + zz * rhs%x
457 if (allocated(self%q).and.allocated(rhs%q)) self%q = self%q + zz * rhs%q
458 if (allocated(self%u).and.allocated(rhs%u)) self%u = self%u + zz * rhs%u
459 if (allocated(self%v).and.allocated(rhs%v)) self%v = self%v + zz * rhs%v
460 if (self%lbc.and.rhs%lbc) then
461  self%x_north = self%x_north+zz*rhs%x_north
462  self%x_south = self%x_south+zz*rhs%x_south
463  self%q_north = self%q_north+zz*rhs%q_north
464  self%q_south = self%q_south+zz*rhs%q_south
465 endif
466 
467 end subroutine qg_fields_axpy
468 ! ------------------------------------------------------------------------------
469 !> Schur product of fields
470 subroutine qg_fields_self_schur(self,rhs)
471 
472 implicit none
473 
474 ! Passed variables
475 type(qg_fields),intent(inout) :: self !< Fields
476 type(qg_fields),intent(in) :: rhs !< Right-hand side
477 
478 ! Check resolution
479 call qg_fields_check_resolution(self,rhs)
480 
481 ! Schur product
482 if (allocated(self%x).and.allocated(rhs%x)) self%x = self%x * rhs%x
483 if (allocated(self%q).and.allocated(rhs%q)) self%q = self%q * rhs%q
484 if (allocated(self%u).and.allocated(rhs%u)) self%u = self%u * rhs%u
485 if (allocated(self%v).and.allocated(rhs%v)) self%v = self%v * rhs%v
486 if (self%lbc.and.rhs%lbc) then
487  self%x_north = self%x_north*rhs%x_north
488  self%x_south = self%x_south*rhs%x_south
489  self%q_north = self%q_north*rhs%q_north
490  self%q_south = self%q_south*rhs%q_south
491 endif
492 
493 end subroutine qg_fields_self_schur
494 ! ------------------------------------------------------------------------------
495 !> Compute dot product for fields
496 subroutine qg_fields_dot_prod(fld1,fld2,zprod)
497 
498 implicit none
499 
500 ! Passed variables
501 type(qg_fields),intent(in) :: fld1 !< First fields
502 type(qg_fields),intent(in) :: fld2 !< Second fields
503 real(kind_real),intent(out) :: zprod !< Dot product
504 
505 ! Check resolution
506 call qg_fields_check_resolution(fld1,fld2)
507 
508 ! Compute dot product
509 zprod = 0.0_kind_real
510 if (allocated(fld1%x).and.allocated(fld2%x)) zprod = zprod+sum(fld1%x*fld2%x)
511 if (allocated(fld1%q).and.allocated(fld2%q)) zprod = zprod+sum(fld1%q*fld2%q)
512 if (allocated(fld1%u).and.allocated(fld2%u)) zprod = zprod+sum(fld1%u*fld2%u)
513 if (allocated(fld1%v).and.allocated(fld2%v)) zprod = zprod+sum(fld1%v*fld2%v)
514 
515 end subroutine qg_fields_dot_prod
516 ! ------------------------------------------------------------------------------
517 !> Add increment to fields
518 subroutine qg_fields_add_incr(self,rhs)
519 
520 implicit none
521 
522 ! Passed variables
523 type(qg_fields),intent(inout) :: self !< Fields
524 type(qg_fields),intent(in) :: rhs !< Right-hand side
525 
526 ! Check fields
527 call qg_fields_check(self)
528 call qg_fields_check(rhs)
529 
530 if ((self%geom%nx==rhs%geom%nx).and.(self%geom%ny==rhs%geom%ny).and.(self%geom%nz==rhs%geom%nz)) then
531  ! Same resolution
532  if (allocated(self%x).and.allocated(rhs%x)) self%x = self%x + rhs%x
533  if (allocated(self%q).and.allocated(rhs%q)) self%q = self%q + rhs%q
534  if (allocated(self%u).and.allocated(rhs%u)) self%u = self%u + rhs%u
535  if (allocated(self%v).and.allocated(rhs%v)) self%v = self%v + rhs%v
536 else
537  ! Different resolutions
538  call abor1_ftn('qg_fields_add_incr: not coded for low res increment yet')
539 endif
540 
541 end subroutine qg_fields_add_incr
542 ! ------------------------------------------------------------------------------
543 !> Compute increment from the difference of two fields
544 subroutine qg_fields_diff_incr(lhs,fld1,fld2)
545 
546 implicit none
547 
548 ! Passed variables
549 type(qg_fields),intent(inout) :: lhs !< Left-hand side
550 type(qg_fields),intent(in) :: fld1 !< First fields
551 type(qg_fields),intent(in) :: fld2 !< Second fields
552 
553 ! Check resolution
554 call qg_fields_check_resolution(fld1,fld2)
555 call qg_fields_check(lhs)
556 
557 ! Initialization
558 call qg_fields_zero(lhs)
559 
560 if ((fld1%geom%nx==lhs%geom%nx).and.(fld1%geom%ny==lhs%geom%ny).and.(fld1%geom%nz==lhs%geom%nz)) then
561  ! Same resolution
562  if (allocated(lhs%x).and.allocated(fld1%x)) lhs%x = fld1%x - fld2%x
563  if (allocated(lhs%q).and.allocated(fld1%q)) lhs%q = fld1%q - fld2%q
564  if (allocated(lhs%u).and.allocated(fld1%u)) lhs%u = fld1%u - fld2%u
565  if (allocated(lhs%v).and.allocated(fld1%v)) lhs%v = fld1%v - fld2%v
566 else
567  ! Different resolutions
568  call abor1_ftn('qg_fields_diff_incr: not coded for low res increment yet')
569 endif
570 
571 end subroutine qg_fields_diff_incr
572 ! ------------------------------------------------------------------------------
573 !> Change fields resolution
574 subroutine qg_fields_change_resol(fld,rhs)
575 
576 implicit none
577 
578 ! Passed variables
579 type(qg_fields),intent(inout) :: fld !< Fields
580 type(qg_fields),intent(in) :: rhs !< Right-hand side
581 
582 integer :: ix,iy,iz
583 real(kind_real), allocatable, dimension(:,:,:) :: q1, q2
584 
585 if ((fld%geom%nx==rhs%geom%nx).and.(fld%geom%ny==rhs%geom%ny).and.(fld%geom%nz==rhs%geom%nz)) then
586  ! Same resolution
587  call qg_fields_copy(fld,rhs)
588 else
589  ! Trilinear interpolation
590  do ix=1,fld%geom%nx
591  do iy=1,fld%geom%ny
592  do iz=1,fld%geom%nz
593  if (allocated(rhs%x).and.allocated(fld%x)) then
594  call qg_interp_trilinear(rhs%geom,fld%geom%lon(ix,iy),fld%geom%lat(ix,iy),fld%geom%z(iz), &
595  rhs%x,fld%x(ix,iy,iz))
596  endif
597  if (allocated(rhs%q).and.allocated(fld%q)) then
598  call qg_interp_trilinear(rhs%geom,fld%geom%lon(ix,iy),fld%geom%lat(ix,iy),fld%geom%z(iz), &
599  rhs%q,fld%q(ix,iy,iz))
600  endif
601  if (allocated(rhs%u).and.allocated(fld%u)) then
602  call qg_interp_trilinear(rhs%geom,fld%geom%lon(ix,iy),fld%geom%lat(ix,iy),fld%geom%z(iz), &
603  rhs%u,fld%u(ix,iy,iz))
604  endif
605  if (allocated(rhs%v).and.allocated(fld%v)) then
606  call qg_interp_trilinear(rhs%geom,fld%geom%lon(ix,iy),fld%geom%lat(ix,iy),fld%geom%z(iz), &
607  rhs%v,fld%v(ix,iy,iz))
608  endif
609  enddo
610  enddo
611  enddo
612 
613  ! Deal with boundary conditions
614  if (fld%lbc) then
615  if (rhs%lbc) then
616  allocate(q1(rhs%geom%nx,rhs%geom%ny,rhs%geom%nz))
617  allocate(q2(fld%geom%nx,fld%geom%ny,fld%geom%nz))
618  do iy=1,rhs%geom%ny
619  q1(:,iy,:) = rhs%q_south
620  enddo
621  do ix=1,fld%geom%nx
622  do iz=1,fld%geom%nz
623  call qg_interp_trilinear(rhs%geom,fld%geom%lon(ix,1),fld%geom%lat(ix,1),fld%geom%z(iz),q1,q2(ix,1,iz) )
624  enddo
625  enddo
626  fld%q_south = q2(:,1,:)
627  do iy=1,rhs%geom%ny
628  q1(:,iy,:) = rhs%q_north
629  enddo
630  do ix=1,fld%geom%nx
631  do iz=1,fld%geom%nz
632  call qg_interp_trilinear(rhs%geom,fld%geom%lon(ix,1),fld%geom%lat(ix,1),fld%geom%z(iz),q1,q2(ix,1,iz))
633  enddo
634  enddo
635  fld%q_north = q2(:,1,:)
636  deallocate(q1,q2)
637  fld%x_north = rhs%x_north
638  fld%x_south = rhs%x_south
639  else
640  fld%x_north = 0.0_kind_real
641  fld%x_south = 0.0_kind_real
642  fld%q_north = 0.0_kind_real
643  fld%q_south = 0.0_kind_real
644  endif
645  endif
646 endif
647 
648 end subroutine qg_fields_change_resol
649 ! ------------------------------------------------------------------------------
650 !> Read fields from file
651 subroutine qg_fields_read_file(fld,f_conf,vdate)
652 use string_utils
653 
654 implicit none
655 
656 ! Passed variables
657 type(qg_fields),intent(inout) :: fld !< Fields
658 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
659 type(datetime),intent(inout) :: vdate !< Date and time
660 
661 ! Local variables
662 integer :: iread,nx,ny,nz,bc
663 integer :: ncid,nx_id,ny_id,nz_id,x_id,q_id,u_id,v_id,x_north_id,x_south_id,q_north_id,q_south_id
664 logical :: lbc
665 character(len=20) :: sdate
666 character(len=1024) :: record,filename
667 character(len=:),allocatable :: str
668 
669 ! Check field
670 call qg_fields_check(fld)
671 
672 ! Check whether the field should be invented or read from file
673 iread = 1
674 if (f_conf%has("read_from_file")) call f_conf%get_or_die("read_from_file",iread)
675 
676 if (iread==0) then
677  ! Invent field
678  call fckit_log%warning('qg_fields_read_file: inventing field')
679  call qg_fields_analytic_init(fld,f_conf,vdate)
680 else
681  ! Read field from file
682 
683  ! Get filename
684  call f_conf%get_or_die("filename",str)
685  call swap_name_member(f_conf, str)
686  filename = str
687  call fckit_log%info('qg_fields_read_file: opening '//trim(filename))
688 
689  ! Initialize field
690  call qg_fields_zero(fld)
691 
692  ! Open NetCDF file
693  call ncerr(nf90_open(trim(filename),nf90_nowrite,ncid))
694 
695  ! Get dimensions ids
696  call ncerr(nf90_inq_dimid(ncid,'nx',nx_id))
697  call ncerr(nf90_inq_dimid(ncid,'ny',ny_id))
698  call ncerr(nf90_inq_dimid(ncid,'nz',nz_id))
699 
700  ! Get dimensions
701  call ncerr(nf90_inquire_dimension(ncid,nx_id,len=nx))
702  call ncerr(nf90_inquire_dimension(ncid,ny_id,len=ny))
703  call ncerr(nf90_inquire_dimension(ncid,nz_id,len=nz))
704 
705  ! Test dimensions consistency with the field geometry
706  if ((nx/=fld%geom%nx).or.(ny/=fld%geom%ny).or.(nz/=fld%geom%nz)) then
707  write (record,*) 'qg_fields_read_file: input fields have wrong dimensions: ',nx,ny,nz
708  call fckit_log%error(record)
709  write (record,*) 'qg_fields_read_file: expected: ',fld%geom%nx,fld%geom%ny,fld%geom%nz
710  call fckit_log%error(record)
711  call abor1_ftn('qg_fields_read_file: input fields have wrong dimensions')
712  endif
713 
714  ! Get attributes
715  call ncerr(nf90_get_att(ncid,nf90_global,'bc',bc))
716  if (bc==0) then
717  lbc = .false.
718  elseif (bc==1) then
719  lbc = .true.
720  else
721  call abor1_ftn('qg_fields_read_file: wrong bc value')
722  endif
723  call ncerr(nf90_get_att(ncid,nf90_global,'sdate',sdate))
724 
725  ! Test attributes consistency with the field
726  if ((.not.lbc).and.fld%lbc) call abor1_ftn('qg_fields_read_file: LBC are missing in NetCDF file')
727 
728  ! Get variables ids
729  if (allocated(fld%x)) call ncerr(nf90_inq_varid(ncid,'x',x_id))
730  if (allocated(fld%q)) call ncerr(nf90_inq_varid(ncid,'q',q_id))
731  if (allocated(fld%u)) call ncerr(nf90_inq_varid(ncid,'u',u_id))
732  if (allocated(fld%v)) call ncerr(nf90_inq_varid(ncid,'v',v_id))
733  if (fld%lbc) then
734  call ncerr(nf90_inq_varid(ncid,'x_north',x_north_id))
735  call ncerr(nf90_inq_varid(ncid,'x_south',x_south_id))
736  call ncerr(nf90_inq_varid(ncid,'q_north',q_north_id))
737  call ncerr(nf90_inq_varid(ncid,'q_south',q_south_id))
738  endif
739 
740  ! Get variables
741  if (allocated(fld%x)) call ncerr(nf90_get_var(ncid,x_id,fld%x))
742  if (allocated(fld%q)) call ncerr(nf90_get_var(ncid,q_id,fld%q))
743  if (allocated(fld%u)) call ncerr(nf90_get_var(ncid,u_id,fld%u))
744  if (allocated(fld%v)) call ncerr(nf90_get_var(ncid,v_id,fld%v))
745  if (fld%lbc) then
746  call ncerr(nf90_get_var(ncid,x_north_id,fld%x_north))
747  call ncerr(nf90_get_var(ncid,x_south_id,fld%x_south))
748  call ncerr(nf90_get_var(ncid,q_north_id,fld%q_north))
749  call ncerr(nf90_get_var(ncid,q_south_id,fld%q_south))
750  endif
751 
752  ! Close NetCDF file
753  call ncerr(nf90_close(ncid))
754 
755  ! Set date
756  call fckit_log%info('qg_fields_read_file: validity date is '//sdate)
757  call datetime_set(sdate,vdate)
758 endif
759 
760 ! Check field
761 call qg_fields_check(fld)
762 
763 end subroutine qg_fields_read_file
764 ! ------------------------------------------------------------------------------
765 !> Write fields to file
766 subroutine qg_fields_write_file(fld,f_conf,vdate)
767 
768 implicit none
769 
770 ! Passed variables
771 type(qg_fields),intent(in) :: fld !< Fields
772 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
773 type(datetime),intent(in) :: vdate !< Date and time
774 
775 ! Local variables
776 integer :: ncid,nx_id,ny_id,nz_id,lon_id,lat_id,z_id,area_id,heat_id,x_id,q_id,u_id,v_id
777 integer :: x_north_id,x_south_id,q_north_id,q_south_id
778 integer :: info
779 character(len=20) :: sdate
780 character(len=1024) :: filename
781 type(oops_variables) :: vars
782 type(qg_fields) :: fld_io
783 
784 ! Check field
785 call qg_fields_check(fld)
786 
787 ! Get all variables
788 vars = oops_variables()
789 call vars%push_back('x')
790 call vars%push_back('q')
791 call vars%push_back('u')
792 call vars%push_back('v')
793 call qg_fields_create(fld_io,fld%geom,vars,.true.)
794 call qg_fields_copy_lbc(fld_io,fld)
795 if (allocated(fld%x)) then
796  fld_io%x = fld%x
797  call qg_fields_complete(fld_io,'x')
798 elseif (allocated(fld%q)) then
799  fld_io%q = fld%q
800  call qg_fields_complete(fld_io,'q')
801 else
802  call abor1_ftn('qg_fields_write_file: x or q required')
803 endif
804 
805 ! Set filename
806 filename = genfilename(f_conf,800,vdate)
807 call fckit_log%info('qg_fields_write_file: writing '//trim(filename))
808 
809 ! Set date
810 call datetime_to_string(vdate,sdate)
811 
812 ! Create NetCDF file
813 call ncerr(nf90_create(trim(filename),or(nf90_clobber,nf90_64bit_offset),ncid))
814 
815 ! Define dimensions
816 call ncerr(nf90_def_dim(ncid,'nx',fld%geom%nx,nx_id))
817 call ncerr(nf90_def_dim(ncid,'ny',fld%geom%ny,ny_id))
818 call ncerr(nf90_def_dim(ncid,'nz',fld%geom%nz,nz_id))
819 
820 ! Define attributes
821 if (fld%lbc) then
822  call ncerr(nf90_put_att(ncid,nf90_global,'bc',1))
823 else
824  call ncerr(nf90_put_att(ncid,nf90_global,'bc',0))
825 endif
826 call ncerr(nf90_put_att(ncid,nf90_global,'sdate',sdate))
827 
828 ! Define variables
829 call ncerr(nf90_def_var(ncid,'lon',nf90_double,(/nx_id,ny_id/),lon_id))
830 call ncerr(nf90_def_var(ncid,'lat',nf90_double,(/nx_id,ny_id/),lat_id))
831 call ncerr(nf90_def_var(ncid,'z',nf90_double,(/nz_id/),z_id))
832 call ncerr(nf90_def_var(ncid,'area',nf90_double,(/nx_id,ny_id/),area_id))
833 call ncerr(nf90_def_var(ncid,'heat',nf90_double,(/nx_id,ny_id/),heat_id))
834 call ncerr(nf90_def_var(ncid,'x',nf90_double,(/nx_id,ny_id,nz_id/),x_id))
835 call ncerr(nf90_put_att(ncid,x_id,'_FillValue',missing_value(1.0_kind_real)))
836 call ncerr(nf90_def_var(ncid,'q',nf90_double,(/nx_id,ny_id,nz_id/),q_id))
837 call ncerr(nf90_put_att(ncid,q_id,'_FillValue',missing_value(1.0_kind_real)))
838 call ncerr(nf90_def_var(ncid,'u',nf90_double,(/nx_id,ny_id,nz_id/),u_id))
839 call ncerr(nf90_put_att(ncid,u_id,'_FillValue',missing_value(1.0_kind_real)))
840 call ncerr(nf90_def_var(ncid,'v',nf90_double,(/nx_id,ny_id,nz_id/),v_id))
841 call ncerr(nf90_put_att(ncid,v_id,'_FillValue',missing_value(1.0_kind_real)))
842 if (fld%lbc) then
843  call ncerr(nf90_def_var(ncid,'x_north',nf90_double,(/nz_id/),x_north_id))
844  call ncerr(nf90_put_att(ncid,x_north_id,'_FillValue',missing_value(1.0_kind_real)))
845  call ncerr(nf90_def_var(ncid,'x_south',nf90_double,(/nz_id/),x_south_id))
846  call ncerr(nf90_put_att(ncid,x_south_id,'_FillValue',missing_value(1.0_kind_real)))
847  call ncerr(nf90_def_var(ncid,'q_north',nf90_double,(/nx_id,nz_id/),q_north_id))
848  call ncerr(nf90_put_att(ncid,q_north_id,'_FillValue',missing_value(1.0_kind_real)))
849  call ncerr(nf90_def_var(ncid,'q_south',nf90_double,(/nx_id,nz_id/),q_south_id))
850  call ncerr(nf90_put_att(ncid,q_south_id,'_FillValue',missing_value(1.0_kind_real)))
851 endif
852 
853 ! End definitions
854 call ncerr(nf90_enddef(ncid))
855 
856 ! Put variables
857 call ncerr(nf90_put_var(ncid,lon_id,fld%geom%lon))
858 call ncerr(nf90_put_var(ncid,lat_id,fld%geom%lat))
859 call ncerr(nf90_put_var(ncid,z_id,fld%geom%z))
860 call ncerr(nf90_put_var(ncid,area_id,fld%geom%area))
861 call ncerr(nf90_put_var(ncid,heat_id,fld%geom%heat))
862 call ncerr(nf90_put_var(ncid,x_id,fld_io%x))
863 call ncerr(nf90_put_var(ncid,q_id,fld_io%q))
864 call ncerr(nf90_put_var(ncid,u_id,fld_io%u))
865 call ncerr(nf90_put_var(ncid,v_id,fld_io%v))
866 if (fld%lbc) then
867  call ncerr(nf90_put_var(ncid,x_north_id,fld%x_north))
868  call ncerr(nf90_put_var(ncid,x_south_id,fld%x_south))
869  call ncerr(nf90_put_var(ncid,q_north_id,fld%q_north))
870  call ncerr(nf90_put_var(ncid,q_south_id,fld%q_south))
871 endif
872 
873 ! Close NetCDF file
874 call ncerr(nf90_close(ncid))
875 
876 ! Release memory
877 call vars%destruct()
878 
879 end subroutine qg_fields_write_file
880 ! ------------------------------------------------------------------------------
881 !> Analytic initialization of fields
882 subroutine qg_fields_analytic_init(fld,f_conf,vdate)
883 
884 implicit none
885 
886 ! Passed variables
887 type(qg_fields),intent(inout) :: fld !< Fields
888 type(fckit_configuration),intent(in) :: f_conf !< FCKIT configuration
889 type(datetime),intent(inout) :: vdate !< Date and time
890 
891 ! Local variables
892 integer :: ix,iy,iz
893 real(kind_real) :: uval
894 real(kind_real),allocatable :: x(:,:,:),q(:,:,:)
895 character(len=30) :: ic
896 character(len=20) :: sdate
897 character(len=:),allocatable :: str
898 
899 ! Check configuration
900 if (f_conf%has("analytic_init")) then
901  call f_conf%get_or_die("analytic_init",str)
902  ic = str
903 else
904  ic = 'baroclinic-instability'
905 endif
906 call fckit_log%warning('qg_fields_analytic_init: '//trim(ic))
907 
908 ! Set date
909 call f_conf%get_or_die("date",str)
910 sdate = str
911 call fckit_log%info('qg_fields_analytic_init: validity date is '//sdate)
912 call datetime_set(sdate,vdate)
913 
914 ! Check boundaries
915 if (.not.fld%lbc) call abor1_ftn('qg_fields_analytic_init: boundaries required')
916 
917 ! Allocation
918 allocate(x(fld%geom%nx,fld%geom%ny,fld%geom%nz))
919 allocate(q(fld%geom%nx,fld%geom%ny,fld%geom%nz))
920 
921 ! Define state
922 select case (trim(ic))
923 case ('baroclinic-instability')
924  ! Baroclinic instability
925  do iz=1,fld%geom%nz
926  do iy=1,fld%geom%ny
927  do ix=1,fld%geom%nx
928  call baroclinic_instability(fld%geom%x(ix),fld%geom%y(iy),fld%geom%z(iz),'x',x(ix,iy,iz))
929  enddo
930  enddo
931  call baroclinic_instability(0.0_kind_real,domain_meridional,fld%geom%z(iz),'x',fld%x_north(iz))
932  call baroclinic_instability(0.0_kind_real,0.0_kind_real,fld%geom%z(iz),'x',fld%x_south(iz))
933  enddo
934 case ('large-vortices')
935  ! Large vortices
936  do iz=1,fld%geom%nz
937  do iy=1,fld%geom%ny
938  do ix=1,fld%geom%nx
939  call large_vortices(fld%geom%x(ix),fld%geom%y(iy),fld%geom%z(iz),'x',x(ix,iy,iz))
940  enddo
941  enddo
942  call large_vortices(0.0_kind_real,domain_meridional,fld%geom%z(iz),'x',fld%x_north(iz))
943  call large_vortices(0.0_kind_real,0.0_kind_real,fld%geom%z(iz),'x',fld%x_south(iz))
944  enddo
945 case ('uniform_field')
946  ! Uniform field
947  call f_conf%get_or_die("uval",uval)
948  x = uval
949 case default
950  call abor1_ftn('qg_fields_analytic_init: unknown initialization')
951 endselect
952 
953 ! Compute q
954 call convert_x_to_q(fld%geom,x,fld%x_north,fld%x_south,q)
955 do iz=1,fld%geom%nz
956  do ix=1,fld%geom%nx
957  fld%q_south(ix,iz) = 2.0*q(ix,1,iz)-q(ix,2,iz)
958  enddo
959  do ix=1,fld%geom%nx
960  fld%q_north(ix,iz) = 2.0*q(ix,fld%geom%ny,iz)-q(ix,fld%geom%ny-1,iz)
961  enddo
962 enddo
963 
964 ! Copy 3d field and ensure consistency
965 if (allocated(x)) then
966  fld%x = x
967  call qg_fields_complete(fld,'x')
968 elseif (allocated(q)) then
969  fld%q = q
970  call qg_fields_complete(fld,'q')
971 else
972  call abor1_ftn('qg_fields_analytic_init: x or q required')
973 endif
974 
975 ! Check field
976 call qg_fields_check(fld)
977 
978 end subroutine qg_fields_analytic_init
979 ! ------------------------------------------------------------------------------
980 !> Fields statistics
981 subroutine qg_fields_gpnorm(fld,vpresent,vmin,vmax,vrms)
982 
983 implicit none
984 
985 ! Passed variables
986 type(qg_fields),intent(in) :: fld !< Fields
987 integer,intent(inout) :: vpresent(6) !< Variables presence flag
988 real(kind_real),intent(inout) :: vmin(6) !< Variables minimum
989 real(kind_real),intent(inout) :: vmax(6) !< Variables maximum
990 real(kind_real),intent(inout) :: vrms(6) !< Variables RMS
991 
992 ! Check field
993 call qg_fields_check(fld)
994 
995 ! Initialization
996 vpresent = 0
997 vmin = 0.0_kind_real
998 vmax = 0.0_kind_real
999 vrms = 0.0_kind_real
1000 
1001 ! 3d fields
1002 if (allocated(fld%x)) then
1003  vpresent(1) = 1
1004  vmin(1) = minval(fld%x)
1005  vmax(1) = maxval(fld%x)
1006  vrms(1) = sqrt(sum(fld%x**2)/real(fld%geom%nx*fld%geom%ny*fld%geom%nz,kind_real))
1007 endif
1008 if (allocated(fld%q)) then
1009  vpresent(2) = 1
1010  vmin(2) = minval(fld%q)
1011  vmax(2) = maxval(fld%q)
1012  vrms(2) = sqrt(sum(fld%q**2)/real(fld%geom%nx*fld%geom%ny*fld%geom%nz,kind_real))
1013 endif
1014 if (allocated(fld%u)) then
1015  vpresent(3) = 1
1016  vmin(3) = minval(fld%u)
1017  vmax(3) = maxval(fld%u)
1018  vrms(3) = sqrt(sum(fld%u**2)/real(fld%geom%nx*fld%geom%ny*fld%geom%nz,kind_real))
1019 endif
1020 if (allocated(fld%v)) then
1021  vpresent(4) = 1
1022  vmin(4) = minval(fld%v)
1023  vmax(4) = maxval(fld%v)
1024  vrms(4) = sqrt(sum(fld%v**2)/real(fld%geom%nx*fld%geom%ny*fld%geom%nz,kind_real))
1025 endif
1026 
1027 ! Boundaries
1028 if (fld%lbc) then
1029  ! Streamfunction
1030  vpresent(5) = 1
1031  vmin(5) = min(minval(fld%x_north),minval(fld%x_south))
1032  vmax(5) = max(maxval(fld%x_north),maxval(fld%x_south))
1033  vrms(5) = sqrt(sum(fld%x_north**2+fld%x_south**2)/real(2*fld%geom%nz,kind_real))
1034 
1035  ! Potential vorticity
1036  vpresent(6) = 1
1037  vmin(6) = min(minval(fld%q_north),minval(fld%q_south))
1038  vmax(6) = max(maxval(fld%q_north),maxval(fld%q_south))
1039  vrms(6) = sqrt(sum(fld%q_north**2+fld%q_south**2)/real(2*fld%geom%nx*fld%geom%nz,kind_real))
1040 endif
1041 
1042 end subroutine qg_fields_gpnorm
1043 ! ------------------------------------------------------------------------------
1044 !> Fields RMS
1045 subroutine qg_fields_rms(fld,prms)
1046 
1047 implicit none
1048 
1049 ! Passed variables
1050 type(qg_fields),intent(in) :: fld !< Fields
1051 real(kind_real),intent(out) :: prms !< RMS
1052 
1053 ! Local variables
1054 integer :: norm
1055 
1056 ! Check field
1057 call qg_fields_check(fld)
1058 
1059 ! Initialization
1060 prms = 0.0_kind_real
1061 norm = 0.0_kind_real
1062 
1063 ! 3d fields
1064 if (allocated(fld%x)) then
1065  prms = prms+sum(fld%x**2)
1066  norm = norm+real(fld%geom%nx*fld%geom%ny*fld%geom%nz,kind_real)
1067 endif
1068 if (allocated(fld%q)) then
1069  prms = prms+sum(fld%q**2)
1070  norm = norm+real(fld%geom%nx*fld%geom%ny*fld%geom%nz,kind_real)
1071 endif
1072 if (allocated(fld%u)) then
1073  prms = prms+sum(fld%u**2)
1074  norm = norm+real(fld%geom%nx*fld%geom%ny*fld%geom%nz,kind_real)
1075 endif
1076 if (allocated(fld%v)) then
1077  prms = prms+sum(fld%v**2)
1078  norm = norm+real(fld%geom%nx*fld%geom%ny*fld%geom%nz,kind_real)
1079 endif
1080 
1081 ! Boundaries
1082 if (fld%lbc) then
1083  prms = prms+sum(fld%x_north**2+fld%x_south**2)+sum(fld%q_north**2+fld%q_south**2)
1084  norm = norm+real(2*(1+fld%geom%nx)*fld%geom%nz,kind_real)
1085 endif
1086 
1087 ! Normalize and square-root
1088 prms = sqrt(prms/norm)
1089 
1090 end subroutine qg_fields_rms
1091 ! ------------------------------------------------------------------------------
1092 !> Get fields geometry
1093 subroutine qg_fields_sizes(fld,nx,ny,nz)
1094 
1095 implicit none
1096 
1097 ! Passed variables
1098 type(qg_fields),intent(in) :: fld !< Fields
1099 integer,intent(out) :: nx !< X size
1100 integer,intent(out) :: ny !< Y size
1101 integer,intent(out) :: nz !< Z size
1102 
1103 ! Copy sizes
1104 nx = fld%geom%nx
1105 ny = fld%geom%ny
1106 nz = fld%geom%nz
1107 
1108 end subroutine qg_fields_sizes
1109 ! ------------------------------------------------------------------------------
1110 !> Get LBC presence
1111 subroutine qg_fields_lbc(fld,lbc)
1112 
1113 implicit none
1114 
1115 ! Passed variables
1116 type(qg_fields),intent(in) :: fld !< Fields
1117 integer,intent(out) :: lbc !< LBC presence
1118 
1119 ! Check field
1120 call qg_fields_check(fld)
1121 
1122 if (fld%lbc) then
1123  lbc = 1
1124 else
1125  lbc = 0
1126 endif
1127 
1128 end subroutine qg_fields_lbc
1129 ! ------------------------------------------------------------------------------
1130 !> Set ATLAS field
1131 subroutine qg_fields_set_atlas(self,vars,afieldset)
1132 
1133 implicit none
1134 
1135 ! Passed variables
1136 type(qg_fields),intent(in) :: self !< Fields
1137 type(oops_variables),intent(in) :: vars !< List of variables
1138 type(atlas_fieldset),intent(inout) :: afieldset !< ATLAS fieldset
1139 
1140 ! Local variables
1141 integer :: jvar
1142 character(len=1024) :: fieldname
1143 type(atlas_field) :: afield
1144 
1145 ! Get or create field
1146 do jvar=1,vars%nvars()
1147  fieldname = vars%variable(jvar)
1148  if (.not.afieldset%has_field(trim(fieldname))) then
1149  ! Create field
1150  afield = self%geom%afunctionspace%create_field(name=trim(fieldname),kind=atlas_real(kind_real),levels=self%geom%nz)
1151 
1152  ! Add field
1153  call afieldset%add(afield)
1154 
1155  ! Release pointer
1156  call afield%final()
1157  endif
1158 enddo
1159 
1160 end subroutine qg_fields_set_atlas
1161 ! ------------------------------------------------------------------------------
1162 !> Convert fields to ATLAS
1163 subroutine qg_fields_to_atlas(self,vars,afieldset)
1164 
1165 implicit none
1166 
1167 ! Passed variables
1168 type(qg_fields),intent(in) :: self !< Fields
1169 type(oops_variables),intent(in) :: vars !< List of variables
1170 type(atlas_fieldset),intent(inout) :: afieldset !< ATLAS fieldset
1171 
1172 ! Local variables
1173 integer :: jvar,ix,iy,iz,inode
1174 real(kind_real),pointer :: ptr(:,:)
1175 character(len=1024) :: fieldname
1176 type(atlas_field) :: afield
1177 
1178 ! Get variable
1179 do jvar=1,vars%nvars()
1180  fieldname = vars%variable(jvar)
1181  if (afieldset%has_field(trim(fieldname))) then
1182  ! Get afield
1183  afield = afieldset%field(trim(fieldname))
1184  else
1185  ! Create field
1186  afield = self%geom%afunctionspace%create_field(name=trim(fieldname),kind=atlas_real(kind_real),levels=self%geom%nz)
1187 
1188  ! Add field
1189  call afieldset%add(afield)
1190  endif
1191 
1192  ! Copy field
1193  call afield%data(ptr)
1194  do iz=1,self%geom%nz
1195  inode = 0
1196  do iy=1,self%geom%ny
1197  do ix=1,self%geom%nx
1198  inode = inode+1
1199  select case (trim(fieldname))
1200  case ('x')
1201  ptr(iz,inode) = self%x(ix,iy,iz)
1202  case ('q')
1203  ptr(iz,inode) = self%q(ix,iy,iz)
1204  case ('u')
1205  ptr(iz,inode) = self%u(ix,iy,iz)
1206  case ('v')
1207  ptr(iz,inode) = self%v(ix,iy,iz)
1208  case default
1209  call abor1_ftn('qg_fields_to_atlas: wrong variable')
1210  endselect
1211  enddo
1212  enddo
1213  enddo
1214 
1215  ! Release pointer
1216  call afield%final()
1217 enddo
1218 
1219 end subroutine qg_fields_to_atlas
1220 ! ------------------------------------------------------------------------------
1221 !> Get fields from ATLAS
1222 subroutine qg_fields_from_atlas(self,vars,afieldset)
1223 
1224 implicit none
1225 
1226 ! Passed variables
1227 type(qg_fields),intent(inout) :: self !< Fields
1228 type(oops_variables),intent(in) :: vars !< List of variables
1229 type(atlas_fieldset),intent(inout) :: afieldset !< ATLAS fieldset
1230 
1231 ! Local variables
1232 integer :: jvar,ix,iy,iz,inode
1233 real(kind_real),pointer :: ptr(:,:)
1234 character(len=1024) :: fieldname
1235 type(atlas_field) :: afield
1236 
1237 ! Get variable
1238 do jvar=1,vars%nvars()
1239  ! Get afield
1240  fieldname = vars%variable(jvar)
1241  afield = afieldset%field(trim(fieldname))
1242 
1243  ! Copy field
1244  call afield%data(ptr)
1245  do iz=1,self%geom%nz
1246  inode = 0
1247  do iy=1,self%geom%ny
1248  do ix=1,self%geom%nx
1249  inode = inode+1
1250  select case (trim(fieldname))
1251  case ('x')
1252  self%x(ix,iy,iz) = ptr(iz,inode)
1253  case ('q')
1254  self%q(ix,iy,iz) = ptr(iz,inode)
1255  case ('u')
1256  self%u(ix,iy,iz) = ptr(iz,inode)
1257  case ('v')
1258  self%v(ix,iy,iz) = ptr(iz,inode)
1259  case default
1260  call abor1_ftn('qg_fields_to_atlas: wrong variable')
1261  endselect
1262  enddo
1263  enddo
1264  enddo
1265 
1266  ! Release pointer
1267  call afield%final()
1268 enddo
1269 
1270 end subroutine qg_fields_from_atlas
1271 ! ------------------------------------------------------------------------------
1272 !> Get points from fields
1273 subroutine qg_fields_getpoint(fld,iter,nval,vals)
1274 
1275 implicit none
1276 
1277 ! Passed variables
1278 type(qg_fields),intent(in) :: fld !< Fields
1279 type(qg_geom_iter),intent(in) :: iter !< Geometry iterator
1280 integer,intent(in) :: nval !< Number of values
1281 real(kind_real),intent(inout) :: vals(nval) !< Values
1282 
1283 ! Local variables
1284 integer :: ii
1285 character(len=1024) :: record
1286 
1287 ! Check
1288 if ((fld%geom%nx/=iter%geom%nx).or.(fld%geom%ny/=iter%geom%ny).or.(fld%geom%nz/=iter%geom%nz)) then
1289  write(record,*) 'qg_fields_getpoint: resolution inconsistency, ',fld%geom%nx,'/',fld%geom%ny,'/',fld%geom%nz, &
1290 & ' and ',iter%geom%nx,'/',iter%geom%ny,'/',iter%geom%nz
1291  call abor1_ftn(record)
1292 endif
1293 if (fld%geom%nz/=nval) then
1294  write(record,*) 'qg_fields_getpoint: array sizes are different: ',fld%geom%nz,'/',nval
1295  call abor1_ftn(record)
1296 endif
1297 
1298 ! Initialization
1299 ii = 0
1300 
1301 ! Get values
1302 if (allocated(fld%x)) then
1303  vals(ii+1:ii+fld%geom%nz) = fld%x(iter%ilon,iter%ilat,:)
1304  ii = ii+fld%geom%nz
1305 endif
1306 if (allocated(fld%q)) then
1307  vals(ii+1:ii+fld%geom%nz) = fld%q(iter%ilon,iter%ilat,:)
1308  ii = ii+fld%geom%nz
1309 endif
1310 if (allocated(fld%u)) then
1311  vals(ii+1:ii+fld%geom%nz) = fld%u(iter%ilon,iter%ilat,:)
1312  ii = ii+fld%geom%nz
1313 endif
1314 if (allocated(fld%v)) then
1315  vals(ii+1:ii+fld%geom%nz) = fld%v(iter%ilon,iter%ilat,:)
1316  ii = ii+fld%geom%nz
1317 endif
1318 
1319 end subroutine qg_fields_getpoint
1320 ! ------------------------------------------------------------------------------
1321 !> Set fields values at a specified gridpoint
1322 subroutine qg_fields_setpoint(fld,iter,nval,vals)
1323 
1324 implicit none
1325 
1326 ! Passed variables
1327 type(qg_fields),intent(inout) :: fld !< Fields
1328 type(qg_geom_iter),intent(in) :: iter !< Geometry iterator
1329 integer,intent(in) :: nval !< Number of values
1330 real(kind_real),intent(in) :: vals(nval) !< Values
1331 
1332 ! Local variables
1333 integer :: ii
1334 character(len=1024) :: record
1335 
1336 ! Check
1337 if ((fld%geom%nx/=iter%geom%nx).or.(fld%geom%ny/=iter%geom%ny).or.(fld%geom%nz/=iter%geom%nz)) then
1338  write(record,*) 'qg_fields_getpoint: resolution inconsistency,',fld%geom%nx,'/',fld%geom%ny,'/',fld%geom%nz, &
1339 & ' and ',iter%geom%nx,'/',iter%geom%ny,'/',iter%geom%nz
1340  call abor1_ftn(record)
1341 endif
1342 if (fld%geom%nz/=nval) then
1343  write(record,*) 'qg_fields_getpoint: array sizes are different:',fld%geom%nz,'/',nval
1344  call abor1_ftn(record)
1345 endif
1346 
1347 ! Initialization
1348 ii = 0
1349 
1350 ! Get values
1351 if (allocated(fld%x)) then
1352  fld%x(iter%ilon,iter%ilat,:) = vals(ii+1:ii+fld%geom%nz)
1353  ii = ii+fld%geom%nz
1354 endif
1355 if (allocated(fld%q)) then
1356  fld%q(iter%ilon,iter%ilat,:) = vals(ii+1:ii+fld%geom%nz)
1357  ii = ii+fld%geom%nz
1358 endif
1359 if (allocated(fld%u)) then
1360  fld%u(iter%ilon,iter%ilat,:) = vals(ii+1:ii+fld%geom%nz)
1361  ii = ii+fld%geom%nz
1362 endif
1363 if (allocated(fld%v)) then
1364  fld%v(iter%ilon,iter%ilat,:) = vals(ii+1:ii+fld%geom%nz)
1365  ii = ii+fld%geom%nz
1366 endif
1367 
1368 end subroutine qg_fields_setpoint
1369 ! ------------------------------------------------------------------------------
1370 !> Serialize fields
1371 subroutine qg_fields_serialize(fld,vsize,vect_fld)
1372 
1373 implicit none
1374 
1375 ! Passed variables
1376 type(qg_fields),intent(in) :: fld !< Fields
1377 integer,intent(in) :: vsize !< Size
1378 real(kind_real),intent(inout) :: vect_fld(vsize) !< Vector
1379 
1380 ! Local variables
1381 integer :: ix,iy,iz,ind
1382 
1383 ! Initialize
1384 ind = 0
1385 
1386 ! Copy
1387 do iz=1,fld%geom%nz
1388  do iy=1,fld%geom%ny
1389  do ix=1,fld%geom%nx
1390  if (allocated(fld%x)) then
1391  ind = ind + 1
1392  vect_fld(ind) = fld%x(ix,iy,iz)
1393  endif
1394  if (allocated(fld%q)) then
1395  ind = ind + 1
1396  vect_fld(ind) = fld%q(ix,iy,iz)
1397  endif
1398  if (allocated(fld%u)) then
1399  ind = ind + 1
1400  vect_fld(ind) = fld%u(ix,iy,iz)
1401  endif
1402  if (allocated(fld%v)) then
1403  ind = ind + 1
1404  vect_fld(ind) = fld%v(ix,iy,iz)
1405  endif
1406  enddo
1407  enddo
1408 enddo
1409 
1410 if (fld%lbc) then
1411  do iz=1,fld%geom%nz
1412  ind = ind + 1
1413  vect_fld(ind) = fld%x_north(iz)
1414  ind = ind + 1
1415  vect_fld(ind) = fld%x_south(iz)
1416  do ix=1,fld%geom%nx
1417  ind = ind + 1
1418  vect_fld(ind) = fld%q_north(ix,iz)
1419  ind = ind + 1
1420  vect_fld(ind) = fld%q_south(ix,iz)
1421  enddo
1422  enddo
1423 endif
1424 
1425 end subroutine qg_fields_serialize
1426 ! ------------------------------------------------------------------------------
1427 !> Deserialize fields
1428 subroutine qg_fields_deserialize(self,vsize,vect_fld,index)
1429 
1430 implicit none
1431 
1432 ! Passed variables
1433 type(qg_fields),intent(inout) :: self !< Fields
1434 integer,intent(in) :: vsize !< Size
1435 real(kind_real),intent(in) :: vect_fld(vsize) !< Vector
1436 integer,intent(inout) :: index !< Index
1437 
1438 ! Local variables
1439 integer :: ix,iy,iz
1440 
1441 ! 3d field
1442 index = 1 + index
1443 do iz=1,self%geom%nz
1444  do iy=1,self%geom%ny
1445  do ix=1,self%geom%nx
1446  if (allocated(self%x)) then
1447  self%x(ix,iy,iz) = vect_fld(index)
1448  index = index+1
1449  endif
1450  if (allocated(self%q)) then
1451  self%x(ix,iy,iz) = vect_fld(index)
1452  index = index+1
1453  endif
1454  if (allocated(self%u)) then
1455  self%x(ix,iy,iz) = vect_fld(index)
1456  index = index+1
1457  endif
1458  if (allocated(self%v)) then
1459  self%x(ix,iy,iz) = vect_fld(index)
1460  index = index+1
1461  endif
1462  enddo
1463  enddo
1464 enddo
1465 
1466 ! Boundaries
1467 if (self%lbc) then
1468  do iz=1,self%geom%nz
1469  index = index + 1
1470  self%x_north(iz) = vect_fld(index)
1471  index = index + 1
1472  self%x_south(iz) = vect_fld(index)
1473  do ix=1,self%geom%nx
1474  index = index + 1
1475  self%q_north(ix,iz) = vect_fld(index)
1476  index = index + 1
1477  self%q_south(ix,iz) = vect_fld(index)
1478  enddo
1479  enddo
1480 endif
1481 
1482 index = index - 1
1483 
1484 end subroutine qg_fields_deserialize
1485 ! ------------------------------------------------------------------------------
1486 !> Complete missing fields consistently
1487 subroutine qg_fields_complete(self,var)
1488 
1489 implicit none
1490 
1491 ! Passed variables
1492 type(qg_fields),intent(inout) :: self !< Fields
1493 character(len=1),intent(in) :: var !< Reference variable ('x' or 'q')
1494 
1495 ! Local variables
1496 real(kind_real) :: x(self%geom%nx,self%geom%ny,self%geom%nz)
1497 real(kind_real) :: q(self%geom%nx,self%geom%ny,self%geom%nz)
1498 real(kind_real) :: u(self%geom%nx,self%geom%ny,self%geom%nz)
1499 real(kind_real) :: v(self%geom%nx,self%geom%ny,self%geom%nz)
1500 
1501 select case (var)
1502 case ('x')
1503  if (allocated(self%q)) then
1504  if (self%lbc) then
1505  call convert_x_to_q(self%geom,self%x,self%x_north,self%x_south,self%q)
1506  else
1507  call convert_x_to_q_tl(self%geom,self%x,self%q)
1508  endif
1509  endif
1510  if (allocated(self%u)) then
1511  if (self%lbc) then
1512  call convert_x_to_u(self%geom,self%x,self%x_north,self%x_south,self%u)
1513  else
1514  call convert_x_to_u_tl(self%geom,self%x,self%u)
1515  endif
1516  endif
1517  if (allocated(self%v)) then
1518  if (self%lbc) then
1519  call convert_x_to_v(self%geom,self%x,self%v)
1520  else
1521  call convert_x_to_v_tl(self%geom,self%x,self%v)
1522  endif
1523  endif
1524 case ('q')
1525  if (allocated(self%x).or.allocated(self%u).or.allocated(self%v)) then
1526  if (self%lbc) then
1527  call convert_q_to_x(self%geom,self%q,self%x_north,self%x_south,x)
1528  else
1529  call convert_q_to_x_tl(self%geom,self%q,x)
1530  endif
1531  if (allocated(self%x)) self%x = x
1532  if (allocated(self%u)) then
1533  if (self%lbc) then
1534  call convert_x_to_u(self%geom,self%x,self%x_north,self%x_south,self%u)
1535  else
1536  call convert_x_to_u_tl(self%geom,self%x,self%u)
1537  endif
1538  endif
1539  if (allocated(self%v)) then
1540  if (self%lbc) then
1541  call convert_x_to_v(self%geom,self%x,self%v)
1542  else
1543  call convert_x_to_v_tl(self%geom,self%x,self%v)
1544  endif
1545  endif
1546  endif
1547 case default
1548  call abor1_ftn('qg_fields_complete: wrong variable')
1549 end select
1550 
1551 end subroutine qg_fields_complete
1552 ! ------------------------------------------------------------------------------
1553 !> Check fields
1554 subroutine qg_fields_check(self)
1555 
1556 implicit none
1557 
1558 ! Passed variables
1559 type(qg_fields),intent(in) :: self !< Fields
1560 
1561 ! Local variables
1562 logical :: bad
1563 character(len=1024) :: record
1564 
1565 ! Initialization
1566 bad = .false.
1567 
1568 ! Check 3d field
1569 bad = bad.or.(.not.(allocated(self%x).or.allocated(self%q).or.allocated(self%u).or.allocated(self%v)))
1570 if (allocated(self%x)) then
1571  bad = bad.or.(size(self%x,1)/=self%geom%nx)
1572  bad = bad.or.(size(self%x,2)/=self%geom%ny)
1573  bad = bad.or.(size(self%x,3)/=self%geom%nz)
1574 endif
1575 if (allocated(self%q)) then
1576  bad = bad.or.(size(self%q,1)/=self%geom%nx)
1577  bad = bad.or.(size(self%q,2)/=self%geom%ny)
1578  bad = bad.or.(size(self%q,3)/=self%geom%nz)
1579 endif
1580 if (allocated(self%u)) then
1581  bad = bad.or.(size(self%u,1)/=self%geom%nx)
1582  bad = bad.or.(size(self%u,2)/=self%geom%ny)
1583  bad = bad.or.(size(self%u,3)/=self%geom%nz)
1584 endif
1585 if (allocated(self%v)) then
1586  bad = bad.or.(size(self%v,1)/=self%geom%nx)
1587  bad = bad.or.(size(self%v,2)/=self%geom%ny)
1588  bad = bad.or.(size(self%v,3)/=self%geom%nz)
1589 endif
1590 
1591 ! Check boundaries
1592 if (self%lbc) then
1593  bad = bad.or.(.not.allocated(self%x_north))
1594  bad = bad.or.(.not.allocated(self%x_south))
1595  bad = bad.or.(.not.allocated(self%q_north))
1596  bad = bad.or.(.not.allocated(self%q_south))
1597  bad = bad.or.(size(self%x_north)/=self%geom%nz)
1598  bad = bad.or.(size(self%x_south)/=self%geom%nz)
1599  bad = bad.or.(size(self%q_north,1)/=self%geom%nx)
1600  bad = bad.or.(size(self%q_north,2)/=self%geom%nz)
1601  bad = bad.or.(size(self%q_south,1)/=self%geom%nx)
1602  bad = bad.or.(size(self%q_south,2)/=self%geom%nz)
1603 else
1604  bad = bad.or.allocated(self%x_north)
1605  bad = bad.or.allocated(self%x_south)
1606  bad = bad.or.allocated(self%q_north)
1607  bad = bad.or.allocated(self%q_south)
1608 endif
1609 
1610 if (bad) then
1611  call fckit_log%info('qg_fields_check: field not consistent')
1612  write(record,*) ' nx,ny,nz,lbc = ',self%geom%nx,self%geom%ny,self%geom%nz,self%lbc
1613  call fckit_log%info(record)
1614  if (allocated(self%x)) then
1615  write(record,*) ' shape(x) = ',shape(self%x)
1616  call fckit_log%info(record)
1617  endif
1618  if (allocated(self%q)) then
1619  write(record,*) ' shape(q) = ',shape(self%q)
1620  call fckit_log%info(record)
1621  endif
1622  if (allocated(self%u)) then
1623  write(record,*) ' shape(u) = ',shape(self%u)
1624  call fckit_log%info(record)
1625  endif
1626  if (allocated(self%v)) then
1627  write(record,*) ' shape(v) = ',shape(self%v)
1628  call fckit_log%info(record)
1629  endif
1630  if (self%lbc) then
1631  write(record,*) ' shape(x_north) = ',shape(self%x_north)
1632  call fckit_log%info(record)
1633  write(record,*) ' shape(x_south) = ',shape(self%x_south)
1634  call fckit_log%info(record)
1635  write(record,*) ' shape(q_north) = ',shape(self%q_north)
1636  call fckit_log%info(record)
1637  write(record,*) ' shape(q_south) = ',shape(self%q_south)
1638  call fckit_log%info(record)
1639  call abor1_ftn ('qg_fields_check: field not consistent')
1640  endif
1641 endif
1642 
1643 end subroutine qg_fields_check
1644 ! ------------------------------------------------------------------------------
1645 !> Check fields resolution
1646 subroutine qg_fields_check_resolution(fld1,fld2)
1647 
1648 implicit none
1649 
1650 ! Passed variables
1651 type(qg_fields),intent(in) :: fld1 !< First fields
1652 type(qg_fields),intent(in) :: fld2 !< Second fields
1653 
1654 ! Local variables
1655 character(len=1024) :: record
1656 
1657 ! Check fields consistency
1658 if ((fld1%geom%nx/=fld2%geom%nx).or.(fld1%geom%ny/=fld2%geom%ny).or.(fld1%geom%nz/=fld2%geom%nz)) then
1659  write(record,*) 'qg_fields_check_resolution: resolution inconsistency, ',fld1%geom%nx,'/',fld1%geom%ny,'/',fld1%geom%nz, &
1660 & ' and ',fld2%geom%nx,'/',fld2%geom%ny,'/',fld2%geom%nz
1661  call abor1_ftn(record)
1662 endif
1663 
1664 ! Check fields independently
1665 call qg_fields_check(fld1)
1666 call qg_fields_check(fld2)
1667 
1668 end subroutine qg_fields_check_resolution
1669 ! ------------------------------------------------------------------------------
1670 end module qg_fields_mod
Fortran interface to Variables.
real(kind_real), parameter, public domain_meridional
Model domain in meridional direction (m)
subroutine, public convert_q_to_x_tl(geom, q, x)
Convert potential vorticity to streamfunction - tangent linear.
subroutine, public convert_q_to_x(geom, q, x_north, x_south, x)
Convert potential vorticity to streamfunction.
subroutine, public convert_x_to_q_tl(geom, x, q)
Convert streamfunction to potential vorticity - tangent linear.
subroutine, public convert_x_to_q(geom, x, x_north, x_south, q)
Convert streamfunction to potential vorticity.
subroutine, public convert_x_to_u(geom, x, x_north, x_south, u)
Convert streafunction to zonal wind.
subroutine, public convert_x_to_u_tl(geom, x, u)
Convert streafunction to zonal wind - tangent Linear.
subroutine, public convert_x_to_v(geom, x, v)
Convert streafunction to meridional wind.
subroutine, public convert_x_to_v_tl(geom, x, v)
Convert streafunction to meridional wind - tangent Linear.
subroutine, public qg_fields_set_atlas(self, vars, afieldset)
Set ATLAS field.
subroutine, public qg_fields_copy_lbc(self, other)
Copy fields LBC.
subroutine, public qg_fields_deserialize(self, vsize, vect_fld, index)
Deserialize fields.
subroutine, public qg_fields_zero(self)
Set fields to zero.
subroutine, public qg_fields_axpy(self, zz, rhs)
Apply axpy operator to fields.
integer, parameter rseed
Random seed (for reproducibility)
subroutine, public qg_fields_to_atlas(self, vars, afieldset)
Convert fields to ATLAS.
subroutine, public qg_fields_delete(self)
Delete fields.
subroutine, public qg_fields_dirac(self, f_conf)
Set fields to Diracs.
subroutine, public qg_fields_getpoint(fld, iter, nval, vals)
Get points from fields.
subroutine, public qg_fields_gpnorm(fld, vpresent, vmin, vmax, vrms)
Fields statistics.
subroutine, public qg_fields_rms(fld, prms)
Fields RMS.
subroutine, public qg_fields_self_schur(self, rhs)
Schur product of fields.
subroutine, public qg_fields_self_mul(self, zz)
Multiply fields by a scalar.
type(registry_t), public qg_fields_registry
Linked list interface - defines registry_t type.
subroutine, public qg_fields_create(self, geom, vars, lbc)
Linked list implementation.
subroutine, public qg_fields_from_atlas(self, vars, afieldset)
Get fields from ATLAS.
subroutine, public qg_fields_complete(self, var)
Complete missing fields consistently.
subroutine, public qg_fields_analytic_init(fld, f_conf, vdate)
Analytic initialization of fields.
subroutine, public qg_fields_check_resolution(fld1, fld2)
Check fields resolution.
subroutine, public qg_fields_self_sub(self, rhs)
Subtract fields.
subroutine, public qg_fields_check(self)
Check fields.
subroutine, public qg_fields_copy(self, other)
Copy fields.
subroutine, public qg_fields_dot_prod(fld1, fld2, zprod)
Compute dot product for fields.
subroutine, public qg_fields_write_file(fld, f_conf, vdate)
Write fields to file.
subroutine, public qg_fields_add_incr(self, rhs)
Add increment to fields.
subroutine, public qg_fields_self_add(self, rhs)
Add fields.
subroutine, public qg_fields_serialize(fld, vsize, vect_fld)
Serialize fields.
subroutine, public qg_fields_read_file(fld, f_conf, vdate)
Read fields from file.
subroutine, public qg_fields_sizes(fld, nx, ny, nz)
Get fields geometry.
subroutine, public qg_fields_setpoint(fld, iter, nval, vals)
Set fields values at a specified gridpoint.
subroutine, public qg_fields_lbc(fld, lbc)
Get LBC presence.
subroutine, public qg_fields_random(self, var)
Generate random fields.
subroutine, public qg_fields_ones(self)
Set fields to ones.
subroutine, public qg_fields_diff_incr(lhs, fld1, fld2)
Compute increment from the difference of two fields.
subroutine, public qg_fields_create_from_other(self, other, geom)
Create fields from another one.
subroutine, public qg_fields_change_resol(fld, rhs)
Change fields resolution.
subroutine, public qg_interp_trilinear(geom, lon, lat, z, field, val)
Trilinear interpolation.
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.
character(len=2 *length) function, public genfilename(f_conf, length, vdate)
Generate filename.