FV3 Bundle
ncutils.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2018 UCAR
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 module to handle temperature profile observations
7 
9 
10  use ncd_kinds
11 
12  implicit none
13  private
14 
15  ! General type definition for marine obs
16  type, public :: simple_marine_obs
17  integer :: station_id
18  real(r_kind) :: observation_type
19  real(r_kind) :: latitude
20  real(r_kind) :: longitude
21  real(r_kind) :: depth
22  real(r_kind) :: time
23  real(r_kind) :: observation
24  real(r_kind) :: observation_error
25  real(r_kind) :: obs_minus_forecast
26  end type simple_marine_obs
27 
28  type, public :: diag_marine_obs
29  integer :: nobs !< Number of obs
30  character(len=120) :: filename !< Netcdf output filename
31  logical :: append=.false. !< If file exist, append to it
32  type(simple_marine_obs), allocatable :: diag(:) !< Data type to hold obs space diagnostics
33  contains
34  procedure :: init
35  procedure :: write_diag
36  procedure :: write_geoval
37  procedure :: finalize
38  end type diag_marine_obs
39 
40 contains
41 
42  ! ------------------------------------------------------------------------------
43 
44  subroutine init(self, nobs, filename)
45 
46  implicit none
47 
48  class(diag_marine_obs) ,intent(out) :: self !< Obs space diagnostics
49  integer ,intent(in) :: nobs !< Number of obs
50  character(len=120) ,intent(in) :: filename !< Filename for netcdf output
51 
52  self%nobs=nobs
53  self%filename=filename
54  allocate(self%diag(nobs))
55 
56  end subroutine init
57 
58  ! ------------------------------------------------------------------------------
59 
60  subroutine finalize(self)
61 
62  implicit none
63 
64  class(diag_marine_obs), intent(inout) :: self !< Obs space diagnostics
65 
66  self%nobs=0
67  self%filename=''
68  deallocate(self%diag)
69 
70  end subroutine finalize
71 
72  ! ------------------------------------------------------------------------------
73 
74  subroutine write_diag(self)
75 
76  use netcdf
77  use ufo_vars_mod
78 
79  implicit none
80 
81  class(diag_marine_obs), intent(inout) :: self !< Obs space diagnostics
82 
83  !netcdf stuff
84  integer(kind=4) :: iNcid,i,iStationNo
85  integer(kind=4) :: iDimStation_ID,iDimTime_ID,iDimLenStringName_ID, iDimLev_ID
86  integer(kind=4) :: iVarLON_ID, iVarLAT_ID, iVarLev_ID, iVarObs_ID
87  integer(kind=4) :: iVarOMF_ID, iVarOMA_ID, iVarSIGO_ID, iVarOBSID_ID
88 
89  integer, allocatable :: obsid(:)
90  integer :: nlev,nobs,iobs
91 
92  ! Create file.
93  call check(nf90_create(self%filename, nf90_clobber, incid))
94  call check(nf90_def_dim(incid, "nobs", nf90_unlimited, idimstation_id))
95  nobs = self%nobs
96 
97  ! Define of variables.
98  call check( nf90_def_var(incid, "OBSID", nf90_int, (/ idimstation_id /), ivarobsid_id) )
99  call check( nf90_def_var(incid, "lon", nf90_real, (/ idimstation_id /), ivarlon_id) )
100  call check( nf90_def_var(incid, "lat", nf90_real, (/ idimstation_id /), ivarlat_id) )
101  call check( nf90_def_var(incid, "lev", nf90_real, (/ idimstation_id /), ivarlev_id) )
102  call check( nf90_def_var(incid, "obs", nf90_real, (/ idimstation_id /), ivarobs_id) )
103  call check( nf90_def_var(incid, "sigo", nf90_real, (/ idimstation_id /), ivarsigo_id) )
104  call check( nf90_def_var(incid, "omf", nf90_real, (/ idimstation_id /), ivaromf_id) )
105 
106  ! End define mode.
107  call check(nf90_enddef(incid))
108 
109  ! Writing
110  call check(nf90_put_var(incid, ivarlon_id , self%diag(:)%Station_ID))
111  call check(nf90_put_var(incid, ivarlon_id , self%diag(:)%Longitude))
112  call check(nf90_put_var(incid, ivarlat_id , self%diag(:)%Latitude))
113  call check(nf90_put_var(incid, ivarlev_id , self%diag(:)%Depth))
114  call check(nf90_put_var(incid, ivarsigo_id , self%diag(:)%Observation_error))
115  call check(nf90_put_var(incid, ivarobs_id , self%diag(:)%Observation))
116  call check(nf90_put_var(incid, ivaromf_id , self%diag(:)%Obs_Minus_Forecast))
117 
118  ! Close file.
119  call check(nf90_close(incid))
120  self%append=.true.
121  end subroutine write_diag
122 
123  ! ------------------------------------------------------------------------------
124 
125  subroutine write_geoval(self,varname,geoval,arg_dim_name)
127  use netcdf
128  use ufo_vars_mod
129  use ufo_geovals_mod
130 
131  implicit none
132 
133  class(diag_marine_obs) , intent(in) :: self !< Obs space diagnostics
134  character(len=MAXVARLEN) , intent(in) :: varname !< One of var_ from ufo_vars_mod
135  type(ufo_geoval) , pointer, intent(in) :: geoval !< 2D array for 1 geoval
136  character(len=MAXVARLEN), optional , intent(in) :: arg_dim_name !< Name for the second dimension
137  !netcdf stuff
138  integer(kind=4) :: iNcid
139  integer(kind=4) :: iDimStation_ID, iDimLev_ID
140  integer(kind=4) :: iVargeoval_ID
141  integer :: nlev, ndims
142  character(len=MAXVARLEN) :: dim_name
143 
144  dim_name="nlev"
145  if (present(arg_dim_name)) dim_name=arg_dim_name
146 
147  if (self%append) then ! If file exists, append to it
148  call check( nf90_open(self%filename, nf90_write, incid) )
149  call check( nf90_inquire(incid, ndimensions = ndims) )
150  call check( nf90_inq_dimid(incid, "nobs", idimstation_id) )
151  if (ndims.eq.2) then
152  call check( nf90_inq_dimid(incid, dim_name, idimlev_id) )
153  end if
154  call check( nf90_redef(incid) )
155  if (ndims.eq.1) then
156  nlev = geoval%nval
157  call check(nf90_def_dim(incid, dim_name, nlev, idimlev_id))
158  end if
159  else
160  call check(nf90_create(self%filename, nf90_clobber, incid))
161  call check(nf90_def_dim(incid, "nobs", nf90_unlimited, idimstation_id))
162  nlev = geoval%nval
163  call check(nf90_def_dim(incid, dim_name, nlev, idimlev_id))
164  end if
165 
166  ! Define of variables.
167  call check( nf90_def_var(incid, varname, nf90_real, (/ idimlev_id, idimstation_id /), ivargeoval_id) )
168 
169  ! End define mode.
170  call check(nf90_enddef(incid))
171 
172  ! Writing
173  call check(nf90_put_var(incid, ivargeoval_id , geoval%vals))
174 
175  ! Close file.
176  call check(nf90_close(incid))
177 
178  end subroutine write_geoval
179 
180  ! ------------------------------------------------------------------------------
181 
182  subroutine check(status)
184  use netcdf
185  IMPLICIT NONE
186  integer(4), intent ( in) :: status
187  if(status /= nf90_noerr) then
188  print *, trim(nf90_strerror(status))
189  stop "Stopped"
190  end if
191  end subroutine check
192 
193 end module ufo_marine_ncutils
subroutine finalize(self)
Definition: ncutils.F90:61
Fortran module to handle temperature profile observations.
Definition: ncutils.F90:8
subroutine check(status)
Definition: ncutils.F90:183
subroutine write_geoval(self, varname, geoval, arg_dim_name)
Definition: ncutils.F90:126
subroutine write_diag(self)
Definition: ncutils.F90:75
subroutine init(self, nobs, filename)
Definition: ncutils.F90:45