FV3 Bundle
qg_wspeed_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 
9 !> Fortran module to handle wind speed observations for the QG model
11 
12 use iso_c_binding
13 use config_mod
14 use duration_mod
15 use qg_obs_data
18 use qg_vars_mod
19 use qg_locs_mod
20 use qg_goms_mod
21 use kinds
22 
23 implicit none
24 private
25 
26 ! ------------------------------------------------------------------------------
27 contains
28 ! ------------------------------------------------------------------------------
29 
30 subroutine c_qg_wspeed_setup(c_key_self, c_conf) bind(c,name='qg_wspeed_setup_f90')
31 implicit none
32 integer(c_int), intent(inout) :: c_key_self
33 type(c_ptr), intent(in) :: c_conf
34 
35 type(qg_obsoper), pointer :: self
36 character(len=1) :: svars(2) = (/"u","v"/)
37 
38 call qg_obsoper_registry%init()
39 call qg_obsoper_registry%add(c_key_self)
40 call qg_obsoper_registry%get(c_key_self, self)
41 
42 call qg_oper_setup(self, c_conf, svars, 1)
43 
44 end subroutine c_qg_wspeed_setup
45 
46 ! ------------------------------------------------------------------------------
47 
48 subroutine c_qg_wspeed_delete(c_key_self) bind(c,name='qg_wspeed_delete_f90')
49 implicit none
50 integer(c_int), intent(inout) :: c_key_self
51 
52 type(qg_obsoper), pointer :: self
53 
54 call qg_obsoper_registry%get(c_key_self, self)
55 call qg_obsoper_registry%remove(c_key_self)
56 
57 end subroutine c_qg_wspeed_delete
58 
59 ! ------------------------------------------------------------------------------
60 subroutine qg_wspeed_eqv(c_key_gom, c_key_hofx, c_bias) bind(c,name='qg_wspeed_eqv_f90')
61 implicit none
62 integer(c_int), intent(in) :: c_key_gom
63 integer(c_int), intent(in) :: c_key_hofx
64 real(c_double), intent(in) :: c_bias
65 type(qg_goms), pointer :: gom
66 type(obs_vect), pointer :: hofx
67 integer :: io, jo
68 character(len=250) :: record
69 real(kind=kind_real) :: zz
70 
71 if (abs(c_bias) > epsilon(c_bias)) call abor1_ftn ("qg_wspeed: bias not implemented")
72 
73 call qg_goms_registry%get(c_key_gom,gom)
74 call qg_obs_vect_registry%get(c_key_hofx,hofx)
75 
76 do jo=1,gom%nobs
77  io=gom%indx(jo)
78  zz=gom%values(1,jo)*gom%values(1,jo)+gom%values(2,jo)*gom%values(2,jo)
79  hofx%values(1,io)=sqrt(zz)
80 enddo
81 
82 end subroutine qg_wspeed_eqv
83 ! ------------------------------------------------------------------------------
84 subroutine qg_wspeed_equiv_tl(c_key_gom, c_key_hofx, c_key_traj, c_bias) &
85  & bind(c,name='qg_wspeed_equiv_tl_f90')
86 implicit none
87 integer(c_int), intent(in) :: c_key_gom
88 integer(c_int), intent(in) :: c_key_hofx
89 integer(c_int), intent(in) :: c_key_traj
90 real(c_double), intent(in) :: c_bias
91 type(qg_goms), pointer :: gom
92 type(obs_vect), pointer :: hofx
93 type(qg_goms), pointer :: traj
94 integer :: io, jo
95 real(kind=kind_real) :: zz, zu, zv, zt
96 
97 call qg_goms_registry%get(c_key_gom,gom)
98 call qg_obs_vect_registry%get(c_key_hofx,hofx)
99 call qg_goms_registry%get(c_key_traj,traj)
100 
101 do jo=1,gom%nobs
102  io=gom%indx(jo)
103  zu=traj%values(1,io)
104  zv=traj%values(2,io)
105  zt=sqrt(zu*zu+zv*zv)
106  if (zt>epsilon(zt)) then
107  zz=2.0_kind_real*zu*gom%values(1,jo) &
108  +2.0_kind_real*zv*gom%values(2,jo)
109  hofx%values(1,io)=zz/(2.0_kind_real*zt)
110  else
111  hofx%values(1,io)=0.0_kind_real
112  endif
113 enddo
114 
115 end subroutine qg_wspeed_equiv_tl
116 ! ------------------------------------------------------------------------------
117 subroutine qg_wspeed_equiv_ad(c_key_gom, c_key_hofx, c_key_traj, c_bias) &
118  & bind(c,name='qg_wspeed_equiv_ad_f90')
119 implicit none
120 integer(c_int), intent(in) :: c_key_gom
121 integer(c_int), intent(in) :: c_key_hofx
122 integer(c_int), intent(in) :: c_key_traj
123 real(c_double), intent(inout) :: c_bias
124 type(qg_goms), pointer :: gom
125 type(obs_vect), pointer :: hofx
126 type(qg_goms), pointer :: traj
127 integer :: io, jo
128 real(kind=kind_real) :: zz, zu, zv, zt
129 
130 call qg_goms_registry%get(c_key_gom,gom)
131 call qg_obs_vect_registry%get(c_key_hofx,hofx)
132 call qg_goms_registry%get(c_key_traj,traj)
133 
134 do jo=1,gom%nobs
135  io=gom%indx(jo)
136  zu=traj%values(1,io)
137  zv=traj%values(2,io)
138  zt=sqrt(zu*zu+zv*zv)
139  if (zt>epsilon(zt)) then
140  zz=hofx%values(1,io)/(2.0_kind_real*zt)
141  gom%values(1,jo)=2.0_kind_real*zu*zz
142  gom%values(2,jo)=2.0_kind_real*zv*zz
143  else
144  gom%values(1,jo)=0.0_kind_real
145  gom%values(2,jo)=0.0_kind_real
146  endif
147 enddo
148 
149 end subroutine qg_wspeed_equiv_ad
150 ! ------------------------------------------------------------------------------
151 subroutine qg_wspeed_gettraj(c_key_self, c_nobs, c_vars, c_key_traj) bind(c,name='qg_wspeed_gettraj_f90')
152 use fckit_log_module, only : fckit_log
153 implicit none
154 integer(c_int), intent(in) :: c_key_self
155 integer(c_int), intent(in) :: c_nobs
156 integer(c_int), dimension(*), intent(in) :: c_vars !< List of variables
157 integer(c_int), intent(inout) :: c_key_traj
158 
159 type(qg_obsoper), pointer :: self
160 type(qg_goms), pointer :: traj
161 type(qg_vars) :: vars
162 integer, allocatable :: mobs(:)
163 integer :: jj
164 
165 call qg_obsoper_registry%get(c_key_self, self)
166 
167 call qg_goms_registry%init()
168 call qg_goms_registry%add(c_key_traj)
169 call qg_goms_registry%get(c_key_traj,traj)
170 call qg_vars_create(vars, c_vars)
171 allocate(mobs(c_nobs))
172 do jj=1,c_nobs
173  mobs(jj)=jj
174 enddo
175 call gom_setup(traj, vars, mobs)
176 deallocate(mobs)
177 
178 end subroutine qg_wspeed_gettraj
179 ! ------------------------------------------------------------------------------
180 subroutine qg_wspeed_settraj(c_key_gom, c_key_traj) bind(c,name='qg_wspeed_settraj_f90')
181 use fckit_log_module, only : fckit_log
182 implicit none
183 integer(c_int), intent(in) :: c_key_gom
184 integer(c_int), intent(in) :: c_key_traj
185 
186 type(qg_goms), pointer :: gom
187 type(qg_goms), pointer :: traj
188 integer :: jo, io
189 
190 call qg_goms_registry%get(c_key_gom,gom)
191 call qg_goms_registry%get(c_key_traj,traj)
192 
193 do jo=1,gom%nobs
194  io=gom%indx(jo)
195  traj%values(1,io)=gom%values(1,jo)
196  traj%values(2,io)=gom%values(2,jo)
197 enddo
198 
199 end subroutine qg_wspeed_settraj
200 ! ------------------------------------------------------------------------------
201 
202 end module qg_wspeed_mod
type(registry_t), public qg_obs_vect_registry
Linked list interface - defines registry_t type.
subroutine qg_wspeed_equiv_tl(c_key_gom, c_key_hofx, c_key_traj, c_bias)
subroutine, public gom_setup(self, vars, kobs)
Definition: qg_goms_mod.F90:90
type(registry_t), public qg_obsoper_registry
Linked list interface - defines registry_t type.
subroutine qg_wspeed_settraj(c_key_gom, c_key_traj)
subroutine c_qg_wspeed_setup(c_key_self, c_conf)
subroutine qg_wspeed_gettraj(c_key_self, c_nobs, c_vars, c_key_traj)
subroutine qg_wspeed_eqv(c_key_gom, c_key_hofx, c_bias)
subroutine c_qg_wspeed_delete(c_key_self)
subroutine, public qg_vars_create(self, kvars)
Linked list implementation.
Definition: qg_vars_mod.F90:46
Fortran module handling observation locations.
Definition: qg_locs_mod.F90:11
Handle observations for the QG model.
Definition: qg_obs_data.F90:11
type(registry_t), public qg_goms_registry
Linked list interface - defines registry_t type.
Definition: qg_goms_mod.F90:42
Fortran module to handle variables for the QG model.
Definition: qg_vars_mod.F90:11
Fortran module to handle wind speed observations for the QG model.
subroutine, public qg_oper_setup(self, c_conf, svars, ncol)
Linked list implementation.
Fortran module handling interpolated (to obs locations) model variables.
Definition: qg_goms_mod.F90:11
Fortran module for streamfunction observations for the QG model.
subroutine qg_wspeed_equiv_ad(c_key_gom, c_key_hofx, c_key_traj, c_bias)
Fortran module handling observation vectors.