Subversion Repositories lagranto.ecmwf

Rev

Details | Last modification | View Log | RSS feed

Rev Author Line No. Line
13 michaesp 1
      subroutine zlev(z,t,q,zb,sp,ie,je,ke,ak,bk)
2
c     ===========================================
3
 
4
c     argument declaration
5
      integer  ie,je,ke,its,iqs
6
      real     z(ie,je,ke),t(ie,je,ke),getps
7
      real     q(ie,je,ke),zb(ie,je),sp(ie,je)
8
      real     ak(ke),bk(ke)
9
 
10
c     variable declaration
11
      integer  i,j,k
12
      real     pu, po, zo, zu, tv
13
      real     rdg,tzero,eps
14
      data     rdg,tzero,eps /29.271, 273.15, 0.6078/
15
 
16
c     t and q are staggered (stag(3)=-0.5; its=int(2*stag(3))
17
      its=-1
18
      iqs=-1
19
 
20
c     computation of height of main levels
21
      do i=1,ie
22
        do j=1,je
23
          zu = zb(i,j)
24
          pu = sp(i,j)
25
          tv = (t(i,j,1)+tzero)*(1.+eps*q(i,j,1))
26
          po = ak(1)+bk(1)*sp(i,j)
27
          if (po.le.0) po=0.5*pu
28
          zo = zu + rdg*tv*alog(pu/po)
29
          z(i,j,1) = zo
30
          pu = po
31
          zu = zo
32
          do k=2,ke
33
            tv = (0.5*(t(i,j,k)+t(i,j,k-1-its))+tzero)*
34
     >           (1.+eps*0.5*(q(i,j,k)+q(i,j,k-1-iqs)))
35
            po = ak(k)+bk(k)*sp(i,j)
36
            if (po.le.0) po=0.5*pu
37
            zo = zu + rdg*tv*alog(pu/po)
38
            z(i,j,k) = zo
39
            pu = po
40
            zu = zo
41
          enddo
42
        enddo
43
      enddo
44
      end
45
 
46
      subroutine flightlev(fl,p,ie,je,ke)
47
c     ===================================
48
 
49
c     Convert an array of pressures (hPa) into an array of altitudes (m)
50
c     or vice versa using the ICAO standard atmosphere.
51
c
52
c     See http://www.pdas.com/coesa.htm
53
c
54
c     The 7 layers of the US standard atmosphere are:
55
c
56
c     h1   h2     dT/dh    h1,h2 geopotential alt in km
57
c      0   11     -6.5     dT/dh in K/km
58
c     11   20      0.0
59
c     20   32      1.0
60
c     32   47      2.8
61
c     47   51      0.0
62
c     51   71     -2.8
63
c     71   84.852 -2.0
64
 
65
c     based upon IDL routine from Dominik Brunner
66
c     converted to f90: Jan 2002  Heini Wernli
67
 
68
      implicit none
69
 
70
c     argument declaration
71
      integer   ie,je,ke
72
      real      fl(ie,je,ke),p(ie,je,ke)
73
 
74
c     variable declaration
75
      integer   i,j,k,n,il
76
      real      limits(0:7),lrs(7)
77
      integer   iszero(7)
78
      real      pb(0:7),tb(0:7)
79
      real      g,r,gmr,feet
80
      data      g,r,feet /9.80665, 287.053, 0.3048/
81
 
82
c     define layer boundaries in m, lapsre rates (9 means 0) and isothermal flag
83
      limits(0)=0.
84
      limits(1)=11000.
85
      limits(2)=20000.
86
      limits(3)=32000.
87
      limits(4)=47000.
88
      limits(5)=51000.
89
      limits(6)=71000.
90
      limits(7)=84852.
91
 
92
      lrs(0)=-6.5/1000.
93
      lrs(1)=9.0/1000.
94
      lrs(2)=1.0/1000.
95
      lrs(3)=2.8/1000.
96
      lrs(4)=9.0/1000.
97
      lrs(5)=-2.8/1000.
98
      lrs(6)=-2.0/1000.
99
 
100
      iszero(0)=0
101
      iszero(1)=1
102
      iszero(2)=0
103
      iszero(3)=0
104
      iszero(4)=1
105
      iszero(5)=0
106
      iszero(6)=0
107
 
108
      gmr=g/r         ! Hydrostatic constant
109
 
110
c     calculate pressures at layer boundaries
111
      pB(0)=1013.25 ! pressure at surface
112
      TB(0)=288.15  ! Temperature at surface
113
 
114
c     loop over layers and get pressures and temperatures at level tops
115
      do i=0,6
116
        TB(i+1)=TB(i)+(1-iszero(i))*lrs(i)*(limits(i+1)-limits(i))
117
        pB(i+1)=(1-iszero(i))*pB(i)*exp(alog(TB(i)/TB(i+1))*gmr/lrs(i))+
118
     >          iszero(i)*PB(i)*exp(-gmr*(limits(i+1)-limits(i))/TB(i))
119
      enddo
120
 
121
c     now calculate which layer each value belongs to
122
c     and calculate the corresponding altitudes
123
      do i=1,ie
124
      do j=1,je
125
      do k=1,ke
126
        do n=0,6
127
          if ((pb(n).ge.p(i,j,k)).and.(pb(n+1).le.p(i,j,k))) then
128
            il=n
129
            goto 100
130
          endif
131
        enddo
132
  100   continue
133
        fl(i,j,k)=iszero(il)*(-ALOG(p(i,j,k)/pB(il))*TB(il)/gmr+
134
     >            limits(il))+(1-iszero(il))*((TB(il)/((p(i,j,k)/
135
     >            pB(il))**(lrs(il)/gmr))-TB(il))/lrs(il)+limits(il))
136
        fl(i,j,k)=fl(i,j,k)/feet/100.
137
      enddo
138
      enddo
139
      enddo
140
 
141
      end
142
 
143
      subroutine pottemp(pt,t,sp,ie,je,ke,ak,bk)
144
c     ==========================================
145
 
146
c     argument declaration
147
      integer   ie,je,ke
148
      real      pt(ie,je,ke),t(ie,je,ke),sp(ie,je),
149
     >     	ak(ke),bk(ke)
150
 
151
c     variable declaration
152
      integer   i,j,k
153
      real      rdcp,tzero
154
      data      rdcp,tzero /0.286,273.15/
155
 
156
c     statement-functions for the computation of pressure
157
      real      pp,psrf
158
      integer   is
159
      pp(is)=ak(is)+bk(is)*psrf
160
 
161
c     computation of potential temperature
162
      do i=1,ie
163
      do j=1,je
164
        psrf=sp(i,j)
165
        do k=1,ke
166
c     distinction of temperature in K and deg C
167
          if (t(i,j,k).lt.100.) then
168
            pt(i,j,k)=(t(i,j,k)+tzero)*( (1000./pp(k))**rdcp )
169
          else
170
            pt(i,j,k)=t(i,j,k)*( (1000./pp(k))**rdcp )
171
          endif
172
        enddo
173
      enddo
174
      enddo
175
      end
176
 
177
      subroutine tnat(tn,t,p,ie,je,ke)
178
c     ================================
179
 
180
c     argument declaration
181
      integer   ie,je,ke
182
      real      tn(ie,je,ke),t(ie,je,ke),p(ie,je,ke)
183
 
184
c     variable declaration
185
      integer   i,j,k
186
      real      tzero
187
      data      tzero /273.15/
188
 
189
      real      mb_to_mm
190
      real      a1,a2,a3,a4,a5,c0,c1,c2,d
191
      real      chiH2O,chiHNO3,pH2O,pHNO3,log_H2O,log_NAT
192
 
193
      mb_to_mm=760./1000.
194
      a1 = -2.7836
195
      a2 = -0.00088
196
      a3 = 89.7674
197
      a4 = -26242.0
198
      a5 = 0.021135
199
 
200
      chiH2O  = 5.e-6
201
      chiHNO3 = 5.e-9
202
 
203
      do i=1,ie
204
      do j=1,je
205
      do k=1,ke
206
        pH2O  = p(i,j,k)*chiH2O
207
        pHNO3 = p(i,j,k)*chiHNO3
208
        log_H2O = log(pH2O*mb_to_mm)
209
        log_NAT = log(pHNO3*mb_to_mm)
210
        c0 = log_NAT-log_H2O*a1-a3
211
        c1 = a2*log_H2O+a5
212
        c2 = a4
213
        d = c0*c0-4.0*c1*c2
214
        tn(i,j,k)=t(i,j,k)+tzero-(c0+sqrt(d))/(2.*c1)
215
      enddo
216
      enddo
217
      enddo
218
      end
219
 
220
      subroutine read_tice(tice_arr)
221
c     ==============================
222
 
223
      integer	n
224
      real	tice_arr(100)
225
 
226
      open(10,
227
     >file='/home/henry/prog/tnat_tice/t_ice_from_10_every_2_hPa')
228
      do n=1,71
229
        read(10,*)tice_arr(n)
230
      enddo
231
      close(10)
232
      end
233
 
234
      subroutine tice(ti,t,p,tice_arr,ie,je,ke)
235
c     =========================================
236
 
237
c     argument declaration
238
      integer   ie,je,ke
239
      real      ti(ie,je,ke),t(ie,je,ke),p(ie,je,ke)
240
 
241
c     variable declaration
242
      integer   i,j,k
243
      real      tzero
244
      data      tzero /273.15/
245
 
246
      real	tice_arr(100)
247
      real      p0,p1,pinc,rind
248
      integer   ind,np
249
 
250
c     define the lowest p-value and the p-increment for the table
251
      p0=10.
252
      pinc=2.
253
      np=71
254
      p1=p0+real(np-1)*pinc
255
 
256
      do i=1,ie
257
      do j=1,je
258
      do k=1,ke
259
 
260
        if (p(i,j,k).lt.p0) then
261
          print*,'*** error: table is not prepared for p < p0 ***'
262
          print*,ie,je,ke,i,j,k,p(i,j,k),p0
263
          call exit(1)
264
        endif
265
 
266
        if (p(i,j,k).gt.p1) then
267
          print*,'*** error: table is not prepared for p > p1 ***'
268
          call exit(1)
269
        endif
270
 
271
c       calculate the real index of the given pressure level
272
        rind=(p(i,j,k)-p0)/pinc+1.
273
        ind=int(rind)
274
        rind=rind-real(ind)
275
 
276
c       interpolate tice
277
        ti(i,j,k)=t(i,j,k)+tzero-
278
     >           (tice_arr(ind)+rind*(tice_arr(ind+1)-tice_arr(ind)))
279
 
280
      enddo
281
      enddo
282
      enddo
283
      end
284
 
285
      subroutine temp(t,pt,p,ie,je,ke,ak,bk)
286
c     ======================================
287
 
288
c     argument declaration
289
      integer   ie,je,ke
290
      real      pt(ie,je,ke),t(ie,je,ke),p(ie,je,ke),
291
     >     	ak(ke),bk(ke)
292
 
293
c     variable declaration
294
      integer   i,j,k
295
      real      rdcp,tzero
296
      data      rdcp,tzero /0.286,273.15/
297
 
298
c     computation of potential temperature
299
      do i=1,ie
300
         do j=1,je
301
            do k=1,ke
302
               t(i,j,k)=pt(i,j,k)*( (p(i,j,k)/1000.)**rdcp )
303
            enddo
304
         enddo
305
      enddo
306
      end
307
 
308
      subroutine laitpv(lpv,pv,th,ie,je,ke)
309
c     =====================================
310
 
311
c     argument declaration
312
      integer   ie,je,ke
313
      real      lpv(ie,je,ke),pv(ie,je,ke),th(ie,je,ke)
314
 
315
c     variable declaration
316
      integer   i,j,k
317
      real      rdcp,tzero
318
      data      rdcp,tzero /0.286,273.15/
319
 
320
c     computation of Lait PV
321
      do i=1,ie
322
      do j=1,je
323
      do k=1,ke
324
        lpv(i,j,k)=pv(i,j,k)*((th(i,j,k)/420.)**(-9./2.))
325
      enddo
326
      enddo
327
      enddo
328
 
329
      end
330
 
331
      subroutine relhum(rh,q,t,sp,ie,je,ke,ak,bk)
332
c     ===========================================
333
 
334
c     argument declaration
335
      integer   ie,je,ke
336
      real      rh(ie,je,ke),t(ie,je,ke),q(ie,je,ke),
337
     >     sp(ie,je),ak(ke),bk(ke)
338
 
339
c     variable declaration
340
      integer   i,j,k
341
      real      rdcp,tzero
342
      real      b1,b2w,b3,b4w,r,rd,gqd,ge
343
      data      rdcp,tzero /0.286,273.15/
344
      data      b1,b2w,b3,b4w,r,rd /6.1078, 17.2693882, 273.16, 35.86,
345
     &     287.05, 461.51/
346
 
347
c     statement-functions for the computation of pressure
348
      real      pp,psrf
349
      integer   is
350
      pp(is)=ak(is)+bk(is)*psrf
351
 
352
      do i=1,ie
353
         do j=1,je
354
            psrf=sp(i,j)
355
            do k=1,ke
356
               ge = b1*exp(b2w*t(i,j,k)/(t(i,j,k)+b3-b4w))
357
               gqd= r/rd*ge/(pp(k)-(1.-r/rd)*ge)
358
               rh(i,j,k)=100.*q(i,j,k)/gqd
359
            enddo
360
         enddo
361
      enddo
362
      end
363
 
364
      subroutine relhum_i(rh,q,p,ie,je,ke,ak,mdv)
365
c     ===========================================
366
 
367
c     argument declaration
368
      integer   ie,je,ke
369
      real      rh(ie,je,ke),p(ie,je,ke),q(ie,je,ke),
370
     >     	ak(ke),mdv
371
 
372
c     variable declaration
373
      integer   i,j,k
374
      real      p0,kappa
375
      data      p0,kappa /1000.,0.286/
376
      real      rdcp,tzero
377
      data      rdcp,tzero /0.286,273.15/
378
      real      b1,b2w,b3,b4w,r,rd,gqd,ge
379
      data      b1,b2w,b3,b4w,r,rd /6.1078, 17.2693882, 273.16, 35.86,
380
     &     287.05, 461.51/
381
 
382
      do i=1,ie
383
      do j=1,je
384
      do k=1,ke
385
        if (p(i,j,k).eq.mdv) then
386
          rh(i,j,k)=mdv
387
        else
388
          tt=ak(k)*((p(i,j,k)/p0)**kappa)-tzero
389
          ge = b1*exp(b2w*tt/(tt+b3-b4w))
390
          gqd= r/rd*ge/(p(i,j,k)-(1.-r/rd)*ge)
391
          rh(i,j,k)=100.*q(i,j,k)/gqd
392
        endif
393
      enddo
394
      enddo
395
      enddo
396
      end
397
 
398
      subroutine gradth(gth,th,sp,levtyp,cl,ie,je,ke,ak,bk,vmin,vmax)
399
C     ===============================================================
400
 
401
c     argument declaration
402
      integer   ie,je,ke,levtyp
403
      real      gth(ie,je,ke),th(ie,je,ke),sp(ie,je),cl(ie,je)
404
      real      ak(ke),bk(ke),vmin(4),vmax(4)
405
 
406
c     variable declaration
407
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dspdx,dspdy
408
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dthdx,dthdy
409
      integer   i,j,k,ind,ind2,stat
410
 
411
      allocate(dthdx(ie,je,ke),STAT=stat)
412
      if (stat.ne.0) print*,'gradth: error allocating dthdx'
413
      allocate(dthdy(ie,je,ke),STAT=stat)
414
      if (stat.ne.0) print*,'gradth: error allocating dthdy'
415
      allocate(dspdx(ie,je),STAT=stat)
416
      if (stat.ne.0) print*,'gradth: error allocating dspdx'
417
      allocate(dspdy(ie,je),STAT=stat)
418
      if (stat.ne.0) print*,'gradth: error allocating dspdy'
419
 
420
      if (levtyp.eq.1) then
421
         call ddh2(th,dthdx,cl,'X',ie,je,1,vmin,vmax)
422
         call ddh2(th,dthdy,cl,'Y',ie,je,1,vmin,vmax)
423
      else if (levtyp.eq.3) then
424
         call ddh2(sp,dspdx,cl,'X',ie,je,1,vmin,vmax)
425
         call ddh2(sp,dspdy,cl,'Y',ie,je,1,vmin,vmax)
426
         call ddh3(th,dthdx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
427
         call ddh3(th,dthdy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
428
      endif
429
 
430
      gth=sqrt(dthdx**2.+dthdy**2.)
431
 
432
      IF (ALLOCATED(dspdx)) DEALLOCATE(dspdx)
433
      IF (ALLOCATED(dspdy)) DEALLOCATE(dspdy)
434
      IF (ALLOCATED(dthdx)) DEALLOCATE(dthdx)
435
      IF (ALLOCATED(dthdy)) DEALLOCATE(dthdy)
436
 
437
      end
438
 
439
       subroutine calc_gradpv(gpv,pv,cl,ie,je,ke,ak,bk,vmin,vmax)
440
C      ===============================================================
441
 
442
c      argument declaration
443
       integer   ie,je,ke,levtyp
444
       real      gpv(ie,je,ke),pv(ie,je,ke),sp(ie,je),cl(ie,je)
445
       real      ak(ke),bk(ke),vmin(4),vmax(4)
446
 
447
c      variable declaration
448
       REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dpvdx,dpvdy
449
       integer   i,j,k,ind,ind2,stat
450
 
451
       allocate(dpvdx(ie,je,ke),STAT=stat)
452
       if (stat.ne.0) print*,'gradpv: error allocating dpvdx'
453
       allocate(dpvdy(ie,je,ke),STAT=stat)
454
       if (stat.ne.0) print*,'gradpv: error allocating dpvdy'
455
 
456
       call ddh2(pv,dpvdx,cl,'X',ie,je,1,vmin,vmax)
457
       call ddh2(pv,dpvdy,cl,'Y',ie,je,1,vmin,vmax)
458
 
459
       gpv=1e6*sqrt(dpvdx**2.+dpvdy**2.)
460
 
461
       IF (ALLOCATED(dpvdx)) DEALLOCATE(dpvdx)
462
       IF (ALLOCATED(dpvdy)) DEALLOCATE(dpvdy)
463
 
464
       end
465
 
466
      subroutine wetbpt(thw,the,sp,ie,je,ke,ak,bk)
467
c     ============================================
468
 
469
c     argument declaration
470
      integer   ie,je,ke
471
      real      thw(ie,je,ke),the(ie,je,ke),
472
     >     sp(ie,je),ak(ke),bk(ke)
473
 
474
c     variable declaration
475
      real      tsa
476
      integer   i,j,k
477
 
478
      do k=1,ke
479
         do j=1,je
480
            do i=1,ie
481
               thw(i,j,k)=tsa(the(i,j,k)-273.15,1000.)+273.15
482
            enddo
483
         enddo
484
      enddo
485
      end
486
 
487
      real function tsa(os,p)
488
c     =======================
489
 
490
C     This function returns the temperature tsa (celsius) on a
491
c     saturation adiabat at pressure p (millibars). os is the equivalent
492
c     potential temperature of the parcel (celsius). sign(a,b) replaces
493
c     the algebraic sign of a with that of b.
494
C     b is an empirical constant approximately equal to 0.001 of the
495
c     latent heat of vaporization for water divided by the specific
496
c     heat at constant pressure for dry air.
497
 
498
      real      a,b,os,p,tq,d,tqk,x,w
499
      integer   i
500
 
501
      data b/2.6518986/
502
      a=os+273.16
503
 
504
C     tq is the first guess for tsa
505
 
506
      tq=253.16
507
 
508
C     d is an initial value used in the iteration below
509
 
510
      d=120.
511
 
512
C     Iterate to obtain sufficient accuracy....see table 1, p.8
513
C     of Stipanuk (1973) for equation used in iteration
514
      do 1 i=1,12
515
         tqk=tq-273.16
516
         d=d/2.
517
         x=a*exp(-b*w(tqk,p)/tq)-tq*((1000./p)**.286)
518
         if (abs(x).lt.1e-7) go to 2
519
         tq=tq+sign(d,x)
520
 1    continue
521
 2    tsa=tq-273.16
522
      end
523
 
524
      real function w(t,p)
525
c     ====================
526
 
527
C     This function returns the mixing ratio (grams of water vapor per
528
C     kilogram of dry air) given the dew point (celsius) and pressure
529
C     (millibars). If the temperature is input instead of the
530
C     dew point, then saturation mixing ratio (same units) is returned.
531
C     The formula is found in most meteorological texts.
532
 
533
      real      t,p,tkel,x,esat
534
 
535
      tkel=t+273.16
536
      x=esat(tkel)              ! our function esat requires t in Kelvin
537
      w=622.*x/(p-x)
538
      end
539
 
540
 
541
      real function esat(t)
542
c     =====================
543
 
544
C     This function returns the saturation vapor pressure over water
545
c     (mb) given the temperature (Kelvin).
546
C     The algorithm is due to Nordquist, W. S. ,1973: "Numerical
547
C     Approximations of Selected Meteorological Parameters for Cloud
548
C     Physics Problems" ECOM-5475, Atmospheric Sciences Laboratory,
549
c     U. S. Army Electronics Command, White Sands Missile Range,
550
c     New Mexico 88002.
551
 
552
      real p1,p2,c1,t
553
 
554
      p1=11.344-0.0303998*t
555
      p2=3.49149-1302.8844/t
556
      c1=23.832241-5.02808*log10(t)
557
      esat=10.**(c1-1.3816e-7*10.**p1+8.1328e-3*10.**p2-2949.076/t)
558
      end
559
 
560
      subroutine equpot(ap,t,q,sp,ie,je,ke,ak,bk)
561
c     ===========================================
562
 
563
c     calculate equivalent potential temperature
564
 
565
c     argument declaration
566
      integer  ie,je,ke
567
      real     ap(ie,je,ke),t(ie,je,ke),sp(ie,je)
568
      real     q(ie,je,ke),ak(ke),bk(ke)
569
 
570
c     variable declaration
571
      integer  i,j,k
572
      real     rdcp,tzero
573
      data     rdcp,tzero /0.286,273.15/
574
 
575
c     statement-functions for the computation of pressure
576
      real      pp,psrf
577
      integer   is
578
      pp(is)=ak(is)+bk(is)*psrf
579
 
580
c     computation of potential temperature
581
      do i=1,ie
582
         do j=1,je
583
            psrf=sp(i,j)
584
            do k=1,ke
585
               ap(i,j,k) = (t(i,j,k)+tzero)*(1000./pp(k))
586
     +              **(0.2854*(1.0-0.28*q(i,j,k)))*exp(
587
     +              (3.376/(2840.0/(3.5*alog(t(i,j,k)+tzero)-alog(
588
     +              100.*pp(k)*max(1.0E-10,q(i,j,k))/(0.622+0.378*
589
     +              q(i,j,k)))-0.1998)+55.0)-0.00254)*1.0E3*
590
     +              max(1.0E-10,q(i,j,k))*(1.0+0.81*q(i,j,k)))
591
            enddo
592
         enddo
593
      enddo
594
      end
595
 
596
      subroutine diabheat(dhr,t,w,rh,sp,ie,je,ke,ak,bk)
597
c     =================================================
598
 
599
c     argument declaration
600
      integer   ie,je,ke
601
      real      dhr(ie,je,ke),t(ie,je,ke),w(ie,je,ke),
602
     &     rh(ie,je,ke),sp(ie,je),ak(ke),bk(ke)
603
 
604
c     variable declaration
605
      integer   i,j,k
606
      real      p0,kappa,tzero
607
      data      p0,kappa,tzero /1000.,0.286,273.15/
608
      real      blog10,cp,r,lw,eps
609
      data      blog10,cp,r,lw,eps /.08006,1004.,287.,2.5e+6,0.622/
610
      real      esat,c,tt
611
 
612
c     statement-functions for the computation of pressure
613
      real      pp,psrf
614
      integer   is
615
      pp(is)=ak(is)+bk(is)*psrf
616
 
617
c     computation of diabatic heating rate
618
      do i=1,ie
619
         do j=1,je
620
            psrf=sp(i,j)
621
            do k=1,ke
622
               if (rh(i,j,k).lt.80.) then ! only moist air of interest
623
                  dhr(i,j,k)=0. ! cond. heating rate set to zero
624
               else if (w(i,j,k).gt.0.) then ! cond. heating only
625
                                             ! for ascent
626
                  dhr(i,j,k)=0.
627
               else
628
                  tt=t(i,j,k)*((pp(k)/p0)**kappa) ! temp. from pot.temp.
629
                  c=lw/cp*eps*blog10*esat(tt)/pp(k)
630
                  dhr(i,j,k)=21600.* ! in units K per 6 hours
631
     >                 (1.-exp(.2*(80.-rh(i,j,k)))) ! weighting fun.
632
                                                    ! for 80<RH<100
633
     >                 *(-c*kappa*t(i,j,k)*w(i,j,k)/(100.*pp(k)))/(1.+c)
634
               endif
635
            enddo
636
         enddo
637
      enddo
638
      end
639
 
640
      subroutine diabheat2(dhr,t,w,rh,sp,ie,je,ke,ak,bk)
641
c     ==================================================
642
 
643
c     argument declaration
644
      integer   ie,je,ke
645
      real      dhr(ie,je,ke),t(ie,je,ke),w(ie,je,ke),
646
     &     rh(ie,je,ke),sp(ie,je),ak(ke),bk(ke)
647
 
648
c     variable declaration
649
      integer   i,j,k
650
      real      p0,kappa,tzero
651
      data      p0,kappa,tzero /1000.,0.286,273.15/
652
      real      blog10,cp,r,lw,eps
653
      data      blog10,cp,r,lw,eps /.08006,1004.,287.,2.5e+6,0.622/
654
      real      esat,c,tt
655
 
656
c     statement-functions for the computation of pressure
657
      real      pp,psrf
658
      integer   is
659
      pp(is)=ak(is)+bk(is)*psrf
660
 
661
c     computation of diabatic heating rate
662
      do i=1,ie
663
         do j=1,je
664
            psrf=sp(i,j)
665
            do k=1,ke
666
               if (rh(i,j,k).lt.80.) then ! only moist air of interest
667
                  dhr(i,j,k)=0. ! cond. heating rate set to zero
668
               else if (w(i,j,k).gt.0.) then ! cond. heating only
669
                                             ! for ascent
670
                  dhr(i,j,k)=0.
671
               else
672
                  c=lw/cp*eps*blog10*esat((t(i,j,k)+273.15))/pp(k)
673
                  tt=(t(i,j,k)+273.15)*((p0/pp(k))**kappa) ! pot.temp from temp
674
                  dhr(i,j,k)=21600.* ! in units K per 6 hours
675
     >                 (1.-exp(.2*(80.-rh(i,j,k)))) ! weighting fun.
676
                                                    ! for 80<RH<100
677
     >                 *(-c*kappa*tt*w(i,j,k)/(100.*pp(k)))/(1.+c)
678
               endif
679
            enddo
680
         enddo
681
      enddo
682
      end
683
 
684
      subroutine diabpvr(dpvr,uu,vv,dhr,sp,cl,f,ie,je,ke,ak,bk,
685
     >     vmin,vmax)
686
C     =========================================================
687
 
688
c     argument declaration
689
      integer   ie,je,ke
690
      real,intent(OUT) :: dpvr(ie,je,ke)
691
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),
692
     >     dhr(ie,je,ke),sp(ie,je),cl(ie,je),f(ie,je)
693
      real,intent(IN)  :: ak(ke),bk(ke),vmin(4),vmax(4)
694
 
695
c     variable declaration
696
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dspdx,dspdy
697
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dhrdp,dudp,dvdp,dvdx
698
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dhrdx,dudy,dhrdy
699
      integer   k,stat
700
 
701
      allocate(dspdx(ie,je),STAT=stat)
702
      if (stat.ne.0) print*,'diabpvr: error allocating dspdx'
703
      allocate(dspdy(ie,je),STAT=stat)
704
      if (stat.ne.0) print*,'diabpvr: error allocating dspdy'
705
 
706
      allocate(dhrdp(ie,je,ke),STAT=stat)
707
      if (stat.ne.0) print*,'diabpvr: error allocating dhrdp'
708
      allocate(dudp(ie,je,ke),STAT=stat)
709
      if (stat.ne.0) print*,'diabpvr: error allocating dudp'
710
      allocate(dvdp(ie,je,ke),STAT=stat)
711
      if (stat.ne.0) print*,'diabpvr: error allocating dvdp'
712
      allocate(dvdx(ie,je,ke),STAT=stat)
713
      if (stat.ne.0) print*,'diabpvr: error allocating dvdx'
714
      allocate(dhrdx(ie,je,ke),STAT=stat)
715
      if (stat.ne.0) print*,'diabpvr: error allocating dhrdx'
716
      allocate(dudy(ie,je,ke),STAT=stat)
717
      if (stat.ne.0) print*,'diabpvr: error allocating dudy'
718
      allocate(dhrdy(ie,je,ke),STAT=stat)
719
      if (stat.ne.0) print*,'diabpvr: error allocating dhrdy'
720
 
721
      call ddp(dhr,dhrdp,sp,ie,je,ke,ak,bk)
722
      call ddp(uu,dudp,sp,ie,je,ke,ak,bk)
723
      call ddp(vv,dvdp,sp,ie,je,ke,ak,bk)
724
      call ddh2(sp,dspdx,cl,'X',ie,je,1,vmin,vmax)
725
      call ddh2(sp,dspdy,cl,'Y',ie,je,1,vmin,vmax)
726
      call ddh3(dhr,dhrdx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
727
      call ddh3(vv,dvdx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
728
      call ddh3(dhr,dhrdy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
729
      call ddh3(uu,dudy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
730
 
731
      do k=1,ke
732
         dpvr(1:ie,1:je,k)=-1.E6*9.80665*(
733
     >        -dvdp(1:ie,1:je,k)*dhrdx(1:ie,1:je,k)
734
     >        +dudp(1:ie,1:je,k)*dhrdy(1:ie,1:je,k)
735
     >        +(-dudy(1:ie,1:je,k)+dvdx(1:ie,1:je,k)
736
     >        +f(1:ie,1:je))*dhrdp(1:ie,1:je,k))
737
      enddo
738
 
739
      IF (ALLOCATED(dspdx)) DEALLOCATE(dspdx)
740
      IF (ALLOCATED(dspdy)) DEALLOCATE(dspdy)
741
      IF (ALLOCATED(dhrdp)) DEALLOCATE(dhrdp)
742
      IF (ALLOCATED(dudp)) DEALLOCATE(dudp)
743
      IF (ALLOCATED(dvdp)) DEALLOCATE(dvdp)
744
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
745
      IF (ALLOCATED(dhrdx)) DEALLOCATE(dhrdx)
746
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
747
      IF (ALLOCATED(dhrdy)) DEALLOCATE(dhrdy)
748
 
749
      end
750
 
751
      subroutine vort(var,uu,vv,sp,cl,ie,je,ke,ak,bk,vmin,vmax)
752
C     =========================================================
753
C     calculate vorticity VORT=D[V:X]-D[U*COS:Y]/COS
754
 
755
c     argument declaration
756
      integer   ie,je,ke
757
      real,intent(OUT) :: var(ie,je,ke)
758
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),
759
     >     sp(ie,je),cl(ie,je)
760
      real      ak(ke),bk(ke),vmin(4),vmax(4)
761
 
762
c     variable declaration
763
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dspdx,dspdy
764
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dvdx,dudy,uucl
765
      integer   stat
766
      real	pole
767
 
768
      allocate(dspdx(ie,je),STAT=stat)
769
      if (stat.ne.0) print*,'vort: error allocating dspdx(ie,je)'
770
      allocate(dspdy(ie,je),STAT=stat)
771
      if (stat.ne.0) print*,'vort: error allocating dspdy(ie,je)'
772
      allocate(dvdx(ie,je,ke),STAT=stat)
773
      if (stat.ne.0) print*,'vort: error allocating dvdx'
774
      allocate(uucl(ie,je,ke),STAT=stat)
775
      if (stat.ne.0) print*,'vort: error allocating uucl'
776
      allocate(dudy(ie,je,ke),STAT=stat)
777
      if (stat.ne.0) print*,'vort: error allocating dudy'
778
 
779
      call ddh2(sp,dspdx,cl,'X',ie,je,1,vmin,vmax)
780
      call ddh2(sp,dspdy,cl,'Y',ie,je,1,vmin,vmax)
781
      call ddh3(vv,dvdx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
782
      do k=1,ke
783
         uucl(1:ie,1:je,k)=uu(1:ie,1:je,k)*cl(1:ie,1:je)
784
      enddo
785
      call ddh3(uucl,dudy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
786
 
787
      do k=1,ke
788
         var(1:ie,1:je,k)=1.E4*(dvdx(1:ie,1:je,k)-
789
     >                          dudy(1:ie,1:je,k)/cl(1:ie,1:je))
790
*        var(1:ie,1:je,k)=1.E4*(dvdx(1:ie,1:je,k)-
791
*    >                          dudy(1:ie,1:je,k))
792
      enddo
793
 
794
c     interpolate poles from neighbouring points on every level
795
      if (vmax(2).eq.90.) then
796
        do k=1,ke
797
          pole=0.
798
          do i=1,ie
799
            pole=pole+var(i,je-1,k)
800
          enddo
801
          pole=pole/real(ie)
802
          do i=1,ie
803
            var(i,je,k)=pole
804
          enddo
805
        enddo
806
      endif
807
 
808
      if (vmin(2).eq.-90.) then
809
        do k=1,ke
810
          pole=0.
811
          do i=1,ie
812
            pole=pole+var(i,2,k)
813
          enddo
814
          pole=pole/real(ie)
815
          do i=1,ie
816
            var(i,1,k)=pole
817
          enddo
818
        enddo
819
      endif
820
 
821
      IF (ALLOCATED(dspdx)) DEALLOCATE(dspdx)
822
      IF (ALLOCATED(dspdy)) DEALLOCATE(dspdy)
823
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
824
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
825
      IF (ALLOCATED(uucl)) DEALLOCATE(uucl)
826
 
827
      end
828
 
829
      subroutine avort(var,uu,vv,sp,cl,f,ie,je,ke,ak,bk,vmin,vmax)
830
C     ============================================================
831
C     calculate absolute vorticity AVO=F+D[V:X]-D[U*COS:Y]/COS
832
 
833
c     argument declaration
834
      integer   ie,je,ke
835
      real,intent(OUT) :: var(ie,je,ke)
836
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),
837
     >                    sp(ie,je),cl(ie,je),f(ie,je)
838
      real      ak(ke),bk(ke),vmin(4),vmax(4)
839
 
840
c     variable declaration
841
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dspdx,dspdy
842
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dvdx,dudy,uucl
843
      integer   stat
844
      real      pole
845
 
846
      allocate(dspdx(ie,je),STAT=stat)
847
      if (stat.ne.0) print*,'vort: error allocating dspdx(ie,je)'
848
      allocate(dspdy(ie,je),STAT=stat)
849
      if (stat.ne.0) print*,'vort: error allocating dspdy(ie,je)'
850
      allocate(dvdx(ie,je,ke),STAT=stat)
851
      if (stat.ne.0) print*,'vort: error allocating dvdx'
852
      allocate(uucl(ie,je,ke),STAT=stat)
853
      if (stat.ne.0) print*,'vort: error allocating uucl'
854
      allocate(dudy(ie,je,ke),STAT=stat)
855
      if (stat.ne.0) print*,'vort: error allocating dudy'
856
 
857
      call ddh2(sp,dspdx,cl,'X',ie,je,1,vmin,vmax)
858
      call ddh2(sp,dspdy,cl,'Y',ie,je,1,vmin,vmax)
859
      call ddh3(vv,dvdx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
860
      do k=1,ke
861
         uucl(1:ie,1:je,k)=uu(1:ie,1:je,k)*cl(1:ie,1:je)
862
      enddo
863
      call ddh3(uucl,dudy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
864
 
865
      do k=1,ke
866
         var(1:ie,1:je,k)=1.E4*(f(1:ie,1:je)+dvdx(1:ie,1:je,k)-
867
     >                          dudy(1:ie,1:je,k)/cl(1:ie,1:je))
868
      enddo
869
 
870
c     interpolate poles from neighbouring points on every level
871
      if (vmax(2).eq.90.) then
872
        do k=1,ke
873
          pole=0.
874
          do i=1,ie
875
            pole=pole+var(i,je-1,k)
876
          enddo
877
          pole=pole/real(ie)
878
          do i=1,ie
879
            var(i,je,k)=pole
880
          enddo
881
        enddo
882
      endif
883
 
884
      if (vmin(2).eq.-90.) then
885
        do k=1,ke
886
          pole=0.
887
          do i=1,ie
888
            pole=pole+var(i,2,k)
889
          enddo
890
          pole=pole/real(ie)
891
          do i=1,ie
892
            var(i,1,k)=pole
893
          enddo
894
        enddo
895
      endif
896
 
897
      IF (ALLOCATED(dspdx)) DEALLOCATE(dspdx)
898
      IF (ALLOCATED(dspdy)) DEALLOCATE(dspdy)
899
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
900
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
901
      IF (ALLOCATED(uucl)) DEALLOCATE(uucl)
902
 
903
      end
904
 
905
      subroutine vort_i(var,uu,vv,cl,ie,je,ke,vmin,vmax,mdv)
906
C     ======================================================
907
C     calculate vorticity VORT=1e4*D[V:X]-D[U*COS:Y]/COS on isentropic
908
C     levels 
909
C     with misdat (mdv) treatment !
910
C     mark liniger 990501
911
 
912
c     argument declaration
913
      integer   ie,je,ke
914
      real,intent(OUT) :: var(ie,je,ke)
915
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),
916
     >     cl(ie,je)
917
      real      vmin(4),vmax(4),mdv
918
 
919
c     variable declaration
920
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dvdx,dudy,uucl
921
      integer   stat
922
 
923
      allocate(dvdx(ie,je),STAT=stat)
924
      if (stat.ne.0) print*,'vort: error allocating dvdx'
925
      allocate(uucl(ie,je),STAT=stat)
926
      if (stat.ne.0) print*,'vort: error allocating uucl'
927
      allocate(dudy(ie,je),STAT=stat)
928
      if (stat.ne.0) print*,'vort: error allocating dudy'
929
 
930
      do k=1,ke
931
 
932
      call ddh2m(vv(1,1,k),dvdx,cl,'X',ie,je,1,vmin,vmax,mdv)
933
 
934
      do i=1,ie
935
         do j=1,je
936
            if (uu(i,j,k).eq.mdv) then
937
               uucl(i,j)=mdv
938
            else
939
               uucl(i,j)=uu(i,j,k)*cl(i,j)
940
            endif
941
         enddo
942
      enddo
943
 
944
 
945
      call ddh2m(uucl,dudy,cl,'Y',ie,je,1,vmin,vmax,mdv)
946
 
947
      do i=1,ie
948
         do j=1,je
949
            if ((dvdx(i,j).eq.mdv).or.
950
     >           (dudy(i,j).eq.mdv).or.
951
     >           (cl(i,j).lt.1e-3)) then
952
               var(i,j,k)=mdv
953
            else
954
               var(i,j,k)=1.E4*(dvdx(i,j)-dudy(i,j)
955
     >              /cl(i,j))
956
            endif
957
         enddo
958
      enddo
959
 
960
      enddo
961
 
962
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
963
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
964
      IF (ALLOCATED(uucl)) DEALLOCATE(uucl)
965
 
966
      end
967
 
968
      subroutine curvo(var,uu,vv,sp,cl,ie,je,ke,ak,bk,vmin,vmax)
969
C     ==========================================================
970
C     calculate curvature vorticity
971
C     =u^2*d[v:x]-v^2*d[u*cos:y]/cos-u*v*(d[u:x]-d[v*cos:y]/cos))/(u^2+v^2)
972
 
973
c     argument declaration
974
      integer   ie,je,ke
975
      real,intent(OUT) :: var(ie,je,ke)
976
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),
977
     >     		  sp(ie,je),cl(ie,je)
978
      real      ak(ke),bk(ke),vmin(4),vmax(4)
979
 
980
c     variable declaration
981
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dspdx,dspdy
982
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dudx,dvdx,dudy,dvdy,
983
     >					     uucl,vvcl
984
      integer   stat
985
 
986
      allocate(dspdx(ie,je),STAT=stat)
987
      if (stat.ne.0) print*,'vort: error allocating dspdx'
988
      allocate(dspdy(ie,je),STAT=stat)
989
      if (stat.ne.0) print*,'vort: error allocating dspdy'
990
      allocate(dudx(ie,je,ke),STAT=stat)
991
      if (stat.ne.0) print*,'vort: error allocating dudx'
992
      allocate(dvdx(ie,je,ke),STAT=stat)
993
      if (stat.ne.0) print*,'vort: error allocating dvdx'
994
      allocate(uucl(ie,je,ke),STAT=stat)
995
      if (stat.ne.0) print*,'vort: error allocating uucl'
996
      allocate(vvcl(ie,je,ke),STAT=stat)
997
      if (stat.ne.0) print*,'vort: error allocating vvcl'
998
      allocate(dudy(ie,je,ke),STAT=stat)
999
      if (stat.ne.0) print*,'vort: error allocating dudy'
1000
      allocate(dvdy(ie,je,ke),STAT=stat)
1001
      if (stat.ne.0) print*,'vort: error allocating dvdy'
1002
 
1003
      call ddh2(sp,dspdx,cl,'X',ie,je,1,vmin,vmax)
1004
      call ddh2(sp,dspdy,cl,'Y',ie,je,1,vmin,vmax)
1005
      call ddh3(uu,dudx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
1006
      call ddh3(vv,dvdx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
1007
      do k=1,ke
1008
        uucl(1:ie,1:je,k)=uu(1:ie,1:je,k)*cl(1:ie,1:je)
1009
        vvcl(1:ie,1:je,k)=vv(1:ie,1:je,k)*cl(1:ie,1:je)
1010
      enddo
1011
      call ddh3(uucl,dudy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
1012
      call ddh3(vvcl,dvdy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
1013
 
1014
      do k=1,ke
1015
      do j=1,je
1016
      do i=1,ie
1017
        var(i,j,k)=1.e4*(uu(i,j,k)**2.*dvdx(i,j,k)-
1018
     >                   vv(i,j,k)**2.*dudy(i,j,k)/cl(i,j)-
1019
     >         uu(i,j,k)*vv(i,j,k)*(dudx(i,j,k)-dvdy(i,j,k)/cl(i,j)))/
1020
     >         max(uu(i,j,k)**2.+vv(i,j,k)**2.,1.e-3)
1021
      enddo
1022
      enddo
1023
      enddo
1024
 
1025
c     set poles values to mdv
1026
      if (vmax(2).eq.90.) then
1027
        do k=1,ke
1028
        do i=1,ie
1029
          var(i,je,k)=-999.98999
1030
        enddo
1031
        enddo
1032
      endif
1033
 
1034
      if (vmin(2).eq.-90.) then
1035
        do k=1,ke
1036
        do i=1,ie
1037
          var(i,1,k)=-999.98999
1038
        enddo
1039
        enddo
1040
      endif
1041
 
1042
      IF (ALLOCATED(dspdx)) DEALLOCATE(dspdx)
1043
      IF (ALLOCATED(dspdy)) DEALLOCATE(dspdy)
1044
      IF (ALLOCATED(dudx)) DEALLOCATE(dudx)
1045
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
1046
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
1047
      IF (ALLOCATED(dvdy)) DEALLOCATE(dvdy)
1048
      IF (ALLOCATED(uucl)) DEALLOCATE(uucl)
1049
      IF (ALLOCATED(vvcl)) DEALLOCATE(vvcl)
1050
 
1051
      end
1052
 
1053
      subroutine potvort(pv,uu,vv,th,sp,cl,f,ie,je,ke,ak,bk,
1054
     >     vmin,vmax)
1055
C     ======================================================
1056
C     calculate PV
1057
 
1058
c     argument declaration
1059
      integer   ie,je,ke
1060
      real      pv(ie,je,ke),uu(ie,je,ke),vv(ie,je,ke),
1061
     >     th(ie,je,ke),sp(ie,je),
1062
     >     cl(ie,je),f(ie,je)
1063
      real      ak(ke),bk(ke),vmin(4),vmax(4)
1064
      real      pvpole
1065
 
1066
c     variable declaration
1067
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dspdx,dspdy
1068
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dthdp,dudp,dvdp
1069
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: uu2,dvdx,dudy,dthdx,dthdy
1070
      integer   i,j,k,stat
1071
 
1072
      allocate(dspdx(ie,je),STAT=stat)
1073
      if (stat.ne.0) print*,'potvort: error allocating dspdx(ie,je)'
1074
      allocate(dspdy(ie,je),STAT=stat)
1075
      if (stat.ne.0) print*,'potvort: error allocating dspdy(ie,je)'
1076
      allocate(dthdp(ie,je,ke),STAT=stat)
1077
      if (stat.ne.0) print*,'potvort: error allocating dthdp'
1078
      allocate(dudp(ie,je,ke),STAT=stat)
1079
      if (stat.ne.0) print*,'potvort: error allocating dudp'
1080
      allocate(dvdp(ie,je,ke),STAT=stat)
1081
      if (stat.ne.0) print*,'potvort: error allocating dvdp'
1082
      allocate(dvdx(ie,je,ke),STAT=stat)
1083
      if (stat.ne.0) print*,'potvort: error allocating dvdx'
1084
      allocate(dudy(ie,je,ke),STAT=stat)
1085
      if (stat.ne.0) print*,'potvort: error allocating dudy'
1086
      allocate(dthdx(ie,je,ke),STAT=stat)
1087
      if (stat.ne.0) print*,'potvort: error allocating dthdx'
1088
      allocate(dthdy(ie,je,ke),STAT=stat)
1089
      if (stat.ne.0) print*,'potvort: error allocating dthdy'
1090
      allocate(uu2(ie,je,ke),STAT=stat)
1091
      if (stat.ne.0) print*,'potvort: error allocating uu2'
1092
 
1093
      call ddp(th,dthdp,sp,ie,je,ke,ak,bk)
1094
      call ddp(uu,dudp,sp,ie,je,ke,ak,bk)
1095
      call ddp(vv,dvdp,sp,ie,je,ke,ak,bk)
1096
      call ddh2(sp,dspdx,cl,'X',ie,je,1,vmin,vmax)
1097
      call ddh2(sp,dspdy,cl,'Y',ie,je,1,vmin,vmax)
1098
      call ddh3(th,dthdx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
1099
      call ddh3(vv,dvdx,sp,dspdx,cl,'X',ie,je,ke,vmin,vmax,ak,bk)
1100
      call ddh3(th,dthdy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
1101
c     conversion of uu for spheric coordinates (uu*cos(phi))
1102
      do k=1,ke
1103
         uu2(1:ie,1:je,k)=uu(1:ie,1:je,k)*cl(1:ie,1:je)
1104
      enddo
1105
      call ddh3(uu2,dudy,sp,dspdy,cl,'Y',ie,je,ke,vmin,vmax,ak,bk)
1106
c     conversion of dudy for spheric coordinates (dudy/cos(phi))
1107
      do k=1,ke
1108
         dudy(1:ie,1:je,k)=dudy(1:ie,1:je,k)/cl(1:ie,1:je)
1109
      enddo
1110
 
1111
      do k=1,ke
1112
         pv(1:ie,1:je,k)=1.E6*9.80665*(
1113
     >        -(-dudy(1:ie,1:je,k)+dvdx(1:ie,1:je,k)
1114
     >        +f(1:ie,1:je))*dthdp(1:ie,1:je,k)
1115
     >        -(dudp(1:ie,1:je,k)*dthdy(1:ie,1:je,k)
1116
     >        -dvdp(1:ie,1:je,k)*dthdx(1:ie,1:je,k)))
1117
      enddo
1118
 
1119
c     interpolate poles from neighbouring points on every level
1120
      if (vmax(2).gt.89.5) then
1121
        do k=1,ke
1122
          pvpole=0.
1123
          do i=1,ie
1124
            pvpole=pvpole+pv(i,je-1,k)
1125
          enddo
1126
          pvpole=pvpole/real(ie)
1127
          do i=1,ie
1128
            pv(i,je,k)=pvpole
1129
          enddo
1130
        enddo
1131
      endif
1132
 
1133
      if (vmin(2).lt.-89.5) then
1134
        do k=1,ke
1135
          pvpole=0.
1136
          do i=1,ie
1137
            pvpole=pvpole+pv(i,2,k)
1138
          enddo
1139
          pvpole=pvpole/real(ie)
1140
          do i=1,ie
1141
            pv(i,1,k)=pvpole
1142
          enddo
1143
        enddo
1144
      endif
1145
 
1146
      IF (ALLOCATED(dspdx)) DEALLOCATE(dspdx)
1147
      IF (ALLOCATED(dspdy)) DEALLOCATE(dspdy)
1148
      IF (ALLOCATED(dthdp)) DEALLOCATE(dthdp)
1149
      IF (ALLOCATED(dudp)) DEALLOCATE(dudp)
1150
      IF (ALLOCATED(dvdp)) DEALLOCATE(dvdp)
1151
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
1152
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
1153
      IF (ALLOCATED(dthdx)) DEALLOCATE(dthdx)
1154
      IF (ALLOCATED(dthdy)) DEALLOCATE(dthdy)
1155
      IF (ALLOCATED(uu2)) DEALLOCATE(uu2)
1156
 
1157
      end
1158
 
1159
      subroutine potvort_i(pv,uu,vv,p,cl,f,ie,je,ke,ak,
1160
     >     vmin,vmax,mdv)
1161
C     =================================================
1162
C     calculate isentropic PV
1163
C     !! calculating PV makes sense only
1164
C     when theta levels are very close to each others !!
1165
C     with misdat (mdv) treatment 
1166
 
1167
c     argument declaration
1168
      integer   ie,je,ke
1169
      real      pv(ie,je,ke),uu(ie,je,ke),vv(ie,je,ke),
1170
     >     	p(ie,je,ke),cl(ie,je),f(ie,je)
1171
      real      ak(ke),vmin(4),vmax(4),mdv
1172
      real      pvpole
1173
 
1174
c     variable declaration
1175
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dpdth
1176
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: uu2,dvdx,dudy
1177
      integer   i,j,k,stat
1178
 
1179
      allocate(dpdth(ie,je,ke),STAT=stat)
1180
      if (stat.ne.0) print*,'potvort_i: error allocating dpdth'
1181
      allocate(dvdx(ie,je,ke),STAT=stat)
1182
      if (stat.ne.0) print*,'potvort_i: error allocating dvdx'
1183
      allocate(dudy(ie,je,ke),STAT=stat)
1184
      if (stat.ne.0) print*,'potvort_i: error allocating dudy'
1185
      allocate(uu2(ie,je,ke),STAT=stat)
1186
      if (stat.ne.0) print*,'potvort_i: error allocating uu2'
1187
 
1188
c     calculate dp/dth
1189
c     k=1
1190
      do i=1,ie
1191
         do j=1,je
1192
            if ((p(i,j,2).eq.mdv).or.(p(i,j,1).eq.mdv)) then
1193
               dpdth(i,j,1)=mdv
1194
            else
1195
               dpdth(i,j,1)=100.*(p(i,j,2)-p(i,j,1))/(ak(2)-ak(1))
1196
            endif
1197
         enddo
1198
      enddo
1199
 
1200
c     k=2,ke-1
1201
      do i=1,ie
1202
      do j=1,je
1203
      do k=2,ke-1
1204
        if ((p(i,j,k+1).eq.mdv).or.(p(i,j,k-1).eq.mdv)) then
1205
          dpdth(i,j,k)=mdv
1206
        else
1207
          dpdth(i,j,k)=100.*(p(i,j,k+1)-p(i,j,k-1))/(ak(k+1)-ak(k-1))
1208
        endif
1209
      enddo
1210
      enddo
1211
      enddo
1212
 
1213
c     k=ke
1214
      do i=1,ie
1215
      do j=1,je
1216
        if ((p(i,j,ke).eq.mdv).or.(p(i,j,ke-1).eq.mdv)) then
1217
          dpdth(i,j,ke)=mdv
1218
        else
1219
          dpdth(i,j,ke)=100.*(p(i,j,ke)-p(i,j,ke-1))/(ak(ke)-ak(ke-1))
1220
        endif
1221
      enddo
1222
      enddo
1223
 
1224
      call ddh2m(vv,dvdx,cl,'X',ie,je,ke,vmin,vmax,mdv)
1225
c     conversion of uu for spheric coordinates (uu*cos(phi))
1226
 
1227
 
1228
      do i=1,ie
1229
         do j=1,je
1230
            do k=1,ke
1231
               if (uu(i,j,k).ne.mdv) then
1232
                  uu2(i,j,k)=uu(i,j,k)*cl(i,j)
1233
               else
1234
                  uu2(i,j,k)=mdv
1235
               endif
1236
            enddo
1237
         enddo
1238
      enddo
1239
      call ddh2m(uu2,dudy,cl,'Y',ie,je,ke,vmin,vmax,mdv)
1240
c     conversion of dudy for spheric coordinates (dudy/cos(phi))
1241
      do i=1,ie
1242
         do j=1,je
1243
            do k=1,ke
1244
               if (dudy(i,j,k).ne.mdv) then
1245
                  dudy(i,j,k)=dudy(i,j,k)/cl(i,j)
1246
               else
1247
                  dudy(i,j,k)=mdv
1248
               endif
1249
            enddo
1250
         enddo
1251
      enddo
1252
 
1253
      do i=1,ie
1254
      do j=1,je
1255
      do k=1,ke
1256
        if ((dudy(i,j,k).eq.mdv).or.(dvdx(i,j,k).eq.mdv).or.
1257
     >      (dpdth(i,j,k).eq.mdv).or.(dpdth(i,j,k).eq.0.)) then
1258
          pv(i,j,k)=mdv
1259
        else
1260
          pv(i,j,k)=-1.E6*9.80665*
1261
     >      (-dudy(i,j,k)+dvdx(i,j,k)+f(i,j))/dpdth(i,j,k)
1262
        endif
1263
      enddo
1264
      enddo
1265
      enddo
1266
 
1267
c     interpolate poles from neighbouring points on every level
1268
      if (vmax(2).eq.90.) then
1269
        do k=1,ke
1270
          pvpole=0.
1271
          do i=1,ie
1272
            pvpole=pvpole+pv(i,je-1,k)
1273
          enddo
1274
          pvpole=pvpole/real(ie)
1275
          do i=1,ie
1276
            pv(i,je,k)=pvpole
1277
          enddo
1278
        enddo
1279
      endif
1280
 
1281
      if (vmin(2).eq.-90.) then
1282
        do k=1,ke
1283
          pvpole=0.
1284
          do i=1,ie
1285
            pvpole=pvpole+pv(i,2,k)
1286
          enddo
1287
          pvpole=pvpole/real(ie)
1288
          do i=1,ie
1289
            pv(i,1,k)=pvpole
1290
          enddo
1291
        enddo
1292
      endif
1293
 
1294
      IF (ALLOCATED(dpdth)) DEALLOCATE(dpdth)
1295
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
1296
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
1297
      IF (ALLOCATED(uu2)) DEALLOCATE(uu2)
1298
 
1299
      end
1300
 
1301
 
1302
      subroutine gpotvort_i2(dpv,uu,vv,p,cl,f,ie,je,ke,ak,
1303
     >     vmin,vmax,mdv)
1304
C     =================================================
1305
C     calculate isentropic PV gradient length
1306
C     GPV=SQRT(D[PV:X]^2+D[PV:Y]^2) 
1307
C     !! it calculates PV with potvort_i, what makes sense only
1308
C     when theta levels are very close to each others !!
1309
C     with misdat (mdv) treatment
1310
 
1311
c     argument declaration
1312
      integer   ie,je,ke
1313
      real      dpv(ie,je,ke),uu(ie,je,ke),vv(ie,je,ke),
1314
     >     p(ie,je,ke),cl(ie,je),f(ie,je)
1315
      real      ak(ke),vmin(4),vmax(4),mdv
1316
 
1317
c     variable declaration
1318
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: pv
1319
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dpvdx,dpvdy
1320
      integer   stat
1321
      real      dpvpole
1322
 
1323
      allocate(pv(ie,je,ke),STAT=stat)
1324
      if (stat.ne.0) print*,'gpotvort_i2: error allocating pv'
1325
      allocate(dpvdx(ie,je),STAT=stat)
1326
      if (stat.ne.0) print*,'gpotvort_i2: error allocating dpvdx'
1327
      allocate(dpvdy(ie,je),STAT=stat)
1328
      if (stat.ne.0) print*,'gpotvort_i2: error allocating dpvdy'
1329
 
1330
      call potvort_i(pv,uu,vv,p,cl,f,ie,je,ke,ak,
1331
     >        vmin,vmax,mdv)
1332
 
1333
      do k=1,ke
1334
         call ddh2m(pv(1,1,k),dpvdx,cl,'X',ie,je,1,vmin,vmax,mdv)
1335
         call ddh2m(pv(1,1,k),dpvdy,cl,'Y',ie,je,1,vmin,vmax,mdv)
1336
         do j=1,je
1337
            do i=1,ie
1338
               if ((dpvdx(i,j).ne.mdv).and.(dpvdy(i,j).ne.mdv)) then
1339
                  dpv(i,j,k)=1.E6*sqrt(dpvdx(i,j)**2.+dpvdy(i,j)**2.)
1340
               else
1341
                  dpv(i,j,k)=mdv
1342
               endif
1343
            enddo
1344
         enddo
1345
      enddo
1346
 
1347
      IF (ALLOCATED(pv)) DEALLOCATE(pv)
1348
      IF (ALLOCATED(dpvdx)) DEALLOCATE(dpvdx)
1349
      IF (ALLOCATED(dpvdy)) DEALLOCATE(dpvdy)
1350
 
1351
      end
1352
 
1353
      subroutine gpotvort_i(dpv,pv,cl,ie,je,ke,
1354
     >     vmin,vmax,mdv)
1355
C     =================================================
1356
C     calculate isentropic PV gradient length
1357
C     GPV=SQRT(D[PV:X]^2+D[PV:Y]^2) 
1358
C     !! it uses PV as input, best calculated with
1359
C     potvort (with data from P file)  !!
1360
C     with misdat (mdv) treatment
1361
C     990201 mark liniger
1362
 
1363
c     argument declaration
1364
      integer   ie,je,ke
1365
      real      dpv(ie,je,ke),pv(ie,je,ke),cl(ie,je)
1366
      real      vmin(4),vmax(4),mdv
1367
 
1368
c     variable declaration
1369
      REAL,ALLOCATABLE, DIMENSION (:,:) :: dpvdx,dpvdy
1370
      integer   stat
1371
 
1372
      allocate(dpvdx(ie,je),STAT=stat)
1373
      if (stat.ne.0) print*,'gpotvort_i: error allocating dpvdx'
1374
      allocate(dpvdy(ie,je),STAT=stat)
1375
      if (stat.ne.0) print*,'gpotvort_i: error allocating dpvdy'
1376
 
1377
      do k=1,ke
1378
         call ddh2m(pv(1,1,k),dpvdx,cl,'X',ie,je,1,vmin,vmax,mdv)
1379
         call ddh2m(pv(1,1,k),dpvdy,cl,'Y',ie,je,1,vmin,vmax,mdv)
1380
         do j=1,je
1381
            do i=1,ie
1382
               if ((dpvdx(i,j).ne.mdv).and.(dpvdy(i,j).ne.mdv)) then
1383
                  dpv(i,j,k)=1.E6*sqrt(dpvdx(i,j)**2.+dpvdy(i,j)**2.)
1384
               else
1385
                  dpv(i,j,k)=mdv
1386
               endif
1387
            enddo
1388
         enddo
1389
      enddo
1390
 
1391
      IF (ALLOCATED(dpvdx)) DEALLOCATE(dpvdx)
1392
      IF (ALLOCATED(dpvdy)) DEALLOCATE(dpvdy)
1393
 
1394
      end
1395
 
1396
 
1397
      subroutine divqu(var,uu,vv,qq,cl,ie,je,ke,vmin,vmax,mdv)
1398
C     ========================================================
1399
C     calculate divergence DIV=D[U:X]+D[V*cos(phi):Y]/cos(phi)
1400
C     with misdat (mdv) and pole treatment
1401
 
1402
c     argument declaration
1403
      integer   ie,je,ke
1404
      real,intent(OUT) :: var(ie,je,ke)
1405
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),qq(ie,je,ke),
1406
     >			  cl(ie,je)
1407
      real      vmin(4),vmax(4),mdv
1408
 
1409
 
1410
c     variable declaration
1411
      REAL,ALLOCATABLE, DIMENSION (:,:) :: qv2,dqudx,dqvdy
1412
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: qu,qv
1413
      integer   stat
1414
      real      pole
1415
 
1416
      allocate(qv2(ie,je),STAT=stat)
1417
      if (stat.ne.0) print*,'divqu: error allocating qv2'
1418
      allocate(dqudx(ie,je),STAT=stat)
1419
      if (stat.ne.0) print*,'divqu: error allocating dqudx'
1420
      allocate(dqvdy(ie,je),STAT=stat)
1421
      if (stat.ne.0) print*,'divqu: error allocating dqvdy'
1422
      allocate(qu(ie,je,ke),STAT=stat)
1423
      if (stat.ne.0) print*,'divqu: error allocating qu'
1424
      allocate(qv(ie,je,ke),STAT=stat)
1425
      if (stat.ne.0) print*,'divqu: error allocating qv'
1426
 
1427
      do k=1,ke
1428
      do j=1,je
1429
      do i=1,ie
1430
        qu(i,j,k)=qq(i,j,k)*uu(i,j,k)
1431
        qv(i,j,k)=qq(i,j,k)*vv(i,j,k)
1432
      enddo
1433
      enddo
1434
      enddo
1435
 
1436
      do k=1,ke
1437
        call ddh2m(qu(1,1,k),dqudx,cl,'X',ie,je,1,vmin,vmax,mdv)
1438
 
1439
c       conversion of qv for spheric coordinates (qv*cos(phi))
1440
        do j=1,je
1441
           do i=1,ie
1442
              if (qv(i,j,k).ne.mdv) then
1443
                 qv2(i,j)=qv(i,j,k)*cl(i,j)
1444
              else
1445
                 qv2(i,k)=mdv
1446
              endif
1447
           enddo
1448
        enddo
1449
 
1450
 
1451
        call ddh2m(qv2,dqvdy,cl,'Y',ie,je,1,vmin,vmax,mdv)
1452
 
1453
c       conversion of dqvdy for spheric coordinates (dqvdy/cos(phi))
1454
        do j=1,je
1455
           do i=1,ie
1456
              if ((dqudx(i,j).ne.mdv).and.(dqvdy(i,j).ne.mdv)) then
1457
                 var(i,j,k)=1.E6*(dqudx(i,j)+dqvdy(i,j)/cl(i,j))
1458
              else
1459
                 var(i,j,k)=mdv
1460
              endif
1461
           enddo
1462
        enddo
1463
      enddo
1464
 
1465
c     interpolate poles from neighbouring points on every level
1466
      if (vmax(2).eq.90.) then
1467
        do k=1,ke
1468
          pole=0.
1469
          counter=0.
1470
          do i=1,ie
1471
             if (var(i,je-1,k).ne.mdv) then
1472
                counter=counter+1
1473
                pole=pole+var(i,je-1,k)
1474
             endif
1475
          enddo
1476
          if (counter.ne.0) then
1477
             pole=pole/counter
1478
          else
1479
             pole=mdv
1480
          endif
1481
          do i=1,ie
1482
            var(i,je,k)=pole
1483
          enddo
1484
        enddo
1485
      endif
1486
 
1487
      if (vmin(2).eq.-90.) then
1488
        do k=1,ke
1489
          pole=0.
1490
          counter=0.
1491
          do i=1,ie
1492
             if (var(i,2,k).ne.mdv) then
1493
                counter=counter+1
1494
                pole=pole+var(i,2,k)
1495
             endif
1496
          enddo
1497
          if (counter.ne.0) then
1498
             pole=pole/counter
1499
          else
1500
             pole=mdv
1501
          endif
1502
          do i=1,ie
1503
            var(i,1,k)=pole
1504
          enddo
1505
        enddo
1506
      endif
1507
 
1508
      IF (ALLOCATED(dqudx)) DEALLOCATE(dqudx)
1509
      IF (ALLOCATED(dqvdy)) DEALLOCATE(dqvdy)
1510
      IF (ALLOCATED(qu)) DEALLOCATE(qu)
1511
      IF (ALLOCATED(qv)) DEALLOCATE(qv)
1512
      IF (ALLOCATED(qv2)) DEALLOCATE(qv2)
1513
 
1514
      end
1515
      subroutine div_i(var,uu,vv,cl,ie,je,ke,vmin,vmax,mdv)
1516
C     =====================================================
1517
C     calculate divergence DIV=D[U:X]+D[V*cos(phi):Y]/cos(phi)
1518
C     with misdat (mdv) and pole treatment
1519
C     990201 mark liniger
1520
 
1521
c     argument declaration
1522
      integer   ie,je,ke
1523
      real,intent(OUT) :: var(ie,je,ke)
1524
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),cl(ie,je)
1525
      real      vmin(4),vmax(4),mdv
1526
 
1527
 
1528
c     variable declaration
1529
      REAL,ALLOCATABLE, DIMENSION (:,:) :: vv2,dudx,dvdy
1530
      integer   stat
1531
      real      pole
1532
 
1533
      allocate(vv2(ie,je),STAT=stat)
1534
      if (stat.ne.0) print*,'div_i: error allocating vv2'
1535
      allocate(dudx(ie,je),STAT=stat)
1536
      if (stat.ne.0) print*,'div_i: error allocating dudx'
1537
      allocate(dvdy(ie,je),STAT=stat)
1538
      if (stat.ne.0) print*,'div_i: error allocating dvdy'
1539
 
1540
 
1541
      do k=1,ke
1542
        call ddh2m(uu(1,1,k),dudx,cl,'X',ie,je,1,vmin,vmax,mdv)
1543
 
1544
c       conversion of vv for spheric coordinates (vv*cos(phi))
1545
        do j=1,je
1546
           do i=1,ie
1547
              if (vv(i,j,k).ne.mdv) then
1548
                 vv2(i,j)=vv(i,j,k)*cl(i,j)
1549
              else
1550
                 vv2(i,k)=mdv
1551
              endif
1552
           enddo
1553
        enddo
1554
 
1555
 
1556
        call ddh2m(vv2,dvdy,cl,'Y',ie,je,1,vmin,vmax,mdv)
1557
 
1558
c       conversion of dvdy for spheric coordinates (dvdy/cos(phi))
1559
        do j=1,je
1560
           do i=1,ie
1561
              if ((dudx(i,j).ne.mdv).and.(dvdy(i,j).ne.mdv)) then
1562
                 var(i,j,k)=1.E6*(dudx(i,j)+dvdy(i,j)/cl(i,j))
1563
              else
1564
                 var(i,j,k)=mdv
1565
              endif
1566
           enddo
1567
        enddo
1568
      enddo
1569
 
1570
c     interpolate poles from neighbouring points on every level
1571
      if (vmax(2).eq.90.) then
1572
        do k=1,ke
1573
          pole=0.
1574
          counter=0.
1575
          do i=1,ie
1576
             if (var(i,je-1,k).ne.mdv) then
1577
                counter=counter+1
1578
                pole=pole+var(i,je-1,k)
1579
             endif
1580
          enddo
1581
          if (counter.ne.0) then
1582
             pole=pole/counter
1583
          else
1584
             pole=mdv
1585
          endif
1586
          do i=1,ie
1587
            var(i,je,k)=pole
1588
          enddo
1589
        enddo
1590
      endif
1591
 
1592
      if (vmin(2).eq.-90.) then
1593
        do k=1,ke
1594
          pole=0.
1595
          counter=0.
1596
          do i=1,ie
1597
             if (var(i,2,k).ne.mdv) then
1598
                counter=counter+1
1599
                pole=pole+var(i,2,k)
1600
             endif
1601
          enddo
1602
          if (counter.ne.0) then
1603
             pole=pole/counter
1604
          else
1605
             pole=mdv
1606
          endif
1607
          do i=1,ie
1608
            var(i,1,k)=pole
1609
          enddo
1610
        enddo
1611
      endif
1612
 
1613
      IF (ALLOCATED(dudx)) DEALLOCATE(dudx)
1614
      IF (ALLOCATED(vv2)) DEALLOCATE(vv2)
1615
      IF (ALLOCATED(dvdy)) DEALLOCATE(dvdy)
1616
 
1617
      end
1618
 
1619
      subroutine def_i(var,uu,vv,cl,ie,je,ke,vmin,vmax,mdv)
1620
C     =====================================================
1621
C     calculate deformation (strain) DEF= 1e6*SQRT(
1622
C       ( D[V:X]+D[U*cos(phi):Y]/cos(phi) )^2
1623
C       ( D[U:X]-D[V*cos(phi):Y]/cos(phi) )^2 )
1624
C     with misdat (mdv) treatment
1625
C     990501 mark liniger
1626
 
1627
c     argument declaration
1628
      integer   ie,je,ke
1629
      real,intent(OUT) :: var(ie,je,ke)
1630
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),cl(ie,je)
1631
      real      vmin(4),vmax(4),mdv
1632
 
1633
c     variable declaration
1634
      REAL,ALLOCATABLE, DIMENSION (:,:) :: uu2,vv2,dudx,dvdx,dudy,dvdy
1635
      integer   stat
1636
      real      pole
1637
 
1638
      allocate(uu2(ie,je),STAT=stat)
1639
      if (stat.ne.0) print*,'dev_i: error allocating uu2'
1640
      allocate(vv2(ie,je),STAT=stat)
1641
      if (stat.ne.0) print*,'dev_i: error allocating vv2'
1642
      allocate(dudx(ie,je),STAT=stat)
1643
      if (stat.ne.0) print*,'dev_i: error allocating dudx'
1644
      allocate(dvdx(ie,je),STAT=stat)
1645
      if (stat.ne.0) print*,'dev_i: error allocating dvdx'
1646
      allocate(dudy(ie,je),STAT=stat)
1647
      if (stat.ne.0) print*,'dev_i: error allocating dudy'
1648
      allocate(dvdy(ie,je),STAT=stat)
1649
      if (stat.ne.0) print*,'dev_i: error allocating dvdy'
1650
 
1651
 
1652
      do k=1,ke
1653
        call ddh2m(uu(1,1,k),dudx,cl,'X',ie,je,1,vmin,vmax,mdv)
1654
        call ddh2m(vv(1,1,k),dvdx,cl,'X',ie,je,1,vmin,vmax,mdv)
1655
 
1656
c       conversion of uu/vv for spheric coordinates (??*cos(phi))
1657
 
1658
        do j=1,je
1659
           do i=1,ie
1660
              if (uu(i,j,k).ne.mdv) then
1661
                 uu2(i,j)=uu(i,j,k)*cl(i,j)
1662
              else
1663
                 uu2(i,j)=mdv
1664
              endif
1665
              if (vv(i,j,k).ne.mdv) then
1666
                 vv2(i,j)=vv(i,j,k)*cl(i,j)
1667
              else
1668
                 vv2(i,j)=mdv
1669
              endif
1670
           enddo
1671
        enddo
1672
 
1673
        call ddh2m(uu2,dudy,cl,'Y',ie,je,1,vmin,vmax,mdv)
1674
        call ddh2m(vv2,dvdy,cl,'Y',ie,je,1,vmin,vmax,mdv)
1675
 
1676
c       conversion of d?dy for spheric coordinates (d?dy/cos(phi))
1677
        do j=1,je
1678
           do i=1,ie
1679
              if ((dudx(i,j).ne.mdv).and.
1680
     >             (dudy(i,j).ne.mdv).and.
1681
     >             (dvdx(i,j).ne.mdv).and.
1682
     >             (dvdy(i,j).ne.mdv)) then
1683
 
1684
                 var(i,j,k)=1.e6*sqrt(
1685
     >                (dvdx(i,j)+dudy(i,j)/cl(i,j))
1686
     >                *(dvdx(i,j)+dudy(i,j)/cl(i,j))
1687
     >                +(dudx(i,j)-dvdy(i,j)/cl(i,j))
1688
     >                *(dudx(i,j)-dvdy(i,j)/cl(i,j)))
1689
 
1690
              else
1691
                 var(i,j,k)=mdv
1692
              endif
1693
           enddo
1694
        enddo
1695
      enddo
1696
 
1697
c     interpolate poles from neighbouring points on every level
1698
      if (vmax(2).eq.90.) then
1699
        do k=1,ke
1700
          pole=0.
1701
          do i=1,ie
1702
            pole=pole+var(i,je-1,k)
1703
          enddo
1704
          pole=pole/real(ie)
1705
          do i=1,ie
1706
            var(i,je,k)=pole
1707
          enddo
1708
        enddo
1709
      endif
1710
 
1711
      if (vmin(2).eq.-90.) then
1712
        do k=1,ke
1713
          pole=0.
1714
          do i=1,ie
1715
            pole=pole+var(i,2,k)
1716
          enddo
1717
          pole=pole/real(ie)
1718
          do i=1,ie
1719
            var(i,1,k)=pole
1720
          enddo
1721
        enddo
1722
      endif
1723
 
1724
      IF (ALLOCATED(uu2)) DEALLOCATE(uu2)
1725
      IF (ALLOCATED(vv2)) DEALLOCATE(vv2)
1726
      IF (ALLOCATED(dudx)) DEALLOCATE(dudx)
1727
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
1728
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
1729
      IF (ALLOCATED(dvdy)) DEALLOCATE(dvdy)
1730
 
1731
      end
1732
 
1733
 
1734
      subroutine defn_i(var,uu,vv,cl,ie,je,ke,vmin,vmax,mdv)
1735
C     =====================================================
1736
C     calculate normal deformation (normal strain) 
1737
C     DEFN= 1e6*( D[U:X]-D[V*cos(phi):Y]/cos(phi) )
1738
C     with misdat (mdv) treatment
1739
C     990727 mark liniger
1740
 
1741
c     argument declaration
1742
      integer   ie,je,ke
1743
      real,intent(OUT) :: var(ie,je,ke)
1744
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),cl(ie,je)
1745
      real      vmin(4),vmax(4),mdv
1746
 
1747
c     variable declaration
1748
      REAL,ALLOCATABLE, DIMENSION (:,:) :: vv2,dudx,dvdy
1749
      integer   stat
1750
      real      pole
1751
 
1752
      allocate(vv2(ie,je),STAT=stat)
1753
      if (stat.ne.0) print*,'dev_i: error allocating vv2'
1754
      allocate(dudx(ie,je),STAT=stat)
1755
      if (stat.ne.0) print*,'dev_i: error allocating dudx'
1756
      allocate(dvdy(ie,je),STAT=stat)
1757
      if (stat.ne.0) print*,'dev_i: error allocating dvdy'
1758
 
1759
 
1760
      do k=1,ke
1761
        call ddh2m(uu(1,1,k),dudx,cl,'X',ie,je,1,vmin,vmax,mdv)
1762
 
1763
c       conversion of uu/vv for spheric coordinates (??*cos(phi)) 
1764
        do j=1,je
1765
           do i=1,ie
1766
              if (vv(i,j,k).ne.mdv) then
1767
                 vv2(i,j)=vv(i,j,k)*cl(i,j)
1768
              else
1769
                 vv2(i,j)=mdv
1770
              endif
1771
           enddo
1772
        enddo
1773
 
1774
        call ddh2m(vv2,dvdy,cl,'Y',ie,je,1,vmin,vmax,mdv)
1775
 
1776
c       final calculation
1777
        do j=1,je
1778
           do i=1,ie
1779
              if ((dudx(i,j).ne.mdv).and.
1780
     >             (dvdy(i,j).ne.mdv)) then 
1781
                 var(i,j,k)=1.e6*(dudx(i,j)-dvdy(i,j)/cl(i,j))
1782
              else
1783
                 var(i,j,k)=mdv
1784
              endif
1785
           enddo
1786
        enddo
1787
      enddo
1788
 
1789
c     interpolate poles from neighbouring points on every level
1790
      if (vmax(2).eq.90.) then
1791
        do k=1,ke
1792
          pole=0.
1793
          do i=1,ie
1794
            pole=pole+var(i,je-1,k)
1795
          enddo
1796
          pole=pole/real(ie)
1797
          do i=1,ie
1798
            var(i,je,k)=pole
1799
          enddo
1800
        enddo
1801
      endif
1802
 
1803
      if (vmin(2).eq.-90.) then
1804
        do k=1,ke
1805
          pole=0.
1806
          do i=1,ie
1807
            pole=pole+var(i,2,k)
1808
          enddo
1809
          pole=pole/real(ie)
1810
          do i=1,ie
1811
            var(i,1,k)=pole
1812
          enddo
1813
        enddo
1814
      endif
1815
 
1816
      IF (ALLOCATED(vv2)) DEALLOCATE(vv2)
1817
      IF (ALLOCATED(dudx)) DEALLOCATE(dudx)
1818
      IF (ALLOCATED(dvdy)) DEALLOCATE(dvdy)
1819
 
1820
      end
1821
 
1822
 
1823
 
1824
 
1825
      subroutine defs_i(var,uu,vv,cl,ie,je,ke,vmin,vmax,mdv)
1826
C     =====================================================
1827
C     calculate shear deformation (shear strain) 
1828
C     DEFS= 1e6*( D[V:X]+D[U*cos(phi):Y]/cos(phi)) 
1829
C     with misdat (mdv) treatment
1830
C     990727 mark liniger
1831
 
1832
c     argument declaration
1833
      integer   ie,je,ke
1834
      real,intent(OUT) :: var(ie,je,ke)
1835
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),cl(ie,je)
1836
      real      vmin(4),vmax(4),mdv
1837
 
1838
c     variable declaration
1839
      REAL,ALLOCATABLE, DIMENSION (:,:) :: uu2,dvdx,dudy
1840
      integer   stat
1841
      real      pole
1842
 
1843
      allocate(uu2(ie,je),STAT=stat)
1844
      if (stat.ne.0) print*,'dev_i: error allocating uu2'
1845
      allocate(dvdx(ie,je),STAT=stat)
1846
      if (stat.ne.0) print*,'dev_i: error allocating dvdx'
1847
      allocate(dudy(ie,je),STAT=stat)
1848
      if (stat.ne.0) print*,'dev_i: error allocating dudy'
1849
 
1850
 
1851
      do k=1,ke
1852
        call ddh2m(vv(1,1,k),dvdx,cl,'X',ie,je,1,vmin,vmax,mdv)
1853
 
1854
c       conversion of uu/vv for spheric coordinates (??*cos(phi)) 
1855
        do j=1,je
1856
           do i=1,ie
1857
              if (uu(i,j,k).ne.mdv) then
1858
                 uu2(i,j)=uu(i,j,k)*cl(i,j)
1859
              else
1860
                 uu2(i,j)=mdv
1861
              endif
1862
           enddo
1863
        enddo
1864
 
1865
        call ddh2m(uu2,dudy,cl,'Y',ie,je,1,vmin,vmax,mdv)
1866
 
1867
c       conversion of d?dy for spheric coordinates (d?dy/cos(phi))
1868
        do j=1,je
1869
           do i=1,ie
1870
              if ((dudy(i,j).ne.mdv).and.
1871
     >             (dvdx(i,j).ne.mdv)) then 
1872
                 var(i,j,k)=1.e6*(dvdx(i,j)+dudy(i,j)/cl(i,j))
1873
              else
1874
                 var(i,j,k)=mdv
1875
              endif
1876
           enddo
1877
        enddo
1878
      enddo
1879
 
1880
c     interpolate poles from neighbouring points on every level
1881
      if (vmax(2).eq.90.) then
1882
        do k=1,ke
1883
          pole=0.
1884
          do i=1,ie
1885
            pole=pole+var(i,je-1,k)
1886
          enddo
1887
          pole=pole/real(ie)
1888
          do i=1,ie
1889
            var(i,je,k)=pole
1890
          enddo
1891
        enddo
1892
      endif
1893
 
1894
      if (vmin(2).eq.-90.) then
1895
        do k=1,ke
1896
          pole=0.
1897
          do i=1,ie
1898
            pole=pole+var(i,2,k)
1899
          enddo
1900
          pole=pole/real(ie)
1901
          do i=1,ie
1902
            var(i,1,k)=pole
1903
          enddo
1904
        enddo
1905
      endif
1906
 
1907
      IF (ALLOCATED(uu2)) DEALLOCATE(uu2)
1908
      IF (ALLOCATED(dudy)) DEALLOCATE(dudy)
1909
      IF (ALLOCATED(dvdx)) DEALLOCATE(dvdx)
1910
 
1911
      end
1912
 
1913
 
1914
 
1915
 
1916
 
1917
      subroutine okuboweiss_i(var,uu,vv,cl,ie,je,ke,vmin,vmax,mdv)
1918
C     =====================================================
1919
C     calculate okubo-weiss parameter= DEF*DEF-VORT*VORT
1920
C     with misdat (mdv) treatment
1921
C     990720 mark liniger
1922
 
1923
c     argument declaration
1924
      integer   ie,je,ke
1925
      real,intent(OUT) :: var(ie,je,ke)
1926
      real,intent(IN)  :: uu(ie,je,ke),vv(ie,je,ke),cl(ie,je)
1927
      real      vmin(4),vmax(4),mdv
1928
 
1929
c     variable declaration
1930
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: def,vort
1931
      integer   stat
1932
      real      pole
1933
 
1934
      allocate(def(ie,je,ke),STAT=stat)
1935
      if (stat.ne.0) print*,'okuboweiss_i: error allocating def'
1936
      allocate(vort(ie,je,ke),STAT=stat)
1937
      if (stat.ne.0) print*,'okuboweiss_i: error allocating vort'
1938
 
1939
      call def_i(def,uu,vv,cl,ie,je,ke,vmin,vmax,mdv)
1940
      call vort_i(vort,uu,vv,cl,ie,je,ke,vmin,vmax,mdv)
1941
 
1942
      do k=1,ke
1943
         do j=1,je
1944
            do i=1,ie
1945
               if ((def(i,j,k).ne.mdv).and.
1946
     >              (vort(i,j,k).ne.mdv)) then
1947
                  var(i,j,k)=def(i,j,k)*def(i,j,k)
1948
     >                 -vort(i,j,k)*vort(i,j,k)*1e4
1949
               else
1950
                  var(i,j,k)=mdv 
1951
               endif
1952
            enddo
1953
         enddo
1954
      enddo
1955
 
1956
      IF (ALLOCATED(def)) DEALLOCATE(def)
1957
      IF (ALLOCATED(vort)) DEALLOCATE(vort)
1958
 
1959
      end
1960
 
1961
      subroutine rich(ri,sp,uu,vv,th,ie,je,ke,ak,bk)
1962
C     ==============================================
1963
C     calculate Richardson number
1964
 
1965
      real      R,p0,kappa
1966
      data      R,p0,kappa /287.,100000.,0.286/
1967
 
1968
c     argument declaration
1969
      integer   ie,je,ke
1970
      real      ri(ie,je,ke),uu(ie,je,ke),vv(ie,je,ke),
1971
     >          th(ie,je,ke),sp(ie,je)
1972
      real      ak(ke),bk(ke)
1973
 
1974
c     variable declaration
1975
      REAL,ALLOCATABLE, DIMENSION (:,:,:) :: dthdp,dudp,dvdp
1976
      integer   i,j,k,stat
1977
 
1978
      allocate(dthdp(ie,je,ke),STAT=stat)
1979
      if (stat.ne.0) print*,'potvort: error allocating dthdp'
1980
*      allocate(vel(ie,je,ke),STAT=stat)
1981
*      if (stat.ne.0) print*,'potvort: error allocating vel'
1982
      allocate(dudp(ie,je,ke),STAT=stat)
1983
      if (stat.ne.0) print*,'potvort: error allocating dudp'
1984
      allocate(dvdp(ie,je,ke),STAT=stat)
1985
      if (stat.ne.0) print*,'potvort: error allocating dvdp'
1986
 
1987
 
1988
*      vel=sqrt(uu**2+vv**2)
1989
      call ddp(th,dthdp,sp,ie,je,ke,ak,bk)
1990
*      call ddp(vel,dveldp,sp,ie,je,ke,ak,bk)
1991
       call ddp(uu,dudp,sp,ie,je,ke,ak,bk)
1992
       call ddp(vv,dvdp,sp,ie,je,ke,ak,bk)
1993
 
1994
      do k=1,ke
1995
      do j=1,je
1996
      do i=1,ie
1997
         pval=(ak(k)+bk(k)*sp(i,j))*100.
1998
        ri(i,j,k)=-R*(p0**(-kappa))*(pval**(kappa-1.))*
1999
     >        dthdp(i,j,k)/(dudp(i,j,k)**2.+dvdp(i,j,k)**2.)
2000
*         ri(i,j,k)=-R*(p0**(-kappa))*(pval**(kappa-1.))*
2001
*    >             dthdp(i,j,k)/(dveldp(i,j,k)**2.)
2002
      enddo
2003
      enddo
2004
      enddo
2005
 
2006
      IF (ALLOCATED(dthdp)) DEALLOCATE(dthdp)
2007
      IF (ALLOCATED(dudp)) DEALLOCATE(dudp)
2008
      IF (ALLOCATED(dvdp)) DEALLOCATE(dvdp)
2009
 
2010
      end
2011
 
2012
      subroutine ddp(a,d,sp,ie,je,ke,ak,bk)
2013
c-----------------------------------------------------------------------
2014
c     Purpose: VERTICAL DERIVATIVE
2015
c     Compute the vertical derivative without missing data checking.
2016
c     The derivative is taken from array 'a' in the direction of 'P'.
2017
c     The result is stored in array 'd'.
2018
c     3 point weighted centered differencing is used.
2019
c     The vertical level-structure of the data is of the form
2020
c     p=ak(k)+bk(k)*ps.
2021
c-----------------------------------------------------------------------
2022
 
2023
c     declaration of arguments
2024
      integer       ie,je,ke
2025
      real          a(ie,je,ke),d(ie,je,ke),sp(ie,je)
2026
 
2027
c     variable declaration
2028
      integer       i,j,k
2029
      real          dpu,dpl,quot,fac,psrf
2030
      real          ak(ke),bk(ke)
2031
 
2032
c     specify scaling factor associated with derivation with respect
2033
c     to pressure
2034
      fac=0.01
2035
 
2036
c     compute vertical 3 point derivative
2037
c     ---------------------------
2038
c     3-point vertical derivative
2039
c     ---------------------------
2040
      do j=1,je
2041
         do i=1,ie
2042
c     get surface pressure at current grid-point
2043
            psrf=sp(i,j)
2044
c     points at k=1
2045
            dpu=(ak(1)+bk(1)*psrf)-(ak(2)+bk(2)*psrf)
2046
            d(i,j,1)=(a(i,j,1)-a(i,j,2))*fac/dpu
2047
c     points at 1<k<ke
2048
            do k=2,ke-1
2049
               dpu=(ak(k)+bk(k)*psrf)-(ak(k+1)+bk(k+1)*psrf)
2050
               dpl=(ak(k-1)+bk(k-1)*psrf)-(ak(k)+bk(k)*psrf)
2051
               quot=dpu/dpl
2052
               d(i,j,k)=(quot*(a(i,j,k-1)-a(i,j,k))
2053
     &              +1./quot*(a(i,j,k)-a(i,j,k+1)))*fac/(dpu+dpl)
2054
            enddo
2055
c     points at k=ke
2056
            dpl=(ak(ke-1)+bk(ke-1)*psrf)-(ak(ke)+bk(ke)*psrf)
2057
            d(i,j,ke)=(a(i,j,ke-1)-a(i,j,ke))*fac/dpl
2058
         enddo
2059
      enddo
2060
      end
2061
 
2062
 
2063
      subroutine ddh3(a,d,ps,dps,cl,dir,ie,je,ke,datmin,datmax,ak,bk)
2064
c-----------------------------------------------------------------------
2065
c     Purpose: HORIZONTAL DERIVATIVE ON PRESSURE-SURFACES WITHOUT
2066
c     MISSING DATA
2067
c     The derivative is taken from array 'a' in the direction of 'dir',
2068
c     where 'dir' is either 'X','Y'. The result is stored in array 'd'.
2069
c     The routine accounts for derivatives at the pole and periodic
2070
c     boundaries in the longitudinal direction (depending on
2071
c     the value of datmin, datmax). If the data-set does not reach to
2072
c     the pole, a one-sided derivative is taken. Pole-treatment is only
2073
c     carried out if the data-set covers 360 deg in longitude, and it
2074
c     requires that ie=4*ii+1, where ii is an integer.
2075
c     History:
2076
c     Daniel Luethi
2077
c-----------------------------------------------------------------------
2078
 
2079
c     declaration of arguments
2080
      integer       ie,je,ke
2081
      real          a(ie,je,ke),d(ie,je,ke),cl(ie,je)
2082
      real          ps(ie,je),dps(ie,je)
2083
      real          datmin(4),datmax(4)
2084
      character*(*) dir
2085
 
2086
c     variable declaration
2087
      integer       i,j,k
2088
      real          ak(ke),bk(ke),as(500),bs(500)
2089
 
2090
c     compute vertical derivatives of ak's and bk's
2091
      do k=2,ke-1
2092
         as(k)=(ak(k-1)-ak(k+1))/2.
2093
         bs(k)=(bk(k-1)-bk(k+1))/2.
2094
      enddo
2095
      as(1 )=ak(1)-ak(2)
2096
      bs(1 )=bk(1)-bk(2)
2097
      as(ke)=ak(ke-1)-ak(ke)
2098
      bs(ke)=bk(ke-1)-bk(ke)
2099
 
2100
c     compute horizontal derivatives on sigma surfaces
2101
      call ddh2(a,d,cl,dir,ie,je,ke,datmin,datmax)
2102
 
2103
c     apply correction for horizontal derivative on p-surfaces
2104
      do j=1,je
2105
         do i=1,ie
2106
            do k=2,ke-1
2107
               d(i,j,k)=d(i,j,k)+bk(k)*dps(i,j)/2./(as(k)+
2108
     &              bs(k)*ps(i,j))*(a(i,j,k+1)-a(i,j,k-1))
2109
            enddo
2110
            k=1
2111
            d(i,j,k)=d(i,j,k)+bk(k)*dps(i,j)/(as(k)+
2112
     &           bs(k)*ps(i,j))*(a(i,j,k+1)-a(i,j,k))
2113
            k=ke
2114
            d(i,j,k)=d(i,j,k)+bk(k)*dps(i,j)/(as(k)+
2115
     &           bs(k)*ps(i,j))*(a(i,j,k)-a(i,j,k-1))
2116
         enddo
2117
      enddo
2118
      end
2119
 
2120
 
2121
      subroutine ddh2(a,d,cl,dir,ie,je,ke,datmin,datmax)
2122
c-----------------------------------------------------------------------
2123
c     Purpose: HORIZONTAL DERIVATIVE ON DATA-SURFACES WITHOUT MISSING
2124
C     DATA
2125
c     Compute the horizontal derivative without missing data checking.
2126
c     The derivative is taken from array 'a' in the direction of 'dir',
2127
c     where 'dir' is either 'X','Y'. The result is stored in array 'd'.
2128
c     The routine accounts for derivatives at the pole and periodic
2129
c     boundaries in the longitudinal direction (depending on
2130
c     the value of datmin, datmax). If the data-set does not reach to
2131
c     the pole, a one-sided derivative is taken. Pole-treatment is only
2132
c     carried out if the data-set covers 360 deg in longitude, and it
2133
c     requires that ie=4*ii+1, where ii is an integer.
2134
c-----------------------------------------------------------------------
2135
 
2136
c     declaration of arguments
2137
      integer       ie,je,ke
2138
      real          a(ie,je,ke),d(ie,je,ke),cl(ie,je)
2139
      real          datmin(4),datmax(4)
2140
      character*(*) dir
2141
 
2142
c     local variable declaration
2143
      integer       i,j,k,ip1,im1,jp1,jm1,ip,im,j1,j2
2144
      real          dlat,dlon,coslat,dx,dy,dxr,dyr
2145
      integer       northpl,southpl,lonper
2146
 
2147
c     rerd and circ are the mean radius and diameter of the earth in
2148
c     meter
2149
      real          rerd,circ,pi
2150
      data          rerd,circ,pi /6.37e6,4.e7,3.141592654/
2151
 
2152
c     compute flags for pole and periodic treatment
2153
      southpl=0
2154
      northpl=0
2155
      lonper =0
2156
      j1=1
2157
      j2=je
2158
      if (abs(datmax(1)-datmin(1)-360.).lt.1.e-3) then
2159
         lonper=1
2160
         if (abs(datmin(2)+90.).lt.1.e-3) then
2161
            southpl=1
2162
            j1=2
2163
         endif
2164
         if (abs(datmax(2)-90.).lt.1.e-3) then
2165
            northpl=1
2166
            j2=je-1
2167
         endif
2168
      endif
2169
 
2170
      dlon=((datmax(1)-datmin(1))/float(ie-1)) *pi/180.
2171
      dlat=((datmax(2)-datmin(2))/float(je-1)) *pi/180.
2172
 
2173
c     print *,'Computing derivative ',dir(1:1),
2174
c     &        ' of an array dimensioned ',ie,je,ke
2175
 
2176
      if (dir(1:1).eq.'X') then
2177
c     -----------------------------
2178
c     derivation in the x-direction
2179
c     -----------------------------
2180
         do k=1,ke
2181
 
2182
c     do gridpoints at j1<=j<=j2
2183
            do j=j1,j2
2184
               coslat=cl(1,j)
2185
 
2186
c     do regular gridpoints at 1<i<ie, 1<j<je
2187
               dx =rerd*coslat*dlon
2188
               dxr=1./(2.*dx)
2189
               do i=2,ie-1
2190
                  ip1=i+1
2191
                  im1=i-1
2192
                  d(i,j,k)=dxr*(a(ip1,j,k)-a(im1,j,k))
2193
               enddo            ! i-loop
2194
c     completed regular gridpoints at 1<i<ie, 1<j<je
2195
 
2196
c     do gridpoints at i=1, i=ie, 1<j<je
2197
               if (lonper.eq.1) then
2198
c     use periodic boundaries
2199
                  i=1
2200
                  ip1=2
2201
                  im1=ie-1
2202
                  d(i,j,k)=dxr*(a(ip1,j,k)-a(im1,j,k))
2203
                  d(ie,j,k)=d(1,j,k)
2204
               else
2205
c     use one-sided derivatives
2206
                  dxr=1./dx
2207
                  i=1
2208
                  ip1=2
2209
                  im1=1
2210
                  d(i,j,k)=dxr*(a(ip1,j,k)-a(im1,j,k))
2211
                  i=ie
2212
                  ip1=ie
2213
                  im1=ie-1
2214
                  d(i,j,k)=dxr*(a(ip1,j,k)-a(im1,j,k))
2215
               endif
2216
c     completed gridpoints at i=1, i=ie, j1<=j<=j2
2217
 
2218
            enddo               ! j-loop
2219
c     completed gridpoints at 1<j<je
2220
 
2221
c     do gridpoints at j=je
2222
            if (northpl.eq.1) then
2223
c     for these gridpoints, the derivative in the x-direction is a
2224
c     derivative in the y-direction at another pole-gridpoint
2225
               dy =rerd*dlat
2226
               dyr=1./(2.*dy)
2227
               j=je
2228
               jp1=je-1
2229
               jm1=je-1
2230
               do i=1,ie
2231
                  ip=mod(i-1+  (ie-1)/4,ie)+1
2232
                  im=mod(i-1+3*(ie-1)/4,ie)+1
2233
                  d(i,j,k)=dyr*(a(ip,jp1,k)-a(im,jm1,k))
2234
               enddo            ! i-loop
2235
c     completed gridpoints at j=je
2236
            endif
2237
c     do gridpoints at j=1
2238
            if (southpl.eq.1) then
2239
               dy =rerd*dlat
2240
               dyr=1./(2.*dy)
2241
               j=1
2242
               jp1=2
2243
               jm1=2
2244
               do i=1,ie
2245
                  ip=mod(i-1+  (ie-1)/4,ie)+1
2246
                  im=mod(i-1+3*(ie-1)/4,ie)+1
2247
                  d(i,j,k)=dyr*(a(ip,jp1,k)-a(im,jm1,k))
2248
               enddo            ! i-loop
2249
            endif
2250
c     completed gridpoints at j=1
2251
 
2252
         enddo                  ! k-loop
2253
 
2254
      else if (dir(1:1).eq.'Y') then
2255
c     -----------------------------
2256
c     derivation in the y-direction
2257
c     -----------------------------
2258
         dy =dlat*rerd
2259
         dyr=1./(2.*dy)
2260
         do k=1,ke
2261
            do i=1,ie
2262
 
2263
c     do regular gridpoints
2264
               do j=2,je-1
2265
                  jp1=j+1
2266
                  jm1=j-1
2267
                  d(i,j,k)=dyr*(a(i,jp1,k)-a(i,jm1,k))
2268
               enddo
2269
 
2270
c     do gridpoints at j=je
2271
               if (northpl.eq.1) then
2272
c     pole-treatment
2273
                  j=je
2274
                  jm1=j-1
2275
                  jp1=j-1
2276
                  ip=mod(i-1+(ie-1)/2,ie)+1
2277
                  im=i
2278
                  d(i,j,k)=dyr*(a(ip,jp1,k)-a(im,jm1,k))
2279
               else
2280
c     one-sided derivative
2281
                  j=je
2282
                  jm1=j-1
2283
                  jp1=j
2284
                  d(i,j,k)=2.*dyr*(a(i,jp1,k)-a(i,jm1,k))
2285
               endif
2286
c     completed gridpoints at j=je
2287
 
2288
c     do gridpoints at j=1
2289
               if (southpl.eq.1) then
2290
c     pole-treatment
2291
                  j=1
2292
                  jm1=2
2293
                  jp1=2
2294
                  ip=i
2295
                  im=mod(i-1+(ie-1)/2,ie)+1
2296
                  d(i,j,k)=dyr*(a(ip,jp1,k)-a(im,jm1,k))
2297
               else
2298
c     one-sided derivative
2299
                  j=1
2300
                  jm1=1
2301
                  jp1=2
2302
                  d(i,j,k)=2.*dyr*(a(i,jp1,k)-a(i,jm1,k))
2303
               endif
2304
c     completed gridpoints at j=1
2305
 
2306
            enddo
2307
         enddo
2308
 
2309
      endif
2310
      end
2311
 
2312
      subroutine ddh2m(a,d,cl,dir,ie,je,ke,datmin,datmax,mdv)
2313
c-----------------------------------------------------------------------
2314
c     Purpose: HORIZONTAL DERIVATIVE ON DATA-SURFACES WITH MISSING
2315
C     DATA
2316
c     Compute the horizontal derivative with missing data checking.
2317
c     The derivative is taken from array 'a' in the direction of 'dir',
2318
c     where 'dir' is either 'X','Y'. The result is stored in array 'd'.
2319
c     The routine accounts for derivatives at the pole and periodic
2320
c     boundaries in the longitudinal direction (depending on
2321
c     the value of datmin, datmax). If the data-set does not reach to
2322
c     the pole, a one-sided derivative is taken. Pole-treatment is only
2323
c     carried out if the data-set covers 360 deg in longitude, and it
2324
c     requires that ie=4*ii+1, where ii is an integer.
2325
c     mdv is the misdat parameter (real)
2326
c-----------------------------------------------------------------------
2327
 
2328
c     declaration of arguments
2329
      integer       ie,je,ke
2330
      real          a(ie,je,ke),d(ie,je,ke),cl(ie,je)
2331
      real          datmin(4),datmax(4),mdv
2332
      character*(*) dir
2333
 
2334
c     local variable declaration
2335
      integer       i,j,k,ip1,im1,jp1,jm1,ip,im,j1,j2
2336
      real          dlat,dlon,coslat,dx,dy,dxr,dyr
2337
      integer       northpl,southpl,lonper
2338
 
2339
c     rerd and circ are the mean radius and diameter of the earth in
2340
c     meter
2341
      real          rerd,circ,pi
2342
      data          rerd,circ,pi /6.37e6,4.e7,3.141592654/
2343
 
2344
c     compute flags for pole and periodic treatment
2345
      southpl=0
2346
      northpl=0
2347
      lonper =0
2348
      j1=1
2349
      j2=je
2350
      if (abs(datmax(1)-datmin(1)-360.).lt.1.e-3) then
2351
         lonper=1
2352
         if (abs(datmin(2)+90.).lt.1.e-3) then
2353
            southpl=1
2354
            j1=2
2355
         endif
2356
         if (abs(datmax(2)-90.).lt.1.e-3) then
2357
            northpl=1
2358
            j2=je-1
2359
         endif
2360
      endif
2361
 
2362
      dlon=((datmax(1)-datmin(1))/float(ie-1)) *pi/180.
2363
      dlat=((datmax(2)-datmin(2))/float(je-1)) *pi/180.
2364
 
2365
c     print *,'Computing derivative ',dir(1:1),
2366
c     &        ' of an array dimensioned ',ie,je,ke
2367
 
2368
      if (dir(1:1).eq.'X') then
2369
c     -----------------------------
2370
c     derivation in the x-direction
2371
c     -----------------------------
2372
         do k=1,ke
2373
 
2374
c     do gridpoints at j1<=j<=j2
2375
            do j=j1,j2
2376
               coslat=cl(1,j)
2377
 
2378
c     do regular gridpoints at 1<i<ie, 1<j<je
2379
               dx =rerd*coslat*dlon
2380
               dxr=1./(2.*dx)
2381
               do i=2,ie-1
2382
                  ip1=i+1
2383
                  im1=i-1
2384
                  if ((a(ip1,j,k).ne.mdv).and.(a(im1,j,k).ne.mdv)) then
2385
                     d(i,j,k)=dxr*(a(ip1,j,k)-a(im1,j,k))
2386
                  else
2387
                     d(i,j,k)=mdv
2388
                  endif
2389
               enddo            ! i-loop
2390
c     completed regular gridpoints at 1<i<ie, 1<j<je
2391
 
2392
c     do gridpoints at i=1, i=ie, 1<j<je
2393
               if (lonper.eq.1) then
2394
c     use periodic boundaries
2395
                  i=1
2396
                  ip1=2
2397
                  im1=ie-1
2398
                  if ((a(ip1,j,k).ne.mdv).and.(a(im1,j,k).ne.mdv)) then
2399
                     d(i,j,k)=dxr*(a(ip1,j,k)-a(im1,j,k))
2400
                  else
2401
                     d(i,j,k)=mdv
2402
                  endif
2403
                  d(ie,j,k)=d(1,j,k)
2404
               else
2405
c     use one-sided derivatives
2406
                  dxr=1./dx
2407
                  i=1
2408
                  ip1=2
2409
                  im1=1
2410
                  if ((a(ip1,j,k).ne.mdv).and.(a(im1,j,k).ne.mdv)) then
2411
                     d(i,j,k)=dxr*(a(ip1,j,k)-a(im1,j,k))
2412
                  else
2413
                     d(i,j,k)=mdv
2414
                  endif
2415
                  i=ie
2416
                  ip1=ie
2417
                  im1=ie-1
2418
                  if ((a(ip1,j,k).ne.mdv).and.(a(im1,j,k).ne.mdv)) then
2419
                     d(i,j,k)=dxr*(a(ip1,j,k)-a(im1,j,k))
2420
                  else
2421
                     d(i,j,k)=mdv
2422
                  endif
2423
               endif
2424
c     completed gridpoints at i=1, i=ie, j1<=j<=j2
2425
 
2426
            enddo               ! j-loop
2427
c     completed gridpoints at 1<j<je
2428
 
2429
c     do gridpoints at j=je
2430
            if (northpl.eq.1) then
2431
c     for these gridpoints, the derivative in the x-direction is a
2432
c     derivative in the y-direction at another pole-gridpoint
2433
               dy =rerd*dlat
2434
               dyr=1./(2.*dy)
2435
               j=je
2436
               jp1=je-1
2437
               jm1=je-1
2438
               do i=1,ie
2439
                  ip=mod(i-1+  (ie-1)/4,ie)+1
2440
                  im=mod(i-1+3*(ie-1)/4,ie)+1
2441
                  if ((a(ip,jp1,k).ne.mdv).and.(a(im,jm1,k).ne.mdv))
2442
     >                 then
2443
                     d(i,j,k)=dyr*(a(ip,jp1,k)-a(im,jm1,k))
2444
                  else
2445
                     d(i,j,k)=mdv
2446
                  endif
2447
               enddo            ! i-loop
2448
c     completed gridpoints at j=je
2449
            endif
2450
c     do gridpoints at j=1
2451
            if (southpl.eq.1) then
2452
               dy =rerd*dlat
2453
               dyr=1./(2.*dy)
2454
               j=1
2455
               jp1=2
2456
               jm1=2
2457
               do i=1,ie
2458
                  ip=mod(i-1+  (ie-1)/4,ie)+1
2459
                  im=mod(i-1+3*(ie-1)/4,ie)+1
2460
                  if ((a(ip,jp1,k).ne.mdv).and.(a(im,jm1,k).ne.mdv))
2461
     >                 then
2462
                     d(i,j,k)=dyr*(a(ip,jp1,k)-a(im,jm1,k))
2463
                  else
2464
                     d(i,j,k)=mdv
2465
                  endif
2466
               enddo            ! i-loop
2467
            endif
2468
c     completed gridpoints at j=1
2469
 
2470
         enddo                  ! k-loop
2471
 
2472
      else if (dir(1:1).eq.'Y') then
2473
c     -----------------------------
2474
c     derivation in the y-direction
2475
c     -----------------------------
2476
         dy =dlat*rerd
2477
         dyr=1./(2.*dy)
2478
         do k=1,ke
2479
            do i=1,ie
2480
 
2481
c     do regular gridpoints
2482
               do j=2,je-1
2483
                  jp1=j+1
2484
                  jm1=j-1
2485
                  if ((a(i,jp1,k).ne.mdv).and.(a(i,jm1,k).ne.mdv))
2486
     >                 then
2487
                     d(i,j,k)=dyr*(a(i,jp1,k)-a(i,jm1,k))
2488
                  else
2489
                     d(i,j,k)=mdv
2490
                  endif
2491
               enddo
2492
 
2493
c     do gridpoints at j=je
2494
               if (northpl.eq.1) then
2495
c     pole-treatment
2496
                  j=je
2497
                  jm1=j-1
2498
                  jp1=j-1
2499
                  ip=mod(i-1+(ie-1)/2,ie)+1
2500
                  im=i
2501
                  if ((a(ip,jp1,k).ne.mdv).and.(a(im,jm1,k).ne.mdv))
2502
     >                 then
2503
                     d(i,j,k)=dyr*(a(ip,jp1,k)-a(im,jm1,k))
2504
                  else
2505
                     d(i,j,k)=mdv
2506
                  endif
2507
               else
2508
c     one-sided derivative
2509
                  j=je
2510
                  jm1=j-1
2511
                  jp1=j
2512
                  if ((a(i,jp1,k).ne.mdv).and.(a(i,jm1,k).ne.mdv))
2513
     >                 then
2514
                     d(i,j,k)=2.*dyr*(a(i,jp1,k)-a(i,jm1,k))
2515
                  else
2516
                     d(i,j,k)=mdv
2517
                  endif
2518
               endif
2519
c     completed gridpoints at j=je
2520
 
2521
c     do gridpoints at j=1
2522
               if (southpl.eq.1) then
2523
c     pole-treatment
2524
                  j=1
2525
                  jm1=2
2526
                  jp1=2
2527
                  ip=i
2528
                  im=mod(i-1+(ie-1)/2,ie)+1
2529
                  if ((a(ip,jp1,k).ne.mdv).and.(a(im,jm1,k).ne.mdv))
2530
     >                 then
2531
                     d(i,j,k)=dyr*(a(ip,jp1,k)-a(im,jm1,k))
2532
                  else
2533
                     d(i,j,k)=mdv
2534
                  endif
2535
               else
2536
c     one-sided derivative
2537
                  j=1
2538
                  jm1=1
2539
                  jp1=2
2540
                  if ((a(i,jp1,k).ne.mdv).and.(a(i,jm1,k).ne.mdv))
2541
     >                 then
2542
                     d(i,j,k)=2.*dyr*(a(i,jp1,k)-a(i,jm1,k))
2543
                  else
2544
                     d(i,j,k)=mdv
2545
                  endif
2546
               endif
2547
c     completed gridpoints at j=1
2548
 
2549
            enddo
2550
         enddo
2551
 
2552
      endif
2553
      end
2554
 
2555
      real function cosd(arg)
2556
c-----------------------------------------------------------------------
2557
c     compute Cos of an argument in Degree instead of Radian
2558
c-----------------------------------------------------------------------
2559
      real,intent(IN) :: arg
2560
      real cos
2561
      real,parameter :: grad2rad=3.1415926/360.
2562
      cosd=cos(arg*grad2rad)
2563
      return
2564
      end
2565
 
2566
      subroutine pressure(pr,sp,stag3,ie,je,ke,aklev,bklev,aklay,bklay)
2567
c     =================================================================
2568
c     argument declaration
2569
      integer  ie,je,ke
2570
      real,intent(OUT) :: pr(ie,je,ke)
2571
      real,intent(IN)  :: sp(ie,je),stag3
2572
      real,intent(IN)  :: aklev(ke),bklev(ke),aklay(ke),bklay(ke)
2573
 
2574
c     variable declaration
2575
      integer  i,j,k
2576
      real     psrf
2577
 
2578
c     statement-functions for the computation of pressure
2579
      real      prlev,prlay
2580
      integer   is
2581
      prlev(is)=aklev(is)+bklev(is)*psrf
2582
      prlay(is)=aklay(is)+bklay(is)*psrf
2583
 
2584
c     computation pressure
2585
      do i=1,ie
2586
         do j=1,je
2587
            psrf=sp(i,j)
2588
            do k=1,ke
2589
               if (stag3.eq.0.) then
2590
                  pr(i,j,k)=prlev(k)
2591
               else
2592
                  pr(i,j,k)=prlay(k)
2593
               endif
2594
            enddo
2595
         enddo
2596
      enddo
2597
      end
2598
 
2599
      subroutine pres(pr,sp,ie,je,ke,ak,bk)
2600
c     =====================================
2601
c     argument declaration
2602
      integer  ie,je,ke
2603
      real,intent(OUT) :: pr(ie,je,ke)
2604
      real,intent(IN)  :: sp(ie,je)
2605
      real,intent(IN)  :: ak(ke),bk(ke)
2606
 
2607
c     variable declaration
2608
      integer  i,j,k
2609
 
2610
c     computation pressure
2611
      do i=1,ie
2612
      do j=1,je
2613
      do k=1,ke
2614
        pr(i,j,k)=ak(k)+bk(k)*sp(i,j)
2615
      enddo
2616
      enddo
2617
      enddo
2618
      end
2619
 
2620
      subroutine coslat(cl,pollon,pollat,vmin,dx,dy,ie,je)
2621
c     ====================================================
2622
c     argument declaration
2623
      integer  ie,je
2624
      real,intent(OUT) :: cl(ie,je)
2625
 
2626
c     variable declaration
2627
      real	pollon,pollat,vmin(3),vmax(3)
2628
      real	lon,lat,rlon,rlat
2629
      integer	i,j
2630
 
2631
      real	phstoph
2632
 
2633
      real	pi
2634
      data	pi /3.141592654/
2635
 
2636
C     Calculate cos(latitude) array and the coriolis parameter
2637
 
2638
      if ((pollon.ne.0.).or.(pollat.ne.90.)) then
2639
         do j=1,je
2640
            rlat=vmin(2)+(j-1)*dy
2641
            do i=1,ie
2642
               rlon=vmin(1)+(i-1)*dx
2643
               yphys=phstoph(rlat,rlon,pollat,pollon)
2644
C              if I use sind(lat in deg): troubles at the N-pole
2645
               cl(i,j)=cos(rlat*pi/180.)
2646
            enddo
2647
         enddo
2648
      else
2649
         do j=1,je
2650
            lat=vmin(2)+(j-1)*dy
2651
            lat=pi*lat/180.
2652
            do i=1,ie
2653
               cl(i,j)=cos(lat)
2654
            enddo
2655
         enddo
2656
      endif
2657
 
2658
      return
2659
      end
2660
 
2661
      subroutine corpar(f,pollon,pollat,vmin,dx,dy,ie,je)
2662
c     ====================================================
2663
c     argument declaration
2664
      integer  ie,je
2665
      real,intent(OUT) :: f(ie,je)
2666
 
2667
c     variable declaration
2668
      real      pollon,pollat,vmin(3),vmax(3)
2669
      real      lon,lat,rlon,rlat
2670
      integer   i,j
2671
 
2672
      real      phstoph
2673
 
2674
      real      pi
2675
      data      pi /3.141592654/
2676
 
2677
C     Calculate cos(latitude) array and the coriolis parameter
2678
 
2679
      if ((pollon.ne.0.).or.(pollat.ne.90.)) then
2680
         do j=1,je
2681
            rlat=vmin(2)+(j-1)*dy
2682
            do i=1,ie
2683
               rlon=vmin(1)+(i-1)*dx
2684
               yphys=phstoph(rlat,rlon,pollat,pollon)
2685
C              if I use sind(lat in deg): troubles at the N-pole
2686
               lat=2.*pi*yphys/360.
2687
               f(i,j)=0.000145444*sin(lat)
2688
            enddo
2689
         enddo
2690
      else
2691
         do j=1,je
2692
            lat=vmin(2)+(j-1)*dy
2693
            lat=pi*lat/180.
2694
            do i=1,ie
2695
               f(i,j)=0.000145444*sin(lat)
2696
            enddo
2697
         enddo
2698
      endif
2699
 
2700
      return
2701
      end
2702
 
2703
      subroutine vint(intfield,field,sp,pmin,pmax,ie,je,ke,ak,bk)
2704
C     ===========================================================
2705
 
2706
c     argument declaration
2707
      integer   ie,je,ke
2708
      real      intfield(ie,je),field(ie,je,ke),sp(ie,je)
2709
      real      ak(ke),bk(ke)
2710
      real      pmin,pmax
2711
 
2712
c     variable declaration
2713
      integer   i,j,k,kmin,kmax
2714
      real      coeff,pval
2715
 
2716
c     statement-functions for the computation of pressure
2717
      real      pp,psrf
2718
      integer   is
2719
      pp(is)=ak(is)+bk(is)*psrf
2720
 
2721
C     ====================================================================
2722
C     ===  begin of main part of this subroutine  ========================
2723
C     ====================================================================
2724
 
2725
      do i=1,ie
2726
      do j=1,je
2727
        intfield(i,j)=0.
2728
        psrf=sp(i,j)
2729
        if (psrf.lt.pmin) then                     
2730
          intfield(i,j)=-999.98999
2731
          goto 556
2732
        endif
2733
        kmin=1
2734
        kmax=1
2735
        do k=2,ke
2736
          if ((pp(k).lt.pmax).and.(pp(k-1).gt.pmax)) kmin=k
2737
          if ((pp(k).lt.pmin).and.(pp(k-1).gt.pmin)) kmax=k-1
2738
        enddo
2739
c       check if field is not equal mdv
2740
        do k=kmin,kmax+1
2741
          if (field(i,j,k).eq.-999.98999) then
2742
            intfield(i,j)=-999.98999
2743
            goto 556
2744
          endif
2745
        enddo
2746
c       pmin and pmax are inbetween same pair of layers
2747
c       interpolate on intermediate value (pmin+pmax)/2.
2748
        if (kmin.eq.kmax+1) then
2749
          pval=(pmin+pmax)/2.
2750
          coeff=(pval-pp(kmin-1))/(pp(kmin)-pp(kmin-1))
2751
          intfield(i,j)=(pmax-pmin)*
2752
     &      (field(i,j,kmin)*coeff+field(i,j,kmax)*(1.-coeff))
2753
          goto 555
2754
        endif
2755
c       one layer is inbetween pmin and pmax
2756
        if (kmin.eq.kmax) then
2757
          if (psrf.lt.pmax) then
2758
            intfield(i,j)=field(i,j,kmin)*(psrf-pmin)
2759
          else
2760
            intfield(i,j)=field(i,j,kmin)*(pmax-pmin)
2761
          endif
2762
          goto 555
2763
        endif
2764
c       loop for the interior levels
2765
        do k=kmin+1,kmax-1
2766
          intfield(i,j)=intfield(i,j)+field(i,j,k)*
2767
     &              0.5*(pp(k-1)-pp(k+1))
2768
        enddo
2769
c       special treatment of the bounding levels
2770
        if (kmin.eq.1) then
2771
          if (psrf.lt.pmax) then
2772
            intfield(i,j)=intfield(i,j)+
2773
     &          field(i,j,1)*(0.5*(pp(1)-pp(2))+psrf-pp(1))
2774
          else
2775
            intfield(i,j)=intfield(i,j)+
2776
     &          field(i,j,1)*(0.5*(pp(1)-pp(2))+pmax-pp(1))
2777
          endif
2778
        else
2779
          coeff=(pmax-pp(kmin-1))/(pp(kmin)-pp(kmin-1))
2780
            intfield(i,j)=intfield(i,j)+
2781
     &      field(i,j,kmin)*0.5*(pmax-pp(kmin+1))+
2782
     &      (field(i,j,kmin)*coeff+field(i,j,kmin-1)*(1.-coeff))*
2783
     &      0.5*(pmax-pp(kmin))
2784
        endif
2785
        if (kmax.eq.ke) then          
2786
          intfield(i,j)=intfield(i,j)+
2787
     &           field(i,j,ke)*0.5*(pp(ke-1)-pp(ke))
2788
        else
2789
          coeff=(pmin-pp(kmax+1))/(pp(kmax)-pp(kmax+1))
2790
          intfield(i,j)=intfield(i,j)+
2791
     &      field(i,j,kmax)*0.5*(pp(kmax-1)-pmin)+
2792
     &      (field(i,j,kmax)*coeff+field(i,j,kmax+1)*(1.-coeff))*
2793
     &      0.5*(pp(kmax)-pmin)
2794
        endif
2795
  555   continue
2796
C       Calculate mean value
2797
        if (psrf.lt.pmax) then
2798
          intfield(i,j)=intfield(i,j)/(psrf-pmin)
2799
        else
2800
          intfield(i,j)=intfield(i,j)/(pmax-pmin)
2801
        endif
2802
  556   continue
2803
      enddo
2804
      enddo
2805
 
2806
      return
2807
      end
2808
 
2809
      subroutine nvint(intfield,field,sp,pmi,pma,
2810
     >                 ie,je,ke,ak,bk,iflag)
2811
C     ===========================================
2812
c
2813
c     iflag     input   =0 average (like SR vint); =1 integration
2814
 
2815
c     argument declaration
2816
      integer   ie,je,ke,iflag
2817
      real      intfield(ie,je),field(ie,je,ke),sp(ie,je)
2818
      real      ak(ke),bk(ke)
2819
      real      pmin,pmax,pmi,pma
2820
 
2821
c     variable declaration
2822
      integer   ilev
2823
      integer   i,j,k,kmin,kmax,kshift
2824
      real      mdv,coeff,pval
2825
 
2826
c     statement-functions for the computation of pressure
2827
      real      pp,psrf
2828
      integer   is
2829
      pp(is)=ak(is)+bk(is)*psrf
2830
 
2831
C     ====================================================================
2832
C     ===  begin of main part of this subroutine  ========================
2833
C     ====================================================================
2834
 
2835
C     set misdat value
2836
      mdv=-999.98999
2837
 
2838
C     determine whether data is on model levels (ilev=0)
2839
C     or pressure levels (ilev=1)
2840
      ilev=1
2841
      do k=1,ke
2842
        if (bk(k).ne.0.) ilev=0
2843
      enddo
2844
 
2845
      do i=1,ie
2846
      do j=1,je
2847
        intfield(i,j)=0.
2848
        psrf=sp(i,j)
2849
        kmin=1
2850
        kmax=1
2851
        pmin=pmi
2852
        pmax=pma
2853
        if (ilev.eq.0) then
2854
          if (psrf.lt.pmin) then
2855
            intfield(i,j)=mdv
2856
            goto 456
2857
          endif
2858
        endif
2859
        do k=2,ke
2860
          if ((pp(k).le.pmax).and.(pp(k-1).gt.pmax)) kmin=k
2861
          if ((pp(k).le.pmin).and.(pp(k-1).gt.pmin)) kmax=k-1
2862
        enddo
2863
c       if model levels: set vint=mdv if one level is equal mdv
2864
        if (ilev.eq.0) then
2865
          do k=kmin,kmax+1
2866
            if (field(i,j,k).eq.mdv) then
2867
              intfield(i,j)=mdv
2868
              goto 456
2869
            endif
2870
          enddo
2871
        endif
2872
c       if pressure levels: check whether lowest levels are not mdv
2873
c       if they are: adjust pmax
2874
        if (ilev.eq.1) then
2875
          kshift=0
2876
          do k=kmin,kmax
2877
            if (field(i,j,k).eq.mdv) then
2878
              kshift=kshift+1
2879
              pmax=0.5*(pp(k)+pp(k+1))
2880
            endif
2881
          enddo
2882
c         if mdv occur go even one level higher
2883
          if (kshift.gt.0) kshift=kshift+1
2884
          kmin=kmin+kshift
2885
          if (kmin.gt.kmax+1) then
2886
            intfield(i,j)=mdv
2887
            goto 456
2888
          endif
2889
        endif
2890
c       pmin and pmax are inbetween same pair of layers
2891
c       interpolate on intermediate value (pmin+pmax)/2.
2892
        if (kmin.eq.kmax+1) then
2893
          pval=(pmin+pmax)/2.
2894
          coeff=(pval-pp(kmin-1))/(pp(kmin)-pp(kmin-1))
2895
          intfield(i,j)=(pmax-pmin)*
2896
     &      (field(i,j,kmin)*coeff+field(i,j,kmax)*(1.-coeff))
2897
          goto 455
2898
        endif
2899
c       one layer is inbetween pmin and pmax
2900
        if (kmin.eq.kmax) then
2901
          if (psrf.lt.pmax) then
2902
            intfield(i,j)=field(i,j,kmin)*(psrf-pmin)
2903
          else
2904
            intfield(i,j)=field(i,j,kmin)*(pmax-pmin)
2905
          endif
2906
          goto 455
2907
        endif
2908
c       loop for the interior levels
2909
        do k=kmin+1,kmax-1
2910
          intfield(i,j)=intfield(i,j)+field(i,j,k)*
2911
     &              0.5*(pp(k-1)-pp(k+1))
2912
        enddo
2913
c       special treatment of the bounding levels
2914
        if (kmin.eq.1) then
2915
          if (psrf.lt.pmax) then
2916
            intfield(i,j)=intfield(i,j)+
2917
     &          field(i,j,1)*(0.5*(pp(1)-pp(2))+psrf-pp(1))
2918
          else
2919
            intfield(i,j)=intfield(i,j)+
2920
     &          field(i,j,1)*(0.5*(pp(1)-pp(2))+pmax-pp(1))
2921
          endif
2922
        else
2923
          coeff=(pmax-pp(kmin-1))/(pp(kmin)-pp(kmin-1))
2924
            intfield(i,j)=intfield(i,j)+
2925
     &      field(i,j,kmin)*0.5*(pmax-pp(kmin+1))+
2926
     &      (field(i,j,kmin)*coeff+field(i,j,kmin-1)*(1.-coeff))*
2927
     &      0.5*(pmax-pp(kmin))
2928
        endif
2929
        if (kmax.eq.ke) then
2930
          intfield(i,j)=intfield(i,j)+
2931
     &           field(i,j,ke)*0.5*(pp(ke-1)-pp(ke))
2932
        else
2933
          coeff=(pmin-pp(kmax+1))/(pp(kmax)-pp(kmax+1))
2934
          intfield(i,j)=intfield(i,j)+
2935
     &      field(i,j,kmax)*0.5*(pp(kmax-1)-pmin)+
2936
     &      (field(i,j,kmax)*coeff+field(i,j,kmax+1)*(1.-coeff))*
2937
     &      0.5*(pp(kmax)-pmin)
2938
        endif
2939
  455   continue
2940
C       Calculate mean or integrated value
2941
        if (iflag.eq.0) then
2942
          if ((ilev.eq.0).and.(psrf.lt.pmax)) then
2943
            intfield(i,j)=intfield(i,j)/(psrf-pmin)
2944
          else
2945
            intfield(i,j)=intfield(i,j)/(pmax-pmin)
2946
          endif
2947
        else
2948
          intfield(i,j)=100./9.8*intfield(i,j)
2949
        endif
2950
  456   continue
2951
      enddo
2952
      enddo
2953
 
2954
      return
2955
      end
2956
 
2957
C     =================================================================
2958
      subroutine z2hno3(z,nmbl,zval,prof,hno3)
2959
C     =================================================================
2960
C     Input is a hno3-profile ("prof") in function of "nmbl" Z-levels
2961
C     stored in zval.
2962
C     Output is the "hno3"-value on this "z"
2963
 
2964
      implicit none
2965
      integer       nmbl
2966
 
2967
C     Argument declaration
2968
      real      zval(nmbl),prof(nmbl),hno3,z
2969
 
2970
C     Further variables declaration
2971
      real      re
2972
      integer   i
2973
 
2974
C     Check if z is inside profile values
2975
      if (z.gt.zval(nmbl)) then
2976
        print*,'*** error: z is over the range of HNO3-profile:'
2977
        print*,z,zval(nmbl)
2978
        stop
2979
      endif
2980
 
2981
      if (z.lt.zval(1)) then
2982
        print*,'*** error: z is below the range of HNO3-profile:'
2983
        print*,z,zval(1)
2984
        stop
2985
      endif
2986
 
2987
C     Interpolate z of table to find the HNO3 value
2988
      do i=1,nmbl-1
2989
        if (z.ge.zval(i) .and. z.lt.zval(i+1)) then
2990
          re=(z-zval(i))/(zval(i+1)-zval(i))
2991
          goto 22
2992
        endif
2993
        if (i.eq.nmbl-1) then
2994
          print*,'***Error: z:',z,' outside HNO3-Profile...'
2995
          stop
2996
        endif
2997
      enddo
2998
 
2999
 22   hno3=prof(i)+re*(prof(i+1)-prof(i))
3000
*     print*,zi,zval(i),zval(i+1),re
3001
      return
3002
 
3003
      end
3004
 
3005
C     =================================================================
3006
      subroutine p2h2o(p,xih2o)
3007
C     =================================================================
3008
 
3009
      implicit none
3010
 
3011
      integer      nn
3012
      parameter    (nn=12)
3013
      real         wprof(nn),pprof(nn)
3014
 
3015
c     variables declaration
3016
      real     p,xih2o,k
3017
      integer  n
3018
      data     pprof/150.,130.,100.,90.,80.,70.,60.,50.,40.,30.,20.,10./
3019
      data     wprof/4.57,4.4,4.25,4.28,4.38,4.53,4.7,4.97,5.4,
3020
     >     5.9,6.6,6.85/
3021
 
3022
      do n=1,nn-1
3023
         if (p.le.pprof(n).and.p.ge.pprof(n+1)) then
3024
            k=(p-pprof(n))/(pprof(n+1)-pprof(n))
3025
            goto 100
3026
         endif
3027
         if (n.eq.nn-1) then
3028
           print*,p,pprof(1),pprof(nn)
3029
           stop'Table too small'
3030
         endif
3031
      enddo
3032
 
3033
 100  continue
3034
 
3035
      xih2o=wprof(n)+k*(wprof(n+1)-wprof(n))
3036
      return
3037
      end
3038
 
3039
C     =================================================================
3040
      subroutine sg(temp,pres,rad,as,xHNO3,xH2O,sedi,growth)
3041
C     =================================================================
3042
C     Input of sg are a Temperature[K], a pressure[hPa], a Radius[um],
3043
C     a particle shape specification stored in "as" with s for
3044
C     spherical particles, a for aspherical a mixing ratio of HNO3.
3045
C     Output is the sedimentation factor and a growth rate for the
3046
C     specified particle.
3047
      implicit none
3048
 
3049
C     Constants declaration
3050
      real       deltarho, g, cunn_a, cunn_b, cunn_c
3051
      real       ff1, ff2
3052
      real       gas_diffusivity_corr, mb2mm
3053
      real       aspect_ratio, d2, k, pi, R
3054
 
3055
      parameter  (g=9.81,  deltarho=1.62e3)   ! deltarho=(NAT Density-Air Density)
3056
      parameter  (cunn_a=1.257, cunn_b=0.4, cunn_c=1.1)!cunn stays for Cunningham...
3057
      parameter  (ff1=1.12, ff2=0.58248) ! Form factors for asphericity 1:3
3058
      parameter  (mb2mm=760./1000.)
3059
C     For aspher. part., changing the next value implies changing others as formfactors
3060
      parameter  (aspect_ratio=1./3.)
3061
      parameter  (d2=3.e-10, k=1.380661e-23, pi=3.1415927, R=8.31441)
3062
 
3063
C     Arguments declaration
3064
      double precision sedi,growth
3065
      real             temp,pres,rad,xHNO3,xH2O
3066
      character*(1)    as
3067
 
3068
C     Functions declaration
3069
      double precision cunn_corr,grown
3070
 
3071
C     Further variables
3072
      double precision cc,wStokes,lamda,knudsen,capacity,deltapress
3073
      real       eta,p_HNO3,p_H2O
3074
 
3075
      eta = 6.45e-8 * temp
3076
 
3077
      wStokes = 2.* rad**2 * g * deltarho / (9. * eta)
3078
 
3079
      lamda= k * temp / ( 2.**0.5 * 100. * pres * pi * d2**2 )
3080
 
3081
      knudsen= lamda / rad
3082
 
3083
      if (as.eq.'s') then
3084
         cc=cunn_corr(knudsen,1.,1.,cunn_a,cunn_b,cunn_c)
3085
         capacity=1.
3086
      else
3087
         cc=cunn_corr(knudsen,ff1,ff2,cunn_a,cunn_b,cunn_c)
3088
         capacity=sqrt(1.-aspect_ratio**2) / (aspect_ratio**(2./3.) *
3089
     >            log( (1. + sqrt(1.-aspect_ratio**2)) / aspect_ratio) )
3090
      endif
3091
 
3092
      sedi= wStokes * cc / rad**2
3093
 
3094
      p_H2O= pres * xH2O * 1.e-6
3095
      p_HNO3= pres * xHNO3 * 1.e-9
3096
      deltapress= p_HNO3 - ( exp( (-2.7836 - 0.00088 * temp) *
3097
     >            log (p_H2O * mb2mm) + 89.7674 - 26242. / temp +
3098
     >            0.021135 * temp ) / mb2mm )
3099
 
3100
      growth= capacity * grown(temp,pres,rad,sedi,deltapress,eta)
3101
C      print*,'T,p,r,shape,xHNO3 -> sedi,growth,',temp,pres,rad,as,
3102
C     >        xHNO3,sedi,growth
3103
      end
3104
 
3105
C     =================================================================
3106
      double precision function grown(temp,pres,rad,sedi,deltapress,eta)
3107
C     =================================================================
3108
      implicit none
3109
 
3110
      real      molmass_g,molmass_s,rho_s
3111
      real      p0,T0,R,pi,gas_diffusivity_corr
3112
 
3113
      parameter (molmass_g=63.e-3 ,molmass_s=117.e-3 ,rho_s=1.62e3 )
3114
      parameter (p0=1013.25,T0=273.15,pi=3.1415927,R=8.3144)
3115
      parameter  (gas_diffusivity_corr=1.)
3116
 
3117
      double precision sedi,deltapress,gasdiffu,gasdiffcorr,prod
3118
      real      temp,pres,rad,eta
3119
      real      meanvel,Reynold,Schmidt,ventil
3120
 
3121
C      if ( temp.lt.233 .or. temp.gt.313 ) print*,'*** Warning:
3122
C     >Temperature',temp,'is out of validity for diffusivity'
3123
 
3124
      gasdiffu= sqrt( 0.018 / molmass_g ) * 0.211e-4 *
3125
     >                   ( temp / T0 )**1.94 * (P0 / pres)
3126
      meanvel= sqrt( 8.* R * temp / (molmass_g * pi) )
3127
      gasdiffcorr= gasdiffu / ( 1. + 4. * gasdiffu /
3128
     >            ( gas_diffusivity_corr * rad * meanvel ) )
3129
      Reynold=  sedi * rad**2 * 2. * rad / eta
3130
      Schmidt= eta / gasdiffu
3131
      prod= Schmidt**(1./3.) * sqrt(Reynold)
3132
 
3133
      if ( prod.lt.1.4) then
3134
         ventil= 1. + 0.108 * prod**2
3135
      else
3136
         ventil= 0.78 + 0.308 * prod
3137
      endif
3138
      grown= gasdiffcorr * ventil * molmass_s * deltapress * 100. /
3139
     >       ( rho_s * R * temp )
3140
      return
3141
      end
3142
 
3143
C     =================================================================
3144
      double precision function cunn_corr(knudsen,ff1,ff2,a,b,c)
3145
C     =================================================================
3146
      double precision knudsen
3147
      real      ff1,ff2,a,b,c
3148
 
3149
      cunn_corr=(1. + (ff2/ff1) * knudsen * (a+b*exp(-c/knudsen)))/ff1
3150
      return
3151
      end