fluxoutput.f90 9.71 KB
Newer Older
1
2
! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
! SPDX-License-Identifier: GPL-3.0-or-later
3

Matthias Langer's avatar
 
Matthias Langer committed
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
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
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
subroutine fluxoutput(itime)
  !                        i
  !*****************************************************************************
  !                                                                            *
  !     Output of the gridded fluxes.                                          *
  !     Eastward, westward, northward, southward, upward and downward gross    *
  !     fluxes are written to output file in either sparse matrix or grid dump *
  !     format, whichever is more efficient.                                   *
  !                                                                            *
  !     Author: A. Stohl                                                       *
  !                                                                            *
  !     04 April 2000                                                          *
  !                                                                            *
  !*****************************************************************************
  !                                                                            *
  ! Variables:                                                                 *
  ! ncellse         number of cells with non-zero values for eastward fluxes   *
  ! sparsee         .true. if in sparse matrix format, else .false.            *
  !                                                                            *
  !*****************************************************************************

  use flux_mod
  use outg_mod
  use par_mod
  use com_mod

  implicit none

  real(kind=dp) :: jul
  integer :: itime,ix,jy,kz,k,nage,jjjjmmdd,ihmmss,kp,i
  integer :: ncellse(maxspec,maxageclass),ncellsw(maxspec,maxageclass)
  integer :: ncellss(maxspec,maxageclass),ncellsn(maxspec,maxageclass)
  integer :: ncellsu(maxspec,maxageclass),ncellsd(maxspec,maxageclass)
  logical :: sparsee(maxspec,maxageclass),sparsew(maxspec,maxageclass)
  logical :: sparses(maxspec,maxageclass),sparsen(maxspec,maxageclass)
  logical :: sparseu(maxspec,maxageclass),sparsed(maxspec,maxageclass)
  character :: adate*8,atime*6


  ! Determine current calendar date, needed for the file name
  !**********************************************************

  jul=bdate+real(itime,kind=dp)/86400._dp
  call caldate(jul,jjjjmmdd,ihmmss)
  write(adate,'(i8.8)') jjjjmmdd
  write(atime,'(i6.6)') ihmmss


  open(unitflux,file=path(2)(1:length(2))//'grid_flux_'//adate// &
       atime,form='unformatted')

  !**************************************************************
  ! Check, whether output of full grid or sparse matrix format is
  ! more efficient in terms of storage space. This is checked for
  ! every species and for every age class
  !**************************************************************

  do k=1,nspec
    do nage=1,nageclass
      ncellse(k,nage)=0
      ncellsw(k,nage)=0
      ncellsn(k,nage)=0
      ncellss(k,nage)=0
      ncellsu(k,nage)=0
      ncellsd(k,nage)=0
    end do
  end do

  do k=1,nspec
  do kp=1,maxpointspec_act
    do nage=1,nageclass
      do jy=0,numygrid-1
        do ix=0,numxgrid-1
          do kz=1,numzgrid
            if (flux(2,ix,jy,kz,k,kp,nage).gt.0) ncellse(k,nage)= &
                 ncellse(k,nage)+1
            if (flux(1,ix,jy,kz,k,kp,nage).gt.0) ncellsw(k,nage)= &
                 ncellsw(k,nage)+1
            if (flux(4,ix,jy,kz,k,kp,nage).gt.0) ncellsn(k,nage)= &
                 ncellsn(k,nage)+1
            if (flux(3,ix,jy,kz,k,kp,nage).gt.0) ncellss(k,nage)= &
                 ncellss(k,nage)+1
            if (flux(5,ix,jy,kz,k,kp,nage).gt.0) ncellsu(k,nage)= &
                 ncellsu(k,nage)+1
            if (flux(6,ix,jy,kz,k,kp,nage).gt.0) ncellsd(k,nage)= &
                 ncellsd(k,nage)+1
          end do
        end do
      end do
    end do
  end do
  end do

  ! Output in sparse matrix format more efficient, if less than
  ! 2/5 of all cells contains concentrations>0
  !************************************************************

  do k=1,nspec
    do nage=1,nageclass
      if (4*ncellse(k,nage).lt.numxgrid*numygrid*numzgrid) then
        sparsee(k,nage)=.true.
      else
        sparsee(k,nage)=.false.
      endif
      if (4*ncellsw(k,nage).lt.numxgrid*numygrid*numzgrid) then
        sparsew(k,nage)=.true.
      else
        sparsew(k,nage)=.false.
      endif
      if (4*ncellsn(k,nage).lt.numxgrid*numygrid*numzgrid) then
        sparsen(k,nage)=.true.
      else
        sparsen(k,nage)=.false.
      endif
      if (4*ncellss(k,nage).lt.numxgrid*numygrid*numzgrid) then
        sparses(k,nage)=.true.
      else
        sparses(k,nage)=.false.
      endif
      if (4*ncellsu(k,nage).lt.numxgrid*numygrid*numzgrid) then
        sparseu(k,nage)=.true.
      else
        sparseu(k,nage)=.false.
      endif
      if (4*ncellsd(k,nage).lt.numxgrid*numygrid*numzgrid) then
        sparsed(k,nage)=.true.
      else
        sparsed(k,nage)=.false.
      endif
    end do
  end do



  ! Flux output: divide by area and time to get flux in ng/m2/s
  !************************************************************

  write(unitflux) itime
  do k=1,nspec
  do kp=1,maxpointspec_act
    do nage=1,nageclass

      if (sparsee(k,nage)) then
        write(unitflux) 1
        do kz=1,numzgrid
          do jy=0,numygrid-1
            do ix=0,numxgrid-1
              if (flux(2,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
                   flux(2,ix,jy,kz,k,kp,nage)/areaeast(ix,jy,kz)/outstep
            end do
          end do
        end do
        write(unitflux) -999,999.
      else
        write(unitflux) 2
        do kz=1,numzgrid
          do ix=0,numxgrid-1
            write(unitflux) (1.e12*flux(2,ix,jy,kz,k,kp,nage)/ &
                 areaeast(ix,jy,kz)/outstep,jy=0,numygrid-1)
          end do
        end do
      endif

      if (sparsew(k,nage)) then
        write(unitflux) 1
        do kz=1,numzgrid
          do jy=0,numygrid-1
            do ix=0,numxgrid-1
              if (flux(1,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
                   flux(1,ix,jy,kz,k,kp,nage)/areaeast(ix,jy,kz)/outstep
            end do
          end do
        end do
        write(unitflux) -999,999.
      else
        write(unitflux) 2
        do kz=1,numzgrid
          do ix=0,numxgrid-1
            write(unitflux) (1.e12*flux(1,ix,jy,kz,k,kp,nage)/ &
                 areaeast(ix,jy,kz)/outstep,jy=0,numygrid-1)
          end do
        end do
      endif

      if (sparses(k,nage)) then
        write(unitflux) 1
        do kz=1,numzgrid
          do jy=0,numygrid-1
            do ix=0,numxgrid-1
              if (flux(3,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
                   flux(3,ix,jy,kz,k,kp,nage)/areanorth(ix,jy,kz)/outstep
            end do
          end do
        end do
        write(unitflux) -999,999.
      else
        write(unitflux) 2
        do kz=1,numzgrid
          do ix=0,numxgrid-1
            write(unitflux) (1.e12*flux(3,ix,jy,kz,k,kp,nage)/ &
                 areanorth(ix,jy,kz)/outstep,jy=0,numygrid-1)
          end do
        end do
      endif

      if (sparsen(k,nage)) then
        write(unitflux) 1
        do kz=1,numzgrid
          do jy=0,numygrid-1
            do ix=0,numxgrid-1 ! north
              if (flux(4,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
                   flux(4,ix,jy,kz,k,kp,nage)/areanorth(ix,jy,kz)/outstep
            end do
          end do
        end do
        write(unitflux) -999,999.
      else
        write(unitflux) 2
        do kz=1,numzgrid
          do ix=0,numxgrid-1
            write(unitflux) (1.e12*flux(4,ix,jy,kz,k,kp,nage)/ &
                 areanorth(ix,jy,kz)/outstep,jy=0,numygrid-1)
          end do
        end do
      endif

      if (sparseu(k,nage)) then
        write(unitflux) 1
        do kz=1,numzgrid
          do jy=0,numygrid-1
            do ix=0,numxgrid-1
              if (flux(5,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
                   flux(5,ix,jy,kz,k,kp,nage)/area(ix,jy)/outstep
            end do
          end do
        end do
        write(unitflux) -999,999.
      else
        write(unitflux) 2
        do kz=1,numzgrid
          do ix=0,numxgrid-1
            write(unitflux) (1.e12*flux(5,ix,jy,kz,k,kp,nage)/ &
                 area(ix,jy)/outstep,jy=0,numygrid-1)
          end do
        end do
      endif

      if (sparsed(k,nage)) then
        write(unitflux) 1
        do kz=1,numzgrid
          do jy=0,numygrid-1
            do ix=0,numxgrid-1
              if (flux(6,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
                   flux(6,ix,jy,kz,k,kp,nage)/area(ix,jy)/outstep
            end do
          end do
        end do
        write(unitflux) -999,999.
      else
        write(unitflux) 2
        do kz=1,numzgrid
          do ix=0,numxgrid-1
            write(unitflux) (1.e12*flux(6,ix,jy,kz,k,kp,nage)/ &
                 area(ix,jy)/outstep,jy=0,numygrid-1)
          end do
        end do
      endif

    end do
  end do
  end do


  close(unitflux)


  ! Reinitialization of grid
  !*************************

  do k=1,nspec
  do kp=1,maxpointspec_act
    do jy=0,numygrid-1
      do ix=0,numxgrid-1
          do kz=1,numzgrid
            do nage=1,nageclass
              do i=1,6
                flux(i,ix,jy,kz,k,kp,nage)=0.
              end do
            end do
          end do
      end do
    end do
  end do
  end do


end subroutine fluxoutput