source: flexpart.git/src/getfields_mpi.f90 @ 61e07ba

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 61e07ba was 78e62dc, checked in by flexpart <>, 9 years ago

New OH parameter in SPECIES files (now 3 instead of 2). New path to OH binariy files.

  • Property mode set to 100644
File size: 10.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,memstat)
23!                       i     o       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!  eso 2014:
45!    MPI version.
46!    If running with number of processes >= mpi_mod::read_grp_min,
47!    only one process (mpi_mod::lmpreader=.true.) does the actual reading, while
48!    the rest call this routine only to update memind, memstat etc.
49!
50!    If mpi_mod::lmp_sync=.true., uses 3 fields instead of 2, to allow reading
51!    the newest in the background.
52!
53!    Return memstat, which is the sum of
54!   
55!    memstat=16:     one new field to be read
56!    memstat=32:     two new fields must be read
57!    memstat=1,2,3:  position(s) in array to read next field
58!    memstat=0:      no new fields to be read
59!         
60!*****************************************************************************
61!                                                                            *
62! Variables:                                                                 *
63! lwindinterval [s]    time difference between the two wind fields read in   *
64! indj                 indicates the number of the wind field to be read in  *
65! indmin               remembers the number of wind fields already treated   *
66! memind(2[3])         pointer, on which place the wind fields are stored    *
67! memtime(2[3]) [s]    times of the wind fields, which are kept in memory    *
68! itime [s]            current time since start date of trajectory calcu-    *
69!                      lation                                                *
70! nstop                > 0, if trajectory has to be terminated               *
71! memstat              keep track of change of status for field pointers     *
72! nx,ny,nuvz,nwz       field dimensions in x,y and z direction               *
73! uu(0:nxmax,0:nymax,nuvzmax,numwfmem)   wind components in x-direction [m/s]*
74! vv(0:nxmax,0:nymax,nuvzmax,numwfmem)   wind components in y-direction [m/s]*
75! ww(0:nxmax,0:nymax,nwzmax,numwfmem)    wind components in z-direction      *
76!                                          [deltaeta/s]                      *
77! tt(0:nxmax,0:nymax,nuvzmax,numwfmem)   temperature [K]                     *
78! ps(0:nxmax,0:nymax,numwfmem)           surface pressure [Pa]               *
79!                                                                            *
80! Constants:                                                                 *
81! idiffmax             maximum allowable time difference between 2 wind      *
82!                      fields                                                *
83!                                                                            *
84!*****************************************************************************
85
86  use par_mod
87  use com_mod
88  use mpi_mod, only: lmpreader,lmp_use_reader,lmp_sync
89
90  implicit none
91
92  integer :: indj,itime,nstop,memaux,mindread
93! mindread: index of where to read 3rd field
94  integer, intent(out) :: memstat
95! keep track of 3rd field index. updated when new fields are read
96  integer, save :: mind3=0
97
98  real :: uuh(0:nxmax-1,0:nymax-1,nuvzmax)
99  real :: vvh(0:nxmax-1,0:nymax-1,nuvzmax)
100  real :: pvh(0:nxmax-1,0:nymax-1,nuvzmax)
101  real :: wwh(0:nxmax-1,0:nymax-1,nwzmax)
102  real :: uuhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
103  real :: vvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
104  real :: pvhn(0:nxmaxn-1,0:nymaxn-1,nuvzmax,maxnests)
105  real :: wwhn(0:nxmaxn-1,0:nymaxn-1,nwzmax,maxnests)
106
107  integer :: indmin = 1
108
109
110! Check, if wind fields are available for the current time step
111!**************************************************************
112
113  nstop=0
114  memstat=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    memstat=0
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! 2 fields in memory
142!*******************
143    if (lmp_sync) then
144      memaux=memind(1)
145      memind(1)=memind(2)
146      memind(2)=memaux
147      memtime(1)=memtime(2)
148      memstat=16+memind(2)
149      mindread=memind(2)
150
151    else
152! 3 fields in memory
153! Note that the reader process never needs to keep 3 fields in memory,
154! (2 is enough) so it is possible to save some memory
155!*********************************************************************
156      if (mind3.eq.32.or.mind3.eq.19) then
157        if (lmpreader) then
158          memind(1)=2
159          memind(2)=3
160          memind(3)=3
161        else
162          memind(1)=2
163          memind(2)=3
164          memind(3)=1
165        end if
166        memstat=17
167      else if (mind3.eq.17) then
168        if (lmpreader) then
169          memind(1)=3
170          memind(2)=1
171          memind(3)=1
172        else
173          memind(1)=3
174          memind(2)=1
175          memind(3)=2
176        end if
177        memstat=18
178      else if (mind3.eq.18) then
179        if (lmpreader) then
180          memind(1)=1
181          memind(2)=2
182          memind(3)=2
183        else
184          memind(1)=1
185          memind(2)=2
186          memind(3)=3
187        end if
188        memstat=19
189      else
190        write(*,*) '#### getfields_mpi.f90> ERROR: ', &
191             & 'unknown mind3, exiting.', mind3,' ####'
192        stop
193      end if
194      mind3=memstat
195      memtime(1)=memtime(2)
196      mindread=memind(3)
197    end if
198
199
200! Read a new wind field and store it on place memind(2)
201! or memind(3) for the 3-field read-ahead version
202!******************************************************
203    do indj=indmin,numbwf-1
204      if (ldirect*wftime(indj+1).gt.ldirect*itime) then
205        if (lmpreader.or..not. lmp_use_reader) then
206          call readwind(indj+1,mindread,uuh,vvh,wwh)
207          call readwind_nests(indj+1,mindread,uuhn,vvhn,wwhn)
208          call calcpar(mindread,uuh,vvh,pvh)
209          call calcpar_nests(mindread,uuhn,vvhn,pvhn)
210          call verttransform(mindread,uuh,vvh,wwh,pvh)
211          call verttransform_nests(mindread,uuhn,vvhn,wwhn,pvhn)
212        end if
213        memtime(2)=wftime(indj+1)
214        nstop = 1
215        goto 40
216      endif
217    end do
21840  indmin=indj
219
220  else
221
222! No wind fields, which can be used, are currently in memory
223! -> read both wind fields
224!***********************************************************
225    memstat=32
226
227    do indj=indmin,numbwf-1
228      if ((ldirect*wftime(indj).le.ldirect*itime).and. &
229           (ldirect*wftime(indj+1).gt.ldirect*itime)) then
230        memind(1)=1
231        if (lmpreader.or..not.lmp_use_reader) then
232          call readwind(indj,memind(1),uuh,vvh,wwh)
233          call readwind_nests(indj,memind(1),uuhn,vvhn,wwhn)
234          call calcpar(memind(1),uuh,vvh,pvh)
235          call calcpar_nests(memind(1),uuhn,vvhn,pvhn)
236          call verttransform(memind(1),uuh,vvh,wwh,pvh)
237          call verttransform_nests(memind(1),uuhn,vvhn,wwhn,pvhn)
238        end if
239        memtime(1)=wftime(indj)
240        memind(2)=2
241        if (lmpreader.or..not.lmp_use_reader) then
242          call readwind(indj+1,memind(2),uuh,vvh,wwh)
243          call readwind_nests(indj+1,memind(2),uuhn,vvhn,wwhn)
244          call calcpar(memind(2),uuh,vvh,pvh)
245          call calcpar_nests(memind(2),uuhn,vvhn,pvhn)
246          call verttransform(memind(2),uuh,vvh,wwh,pvh)
247          call verttransform_nests(memind(2),uuhn,vvhn,wwhn,pvhn)
248        end if
249        memtime(2)=wftime(indj+1)
250! DEV: not used
251        if (.not.lmp_sync) memind(3)=3 ! indicate position of next read
252        nstop = 1
253        goto 60
254      endif
255    end do
25660  indmin=indj
257
258    mind3=memstat
259
260  endif
261
262  lwindinterv=abs(memtime(2)-memtime(1))
263
264  if (lwindinterv.gt.idiffmax) nstop=3
265
266end subroutine getfields
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG