Subversion Repositories lagranto.icon

Rev

Details | Last modification | View Log | RSS feed

Rev Author Line No. Line
3 michaesp 1
c     -------------------------------------------------------------------------
2
c     Potential temperature (TH)
3
c     -------------------------------------------------------------------------
4
 
5
      subroutine calc_TH (pt,t,p)
6
 
7
      implicit none
8
 
9
c     Argument declaration
10
      real  pt     ! Potential temperature [K]
11
      real  t      ! Temperature [either in C or in K]
12
      real  p      ! Pressure [hPa]
13
 
14
c     Physical parameters
15
      real      rdcp,tzero
16
      data      rdcp,tzero /0.286,273.15/
17
 
18
c     Calculation - distinction between temperature in C or in K
19
      if (t.lt.100.) then
20
         pt = (t+tzero) * ( (1000./p)**rdcp ) 
21
      else
22
         pt = t * ( (1000./p)**rdcp )
23
      endif
24
 
25
      end
26
 
27
c     -------------------------------------------------------------------------
28
c     Density (RHO)
29
c     -------------------------------------------------------------------------
30
 
31
      subroutine calc_RHO (rho,t,p)
32
 
33
      implicit none
34
 
35
c     Argument declaration
36
      real  rho    ! Density [kg/m^3]
37
      real  t      ! Temperature [either in C or in K]
38
      real  p      ! Pressure [hPa]
39
 
40
c     Physical parameters
41
      real      rd,tzero
42
      data      rd,tzero /287.05,273.15/
43
 
44
c     Auxiliary variables
45
      real      tk
46
 
47
c     Calculation - distinction between temperature in C or in K
48
      if (t.lt.100.) then
49
         tk = t + tzero
50
      else
51
         tk = t
52
      endif
53
 
54
      rho = 100.*p/( tk * rd ) 
55
 
56
      end
57
 
58
c     -------------------------------------------------------------------------
59
c     Relative humidity (RH)
60
c     -------------------------------------------------------------------------
61
 
62
      subroutine calc_RH (rh,t,p,q)
63
 
64
      implicit none
65
 
66
c     Argument declaration
67
      real  rh     ! Relative humidity [%]
68
      real  t      ! Temperature [either in C or in K]
69
      real  p      ! Pressure [hPa]
70
      real  q      ! Specific humidity [g/kg]
71
 
72
c     Physical parameters
73
      real      rdcp,tzero
74
      data      rdcp,tzero /0.286,273.15/
75
      real      b1,b2w,b3,b4w,r,rd
76
      data      b1,b2w,b3,b4w,r,rd /6.1078, 17.2693882, 273.16, 35.86,
77
     &                              287.05, 461.51/
78
 
79
c     Auxiliary variables
80
      real      ge
81
      real      gqd
82
      real      tc
83
      real      pp,qk
84
 
85
c     Calculation - distinction between temperature in C or in K
86
      if (t.gt.100.) then
87
         tc = t - tzero
88
      else
89
         tc = t
90
      endif
91
      qk = 0.001 * q
92
 
93
      ge  = b1*exp(b2w*tc/(tc+b3-b4w))
94
      gqd = r/rd*ge/(p-(1.-r/rd)*ge)
95
      rh  = 100.*qk/gqd
96
 
97
      end
98
 
99
c     -------------------------------------------------------------------------
100
c     Equivalent potential temperature (THE)
101
c     -------------------------------------------------------------------------
102
 
103
      subroutine calc_THE (the,t,p,q)
104
 
105
      implicit none
106
 
107
c     Argument declaration
108
      real  the    ! Equivalent potential temperature [K]
109
      real  t      ! Temperature [either in C or in K]
110
      real  p      ! Pressure [hPa]
111
      real  q      ! Specific humidity [g/kg]
112
 
113
c     Physical parameters 
114
      real     rdcp,tzero
115
      data     rdcp,tzero /0.286,273.15/
116
 
