Changeset 8a65cb0 in flexpart.git for src/getfields.f90


Ignore:
Timestamp:
Mar 2, 2015, 3:11:55 PM (9 years ago)
Author:
Espen Sollum ATMOS <espen@…>
Branches:
master, 10.4.1_pesei, GFS_025, bugfixes+enhancements, dev, release-10, release-10.4.1, scaling-bug, univie
Children:
1d207bb
Parents:
60403cd
Message:

Added code, makefile for dev branch

File:
1 edited

Legend:

Unmodified
Added
Removed
  • src/getfields.f90

    re200b7a r8a65cb0  
    2121
    2222subroutine 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   !*****************************************************************************
     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!*****************************************************************************
    6666
    6767  use par_mod
     
    8484
    8585
    86   ! Check, if wind fields are available for the current time step
    87   !**************************************************************
     86! Check, if wind fields are available for the current time step
     87!**************************************************************
    8888
    8989  nstop=0
     
    101101       (ldirect*memtime(2).gt.ldirect*itime)) then
    102102
    103   ! The right wind fields are already in memory -> don't do anything
    104   !*****************************************************************
     103! The right wind fields are already in memory -> don't do anything
     104!*****************************************************************
    105105
    106106    continue
     
    110110
    111111
    112   ! Current time is after 2nd wind field
    113   ! -> Resort wind field pointers, so that current time is between 1st and 2nd
    114   !***************************************************************************
     112! Current time is after 2nd wind field
     113! -> Resort wind field pointers, so that current time is between 1st and 2nd
     114!***************************************************************************
    115115
    116116    memaux=memind(1)
     
    120120
    121121
    122   ! Read a new wind field and store it on place memind(2)
    123   !******************************************************
     122! Read a new wind field and store it on place memind(2)
     123!******************************************************
    124124
    125125    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
     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
    137137    end do
    138  40   indmin=indj
     13840  indmin=indj
    139139
    140140  else
    141141
    142   ! No wind fields, which can be used, are currently in memory
    143   ! -> read both wind fields
    144   !***********************************************************
     142! No wind fields, which can be used, are currently in memory
     143! -> read both wind fields
     144!***********************************************************
    145145
    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
    169  60   indmin=indj
     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
    170170
    171171  endif
Note: See TracChangeset for help on using the changeset viewer.
hosted by ZAMG