source: LMDZ6/branches/Amaury_dev/tools/make_sso/make_sso_SpherePack.f90 @ 5193

Last change on this file since 5193 was 5101, checked in by abarral, 4 months ago

Handle DEBUG_IO in lmdz_cppkeys_wrapper.F90
Transform some files .F -> .[fF]90
[ne compile pas à cause de writefield_u non défini - en attente de réponse Laurent]

File size: 48.5 KB
Line 
1!-------------------------------------------------------------------------------
2
3PROGRAM make_sso
4
5!-------------------------------------------------------------------------------
6! Purpose: Project ETOPO file (GMT4 axes conventions) on spherical harmonics.
7!-------------------------------------------------------------------------------
8  USE netcdf, ONLY: nf90_noerr,nf90_strerror,nf90_close,nf90_put_var,nf90_enddef,&
9          nf90_put_att,nf90_global,nf90_real,nf90_def_var,nf90_def_dim,nf90_inq_varid,&
10          nf90_nowrite,nf90_inquire_dimension,nf90_inquire_variable,nf90_open
11!  USE sphpack
12  IMPLICIT NONE
13!-------------------------------------------------------------------------------
14! Local variables:
15  REAL,    PARAMETER :: xpi=4.*ATAN(1.),         &    !--- pi
16                        hp=xpi/2., d2r=xpi/180., &    !--- pi, pi/2, pi/180 
17                        REarth=2.E7/xpi,         &    !--- Earth radius
18                        eps=1.E-6                     !--- Tolerance for angles
19  !=== Transition length for spectral filters (to avoid side lobes)
20!  REAL,    PARAMETER :: del=2000./REarth              !--- On unit sphere
21  !=== .TRUE./.FALSE.: Fast/slow routine with heavy O(N^3)/light O(N^2) storage
22  LOGICAL, PARAMETER :: lfast=.TRUE.                  !--- Fast/slow routine
23
24!--- From the input arguments
25  CHARACTER(LEN=256), ALLOCATABLE :: args(:)          !--- All the arguments
26  CHARACTER(LEN=2048):: call_seq                      !--- Calling sequence
27  INTEGER            :: nargs                         !--- Input args number
28  CHARACTER(LEN=256) :: f_in                          !--- Input file name
29  CHARACTER(LEN=256) :: vnam                          !--- Height variable name
30  CHARACTER(LEN=256) :: res_in                        !--- Input resolution
31  CHARACTER(LEN=256) :: fmsk                          !--- Mask file or keyword
32  !=== Filter on output field, of resolution=ffactor*output resolution
33  LOGICAL            :: lfilter=.TRUE.                !--- Filter activation
34  REAL               :: ffactor                       !--- Filter resol. factor
35  REAL               :: tfactor                       !--- Transition factor
36  INTEGER            :: nlon_ou, nlat_ou              !--- Output grid resol.
37
38!--- Quantities related to input field grid
39  CHARACTER(LEN=256) :: lonn, latn                    !--- Axes vars names
40  CHARACTER(LEN=256) :: lonu, latu                    !--- Angular units
41  REAL, ALLOCATABLE  :: lon_in(:), lat_in(:)          !--- Cells centers
42  INTEGER            :: nlon_in,   nlat_in            !--- Input  grid dims
43  INTEGER            :: igrid_in(2)                   !--- Grid identification
44  REAL               :: dlon, dlat                    !--- Domain sizes
45  REAL               :: D1, del1                      !--- Resol. at equator
46  REAL               :: al                            !--- Longitude ratio
47  LOGICAL            :: lcyc                          !--- T: redundant longitude
48  LOGICAL            :: lXdec, lYinc                  !--- T: dec/increasing lon/lat
49  REAL, ALLOCATABLE, DIMENSION(:,:) ::v,h , a , b,mski!--- Heights, coeffs, mask
50  REAL, ALLOCATABLE, DIMENSION(:,:,:) ::hs, as, bs, & !--- Same, filtered
51    dxZ, dyZ, dxZ2, dyZ2, dxyZ, Xk, Xl, Xm, Xp, Xq, dH!--- Grads, correlations
52
53!--- Quantities related to output grid
54  CHARACTER(LEN=256) :: f_ou                          !--- Output file name
55  CHARACTER(LEN=256) :: res_ou                        !--- Output resolution
56  INTEGER(KIND=1), ALLOCATABLE :: mask(:,:)           !--- Boolean out land mask
57  REAL, ALLOCATABLE  :: lon_ou(:), lat_ou(:), tmp1(:) !--- Output grid coords
58  INTEGER            :: igrid_ou(2)                   !--- Grid identification
59  REAL               :: D2, del2                      !--- Resol. at equator
60  REAL, ALLOCATABLE, DIMENSION(:,:) ::   msko, &      !--- Fractional land mask
61                                           h0, &      !--- Large scale orography
62                   Zphi, Zsig, Zthe, Zgam, mu, tmp2   !--- Geopo, SSO parameters
63
64!--- Miscellanous Netcdf-related variables
65  INTEGER :: fID, xID, loID, phiID, sigID, theID, picID, dIDs(4), nn
66  INTEGER :: vID, yID, laID, meaID, gamID, mskID, valID, muID, ix180
67  LOGICAL :: ll
68
69!--- Miscellanous SPHEREPACK-specific quantities
70  INTEGER :: l1, idh, mdab, lsh, la1, la2, lwa, lsvmin, lsave, ierr, intl
71  INTEGER :: l2, jdh, ndab, lvh, lb1, lb2, lwb, lwkmin, lwork, ldwork
72  REAL, ALLOCATABLE :: wsha(:), wsave(:), wvh (:)
73  REAL, ALLOCATABLE :: wshs(:), dwork(:), work(:)
74  INTEGER :: k, n, nt, isym
75
76  CHARACTER(LEN=256) :: sub, fnam, msg, arg
77
78!-------------------------------------------------------------------------------
79  sub='make_sso'
80
81!=== DEAL WITH INPUT ARGUMENTS =================================================
82  call_seq='./'//TRIM(sub)//'.e'
83  nargs=COMMAND_ARGUMENT_COUNT(); ALLOCATE(args(nargs))
84!--- MANUAL DISPLAYED IF NO ARGUMENTS ARE PROVIDED
85  IF(nargs==0) THEN
86    WRITE(*,*)''
87    WRITE(*,*)'General calling sequence:'
88    WRITE(*,*)''
89    WRITE(*,*)TRIM(call_seq)//'  -i <file>  -v <var>  -res <nlon> <nlat>'      &
90      //'  -m <mask>  -f <ffac>  -t <tfac>'
91    WRITE(*,*)''
92    WRITE(*,*)'Where:'
93    WRITE(*,*)''
94    WRITE(*,*)' <file>:   input file name'
95    WRITE(*,*)' <var> :   height variable name'
96    WRITE(*,*)' <nlon>:   number of distinct longitude points of output grid'
97    WRITE(*,*)' <nlat>:   number of distinct latitude  points of output grid'
98    WRITE(*,*)' <mask>:   can be:  1. a mask file (o2a.nc)'
99    WRITE(*,*)'                    2. "noro": computed with grid_noro'
100    WRITE(*,*)'           default: 3. "spec": computed internally'
101    WRITE(*,*)' <ffac>:   filtering  width: ffac*Do         (Do:output resol.)'
102    WRITE(*,*)' <tfac>:   transition width: Di+tfac*(Do-Di) (Di: input resol.)'
103    WRITE(*,*)
104    WRITE(*,*)'Note that latitudes grid contains both poles.'
105    WRITE(*,*)
106    STOP
107  END IF
108!--- ARGUMENTS PARSING
109  DO k=1,nargs
110    CALL GET_COMMAND_ARGUMENT(NUMBER=k,VALUE=args(k))
111    call_seq=TRIM(call_seq)//' '//TRIM(args(k))
112  END DO
113  k=1
114  f_in=''; nlon_ou=0; nlat_ou=0; vnam=''; fmsk=''; ffactor=2.0; tfactor=0.
115  DO
116    k=k+1; ll=.FALSE.; arg=args(k-1)
117    SELECT CASE(arg)
118      CASE('-i');   f_in=args(k); msg='Missing file "'//TRIM(f_in)//'".'
119                    ll=nf90_open(f_in,nf90_nowrite,fID)/=nf90_noerr
120                    IF(.NOT.ll) n=nf90_close(fID)
121      CASE('-res'); nlon_ou=str2int(args(k)); k=k+1
122                    nlat_ou=str2int(args(k))
123      CASE('-v');   vnam=args(k)
124      CASE('-m');   fmsk=args(k)
125      CASE('-f');   ffactor=str2real(args(k))
126      CASE('-t');   tfactor=str2real(args(k))
127      CASE DEFAULT; msg='Unrecognized argument "'//TRIM(arg)//'".'; ll=.TRUE.
128    END SELECT
129    IF(ll) THEN; WRITE(*,*)msg; STOP; END IF
130    k=k+1
131    IF(k>nargs) EXIT
132  END DO
133!--- CHECK OUTPUT GRID IS DEFINED
134  IF(TRIM(fmsk)=='') fmsk='spec'
135  IF(ALL(['noro','spec']/=fmsk)) THEN
136    msg='Missing or wrong "-m" option ; can be "noro", "spec" or a mask file'
137    CALL err(nf90_open(fmsk,nf90_nowrite,fID)/=nf90_noerr,msg)
138    CALL nc(nf90_inq_varid(fID,"MaskOcean",vID),"MaskOcean")    !--- MASK ID
139    CALL nc(nf90_inquire_variable(fID,vID,dimids=dIDs))         !--- DIMS IDS
140    CALL nc(nf90_inquire_dimension(fID,dIDs(1),len=nlon_ou),'x')!--- NB LONG
141    CALL nc(nf90_inquire_dimension(fID,dIDs(2),len=nlat_ou),'y')!--- NB LAT
142    CALL nc(nf90_close(fID))
143  END IF
144  IF(nlon_ou<=0.OR.nlat_ou<=0) THEN
145    WRITE(*,*)'Missing or wrong "-res" arguments for output grid'; STOP
146  END IF
147
148!=== READ THE INPUT FIELD =======================================================
149  CALL nc(nf90_open(f_in,nf90_nowrite,fID))
150  WRITE(*,*)'>> Reading variable "'//TRIM(vnam)//'" from "'//TRIM(f_in)//'"...'
151
152!--- GET HEIGHT VARIABLE AND ITS DIMENSIONS
153  CALL nc(nf90_inq_varid(fID,vnam,vID))
154  CALL nc(nf90_inquire_variable (fID,vID,dimids=dIDs))
155  CALL nc(nf90_inquire_dimension(fID,dIDs(1),len=nlon_in,name=lonn))
156  CALL nc(nf90_inquire_dimension(fID,dIDs(2),len=nlat_in,name=latn))
157  WRITE(*,*)TRIM(vnam)//' is '//TRIM(lonn)//'('//TRIM(int2str(nlon_in))//')*'   &
158                              //TRIM(latn)//'('//TRIM(int2str(nlat_in))//')'
159  nlon_in=nlon_in-1                                   !--- DISTINCT POINTS NUMBER
160  ALLOCATE(lon_in(nlon_in+1),lat_in(nlat_in),h(nlon_in+1,nlat_in))
161  CALL nc(nf90_inq_varid(fID,lonn,loID)        ,lonn)
162  CALL nc(nf90_inq_varid(fID,latn,laID)        ,latn)
163  CALL nc(nf90_get_var  (fID,loID,lon_in)      ,lonn)
164  CALL nc(nf90_get_var  (fID,laID,lat_in)      ,latn)
165  CALL nc(nf90_get_att  (fID,loID,'units',lonu),lonn)
166  CALL nc(nf90_get_att  (fID,laID,'units',latu),latn)
167  CALL nc(nf90_get_var  (fID,vID,h(:,:),[1,1],[nlon_in+1,nlat_in]),vnam)
168  CALL nc(nf90_close(fID))
169
170!--- CHECK WETHER GRID IS CORRECT (GLOBAL DOMAIN, IDENTIFIED UNITS...)
171  dlon=ABS(lon_in(nlon_in+1)-lon_in(1)); dlat=ABS(lat_in(nlat_in)-lat_in(1))
172  CALL err(lonu(1:6)/='degree','Invalid longitudes unit: "'//TRIM(lonu)//'"')
173  CALL err(latu(1:6)/='degree','Invalid latitudes unit: "' //TRIM(latu)//'"')
174  CALL err(ABS(dlon-360.)>eps,'Longitudes domain is not global.')
175  CALL err(ABS(dlat-180.)>eps, 'Latitudes domain is not global.')
176
177!--- ENSURE LONGITUDES ARE INCREASING AND LATITUDES ARE DECREASING
178  lXdec=lon_in(1)>lon_in(nlon_in)
179  lYinc=lat_in(1)<lat_in(nlat_in)
180  IF(lXdec) THEN; lon_in=lon_in(nlon_in+1:1:-1); h(:,:)=h(nlon_in+1:1:-1,:); END IF
181  IF(lYinc) THEN; lat_in=lat_in(nlat_in  :1:-1); h(:,:)=h(:,nlat_in  :1:-1); END IF
182
183!--- ENSURE LONGITUDES ARE BETWEEN -180 and 180deg WITHOUT REDUNDANT END POINT
184  DO ix180=1,nlon_in+1; IF(lon_in(ix180)>=180.) EXIT; END DO
185  IF(ix180>=nlon_in+1) THEN
186    lon_in(:)=[lon_in(ix180:nlon_in)-360.,lon_in(1:ix180)]
187    ALLOCATE(tmp2(nlon_in+1,nlat_in)); nn=nlon_in-ix180+1
188    tmp2(1:nn,:)=h(ix180:nlon_in,:); tmp2(nn+1:nlon_in+1,:)=h(1:ix180,:)
189    DEALLOCATE(tmp2)
190  END IF
191
192!--- REMOVE OCEAN ; FLIP FIELD (SPHEREPACK: LAT FIRST) ; INPUT BOOLEAN MASK
193  WHERE(h(:,:)<0.) h(:,:)=0
194  idh=nlat_in; jdh=nlon_in; CALL flip(h)
195  ALLOCATE(mski(idh,jdh)); mski(:,:)=0.; WHERE(h(:,:)>0.) mski(:,:)=1.
196
197!=== SPHERICAL HARMONICS ANALYSIS ===============================================
198!--- Allocate work vector (min. length depends on used routines)
199  nt=2                                                        ! MAX nt
200  l1=MIN(nlat_in,(nlon_in+2)/2)
201  IF(MODULO(nlon_in,2)==1) l1=MIN(nlat_in,(nlon_in+1)/2)
202  l2=(nlat_in+1)/2
203  IF(lfast) THEN
204    lwork=2*(nlat_in+1)                                       ! shaeci/shseci
205    lwork=MAX(lwork,nlat_in*(nt*nlon_in+MAX(3*l2,nlon_in)))   ! shaec /shaes
206    lwork=MAX(lwork,4*(nlat_in+1))                            ! vhseci
207    lwork=MAX(lwork,nlat_in*(2*nt*nlon_in+MAX(6*l2,nlon_in))+nlat_in*(2*l1*nt+1)) ! gradec
208    lsh=2*nlat_in*l2+3*(  (l1-2  )*(nlat_in+nlat_in-l1-1))/2+nlon_in+15 ! shaec /shsec
209    lvh=4*nlat_in*l2+3*MAX(l1-2,0)*(nlat_in+nlat_in-l1-1)   +nlon_in+15 ! vhseci/gradec
210  ELSE
211    lwork=5*nlat_in*l2+3*((l1-2)*(nlat_in+nlat_in-l1-1))/2    ! shaesi/shsesi
212    lwork=MAX(lwork,(nt+1)*nlat_in*nlon_in)                   ! shaes /shses
213    lwork=MAX(lwork,3*(MAX(l1-2,0)*(nlat_in+nlat_in-l1-1))/2+5*l2*nlat_in)     ! vhsesi
214    lwork=MAX(lwork,nlat_in*((2*nt+1)*nlon_in+2*l1*nt+1))     ! grades
215    lsh=(l1*l2*(nlat_in+nlat_in-l1+1))/2+nlon_in+15           ! shaes /shses
216    lvh= l1*l2*(nlat_in+nlat_in-l1+1)   +nlon_in+15           ! vhsesi/grades
217  END IF
218  ALLOCATE(work(lwork))
219  ALLOCATE(wsha(lsh),wshs(lsh),wvh(lvh))
220
221!--- Allocate spectral coefficients vectors
222  mdab=l1; ndab=nlat_in; IF(.NOT.lfast) mdab=MIN(nlat_in,(nlon_in+2)/2)
223  ALLOCATE(a(mdab,ndab),b(mdab,ndab))
224  isym=0; nt=1                                                !--- Symetries, fields nb
225
226!--- FIELD PROJECTION ON THE SPHERICAL HARMONICS BASE.
227  WRITE(*,*)'>> Projecting fields on the spherical harmonics...'
228  IF(lfast) THEN
229    CALL shaeci(nlat_in,nlon_in,wsha,lsh,work,lwork,ierr)
230    CALL error('shaeci',ierr)
231    CALL shaec(nlat_in,nlon_in,isym,nt,h,idh,jdh,a,b,mdab,ndab,wsha,lsh,work,lwork,ierr)
232    CALL error('shaec',ierr)
233  ELSE
234    CALL shaesi(nlat_in,nlon_in,wsha,lsh,work,lwork,ierr)
235    CALL error('shaesi',ierr)
236    CALL shaes(nlat_in,nlon_in,isym,nt,h,idh,jdh,a,b,mdab,ndab,wsha,lsh,work,lwork,ierr)
237    CALL error('shaes',ierr)
238  END IF
239
240!--- ALIASING DAMPLING FILTER
241!  WRITE(*,*)'>> Bi-Laplacian filtering to avoid Gibbs phenomenon...'
242!  CALL diffusion_filter(a,b)
243
244!=== INPUT/OUTPUT RESOLUTIONS AT EQUATOR (normalized with Earth radius)
245!--- The area of a cell at (lam,phi) of size(dl,dp) is: dl.SIN(dp).COS(phi)
246!--- At equator, Res=SQRT(dl.SIN(dp)) ; with dl=2.pi/nlon and dp=pi/nlat:
247!---             Res=SQRT(2.pi.SIN(pi/nlat)/nlon)~pi*SQRT(2/(nlon*nlat))
248!--- Note: nlat is the slices number = latitude points number -1
249  D1=SQRT(2.*xpi*SIN(xpi/(nlat_in-1))/nlon_in)   !--- Input  grid
250  D2=SQRT(2.*xpi*SIN(xpi/(nlat_ou-1))/nlon_ou)   !--- Output grid
251
252!=== FILTER FACTOR EQUAL TO ZERO => NO FILTERING
253  lfilter=ffactor/=0.
254
255!=== DEFAULT TRANSITION SCALE: INPUT GRID RESOLUTION
256  WRITE(*,*)'tfactor=',tfactor
257  del1=D1; del2=D1+(D2-D1)*tfactor
258  WRITE(*,*)'2*D1,del1,D2,del2=',2*D1,del1,D2,del2
259
260!=== SCALES SEPARATION AT THE OUTPUT GRID RESOLUTION ============================
261  nt=2
262  WRITE(*,*)'>> Scales separation at '//TRIM(real2str(D2/d2r))//' degrees (deta'&
263      //'ils below '//TRIM(real2str(60*D1/d2r))//' minutes removed)...'
264  ALLOCATE(as(mdab,ndab,nt),bs(mdab,ndab,nt))
265  DO n=1,nt; as(:,:,n)=a; bs(:,:,n)=b; END DO
266  CALL  lowpass_filter(as(:,:,1),bs(:,:,1),          D2,del2)   !--- LARGE SCALES
267  CALL bandpass_filter(as(:,:,2),bs(:,:,2),2*D1,del1,D2,del2)   !--- SMALL SCALES
268
269!--- BACKWARD TRANSFORM FOR LARGE AND SUB-CELL SCALES OROGRAPHYS
270  WRITE(*,*)'>> Backward transform of large and sub-cell scales orographies...'
271  ALLOCATE(hs(idh,jdh,nt))
272  IF(lfast) THEN
273    CALL shseci(nlat_in,nlon_in,wshs,lsh,work,lwork,ierr)
274    CALL error('shseci',ierr)
275    CALL shsec(nlat_in,nlon_in,isym,nt,hs,idh,jdh,as,bs,mdab,ndab,wshs,lsh,work,lwork,ierr)
276    CALL error('shsec',ierr)
277  ELSE
278    CALL shsesi(nlat_in,nlon_in,wshs,lsh,work,lwork,ierr)
279    CALL error('shsesi',ierr)
280    CALL shses(nlat_in,nlon_in,isym,nt,hs,idh,jdh,as,bs,mdab,ndab,wshs,lsh,work,lwork,ierr)
281    CALL error('shses',ierr)
282  END IF
283  WHERE(mski==0.); hs(:,:,1)=0.; hs(:,:,2)=0.; END WHERE
284  hs(:,:,2)=hs(:,:,2)**2
285
286!=== COMPUTE GRADIENT ; n=1: large scale ; n=2: sub-cell grid ===================
287  WRITE(*,*)'>> Computing large and sub-cell scales gradients...'
288  ALLOCATE(dxZ(idh,jdh,nt),dyZ(idh,jdh,nt))
289  IF(lfast) THEN
290    CALL vhseci(nlat_in,nlon_in,wvh,lvh,work,lwork,ierr)
291    CALL error('vhseci',ierr)
292    DO n=1,2
293      CALL gradec(nlat_in,nlon_in,isym,1,dyZ(:,:,n),dxZ(:,:,n),idh,jdh,         &
294                           as(:,:,n),bs(:,:,n),mdab,ndab,wvh,lvh,work,lwork,ierr)
295    END DO
296    CALL error('gradec',ierr)
297  ELSE
298    CALL vhsesi(nlat_in,nlon_in,wvh,lvh,work,lwork,ierr)
299    CALL error('vhsesi',ierr)
300    DO n=1,2
301      CALL gradec(nlat_in,nlon_in,isym,1,dyZ(:,:,n),dxZ(:,:,n),idh,jdh,as(:,:,n),&
302                                   bs(:,:,n),mdab,ndab,wvh,lvh,work,lwork,ierr)
303    END DO
304    CALL error('grades',ierr)
305  END IF
306  DEALLOCATE(as,bs)
307  dxZ=dxZ/REarth; dyZ=dyZ/REarth
308  IF(lat_in(nlat_in)>lat_in(1)) dyZ=-dyZ
309  WHERE(mski==0.)
310    dxZ(:,:,1)=0.; dyZ(:,:,1)=0.
311    dxZ(:,:,2)=0.; dyZ(:,:,2)=0.
312  END WHERE
313  ALLOCATE(dH(idh,jdh,6))
314  dH(:,:,1)=dxZ(:,:,1)*dxZ(:,:,1); dH(:,:,4)=dxZ(:,:,2)*dxZ(:,:,2)
315  dH(:,:,2)=dyZ(:,:,1)*dyZ(:,:,1); dH(:,:,5)=dyZ(:,:,2)*dyZ(:,:,2)
316  dH(:,:,3)=dxZ(:,:,1)*dyZ(:,:,1); dH(:,:,6)=dxZ(:,:,2)*dyZ(:,:,2)
317  DEALLOCATE(dxZ,dyZ)
318
319!===============================================================================
320  IF(lfilter) THEN    !=== FILTER TO REMOVE ALIASING ===========================
321!===============================================================================
322    WRITE(*,*)'>> Smoothing large and sub-cell scales fields at '//            &
323              TRIM(real2str(ffactor*D2/d2r))//' degrees...'
324    CALL spatial_filter_m(dH       ,ffactor*D2,del2)
325    CALL spatial_filter_1(hs(:,:,2),ffactor*D2,del2)
326  END IF
327  DEALLOCATE(work)
328!===============================================================================
329
330!=== DEAL WITH OUTPUT GRID =====================================================
331  intl=0
332!--- OUTPUT GRID (CONTAINS A LONGITUDE REDUNDANT POINT FOR NOW)
333  lon_ou = [(360.*DBLE(n)/DBLE(nlon_ou  ), n = 0, nlon_ou)] + lon_in(1)
334  lat_ou = [(180.*DBLE(n)/DBLE(nlat_ou-1), n = nlat_ou-1, 0, -1)] - 90.
335  la1 = MIN(nlat_in,(nlon_in+1)/2); la2 = (nlat_in+1)/2
336  lb1 = MIN(nlat_ou,(nlon_ou+1)/2); lb2 = (nlat_ou+1)/2
337  lwa = 4*nlat_in*la2+3*MAX(la1-2,0)*(2*nlat_in-la1-1)+la2+nlon_in+15
338  lwb = 4*nlat_ou*lb2+3*MAX(lb1-2,0)*(2*nlat_ou-lb1-1)+    nlon_ou+15
339  lsave=lwa+lwb
340  l1 = MAX(nlat_in,nlat_ou)
341  l2 = MAX(nlon_in,nlon_ou)
342  lwork  = 2*l1*(8*MIN(l1,(l2+2)/2)+4*l2+3)
343  ldwork = 2*l1*(l1+1)+1
344  ALLOCATE(wsave(lsave),work(lwork),dwork(ldwork))
345
346!--- Latitude grid: equally spaced partition of [-pi/2,pi/2]
347  igrid_in(:)=[-1,1]                                  !--- NP->SP, nlat*nlon
348  igrid_ou(:)=[-1,0]                                  !--- NP->SP, nlon*nlat
349
350!--- BUILD FRACTIONAL LAND MASK msko AND GEOPOTENTIAL ('noro' CASE)
351  SELECT CASE(fmsk)
352    CASE('noro')
353      WRITE(*,*)'>> Computing mask & geop (grid_noro compatible)...'
354      ALLOCATE(msko(nlon_ou+1,nlat_ou), Zphi(nlon_ou+1,nlat_ou))
355      CALL flip(h)
356      CALL grid_noro0(lon_in*d2r, lat_in*d2r, h,lon_ou*d2r, lat_ou*d2r, Zphi, msko)
357      msko = msko(1:nlon_ou,:)
358      Zphi = Zphi(1:nlon_ou,:)
359    CASE('spec')
360      WRITE(*,*)'>> Computing mask in spectral space...'
361      ALLOCATE(msko(nlon_ou,nlat_ou))
362      CALL trssph(intl, igrid_in, nlon_in, nlat_in, mski,                      &
363                        igrid_ou, nlon_ou, nlat_ou, msko, wsave, lsave, lsvmin,&
364                                     work, lwork, lwkmin, dwork, ldwork, ierr)
365      intl=1
366    CASE DEFAULT
367      WRITE(*,*)'>> Trying to read mask from "'//TRIM(fmsk)//'" file'
368      ALLOCATE(msko(nlon_ou,nlat_ou))
369      msg='Missing or wrong "-m" option ; can be "noro", "spec" or a mask file'
370      CALL err(nf90_open(fmsk,nf90_nowrite,fID)/=nf90_noerr,msg)
371      CALL nc(nf90_inq_varid(fID,"MaskOcean",vID),"MaskOcean")    !--- MASK ID
372      CALL nc(nf90_get_var(fID,vID,msko(:,:)))                    !--- MASK
373      CALL nc(nf90_close(fID))
374      msko(:,:)=1.0-msko(:,:)
375  END SELECT
376  DEALLOCATE(mski,h)
377
378!--- GET RID OF REDUNDANT LONGITUDE (USELESS FOR SPHEREPACK)
379  lon_ou=lon_ou(1:nlon_ou)
380
381!=== REGRIDDING ON THE TARGET OUTPUT GRID ======================================
382  WRITE(*,*)'>> Regridding on the target grid...'
383
384!--- COMPUTE BOOLEAN OUTPUT MASK
385  ALLOCATE(mask(nlon_ou,nlat_ou)); mask=0_1; WHERE(msko>=0.5) mask=1_1
386
387!--- MEAN OROGRAPHY
388  ALLOCATE(h0(nlon_ou,nlat_ou))
389  CALL trssph(intl, igrid_in, nlon_in, nlat_in, hs(:,:,1),                     &
390                    igrid_ou, nlon_ou, nlat_ou, h0,                            &
391    wsave, lsave, lsvmin, work, lwork, lwkmin, dwork, ldwork, ierr)
392  WHERE(h0<0..OR.mask==0_1) h0(:,:)=0.
393
394!--- STANDARD DEVIATION OF SUB-CELL SCALES OROGRAPHY
395  ALLOCATE(mu(nlon_ou,nlat_ou))
396  CALL trssph(1, igrid_in, nlon_in, nlat_in, hs(:,:,2),                        &
397                 igrid_ou, nlon_ou, nlat_ou, mu,                               &
398    wsave, lsave, lsvmin, work, lwork, lwkmin, dwork, ldwork, ierr)
399  DEALLOCATE(hs)
400  WHERE(mu(:,:)<0..OR.mask==0_1) mu(:,:)=0.; mu=SQRT(mu)
401
402!--- SQUARE OF LARGE AND SUB-CELL SCALES ZONAL SLOPE
403  ALLOCATE(dxZ2(nlon_ou,nlat_ou,2))
404  CALL trssph(1, igrid_in, nlon_in, nlat_in, dH  (:,:,1),                      &
405                 igrid_ou, nlon_ou, nlat_ou, dxZ2(:,:,1),                      &
406    wsave, lsave, lsvmin, work, lwork, lwkmin, dwork, ldwork, ierr)
407  CALL trssph(1, igrid_in, nlon_in, nlat_in, dH  (:,:,4),                      &
408                 igrid_ou, nlon_ou, nlat_ou, dxZ2(:,:,2),                      &
409    wsave, lsave, lsvmin, work, lwork, lwkmin, dwork, ldwork, ierr)
410  WHERE(mask==0_1); dxZ2(:,:,1)=0.; dxZ2(:,:,2)=0.; END WHERE
411
412!--- SQUARE OF LARGE AND SUB-CELL SCALES MERIDIONAL SLOPE
413  ALLOCATE(dyZ2(nlon_ou,nlat_ou,2))
414  CALL trssph(1, igrid_in, nlon_in, nlat_in, dH  (:,:,2),                      &
415                 igrid_ou, nlon_ou, nlat_ou, dyZ2(:,:,1),                      &
416    wsave, lsave, lsvmin, work, lwork, lwkmin, dwork, ldwork, ierr)
417  CALL trssph(1, igrid_in, nlon_in, nlat_in, dH  (:,:,5),                      &
418                 igrid_ou, nlon_ou, nlat_ou, dyZ2(:,:,2),                      &
419    wsave, lsave, lsvmin, work, lwork, lwkmin, dwork, ldwork, ierr)
420  WHERE(mask==0_1); dyZ2(:,:,1)=0.; dyZ2(:,:,2)=0.; END WHERE
421
422!--- PRODUCT OF LARGE AND SUB-CELL SCALES MERIDIONAL AND ZONAL SLOPES
423  ALLOCATE(dxyZ(nlon_ou,nlat_ou,2))
424  CALL trssph(1, igrid_in, nlon_in, nlat_in, dH  (:,:,3),                      &
425                 igrid_ou, nlon_ou, nlat_ou, dxyZ(:,:,1),                      &
426    wsave, lsave, lsvmin, work, lwork, lwkmin, dwork, ldwork, ierr)
427  CALL trssph(1, igrid_in, nlon_in, nlat_in, dH  (:,:,6),                      &
428                 igrid_ou, nlon_ou, nlat_ou, dxyZ(:,:,2),                      &
429    wsave, lsave, lsvmin, work, lwork, lwkmin, dwork, ldwork, ierr)
430  WHERE(mask==0_1); dxyZ(:,:,1)=0.; dxyZ(:,:,2)=0.; END WHERE
431  DEALLOCATE(dH)
432
433!--- COMPUTE LARGE AND SUB-CELL SCALES DERIVED QUANTITIES
434  ALLOCATE(Xk(nlon_ou,nlat_ou,2)); Xk=(dxZ2+dyZ2)/2.
435  ALLOCATE(Xl(nlon_ou,nlat_ou,2)); Xl=(dxZ2-dyZ2)/2.
436  DEALLOCATE(dxZ2,dyZ2)
437  ALLOCATE(Xm(nlon_ou,nlat_ou,2)); Xm=dxyZ
438  DEALLOCATE(dxyZ)
439  ALLOCATE(Xp(nlon_ou,nlat_ou,2),Xq(nlon_ou,nlat_ou,2))
440  Xp=Xk; Xq=Xk; Xk=SQRT(Xl**2+Xm**2); Xp=Xp-Xk; Xq=Xq+Xk
441  DEALLOCATE(Xk)
442  WHERE(ABS(Xm)<=1e-16) Xm=1e-16*SIGN(1.,Xm)
443  ALLOCATE(Zthe(nlon_ou,nlat_ou)); Zthe=0.5*ATAN2(Xm(:,:,1),Xl(:,:,1))/d2r
444  DEALLOCATE(Xm,Xl)
445  ALLOCATE(Zgam(nlon_ou,nlat_ou)); Zgam=0.
446  WHERE(Xq(:,:,1)>0..AND.Xp(:,:,1)>1.E-8) Zgam(:,:)=Xp(:,:,1)/Xq(:,:,1)
447  DEALLOCATE(Xp)
448  ALLOCATE(Zsig(nlon_ou,nlat_ou)); Zsig=SQRT(Xq(:,:,2))
449  WHERE(Xq(:,:,2)<=1.E-16) Zsig(:,:)=1.E-8
450  DEALLOCATE(Xq)
451  WHERE(mask==0_1); Zthe=0.; Zgam=0.; Zsig=0.; END WHERE
452
453!--- Save the data
454  WRITE(*,*)'>> Saving fields in a netcdf file...'
455  res_in=TRIM(int2str(nlon_in))//'x'//TRIM(int2str(nlat_in))
456  res_ou=TRIM(int2str(nlon_ou))//'x'//TRIM(int2str(nlat_ou))
457  f_ou='make_sso_'//TRIM(res_ou)//'_'//TRIM(f_in)
458  fnam=f_ou
459  CALL nc(nf90_create(f_ou,nf90_clobber,fID))
460
461  CALL nc(nf90_def_dim(fID,'x',nlon_ou,xID))
462  CALL nc(nf90_def_var(fID,'x',nf90_real,xID,loID)      ,'x')
463  CALL nc(nf90_put_att(fID,loID,'long_name','Longitude'),'x')
464  CALL nc(nf90_put_att(fID,loID,'units','degrees_east') ,'x')
465
466  CALL nc(nf90_def_dim(fID,'y',nlat_ou,yID))
467  CALL nc(nf90_def_var(fID,'y',nf90_real,yID,laID)      ,'y')
468  CALL nc(nf90_put_att(fID,laID,'long_name','Latitude') ,'y')
469  CALL nc(nf90_put_att(fID,laID,'units','degrees_north'),'y')
470
471  CALL nc(nf90_def_var(fID,'mask',nf90_real,[xID,yID],mskID),'mask')
472  IF(fmsk=='noro') &
473  CALL nc(nf90_def_var(fID,'Zphi',nf90_real,[xID,yID],phiID),'Zphi')
474  CALL nc(nf90_def_var(fID,'Zmea',nf90_real,[xID,yID],meaID),'Zmea')
475  CALL nc(nf90_def_var(fID,'mu'  ,nf90_real,[xID,yID], muID),'mu'  )
476  CALL nc(nf90_def_var(fID,'Zsig',nf90_real,[xID,yID],sigID),'Zsig')
477  CALL nc(nf90_def_var(fID,'Zgam',nf90_real,[xID,yID],gamID),'Zgam')
478  CALL nc(nf90_def_var(fID,'Zthe',nf90_real,[xID,yID],theID),'Zthe')
479  CALL nc(nf90_def_var(fID,'Zpic',nf90_real,[xID,yID],picID),'Zpic')
480  CALL nc(nf90_def_var(fID,'Zval',nf90_real,[xID,yID],valID),'Zval')
481
482  CALL nc(nf90_put_att(fID,mskID,'long_name','Fractional land mask'                      ),'mask')
483  IF(fmsk=='noro') &
484  CALL nc(nf90_put_att(fID,phiID,'long_name','Geopotential'                              ),'Zphi')
485  CALL nc(nf90_put_att(fID,meaID,'long_name','Mean orography'                            ),'Zmea')
486  CALL nc(nf90_put_att(fID, muID,'long_name','Std deviation of sub-cell scales orography'),'mu'  )
487  CALL nc(nf90_put_att(fID,sigID,'long_name','Slope along principal axis'                ),'Zsig')
488  CALL nc(nf90_put_att(fID,gamID,'long_name','Anisotropy (aspect ratio)'                 ),'Zgam')
489  CALL nc(nf90_put_att(fID,theID,'long_name','Orientation (principal axis)'              ),'Zthe')
490  CALL nc(nf90_put_att(fID,picID,'long_name','Maximum height'                            ),'Zpic')
491  CALL nc(nf90_put_att(fID,valID,'long_name','Minimum height'                            ),'Zval')
492
493  CALL nc(nf90_put_att(fID,mskID,'units','none'   ),'mask')
494  IF(fmsk=='noro') &
495  CALL nc(nf90_put_att(fID,phiID,'units','m'      ),'Zphi')
496  CALL nc(nf90_put_att(fID,meaID,'units','m'      ),'Zmea')
497  CALL nc(nf90_put_att(fID, muID,'units','m'      ),'mu'  )
498  CALL nc(nf90_put_att(fID,sigID,'units','m/m'    ),'Zsig')
499  CALL nc(nf90_put_att(fID,gamID,'units','none'   ),'Zgam')
500  CALL nc(nf90_put_att(fID,theID,'units','degrees'),'Zthe')
501  CALL nc(nf90_put_att(fID,picID,'units','m'      ),'Zpic')
502  CALL nc(nf90_put_att(fID,valID,'units','m'      ),'Zval')
503
504  CALL nc(nf90_put_att(fID,nf90_global,'Conventions','COARDS/CF-1.0'))
505  CALL nc(nf90_put_att(fID,nf90_global,'Initial_Grid',TRIM(res_in)))
506  CALL nc(nf90_put_att(fID,nf90_global,'history',TRIM(call_seq)))
507  CALL nc(nf90_enddef(fID))
508
509  CALL nc(nf90_put_var(fID, loID,lon_ou),'x' )
510  CALL nc(nf90_put_var(fID, laID,lat_ou),'y' )
511  CALL nc(nf90_put_var(fID,mskID,msko),'mask')
512  IF(fmsk=='noro') &
513  CALL nc(nf90_put_var(fID,phiID,Zphi),'Zphi')
514  CALL nc(nf90_put_var(fID,meaID,h0  ),'Zmea')
515  CALL nc(nf90_put_var(fID, muID,mu  ),'mu'  )
516  CALL nc(nf90_put_var(fID,sigID,Zsig),'Zsig')
517  CALL nc(nf90_put_var(fID,gamID,Zgam),'Zgam')
518  CALL nc(nf90_put_var(fID,theID,Zthe),'Zthe')
519  CALL nc(nf90_put_var(fID,picID,       h0+2*mu ),'Zpic')
520  CALL nc(nf90_put_var(fID,valID,MAX(0.,h0-2*mu)),'Zval')
521  CALL nc(nf90_close(fID))
522  WRITE(*,*)'Finished.'
523
524CONTAINS
525
526
527!-------------------------------------------------------------------------------
528
529SUBROUTINE diffusion_filter(a,b)
530
531!-------------------------------------------------------------------------------
532! Purpose: ECMWF bi-laplacian diffusion filter to limit aliasing.
533! Note that one third of the frequencies is skipped to allow some accuracy for
534! derived quantities computation.
535!-------------------------------------------------------------------------------
536  IMPLICIT NONE
537!-------------------------------------------------------------------------------
538! Arguments:
539  REAL, INTENT(INOUT) :: a(:,:), b(:,:)
540!-------------------------------------------------------------------------------
541! Local variables:
542  INTEGER :: m
543  REAL, ALLOCATABLE :: fn(:)
544!-------------------------------------------------------------------------------
545  ALLOCATE(fn(ndab)); fn=[(REAL(m*(m-1)),m=1,ndab)]
546  a(1,:)=a(1,:)/(1+4.*fn(:)/fn(ndab))**2
547  b(1,:)=b(1,:)/(1+4.*fn(:)/fn(ndab))**2
548  DO m=1,ndab
549    IF(m<=2*ndab/3) THEN
550      a(m,:)=a(m,:)/(1+5.*fn(:)/fn(ndab))**2
551      b(m,:)=b(m,:)/(1+5.*fn(:)/fn(ndab))**2
552    ELSE
553      a(m,:)=0.
554      b(m,:)=0.
555    END IF
556  END DO
557
558END SUBROUTINE diffusion_filter
559
560!-------------------------------------------------------------------------------
561
562
563!-------------------------------------------------------------------------------
564
565SUBROUTINE lowpass_filter(a,b,d,dd)
566
567!-------------------------------------------------------------------------------
568! Purpose: ECMWF lowpass filter.
569!-------------------------------------------------------------------------------
570  IMPLICIT NONE
571!-------------------------------------------------------------------------------
572! Arguments:
573  REAL, INTENT(INOUT) :: a(:,:), b(:,:)
574  REAL, INTENT(IN)    :: d, dd
575!-------------------------------------------------------------------------------
576! Local variables:
577  INTEGER :: m, na
578  REAL, ALLOCATABLE  :: h(:)
579!-------------------------------------------------------------------------------
580  na=SIZE(a,1)
581  IF(ALLOCATED(h)) THEN; IF(SIZE(h)/=na) DEALLOCATE(h); END IF
582  IF(.NOT.ALLOCATED(h)) ALLOCATE(h(na))
583  CALL lowpass(d,dd,h)
584  DO m=1,SIZE(h); a(m,:)=a(m,:)*h; b(m,:)=b(m,:)*h; END DO
585
586END SUBROUTINE lowpass_filter
587
588!-------------------------------------------------------------------------------
589
590
591!-------------------------------------------------------------------------------
592
593SUBROUTINE bandpass_filter(a,b,d1,dd1,d2,dd2)
594
595!-------------------------------------------------------------------------------
596! Purpose: ECMWF bandpass filter. Calls lowpass filter with 2 different scales.
597!-------------------------------------------------------------------------------
598  IMPLICIT NONE
599!-------------------------------------------------------------------------------
600! Arguments:
601  REAL, INTENT(INOUT) :: a(:,:), b(:,:)
602  REAL, INTENT(IN)    :: d1, dd1, d2, dd2
603!-------------------------------------------------------------------------------
604! Local variables:
605  INTEGER :: m, na
606  REAL, ALLOCATABLE  :: h(:), hh(:)
607!-------------------------------------------------------------------------------
608  na=SIZE(a,1)
609  IF(ALLOCATED(hh)) THEN; IF(SIZE(hh)/=na) DEALLOCATE(hh); END IF
610  IF(ALLOCATED(h )) THEN; IF(SIZE(h )/=na) DEALLOCATE(h ); END IF
611  IF(.NOT.ALLOCATED(hh)) ALLOCATE(hh(na))
612  IF(.NOT.ALLOCATED(h )) ALLOCATE(h (na))
613  CALL lowpass(d1,dd1,hh)
614  CALL lowpass(d2,dd2,h ); hh=hh-h
615  DO m=1,SIZE(h); a(m,:)=a(m,:)*hh; b(m,:)=b(m,:)*hh; END DO
616
617END SUBROUTINE bandpass_filter
618
619!-------------------------------------------------------------------------------
620
621
622!-------------------------------------------------------------------------------
623
624SUBROUTINE lowpass(d,dd,h)
625
626!-------------------------------------------------------------------------------
627  IMPLICIT NONE
628!-------------------------------------------------------------------------------
629! Arguments:
630  REAL, INTENT(IN)    :: d, dd   ! d is DELTA/REarth ; dd is delta/REarth
631  REAL, INTENT(INOUT) :: h(:)
632!-------------------------------------------------------------------------------
633  INTEGER :: m, n
634  REAL    :: hd
635!-------------------------------------------------------------------------------
636  n=SIZE(h); hd=0.5*d
637  h(:)=[((SIN(m*(hd-dd))/m+(COS(hp+m*hd)*SIN(hp+m*dd))/(hp/dd+m)               &
638        + SIN(m*(hd+dd))/m+(COS(hp-m*hd)*SIN(hp-m*dd))/(hp/dd-m))**2/d/d,m=1,n)]
639
640END SUBROUTINE lowpass
641
642!-------------------------------------------------------------------------------
643
644
645!-------------------------------------------------------------------------------
646
647SUBROUTINE spatial_filter_1(u,D,dd)
648
649!-------------------------------------------------------------------------------
650! Purpose: Spatial filter ; using spectral lowpass filter.
651!-------------------------------------------------------------------------------
652  IMPLICIT NONE
653!-------------------------------------------------------------------------------
654! Arguments:
655  REAL, INTENT(INOUT) :: u(:,:)
656  REAL, INTENT(IN)    :: D, dd
657!-------------------------------------------------------------------------------
658! Local variables:
659  REAL, ALLOCATABLE  :: uu(:,:,:)
660!-------------------------------------------------------------------------------
661  ALLOCATE(uu(SIZE(u,1),SIZE(u,2),1)); uu(:,:,1)=u
662  CALL spatial_filter_m(uu,D,dd)
663  u=uu(:,:,1)
664  DEALLOCATE(uu)
665
666END SUBROUTINE spatial_filter_1
667
668!-------------------------------------------------------------------------------
669
670
671!-------------------------------------------------------------------------------
672
673SUBROUTINE spatial_filter_m(u,D,dd)
674
675!-------------------------------------------------------------------------------
676! Purpose: Spatial filter ; using spectral lowpass filter.
677!-------------------------------------------------------------------------------
678  IMPLICIT NONE
679!-------------------------------------------------------------------------------
680! Arguments:
681  REAL, INTENT(INOUT) :: u(:,:,:)
682  REAL, INTENT(IN)    :: D, dd
683!-------------------------------------------------------------------------------
684! Local variables:
685  INTEGER :: ku
686  REAL, ALLOCATABLE  :: aa(:,:,:),bb(:,:,:)
687!-------------------------------------------------------------------------------
688  ku=SIZE(u,3)                                   !--- NUMBER OF FIELDS TO SMOOTH
689  ALLOCATE(aa(mdab,ndab,ku),bb(mdab,ndab,ku))
690
691!--- PROJECTION OF THE SUB-CELL SCALES FIELDS ON THE SPHERICAL HARMONICS BASE.
692  IF(lfast) THEN
693    CALL shaec(nlat_in,nlon_in,0,ku,u,idh,jdh,aa,bb,mdab,ndab,wsha,lsh,work,lwork,ierr)
694    CALL error('shaec',ierr)
695  ELSE
696    CALL shaes(nlat_in,nlon_in,0,ku,u,idh,jdh,aa,bb,mdab,ndab,wsha,lsh,work,lwork,ierr)
697    CALL error('shaes',ierr)
698  END IF
699
700!--- LOW PASS FILTER AT THE OUTPUT GRID RESOLUTION
701  DO n=1,ku; CALL lowpass_filter(aa(:,:,n),bb(:,:,n),D,dd); END DO
702
703!--- BACKWARD TRANSFORM FOR MEAN SUB-CELL SCALES PARAMETERS
704  IF(lfast) THEN
705    CALL shsec(nlat_in,nlon_in,0,ku,u,idh,jdh,aa,bb,mdab,ndab,wshs,lsh,work,lwork,ierr)
706    CALL error('shsec',ierr)
707  ELSE
708    CALL shses(nlat_in,nlon_in,0,ku,u,idh,jdh,aa,bb,mdab,ndab,wshs,lsh,work,lwork,ierr)
709    CALL error('shses',ierr)
710  END IF
711  DEALLOCATE(aa,bb)
712
713END SUBROUTINE spatial_filter_m
714
715!-------------------------------------------------------------------------------
716
717
718!-------------------------------------------------------------------------------
719
720SUBROUTINE grid_noro0(xin,yin,zin,xou,you,zphi,mask)
721
722!-------------------------------------------------------------------------------
723! Purpose: Sub-cell scales orographic parameters coomputation. Angles in radians
724! Notes: Input coordinates contain a longitude periodic redundant point.
725!-------------------------------------------------------------------------------
726  IMPLICIT NONE
727!-------------------------------------------------------------------------------
728! Arguments:
729  REAL, INTENT(IN)  :: xin(:), yin(:)   !--- INPUT  COORDINATES (nxi)   (nyi)
730  REAL, INTENT(IN)  :: zin(:,:)         !--- INPUT  FIELD       (nxi,    nyi)
731  REAL, INTENT(IN)  :: xou(:), you(:)   !--- OUTPUT COORDINATES (nxo+1) (nyo)
732  REAL, INTENT(OUT) :: zphi(:,:)        !--- GEOPOTENTIAL       (nxo+1,  nyo)
733  REAL, INTENT(OUT) :: mask(:,:)        !--- FRACTIONAL MASK    (nxo+1,  nyo)
734!-------------------------------------------------------------------------------
735! Local variables:
736  REAL, ALLOCATABLE :: xi(:), xo_w(:), yo_s(:), zi(:,:), num_tot(:,:)
737  REAL, ALLOCATABLE :: yi(:), xo_e(:), yo_n(:),  w(:,:), num_lan(:,:)
738  REAL    :: dx, dy, wx, wy, zx, zy, ex, ey
739  INTEGER :: nxi, nxo, i, ii, nn
740  INTEGER :: nyi, nyo, j, jj, di
741!-------------------------------------------------------------------------------
742  nxi=SIZE(xin); nxo=SIZE(xou)-1; ex=SIGN(1.,xin(nxi)-xin(1))
743  nyi=SIZE(yin); nyo=SIZE(you);   ey=SIGN(1.,yin(nyi)-yin(1))
744
745!--- MARGIN TO ENSURE OUTPUT CELLS HAVE MATCHING INPUT CELLS
746  di=CEILING(1.+MAXVAL([xou(2:nxo+1)-xou(1:nxo)])/MINVAL([xin(2:nxi)-xin(1:nxi-1)]))
747
748!--- EXTENSION OF THE INPUT DATABASE TO PROCEED COMPUTATIONS AT BOUNDARIES
749  ALLOCATE(xi(nxi+2*di)); xi(:)=[xin(nxi-di:nxi-1)-2.*xpi,xin,xin(2:di+1)+2.*xpi]
750  ALLOCATE(yi(nyi+2));    yi(:)=[2.*yin(1)-yin(2),yin,2.*yin(nyi)-yin(nyi-1)]
751
752  ALLOCATE(zi(nxi+2*di,nyi+2))
753  zi(1      +di:nxi+  di,2:nyi+1)=zin(:,:)
754  zi(1         :      di,2:nyi+1)=zin(nxi-di:nxi-1,:)
755  zi(1+nxi  +di:nxi+2*di,2:nyi+1)=zin(2     :di+1 ,:)
756  zi(1         :nxi/2+di,      1)=zi (1+nxi/2:nxi  +di,    2)
757  zi(1+nxi/2+di:nxi+2*di,      1)=zi (1      :nxi/2+di,    2)
758  zi(1         :nxi/2+di,  nyi+2)=zi (1+nxi/2:nxi  +di,nyi+1)
759  zi(1+nxi/2+di:nxi+2*di,  nyi+2)=zi (1      :nxi/2+di,nyi+1)
760
761!  write(*,'(a,2f10.5,a,2f10.5)')'xin=',xin(1:2),'...',xin(nxi-1:nxi)
762!  write(*,'(a,2f10.5,a,2f10.5)')'yin=',yin(1:2),'...',yin(nyi-1:nyi)
763!  write(*,'(a,2f10.5,a,2f10.5)')'xou=',xou(1:2),'...',xou(nxo:nxo+1)
764!  write(*,'(a,2f10.5,a,2f10.5)')'you=',you(1:2),'...',you(nyo-1:nyo)
765!  write(*,*)'xi=',xi(1:di+1),'...',xi(nxi+di:nxi+2*di)
766!  write(*,'(a,3f10.5,a,3f10.5)')'yi=',yi(1:3),'...',yi(nyi:nyi+2)
767
768!--- COMPUTE LIMITS OF MODEL GRIDPOINT AREA (REGULAR GRID) bdac
769  ALLOCATE(xo_w(nxo+1)); xo_w=(xou+[2.*xou(    1)-xou(  2),xou(1:nxo  )])/2.0
770  ALLOCATE(yo_s(nyo  )); yo_s=(you+[2.*you(    1)-you(  2),you(1:nyo-1)])/2.0
771  ALLOCATE(xo_e(nxo+1)); xo_e=(xou+[xou(2:nxo+1),2.*xou(nxo+1)-xou(nxo)])/2.0
772  ALLOCATE(yo_n(nyo  )); yo_n=(you+[you(2:nyo  ),2.*you(nyo)-you(nyo-1)])/2.0
773
774!--- INITIALIZATIONS:
775  ALLOCATE(w(nxo+1,nyo)); w(:,:)=0.0; zphi(:,:)=0.0
776
777!--- SUMMATION OVER GRIDPOINT AREA
778  zx=2.*xpi/REAL(nxi)                            !--- LON CELL LENGTH (1-SPHERE)
779  dx=zx/2.                                       !--- HALF "
780  ALLOCATE(num_tot(nxo+1,nyo)); num_tot(:,:)=0.
781  ALLOCATE(num_lan(nxo+1,nyo)); num_lan(:,:)=0.
782  DO ii=1,nxo+1
783    DO jj=1,nyo
784      DO j=2,nyi+1
785        zy=xpi/REAL(nyi)*COS(yi(j))              !--- LAT CELL LENGTH (1-SPHERE)
786        dy=zy/2.
787        wy=MAX(0.,MIN(zy,MINVAL(dy+ey*[yo_n(jj)-yi(j),yi(j)-yo_s(jj)])))
788        IF(wy==0.) CYCLE
789        DO i=2,nxi+2*di-1
790          wx=MAX(0.,MIN(zx,MINVAL(dx+ex*[xo_e(ii)-xi(i),xi(i)-xo_w(ii)])*COS(yi(j))))
791          IF(wx==0.) CYCLE
792          num_tot(ii,jj)=num_tot(ii,jj)+1.0
793          IF(zi(i,j)>=1.) num_lan(ii,jj)=num_lan(ii,jj)+1.0
794          w   (ii,jj)=   w(ii,jj)+wx*wy
795          zphi(ii,jj)=zphi(ii,jj)+wx*wy*zi(i,j)
796        END DO
797      END DO
798    END DO
799  END DO
800
801!--- COMPUTE FRACTIONAL MASK AND GEOPOTENTIAL
802  WHERE(w(:,:)/=0.0) mask=num_lan(:,:)/num_tot(:,:)
803  nn=COUNT(w(:,:)==0.0)
804  IF(nn/=0) WRITE(*,*)'Problem with weight ; vanishing occurrences: ',nn
805  WHERE(w/=0.0) zphi=zphi/w
806  WHERE(w==0.0.OR.mask<0.1) zphi=0.0
807  zphi(nxo+1,:)=zphi(1,:)
808  zphi(:,  1)=SUM(w(1:nxo,  1)*zphi(1:nxo,  1))/SUM(w(1:nxo,  1))
809  zphi(:,nyo)=SUM(w(1:nxo,nyo)*zphi(1:nxo,nyo))/SUM(w(1:nxo,nyo))
810
811END SUBROUTINE grid_noro0
812
813!-------------------------------------------------------------------------------
814
815
816!-------------------------------------------------------------------------------
817
818SUBROUTINE nc(ncres,var)
819
820!-------------------------------------------------------------------------------
821! Purpose: NetCDF errors handling.
822!-------------------------------------------------------------------------------
823  IMPLICIT NONE
824!-------------------------------------------------------------------------------
825! Arguments:
826  INTEGER,          INTENT(IN)           :: ncres
827  CHARACTER(LEN=*), OPTIONAL, INTENT(IN) :: var
828!-------------------------------------------------------------------------------
829! Local variables:
830  CHARACTER(LEN=256) :: msg
831!-------------------------------------------------------------------------------
832  IF(ncres/=nf90_noerr) THEN
833    msg='Error in routine '//TRIM(sub)
834    IF(fnam/='')     msg=TRIM(msg)//' for file "'//TRIM(fnam)//'"'
835    IF(PRESENT(var)) msg=TRIM(msg)//'" and variable "'//TRIM(var)//'"'
836    WRITE(*,*)TRIM(msg)//': '//nf90_strerror(ncres); STOP
837  END IF
838
839END SUBROUTINE nc
840
841!-------------------------------------------------------------------------------
842
843
844!-------------------------------------------------------------------------------
845
846ELEMENTAL FUNCTION int2str(val)
847
848!-------------------------------------------------------------------------------
849  IMPLICIT NONE
850!-------------------------------------------------------------------------------
851! Purpose: convert an integer into a string.
852!-------------------------------------------------------------------------------
853! Arguments:
854  CHARACTER(LEN=128)  :: int2str
855  INTEGER, INTENT(IN) :: val
856!-------------------------------------------------------------------------------
857! Local variable:
858  CHARACTER(LEN=128) :: s1
859!-------------------------------------------------------------------------------
860  WRITE(s1,*)val; int2str=TRIM(ADJUSTL(s1))
861
862END FUNCTION int2str
863!-------------------------------------------------------------------------------
864
865
866!-------------------------------------------------------------------------------
867
868SUBROUTINE err(ierr,str)
869
870!-------------------------------------------------------------------------------
871  IMPLICIT NONE
872!-------------------------------------------------------------------------------
873! Argument:
874  LOGICAL,          INTENT(IN) :: ierr
875  CHARACTER(LEN=*), INTENT(IN) :: str
876!-------------------------------------------------------------------------------
877  IF(.NOT.ierr) RETURN
878  WRITE(*,*)TRIM(str)
879  STOP
880
881END SUBROUTINE err
882
883!-------------------------------------------------------------------------------
884
885
886!-------------------------------------------------------------------------------
887
888SUBROUTINE error(sub,ierr)
889
890!-------------------------------------------------------------------------------
891  IMPLICIT NONE
892!-------------------------------------------------------------------------------
893! Arguments:
894  CHARACTER(LEN=*), INTENT(IN) :: sub
895  INTEGER,          INTENT(IN) :: ierr
896!-------------------------------------------------------------------------------
897  IF(ierr==0) RETURN
898  WRITE(*,*)'In routine "'//TRIM(sub)//'": error in the specification of argume&
899       nt nr. '//TRIM(int2str(ierr))
900  STOP
901END SUBROUTINE error
902
903!-------------------------------------------------------------------------------
904
905
906!-------------------------------------------------------------------------------
907
908ELEMENTAL FUNCTION real2str(val)
909
910!-------------------------------------------------------------------------------
911  IMPLICIT NONE
912!-------------------------------------------------------------------------------
913! Purpose: convert an integer into a string.
914!-------------------------------------------------------------------------------
915! Arguments:
916  CHARACTER(LEN=128) :: real2str
917  REAL,   INTENT(IN) :: val
918  INTEGER :: r, s
919!-------------------------------------------------------------------------------
920                s=NINT(val); r=NINT(1000.*(val-REAL(s)))
921  IF(r<0) THEN; s= INT(val); r= INT(1000.*(val-REAL(s))); END IF
922  real2str=TRIM(ADJUSTL(int2str(s)))//'.'//TRIM(ADJUSTL(int2str(r)))
923
924END FUNCTION real2str
925!-------------------------------------------------------------------------------
926
927
928!-------------------------------------------------------------------------------
929
930SUBROUTINE flip(a)
931
932!-------------------------------------------------------------------------------
933  IMPLICIT NONE
934!-------------------------------------------------------------------------------
935! Arguments:
936  REAL, ALLOCATABLE, INTENT(INOUT) :: a(:,:)
937!-------------------------------------------------------------------------------
938! Local variables:
939  INTEGER :: k, n
940  REAL, ALLOCATABLE :: tmp(:,:)
941!-------------------------------------------------------------------------------
942  n=SIZE(a,2)
943  ALLOCATE(tmp(n,SIZE(a,1))); DO k=1,n; tmp(k,:)=a(:,k); END DO
944  CALL MOVE_ALLOC(FROM=tmp,TO=a)
945
946END SUBROUTINE flip
947
948!-------------------------------------------------------------------------------
949
950
951!-------------------------------------------------------------------------------
952ELEMENTAL LOGICAL FUNCTION is_int(str) RESULT(out)
953!-------------------------------------------------------------------------------
954  IMPLICIT NONE
955!-------------------------------------------------------------------------------
956! Purpose: tests if str is a valid integer number
957!-------------------------------------------------------------------------------
958  CHARACTER(LEN=*), INTENT(IN) :: str
959  INTEGER :: k, k0
960!-------------------------------------------------------------------------------
961  out=LEN_TRIM(str)/=0; IF(.NOT.out) RETURN
962  k0=1; IF(INDEX('+-',str(1:1))/=0) k0=2
963  DO k=k0,LEN_TRIM(str); IF(INDEX('0123456789',str(k:k))==0) out=.FALSE.; END DO
964END FUNCTION is_int
965!-------------------------------------------------------------------------------
966
967!-------------------------------------------------------------------------------
968ELEMENTAL LOGICAL FUNCTION is_dec(str) RESULT(out)
969!-------------------------------------------------------------------------------
970  IMPLICIT NONE
971!-------------------------------------------------------------------------------
972! Purpose: tests if str is a valid decimal number or an integer
973!-------------------------------------------------------------------------------
974  CHARACTER(LEN=*), INTENT(IN) :: str
975  INTEGER :: ix
976!-------------------------------------------------------------------------------
977  out=LEN_TRIM(str)/=0; IF(.NOT.out) RETURN
978
979  !--- TEST IF THE NUMBER IS A PURE INTEGER
980  IF(is_int(str)) RETURN
981
982  !--- CHECK FOR POINT SEPARATOR (SHOULD HAVE A STRING INT BEFORE AND/OR AFTER
983  ix=INDEX(str,'.'); out=ix>0.AND.LEN_TRIM(str)/=1; IF(.NOT.out) RETURN
984
985  !--- CHECK SURROUNDING STRINGS ARE INTEGERS
986  out=is_int(str(1:ix-1)); IF(ix==LEN_TRIM(str)) RETURN
987  out=out.AND.is_int(str(ix+1:LEN_TRIM(str)))
988END FUNCTION is_dec
989!-------------------------------------------------------------------------------
990
991!-------------------------------------------------------------------------------
992ELEMENTAL LOGICAL FUNCTION is_num(str) RESULT(out)
993!-------------------------------------------------------------------------------
994  IMPLICIT NONE
995!-------------------------------------------------------------------------------
996! Purpose: test if str is a valid integer/decimal/scientific notation number
997!-------------------------------------------------------------------------------
998  CHARACTER(LEN=*),  INTENT(IN) :: str
999  INTEGER :: ix, ls
1000!-------------------------------------------------------------------------------
1001  out=LEN_TRIM(str)/=0; IF(.NOT.out) RETURN; ls=LEN_TRIM(str)
1002
1003  !--- TEST IF THE NUMBER IS A PURE INTEGER OR A PURE DECIMAL NUMBER
1004  IF(is_int(str)) RETURN
1005  IF(is_dec(str)) RETURN
1006
1007  !--- CHECK FOR "D" OR "E" SEPARATOR
1008  ix=MAX(INDEX(strLow(str),'d'),INDEX(strLow(str),'e')); out=ix>0
1009  IF(.NOT.out) RETURN
1010
1011  !--- CHECK SURROUNDING STRINGS ARE NUMBERS
1012  out=is_dec(str(1:ix-1)).AND.is_int(str(ix+1:ls))
1013END FUNCTION is_num
1014!-------------------------------------------------------------------------------
1015
1016!-------------------------------------------------------------------------------
1017INTEGER FUNCTION str2int(str) RESULT(out)
1018!-------------------------------------------------------------------------------
1019  IMPLICIT NONE
1020!-------------------------------------------------------------------------------
1021! Purpose: convert a string into an integer
1022!-------------------------------------------------------------------------------
1023  CHARACTER(LEN=*), INTENT(IN) :: str
1024  INTEGER :: k
1025!-------------------------------------------------------------------------------
1026  CALL err(.NOT.is_int(str),'Invalid integer "'//TRIM(str)//'"')
1027  out=0
1028  DO k=1,LEN_TRIM(str); out=10*out+INT(INDEX('0123456789',str(k:k))-1); END DO
1029END FUNCTION str2int
1030!-------------------------------------------------------------------------------
1031
1032!-------------------------------------------------------------------------------
1033REAL FUNCTION str2real(str) RESULT(out)
1034!-------------------------------------------------------------------------------
1035  IMPLICIT NONE
1036!-------------------------------------------------------------------------------
1037! Purpose: Convert a string into an real
1038!-------------------------------------------------------------------------------
1039  CHARACTER(LEN=*), INTENT(IN) :: str
1040  INTEGER :: ix, l
1041!-------------------------------------------------------------------------------
1042  CALL err(.NOT.is_num(str),'Invalid number "'//TRIM(str)//'"')
1043  ix=MAX(INDEX(strLow(str),'d'),INDEX(strLow(str),'e')); l=LEN_TRIM(str)
1044  IF(ix==0) out=dec2real(str)
1045  IF(ix/=0) out=dec2real(str(1:ix-1))*10.**str2int(str(ix+1:l))
1046END FUNCTION str2real
1047!-------------------------------------------------------------------------------
1048
1049!-------------------------------------------------------------------------------
1050ELEMENTAL REAL FUNCTION dec2real(str) RESULT(out)
1051!-------------------------------------------------------------------------------
1052  IMPLICIT NONE
1053!-------------------------------------------------------------------------------
1054! Purpose: convert a valid (already checked) string into an real
1055!-------------------------------------------------------------------------------
1056  CHARACTER(LEN=*), INTENT(IN) :: str
1057  INTEGER :: k, k0, ix
1058!-------------------------------------------------------------------------------
1059  out=0.
1060  k0=1; IF(INDEX('+-',str(1:1))/=0) k0=2; ix=INDEX(str,'.')     !--- CHECK SIGN
1061  DO k=k0,LEN_TRIM(str); IF(k==ix) CYCLE                        !--- SKIP .+-
1062    out=out*10.+REAL(INDEX('0123456789',str(k:k))-1)
1063  END DO
1064  IF(ix/=0) out=out*10.**(ix-LEN_TRIM(str))                     !--- COMA
1065  IF(k0==2) out=-out                                            !--- SIGN
1066END FUNCTION dec2real
1067!-------------------------------------------------------------------------------
1068
1069!-------------------------------------------------------------------------------
1070
1071ELEMENTAL CHARACTER(LEN=256) FUNCTION strLow(str) RESULT(out)
1072
1073!-------------------------------------------------------------------------------
1074  IMPLICIT NONE
1075!-------------------------------------------------------------------------------
1076! Purpose: Lower case conversion.
1077!-------------------------------------------------------------------------------
1078! Parameters:
1079  CHARACTER(LEN=*), INTENT(IN) :: str
1080!-------------------------------------------------------------------------------
1081! Local variable:
1082  INTEGER :: k, ix
1083!-------------------------------------------------------------------------------
1084  out=str
1085  DO k=1,LEN(str)
1086    ix=IACHAR(str(k:k)); IF(64<ix.AND.ix<91) out(k:k)=ACHAR(ix+32)
1087  END DO
1088
1089END FUNCTION strLow
1090
1091!-------------------------------------------------------------------------------
1092
1093END PROGRAM make_sso
1094
1095!-------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.