117
c     Auxiliary variables 
118
      real     tk,qk
119
 
120
c     Calculation - distinction between temperature in C or in K
121
      if (t.lt.100.) then
122
         tk = t + tzero
123
      else
124
         tk = t
125
      endif
126
      qk = 0.001 * q
127
 
128
      the = tk*(1000./p)
129
     +      **(0.2854*(1.0-0.28*qk))*exp(
130
     +      (3.376/(2840.0/(3.5*alog(tk)-alog(
131
     +      100.*p*max(1.0E-10,qk)/(0.622+0.378*
132
     +      q))-0.1998)+55.0)-0.00254)*1.0E3*
133
     +      max(1.0E-10,qk)*(1.0+0.81*qk))
134
 
135
      end
136
 
137
c     -------------------------------------------------------------------------
138
c     Latent heating rate (LHR)
139
c     -------------------------------------------------------------------------
140
 
141
      subroutine calc_LHR (lhr,t,p,q,w,rh)
142
 
143
      implicit none 
144
 
145
c     Argument declaration
146
      real  lhr    ! Latent heating rate [K/h]
147
      real  t      ! Temperature [either in C or in K]
148
      real  p      ! Pressure [hPa]
149
      real  q      ! Specific humidity [kg/kg]
150
      real  w      ! Vertical velocity [m/s]
151
      real  rh     ! Relative humidity [%] 
152
 
153
c     Physical parameters 
154
      real      p0,kappa,tzero,g
155
      data      p0,kappa,tzero,g /1000.,0.286,273.15,9.81/
156
      real      blog10,cp,r,lw,eps
157
      data      blog10,cp,r,lw,eps /.08006,1004.,287.,2.5e+6,0.622/
158
 
159
c     Auxiliary variables
160
      real      tk
161
      real      qk
162
      real      tt
163
      real      esat,c
164
	  real      omega
165
 
166
c     Calculation - distinction between temperature in C or in K
167
      if (t.lt.100.) then
168
         tk = t + tzero
169
      else
170
         tk = t
171
      endif
172
      qk = 0.001 * q
173
 
174
      omega = -p/(r*tk) * g * w
175
 
176
      if (rh.lt.80.) then 
177
         lhr = 0.    
178
      else if (omega.gt.0.) then
179
         lhr = 0.
180
      else
181
         c   = lw/cp*eps*blog10*esat(tk)/p
182
         tt  = (tk*(p0/p)**kappa) 
183
         lhr = 3600.*
184
     >        (1.-exp(.2*(80.-rh))) 
185
     >        *(-c*kappa*tt*omega/(100.*p))/(1.+c)   
186
 
187
      endif
188
 
189
      end
190
 
191
c     -------------------------------------------------------------------------
192
c     Wind speed (VEL)
193
c     -------------------------------------------------------------------------
194
 
195
      subroutine calc_VEL (vel,u,v)
196
 
197
      implicit none
198
 
199
c     Argument declaration
200
      real  vel    ! Wind speed [m/s]
201
      real  u      ! Zonal wind component [m/s]
202
      real  v      ! Meridional wind component [m/s]
203
 
204
      vel = sqrt ( u*u + v*v )
205
 
206
      end
207
 
208
c     -------------------------------------------------------------------------
209
c     Wind direction (DIR)
210
c     -------------------------------------------------------------------------
211
 
212
      subroutine calc_DIR (dir,u,v)
213
 
214
      implicit none
215
 
216
c     Argument declaration
217
      real  dir    ! Wind direction [deg]
218
      real  u      ! Zonal wind component [m/s]
219
      real  v      ! Meridional wind component [m/s]
220
 
221
      call getangle(1.,0.,u,v,dir)
222
 
223
      end
224
 
225
c     -------------------------------------------------------------------------
226
c     Zonal derivative of U (DUDX)
227
c     -------------------------------------------------------------------------
228
 
229
      subroutine calc_DUDX (dudx,u1,u0,lat)
