source: flexpart.git/src/outgrid_init_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: 8.0 KB
Line 
1! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt
2! SPDX-License-Identifier: GPL-3.0-or-later
3
4subroutine outgrid_init_nest
5
6!*****************************************************************************
7!                                                                            *
8!  This routine calculates, for each grid cell of the output nest, the       *
9!  volume and the surface area.                                              *
10!                                                                            *
11!     Author: A. Stohl                                                       *
12!                                                                            *
13!    30 August 2004                                                          *
14!                                                                            *
15!*****************************************************************************
16!                                                                            *
17! Variables:                                                                 *
18!                                                                            *
19! arean              surface area of all output nest cells                   *
20! volumen            volumes of all output nest cells                        *
21!                                                                            *
22!*****************************************************************************
23
24  use unc_mod
25  use outg_mod
26  use par_mod
27  use com_mod
28
29  implicit none
30
31  integer :: ix,jy,kz,ks,kp,nage,l,iix,jjy,ixp,jyp,i1,j1,j,ngrid
32  integer :: stat
33  real :: ylat,gridarea,ylatp,ylatm,hzone,cosfactm,cosfactp
34  real :: xlon,xl,yl,ddx,ddy,rddx,rddy,p1,p2,p3,p4,xtn,ytn,oroh
35  real,parameter :: eps=nxmax/3.e5
36
37
38
39! gridunc,griduncn        uncertainty of outputted concentrations
40  allocate(griduncn(0:numxgridn-1,0:numygridn-1,numzgrid,maxspec, &
41       maxpointspec_act,nclassunc,maxageclass),stat=stat)
42  if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
43
44  if (ldirect.gt.0) then
45    allocate(wetgriduncn(0:numxgridn-1,0:numygridn-1,maxspec, &
46         maxpointspec_act,nclassunc,maxageclass),stat=stat)
47    if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
48    allocate(drygriduncn(0:numxgridn-1,0:numygridn-1,maxspec, &
49         maxpointspec_act,nclassunc,maxageclass),stat=stat)
50    if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
51  endif
52
53#ifdef USE_MPIINPLACE
54#else
55! Extra field for totals at MPI root process
56  if (lroot.and.mpi_mode.gt.0) then
57! If MPI_IN_PLACE option is not used in mpi_mod.f90::mpif_tm_reduce_grid_nest(),
58! then an aux array is needed for parallel grid reduction
59    allocate(griduncn0(0:numxgridn-1,0:numygridn-1,numzgrid,maxspec, &
60         maxpointspec_act,nclassunc,maxageclass),stat=stat)
61    if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
62! allocate a dummy to avoid compilator complaints
63  else if (.not.lroot.and.mpi_mode.gt.0) then
64    allocate(griduncn0(1,1,1,1,1,1,1),stat=stat)
65  end if
66#endif
67  if (ldirect.gt.0) then
68    if (lroot.and.mpi_mode.gt.0) then
69      allocate(wetgriduncn0(0:numxgridn-1,0:numygridn-1,maxspec, &
70           maxpointspec_act,nclassunc,maxageclass),stat=stat)
71      if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
72      allocate(drygriduncn0(0:numxgridn-1,0:numygridn-1,maxspec, &
73           maxpointspec_act,nclassunc,maxageclass),stat=stat)
74      if (stat.ne.0) write(*,*)'ERROR:could not allocate nested gridunc'
75!  endif
76! allocate a dummy to avoid compilator complaints
77    else if (.not.lroot.and.mpi_mode.gt.0) then
78      allocate(wetgriduncn0(1,1,1,1,1,1),stat=stat)
79      allocate(drygriduncn0(1,1,1,1,1,1),stat=stat)
80    end if
81  end if
82
83! Compute surface area and volume of each grid cell: area, volume;
84! and the areas of the northward and eastward facing walls: areaeast, areanorth
85!***********************************************************************
86
87  do jy=0,numygridn-1
88    ylat=outlat0n+(real(jy)+0.5)*dyoutn
89    ylatp=ylat+0.5*dyoutn
90    ylatm=ylat-0.5*dyoutn
91    if ((ylatm.lt.0).and.(ylatp.gt.0.)) then
92      hzone=dyoutn*r_earth*pi180
93    else
94
95! Calculate area of grid cell with formula M=2*pi*R*h*dx/360,
96! see Netz, Formeln der Mathematik, 5. Auflage (1983), p.90
97!************************************************************
98
99      cosfactp=cos(ylatp*pi180)
100      cosfactm=cos(ylatm*pi180)
101      if (cosfactp.lt.cosfactm) then
102        hzone=sqrt(1-cosfactp**2)- &
103             sqrt(1-cosfactm**2)
104        hzone=hzone*r_earth
105      else
106        hzone=sqrt(1-cosfactm**2)- &
107             sqrt(1-cosfactp**2)
108        hzone=hzone*r_earth
109      endif
110    endif
111
112
113
114! Surface are of a grid cell at a latitude ylat
115!**********************************************
116
117    gridarea=2.*pi*r_earth*hzone*dxoutn/360.
118
119    do ix=0,numxgridn-1
120      arean(ix,jy)=gridarea
121
122! Volume = area x box height
123!***************************
124
125      volumen(ix,jy,1)=arean(ix,jy)*outheight(1)
126      do kz=2,numzgrid
127        volumen(ix,jy,kz)=arean(ix,jy)*(outheight(kz)-outheight(kz-1))
128      end do
129    end do
130  end do
131
132
133!**************************************************************************
134! Determine average height of model topography in nesteed output grid cells
135!**************************************************************************
136
137! Loop over all output grid cells
138!********************************
139
140  do jjy=0,numygridn-1
141    do iix=0,numxgridn-1
142      oroh=0.
143
144! Take 100 samples of the topography in every grid cell
145!******************************************************
146
147      do j1=1,10
148        ylat=outlat0n+(real(jjy)+real(j1)/10.-0.05)*dyoutn
149        yl=(ylat-ylat0)/dy
150        do i1=1,10
151          xlon=outlon0n+(real(iix)+real(i1)/10.-0.05)*dxoutn
152          xl=(xlon-xlon0)/dx
153
154! Determine the nest we are in
155!*****************************
156
157          ngrid=0
158          do j=numbnests,1,-1
159            if ((xl.gt.xln(j)+eps).and.(xl.lt.xrn(j)-eps).and. &
160                 (yl.gt.yln(j)+eps).and.(yl.lt.yrn(j)-eps)) then
161              ngrid=j
162              goto 43
163            endif
164          end do
16543        continue
166
167! Determine (nested) grid coordinates and auxiliary parameters used for interpolation
168!*****************************************************************************
169
170          if (ngrid.gt.0) then
171            xtn=(xl-xln(ngrid))*xresoln(ngrid)
172            ytn=(yl-yln(ngrid))*yresoln(ngrid)
173            ix=int(xtn)
174            jy=int(ytn)
175            ddy=ytn-real(jy)
176            ddx=xtn-real(ix)
177          else
178            ix=int(xl)
179            jy=int(yl)
180            ddy=yl-real(jy)
181            ddx=xl-real(ix)
182          endif
183          ixp=ix+1
184          jyp=jy+1
185          rddx=1.-ddx
186          rddy=1.-ddy
187          p1=rddx*rddy
188          p2=ddx*rddy
189          p3=rddx*ddy
190          p4=ddx*ddy
191
192          if (ngrid.gt.0) then
193            oroh=oroh+p1*oron(ix ,jy ,ngrid) &
194                 + p2*oron(ixp,jy ,ngrid) &
195                 + p3*oron(ix ,jyp,ngrid) &
196                 + p4*oron(ixp,jyp,ngrid)
197          else
198            oroh=oroh+p1*oro(ix ,jy) &
199                 + p2*oro(ixp,jy) &
200                 + p3*oro(ix ,jyp) &
201                 + p4*oro(ixp,jyp)
202          endif
203        end do
204      end do
205
206! Divide by the number of samples taken
207!**************************************
208
209      orooutn(iix,jjy)=oroh/100.
210    end do
211  end do
212
213
214
215!*******************************
216! Initialization of output grids
217!*******************************
218
219  do kp=1,maxpointspec_act
220    do ks=1,nspec
221      do nage=1,nageclass
222        do jy=0,numygridn-1
223          do ix=0,numxgridn-1
224            do l=1,nclassunc
225! Deposition fields
226              if (ldirect.gt.0) then
227                wetgriduncn(ix,jy,ks,kp,l,nage)=0.
228                drygriduncn(ix,jy,ks,kp,l,nage)=0.
229              endif
230! Concentration fields
231              do kz=1,numzgrid
232                griduncn(ix,jy,kz,ks,kp,l,nage)=0.
233              end do
234            end do
235          end do
236        end do
237      end do
238    end do
239  end do
240
241
242end subroutine outgrid_init_nest
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG