11 use fckit_configuration_module,
only: fckit_configuration
12 use fckit_log_module,
only: fckit_log
35 real(kind_real) :: sigma
36 real(kind_real),
allocatable :: sqrt_zonal(:)
37 real(kind_real),
allocatable :: sqrt_merid(:,:)
38 real(kind_real),
allocatable :: sqrt_vert(:,:)
39 real(kind_real),
allocatable :: norm(:,:)
42 real(kind_real),
parameter ::
eps_ad = 1.0e-10
44 #define LISTED_TYPE qg_error_covariance_config
47 #include "oops/util/linkedList_i.f"
57 #include "oops/util/linkedList_c.f"
66 type(fckit_configuration),
intent(in) :: f_conf
67 type(
qg_geom),
intent(in) :: geom
70 integer :: ix,iy,jy,ky,iz,jz,kz,info
71 real(kind_real) :: horizontal_length_scale,vertical_length_scale
72 real(kind_real) :: distx,condition_number,threshold
73 real(kind_real),
allocatable :: struct_fn(:),workx(:)
74 real(kind_real),
allocatable :: evalsy(:),worky(:),evectsy(:,:),revalsy(:)
75 real(kind_real),
allocatable :: evalsz(:),workz(:),evectsz(:,:),revalsz(:)
76 real(kind_real),
allocatable :: norm(:,:)
77 character(len=160) :: record
82 call f_conf%get_or_die(
"standard_deviation",self%sigma)
83 call f_conf%get_or_die(
"horizontal_length_scale",horizontal_length_scale)
84 call f_conf%get_or_die(
"vertical_length_scale",vertical_length_scale)
85 call f_conf%get_or_die(
"maximum_condition_number",condition_number)
88 if (mod(geom%nx,2)/=0)
then
89 write(record,*)
'qg_error_covariance_setup: number of zonal gridpoints nx=',geom%nx
90 call fckit_log%error(record)
91 call abor1_ftn(
'qg_error_covariance_setup: odd number of zonal grid points')
100 allocate(self%sqrt_merid(geom%ny,geom%ny))
101 allocate(self%sqrt_vert(geom%nz,geom%nz))
102 allocate(self%norm(geom%ny,geom%nz))
103 allocate(struct_fn(geom%nx))
104 allocate(workx(geom%nx+2))
105 allocate(self%sqrt_zonal(0:geom%nx/2))
106 allocate(evectsy(geom%ny,geom%ny))
107 allocate(evalsy(geom%ny))
108 allocate(worky((geom%ny+3)*geom%ny))
109 allocate(revalsy(geom%ny))
110 allocate(evectsz(geom%nz,geom%nz))
111 allocate(evalsz(geom%nz))
112 allocate(workz((geom%nz+3)*geom%nz))
113 allocate(revalsz(geom%nz))
114 allocate(norm(geom%ny,geom%nz))
119 distx = geom%x(ix)-geom%x(1)
121 struct_fn(ix) = exp(-0.5*(distx/horizontal_length_scale)**2)
123 call fft_fwd(geom%nx,struct_fn,workx)
124 workx = workx/real(geom%nx,kind_real)
131 threshold = max(threshold,workx(1+2*ix))
133 threshold = threshold/condition_number
135 self%sqrt_zonal(ix) = sqrt(max(threshold,real(geom%nx,kind_real)*workx(1+2*ix)))
141 evectsy(iy,jy) = exp(-0.5*((geom%y(iy)-geom%y(jy))/horizontal_length_scale)**2)
142 evectsy(jy,iy) = evectsy(iy,jy)
147 call dsyev(
'V',
'L',geom%ny,evectsy,geom%ny,evalsy,worky,
size(worky),info)
149 write(record,*)
'qg_error_covariance_setup: dsyev returns info=',info
150 call fckit_log%error(record)
151 call abor1_ftn (
'qg_error_covariance_setup: dsyev failed to sqrt matrix')
155 threshold = maxval(evalsy)/condition_number
157 evalsy(iy) = sqrt(max(threshold,evalsy(iy)))
158 revalsy(iy) = 1.0/evalsy(iy)
162 self%sqrt_merid(iy,jy) = 0.0
164 self%sqrt_merid(iy,jy) = self%sqrt_merid(iy,jy)+evectsy(iy,ky)*evectsy(jy,ky)*evalsy(ky)
172 evectsz(iz,jz) = exp(-0.5*((geom%z(iz)-geom%z(jz))/vertical_length_scale)**2)
173 evectsz(jz,iz) = evectsz(iz,jz)
178 call dsyev(
'V',
'L',geom%nz,evectsz,geom%nz,evalsz,workz,
size(workz),info)
180 write(record,*)
'qg_error_covariance_setup: dsyev returns info=',info
181 call fckit_log%error(record)
182 call abor1_ftn (
'qg_error_covariance_setup: dsyev failed to sqrt matrix')
186 threshold = maxval(evalsz)/condition_number
188 evalsz(iz) = sqrt(max(threshold,evalsz(iz)))
189 revalsz(iz) = 1.0/evalsz(iz)
193 self%sqrt_vert(iz,jz) = 0.0
195 self%sqrt_vert(iz,jz) = self%sqrt_vert(iz,jz)+evectsz(iz,kz)*evectsz(jz,kz)*evalsz(kz)
202 call vars%push_back(
'x')
209 fld_in%x(1,iy,iz) = 1.0
211 norm(iy,iz) = 1.0/sqrt(fld_out%x(1,iy,iz))
214 self%norm = norm*self%sigma
217 deallocate(struct_fn)
242 deallocate(self%sqrt_zonal)
243 deallocate(self%sqrt_merid)
244 deallocate(self%sqrt_vert)
245 deallocate(self%norm)
303 real(kind_real),
intent(inout) :: fld(self%nx,self%ny,self%nz)
306 integer :: iy,iz,m,iri
307 real(kind_real) :: zfour(self%nx+2)
311 call fft_fwd(self%nx,fld(:,iy,iz),zfour)
315 zfour(2*m+iri) = zfour(2*m+iri)*self%sqrt_zonal(m)
319 call fft_inv(self%nx,zfour,fld(:,iy,iz))
332 real(kind_real),
intent(inout) :: fld(self%nx,self%ny,self%nz)
336 real(kind_real),
allocatable :: arr_in(:),arr_out(:)
342 allocate(arr_in(self%ny))
343 allocate(arr_out(self%ny))
346 arr_in = fld(ix,:,iz)
349 call dsymv(
'L',self%ny,1.0_kind_real,self%sqrt_merid,self%ny,arr_in,1,0.0_kind_real,arr_out,1)
352 fld(ix,:,iz) = arr_out
370 real(kind_real),
intent(inout) :: fld(self%nx,self%ny,self%nz)
374 real(kind_real),
allocatable :: arr_in(:),arr_out(:)
380 allocate(arr_in(self%nz))
381 allocate(arr_out(self%nz))
384 arr_in = fld(ix,iy,:)
387 call dsymv(
'L',self%nz,1.0_kind_real,self%sqrt_vert,self%nz,arr_in,1,0.0_kind_real,arr_out,1)
390 fld(ix,iy,:) = arr_out
415 if (.not.
allocated(fld_in%x))
call abor1_ftn(
"qg_error_covariance_sqrt_mult: x required as input")
416 if (.not.
allocated(fld_out%x))
call abor1_ftn(
"qg_error_covariance_sqrt_mult: x required as output")
432 do ix=1,fld_out%geom%nx
433 fld_out%x(ix,:,:) = fld_out%x(ix,:,:)*self%norm
438 fld_out%x = fld_out%x*self%sigma
456 if (.not.
allocated(fld_in%x))
call abor1_ftn(
"qg_error_covariance_sqrt_mult: x required as input")
457 if (.not.
allocated(fld_out%x))
call abor1_ftn(
"qg_error_covariance_sqrt_mult: x required as output")
463 fld_out%x = fld_out%x*self%sigma
467 do ix=1,fld_out%geom%nx
468 fld_out%x(ix,:,:) = fld_out%x(ix,:,:)*self%norm
Fortran module for Eigen FFTs.
subroutine, public fft_inv(kk, pf, pg)
subroutine, public fft_fwd(kk, pg, pf)
Fortran interface to Variables.
real(kind_real), parameter, public domain_zonal
Model domain in zonal direction (m)
type(registry_t), public qg_error_covariance_registry
Linked list interface - defines registry_t type.
subroutine qg_error_covariance_sqrt_mult_vertical(self, fld)
Multiply by error covariance matrix square-root - vertical part.
subroutine qg_error_covariance_sqrt_mult_ad(self, fld_in, fld_out)
Multiply by error covariance matrix square-root - adjoint.
subroutine, public qg_error_covariance_setup(self, f_conf, geom)
Linked list implementation.
subroutine qg_error_covariance_sqrt_mult_zonal(self, fld)
Multiply by error covariance matrix square-root, zonal part.
subroutine, public qg_error_covariance_delete(self)
Delete error covariance matrix.
subroutine, public qg_error_covariance_randomize(self, fld_out)
Randomize error covariance.
subroutine qg_error_covariance_sqrt_mult(self, fld_in, fld_out)
Multiply by error covariance matrix square-root.
real(kind_real), parameter eps_ad
Epsilon value for adjoint tests.
subroutine qg_error_covariance_sqrt_mult_meridional(self, fld)
Multiply by error covariance matrix square-root - meridional part.
subroutine, public qg_error_covariance_mult(self, fld_in, fld_out)
Multiply by error covariance matrix.
subroutine, public qg_fields_zero(self)
Set fields to zero.
subroutine, public qg_fields_delete(self)
Delete fields.
subroutine, public qg_fields_create(self, geom, vars, lbc)
Linked list implementation.
subroutine, public qg_fields_copy(self, other)
Copy fields.
subroutine, public qg_fields_random(self, var)
Generate random fields.
subroutine, public qg_fields_create_from_other(self, other, geom)
Create fields from another one.