230
 
231
      implicit none
232
 
233
c     Argument declaration
234
      real  dudx   ! Derivative of U in zonal direction [s^-1]
235
      real  u1     ! U @ LON + 1 DLON [m/s]
236
      real  u0     ! U @ LON - 1 DLON [m/s]
237
      real  lat    ! Latitude [deg]
238
 
239
c     Physical parameters
240
      real      pi180  
241
      parameter (pi180=31.14159/180.)
242
      real      deltay
243
      parameter (deltay =1.11e5)
244
 
245
      dudx = (u1-u0) / ( 2. * deltay * cos(pi180 * lat) )
246
 
247
      end
248
 
249
c     -------------------------------------------------------------------------
250
c     Zonal derivative of V (DVDX)
251
c     -------------------------------------------------------------------------
252
 
253
      subroutine calc_DVDX (dvdx,v1,v0,lat)
254
 
255
c     Argument declaration
256
      real  dvdx   ! Derivative of V in zonal direction [s^-1]
257
      real  v1     ! V @ LON + 1 DLON [m/s]
258
      real  v0     ! V @ LON - 1 DLON [m/s]
259
      real  lat    ! Latitude [deg]
260
 
261
c     Physical parameters
262
      real      pi180  
263
      parameter (pi180=3.14159/180.)
264
      real      deltay
265
      parameter (deltay =1.11e5)
266
 
267
      dvdx = (v1-v0) / ( 2. * deltay * cos(pi180 * lat) )
268
 
269
      end
270
 
271
c     -------------------------------------------------------------------------
272
c     Zonal derivative of T (DTDX)
273
c     -------------------------------------------------------------------------
274
 
275
      subroutine calc_DTDX (dtdx,t1,t0,lat)
276
 
277
      implicit none
278
 
279
c     Argument declaration
280
      real  dtdx   ! Derivative of T in zonal direction [K/m]
281
      real  t1     ! T @ LON + 1 DLON [K]
282
      real  t0     ! T @ LON - 1 DLON [K]
283
      real  lat    ! Latitude [deg]
284
 
285
c     Physical parameters
286
      real      pi180  
287
      parameter (pi180=3.14159/180.)
288
      real      deltay
289
      parameter (deltay =1.11e5)
290
 
291
      dtdx = (t1-t0) / ( 2. * deltay * cos(pi180 * lat) )
292
 
293
      end
294
 
295
c     -------------------------------------------------------------------------
296
c     Zonal derivative of TH (DTHDX)
297
c     -------------------------------------------------------------------------
298
 
299
      subroutine calc_DTHDX (dthdx,t1,t0,p,lat)
300
 
301
      implicit none
302
 
303
c     Argument declaration
304
      real  dthdx  ! Derivative of TH in zonal direction [K/m]
305
      real  t1     ! T @ LON + 1 DLON [K]
306
      real  t0     ! T @ LON - 1 DLON [K]
307
      real  p      ! P [hPa]
308
      real  lat    ! Latitude [deg]
309
 
310
c     Physical parameters
311
      real      pi180  
312
      parameter (pi180=3.14159/180.)
313
      real      deltay
314
      parameter (deltay =1.11e5)
315
      real      rdcp,pref
316
      data      rdcp,pref /0.286,1000./
317
 
318
      dthdx = (pref/p)**rdcp * 
319
     >        (t1-t0) / ( 2. * deltay * cos(pi180 * lat) )
320
 
321
      end
322
 
323
c     -------------------------------------------------------------------------
324
c     Meridional derivative of U (DUDY)
325
c     -------------------------------------------------------------------------
326
 
327
      subroutine calc_DUDY (dudy,u1,u0)
328
 
329
      implicit none
330
 
331
c     Argument declaration
332
      real  dudy   ! Derivative of U in meridional direction [s^-1]
333
      real  u1     ! U @ LAT + 1 DLAT [m/s]
334
      real  u0     ! U @ LAT - 1 DLAT [m/s]
335
 
336
c     Physical parameters
337
      real      deltay
338
      parameter (deltay =1.11e5)
339
 
340
      dudy = (u1-u0) / ( 2. * deltay )
341
 
342
      end
343
 
344
c     -------------------------------------------------------------------------
345
c     Meridional derivative of V (DVDY)
346
c     -------------------------------------------------------------------------
347
 
348
      subroutine calc_DVDY (dvdy,v1,v0)
349
 
350
      implicit none
351
 
352
c     Argument declaration
353
      real  dvdy   ! Derivative of V in meridional direction [s^-1]
354
      real  v1     ! V @ LAT + 1 DLAT [m/s]
355
      real  v0     ! V @ LAT - 1 DLAT [m/s]
356
 
357
c     Physical parameters
358
      real      deltay
359
      parameter (deltay =1.11e5)
360
 
361
      dvdy = (v1-v0) / ( 2. * deltay )
362
 
363
      end
364
 
365
c     -------------------------------------------------------------------------
366
c     Meridional derivative of T (DTDY)
367
c     -------------------------------------------------------------------------
368
 
369
      subroutine calc_DTDY (dtdy,t1,t0)
370
 
371
      implicit none
372
 
373
c     Argument declaration
374
      real  dtdy   ! Derivative of T in meridional direction [K/m]
375
      real  t1     ! T @ LAT + 1 DLAT [K]
376
      real  t0     ! T @ LAT - 1 DLAT [K]
377
 
378
c     Physical parameters
379
      real      deltay
380
      parameter (deltay =1.11e5)
381
 
382
      dtdy = (t1-t0) / ( 2. * deltay )
383
 
384
      end
385
 
386
c     -------------------------------------------------------------------------
387
c     Meridional derivative of TH (DTHDY)
388
c     -------------------------------------------------------------------------
389
 
390
      subroutine calc_DTHDY (dthdy,t1,t0,p)
391
 
392
      implicit none
393
 
394
c     Argument declaration
395
      real  dthdy  ! Derivative of TH in meridional direction [K/m]
396
      real  t1     ! TH @ LAT + 1 DLAT [K]
397
      real  t0     ! TH @ LAT - 1 DLAT [K]
398
      real  p      ! P [hPa]
399
 
400
c     Physical parameters
401
      real      deltay
402
      parameter (deltay =1.11e5)
403
      real      rdcp,pref
404
      data      rdcp,pref /0.286,1000./
405
 
406
      dthdy = (pref/p)**rdcp * (t1-t0) / ( 2. * deltay )
407
 
408
      end
409
 
410
c     -------------------------------------------------------------------------
411
c     Wind shear of U (DUDZ)
412
c     -------------------------------------------------------------------------
413
 
414
      subroutine calc_DUDZ (dudz,u1,u0,z1,z0)
415
 
416
      implicit none
417
 
418
c     Argument declaration
419
      real  dudz   ! Wind shear [m/s per m = 1/m]
420
      real  u1     ! U @ Z + 1 DZ [m/s]
421
      real  u0     ! U @ Z - 1 DZ [m/s]
422
      real  z1     ! Z + 1 DZ [m]
423
      real  z0     ! Z - 1 DZ [m]
424
 
425
      dudz = (u1-u0) / (z1-z0)
426
 
427
      end
428
 
429
c     -------------------------------------------------------------------------
430
c     Wind shear of V (DVDZ)
431
c     -------------------------------------------------------------------------
432
 
433
      subroutine calc_DVDZ (dvdz,v1,v0,z1,z0)
434
 
435
      implicit none
436
 
437
c     Argument declaration
438
      real  dvdz   ! Wind shear [m/s per m = 1/m]
439
      real  v1     ! V @ Z + 1 DZ [m/s]
440
      real  v0     ! V @ Z - 1 DZ [m/s]
441
      real  z1     ! Z + 1 DZ [m]
