/[lmdze]/trunk/dyn3d/dist_sphe.f
ViewVC logotype

Contents of /trunk/dyn3d/dist_sphe.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 212 - (show annotations)
Thu Jan 12 12:31:31 2017 UTC (7 years, 4 months ago) by guez
Original Path: trunk/Sources/dyn3d/dist_sphe.f
File size: 1466 byte(s)
Moved variables from module com_io_dyn to module inithist_m, where
they are defined.

Split grid_atob.f into grille_m.f and dist_sphe.f. Extracted ASCCI art
to documentation. In grille_m, use automatic arrays instead of maximum
size. In grille_m, instead of printing data for every problematic
point, print a single diagnostic message.

Removed variables top_height, overlap, lev_histhf, lev_histday,
lev_histmth, type_run, ok_isccp, ok_regdyn, lonmin_ins, lonmax_ins,
latmin_ins, latmax_ins of module clesphys, not used.

Removed variable itap of module histwrite_phy_m, not used. There is a
variable itap in module time_phylmdz.

Added output of tro3.

In physiq, no need to compute wo at every time-step, since we only use
it in radlwsw.

1 module dist_sphe_m
2
3 ! From grid_atob.F, v 1.1.1.1 2004/05/19 12:53:05
4
5 IMPLICIT none
6
7 contains
8
9 SUBROUTINE dist_sphe(rf_lon, rf_lat, rlon, rlat, im, jm, distance)
10
11 ! Auteur: Laurent Li (le 30 decembre 1996)
12
13 ! Ce programme calcule la distance minimale (selon le grand cercle)
14 ! entre deux points sur la terre
15
16 INTEGER, intent(in):: im, jm ! dimensions
17 REAL, intent(in):: rf_lon ! longitude du point de reference (degres)
18 REAL, intent(in):: rf_lat ! latitude du point de reference (degres)
19 REAL, intent(in):: rlon(im), rlat(jm) ! longitude et latitude des points
20
21 REAL, intent(out):: distance(im, jm) ! distances en metre
22
23 REAL rlon1, rlat1
24 REAL rlon2, rlat2
25 REAL dist
26 REAL pa, pb, p, pi
27
28 REAL radius
29 PARAMETER (radius=6371229.)
30 integer i, j
31
32 !---------------------------------------------------------------------
33
34 pi = 4.0 * ATAN(1.0)
35
36 DO j = 1, jm
37 DO i = 1, im
38
39 rlon1=rf_lon
40 rlat1=rf_lat
41 rlon2=rlon(i)
42 rlat2=rlat(j)
43 pa = pi/2.0 - rlat1*pi/180.0 ! dist. entre pole n et point a
44 pb = pi/2.0 - rlat2*pi/180.0 ! dist. entre pole n et point b
45 p = (rlon1-rlon2)*pi/180.0 ! angle entre a et b (leurs meridiens)
46
47 dist = ACOS(COS(pa)*COS(pb) + SIN(pa)*SIN(pb)*COS(p))
48 dist = radius * dist
49 distance(i, j) = dist
50
51 end DO
52 end DO
53
54 END SUBROUTINE dist_sphe
55
56 end module dist_sphe_m

  ViewVC Help
Powered by ViewVC 1.1.21