source: flexpart.git/src/outgrid_init.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: 12.1 KB
Line 
1! DJM - 2017-05-09 - added #ifdef USE_MPIINPLACE cpp directive to     *
2! enable allocation of a gridunc0 array if required by MPI code in    *
3! mpi_mod.f90                                                         *
4!                                                                     *
5!**********************************************************************
6
7
8subroutine outgrid_init
9  !
10  !*****************************************************************************
11  !                                                                            *
12  !  This routine initializes the output grids                                 *
13  !                                                                            *
14  !     Author: A. Stohl                                                       *
15  !                                                                            *
16  !     7 August 2002                                                          *
17  !                                                                            *
18  !*****************************************************************************
19  !                                                                            *
20  ! Variables:                                                                 *
21  !                                                                            *
22  ! area               surface area of all output grid cells                   *
23  ! areaeast           eastward facing wall area of all output grid cells      *
24  ! areanorth          northward facing wall area of all output grid cells     *
25  ! volume             volumes of all output grid cells                        *
26  !                                                                            *
27  !*****************************************************************************
28
29  use flux_mod
30  use oh_mod
31  use unc_mod
32  use outg_mod
33  use par_mod
34  use com_mod
35
36  implicit none
37
38  integer :: ix,jy,kz,i,nage,l,iix,jjy,ixp,jyp,i1,j1,j,ngrid
39  integer :: ks,kp,stat
40  real :: ylat,gridarea,ylatp,ylatm,hzone,cosfactm,cosfactp
41  real :: xlon,xl,yl,ddx,ddy,rddx,rddy,p1,p2,p3,p4,xtn,ytn,oroh
42  real,parameter :: eps=nxmax/3.e5
43
44
45  ! Compute surface area and volume of each grid cell: area, volume;
46  ! and the areas of the northward and eastward facing walls: areaeast, areanorth
47  !***********************************************************************
48  do jy=0,numygrid-1
49    ylat=outlat0+(real(jy)+0.5)*dyout
50    ylatp=ylat+0.5*dyout
51    ylatm=ylat-0.5*dyout
52    if ((ylatm.lt.0).and.(ylatp.gt.0.)) then
53      hzone=dyout*r_earth*pi180
54    else
55
56  ! Calculate area of grid cell with formula M=2*pi*R*h*dx/360,
57  ! see Netz, Formeln der Mathematik, 5. Auflage (1983), p.90
58  !************************************************************
59
60      cosfactp=cos(ylatp*pi180)
61      cosfactm=cos(ylatm*pi180)
62      if (cosfactp.lt.cosfactm) then
63        hzone=sqrt(1-cosfactp**2)- &
64             sqrt(1-cosfactm**2)
65        hzone=hzone*r_earth
66      else
67        hzone=sqrt(1-cosfactm**2)- &
68             sqrt(1-cosfactp**2)
69        hzone=hzone*r_earth
70      endif
71    endif
72
73  ! Surface are of a grid cell at a latitude ylat
74  !**********************************************
75
76    gridarea=2.*pi*r_earth*hzone*dxout/360.
77
78    do ix=0,numxgrid-1
79      area(ix,jy)=gridarea
80
81  ! Volume = area x box height
82  !***************************
83
84      volume(ix,jy,1)=area(ix,jy)*outheight(1)
85      areaeast(ix,jy,1)=dyout*r_earth*pi180*outheight(1)
86      areanorth(ix,jy,1)=cos(ylat*pi180)*dxout*r_earth*pi180* &
87           outheight(1)
88      do kz=2,numzgrid
89        areaeast(ix,jy,kz)=dyout*r_earth*pi180* &
90             (outheight(kz)-outheight(kz-1))
91        areanorth(ix,jy,kz)=cos(ylat*pi180)*dxout*r_earth*pi180* &
92             (outheight(kz)-outheight(kz-1))
93        volume(ix,jy,kz)=area(ix,jy)*(outheight(kz)-outheight(kz-1))
94      end do
95    end do
96  end do
97
98
99
100
101  !******************************************************************
102  ! Determine average height of model topography in output grid cells
103  !******************************************************************
104
105  ! Loop over all output grid cells
106  !********************************
107
108  do jjy=0,numygrid-1
109    do iix=0,numxgrid-1
110      oroh=0.
111
112  ! Take 100 samples of the topography in every grid cell
113  !******************************************************
114
115      do j1=1,10
116        ylat=outlat0+(real(jjy)+real(j1)/10.-0.05)*dyout
117        yl=(ylat-ylat0)/dy
118        do i1=1,10
119          xlon=outlon0+(real(iix)+real(i1)/10.-0.05)*dxout
120          xl=(xlon-xlon0)/dx
121
122  ! Determine the nest we are in
123  !*****************************
124
125          ngrid=0
126          do j=numbnests,1,-1
127            if ((xl.gt.xln(j)+eps).and.(xl.lt.xrn(j)-eps).and. &
128                 (yl.gt.yln(j)+eps).and.(yl.lt.yrn(j)-eps)) then
129              ngrid=j
130              goto 43
131            endif
132          end do
13343        continue
134
135  ! Determine (nested) grid coordinates and auxiliary parameters used for interpolation
136  !*****************************************************************************
137
138          if (ngrid.gt.0) then
139            xtn=(xl-xln(ngrid))*xresoln(ngrid)
140            ytn=(yl-yln(ngrid))*yresoln(ngrid)
141            ix=int(xtn)
142            jy=int(ytn)
143            ddy=ytn-real(jy)
144            ddx=xtn-real(ix)
145          else
146            ix=int(xl)
147            jy=int(yl)
148            ddy=yl-real(jy)
149            ddx=xl-real(ix)
150          endif
151          ixp=ix+1
152          jyp=jy+1
153          rddx=1.-ddx
154          rddy=1.-ddy
155          p1=rddx*rddy
156          p2=ddx*rddy
157          p3=rddx*ddy
158          p4=ddx*ddy
159
160          if (ngrid.gt.0) then
161            oroh=oroh+p1*oron(ix ,jy ,ngrid) &
162                 + p2*oron(ixp,jy ,ngrid) &
163                 + p3*oron(ix ,jyp,ngrid) &
164                 + p4*oron(ixp,jyp,ngrid)
165          else
166            oroh=oroh+p1*oro(ix ,jy) &
167                 + p2*oro(ixp,jy) &
168                 + p3*oro(ix ,jyp) &
169                 + p4*oro(ixp,jyp)
170          endif
171        end do
172      end do
173
174  ! Divide by the number of samples taken
175  !**************************************
176
177      oroout(iix,jjy)=oroh/100.
178    end do
179  end do
180
181  ! if necessary allocate flux fields
182  if (iflux.eq.1) then
183    allocate(flux(6,0:numxgrid-1,0:numygrid-1,numzgrid, &
184         1:nspec,1:maxpointspec_act,1:nageclass),stat=stat)
185    if (stat.ne.0) write(*,*)'ERROR: could not allocate flux array '
186  endif
187
188  ! gridunc,griduncn        uncertainty of outputted concentrations
189  allocate(gridunc(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, &
190       maxpointspec_act,nclassunc,maxageclass),stat=stat)
191  if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
192  if (ldirect.gt.0) then
193    allocate(wetgridunc(0:numxgrid-1,0:numygrid-1,maxspec, &
194         maxpointspec_act,nclassunc,maxageclass),stat=stat)
195    if (stat.ne.0) write(*,*)'ERROR: could not allocate wetgridunc'
196    allocate(drygridunc(0:numxgrid-1,0:numygrid-1,maxspec, &
197         maxpointspec_act,nclassunc,maxageclass),stat=stat)
198    if (stat.ne.0) write(*,*)'ERROR: could not allocate drygridunc'
199  endif
200
201#ifdef USE_MPIINPLACE
202#else
203! Extra field for totals at MPI root process
204  if (lroot.and.mpi_mode.gt.0) then
205! If MPI_IN_PLACE option is not used in mpi_mod.f90::mpif_tm_reduce_grid(),
206! then an aux array is needed for parallel grid reduction
207    allocate(gridunc0(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, &
208         maxpointspec_act,nclassunc,maxageclass),stat=stat)
209    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc0'
210  else if (.not.lroot.and.mpi_mode.gt.0) then
211    allocate(gridunc0(1,1,1,1,1,1,1),stat=stat)
212    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc0'
213  end if
214#endif
215  if (ldirect.gt.0) then
216    if (lroot.and.mpi_mode.gt.0) then
217      allocate(wetgridunc0(0:numxgrid-1,0:numygrid-1,maxspec, &
218           maxpointspec_act,nclassunc,maxageclass),stat=stat)
219      if (stat.ne.0) write(*,*)'ERROR: could not allocate wetgridunc0'
220      allocate(drygridunc0(0:numxgrid-1,0:numygrid-1,maxspec, &
221           maxpointspec_act,nclassunc,maxageclass),stat=stat)
222      if (stat.ne.0) write(*,*)'ERROR: could not allocate drygridunc0'
223
224! allocate a dummy to avoid compilator complaints
225    else if (.not.lroot.and.mpi_mode.gt.0) then
226      allocate(wetgridunc0(1,1,1,1,1,1),stat=stat)
227      allocate(drygridunc0(1,1,1,1,1,1),stat=stat)
228    end if
229  end if
230
231  !write (*,*) 'Dimensions for fields', numxgrid,numygrid, &
232  !     maxspec,maxpointspec_act,nclassunc,maxageclass
233
234  if (lroot) then
235    write (*,*) 'Allocating fields for global output (x,y): ', numxgrid,numygrid
236    write (*,*) 'Allocating fields for nested output (x,y): ', numxgridn,numygridn
237  end if
238
239  ! allocate fields for concoutput with maximum dimension of outgrid
240  ! and outgrid_nest
241
242  allocate(gridsigma(0:max(numxgrid,numxgridn)-1, &
243       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
244    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
245  allocate(grid(0:max(numxgrid,numxgridn)-1, &
246       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
247    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
248  allocate(densityoutgrid(0:max(numxgrid,numxgridn)-1, &
249       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
250    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
251! RLT
252  allocate(densitydrygrid(0:max(numxgrid,numxgridn)-1, &
253       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
254    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
255  allocate(factor_drygrid(0:max(numxgrid,numxgridn)-1, &
256       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
257    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
258
259  allocate(factor3d(0:max(numxgrid,numxgridn)-1, &
260       0:max(numygrid,numygridn)-1,numzgrid),stat=stat)
261    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
262  allocate(sparse_dump_r(max(numxgrid,numxgridn)* &
263       max(numygrid,numygridn)*numzgrid),stat=stat)
264    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
265
266   allocate(sparse_dump_u(max(numxgrid,numxgridn)* &
267       max(numygrid,numygridn)*numzgrid),stat=stat)
268        if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
269
270  allocate(sparse_dump_i(max(numxgrid,numxgridn)* &
271       max(numygrid,numygridn)*numzgrid),stat=stat)
272    if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
273
274  ! deposition fields are only allocated for forward runs
275  if (ldirect.gt.0) then
276     allocate(wetgridsigma(0:max(numxgrid,numxgridn)-1, &
277          0:max(numygrid,numygridn)-1),stat=stat)
278     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
279     allocate(drygridsigma(0:max(numxgrid,numxgridn)-1, &
280          0:max(numygrid,numygridn)-1),stat=stat)
281     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
282     allocate(wetgrid(0:max(numxgrid,numxgridn)-1, &
283          0:max(numygrid,numygridn)-1),stat=stat)
284     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
285     allocate(drygrid(0:max(numxgrid,numxgridn)-1, &
286          0:max(numygrid,numygridn)-1),stat=stat)
287     if (stat.ne.0) write(*,*)'ERROR: could not allocate gridunc'
288  endif
289
290  ! Initial condition field
291
292  if (linit_cond.gt.0) then
293    allocate(init_cond(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, &
294         maxpointspec_act),stat=stat)
295    if (stat.ne.0) write(*,*)'ERROR: could not allocate init_cond'
296  endif
297
298  !************************
299  ! Initialize output grids
300  !************************
301
302  do ks=1,nspec
303  do kp=1,maxpointspec_act
304    do i=1,numreceptor
305  ! Receptor points
306      creceptor(i,ks)=0.
307    end do
308    do nage=1,nageclass
309      do jy=0,numygrid-1
310        do ix=0,numxgrid-1
311          do l=1,nclassunc
312  ! Deposition fields
313            if (ldirect.gt.0) then
314            wetgridunc(ix,jy,ks,kp,l,nage)=0.
315            drygridunc(ix,jy,ks,kp,l,nage)=0.
316            endif
317            do kz=1,numzgrid
318              if (iflux.eq.1) then
319  ! Flux fields
320                 do i=1,5
321                   flux(i,ix,jy,kz,ks,kp,nage)=0.
322                 end do
323              endif
324  ! Initial condition field
325              if ((l.eq.1).and.(nage.eq.1).and.(linit_cond.gt.0)) &
326                   init_cond(ix,jy,kz,ks,kp)=0.
327  ! Concentration fields
328              gridunc(ix,jy,kz,ks,kp,l,nage)=0.
329            end do
330          end do
331        end do
332      end do
333    end do
334  end do
335  end do
336
337
338
339end subroutine outgrid_init
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG