!********************************************************************** ! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010 * ! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa, * ! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann * ! * ! This file is part of FLEXPART. * ! * ! FLEXPART is free software: you can redistribute it and/or modify * ! it under the terms of the GNU General Public License as published by* ! the Free Software Foundation, either version 3 of the License, or * ! (at your option) any later version. * ! * ! FLEXPART is distributed in the hope that it will be useful, * ! but WITHOUT ANY WARRANTY; without even the implied warranty of * ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * ! GNU General Public License for more details. * ! * ! You should have received a copy of the GNU General Public License * ! along with FLEXPART. If not, see . * !********************************************************************** subroutine gridcheck_nests !***************************************************************************** ! * ! This routine checks the grid specification for the nested model * ! domains. It is similar to subroutine gridcheck, which checks the * ! mother domain. * ! * ! Authors: A. Stohl, G. Wotawa * ! * ! 8 February 1999 * ! * !***************************************************************************** use par_mod use com_mod implicit none integer :: i,j,k,l,ifn,ifield,iumax,iwmax,numskip,nlev_ecn integer :: nuvzn,nwzn real :: akmn(nwzmax),bkmn(nwzmax),akzn(nuvzmax),bkzn(nuvzmax) real :: xaux1,xaux2,yaux1,yaux2 ! VARIABLES AND ARRAYS NEEDED FOR GRIB DECODING ! dimension of isec2 at least (22+n), where n is the number of parallels or ! meridians in a quasi-regular (reduced) Gaussian or lat/long grid ! dimension of zsec2 at least (10+nn), where nn is the number of vertical ! coordinate parameters integer :: isec0(2),isec1(56),isec2(22+nxmaxn+nymaxn),isec3(2) integer :: isec4(64),inbuff(jpack),ilen,ierr,iword,lunit !integer iswap real :: zsec2(60+2*nuvzmax),zsec3(2),zsec4(jpunp) character(len=1) :: yoper = 'D' xresoln(0)=1. ! resolution enhancement for mother grid yresoln(0)=1. ! resolution enhancement for mother grid ! Loop about all nesting levels !****************************** do l=1,numbnests iumax=0 iwmax=0 if(ideltas.gt.0) then ifn=1 else ifn=numbwf endif ! ! OPENING OF DATA FILE (GRIB CODE) ! 5 call pbopen(lunit,path(numpath+2*(l-1)+1) & (1:length(numpath+2*(l-1)+1))//wfnamen(l,ifn),'r',ierr) if(ierr.lt.0) goto 999 ifield=0 10 ifield=ifield+1 ! ! GET NEXT FIELDS ! call pbgrib(lunit,inbuff,jpack,ilen,ierr) if(ierr.eq.-1) goto 30 ! EOF DETECTED if(ierr.lt.-1) goto 999 ! ERROR DETECTED ierr=1 ! Check whether we are on a little endian or on a big endian computer !******************************************************************** ! if (inbuff(1).eq.1112101447) then ! little endian, swap bytes ! iswap=1+ilen/4 ! call swap32(inbuff,iswap) ! else if (inbuff(1).ne.1196575042) then ! big endian ! stop 'subroutine gridcheck: corrupt GRIB data' ! endif call gribex(isec0,isec1,isec2,zsec2,isec3,zsec3,isec4, & zsec4,jpunp,inbuff,jpack,iword,yoper,ierr) if (ierr.ne.0) goto 999 ! ERROR DETECTED if(ifield.eq.1) then nxn(l)=isec2(2) nyn(l)=isec2(3) xaux1=real(isec2(5))/1000. xaux2=real(isec2(8))/1000. if(xaux1.gt.180) xaux1=xaux1-360.0 if(xaux2.gt.180) xaux2=xaux2-360.0 if(xaux1.lt.-180) xaux1=xaux1+360.0 if(xaux2.lt.-180) xaux2=xaux2+360.0 if (xaux2.lt.xaux1) xaux2=xaux2+360. yaux1=real(isec2(7))/1000. yaux2=real(isec2(4))/1000. xlon0n(l)=xaux1 ylat0n(l)=yaux1 dxn(l)=(xaux2-xaux1)/real(nxn(l)-1) dyn(l)=(yaux2-yaux1)/real(nyn(l)-1) nlev_ecn=isec2(12)/2-1 endif k=isec1(8) if(isec1(6).eq.131) iumax=max(iumax,nlev_ec-k+1) if(isec1(6).eq.135) iwmax=max(iwmax,nlev_ec-k+1) if(isec1(6).eq.129) then do j=0,nyn(l)-1 do i=0,nxn(l)-1 oron(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga end do end do endif if(isec1(6).eq.172) then do j=0,nyn(l)-1 do i=0,nxn(l)-1 lsmn(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga end do end do endif if(isec1(6).eq.160) then do j=0,nyn(l)-1 do i=0,nxn(l)-1 excessoron(i,j,l)=zsec4(nxn(l)*(nyn(l)-j-1)+i+1)/ga end do end do endif goto 10 !! READ NEXT LEVEL OR PARAMETER ! ! CLOSING OF INPUT DATA FILE ! 30 call pbclose(lunit,ierr) !! FINISHED READING / CLOSING GRIB FILE nuvzn=iumax nwzn=iwmax if(nuvzn.eq.nlev_ec) nwzn=nlev_ecn+1 if (nxn(l).gt.nxmaxn) then write(*,*) 'FLEXPART error: Too many grid points in x direction.' write(*,*) 'Reduce resolution of wind fields (file GRIDSPEC)' write(*,*) 'for nesting level ',l write(*,*) 'Or change parameter settings in file par_mod.' write(*,*) nxn(l),nxmaxn stop endif if (nyn(l).gt.nymaxn) then write(*,*) 'FLEXPART error: Too many grid points in y direction.' write(*,*) 'Reduce resolution of wind fields (file GRIDSPEC)' write(*,*) 'for nesting level ',l write(*,*) 'Or change parameter settings in file par_mod.' write(*,*) nyn(l),nymaxn stop endif if ((nuvzn.gt.nuvzmax).or.(nwzn.gt.nwzmax)) then write(*,*) 'FLEXPART error: Nested wind fields have too many'// & 'vertical levels.' write(*,*) 'Problem was encountered for nesting level ',l stop endif ! Output of grid info !******************** write(*,'(a,i2)') 'Nested domain #: ',l write(*,'(a,f10.2,a1,f10.2,a,f10.2)') ' Longitude range: ', & xlon0n(l),' to ',xlon0n(l)+(nxn(l)-1)*dxn(l), & ' Grid distance: ',dxn(l) write(*,'(a,f10.2,a1,f10.2,a,f10.2)') ' Latitude range: ', & ylat0n(l),' to ',ylat0n(l)+(nyn(l)-1)*dyn(l), & ' Grid distance: ',dyn(l) write(*,*) ! Determine, how much the resolutions in the nests are enhanced as ! compared to the mother grid !***************************************************************** xresoln(l)=dx/dxn(l) yresoln(l)=dy/dyn(l) ! Determine the mother grid coordinates of the corner points of the ! nested grids ! Convert first to geographical coordinates, then to grid coordinates !******************************************************************** xaux1=xlon0n(l) xaux2=xlon0n(l)+real(nxn(l)-1)*dxn(l) yaux1=ylat0n(l) yaux2=ylat0n(l)+real(nyn(l)-1)*dyn(l) xln(l)=(xaux1-xlon0)/dx xrn(l)=(xaux2-xlon0)/dx yln(l)=(yaux1-ylat0)/dy yrn(l)=(yaux2-ylat0)/dy if ((xln(l).lt.0.).or.(yln(l).lt.0.).or. & (xrn(l).gt.real(nxmin1)).or.(yrn(l).gt.real(nymin1))) then write(*,*) 'Nested domain does not fit into mother domain' write(*,*) 'For global mother domain fields, you can shift' write(*,*) 'shift the mother domain into x-direction' write(*,*) 'by setting nxshift (file par_mod) to a' write(*,*) 'positive value. Execution is terminated.' stop endif ! CALCULATE VERTICAL DISCRETIZATION OF ECMWF MODEL ! PARAMETER akm,bkm DESCRIBE THE HYBRID "ETA" COORDINATE SYSTEM numskip=nlev_ecn-nuvzn ! number of ecmwf model layers not used by FLEXPART do i=1,nwzn j=numskip+i k=nlev_ecn+1+numskip+i akmn(nwzn-i+1)=zsec2(10+j) bkmn(nwzn-i+1)=zsec2(10+k) end do ! ! CALCULATION OF AKZ, BKZ ! AKZ,BKZ: model discretization parameters at the center of each model ! layer ! ! Assign the 10 m winds to an artificial model level with akz=0 and bkz=1.0, ! i.e. ground level !***************************************************************************** akzn(1)=0. bkzn(1)=1.0 do i=1,nuvzn akzn(i+1)=0.5*(akmn(i+1)+akmn(i)) bkzn(i+1)=0.5*(bkmn(i+1)+bkmn(i)) end do nuvzn=nuvzn+1 ! Check, whether the heights of the model levels of the nested ! wind fields are consistent with those of the mother domain. ! If not, terminate model run. !************************************************************* do i=1,nuvz if ((akzn(i).ne.akz(i)).or.(bkzn(i).ne.bkz(i))) then write(*,*) 'FLEXPART error: The wind fields of nesting level',l write(*,*) 'are not consistent with the mother domain:' write(*,*) 'Differences in vertical levels detected.' stop endif end do do i=1,nwz if ((akmn(i).ne.akm(i)).or.(bkmn(i).ne.bkm(i))) then write(*,*) 'FLEXPART error: The wind fields of nesting level',l write(*,*) 'are not consistent with the mother domain:' write(*,*) 'Differences in vertical levels detected.' stop endif end do end do return 999 write(*,*) write(*,*) ' ###########################################'// & '###### ' write(*,*) ' FLEXPART SUBROUTINE GRIDCHECK:' write(*,*) ' CAN NOT OPEN INPUT DATA FILE '//wfnamen(l,ifn) write(*,*) ' FOR NESTING LEVEL ',k write(*,*) ' ###########################################'// & '###### ' stop end subroutine gridcheck_nests