outgrid_init.f90 13.5 KB
Newer Older
Matthias Langer's avatar
 
Matthias Langer committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
!**********************************************************************
! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
!                                                                     *
! This file is part of FLEXPART.                                      *
!                                                                     *
! FLEXPART is free software: you can redistribute it and/or modify    *
! it under the terms of the GNU General Public License as published by*
! the Free Software Foundation, either version 3 of the License, or   *
! (at your option) any later version.                                 *
!                                                                     *
! FLEXPART is distributed in the hope that it will be useful,         *
! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
! GNU General Public License for more details.                        *
!                                                                     *
! You should have received a copy of the GNU General Public License   *
! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
!**********************************************************************
21
22
23
24
25
26
27
!                                                                     *
! DJM - 2017-05-09 - added #ifdef USE_MPIINPLACE cpp directive to     *
! enable allocation of a gridunc0 array if required by MPI code in    *
! mpi_mod.f90                                                         *
!                                                                     *
!**********************************************************************

Matthias Langer's avatar
 
Matthias Langer committed
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212

subroutine outgrid_init
  !
  !*****************************************************************************
  !                                                                            *
  !  This routine initializes the output grids                                 *
  !                                                                            *
  !     Author: A. Stohl                                                       *
  !                                                                            *
  !     7 August 2002                                                          *
  !                                                                            *
  !*****************************************************************************
  !                                                                            *
  ! Variables:                                                                 *
  !                                                                            *
  ! area               surface area of all output grid cells                   *
  ! areaeast           eastward facing wall area of all output grid cells      *
  ! areanorth          northward facing wall area of all output grid cells     *
  ! volume             volumes of all output grid cells                        *
  !                                                                            *
  !*****************************************************************************

  use flux_mod
  use oh_mod
  use unc_mod
  use outg_mod
  use par_mod
  use com_mod

  implicit none

  integer :: ix,jy,kz,i,nage,l,iix,jjy,ixp,jyp,i1,j1,j,ngrid
  integer :: ks,kp,stat
  real :: ylat,gridarea,ylatp,ylatm,hzone,cosfactm,cosfactp
  real :: xlon,xl,yl,ddx,ddy,rddx,rddy,p1,p2,p3,p4,xtn,ytn,oroh
  real,parameter :: eps=nxmax/3.e5


  ! Compute surface area and volume of each grid cell: area, volume;
  ! and the areas of the northward and eastward facing walls: areaeast, areanorth
  !***********************************************************************
  do jy=0,numygrid-1
    ylat=outlat0+(real(jy)+0.5)*dyout
    ylatp=ylat+0.5*dyout
    ylatm=ylat-0.5*dyout
    if ((ylatm.lt.0).and.(ylatp.gt.0.)) then
      hzone=dyout*r_earth*pi180
    else

  ! Calculate area of grid cell with formula M=2*pi*R*h*dx/360,
  ! see Netz, Formeln der Mathematik, 5. Auflage (1983), p.90
  !************************************************************

      cosfactp=cos(ylatp*pi180)
      cosfactm=cos(ylatm*pi180)
      if (cosfactp.lt.cosfactm) then
        hzone=sqrt(1-cosfactp**2)- &
             sqrt(1-cosfactm**2)
        hzone=hzone*r_earth
      else
        hzone=sqrt(1-cosfactm**2)- &
             sqrt(1-cosfactp**2)
        hzone=hzone*r_earth
      endif
    endif

  ! Surface are of a grid cell at a latitude ylat
  !**********************************************

    gridarea=2.*pi*r_earth*hzone*dxout/360.

    do ix=0,numxgrid-1
      area(ix,jy)=gridarea

  ! Volume = area x box height
  !***************************

      volume(ix,jy,1)=area(ix,jy)*outheight(1)
      areaeast(ix,jy,1)=dyout*r_earth*pi180*outheight(1)
      areanorth(ix,jy,1)=cos(ylat*pi180)*dxout*r_earth*pi180* &
           outheight(1)
      do kz=2,numzgrid
        areaeast(ix,jy,kz)=dyout*r_earth*pi180* &
             (outheight(kz)-outheight(kz-1))
        areanorth(ix,jy,kz)=cos(ylat*pi180)*dxout*r_earth*pi180* &
             (outheight(kz)-outheight(kz-1))
        volume(ix,jy,kz)=area(ix,jy)*(outheight(kz)-outheight(kz-1))
      end do
    end do
  end do




  !******************************************************************
  ! Determine average height of model topography in output grid cells
  !******************************************************************

  ! Loop over all output grid cells
  !********************************

  do jjy=0,numygrid-1
    do iix=0,numxgrid-1
      oroh=0.

  ! Take 100 samples of the topography in every grid cell
  !******************************************************

      do j1=1,10
        ylat=outlat0+(real(jjy)+real(j1)/10.-0.05)*dyout
        yl=(ylat-ylat0)/dy
        do i1=1,10
          xlon=outlon0+(real(iix)+real(i1)/10.-0.05)*dxout
          xl=(xlon-xlon0)/dx

  ! Determine the nest we are in
  !*****************************

          ngrid=0
          do j=numbnests,1,-1
            if ((xl.gt.xln(j)+eps).and.(xl.lt.xrn(j)-eps).and. &
                 (yl.gt.yln(j)+eps).and.(yl.lt.yrn(j)-eps)) then
              ngrid=j
              goto 43
            endif
          end do
