source: flexpart.git/flexpart_code/getfields.F90 @ b398fb6

FPv9.3.1FPv9.3.1b_testingFPv9.3.2fp9.3.1-20161214-nc4grib2nc4_repair
Last change on this file since b398fb6 was 496c607, checked in by Don Morton <Don.Morton@…>, 8 years ago

Initial commit of FPv9.3.1

Currently, this is a clone of snapshot FPv9.3.0

  • Property mode set to 100644
File size: 12.7 KB
Line 
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
22subroutine getfields(itime,nstop,metdata_format)
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
53  implicit none
54
55  !****************************************************************************
56  ! Subroutine Parameters:                                                    *
57  !    input                                                                  *
58  ! itime [s]            current time since start date of trajectory calcu-   *
59  !                      lation                                               *
60  ! nstop                > 0, if trajectory has to be terminated              *
61  !    output                                                                 *
62  ! metdata_format       0 = unknown, 1 = ecmwf, 2 = gfs meteorological data  *
63  ! 
64  integer :: itime, nstop, metdata_format
65
66  !****************************************************************************
67
68  !****************************************************************************
69  ! Local variables                                                           *
70  !
71  ! indj                 indicates the number of the wind field to be read in *
72  ! memaux               auxiliary variable to shuffle winds in memory        *
73  ! indmin               remembers the number of wind fields already treated  *
74  ! uu(0:nxmax,0:nymax,nuvzmax,2)  wind components in x-direction [m/s]       *
75  ! vv(0:nxmax,0:nymax,nuvzmax,2)  wind components in y-direction [m/s]       *
76  ! ww(0:nxmax,0:nymax,nwzmax,2)   wind components in z-direction [deltaeta/s]*
77  ! tt(0:nxmax,0:nymax,nuvzmax,2)  temperature [K]                            *
78  ! ps(0:nxmax,0:nymax,2)          surface pressure [Pa]                      *
79  !
80  integer :: indj,memaux
81  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
82  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
83  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
84  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
85  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
86  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
87  real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
88  real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
89  integer :: indmin = 1
90  !****************************************************************************
91
92
93  !****************************************************************************
94  ! Global variables                                                          *
95  !     from par_mod and com_mod                                              *
96  ! nx,ny,nuvz,nwz       field dimensions in x,y and z direction              *
97  ! indmin               remembers the number of wind fields already treated  *
98  ! memind(2)            pointer, on which place the wind fields are stored   *
99  ! memtime(2) [s]       times of the wind fields, which are kept in memory   *
100  ! ndinterval [s]    time difference between the two wind fields read in     *
101  ! ldirect                 1 for forward, -1 for backward simulation         *
102  ! numbwf                  actual number of wind fields                      *
103  ! wftime(maxwf) [s]       times relative to beginning time of wind fields   *
104  !****************************************************************************
105
106#ifdef PERFTIMER
107  INTEGER millisecs_start, millisecs_stop, count_rate, count_max
108#endif
109!-----------------------------------------------------------------------------
110
111
112  ! Check, if wind fields are available for the current time step
113  !**************************************************************
114  nstop=0
115
116  if ((ldirect*wftime(1).gt.ldirect*itime).or. &
117       (ldirect*wftime(numbwf).lt.ldirect*itime)) then
118    write(*,*) 'FLEXPART WARNING: NO WIND FIELDS ARE AVAILABLE.'
119    write(*,*) 'A TRAJECTORY HAS TO BE TERMINATED.'
120    nstop=4
121    return
122  endif
123
124
125  if ((ldirect*memtime(1).le.ldirect*itime).and. &
126       (ldirect*memtime(2).gt.ldirect*itime)) then
127
128  ! The right wind fields are already in memory -> don't do anything
129  !*****************************************************************
130
131    continue
132
133  else if ((ldirect*memtime(2).le.ldirect*itime).and. &
134       (memtime(2).ne.999999999)) then
135
136
137  ! Current time is after 2nd wind field
138  ! -> Resort wind field pointers, so that current time is between 1st and 2nd
139  !***************************************************************************
140
141    memaux=memind(1)
142    memind(1)=memind(2)
143    memind(2)=memaux
144    memtime(1)=memtime(2)
145
146
147  ! Read a new wind field and store it on place memind(2)
148  !******************************************************
149
150    do indj=indmin,numbwf-1
151       if (ldirect*wftime(indj+1).gt.ldirect*itime) then
152
153#ifdef PERFTIMER
154         CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max)
155#endif
156          if (metdata_format.eq.ecmwf_metdata) then
157             call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh)
158             call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
159             call calcpar_ecmwf(memind(2),uuh,vvh,pvh)
160             call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format)
161             call verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh)
162             call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
163             memtime(2)=wftime(indj+1)
164             nstop = 1
165          endif
166          if (metdata_format.eq.gfs_metdata) then
167             call readwind_gfs(indj+1,memind(2),uuh,vvh,wwh)
168             call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
169             call calcpar_gfs(memind(2),uuh,vvh,pvh)
170             call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format)
171             call verttransform_gfs(memind(2),uuh,vvh,wwh,pvh)
172             call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
173             memtime(2)=wftime(indj+1)
174             nstop = 1
175          endif
176
177#ifdef PERFTIMER
178         CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max)
179         PRINT *, 'Wall time to process: ', TRIM(wfname(indj+1)), &
180                  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
181#endif
182
183          goto 40
184       endif
185    end do
186 40   indmin=indj
187
188  else
189
190  ! No wind fields, which can be used, are currently in memory
191  ! -> read both wind fields
192  !***********************************************************
193
194     do indj=indmin,numbwf-1
195        if ((ldirect*wftime(indj).le.ldirect*itime).and. &
196             (ldirect*wftime(indj+1).gt.ldirect*itime)) then
197           memind(1)=1
198           if (metdata_format.eq.ecmwf_metdata) then
199#ifdef PERFTIMER
200         CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max)
201#endif
202              call readwind_ecmwf(indj,memind(1),uuh,vvh,wwh)
203              call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn)
204              call calcpar_ecmwf(memind(1),uuh,vvh,pvh)
205              call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format)
206              call verttransform_ecmwf(memind(1),uuh,vvh,wwh,pvh)
207              call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
208              memtime(1)=wftime(indj)
209              memind(2)=2
210#ifdef PERFTIMER
211         CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max)
212         PRINT *, 'Wall time to process: ', TRIM(wfname(indj)), &
213                  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
214#endif
215
216#ifdef PERFTIMER
217         CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max)
218#endif
219              call readwind_ecmwf(indj+1,memind(2),uuh,vvh,wwh)
220              call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
221              call calcpar_ecmwf(memind(2),uuh,vvh,pvh)
222              call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format)
223              call verttransform_ecmwf(memind(2),uuh,vvh,wwh,pvh)
224              call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
225              memtime(2)=wftime(indj+1)
226              nstop = 1
227#ifdef PERFTIMER
228         CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max)
229         PRINT *, 'Wall time to process: ', TRIM(wfname(indj+1)), &
230                  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
231#endif
232           endif
233           if (metdata_format.eq.gfs_metdata) then
234#ifdef PERFTIMER
235         CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max)
236#endif
237              call readwind_gfs(indj,memind(1),uuh,vvh,wwh)
238              call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn)
239              call calcpar_gfs(memind(1),uuh,vvh,pvh)
240              call calcpar_nests(memind(1),uuhn,vvhn,pvhn,metdata_format)
241              call verttransform_gfs(memind(1),uuh,vvh,wwh,pvh)
242              call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
243              memtime(1)=wftime(indj)
244              memind(2)=2
245#ifdef PERFTIMER
246         CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max)
247         PRINT *, 'Wall time to process: ', TRIM(wfname(indj)), &
248                  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
249#endif
250
251#ifdef PERFTIMER
252         CALL SYSTEM_CLOCK(millisecs_start, count_rate, count_max)
253#endif
254              call readwind_gfs(indj+1,memind(2),uuh,vvh,wwh)
255              call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
256              call calcpar_gfs(memind(2),uuh,vvh,pvh)
257              call calcpar_nests(memind(2),uuhn,vvhn,pvhn,metdata_format)
258              call verttransform_gfs(memind(2),uuh,vvh,wwh,pvh)
259              call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
260              memtime(2)=wftime(indj+1)
261              nstop = 1
262#ifdef PERFTIMER
263         CALL SYSTEM_CLOCK(millisecs_stop, count_rate, count_max)
264         PRINT *, 'Wall time to process: ', TRIM(wfname(indj+1)), &
265                  ': ', (millisecs_stop-millisecs_start)/1000.0, ' seconds'
266#endif
267           endif
268           goto 60
269        endif
270     end do
271 60   indmin=indj
272
273  endif
274
275  lwindinterv=abs(memtime(2)-memtime(1))
276
277  if (lwindinterv.gt.idiffmax) nstop=3
278
279end subroutine getfields
280
281
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG