source: flexpart.git/src/concoutput_nest.f90

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

add SPDX-License-Identifier to all .f90 files

  • Property mode set to 100644
File size: 21.6 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine concoutput_nest(itime,outnum)
5  !                        i     i
6  !*****************************************************************************
7  !                                                                            *
8  !     Output of the concentration grid and the receptor concentrations.      *
9  !                                                                            *
10  !     Author: A. Stohl                                                       *
11  !                                                                            *
12  !     24 May 1995                                                            *
13  !                                                                            *
14  !     13 April 1999, Major update: if output size is smaller, dump output    *
15  !                    in sparse matrix format; additional output of           *
16  !                    uncertainty                                             *
17  !                                                                            *
18  !     05 April 2000, Major update: output of age classes; output for backward*
19  !                    runs is time spent in grid cell times total mass of     *
20  !                    species.                                                *
21  !                                                                            *
22  !     17 February 2002, Appropriate dimensions for backward and forward runs *
23  !                       are now specified in file par_mod                    *
24  !                                                                            *
25  !     June 2006, write grid in sparse matrix with a single write command     *
26  !                in order to save disk space                                 *
27  !                                                                            *
28  !     2008 new sparse matrix format                                          *
29  !                                                                            *
30  !*****************************************************************************
31  !                                                                            *
32  ! Variables:                                                                 *
33  ! outnum          number of samples                                          *
34  ! ncells          number of cells with non-zero concentrations               *
35  ! sparse          .true. if in sparse matrix format, else .false.            *
36  ! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
37  !                                                                            *
38  !*****************************************************************************
39
40  use unc_mod
41  use point_mod
42  use outg_mod
43  use par_mod
44  use com_mod
45  use mean_mod
46
47  implicit none
48
49  real(kind=dp) :: jul
50  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
51  integer :: sp_count_i,sp_count_r
52  real :: sp_fact
53  real :: outnum,densityoutrecept(maxreceptor),xl,yl
54! RLT
55  real :: densitydryrecept(maxreceptor)
56  real :: factor_dryrecept(maxreceptor)
57
58  !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
59  !    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
60  !    +    maxageclass)
61  !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
62  !    +       maxageclass)
63  !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
64  !    +       maxpointspec_act,maxageclass)
65  !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
66  !    +       maxpointspec_act,maxageclass),
67  !    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
68  !    +     maxpointspec_act,maxageclass),
69  !    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
70  !    +     maxpointspec_act,maxageclass)
71  !real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
72  !real sparse_dump_r(numxgrid*numygrid*numzgrid)
73  !integer sparse_dump_i(numxgrid*numygrid*numzgrid)
74
75  !real sparse_dump_u(numxgrid*numygrid*numzgrid)
76  real(dep_prec) :: auxgrid(nclassunc)
77  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
78  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
79  real,parameter :: weightair=28.97
80  logical :: sp_zer
81  character :: adate*8,atime*6
82  character(len=3) :: anspec
83  logical :: lexist
84  integer :: mind
85! mind        eso:added to ensure identical results between 2&3-fields versions
86
87
88  ! Determine current calendar date, needed for the file name
89  !**********************************************************
90
91  jul=bdate+real(itime,kind=dp)/86400._dp
92  call caldate(jul,jjjjmmdd,ihmmss)
93  write(adate,'(i8.8)') jjjjmmdd
94  write(atime,'(i6.6)') ihmmss
95
96
97  ! For forward simulations, output fields have dimension MAXSPEC,
98  ! for backward simulations, output fields have dimension MAXPOINT.
99  ! Thus, make loops either about nspec, or about numpoint
100  !*****************************************************************
101
102
103    if (ldirect.eq.1) then
104       do ks=1,nspec
105         do kp=1,maxpointspec_act
106           tot_mu(ks,kp)=1
107         end do
108       end do
109   else
110      do ks=1,nspec
111             do kp=1,maxpointspec_act
112               tot_mu(ks,kp)=xmass(kp,ks)
113             end do
114      end do
115    endif
116
117
118  !*******************************************************************
119  ! Compute air density: sufficiently accurate to take it
120  ! from coarse grid at some time
121  ! Determine center altitude of output layer, and interpolate density
122  ! data to that altitude
123  !*******************************************************************
124
125  mind=memind(2)
126  do kz=1,numzgrid
127    if (kz.eq.1) then
128      halfheight=outheight(1)/2.
129    else
130      halfheight=(outheight(kz)+outheight(kz-1))/2.
131    endif
132    do kzz=2,nz
133      if ((height(kzz-1).lt.halfheight).and. &
134           (height(kzz).gt.halfheight)) goto 46
135    end do
13646   kzz=max(min(kzz,nz),2)
137    dz1=halfheight-height(kzz-1)
138    dz2=height(kzz)-halfheight
139    dz=dz1+dz2
140    do jy=0,numygridn-1
141      do ix=0,numxgridn-1
142        xl=outlon0n+real(ix)*dxoutn
143        yl=outlat0n+real(jy)*dyoutn
144        xl=(xl-xlon0)/dx
145        yl=(yl-ylat0)/dy
146        iix=max(min(nint(xl),nxmin1),0)
147        jjy=max(min(nint(yl),nymin1),0)
148        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
149        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
150        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
151             rho(iix,jjy,kzz-1,mind)*dz2)/dz
152! RLT
153        densitydrygrid(ix,jy,kz)=(rho_dry(iix,jjy,kzz,mind)*dz1+ &
154             rho_dry(iix,jjy,kzz-1,mind)*dz2)/dz
155      end do
156    end do
157  end do
158
159  do i=1,numreceptor
160    xl=xreceptor(i)
161    yl=yreceptor(i)
162    iix=max(min(nint(xl),nxmin1),0)
163    jjy=max(min(nint(yl),nymin1),0)
164    !densityoutrecept(i)=rho(iix,jjy,1,2)
165    densityoutrecept(i)=rho(iix,jjy,1,mind)
166! RLT
167    densitydryrecept(i)=rho_dry(iix,jjy,1,mind)
168  end do
169
170! RLT
171! conversion factor for output relative to dry air
172  factor_drygrid=densityoutgrid/densitydrygrid
173  factor_dryrecept=densityoutrecept/densitydryrecept
174
175  ! Output is different for forward and backward simulations
176    do kz=1,numzgrid
177      do jy=0,numygridn-1
178        do ix=0,numxgridn-1
179          if (ldirect.eq.1) then
180            factor3d(ix,jy,kz)=1.e12/volumen(ix,jy,kz)/outnum
181          else
182            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
183          endif
184        end do
185      end do
186    end do
187
188  !*********************************************************************
189  ! Determine the standard deviation of the mean concentration or mixing
190  ! ratio (uncertainty of the output) and the dry and wet deposition
191  !*********************************************************************
192
193  do ks=1,nspec
194
195  write(anspec,'(i3.3)') ks
196 
197  if (DRYBKDEP.or.WETBKDEP) then !scavdep output
198      if (DRYBKDEP) &
199      open(unitoutgrid,file=path(2)(1:length(2))//'grid_drydep_nest_'//adate// &
200           atime//'_'//anspec,form='unformatted')
201      if (WETBKDEP) &
202      open(unitoutgrid,file=path(2)(1:length(2))//'grid_wetdep_nest_'//adate// &
203           atime//'_'//anspec,form='unformatted')
204      write(unitoutgrid) itime
205  else
206     if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
207       if (ldirect.eq.1) then
208         open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_nest_' &
209           //adate// &
210           atime//'_'//anspec,form='unformatted')
211       else
212         open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_nest_' &
213           //adate// &
214           atime//'_'//anspec,form='unformatted')
215       endif
216       write(unitoutgrid) itime
217     endif
218  endif
219
220  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
221   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_nest_' &
222        //adate// &
223        atime//'_'//anspec,form='unformatted')
224
225    write(unitoutgridppt) itime
226  endif
227
228  do kp=1,maxpointspec_act
229  do nage=1,nageclass
230
231    do jy=0,numygridn-1
232      do ix=0,numxgridn-1
233
234  ! WET DEPOSITION
235        if ((WETDEP).and.(ldirect.gt.0)) then
236            do l=1,nclassunc
237              auxgrid(l)=wetgriduncn(ix,jy,ks,kp,l,nage)
238            end do
239            call mean(auxgrid,wetgrid(ix,jy), &
240                 wetgridsigma(ix,jy),nclassunc)
241  ! Multiply by number of classes to get total concentration
242            wetgrid(ix,jy)=wetgrid(ix,jy) &
243                 *nclassunc
244  ! Calculate standard deviation of the mean
245            wetgridsigma(ix,jy)= &
246                 wetgridsigma(ix,jy)* &
247                 sqrt(real(nclassunc))
248        endif
249
250  ! DRY DEPOSITION
251        if ((DRYDEP).and.(ldirect.gt.0)) then
252            do l=1,nclassunc
253              auxgrid(l)=drygriduncn(ix,jy,ks,kp,l,nage)
254            end do
255            call mean(auxgrid,drygrid(ix,jy), &
256                 drygridsigma(ix,jy),nclassunc)
257  ! Multiply by number of classes to get total concentration
258            drygrid(ix,jy)=drygrid(ix,jy)* &
259                 nclassunc
260  ! Calculate standard deviation of the mean
261            drygridsigma(ix,jy)= &
262                 drygridsigma(ix,jy)* &
263                 sqrt(real(nclassunc))
264        endif
265
266  ! CONCENTRATION OR MIXING RATIO
267        do kz=1,numzgrid
268            do l=1,nclassunc
269              auxgrid(l)=griduncn(ix,jy,kz,ks,kp,l,nage)
270            end do
271            call mean(auxgrid,grid(ix,jy,kz), &
272                 gridsigma(ix,jy,kz),nclassunc)
273  ! Multiply by number of classes to get total concentration
274            grid(ix,jy,kz)= &
275                 grid(ix,jy,kz)*nclassunc
276  ! Calculate standard deviation of the mean
277            gridsigma(ix,jy,kz)= &
278                 gridsigma(ix,jy,kz)* &
279                 sqrt(real(nclassunc))
280        end do
281      end do
282    end do
283
284
285  !*******************************************************************
286  ! Generate output: may be in concentration (ng/m3) or in mixing
287  ! ratio (ppt) or both
288  ! Output the position and the values alternated multiplied by
289  ! 1 or -1, first line is number of values, number of positions
290  ! For backward simulations, the unit is seconds, stored in grid_time
291  !*******************************************************************
292
293  ! Concentration output
294  !*********************
295  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
296
297  ! Wet deposition
298         sp_count_i=0
299         sp_count_r=0
300         sp_fact=-1.
301         sp_zer=.true.
302         if ((ldirect.eq.1).and.(WETDEP)) then
303         do jy=0,numygridn-1
304            do ix=0,numxgridn-1
305  !oncentraion greater zero
306              if (wetgrid(ix,jy).gt.smallnum) then
307                 if (sp_zer.eqv..true.) then ! first non zero value
308                    sp_count_i=sp_count_i+1
309                    sparse_dump_i(sp_count_i)=ix+jy*numxgridn
310                    sp_zer=.false.
311                    sp_fact=sp_fact*(-1.)
312                 endif
313                 sp_count_r=sp_count_r+1
314                 sparse_dump_r(sp_count_r)= &
315                      sp_fact*1.e12*wetgrid(ix,jy)/arean(ix,jy)
316  !                sparse_dump_u(sp_count_r)=
317  !+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
318              else ! concentration is zero
319                  sp_zer=.true.
320              endif
321            end do
322         end do
323         else
324            sp_count_i=0
325            sp_count_r=0
326         endif
327         write(unitoutgrid) sp_count_i
328         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
329         write(unitoutgrid) sp_count_r
330         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
331  !       write(unitoutgrid) sp_count_u
332  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
333
334  ! Dry deposition
335         sp_count_i=0
336         sp_count_r=0
337         sp_fact=-1.
338         sp_zer=.true.
339         if ((ldirect.eq.1).and.(DRYDEP)) then
340          do jy=0,numygridn-1
341            do ix=0,numxgridn-1
342              if (drygrid(ix,jy).gt.smallnum) then
343                 if (sp_zer.eqv..true.) then ! first non zero value
344                    sp_count_i=sp_count_i+1
345                    sparse_dump_i(sp_count_i)=ix+jy*numxgridn
346                    sp_zer=.false.
347                    sp_fact=sp_fact*(-1.)
348                 endif
349                 sp_count_r=sp_count_r+1
350                 sparse_dump_r(sp_count_r)= &
351                      sp_fact* &
352                      1.e12*drygrid(ix,jy)/arean(ix,jy)
353  !                sparse_dump_u(sp_count_r)=
354  !+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
355              else ! concentration is zero
356                  sp_zer=.true.
357              endif
358            end do
359          end do
360         else
361            sp_count_i=0
362            sp_count_r=0
363         endif
364         write(unitoutgrid) sp_count_i
365         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
366         write(unitoutgrid) sp_count_r
367         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
368  !       write(*,*) sp_count_u
369  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
370
371
372
373  ! Concentrations
374         sp_count_i=0
375         sp_count_r=0
376         sp_fact=-1.
377         sp_zer=.true.
378          do kz=1,numzgrid
379            do jy=0,numygridn-1
380              do ix=0,numxgridn-1
381                if (grid(ix,jy,kz).gt.smallnum) then
382                  if (sp_zer.eqv..true.) then ! first non zero value
383                    sp_count_i=sp_count_i+1
384                    sparse_dump_i(sp_count_i)= &
385                         ix+jy*numxgridn+kz*numxgridn*numygridn
386                    sp_zer=.false.
387                    sp_fact=sp_fact*(-1.)
388                   endif
389                   sp_count_r=sp_count_r+1
390                   sparse_dump_r(sp_count_r)= &
391                        sp_fact* &
392                        grid(ix,jy,kz)* &
393                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
394  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
395  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
396  !                sparse_dump_u(sp_count_r)=
397  !+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
398  !+               factor(ix,jy,kz)/tot_mu(ks,kp)
399              else ! concentration is zero
400                  sp_zer=.true.
401              endif
402              end do
403            end do
404          end do
405         write(unitoutgrid) sp_count_i
406         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
407         write(unitoutgrid) sp_count_r
408         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
409  !       write(unitoutgrid) sp_count_u
410  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
411
412
413
414    endif !  concentration output
415
416  ! Mixing ratio output
417  !********************
418
419  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
420
421  ! Wet deposition
422         sp_count_i=0
423         sp_count_r=0
424         sp_fact=-1.
425         sp_zer=.true.
426         if ((ldirect.eq.1).and.(WETDEP)) then
427          do jy=0,numygridn-1
428            do ix=0,numxgridn-1
429                if (wetgrid(ix,jy).gt.smallnum) then
430                  if (sp_zer.eqv..true.) then ! first non zero value
431                    sp_count_i=sp_count_i+1
432                    sparse_dump_i(sp_count_i)= &
433                         ix+jy*numxgridn
434                    sp_zer=.false.
435                    sp_fact=sp_fact*(-1.)
436                 endif
437                 sp_count_r=sp_count_r+1
438                 sparse_dump_r(sp_count_r)= &
439                      sp_fact* &
440                      1.e12*wetgrid(ix,jy)/arean(ix,jy)
441  !                sparse_dump_u(sp_count_r)=
442  !    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
443              else ! concentration is zero
444                  sp_zer=.true.
445              endif
446            end do
447          end do
448         else
449           sp_count_i=0
450           sp_count_r=0
451         endif
452         write(unitoutgridppt) sp_count_i
453         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
454         write(unitoutgridppt) sp_count_r
455         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
456  !       write(unitoutgridppt) sp_count_u
457  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
458
459
460  ! Dry deposition
461         sp_count_i=0
462         sp_count_r=0
463         sp_fact=-1.
464         sp_zer=.true.
465         if ((ldirect.eq.1).and.(DRYDEP)) then
466          do jy=0,numygridn-1
467            do ix=0,numxgridn-1
468                if (drygrid(ix,jy).gt.smallnum) then
469                  if (sp_zer.eqv..true.) then ! first non zero value
470                    sp_count_i=sp_count_i+1
471                    sparse_dump_i(sp_count_i)= &
472                         ix+jy*numxgridn
473                    sp_zer=.false.
474                    sp_fact=sp_fact*(-1)
475                 endif
476                 sp_count_r=sp_count_r+1
477                 sparse_dump_r(sp_count_r)= &
478                      sp_fact* &
479                      1.e12*drygrid(ix,jy)/arean(ix,jy)
480  !                sparse_dump_u(sp_count_r)=
481  !    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
482              else ! concentration is zero
483                  sp_zer=.true.
484              endif
485            end do
486          end do
487         else
488           sp_count_i=0
489           sp_count_r=0
490         endif
491         write(unitoutgridppt) sp_count_i
492         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
493         write(unitoutgridppt) sp_count_r
494         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
495  !       write(unitoutgridppt) sp_count_u
496  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
497
498
499  ! Mixing ratios
500         sp_count_i=0
501         sp_count_r=0
502         sp_fact=-1.
503         sp_zer=.true.
504          do kz=1,numzgrid
505            do jy=0,numygridn-1
506              do ix=0,numxgridn-1
507                if (grid(ix,jy,kz).gt.smallnum) then
508                  if (sp_zer.eqv..true.) then ! first non zero value
509                    sp_count_i=sp_count_i+1
510                    sparse_dump_i(sp_count_i)= &
511                         ix+jy*numxgridn+kz*numxgridn*numygridn
512                    sp_zer=.false.
513                    sp_fact=sp_fact*(-1.)
514                 endif
515                 sp_count_r=sp_count_r+1
516                 sparse_dump_r(sp_count_r)= &
517                      sp_fact* &
518                      1.e12*grid(ix,jy,kz) &
519                      /volumen(ix,jy,kz)/outnum* &
520                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
521  !                sparse_dump_u(sp_count_r)=
522  !+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
523  !+              outnum*weightair/weightmolar(ks)/
524  !+              densityoutgrid(ix,jy,kz)
525              else ! concentration is zero
526                  sp_zer=.true.
527              endif
528              end do
529            end do
530          end do
531         write(unitoutgridppt) sp_count_i
532         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
533         write(unitoutgridppt) sp_count_r
534         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
535  !       write(unitoutgridppt) sp_count_u
536  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
537
538      endif ! output for ppt
539
540  end do
541  end do
542
543    close(unitoutgridppt)
544    close(unitoutgrid)
545
546  end do
547
548! RLT Aug 2017
549! Write out conversion factor for dry air
550  inquire(file=path(2)(1:length(2))//'factor_drygrid_nest',exist=lexist)
551  if (lexist) then
552    ! open and append
553    open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid_nest',form='unformatted',&
554            status='old',action='write',access='append')
555  else
556    ! create new
557    open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid_nest',form='unformatted',&
558            status='new',action='write')
559  endif
560  sp_count_i=0
561  sp_count_r=0
562  sp_fact=-1.
563  sp_zer=.true.
564  do kz=1,numzgrid
565    do jy=0,numygridn-1
566      do ix=0,numxgridn-1
567        if (factor_drygrid(ix,jy,kz).gt.(1.+smallnum).or.factor_drygrid(ix,jy,kz).lt.(1.-smallnum)) then
568          if (sp_zer.eqv..true.) then ! first value not equal to one
569            sp_count_i=sp_count_i+1
570            sparse_dump_i(sp_count_i)= &
571                  ix+jy*numxgridn+kz*numxgridn*numygridn
572            sp_zer=.false.
573            sp_fact=sp_fact*(-1.)
574          endif
575          sp_count_r=sp_count_r+1
576          sparse_dump_r(sp_count_r)= &
577               sp_fact*factor_drygrid(ix,jy,kz)
578        else ! factor is one
579          sp_zer=.true.
580        endif
581      end do
582    end do
583  end do
584  write(unitoutfactor) sp_count_i
585  write(unitoutfactor) (sparse_dump_i(i),i=1,sp_count_i)
586  write(unitoutfactor) sp_count_r
587  write(unitoutfactor) (sparse_dump_r(i),i=1,sp_count_r)
588  close(unitoutfactor)
589
590
591  ! Reinitialization of grid
592  !*************************
593
594  ! do ks=1,nspec
595  ! do kp=1,maxpointspec_act
596  !   do i=1,numreceptor
597  !     creceptor(i,ks)=0.
598  !   end do
599  !   do jy=0,numygridn-1
600  !     do ix=0,numxgridn-1
601  !       do l=1,nclassunc
602  !         do nage=1,nageclass
603  !           do kz=1,numzgrid
604  !             griduncn(ix,jy,kz,ks,kp,l,nage)=0.
605  !           end do
606  !         end do
607  !       end do
608  !     end do
609  !   end do
610  ! end do
611  ! end do
612  creceptor(:,:)=0.
613  griduncn(:,:,:,:,:,:,:)=0.
614
615
616end subroutine concoutput_nest
617
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG