source: flexpart.git/src/getfields.f90 @ 8a65cb0

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 8a65cb0 was 8a65cb0, checked in by Espen Sollum ATMOS <espen@…>, 9 years ago

Added code, makefile for dev branch

  • Property mode set to 100644
File size: 8.0 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)
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! Variables:                                                                 *
46! lwindinterval [s]    time difference between the two wind fields read in   *
47! indj                 indicates the number of the wind field to be read in  *
48! indmin               remembers the number of wind fields already treated   *
49! memind(2)            pointer, on which place the wind fields are stored    *
50! memtime(2) [s]       times of the wind fields, which are kept in memory    *
51! itime [s]            current time since start date of trajectory calcu-    *
52!                      lation                                                *
53! nstop                > 0, if trajectory has to be terminated               *
54! nx,ny,nuvz,nwz       field dimensions in x,y and z direction               *
55! uu(0:nxmax,0:nymax,nuvzmax,2)   wind components in x-direction [m/s]       *
56! vv(0:nxmax,0:nymax,nuvzmax,2)   wind components in y-direction [m/s]       *
57! ww(0:nxmax,0:nymax,nwzmax,2)    wind components in z-direction [deltaeta/s]*
58! tt(0:nxmax,0:nymax,nuvzmax,2)   temperature [K]                            *
59! ps(0:nxmax,0:nymax,2)           surface pressure [Pa]                      *
60!                                                                            *
61! Constants:                                                                 *
62! idiffmax             maximum allowable time difference between 2 wind      *
63!                      fields                                                *
64!                                                                            *
65!*****************************************************************************
66
67  use par_mod
68  use com_mod
69
70  implicit none
71
72  integer :: indj,itime,nstop,memaux
73
74  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
75  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
76  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
77  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
78  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
79  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
80  real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
81  real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
82
83  integer :: indmin = 1
84
85
86! Check, if wind fields are available for the current time step
87!**************************************************************
88
89  nstop=0
90
91  if ((ldirect*wftime(1).gt.ldirect*itime).or. &
92       (ldirect*wftime(numbwf).lt.ldirect*itime)) then
93    write(*,*) 'FLEXPART WARNING: NO WIND FIELDS ARE AVAILABLE.'
94    write(*,*) 'A TRAJECTORY HAS TO BE TERMINATED.'
95    nstop=4
96    return
97  endif
98
99
100  if ((ldirect*memtime(1).le.ldirect*itime).and. &
101       (ldirect*memtime(2).gt.ldirect*itime)) then
102
103! The right wind fields are already in memory -> don't do anything
104!*****************************************************************
105
106    continue
107
108  else if ((ldirect*memtime(2).le.ldirect*itime).and. &
109       (memtime(2).ne.999999999)) then
110
111
112! Current time is after 2nd wind field
113! -> Resort wind field pointers, so that current time is between 1st and 2nd
114!***************************************************************************
115
116    memaux=memind(1)
117    memind(1)=memind(2)
118    memind(2)=memaux
119    memtime(1)=memtime(2)
120
121
122! Read a new wind field and store it on place memind(2)
123!******************************************************
124
125    do indj=indmin,numbwf-1
126      if (ldirect*wftime(indj+1).gt.ldirect*itime) then
127        call readwind(indj+1,memind(2),uuh,vvh,wwh)
128        call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
129        call calcpar(memind(2),uuh,vvh,pvh)
130        call calcpar_nests(memind(2),uuhn,vvhn,pvhn)
131        call verttransform(memind(2),uuh,vvh,wwh,pvh)
132        call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
133        memtime(2)=wftime(indj+1)
134        nstop = 1
135        goto 40
136      endif
137    end do
13840  indmin=indj
139
140  else
141
142! No wind fields, which can be used, are currently in memory
143! -> read both wind fields
144!***********************************************************
145
146    do indj=indmin,numbwf-1
147      if ((ldirect*wftime(indj).le.ldirect*itime).and. &
148           (ldirect*wftime(indj+1).gt.ldirect*itime)) then
149        memind(1)=1
150        call readwind(indj,memind(1),uuh,vvh,wwh)
151        call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn)
152        call calcpar(memind(1),uuh,vvh,pvh)
153        call calcpar_nests(memind(1),uuhn,vvhn,pvhn)
154        call verttransform(memind(1),uuh,vvh,wwh,pvh)
155        call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
156        memtime(1)=wftime(indj)
157        memind(2)=2
158        call readwind(indj+1,memind(2),uuh,vvh,wwh)
159        call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
160        call calcpar(memind(2),uuh,vvh,pvh)
161        call calcpar_nests(memind(2),uuhn,vvhn,pvhn)
162        call verttransform(memind(2),uuh,vvh,wwh,pvh)
163        call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
164        memtime(2)=wftime(indj+1)
165        nstop = 1
166        goto 60
167      endif
168    end do
16960  indmin=indj
170
171  endif
172
173  lwindinterv=abs(memtime(2)-memtime(1))
174
175  if (lwindinterv.gt.idiffmax) nstop=3
176
177end subroutine getfields
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG