FV3 Bundle
c_qg_bmatrix.f90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 ! In applying this licence, ECMWF does not waive the privileges and immunities
6 ! granted to it by virtue of its status as an intergovernmental organisation nor
7 ! does it submit to any jurisdiction.
8 
9 ! ------------------------------------------------------------------------------
10 
11 !> Setup for the QG model's background error covariance matrix
12 
13 subroutine c_qg_b_setup(c_key_self, c_conf, c_key_geom) &
14  & bind(c,name='qg_b_setup_f90')
15 
16 use iso_c_binding
18 use qg_geom_mod
19 
20 implicit none
21 integer(c_int), intent(inout) :: c_key_self !< The background covariance structure
22 type(c_ptr), intent(in) :: c_conf !< The configuration
23 integer(c_int), intent(in) :: c_key_geom !< Geometry
24 type(qg_3d_covar_config), pointer :: self
25 type(qg_geom), pointer :: geom
26 
27 call qg_geom_registry%get(c_key_geom, geom)
28 call qg_3d_cov_registry%init()
29 call qg_3d_cov_registry%add(c_key_self)
30 call qg_3d_cov_registry%get(c_key_self, self)
31 
32 call qg_3d_covar_setup(c_conf, geom, self)
33 
34 end subroutine c_qg_b_setup
35 
36 ! ------------------------------------------------------------------------------
37 !> Delete for the QG model's background error covariance matrix
38 
39 subroutine c_qg_b_delete(c_key_self) bind (c,name='qg_b_delete_f90')
40 
41 use iso_c_binding
43 
44 implicit none
45 integer(c_int), intent(inout) :: c_key_self !< The background covariance structure
46 type(qg_3d_covar_config), pointer :: self
47 
48 call qg_3d_cov_registry%get(c_key_self,self)
49 call qg_3d_covar_delete(self)
50 call qg_3d_cov_registry%remove(c_key_self)
51 
52 end subroutine c_qg_b_delete
53 
54 ! ------------------------------------------------------------------------------
55 
56 !> Multiply streamfunction by inverse of covariance
57 
58 subroutine c_qg_b_inv_mult(c_key_conf, c_key_in, c_key_out) bind(c,name='qg_b_invmult_f90')
59 
60 use iso_c_binding
62 use qg_fields
63 use kinds
64 
65 implicit none
66 integer(c_int), intent(in) :: c_key_conf !< covar config structure
67 integer(c_int), intent(in) :: c_key_in !< Streamfunction: psi
68 integer(c_int), intent(in) :: c_key_out !< Streamfunction: psi
69 type(qg_3d_covar_config), pointer :: conf
70 type(qg_field), pointer :: xin
71 type(qg_field), pointer :: xout
72 real(kind=kind_real), allocatable :: xctl(:,:,:) ! Control vector
73 
74 call qg_3d_cov_registry%get(c_key_conf,conf)
75 call qg_field_registry%get(c_key_in,xin)
76 call qg_field_registry%get(c_key_out,xout)
77 
78 allocate(xctl(conf%nx, conf%ny, 2))
79 xctl(:,:,:)=0.0_kind_real
80 
81 call qg_3d_covar_sqrt_inv_mult(conf%nx,conf%ny,xctl,xin,conf)
82 call zeros(xout)
83 call qg_3d_covar_sqrt_inv_mult_ad(conf%nx,conf%ny,xctl,xout,conf)
84 
85 deallocate(xctl)
86 
87 end subroutine c_qg_b_inv_mult
88 
89 ! ------------------------------------------------------------------------------
90 
91 !> Multiply streamfunction by covariance
92 
93 subroutine c_qg_b_mult(c_key_conf, c_key_in, c_key_out) bind(c,name='qg_b_mult_f90')
94 
95 use iso_c_binding
97 use qg_fields
98 use kinds
99 
100 implicit none
101 integer(c_int), intent(in) :: c_key_conf !< covar config structure
102 integer(c_int), intent(in) :: c_key_in !< Streamfunction: psi
103 integer(c_int), intent(in) :: c_key_out !< Streamfunction: psi
104 type(qg_3d_covar_config), pointer :: conf
105 type(qg_field), pointer :: xin
106 type(qg_field), pointer :: xout
107 real(kind=kind_real), allocatable :: xctl(:,:,:) ! Control vector
108 
109 call qg_3d_cov_registry%get(c_key_conf,conf)
110 call qg_field_registry%get(c_key_in,xin)
111 call qg_field_registry%get(c_key_out,xout)
112 
113 allocate(xctl(conf%nx, conf%ny, 2))
114 
115 xctl(:,:,:)=0.0_kind_real
116 call qg_3d_covar_sqrt_mult_ad(conf%nx,conf%ny,xin,xctl,conf)
117 call zeros(xout)
118 call qg_3d_covar_sqrt_mult(conf%nx,conf%ny,xout,xctl,conf)
119 
120 deallocate(xctl)
121 
122 end subroutine c_qg_b_mult
123 
124 ! ------------------------------------------------------------------------------
125 
126 !> Generate randomized increment
127 
128 subroutine c_qg_b_randomize(c_key_conf, c_key_out) bind(c,name='qg_b_randomize_f90')
130 use iso_c_binding
132 use qg_fields
134 use kinds
135 
136 implicit none
137 integer(c_int), intent(in) :: c_key_conf !< covar config structure
138 integer(c_int), intent(in) :: c_key_out !< Randomized increment
139 type(qg_3d_covar_config), pointer :: conf
140 type(qg_field), pointer :: xout
141 real(kind=kind_real), allocatable :: xctl(:,:,:) ! Control vector
142 
143 call qg_3d_cov_registry%get(c_key_conf,conf)
144 call qg_field_registry%get(c_key_out,xout)
145 
146 allocate(xctl(conf%nx, conf%ny, 2))
147 
148 call random_vector(xctl(:,:,:))
149 call zeros(xout)
150 call qg_3d_covar_sqrt_mult(conf%nx,conf%ny,xout,xctl,conf)
151 
152 deallocate(xctl)
153 
154 end subroutine c_qg_b_randomize
155 
156 ! ------------------------------------------------------------------------------
Fortran generic for generating random 1d, 2d and 3d arrays.
type(registry_t), public qg_field_registry
Linked list interface - defines registry_t type.
Definition: qg_fields.F90:65
subroutine qg_3d_covar_sqrt_inv_mult(kx, ky, xctl, xincr, config)
Multiply streamfunction by inverse(sqrt(C)), where C is 3d covariance matrix.
Structure holding configuration variables for the 3d error covariance matrices of the QG analysis...
Definition: conf.py:1
subroutine c_qg_b_inv_mult(c_key_conf, c_key_in, c_key_out)
Multiply streamfunction by inverse of covariance.
subroutine qg_3d_covar_delete(self)
Delete for the QG model&#39;s 3d error covariance matrices.
subroutine c_qg_b_setup(c_key_self, c_conf, c_key_geom)
Setup for the QG model&#39;s background error covariance matrix.
subroutine c_qg_b_mult(c_key_conf, c_key_in, c_key_out)
Multiply streamfunction by covariance.
type(registry_t), public qg_geom_registry
Linked list interface - defines registry_t type.
Definition: qg_geom_mod.F90:40
Fortran module handling geometry for the QG model.
Definition: qg_geom_mod.F90:11
subroutine qg_3d_covar_setup(c_model, geom, config)
Linked list implementation.
subroutine qg_3d_covar_sqrt_mult(kx, ky, xincr, xctl, config)
Multiply streamfunction by sqrt(C), where C is a 3d covariance matrix.
subroutine c_qg_b_delete(c_key_self)
Delete for the QG model&#39;s background error covariance matrix.
subroutine c_qg_b_randomize(c_key_conf, c_key_out)
Generate randomized increment.
type(registry_t) qg_3d_cov_registry
Linked list interface - defines registry_t type.
subroutine qg_3d_covar_sqrt_inv_mult_ad(kx, ky, xctl, xincr, config)
Multiply streamfunction by inverse(sqrt(C)) - Adjoint.
Handle fields for the QG model.
Definition: qg_fields.F90:11
Fortran module for generating random vectors.
subroutine, public zeros(self)
Definition: qg_fields.F90:151
subroutine qg_3d_covar_sqrt_mult_ad(kx, ky, xincr, xctl, config)
Multiply streamfunction by sqrt(C) - Adjoint.