source: flexpart.git/src/fluxoutput.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 9.6 KB
Line 
1subroutine fluxoutput(itime)
2  !                        i
3  !*****************************************************************************
4  !                                                                            *
5  !     Output of the gridded fluxes.                                          *
6  !     Eastward, westward, northward, southward, upward and downward gross    *
7  !     fluxes are written to output file in either sparse matrix or grid dump *
8  !     format, whichever is more efficient.                                   *
9  !                                                                            *
10  !     Author: A. Stohl                                                       *
11  !                                                                            *
12  !     04 April 2000                                                          *
13  !                                                                            *
14  !*****************************************************************************
15  !                                                                            *
16  ! Variables:                                                                 *
17  ! ncellse         number of cells with non-zero values for eastward fluxes   *
18  ! sparsee         .true. if in sparse matrix format, else .false.            *
19  !                                                                            *
20  !*****************************************************************************
21
22  use flux_mod
23  use outg_mod
24  use par_mod
25  use com_mod
26
27  implicit none
28
29  real(kind=dp) :: jul
30  integer :: itime,ix,jy,kz,k,nage,jjjjmmdd,ihmmss,kp,i
31  integer :: ncellse(maxspec,maxageclass),ncellsw(maxspec,maxageclass)
32  integer :: ncellss(maxspec,maxageclass),ncellsn(maxspec,maxageclass)
33  integer :: ncellsu(maxspec,maxageclass),ncellsd(maxspec,maxageclass)
34  logical :: sparsee(maxspec,maxageclass),sparsew(maxspec,maxageclass)
35  logical :: sparses(maxspec,maxageclass),sparsen(maxspec,maxageclass)
36  logical :: sparseu(maxspec,maxageclass),sparsed(maxspec,maxageclass)
37  character :: adate*8,atime*6
38
39
40  ! Determine current calendar date, needed for the file name
41  !**********************************************************
42
43  jul=bdate+real(itime,kind=dp)/86400._dp
44  call caldate(jul,jjjjmmdd,ihmmss)
45  write(adate,'(i8.8)') jjjjmmdd
46  write(atime,'(i6.6)') ihmmss
47
48
49  open(unitflux,file=path(2)(1:length(2))//'grid_flux_'//adate// &
50       atime,form='unformatted')
51
52  !**************************************************************
53  ! Check, whether output of full grid or sparse matrix format is
54  ! more efficient in terms of storage space. This is checked for
55  ! every species and for every age class
56  !**************************************************************
57
58  do k=1,nspec
59    do nage=1,nageclass
60      ncellse(k,nage)=0
61      ncellsw(k,nage)=0
62      ncellsn(k,nage)=0
63      ncellss(k,nage)=0
64      ncellsu(k,nage)=0
65      ncellsd(k,nage)=0
66    end do
67  end do
68
69  do k=1,nspec
70  do kp=1,maxpointspec_act
71    do nage=1,nageclass
72      do jy=0,numygrid-1
73        do ix=0,numxgrid-1
74          do kz=1,numzgrid
75            if (flux(2,ix,jy,kz,k,kp,nage).gt.0) ncellse(k,nage)= &
76                 ncellse(k,nage)+1
77            if (flux(1,ix,jy,kz,k,kp,nage).gt.0) ncellsw(k,nage)= &
78                 ncellsw(k,nage)+1
79            if (flux(4,ix,jy,kz,k,kp,nage).gt.0) ncellsn(k,nage)= &
80                 ncellsn(k,nage)+1
81            if (flux(3,ix,jy,kz,k,kp,nage).gt.0) ncellss(k,nage)= &
82                 ncellss(k,nage)+1
83            if (flux(5,ix,jy,kz,k,kp,nage).gt.0) ncellsu(k,nage)= &
84                 ncellsu(k,nage)+1
85            if (flux(6,ix,jy,kz,k,kp,nage).gt.0) ncellsd(k,nage)= &
86                 ncellsd(k,nage)+1
87          end do
88        end do
89      end do
90    end do
91  end do
92  end do
93
94  ! Output in sparse matrix format more efficient, if less than
95  ! 2/5 of all cells contains concentrations>0
96  !************************************************************
97
98  do k=1,nspec
99    do nage=1,nageclass
100      if (4*ncellse(k,nage).lt.numxgrid*numygrid*numzgrid) then
101        sparsee(k,nage)=.true.
102      else
103        sparsee(k,nage)=.false.
104      endif
105      if (4*ncellsw(k,nage).lt.numxgrid*numygrid*numzgrid) then
106        sparsew(k,nage)=.true.
107      else
108        sparsew(k,nage)=.false.
109      endif
110      if (4*ncellsn(k,nage).lt.numxgrid*numygrid*numzgrid) then
111        sparsen(k,nage)=.true.
112      else
113        sparsen(k,nage)=.false.
114      endif
115      if (4*ncellss(k,nage).lt.numxgrid*numygrid*numzgrid) then
116        sparses(k,nage)=.true.
117      else
118        sparses(k,nage)=.false.
119      endif
120      if (4*ncellsu(k,nage).lt.numxgrid*numygrid*numzgrid) then
121        sparseu(k,nage)=.true.
122      else
123        sparseu(k,nage)=.false.
124      endif
125      if (4*ncellsd(k,nage).lt.numxgrid*numygrid*numzgrid) then
126        sparsed(k,nage)=.true.
127      else
128        sparsed(k,nage)=.false.
129      endif
130    end do
131  end do
132
133
134
135  ! Flux output: divide by area and time to get flux in ng/m2/s
136  !************************************************************
137
138  write(unitflux) itime
139  do k=1,nspec
140  do kp=1,maxpointspec_act
141    do nage=1,nageclass
142
143      if (sparsee(k,nage)) then
144        write(unitflux) 1
145        do kz=1,numzgrid
146          do jy=0,numygrid-1
147            do ix=0,numxgrid-1
148              if (flux(2,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
149                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
150                   flux(2,ix,jy,kz,k,kp,nage)/areaeast(ix,jy,kz)/outstep
151            end do
152          end do
153        end do
154        write(unitflux) -999,999.
155      else
156        write(unitflux) 2
157        do kz=1,numzgrid
158          do ix=0,numxgrid-1
159            write(unitflux) (1.e12*flux(2,ix,jy,kz,k,kp,nage)/ &
160                 areaeast(ix,jy,kz)/outstep,jy=0,numygrid-1)
161          end do
162        end do
163      endif
164
165      if (sparsew(k,nage)) then
166        write(unitflux) 1
167        do kz=1,numzgrid
168          do jy=0,numygrid-1
169            do ix=0,numxgrid-1
170              if (flux(1,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
171                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
172                   flux(1,ix,jy,kz,k,kp,nage)/areaeast(ix,jy,kz)/outstep
173            end do
174          end do
175        end do
176        write(unitflux) -999,999.
177      else
178        write(unitflux) 2
179        do kz=1,numzgrid
180          do ix=0,numxgrid-1
181            write(unitflux) (1.e12*flux(1,ix,jy,kz,k,kp,nage)/ &
182                 areaeast(ix,jy,kz)/outstep,jy=0,numygrid-1)
183          end do
184        end do
185      endif
186
187      if (sparses(k,nage)) then
188        write(unitflux) 1
189        do kz=1,numzgrid
190          do jy=0,numygrid-1
191            do ix=0,numxgrid-1
192              if (flux(3,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
193                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
194                   flux(3,ix,jy,kz,k,kp,nage)/areanorth(ix,jy,kz)/outstep
195            end do
196          end do
197        end do
198        write(unitflux) -999,999.
199      else
200        write(unitflux) 2
201        do kz=1,numzgrid
202          do ix=0,numxgrid-1
203            write(unitflux) (1.e12*flux(3,ix,jy,kz,k,kp,nage)/ &
204                 areanorth(ix,jy,kz)/outstep,jy=0,numygrid-1)
205          end do
206        end do
207      endif
208
209      if (sparsen(k,nage)) then
210        write(unitflux) 1
211        do kz=1,numzgrid
212          do jy=0,numygrid-1
213            do ix=0,numxgrid-1 ! north
214              if (flux(4,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
215                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
216                   flux(4,ix,jy,kz,k,kp,nage)/areanorth(ix,jy,kz)/outstep
217            end do
218          end do
219        end do
220        write(unitflux) -999,999.
221      else
222        write(unitflux) 2
223        do kz=1,numzgrid
224          do ix=0,numxgrid-1
225            write(unitflux) (1.e12*flux(4,ix,jy,kz,k,kp,nage)/ &
226                 areanorth(ix,jy,kz)/outstep,jy=0,numygrid-1)
227          end do
228        end do
229      endif
230
231      if (sparseu(k,nage)) then
232        write(unitflux) 1
233        do kz=1,numzgrid
234          do jy=0,numygrid-1
235            do ix=0,numxgrid-1
236              if (flux(5,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
237                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
238                   flux(5,ix,jy,kz,k,kp,nage)/area(ix,jy)/outstep
239            end do
240          end do
241        end do
242        write(unitflux) -999,999.
243      else
244        write(unitflux) 2
245        do kz=1,numzgrid
246          do ix=0,numxgrid-1
247            write(unitflux) (1.e12*flux(5,ix,jy,kz,k,kp,nage)/ &
248                 area(ix,jy)/outstep,jy=0,numygrid-1)
249          end do
250        end do
251      endif
252
253      if (sparsed(k,nage)) then
254        write(unitflux) 1
255        do kz=1,numzgrid
256          do jy=0,numygrid-1
257            do ix=0,numxgrid-1
258              if (flux(6,ix,jy,kz,k,kp,nage).gt.0.) write(unitflux) &
259                   ix+jy*numxgrid+kz*numxgrid*numygrid,1.e12* &
260                   flux(6,ix,jy,kz,k,kp,nage)/area(ix,jy)/outstep
261            end do
262          end do
263        end do
264        write(unitflux) -999,999.
265      else
266        write(unitflux) 2
267        do kz=1,numzgrid
268          do ix=0,numxgrid-1
269            write(unitflux) (1.e12*flux(6,ix,jy,kz,k,kp,nage)/ &
270                 area(ix,jy)/outstep,jy=0,numygrid-1)
271          end do
272        end do
273      endif
274
275    end do
276  end do
277  end do
278
279
280  close(unitflux)
281
282
283  ! Reinitialization of grid
284  !*************************
285
286  do k=1,nspec
287  do kp=1,maxpointspec_act
288    do jy=0,numygrid-1
289      do ix=0,numxgrid-1
290          do kz=1,numzgrid
291            do nage=1,nageclass
292              do i=1,6
293                flux(i,ix,jy,kz,k,kp,nage)=0.
294              end do
295            end do
296          end do
297      end do
298    end do
299  end do
300  end do
301
302
303end subroutine fluxoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG