OOPS
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 kinds
13 use iso_c_binding
14 use qg_gom_mod
15 use qg_obsvec_mod
17 
18 implicit none
19 
20 private
22 ! ------------------------------------------------------------------------------
23 contains
24 ! ------------------------------------------------------------------------------
25 ! Public
26 ! ------------------------------------------------------------------------------
27 !> Get equivalent for wind speed
28 subroutine qg_wspeed_equiv(gom,hofx,bias)
29 
30 implicit none
31 
32 ! Passed variables
33 type(qg_gom),intent(in) :: gom !< GOM
34 type(qg_obsvec),intent(inout) :: hofx !< Observation vector
35 real(kind_real),intent(in) :: bias !< Bias
36 
37 ! Local variables
38 integer :: iobs
39 
40 ! Check bias
41 if (abs(bias)>epsilon(bias)) call abor1_ftn ('qg_wspeed_equiv: bias not implemented')
42 
43 ! Loop over observations
44 do iobs=1,gom%nobs
45  hofx%values(1,iobs) = sqrt(gom%u(iobs)*gom%u(iobs)+gom%v(iobs)*gom%v(iobs))
46 enddo
47 
48 end subroutine qg_wspeed_equiv
49 ! ------------------------------------------------------------------------------
50 !> Get equivalent for wind speed - tangent linear
51 subroutine qg_wspeed_equiv_tl(gom,hofx,traj,bias)
52 
53 implicit none
54 
55 ! Passed variables
56 type(qg_gom),intent(in) :: gom !< GOM
57 type(qg_obsvec),intent(inout) :: hofx !< Observation vector
58 type(qg_gom),intent(in) :: traj !< GOM trajectory
59 real(kind_real),intent(in) :: bias !< Bias
60 
61 ! Local variables
62 integer :: iobs
63 real(kind_real) :: zu,zv,zt
64 
65 ! Loop over observations
66 do iobs=1,gom%nobs
67  zu = traj%u(iobs)
68  zv = traj%v(iobs)
69  zt = sqrt(zu**2+zv**2)
70  if (zt>epsilon(zt)) then
71  hofx%values(1,iobs) = (zu*gom%u(iobs)+zv*gom%v(iobs))/zt
72  else
73  hofx%values(1,iobs) = 0.0
74  endif
75 enddo
76 
77 end subroutine qg_wspeed_equiv_tl
78 ! ------------------------------------------------------------------------------
79 !> Get equivalent for wind speed - adjoint
80 subroutine qg_wspeed_equiv_ad(gom,hofx,traj,bias)
81 
82 implicit none
83 
84 ! Passed variables
85 type(qg_gom),intent(inout) :: gom !< GOM
86 type(qg_obsvec),intent(in) :: hofx !< Observation vector
87 type(qg_gom),intent(in) :: traj !< GOM trajectory
88 real(kind_real),intent(inout) :: bias !< Bias
89 
90 ! Local variables
91 integer :: iobs
92 real(kind_real) :: zu,zv,zt
93 
94 ! Loop over observations
95 do iobs=1,gom%nobs
96  zu = traj%u(iobs)
97  zv = traj%v(iobs)
98  zt = sqrt(zu**2+zv**2)
99  if (zt>epsilon(zt)) then
100  gom%u(iobs) = zu*hofx%values(1,iobs)/zt
101  gom%v(iobs) = zv*hofx%values(1,iobs)/zt
102  else
103  gom%u(iobs) = 0.0
104  gom%v(iobs) = 0.0
105  endif
106 enddo
107 
108 end subroutine qg_wspeed_equiv_ad
109 ! ------------------------------------------------------------------------------
110 !> Set wind speed trajectory
111 subroutine qg_wspeed_settraj(gom,traj)
112 
113 implicit none
114 
115 ! Passed variables
116 type(qg_gom),intent(in) :: gom !< GOM
117 type(qg_gom),intent(inout) :: traj !< GOM trajectory
118 
119 ! Local variables
120 integer :: iobs
121 
122 ! Loop over observations
123 do iobs=1,gom%nobs
124  traj%u(iobs) = gom%u(iobs)
125  traj%v(iobs) = gom%v(iobs)
126 end do
127 
128 end subroutine qg_wspeed_settraj
129 ! ------------------------------------------------------------------------------
130 end module qg_wspeed_mod
Fortran interface to Variables.
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.