FLEXPART CTBTO WO8
 All Classes Files Functions Variables
convertfields.F90
Go to the documentation of this file.
1 subroutine convertfields(ind,metdata_format,dumpPath)
2  ! i o
3  !*****************************************************************************
4  ! *
5 
6  ! This subrotuine produces all the meteorological data in binary (fp) *
7  ! format. It includes all the subroutines that the traditional getfields *
8  ! would use to generate the meteorological data *
9  ! *
10  ! This routine is essential part of the grib2flexpart utility *
11  !
12  ! Author: M. Harustak (AWST) *
13  ! D. Morton (Boreal Scienctific Computing) *
14  ! D. Arnold (ZAMG) *
15  ! 15 October 2015 *
16  ! *
17  !*****************************************************************************
18  !*****************************************************************************
19 
20  ! Changes Arnold, D. and Morton, D. (2015): *
21  ! -- description of local and common variables *
22  !*****************************************************************************
23 
24 
25  use par_mod
26  use com_mod
27  use fpmetbinary_mod
28 
29  implicit none
30 
31  !****************************************************************************
32  ! Subroutine Parameters: *
33  ! input *
34  ! nstop > 0, if trajectory has to be terminated *
35  ! output *
36  ! metdata_format 0 = unknown, 1 = ecmwf, 2 = gfs meteorological data *
37  !
38  integer :: ind, metdata_format, i, lastslash
39  integer :: dumpdata
40  character(len=512):: fpfname, dumppath, filename
41 
42 
43  !****************************************************************************
44 
45  !****************************************************************************
46  ! Local variables *
47  !
48  ! indj indicates the number of the wind field to be read in *
49  ! memaux auxiliary variable to shuffle winds in memory *
50  ! indmin remembers the number of wind fields already treated *
51  ! uu(0:nxmax,0:nymax,nuvzmax,2) wind components in x-direction [m/s] *
52  ! vv(0:nxmax,0:nymax,nuvzmax,2) wind components in y-direction [m/s] *
53  ! ww(0:nxmax,0:nymax,nwzmax,2) wind components in z-direction [deltaeta/s]*
54  ! tt(0:nxmax,0:nymax,nuvzmax,2) temperature [K] *
55  ! ps(0:nxmax,0:nymax,2) surface pressure [Pa] *
56  !
57  integer :: memaux
58  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
59  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
60  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
61  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
62  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
63  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
64  real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
65  real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
66  !****************************************************************************
67 
68 
69  !****************************************************************************
70  ! Global variables *
71  ! from par_mod and com_mod *
72  ! nx,ny,nuvz,nwz field dimensions in x,y and z direction *
73  ! indmin remembers the number of wind fields already treated *
74  ! memind(2) pointer, on which place the wind fields are stored *
75  ! memtime(2) [s] times of the wind fields, which are kept in memory *
76  ! ndinterval [s] time difference between the two wind fields read in *
77  ! ldirect 1 for forward, -1 for backward simulation *
78  ! numbwf actual number of wind fields *
79  ! wftime(maxwf) [s] times relative to beginning time of wind fields *
80  ! wfname(maxwf) name of met file (used for performance timing out)*
81  !****************************************************************************
82 
83 
84 #ifdef PERFTIMER
85  INTEGER millisecs_start, millisecs_stop, count_rate, count_max
86 #endif
87 
88 !-----------------------------------------------------------------------------
89 
90 
91  ! Current time is after 2nd wind field
92  ! -> Resort wind field pointers, so that current time is between 1st and 2nd
93  !***************************************************************************
94 
95  !memaux=memind(1)
96  !memind(1)=memind(2)
97  !memind(2)=memaux
98  !memtime(1)=memtime(2)
99 
100 
101  ! Read a new wind field and store it on place memind(2)
102  !******************************************************
103 #ifdef PERFTIMER
104  CALL system_clock(millisecs_start, count_rate, count_max)
105 #endif
106  if (metdata_format.eq.ecmwf_metdata) then
107  call readwind_ecmwf(ind,memind(1),uuh,vvh,wwh)
108  call readwind_nests(ind,memind(1),uuhn,vvhn,wwhn)
109  call calcpar_ecmwf(memind(1),uuh,vvh,pvh)
110  call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format)
111  call verttransform_ecmwf(memind(1),uuh,vvh,wwh,pvh)
112  call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
113  memtime(1)=wftime(ind)
114  endif
115  if (metdata_format.eq.gfs_metdata) then
116  call readwind_gfs(ind,memind(1),uuh,vvh,wwh)
117  call readwind_nests(ind,memind(1),uuhn,vvhn,wwhn)
118  call calcpar_gfs(memind(1),uuh,vvh,pvh)
119  call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format)
120  call verttransform_gfs(memind(1),uuh,vvh,wwh,pvh)
121  call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
122  memtime(1)=wftime(ind)
123  endif
124 
125 #ifdef PERFTIMER
126  CALL system_clock(millisecs_stop, count_rate, count_max)
127  print *, 'Wall time to process: ', trim(wfname(ind)), &
128  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
129 #endif
130 
131 
132  lastslash = 0
133  do i=1,len(wfname(ind))
134  if (wfname(ind)(i:i).eq.'/') then
135  lastslash = i
136  end if
137  end do
138  filename = wfname(ind)(lastslash+1:len(wfname(ind)))
139 
140  if ( ldirect.eq.1 ) then
141  fpfname = trim(filename) // '_fwd.fp'
142  else
143  fpfname = trim(filename) // '_bwd.fp'
144  endif
145  print *, 'writing ', trim(dumppath) // '/' // trim(fpfname)
146 
147 #ifdef PERFTIMER
148  CALL system_clock(millisecs_start, count_rate, count_max)
149 #endif
150 
151  CALL fpmetbinary_dump( trim(dumppath) // '/' // trim(fpfname), memind(1))
152 
153 #ifdef PERFTIMER
154  CALL system_clock(millisecs_stop, count_rate, count_max)
155  print *, 'Wall time to dump to: ', trim(fpfname), &
156  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
157 #endif
158 
159 end subroutine convertfields
160 
161 
subroutine convertfields(ind, metdata_format, dumpPath)
subroutine fpmetbinary_dump(filename, cm_index)
subroutine readwind_nests(indj, n, uuhn, vvhn, wwhn)
subroutine verttransform_ecmwf(n, uuh, vvh, wwh, pvh)
subroutine verttransform_nests(n, uuhn, vvhn, wwhn, pvhn)
subroutine calcpar_nests(n, uuhn, vvhn, pvhn, metdata_format)
subroutine readwind_gfs(indj, n, uuh, vvh, wwh)
subroutine readwind_ecmwf(indj, n, uuh, vvh, wwh)
Definition: readwind.F90:22
subroutine calcpar_gfs(n, uuh, vvh, pvh)
Definition: calcpar_gfs.f90:22
subroutine calcpar_ecmwf(n, uuh, vvh, pvh)
Definition: calcpar.f90:22
subroutine verttransform_gfs(n, uuh, vvh, wwh, pvh)