OOPS
qg_wspeed_interface.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 fckit_configuration_module, only: fckit_configuration
12 use iso_c_binding
13 use qg_gom_mod
14 use qg_obsvec_mod
15 use qg_wspeed_mod
17 
18 implicit none
19 
20 private
21 ! ------------------------------------------------------------------------------
22 contains
23 ! ------------------------------------------------------------------------------
24 !> Get equivalent for wind speed
25 subroutine qg_wspeed_equiv_c(c_key_gom,c_key_hofx,c_bias) bind(c,name='qg_wspeed_equiv_f90')
26 
27 implicit none
28 
29 ! Passed variables
30 integer(c_int),intent(in) :: c_key_gom !< GOM
31 integer(c_int),intent(in) :: c_key_hofx !< Observation vector
32 real(c_double),intent(in) :: c_bias !< Bias
33 
34 ! Local variables
35 type(qg_gom),pointer :: gom
36 type(qg_obsvec),pointer :: hofx
37 
38 ! Interface
39 call qg_gom_registry%get(c_key_gom,gom)
40 call qg_obsvec_registry%get(c_key_hofx,hofx)
41 
42 ! Call Fortran
43 call qg_wspeed_equiv(gom,hofx,c_bias)
44 
45 end subroutine qg_wspeed_equiv_c
46 ! ------------------------------------------------------------------------------
47 !> Get equivalent for wind speed - tangent linear
48 subroutine qg_wspeed_equiv_tl_c(c_key_gom,c_key_hofx,c_key_traj,c_bias) bind(c,name='qg_wspeed_equiv_tl_f90')
49 
50 implicit none
51 
52 ! Passed variables
53 integer(c_int),intent(in) :: c_key_gom !< GOM
54 integer(c_int),intent(in) :: c_key_hofx !< Observation vector
55 integer(c_int),intent(in) :: c_key_traj !< GOM trajectory
56 real(c_double),intent(in) :: c_bias !< Bias
57 
58 ! Local variables
59 type(qg_gom),pointer :: gom,traj
60 type(qg_obsvec),pointer :: hofx
61 
62 ! Interface
63 call qg_gom_registry%get(c_key_gom,gom)
64 call qg_obsvec_registry%get(c_key_hofx,hofx)
65 call qg_gom_registry%get(c_key_traj,traj)
66 
67 ! Call Fortran
68 call qg_wspeed_equiv_tl(gom,hofx,traj,c_bias)
69 
70 end subroutine qg_wspeed_equiv_tl_c
71 ! ------------------------------------------------------------------------------
72 !> Get equivalent for wind speed - adjoint
73 subroutine qg_wspeed_equiv_ad_c(c_key_gom,c_key_hofx,c_key_traj,c_bias) bind(c,name='qg_wspeed_equiv_ad_f90')
74 
75 implicit none
76 
77 ! Passed variables
78 integer(c_int),intent(in) :: c_key_gom !< GOM
79 integer(c_int),intent(in) :: c_key_hofx !< Observation vector
80 integer(c_int),intent(in) :: c_key_traj !< GOM trajectory
81 real(c_double),intent(inout) :: c_bias !< Bias
82 
83 ! Local variables
84 type(qg_gom),pointer :: gom,traj
85 type(qg_obsvec),pointer :: hofx
86 
87 ! Interface
88 call qg_gom_registry%get(c_key_gom,gom)
89 call qg_obsvec_registry%get(c_key_hofx,hofx)
90 call qg_gom_registry%get(c_key_traj,traj)
91 
92 ! Call Fortran
93 call qg_wspeed_equiv_ad(gom,hofx,traj,c_bias)
94 
95 end subroutine qg_wspeed_equiv_ad_c
96 ! ------------------------------------------------------------------------------
97 !> Get wind speed trajectory
98 subroutine qg_wspeed_gettraj_c(c_nobs,c_vars,c_key_traj) bind(c,name='qg_wspeed_gettraj_f90')
99 
100 implicit none
101 
102 ! Passed variables
103 integer(c_int),intent(in) :: c_nobs !< Number of observations
104 type(c_ptr),value,intent(in) :: c_vars !< Variables
105 integer(c_int),intent(inout) :: c_key_traj !< GOM trajectory
106 
107 ! Local variables
108 type(qg_gom),pointer :: traj
109 
110 ! Interface
111 call qg_gom_registry%init()
112 call qg_gom_registry%add(c_key_traj)
113 call qg_gom_registry%get(c_key_traj,traj)
114 traj%vars = oops_variables(c_vars)
115 
116 ! Call Fortran
117 call qg_gom_setup(traj,c_nobs)
118 
119 end subroutine qg_wspeed_gettraj_c
120 ! ------------------------------------------------------------------------------
121 !> Set wind speed trajectory
122 subroutine qg_wspeed_settraj_c(c_key_gom,c_key_traj) bind(c,name='qg_wspeed_settraj_f90')
123 
124 implicit none
125 
126 ! Passed variables
127 integer(c_int),intent(in) :: c_key_gom !< GOM
128 integer(c_int),intent(in) :: c_key_traj !< GOM trajectory
129 
130 ! Local variables
131 type(qg_gom),pointer :: gom
132 type(qg_gom),pointer :: traj
133 
134 ! Interface
135 call qg_gom_registry%get(c_key_gom,gom)
136 call qg_gom_registry%get(c_key_traj,traj)
137 
138 ! Call Fortran
139 call qg_wspeed_settraj(gom,traj)
140 
141 end subroutine qg_wspeed_settraj_c
142 ! ------------------------------------------------------------------------------
143 end module qg_wspeed_interface
Fortran interface to Variables.
type(registry_t), public qg_gom_registry
Linked list interface - defines registry_t type.
Definition: qg_gom_mod.F90:49
subroutine, public qg_gom_setup(self, nobs)
Linked list implementation.
Definition: qg_gom_mod.F90:60
type(registry_t), public qg_obsvec_registry
Linked list interface - defines registry_t type.
subroutine qg_wspeed_equiv_c(c_key_gom, c_key_hofx, c_bias)
Get equivalent for wind speed.
subroutine qg_wspeed_equiv_tl_c(c_key_gom, c_key_hofx, c_key_traj, c_bias)
Get equivalent for wind speed - tangent linear.
subroutine qg_wspeed_gettraj_c(c_nobs, c_vars, c_key_traj)
Get wind speed trajectory.
subroutine qg_wspeed_settraj_c(c_key_gom, c_key_traj)
Set wind speed trajectory.
subroutine qg_wspeed_equiv_ad_c(c_key_gom, c_key_hofx, c_key_traj, c_bias)
Get equivalent for wind speed - adjoint.
Fortran module to handle wind speed observations for the QG model.
subroutine, public qg_wspeed_equiv_ad(gom, hofx, traj, bias)
Get equivalent for wind speed - adjoint.
subroutine, public qg_wspeed_equiv(gom, hofx, bias)
Get equivalent for wind speed.
subroutine, public qg_wspeed_settraj(gom, traj)
Set wind speed trajectory.
subroutine, public qg_wspeed_equiv_tl(gom, hofx, traj, bias)
Get equivalent for wind speed - tangent linear.