442
      real  z0     ! Z - 1 DZ [m]
443
 
444
      dvdz = (v1-v0) / (z1-z0)
445
 
446
      end
447
 
448
c     -------------------------------------------------------------------------
449
c     Vertical derivative of T (DTDZ)
450
c     -------------------------------------------------------------------------
451
 
452
      subroutine calc_DTDZ (dtdz,t1,t0,z1,z0)
453
 
454
      implicit none
455
 
456
c     Argument declaration
457
      real  dtdz   ! Vertical derivative of T [K/m]
458
      real  t1     ! T @ Z + 1 DZ [K]
459
      real  t0     ! T @ Z - 1 DZ [K]
460
      real  z1     ! Z + 1 DZ [m]
461
      real  z0     ! Z - 1 DZ [m]
462
 
463
      dtdz = (t1-t0) / (z1-z0)
464
 
465
      end
466
 
467
c     -------------------------------------------------------------------------
468
c     Vertical derivative of TH (DTHDẐ)
469
c     -------------------------------------------------------------------------
470
 
471
      subroutine calc_DTHDZ (dthdz,t1,t0,z1,z0,p,t)
472
 
473
      implicit none
474
 
475
c     Argument declaration
476
      real  dthdz  ! Vertical derivative of TH [K/m]
477
      real  t1     ! T @ Z + 1 DZ [K]
478
      real  t0     ! T @ Z - 1 DZ [K]
479
      real  t      ! T [K]
480
      real  z1     ! Z + 1 DZ [m]
481
      real  z0     ! Z - 1 DZ [m]
482
      real  p      ! P [hPa]
483
 
484
c     Physical parameters
485
      real      rdcp,tzero,pref,g,r
486
      data      rdcp,tzero,pref,g,r /0.286,273.15,1000.,287.,9.81/
487
 
488
c     Auxiliary variables 
489
      real      tk1,tk0,tk
490
 
491
      if (t0.lt.100.) then
492
         tk0 = t0 + tzero
493
      endif
494
      if (t1.lt.100.) then
495
         tk1 = t1 + tzero
496
      endif      
497
      if (t.lt.100.) then
498
         tk  = t  + tzero
499
      endif      
500
 
501
      dthdz = (pref/p)**rdcp *
502
     >        ( (tk1-tk0)/(z1-z0) + rdcp * tk * g / r )
503
 
504
      end
505
 
506
c     -------------------------------------------------------------------------
507
c     Squared Brunt-Vaisäla frequency (NSQ)
508
c     -------------------------------------------------------------------------
509
 
510
      subroutine calc_NSQ (nsq,dthdz,th)
511
 
512
      implicit none
513
 
514
c     Argument declaration
515
      real  nsq    ! Squared Brunt-Vaisäla frequency [s^-1]
516
      real  dthdz  ! D(TH)/DZ [K/m]
517
      real  th     ! K
518
 
519
c     Physical parameters
520
      real      g 
521
      parameter (g=9.81)
522
 
523
      nsq = g/th * dthdz
524
 
525
      end
526
 
527
c     -------------------------------------------------------------------------
528
c     Relative vorticity (RELVORT)
529
c     -------------------------------------------------------------------------
530
 
531
      subroutine calc_RELVORT (relvort,dudy,dvdx,u,lat)
532
 
533
      implicit none
534
 
535
c     Argument declaration
536
      real  relvort    ! Relative vorticity [s^-1]
537
      real  u          ! Zonal wind component [m/s]
538
      real  dudy       ! du/dy [s^-1]
539
      real  dvdx       ! dv/dx [s^-1]
540
      real  lat        ! Latitude [deg]
541
 
542
c     Physical parameters
543
      real      pi180
544
      parameter (pi180=3.14159/180.)
545
      real      deltay
546
      parameter (deltay =1.11e5)
547
 
548
      relvort = dvdx - dudy + u * pi180/deltay * tan(pi180 * lat) 
549
 
550
      end
551
 
552
c     -------------------------------------------------------------------------
553
c     Absolute vorticity (ABSVORT)
554
c     -------------------------------------------------------------------------
555
 
556
      subroutine calc_ABSVORT (absvort,dudy,dvdx,u,lat,lon,
557
     >                                 pollon,pollat)
558
 
559
      implicit none
560
 
561
c     Argument declaration
562
      real  absvort       ! Absolute vorticity [s^-1]
563
      real  u             ! Zonal wind component [m/s]
564
      real  dudy          ! du/dy [s^-1]
565
      real  dvdx          ! dv/dx [s^-1]
566
      real  lat           ! Latitude [deg]
567
      real  lon           ! Longitude [deg]
568
      real  pollon,pollat ! Rotated north pole [deg]
569
 
570
c     Physical parameters
571
      real      pi180
572
      parameter (pi180=3.14159/180.)
573
      real      deltay
574
      parameter (deltay =1.11e5)
575
      real      omega
576
      parameter (omega=7.292e-5)
577
 
578
c     Externals      
579
      real      lmtolms 
580
      real      phtophs    
581
      real      lmstolm
582
      real      phstoph        
583
      external  lmtolms,phtophs
584
      external  lmstolm,phstoph
585
 
586
c     Auxiliary varaibles
587
      real      truelat
588
 
589
      truelat = phstoph(lat,lon,pollat,pollon)
590
 
591
      absvort = dvdx - dudy + u * pi180/deltay * tan(pi180 * lat) +
592
     >          2. * omega * sin(pi180 * truelat) 
593
 
594
      end
595
 
596
c     -------------------------------------------------------------------------
597
c     Divergence (DIV)
598
c     -------------------------------------------------------------------------
599
 
600
      subroutine calc_DIV (div,dudx,dvdy,v,lat)
601
 
602
      implicit none
603
 
604
c     Argument declaration
605
      real  div        ! Divergence [s^-1]
606
      real  v          ! Meridional wind component [m/s]
607
      real  dudx       ! du/dx [s^-1]
608
      real  dvdy       ! dv/dy [s^-1]
609
      real  lat        ! Latitude [deg]
610
 
611
c     Physical parameters
612
      real      pi180
613
      parameter (pi180=3.14159/180.)
614
      real      deltay
615
      parameter (deltay =1.11e5)
616
 
617
      div = dudx + dvdy - v * pi180/deltay * tan(pi180 * lat) 
618
 
619
      end
620
 
621
c     -------------------------------------------------------------------------
622
c     Deformation (DEF)
623
c     -------------------------------------------------------------------------
624
 
625
      subroutine calc_DEF (def,dudx,dvdx,dudy,dvdy)
626
 
627
      implicit none
628
 
629
c     Argument declaration
630
      real  def        ! Deformation [s^-1]
631
      real  dudx       ! du/dx [s^-1]
632
      real  dvdx       ! dv/dy [s^-1]
633
      real  dudy       ! du/dx [s^-1]
634
      real  dvdy       ! dv/dy [s^-1]
635
 
636
c     Physical parameters
637
      real      pi180
638
      parameter (pi180=3.14159/180.)
639
 
640
      def = sqrt( (dvdx+dudy)**2 + (dudx-dvdy)**2 )
641
 
642
      end
643
 
644
c     -------------------------------------------------------------------------
645
c     Potential Vorticity (PV)
646
c     -------------------------------------------------------------------------
647
 
648
      subroutine calc_PV (pv,absvort,dthdz,dudz,dvdz,dthdx,dthdy,rho)
649
 
650
      implicit none
651
 
652
c     Argument declaration
653
      real  pv         ! Ertel-PV [PVU]
654
      real  absvort    ! Absolute vorticity [s^-1]
655
      real  dthdz      ! dth/dz [K/m]
656
      real  dudz       ! du/dz [m/s per m = 1/s]
657
      real  dvdz       ! dv/dz [m/s per m = 1/s]
658
      real  dthdx      ! dth/dx [K/m]
659
      real  dthdy      ! dth/dy [K/m]
660
      real  rho        ! Density [kg/m^3]
661
 
662
c     Physical and numerical parameters
663
      real      scale
664
      parameter (scale=1.E6)
665
      real      g
666
      parameter (g=9.80665)
667
 
668
      pv = scale * 1./rho *
669
     >      ( absvort * dthdz + dudz * dthdy - dvdz * dthdx)
670
 
671
      end
672
 
673
c     -------------------------------------------------------------------------
674
c     Richardson number (RI)
675
c     -------------------------------------------------------------------------
676
 
677
      subroutine calc_RI (ri,dudz,dvdz,nsq)
678
 
679
      implicit none
680
 
681
c     Argument declaration
682
      real  ri         ! Richardson number
683
      real  dudz       ! Du/Dp [m/s per m = 1/s]
684
      real  dvdz       ! Dv/Dp [m/s per m = 1/s]
685
      real  nsq        ! Squared Brunt-Vailälä frequency [s^-1]
686
      real  rho        ! Density [kg/m^3]
687
 
688
c     Physical and numerical parameters
689
      real      g
690
      parameter (g=9.80665)
691
 
692
      ri = nsq / ( dudz**2 + dvdz**2 )
693
 
694
      end
695
 
696
c     -------------------------------------------------------------------------
697
c     Ellrod and Knapp's turbulence indicator (TI)
698
c     -------------------------------------------------------------------------
699
 
700
      subroutine calc_TI (ti,def,dudz,dvdz)
701
 
702
      implicit none
703
 
704
c     Argument declaration
705
      real  ti         ! Turbulence idicator
706
      real  def        ! Deformation [s^-1]
707
      real  dudz       ! Du/Dp [m/s per m = 1/s]
708
      real  dvdz       ! Dv/Dp [m/s per m = 1/s]
709
 
710
      ti = def * sqrt ( dudz**2 + dvdz**2 )
711
 
712
      end
713
 
714
c     -------------------------------------------------------------------------
715
c     Distance from starting position
716
c     -------------------------------------------------------------------------
717
 
718
      subroutine calc_DIST0 (dist0,lon0,lat0,lon1,lat1)
719
 
720
      implicit none
721
 
722
c     Argument declaration
723
      real  dist0      ! Distance from starting position [km]
724
      real  lon0,lat0  ! Starting position
725
      real  lon1,lat1  ! New position 
726
 
727
c     Externals 
728
      real     sdis
729
      external sdis
730
 
731
      dist0 = sdis(lon0,lat0,lon1,lat1)
732
 
733
      end
734
 
735
c     -------------------------------------------------------------------------
736
c     Heading of the trajectory (HEAD)
737
c     -------------------------------------------------------------------------
738
 
739
      subroutine calc_HEAD (head,lon0,lat0,lon1,lat1)
740
 
741
      implicit none
742
 
743
c     Argument declaration
744
      real  head       ! Heading angle (in deg) relativ to zonal direction
745
      real  lon0,lat0  ! Starting position
746
      real  lon1,lat1  ! New position 
747
 
748
c     Physical parameters
749
      real      pi180
750
      parameter (pi180=3.14159/180.)
751
 
752
c     Auixiliary variables
753
      real     dx,dy
754
 
755
      dx = (lon1-lon0) * cos(pi180*0.5*(lat0+lat1))
756
      dy = lat1-lat0
757
 
758
      call getangle(1.,0.,dx,dy,head)
759
 
760
      end
761
 
762
c 
763
c     *************************************************************************
764
c     Auxiliary subroutines and functions
765
c     *************************************************************************
766
 
767
c     -------------------------------------------------------------------------
768
c     Saturation vapor pressure over water
769
c     -------------------------------------------------------------------------
770
 
771
      real function esat(t)
