FV3 Bundle
m_diag_marine_conv_mod.F90
Go to the documentation of this file.
2  use fckit_log_module, only: fckit_log
4  use nc_diag_write_mod, only: nc_diag_init, nc_diag_metadata, nc_diag_write
5 
6  use nc_diag_read_mod, only: nc_diag_read_get_dim
8  use nc_diag_read_mod, only: nc_diag_read_get_var
9  use nc_diag_read_mod, only: nc_diag_read_get_global_attr
10 
11  implicit none
12 
13  private
14  save
15 
16  public :: diag_marine_conv_header
17  public :: diag_marine_conv_tracer !generic name for tao T,S obs - obs w/ single obs associated
18 
20 
23 
25  character(3), dimension(:), allocatable :: obstype
26  integer(i_kind) :: n_obstype
27  integer(i_kind) :: n_observations_tracer
28  integer(i_kind) :: date
30 
32  integer :: station_id
33  real(r_kind) :: observation_type
34  real(r_kind) :: latitude
35  real(r_kind) :: longitude
36  real(r_kind) :: station_depth
37  real(r_kind) :: time
38  real(r_kind) :: observation
39  real(r_kind) :: obs_minus_forecast_adjusted
40  real(r_kind) :: obs_minus_forecast_unadjusted
42 
43  integer,parameter :: maxobstype=30
44  integer :: nobstype
45  integer :: tao_tracer_type = 120 ! type for TAO T and S
46  integer :: t_qcmark = 0 ! 0=tv; 1=tdry
47 
48 contains
49 
50  subroutine write_split_marine_conv_diag_nc(infn,tao_header, tao_tracer, append_suffix)
51  character(120), intent(in) :: infn
52  type(diag_marine_conv_header), intent(in) :: tao_header
53  type(diag_marine_conv_tracer),dimension(tao_header%n_Observations_Tracer),intent(in) :: tao_tracer
54  logical, intent(in) :: append_suffix
55 
56  end subroutine write_split_marine_conv_diag_nc
57 
58  subroutine read_marine_conv_diag_nc_header(infn,tao_header)
59  character(len=*), intent(in) :: infn
60  type(diag_marine_conv_header), intent(inout) :: tao_header
61 
62  end subroutine read_marine_conv_diag_nc_header
63 
64  subroutine read_marine_conv_diag_nc_tracer(infn, tao_tracer, ierr)
65  use netcdf
66  implicit none
67  character(len=*), intent(in) :: infn
68  type(diag_marine_conv_tracer),pointer, intent(inout) :: tao_tracer(:)
69  integer, intent(out) :: ierr
70  character(20) :: str, str2
71  integer(i_kind) :: i, itype, fid, nobs
72  real, allocatable :: obs(:,:,:,:), lon(:), lat(:), depth(:), time(:)
73  integer :: nlon, nlat, ndepth, ntime
74  integer :: ilon, ilat, idepth, itime, cnt, varid
75 
76  call nc_diag_read_init(infn,fid)
77 
78  nlon = nc_diag_read_get_dim(fid,'lon')
79  nlat = nc_diag_read_get_dim(fid,'lat')
80  ndepth = nc_diag_read_get_dim(fid,'depth')
81  ntime = nc_diag_read_get_dim(fid,'time')
82 
83  nobs=nlon*nlat*ndepth*ntime
84 
85  allocate(tao_tracer(nobs))
86 
87  !allocate(obs(ntime, ndepth, nlat, nlon))
88  allocate(obs(nlon,nlat,ndepth,ntime))
89  allocate(time(ntime),lon(nlon),lat(nlat),depth(ndepth))
90 
91  call nc_diag_read_get_var(fid,"lat", lat )
92  call nc_diag_read_get_var(fid,"lon", lon )
93  call nc_diag_read_get_var(fid,"depth", depth )
94  call nc_diag_read_get_var(fid,"time", time )
95  !call nc_diag_read_get_var(fid,"T_20", obs ) ! reading of 4d array not im[lemented yet
96  call nc_diag_read_close(infn)
97 
98  call check( nf90_open(infn, nf90_nowrite, fid) )
99  call check( nf90_inq_varid(fid, "T_20", varid) )
100  call check( nf90_get_var(fid, varid, obs) )
101  call check( nf90_close(fid) )
102  cnt=1
103 
104  do ilon=1,nlon
105  do ilat=1,nlat
106  do idepth=1,ndepth
107  tao_tracer(cnt)%Latitude = lat(ilat)
108  tao_tracer(cnt)%Longitude = lon(ilon)
109  tao_tracer(cnt)%Station_Depth = depth(idepth)
110  tao_tracer(cnt)%Time = time(1)
111  !tao_tracer(cnt)%Observation = obs(1,idepth,ilat,ilon)!obs(ilon,ilat,idepth,1)
112  tao_tracer(cnt)%Observation = obs(ilon,ilat,idepth,1)
113  cnt=cnt+1
114  end do
115  end do
116  end do
117 
118  deallocate(lat,lon,depth,time,obs)
119 
120  end subroutine read_marine_conv_diag_nc_tracer
121 
122  subroutine check(status)
123  use netcdf
124  implicit none
125  integer, intent ( in) :: status
126 
127  if(status /= nf90_noerr) then
128  print *, trim(nf90_strerror(status))
129  stop "Stopped"
130  end if
131  end subroutine check
132 end module m_diag_marine_conv
subroutine, public read_marine_conv_diag_nc_tracer(infn, tao_tracer, ierr)
integer, parameter, public i_kind
Definition: ncd_kinds.F90:71
subroutine nc_diag_init(filename, append)
integer, parameter maxobstype
subroutine, public write_split_marine_conv_diag_nc(infn, tao_header, tao_tracer, append_suffix)
subroutine nc_diag_read_close(filename, file_ncdr_id, from_pop)
subroutine, public read_marine_conv_diag_nc_header(infn, tao_header)
integer, parameter, public r_double
Definition: ncd_kinds.F90:80
integer, parameter, public r_single
Definition: ncd_kinds.F90:79
integer, parameter, public r_kind
Definition: ncd_kinds.F90:108
subroutine nc_diag_read_init(filename, file_ncdr_id, from_push)