11 use atlas_module,
only: atlas_field, atlas_fieldset, atlas_real
12 use fckit_configuration_module,
only: fckit_configuration
15 use fckit_log_module,
only: fckit_log
18 use missing_values_mod
53 real(kind_real),
allocatable :: x(:,:,:)
54 real(kind_real),
allocatable :: q(:,:,:)
55 real(kind_real),
allocatable :: u(:,:,:)
56 real(kind_real),
allocatable :: v(:,:,:)
57 real(kind_real),
allocatable :: x_north(:)
58 real(kind_real),
allocatable :: x_south(:)
59 real(kind_real),
allocatable :: q_north(:,:)
60 real(kind_real),
allocatable :: q_south(:,:)
63 #define LISTED_TYPE qg_fields
66 #include "oops/util/linkedList_i.f"
76 #include "oops/util/linkedList_c.f"
85 type(
qg_geom),
target,
intent(in) :: geom
87 logical,
intent(in) :: lbc
90 character(len=1024) :: record
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))
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))
126 type(
qg_geom),
target,
intent(in) :: geom
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))
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))
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)
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
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
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
231 type(fckit_configuration),
intent(in) :: f_conf
235 integer,
allocatable :: ixdir(:),iydir(:),izdir(:)
236 character(len=1) :: var
237 character(len=:),
allocatable :: str
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')
248 allocate(ixdir(ndir))
249 allocate(iydir(ndir))
250 allocate(izdir(ndir))
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)
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')
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
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
275 call abor1_ftn(
'qg_fields_dirac: wrong variable')
291 character(len=1),
intent(in) :: var
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)
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)
305 call abor1_ftn(
'qg_fields_random: wrong variable')
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
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
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
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
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
423 real(kind_real),
intent(in) :: zz
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
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
449 real(kind_real),
intent(in) :: zz
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
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
503 real(kind_real),
intent(out) :: zprod
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)
530 if ((self%geom%nx==rhs%geom%nx).and.(self%geom%ny==rhs%geom%ny).and.(self%geom%nz==rhs%geom%nz))
then
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
538 call abor1_ftn(
'qg_fields_add_incr: not coded for low res increment yet')
560 if ((fld1%geom%nx==lhs%geom%nx).and.(fld1%geom%ny==lhs%geom%ny).and.(fld1%geom%nz==lhs%geom%nz))
then
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
568 call abor1_ftn(
'qg_fields_diff_incr: not coded for low res increment yet')
583 real(kind_real),
allocatable,
dimension(:,:,:) :: q1, q2
585 if ((fld%geom%nx==rhs%geom%nx).and.(fld%geom%ny==rhs%geom%ny).and.(fld%geom%nz==rhs%geom%nz))
then
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))
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))
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))
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))
616 allocate(q1(rhs%geom%nx,rhs%geom%ny,rhs%geom%nz))
617 allocate(q2(fld%geom%nx,fld%geom%ny,fld%geom%nz))
619 q1(:,iy,:) = rhs%q_south
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) )
626 fld%q_south = q2(:,1,:)
628 q1(:,iy,:) = rhs%q_north
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))
635 fld%q_north = q2(:,1,:)
637 fld%x_north = rhs%x_north
638 fld%x_south = rhs%x_south
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
658 type(fckit_configuration),
intent(in) :: f_conf
659 type(datetime),
intent(inout) :: vdate
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
665 character(len=20) :: sdate
666 character(len=1024) :: record,filename
667 character(len=:),
allocatable :: str
674 if (f_conf%has(
"read_from_file"))
call f_conf%get_or_die(
"read_from_file",iread)
678 call fckit_log%warning(
'qg_fields_read_file: inventing field')
684 call f_conf%get_or_die(
"filename",str)
685 call swap_name_member(f_conf, str)
687 call fckit_log%info(
'qg_fields_read_file: opening '//trim(filename))
693 call ncerr(nf90_open(trim(filename),nf90_nowrite,ncid))
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))
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))
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')
715 call ncerr(nf90_get_att(ncid,nf90_global,
'bc',bc))
721 call abor1_ftn(
'qg_fields_read_file: wrong bc value')
723 call ncerr(nf90_get_att(ncid,nf90_global,
'sdate',sdate))
726 if ((.not.lbc).and.fld%lbc)
call abor1_ftn(
'qg_fields_read_file: LBC are missing in NetCDF file')
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))
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))
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))
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))
753 call ncerr(nf90_close(ncid))
756 call fckit_log%info(
'qg_fields_read_file: validity date is '//sdate)
757 call datetime_set(sdate,vdate)
772 type(fckit_configuration),
intent(in) :: f_conf
773 type(datetime),
intent(in) :: vdate
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
779 character(len=20) :: sdate
780 character(len=1024) :: filename
789 call vars%push_back(
'x')
790 call vars%push_back(
'q')
791 call vars%push_back(
'u')
792 call vars%push_back(
'v')
795 if (
allocated(fld%x))
then
798 elseif (
allocated(fld%q))
then
802 call abor1_ftn(
'qg_fields_write_file: x or q required')
807 call fckit_log%info(
'qg_fields_write_file: writing '//trim(filename))
810 call datetime_to_string(vdate,sdate)
813 call ncerr(nf90_create(trim(filename),or(nf90_clobber,nf90_64bit_offset),ncid))
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))
822 call ncerr(nf90_put_att(ncid,nf90_global,
'bc',1))
824 call ncerr(nf90_put_att(ncid,nf90_global,
'bc',0))
826 call ncerr(nf90_put_att(ncid,nf90_global,
'sdate',sdate))
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)))
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)))
854 call ncerr(nf90_enddef(ncid))
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))
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))
874 call ncerr(nf90_close(ncid))
888 type(fckit_configuration),
intent(in) :: f_conf
889 type(datetime),
intent(inout) :: vdate
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
900 if (f_conf%has(
"analytic_init"))
then
901 call f_conf%get_or_die(
"analytic_init",str)
904 ic =
'baroclinic-instability'
906 call fckit_log%warning(
'qg_fields_analytic_init: '//trim(ic))
909 call f_conf%get_or_die(
"date",str)
911 call fckit_log%info(
'qg_fields_analytic_init: validity date is '//sdate)
912 call datetime_set(sdate,vdate)
915 if (.not.fld%lbc)
call abor1_ftn(
'qg_fields_analytic_init: boundaries required')
918 allocate(x(fld%geom%nx,fld%geom%ny,fld%geom%nz))
919 allocate(q(fld%geom%nx,fld%geom%ny,fld%geom%nz))
922 select case (trim(ic))
923 case (
'baroclinic-instability')
934 case (
'large-vortices')
939 call large_vortices(fld%geom%x(ix),fld%geom%y(iy),fld%geom%z(iz),
'x',x(ix,iy,iz))
943 call large_vortices(0.0_kind_real,0.0_kind_real,fld%geom%z(iz),
'x',fld%x_south(iz))
945 case (
'uniform_field')
947 call f_conf%get_or_die(
"uval",uval)
950 call abor1_ftn(
'qg_fields_analytic_init: unknown initialization')
957 fld%q_south(ix,iz) = 2.0*q(ix,1,iz)-q(ix,2,iz)
960 fld%q_north(ix,iz) = 2.0*q(ix,fld%geom%ny,iz)-q(ix,fld%geom%ny-1,iz)
965 if (
allocated(x))
then
968 elseif (
allocated(q))
then
972 call abor1_ftn(
'qg_fields_analytic_init: x or q required')
987 integer,
intent(inout) :: vpresent(6)
988 real(kind_real),
intent(inout) :: vmin(6)
989 real(kind_real),
intent(inout) :: vmax(6)
990 real(kind_real),
intent(inout) :: vrms(6)
1002 if (
allocated(fld%x))
then
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))
1008 if (
allocated(fld%q))
then
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))
1014 if (
allocated(fld%u))
then
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))
1020 if (
allocated(fld%v))
then
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))
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))
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))
1051 real(kind_real),
intent(out) :: prms
1060 prms = 0.0_kind_real
1061 norm = 0.0_kind_real
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)
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)
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)
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)
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)
1088 prms = sqrt(prms/norm)
1099 integer,
intent(out) :: nx
1100 integer,
intent(out) :: ny
1101 integer,
intent(out) :: nz
1117 integer,
intent(out) :: lbc
1138 type(atlas_fieldset),
intent(inout) :: afieldset
1142 character(len=1024) :: fieldname
1143 type(atlas_field) :: afield
1146 do jvar=1,vars%nvars()
1147 fieldname = vars%variable(jvar)
1148 if (.not.afieldset%has_field(trim(fieldname)))
then
1150 afield = self%geom%afunctionspace%create_field(name=trim(fieldname),kind=atlas_real(kind_real),levels=self%geom%nz)
1153 call afieldset%add(afield)
1170 type(atlas_fieldset),
intent(inout) :: afieldset
1173 integer :: jvar,ix,iy,iz,inode
1174 real(kind_real),
pointer :: ptr(:,:)
1175 character(len=1024) :: fieldname
1176 type(atlas_field) :: afield
1179 do jvar=1,vars%nvars()
1180 fieldname = vars%variable(jvar)
1181 if (afieldset%has_field(trim(fieldname)))
then
1183 afield = afieldset%field(trim(fieldname))
1186 afield = self%geom%afunctionspace%create_field(name=trim(fieldname),kind=atlas_real(kind_real),levels=self%geom%nz)
1189 call afieldset%add(afield)
1193 call afield%data(ptr)
1194 do iz=1,self%geom%nz
1196 do iy=1,self%geom%ny
1197 do ix=1,self%geom%nx
1199 select case (trim(fieldname))
1201 ptr(iz,inode) = self%x(ix,iy,iz)
1203 ptr(iz,inode) = self%q(ix,iy,iz)
1205 ptr(iz,inode) = self%u(ix,iy,iz)
1207 ptr(iz,inode) = self%v(ix,iy,iz)
1209 call abor1_ftn(
'qg_fields_to_atlas: wrong variable')
1229 type(atlas_fieldset),
intent(inout) :: afieldset
1232 integer :: jvar,ix,iy,iz,inode
1233 real(kind_real),
pointer :: ptr(:,:)
1234 character(len=1024) :: fieldname
1235 type(atlas_field) :: afield
1238 do jvar=1,vars%nvars()
1240 fieldname = vars%variable(jvar)
1241 afield = afieldset%field(trim(fieldname))
1244 call afield%data(ptr)
1245 do iz=1,self%geom%nz
1247 do iy=1,self%geom%ny
1248 do ix=1,self%geom%nx
1250 select case (trim(fieldname))
1252 self%x(ix,iy,iz) = ptr(iz,inode)
1254 self%q(ix,iy,iz) = ptr(iz,inode)
1256 self%u(ix,iy,iz) = ptr(iz,inode)
1258 self%v(ix,iy,iz) = ptr(iz,inode)
1260 call abor1_ftn(
'qg_fields_to_atlas: wrong variable')
1280 integer,
intent(in) :: nval
1281 real(kind_real),
intent(inout) :: vals(nval)
1285 character(len=1024) :: record
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)
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)
1302 if (
allocated(fld%x))
then
1303 vals(ii+1:ii+fld%geom%nz) = fld%x(iter%ilon,iter%ilat,:)
1306 if (
allocated(fld%q))
then
1307 vals(ii+1:ii+fld%geom%nz) = fld%q(iter%ilon,iter%ilat,:)
1310 if (
allocated(fld%u))
then
1311 vals(ii+1:ii+fld%geom%nz) = fld%u(iter%ilon,iter%ilat,:)
1314 if (
allocated(fld%v))
then
1315 vals(ii+1:ii+fld%geom%nz) = fld%v(iter%ilon,iter%ilat,:)
1329 integer,
intent(in) :: nval
1330 real(kind_real),
intent(in) :: vals(nval)
1334 character(len=1024) :: record
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)
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)
1351 if (
allocated(fld%x))
then
1352 fld%x(iter%ilon,iter%ilat,:) = vals(ii+1:ii+fld%geom%nz)
1355 if (
allocated(fld%q))
then
1356 fld%q(iter%ilon,iter%ilat,:) = vals(ii+1:ii+fld%geom%nz)
1359 if (
allocated(fld%u))
then
1360 fld%u(iter%ilon,iter%ilat,:) = vals(ii+1:ii+fld%geom%nz)
1363 if (
allocated(fld%v))
then
1364 fld%v(iter%ilon,iter%ilat,:) = vals(ii+1:ii+fld%geom%nz)
1377 integer,
intent(in) :: vsize
1378 real(kind_real),
intent(inout) :: vect_fld(vsize)
1381 integer :: ix,iy,iz,ind
1390 if (
allocated(fld%x))
then
1392 vect_fld(ind) = fld%x(ix,iy,iz)
1394 if (
allocated(fld%q))
then
1396 vect_fld(ind) = fld%q(ix,iy,iz)
1398 if (
allocated(fld%u))
then
1400 vect_fld(ind) = fld%u(ix,iy,iz)
1402 if (
allocated(fld%v))
then
1404 vect_fld(ind) = fld%v(ix,iy,iz)
1413 vect_fld(ind) = fld%x_north(iz)
1415 vect_fld(ind) = fld%x_south(iz)
1418 vect_fld(ind) = fld%q_north(ix,iz)
1420 vect_fld(ind) = fld%q_south(ix,iz)
1434 integer,
intent(in) :: vsize
1435 real(kind_real),
intent(in) :: vect_fld(vsize)
1436 integer,
intent(inout) :: 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)
1450 if (
allocated(self%q))
then
1451 self%x(ix,iy,iz) = vect_fld(index)
1454 if (
allocated(self%u))
then
1455 self%x(ix,iy,iz) = vect_fld(index)
1458 if (
allocated(self%v))
then
1459 self%x(ix,iy,iz) = vect_fld(index)
1468 do iz=1,self%geom%nz
1470 self%x_north(iz) = vect_fld(index)
1472 self%x_south(iz) = vect_fld(index)
1473 do ix=1,self%geom%nx
1475 self%q_north(ix,iz) = vect_fld(index)
1477 self%q_south(ix,iz) = vect_fld(index)
1493 character(len=1),
intent(in) :: var
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)
1503 if (
allocated(self%q))
then
1505 call convert_x_to_q(self%geom,self%x,self%x_north,self%x_south,self%q)
1510 if (
allocated(self%u))
then
1512 call convert_x_to_u(self%geom,self%x,self%x_north,self%x_south,self%u)
1517 if (
allocated(self%v))
then
1525 if (
allocated(self%x).or.
allocated(self%u).or.
allocated(self%v))
then
1527 call convert_q_to_x(self%geom,self%q,self%x_north,self%x_south,x)
1531 if (
allocated(self%x)) self%x = x
1532 if (
allocated(self%u))
then
1534 call convert_x_to_u(self%geom,self%x,self%x_north,self%x_south,self%u)
1539 if (
allocated(self%v))
then
1548 call abor1_ftn(
'qg_fields_complete: wrong variable')
1563 character(len=1024) :: record
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)
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)
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)
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)
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)
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)
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)
1618 if (
allocated(self%q))
then
1619 write(record,*)
' shape(q) = ',shape(self%q)
1620 call fckit_log%info(record)
1622 if (
allocated(self%u))
then
1623 write(record,*)
' shape(u) = ',shape(self%u)
1624 call fckit_log%info(record)
1626 if (
allocated(self%v))
then
1627 write(record,*)
' shape(v) = ',shape(self%v)
1628 call fckit_log%info(record)
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')
1655 character(len=1024) :: record
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)
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.