CTBTO FLEXPART WO4 (2015-10-15)
 All Classes Files Functions Variables
convertfields.F90
Go to the documentation of this file.
1 !**********************************************************************
2 ! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 *
3 ! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, *
4 ! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann *
5 ! *
6 ! This file is part of FLEXPART. *
7 ! *
8 ! FLEXPART is free software: you can redistribute it and/or modify *
9 ! it under the terms of the GNU General Public License as published by*
10 ! the Free Software Foundation, either version 3 of the License, or *
11 ! (at your option) any later version. *
12 ! *
13 ! FLEXPART is distributed in the hope that it will be useful, *
14 ! but WITHOUT ANY WARRANTY; without even the implied warranty of *
15 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
16 ! GNU General Public License for more details. *
17 ! *
18 ! You should have received a copy of the GNU General Public License *
19 ! along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. *
20 !**********************************************************************
21 
22 subroutine convertfields(ind,metdata_format,dumpPath)
23  ! i o
24  !*****************************************************************************
25  ! *
26  ! This subroutine manages the 3 data fields to be kept in memory. *
27  ! During the first time step of petterssen it has to be fulfilled that the *
28  ! first data field must have |wftime|<itime, i.e. the absolute value of *
29  ! wftime must be smaller than the absolute value of the current time in [s].*
30  ! The other 2 fields are the next in time after the first one. *
31  ! Pointers (memind) are used, because otherwise one would have to resort the*
32  ! wind fields, which costs a lot of computing time. Here only the pointers *
33  ! are resorted. *
34  ! *
35  ! Author: A. Stohl *
36  ! *
37  ! 29 April 1994 *
38  ! *
39  !*****************************************************************************
40  ! Changes, Bernd C. Krueger, Feb. 2001:
41  ! Variables tth,qvh,tthn,qvhn (on eta coordinates) in common block.
42  ! Function of nstop extended.
43  !*****************************************************************************
44 
45  ! Changes Arnold, D. and Morton, D. (2015): *
46  ! -- description of local and common variables *
47  !*****************************************************************************
48 
49 
50  use par_mod
51  use com_mod
52  use fpmetbinary_mod
53 
54  implicit none
55 
56  !****************************************************************************
57  ! Subroutine Parameters: *
58  ! input *
59  ! nstop > 0, if trajectory has to be terminated *
60  ! output *
61  ! metdata_format 0 = unknown, 1 = ecmwf, 2 = gfs meteorological data *
62  !
63  integer :: ind, metdata_format, i, lastslash
64  integer :: dumpdata
65  character(len=512):: fpfname, dumppath, filename
66 
67 
68  !****************************************************************************
69 
70  !****************************************************************************
71  ! Local variables *
72  !
73  ! indj indicates the number of the wind field to be read in *
74  ! memaux auxiliary variable to shuffle winds in memory *
75  ! indmin remembers the number of wind fields already treated *
76  ! uu(0:nxmax,0:nymax,nuvzmax,2) wind components in x-direction [m/s] *
77  ! vv(0:nxmax,0:nymax,nuvzmax,2) wind components in y-direction [m/s] *
78  ! ww(0:nxmax,0:nymax,nwzmax,2) wind components in z-direction [deltaeta/s]*
79  ! tt(0:nxmax,0:nymax,nuvzmax,2) temperature [K] *
80  ! ps(0:nxmax,0:nymax,2) surface pressure [Pa] *
81  !
82  integer :: memaux
83  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
84  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
85  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
86  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
87  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
88  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
89  real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
90  real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
91  !****************************************************************************
92 
93 
94  !****************************************************************************
95  ! Global variables *
96  ! from par_mod and com_mod *
97  ! nx,ny,nuvz,nwz field dimensions in x,y and z direction *
98  ! indmin remembers the number of wind fields already treated *
99  ! memind(2) pointer, on which place the wind fields are stored *
100  ! memtime(2) [s] times of the wind fields, which are kept in memory *
101  ! ndinterval [s] time difference between the two wind fields read in *
102  ! ldirect 1 for forward, -1 for backward simulation *
103  ! numbwf actual number of wind fields *
104  ! wftime(maxwf) [s] times relative to beginning time of wind fields *
105  ! wfname(maxwf) name of met file (used for performance timing out)*
106  !****************************************************************************
107 
108 
109 #ifdef PERFTIMER
110  INTEGER millisecs_start, millisecs_stop, count_rate, count_max
111 #endif
112 
113 !-----------------------------------------------------------------------------
114 
115 
116  ! Current time is after 2nd wind field
117  ! -> Resort wind field pointers, so that current time is between 1st and 2nd
118  !***************************************************************************
119 
120  !memaux=memind(1)
121  !memind(1)=memind(2)
122  !memind(2)=memaux
123  !memtime(1)=memtime(2)
124 
125 
126  ! Read a new wind field and store it on place memind(2)
127  !******************************************************
128 #ifdef PERFTIMER
129  CALL system_clock(millisecs_start, count_rate, count_max)
130 #endif
131  if (metdata_format.eq.ecmwf_metdata) then
132  call readwind_ecmwf(ind,memind(1),uuh,vvh,wwh)
133  call readwind_nests(ind,memind(1),uuhn,vvhn,wwhn)
134  call calcpar_ecmwf(memind(1),uuh,vvh,pvh)
135  call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format)
136  call verttransform_ecmwf(memind(1),uuh,vvh,wwh,pvh)
137  call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
138  memtime(1)=wftime(ind)
139  endif
140  if (metdata_format.eq.gfs_metdata) then
141  call readwind_gfs(ind,memind(1),uuh,vvh,wwh)
142  call readwind_nests(ind,memind(1),uuhn,vvhn,wwhn)
143  call calcpar_gfs(memind(1),uuh,vvh,pvh)
144  call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format)
145  call verttransform_gfs(memind(1),uuh,vvh,wwh,pvh)
146  call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
147  memtime(1)=wftime(ind)
148  endif
149 
150 #ifdef PERFTIMER
151  CALL system_clock(millisecs_stop, count_rate, count_max)
152  print *, 'Wall time to process: ', trim(wfname(ind)), &
153  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
154 #endif
155 
156 
157  lastslash = 0
158  do i=1,len(wfname(ind))
159  if (wfname(ind)(i:i).eq.'/') then
160  lastslash = i
161  end if
162  end do
163  filename = wfname(ind)(lastslash+1:len(wfname(ind)))
164 
165  if ( ldirect.eq.1 ) then
166  fpfname = trim(filename) // '_fwd.fp'
167  else
168  fpfname = trim(filename) // '_bwd.fp'
169  endif
170  print *, 'writing ', trim(dumppath) // '/' // trim(fpfname)
171 
172 #ifdef PERFTIMER
173  CALL system_clock(millisecs_start, count_rate, count_max)
174 #endif
175 
176  CALL fpmetbinary_dump( trim(dumppath) // '/' // trim(fpfname), memind(1))
177 
178 #ifdef PERFTIMER
179  CALL system_clock(millisecs_stop, count_rate, count_max)
180  print *, 'Wall time to dump to: ', trim(fpfname), &
181  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
182 #endif
183 
184 end subroutine convertfields
185 
186 
subroutine convertfields(ind, metdata_format, dumpPath)
subroutine fpmetbinary_dump(filename, cm_index)
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_ecmwf(indj, n, uuh, vvh, wwh)
Definition: readwind.F90:22
subroutine readwind_gfs(indj, n, uuh, vvh, wwh)
subroutine calcpar_gfs(n, uuh, vvh, pvh)
Definition: calcpar_gfs.f90:22
subroutine readwind_nests(indj, n, uuhn, vvhn, wwhn)
subroutine calcpar_ecmwf(n, uuh, vvh, pvh)
Definition: calcpar.f90:22
subroutine verttransform_gfs(n, uuh, vvh, wwh, pvh)