FV3 Bundle
qg_wind_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 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_wind_setup(c_key_self, c_conf) bind(c,name='qg_wind_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, 2)
43 
44 end subroutine c_qg_wind_setup
45 
46 ! ------------------------------------------------------------------------------
47 
48 subroutine c_qg_wind_delete(c_key_self) bind(c,name='qg_wind_delete_f90')
49 implicit none
50 integer(c_int), intent(inout) :: c_key_self
51 
52 type(qg_obsoper), pointer :: self
53 call qg_obsoper_registry%get(c_key_self, self)
54 call qg_obsoper_registry%remove(c_key_self)
55 
56 end subroutine c_qg_wind_delete
57 
58 ! ------------------------------------------------------------------------------
59 subroutine qg_wind_equiv(c_key_gom, c_key_hofx, c_bias) bind(c,name='qg_wind_equiv_f90')
60 implicit none
61 integer(c_int), intent(in) :: c_key_gom
62 integer(c_int), intent(in) :: c_key_hofx
63 real(c_double), intent(in) :: c_bias(2)
64 type(qg_goms), pointer :: gom
65 type(obs_vect), pointer :: hofx
66 integer :: io, jo
67 
68 call qg_goms_registry%get(c_key_gom,gom)
69 call qg_obs_vect_registry%get(c_key_hofx,hofx)
70 
71 do jo=1,gom%nobs
72  io=gom%indx(jo)
73  hofx%values(1,io)=gom%values(1,jo) + c_bias(1)
74  hofx%values(2,io)=gom%values(2,jo) + c_bias(2)
75 enddo
76 
77 end subroutine qg_wind_equiv
78 ! ------------------------------------------------------------------------------
79 subroutine qg_wind_equiv_tl(c_key_gom, c_key_hofx, c_bias) bind(c,name='qg_wind_equiv_tl_f90')
80 implicit none
81 integer(c_int), intent(in) :: c_key_gom
82 integer(c_int), intent(in) :: c_key_hofx
83 real(c_double), intent(in) :: c_bias(2)
84 type(qg_goms), pointer :: gom
85 type(obs_vect), pointer :: hofx
86 integer :: io, jo
87 
88 call qg_goms_registry%get(c_key_gom,gom)
89 call qg_obs_vect_registry%get(c_key_hofx,hofx)
90 
91 do jo=1,gom%nobs
92  io=gom%indx(jo)
93  hofx%values(1,io)=gom%values(1,jo) + c_bias(1)
94  hofx%values(2,io)=gom%values(2,jo) + c_bias(2)
95 enddo
96 
97 end subroutine qg_wind_equiv_tl
98 ! ------------------------------------------------------------------------------
99 subroutine qg_wind_equiv_ad(c_key_gom, c_key_hofx, c_bias) bind(c,name='qg_wind_equiv_ad_f90')
100 implicit none
101 integer(c_int), intent(in) :: c_key_gom
102 integer(c_int), intent(in) :: c_key_hofx
103 real(c_double), intent(inout) :: c_bias(2)
104 type(qg_goms), pointer :: gom
105 type(obs_vect), pointer :: hofx
106 integer :: io, jo
107 
108 call qg_goms_registry%get(c_key_gom,gom)
109 call qg_obs_vect_registry%get(c_key_hofx,hofx)
110 
111 do jo=1,gom%nobs
112  io=gom%indx(jo)
113  gom%values(1,jo)=hofx%values(1,io)
114  gom%values(2,jo)=hofx%values(2,io)
115  c_bias(1) = c_bias(1) + hofx%values(1,io)
116  c_bias(2) = c_bias(2) + hofx%values(2,io)
117 enddo
118 
119 end subroutine qg_wind_equiv_ad
120 ! ------------------------------------------------------------------------------
121 
122 end module qg_wind_mod
type(registry_t), public qg_obs_vect_registry
Linked list interface - defines registry_t type.
subroutine c_qg_wind_delete(c_key_self)
Definition: qg_wind_mod.F90:49
type(registry_t), public qg_obsoper_registry
Linked list interface - defines registry_t type.
subroutine qg_wind_equiv_ad(c_key_gom, c_key_hofx, c_bias)
Fortran module to handle wind observations for the QG model.
Definition: qg_wind_mod.F90:10
subroutine qg_wind_equiv(c_key_gom, c_key_hofx, c_bias)
Definition: qg_wind_mod.F90:60
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
subroutine c_qg_wind_setup(c_key_self, c_conf)
Definition: qg_wind_mod.F90:31
Fortran module to handle variables for the QG model.
Definition: qg_vars_mod.F90:11
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.
Fortran module handling observation vectors.
subroutine qg_wind_equiv_tl(c_key_gom, c_key_hofx, c_bias)
Definition: qg_wind_mod.F90:80