UFO
ufo_metoffice_rmatrixradiance_mod.f90
Go to the documentation of this file.
1 ! (C) British Crown Copyright 2017-2018 Met Office
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 
6 !> Fortran derived type to hold data for the observation covariance
7 
9 
10 use kinds
11 use netcdf
12 
13 implicit none
14 private
15 
17  integer, allocatable :: channels(:) !< array with instruemnt channel numbers
18  integer :: wmo_id !< wmo id for satellite
19  integer :: rtype !< type of r-matrix (1=full; 2=diagonal)
20  integer :: nchans !< number of channels in rmatrix
21  real(kind_real), allocatable :: errors(:) !< vector of errors
22 
23 contains
24  procedure :: setup => rmatrix_setup
25  procedure :: delete => rmatrix_delete
26  procedure :: print => rmatrix_print
27 
29 
30 ! ------------------------------------------------------------------------------
31 contains
32 ! ------------------------------------------------------------------------------
33 !> Setup for the full r_matrix
34 !!
35 !! \author Met Office
36 !!
37 !! \date 09/06/2020: Created
38 !!
39 subroutine rmatrix_setup(self, filename)
40 
41 implicit none
42 class(ufo_metoffice_rmatrixradiance), intent(inout) :: self !< Full R matrix structure
43 character(len=*), intent(in) :: filename !< Path to input filename
44 
45 integer :: error ! error code for read
46 integer :: lun ! value for identifying file
47 integer :: dimid ! value for dimension id
48 integer :: varid ! value for variable id
49 integer :: nchans ! number of channels
50 
51 ! defaults
52 self % nchans = 0
53 
54 ! Open netcdf file
55 error = nf90_open(trim(filename), nf90_nowrite, lun)
56 if (error /= 0) call abor1_ftn("error in opening file")
57 
58 ! Read wmo_id
59 error = nf90_inq_varid(lun, "wmo_id", varid)
60 error = nf90_get_var(lun, varid, self % wmo_id)
61 if (error /= 0) call abor1_ftn("error in reading the WMO ID")
62 
63 ! Read rtype
64 error = nf90_inq_varid(lun, "r_type", varid)
65 error = nf90_get_var(lun, varid, self % rtype)
66 if (error /= 0) call abor1_ftn("error in reading the rmatrix type")
67 
68 ! Get dimensions of arrays
69 error = nf90_inq_dimid(lun, "nchans", dimid)
70 error = nf90_inquire_dimension(lun, dimid, len=self % nchans)
71 if (error /= 0) call abor1_ftn("error in reading the dimensions from numchans")
72 
73 ! Allocate arrays
74 allocate(self % channels(self % nchans))
75 allocate(self % errors(self % nchans))
76 
77 ! Read channels from files
78 error = nf90_inq_varid(lun, "channels", varid)
79 error = nf90_get_var(lun, varid, self % channels)
80 if (error /= 0) call abor1_ftn("error in reading the rmatrix channels")
81 
82 ! Read channels from files
83 error = nf90_inq_varid(lun, "obs_error", varid)
84 error = nf90_get_var(lun, varid, self % errors)
85 if (error /= 0) call abor1_ftn("error in reading the rmatrix errors")
86 
87 ! Close netcdf file
88 error = nf90_close(lun)
89 if (error /= 0) call abor1_ftn("error in closing the rmatrix")
90 
91 write(*,*) "Successfully opened, read and closed r matrix netcdf file"
92 
93 end subroutine rmatrix_setup
94 
95 ! ------------------------------------------------------------------------------
96 !> Delete method for the full r_matrix
97 !!
98 !! \author Met Office
99 !!
100 !! \date 09/06/2020: Created
101 !!
102 subroutine rmatrix_delete(self)
103 
104 implicit none
105 class(ufo_metoffice_rmatrixradiance), intent(inout) :: self !< Full R matrix structure
106 
107 if (allocated(self % channels)) deallocate(self % channels)
108 if (allocated(self % errors)) deallocate(self % errors)
109 
110 end subroutine rmatrix_delete
111 
112 ! ------------------------------------------------------------------------------
113 !> Print method for the full r_matrix
114 !!
115 !! \author Met Office
116 !!
117 !! \date 09/06/2020: Created
118 !!
119 subroutine rmatrix_print(self)
120 
121 implicit none
122 class(ufo_metoffice_rmatrixradiance), intent(inout) :: self !< Full R matrix structure
123 
124 write(*,*) "wmo_id = ", self % wmo_id
125 write(*,*) "rtype = ", self % rtype
126 write(*,*) "channels = ", self % channels(:)
127 write(*,*) "errors = ", self % errors(:)
128 
129 end subroutine rmatrix_print
130 
Fortran derived type to hold data for the observation covariance.
subroutine rmatrix_print(self)
Print method for the full r_matrix.
subroutine rmatrix_setup(self, filename)
Setup for the full r_matrix.
subroutine rmatrix_delete(self)
Delete method for the full r_matrix.