772
 
773
C     This function returns the saturation vapor pressure over water
774
c     (mb) given the temperature (Kelvin).
775
C     The algorithm is due to Nordquist, W. S. ,1973: "Numerical
776
C     Approximations of Selected Meteorological Parameters for Cloud
777
C     Physics Problems" ECOM-5475, Atmospheric Sciences Laboratory,
778
c     U. S. Army Electronics Command, White Sands Missile Range,
779
c     New Mexico 88002.
780
 
781
      real p1,p2,c1,t
782
 
783
      p1=11.344-0.0303998*t
784
      p2=3.49149-1302.8844/t
785
      c1=23.832241-5.02808*log10(t)
786
      esat=10.**(c1-1.3816e-7*10.**p1+8.1328e-3*10.**p2-2949.076/t)
787
 
788
      end
789
 
790
c     --------------------------------------------------------------------------
791
c     Angle between two vectors
792
c     --------------------------------------------------------------------------
793
 
794
      SUBROUTINE getangle (ux1,uy1,ux2,uy2,angle)
795
 
796
c     Given two vectors <ux1,uy1> and <ux2,uy2>, determine the angle (in deg)
797
c     between the two vectors.
798
 
799
      implicit none
800
 
801
c     Declaration of subroutine parameters
802
      real ux1,uy1
803
      real ux2,uy2
804
      real angle
805
 
806
c     Auxiliary variables and parameters
807
      real len1,len2,len3
808
      real val1,val2,val3
809
      real vx1,vy1
810
      real vx2,vy2
811
      real pi
812
      parameter (pi=3.14159265359)
813
 
814
      vx1 = ux1
815
      vx2 = ux2
816
      vy1 = uy1
817
      vy2 = uy2
818
 
819
      len1=sqrt(vx1*vx1+vy1*vy1)
820
      len2=sqrt(vx2*vx2+vy2*vy2)
821
 
822
      if ((len1.gt.0.).and.(len2.gt.0.)) then
823
         vx1=vx1/len1
824
         vy1=vy1/len1
825
         vx2=vx2/len2
826
         vy2=vy2/len2
827
 
828
         val1=vx1*vx2+vy1*vy2
829
         val2=-vy1*vx2+vx1*vy2
830
 
831
         len3=sqrt(val1*val1+val2*val2)
832
 
833
         if ( (val1.ge.0.).and.(val2.ge.0.) ) then
834
            val3=acos(val1/len3)
835
         else if ( (val1.lt.0.).and.(val2.ge.0.) ) then
836
            val3=pi-acos(abs(val1)/len3)
837
         else if ( (val1.ge.0.).and.(val2.le.0.) ) then
838
            val3=-acos(val1/len3)
839
         else if ( (val1.lt.0.).and.(val2.le.0.) ) then
840
            val3=-pi+acos(abs(val1)/len3)
841
         endif
842
      else
843
         val3=0.
844
      endif
845
 
846
      angle=180./pi*val3
847
 
848
      END
849
 
850
 
851
c     --------------------------------------------------------------------------
852
c     Spherical distance between lat/lon points                                                          
853
c     --------------------------------------------------------------------------
854
 
855
      real function sdis(xp,yp,xq,yq)
856
c
857
c     calculates spherical distance (in km) between two points given
858
c     by their spherical coordinates (xp,yp) and (xq,yq), respectively.
859
c
860
      real      re
861
      parameter (re=6370.)
862
      real      pi180
863
      parameter (pi180=3.14159/180.)
864
      real      xp,yp,xq,yq,arg
865
 
866
      arg=sin(pi180*yp)*sin(pi180*yq)+
867
     >    cos(pi180*yp)*cos(pi180*yq)*cos(pi180*(xp-xq))
868
      if (arg.lt.-1.) arg=-1.
869
      if (arg.gt.1.) arg=1.
870
 
871
      sdis=re*acos(arg)
872
 
873
      end
874