43        continue

  ! Determine (nested) grid coordinates and auxiliary parameters used for interpolation
  !*****************************************************************************

          if (ngrid.gt.0) then
            xtn=(xl-xln(ngrid))*xresoln(ngrid)
            ytn=(yl-yln(ngrid))*yresoln(ngrid)
            ix=int(xtn)
            jy=int(ytn)
            ddy=ytn-real(jy)
            ddx=xtn-real(ix)
          else
            ix=int(xl)
            jy=int(yl)
            ddy=yl-real(jy)
            ddx=xl-real(ix)
          endif
          ixp=ix+1
          jyp=jy+1
          rddx=1.-ddx
          rddy=1.-ddy
          p1=rddx*rddy
          p2=ddx*rddy
          p3=rddx*ddy
          p4=ddx*ddy

          if (ngrid.gt.0) then
            oroh=oroh+p1*oron(ix ,jy ,ngrid) &
                 + p2*oron(ixp,jy ,ngrid) &
                 + p3*oron(ix ,jyp,ngrid) &
                 + p4*oron(ixp,jyp,ngrid)
          else
            oroh=oroh+p1*oro(ix ,jy) &
                 + p2*oro(ixp,jy) &
                 + p3*oro(ix ,jyp) &
                 + p4*oro(ixp,jyp)
          endif
        end do
      end do

  ! Divide by the number of samples taken
  !**************************************

      oroout(iix,jjy)=oroh/100.
    end do
  end do

  ! if necessary allocate flux fields
  if (iflux.eq.1) then
    allocate(flux(6,0:numxgrid-1,0:numygrid-1,numzgrid, &
         1:nspec,1:maxpointspec_act,1:nageclass),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate flux array '
  endif

  ! gridunc,griduncn        uncertainty of outputted concentrations
  allocate(gridunc(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, &
       maxpointspec_act,nclassunc,maxageclass),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
213
  if ((ldirect.gt.0).or.(WATERCYCLE)) then
214
    allocate(wetgridunc(0:numxgrid-1,0:numygrid-1,maxspec, &
Matthias Langer's avatar
 
Matthias Langer committed
215
       maxpointspec_act,nclassunc,maxageclass),stat=stat)
216
    if (stat.ne.0) write(*,*)'ERROR: could not allocate wetgridunc'
Matthias Langer's avatar
 
Matthias Langer committed
217
218
  allocate(drygridunc(0:numxgrid-1,0:numygrid-1,maxspec, &
       maxpointspec_act,nclassunc,maxageclass),stat=stat)
219
    if (stat.ne.0) write(*,*)'ERROR: could not allocate drygridunc'
Matthias Langer's avatar
 
Matthias Langer committed
220
  endif
221
222
223
  if (WATERCYCLE) then
      allocate(uul(0:nxmax-1,0:nymax-1,nuvzmax,2))
      allocate(vvl(0:nxmax-1,0:nymax-1,nuvzmax,2))
Sabine's avatar
Sabine committed
224
225
      allocate(uuln(0:nxmaxn-1,0:nymaxn-1,nuvzmax,2,maxnests))
      allocate(vvln(0:nxmaxn-1,0:nymaxn-1,nuvzmax,2,maxnests))
226
  endif
227
228
229

! Extra field for totals at MPI root process
  if (lroot.and.mpi_mode.gt.0) then
230
231
232
233
234
235
236
237
238

#ifdef USE_MPIINPLACE
#else
    ! If MPI_IN_PLACE option is not used in mpi_mod.f90::mpif_tm_reduce_grid(),
    ! then an aux array is needed for parallel grid reduction
    allocate(gridunc0(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, &
         maxpointspec_act,nclassunc,maxageclass),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc0'
#endif
239
    if ((ldirect.gt.0).or.(WATERCYCLE)) then
240
241
242
243
244
245
246
247
248
249
250
251
252
      allocate(wetgridunc0(0:numxgrid-1,0:numygrid-1,maxspec, &
           maxpointspec_act,nclassunc,maxageclass),stat=stat)
      if (stat.ne.0) write(*,*)'ERROR: could not allocate wetgridunc0'
      allocate(drygridunc0(0:numxgrid-1,0:numygrid-1,maxspec, &
           maxpointspec_act,nclassunc,maxageclass),stat=stat)
      if (stat.ne.0) write(*,*)'ERROR: could not allocate drygridunc0'
    endif
! allocate a dummy to avoid compilator complaints
  else if (.not.lroot.and.mpi_mode.gt.0) then
    allocate(wetgridunc0(1,1,1,1,1,1),stat=stat)
    allocate(drygridunc0(1,1,1,1,1,1),stat=stat)
  end if

Matthias Langer's avatar
 
Matthias Langer committed
253
254
255
  !write (*,*) 'Dimensions for fields', numxgrid,numygrid, &
  !     maxspec,maxpointspec_act,nclassunc,maxageclass

256
257
258
259
  if (lroot) then
    write (*,*) 'Allocating fields for global output (x,y): ', numxgrid,numygrid
    write (*,*) 'Allocating fields for nested output (x,y): ', numxgridn,numygridn 
  end if
Matthias Langer's avatar
 
Matthias Langer committed
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279

  ! allocate fields for concoutput with maximum dimension of outgrid
  ! and outgrid_nest

  allocate(gridsigma(0:max(numxgrid,numxgridn)-1, &
       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
  allocate(grid(0:max(numxgrid,numxgridn)-1, &
       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
  allocate(densityoutgrid(0:max(numxgrid,numxgridn)-1, &
       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'

  allocate(factor3d(0:max(numxgrid,numxgridn)-1, &
       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
  allocate(sparse_dump_r(max(numxgrid,numxgridn)* &
       max(numygrid,numygridn)*numzgrid),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
280
281
282
283
284

   allocate(sparse_dump_u(max(numxgrid,numxgridn)* &
       max(numygrid,numygridn)*numzgrid),stat=stat)
        if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'

Matthias Langer's avatar
 
Matthias Langer committed
285
286
287
288
289
  allocate(sparse_dump_i(max(numxgrid,numxgridn)* &
       max(numygrid,numygridn)*numzgrid),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'

  ! deposition fields are only allocated for forward runs
290
  if ((ldirect.gt.0).or.(WATERCYCLE)) then
Matthias Langer's avatar
 
Matthias Langer committed
291
292
293
294
295
296
297
298
299
300
301
302
     allocate(wetgridsigma(0:max(numxgrid,numxgridn)-1, &
          0:max(numygrid,numygridn)-1),stat=stat)
     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
     allocate(drygridsigma(0:max(numxgrid,numxgridn)-1, &
          0:max(numygrid,numygridn)-1),stat=stat)
     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
     allocate(wetgrid(0:max(numxgrid,numxgridn)-1, &
          0:max(numygrid,numygridn)-1),stat=stat)
     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
     allocate(drygrid(0:max(numxgrid,numxgridn)-1, &
          0:max(numygrid,numygridn)-1),stat=stat)
     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
303
     if (WATERCYCLE) then
Sabine's avatar
Sabine committed
304
       allocate(e_minus_p(0:nxmax-1,0:nymax-1,2))
Sabine's avatar
Sabine committed
305
       allocate(e_minus_p_nest(0:nxmaxn-1,0:nymaxn-1,2,maxnests))
306
     endif
Matthias Langer's avatar
 
Matthias Langer committed
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
  endif

  ! Initial condition field

  if (linit_cond.gt.0) then
    allocate(init_cond(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, &
         maxpointspec_act),stat=stat)
    if (stat.ne.0) write(*,*)'ERROR: could not allocate init_cond'
  endif

  !************************
  ! Initialize output grids
  !************************

  do ks=1,nspec
  do kp=1,maxpointspec_act
    do i=1,numreceptor
  ! Receptor points
      creceptor(i,ks)=0.
    end do
    do nage=1,nageclass
      do jy=0,numygrid-1
        do ix=0,numxgrid-1
          do l=1,nclassunc
  ! Deposition fields
332
            if ((ldirect.gt.0).or.(WATERCYCLE)) then
Matthias Langer's avatar
 
Matthias Langer committed
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
            wetgridunc(ix,jy,ks,kp,l,nage)=0.
            drygridunc(ix,jy,ks,kp,l,nage)=0.
            endif
            do kz=1,numzgrid
              if (iflux.eq.1) then
  ! Flux fields
                 do i=1,5
                   flux(i,ix,jy,kz,ks,kp,nage)=0.
                 end do
              endif
  ! Initial condition field
              if ((l.eq.1).and.(nage.eq.1).and.(linit_cond.gt.0)) &
                   init_cond(ix,jy,kz,ks,kp)=0.
  ! Concentration fields
              gridunc(ix,jy,kz,ks,kp,l,nage)=0.
            end do
          end do
        end do
      end do
    end do
  end do
  end do



end subroutine outgrid_init