source: LMDZ6/branches/Amaury_posttrusting/libf/phylmd/orbite.f90 @ 5333

Last change on this file since 5333 was 5285, checked in by abarral, 3 weeks ago

As discussed internally, remove generic ONLY: ... for new _mod_h modules

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 10.2 KB
Line 
1
2! $Header$
3
4! ======================================================================
5SUBROUTINE orbite(xjour, longi, dist)
6  USE yomcst_mod_h
7IMPLICIT NONE
8  ! ======================================================================
9  ! Auteur(s): Z.X. Li (LMD/CNRS) (adapte du GCM du LMD) date: 19930818
10  ! Objet: pour un jour donne, calculer la longitude vraie de la terre
11  ! (par rapport au point vernal-21 mars) dans son orbite solaire
12  ! calculer aussi la distance terre-soleil (unite astronomique)
13  ! ======================================================================
14  ! Arguments:
15  ! xjour--INPUT--R- jour de l'annee a compter du 1er janvier
16  ! longi--OUTPUT-R- longitude vraie en degres par rapport au point
17  ! vernal (21 mars) en degres
18  ! dist---OUTPUT-R- distance terre-soleil (par rapport a la moyenne)
19  REAL xjour, longi, dist
20  ! ======================================================================
21
22
23  ! -- Variables dynamiques locales
24  REAL pir, xl, xllp, xee, xse, xlam, dlamm, anm, ranm, anv, ranv
25
26  pir = 4.0*atan(1.0)/180.0
27  xl = r_peri + 180.0
28  xllp = xl*pir
29  xee = r_ecc*r_ecc
30  xse = sqrt(1.0-xee)
31  xlam = (r_ecc/2.0+r_ecc*xee/8.0)*(1.0+xse)*sin(xllp) - &
32    xee/4.0*(0.5+xse)*sin(2.0*xllp) + r_ecc*xee/8.0*(1.0/3.0+xse)*sin(3.0* &
33    xllp)
34  xlam = 2.0*xlam/pir
35  dlamm = xlam + (xjour-81.0)
36  anm = dlamm - xl
37  ranm = anm*pir
38  xee = xee*r_ecc
39  ranv = ranm + (2.0*r_ecc-xee/4.0)*sin(ranm) + 5.0/4.0*r_ecc*r_ecc*sin(2.0* &
40    ranm) + 13.0/12.0*xee*sin(3.0*ranm)
41
42  anv = ranv/pir
43  longi = anv + xl
44
45  dist = (1-r_ecc*r_ecc)/(1+r_ecc*cos(pir*(longi-(r_peri+180.0))))
46  RETURN
47END SUBROUTINE orbite
48! ======================================================================
49SUBROUTINE angle(longi, lati, frac, muzero)
50  USE dimphy
51  USE yomcst_mod_h
52IMPLICIT NONE
53  ! ======================================================================
54  ! Auteur(s): Z.X. Li (LMD/CNRS) date: 19930818
55  ! Objet: Calculer la duree d'ensoleillement pour un jour et la hauteur
56  ! du soleil (cosinus de l'angle zinithal) moyenne sur la journee
57  ! ======================================================================
58  ! Arguments:
59  ! longi----INPUT-R- la longitude vraie de la terre dans son plan
60  ! solaire a partir de l'equinoxe de printemps (degre)
61  ! lati-----INPUT-R- la latitude d'un point sur la terre (degre)
62  ! frac-----OUTPUT-R la duree d'ensoleillement dans la journee divisee
63  ! par 24 heures (unite en fraction de 0 a 1)
64  ! muzero---OUTPUT-R la moyenne du cosinus de l'angle zinithal sur
65  ! la journee (0 a 1)
66  ! ======================================================================
67  REAL longi
68  REAL lati(klon), frac(klon), muzero(klon)
69
70  REAL lat, omega, lon_sun, lat_sun
71  REAL pi_local, incl
72  INTEGER i
73
74  pi_local = 4.0*atan(1.0)
75  incl = r_incl*pi_local/180.
76
77  lon_sun = longi*pi_local/180.0
78  lat_sun = asin(sin(lon_sun)*sin(incl))
79
80  DO i = 1, klon
81    lat = lati(i)*pi_local/180.0
82
83    IF (lat>=(pi_local/2.+lat_sun) .OR. lat<=(-pi_local/2.+lat_sun)) THEN
84      omega = 0.0 ! nuit polaire
85    ELSE IF (lat>=(pi_local/2.-lat_sun) .OR. lat<=(-pi_local/2.-lat_sun)) &
86        THEN
87      omega = pi_local ! journee polaire
88    ELSE
89      omega = -tan(lat)*tan(lat_sun)
90      omega = acos(omega)
91    END IF
92
93    frac(i) = omega/pi_local
94
95    IF (omega>0.0) THEN
96      muzero(i) = sin(lat)*sin(lat_sun) + cos(lat)*cos(lat_sun)*sin(omega)/ &
97        omega
98    ELSE
99      muzero(i) = 0.0
100    END IF
101  END DO
102
103  RETURN
104END SUBROUTINE angle
105! ====================================================================
106SUBROUTINE zenang(longi, gmtime, pdtrad1, pdtrad2, lat, long, pmu0, frac)
107  USE dimphy
108  USE yomcst_mod_h
109IMPLICIT NONE
110  ! =============================================================
111  ! Auteur : O. Boucher (LMD/CNRS)
112  ! d'apres les routines zenith et angle de Z.X. Li
113  ! Objet  : calculer les valeurs moyennes du cos de l'angle zenithal
114  ! et l'ensoleillement moyen entre gmtime1 et gmtime2
115  ! connaissant la declinaison, la latitude et la longitude.
116  ! Rque   : Different de la routine angle en ce sens que zenang
117  ! fournit des moyennes de pmu0 et non des valeurs
118  ! instantanees, du coup frac prend toutes les valeurs
119  ! entre 0 et 1. La routine integre entre gmtime+pdtrad1 et
120  ! gmtime+pdtrad2 avec pdtrad1 et pdtrad2 exprimes en secondes.
121  ! Date   : premiere version le 13 decembre 1994
122  ! revu pour  GCM  le 30 septembre 1996
123  ! revu le 3 septembre 2015 pour les bornes de l'integrale
124  ! ===============================================================
125  ! longi : la longitude vraie de la terre dans son plan
126  ! solaire a partir de l'equinoxe de printemps (degre)
127  ! gmtime : temps universel en fraction de jour
128  ! pdtrad1 : borne inferieure du pas de temps du rayonnement (secondes)
129  ! pdtrad2 : borne inferieure du pas de temps du rayonnement (secondes)
130  ! pdtrad2-pdtrad1 correspond a pdtrad, le pas de temps du rayonnement (secondes)
131  ! lat------INPUT : latitude en degres
132  ! long-----INPUT : longitude en degres
133  ! pmu0-----OUTPUT: angle zenithal moyen entre gmtime+pdtrad1 et gmtime+pdtrad2
134  ! frac-----OUTPUT: ensoleillement moyen entre gmtime+pdtrad1 et gmtime+pdtrad2
135  ! ================================================================
136
137  ! ================================================================
138  REAL, INTENT (IN) :: longi, gmtime, pdtrad1, pdtrad2
139  REAL lat(klon), long(klon), pmu0(klon), frac(klon)
140  ! ================================================================
141  INTEGER i
142  REAL gmtime1, gmtime2
143  REAL pi_local, deux_pi_local, incl
144  REAL omega1, omega2, omega
145  ! omega1, omega2 : temps 1 et 2 exprime en radian avec 0 a midi.
146  ! omega : heure en radian du coucher de soleil
147  ! -omega est donc l'heure en radian de lever du soleil
148  REAL omegadeb, omegafin
149  REAL zfrac1, zfrac2, z1_mu, z2_mu
150  REAL lat_sun ! declinaison en radian
151  REAL lon_sun ! longitude solaire en radian
152  REAL latr    ! latitude du pt de grille en radian
153  ! ================================================================
154
155  pi_local = 4.0*atan(1.0)
156  deux_pi_local = 2.0*pi_local
157  incl = r_incl*pi_local/180.
158
159  lon_sun = longi*pi_local/180.0
160  lat_sun = asin(sin(lon_sun)*sin(incl))
161
162  gmtime1 = gmtime*86400. + pdtrad1
163  gmtime2 = gmtime*86400. + pdtrad2
164
165  DO i = 1, klon
166
167    latr = lat(i)*pi_local/180.
168
169    omega = 0.0 !--nuit polaire
170
171    IF (latr>=(pi_local/2.-lat_sun) .OR. latr<=(-pi_local/2.-lat_sun)) THEN
172      omega = pi_local ! journee polaire
173    END IF
174
175    IF (latr<(pi_local/2.+lat_sun) .AND. latr>(-pi_local/2.+lat_sun) .AND. &
176        latr<(pi_local/2.-lat_sun) .AND. latr>(-pi_local/2.-lat_sun)) THEN
177      omega = -tan(latr)*tan(lat_sun)
178      omega = acos(omega)
179    END IF
180
181    omega1 = gmtime1 + long(i)*86400.0/360.0
182    omega1 = omega1/86400.0*deux_pi_local
183    omega1 = mod(omega1+deux_pi_local, deux_pi_local)
184    omega1 = omega1 - pi_local
185
186    omega2 = gmtime2 + long(i)*86400.0/360.0
187    omega2 = omega2/86400.0*deux_pi_local
188    omega2 = mod(omega2+deux_pi_local, deux_pi_local)
189    omega2 = omega2 - pi_local
190
191    IF (omega1<=omega2) THEN !--on est dans la meme journee locale
192
193      IF (omega2<=-omega .OR. omega1>=omega .OR. omega<1E-5) THEN !--nuit
194        frac(i) = 0.0
195        pmu0(i) = 0.0
196      ELSE !--jour+nuit/jour
197        omegadeb = max(-omega, omega1)
198        omegafin = min(omega, omega2)
199        frac(i) = (omegafin-omegadeb)/(omega2-omega1)
200        pmu0(i) = sin(latr)*sin(lat_sun) + cos(latr)*cos(lat_sun)*(sin( &
201          omegafin)-sin(omegadeb))/(omegafin-omegadeb)
202      END IF
203
204    ELSE !---omega1 GT omega2 -- a cheval sur deux journees
205
206      ! -------------------entre omega1 et pi
207      IF (omega1>=omega) THEN !--nuit
208        zfrac1 = 0.0
209        z1_mu = 0.0
210      ELSE !--jour+nuit
211        omegadeb = max(-omega, omega1)
212        omegafin = omega
213        zfrac1 = omegafin - omegadeb
214        z1_mu = sin(latr)*sin(lat_sun) + cos(latr)*cos(lat_sun)*(sin(omegafin &
215          )-sin(omegadeb))/(omegafin-omegadeb)
216      END IF
217      ! ---------------------entre -pi et omega2
218      IF (omega2<=-omega) THEN !--nuit
219        zfrac2 = 0.0
220        z2_mu = 0.0
221      ELSE !--jour+nuit
222        omegadeb = -omega
223        omegafin = min(omega, omega2)
224        zfrac2 = omegafin - omegadeb
225        z2_mu = sin(latr)*sin(lat_sun) + cos(latr)*cos(lat_sun)*(sin(omegafin &
226          )-sin(omegadeb))/(omegafin-omegadeb)
227
228      END IF
229      ! -----------------------moyenne
230      frac(i) = (zfrac1+zfrac2)/(omega2+deux_pi_local-omega1)
231      pmu0(i) = (zfrac1*z1_mu+zfrac2*z2_mu)/max(zfrac1+zfrac2, 1.E-10)
232
233    END IF !---comparaison omega1 et omega2
234
235  END DO
236
237END SUBROUTINE zenang
238! ===================================================================
239SUBROUTINE zenith(longi, gmtime, lat, long, pmu0, fract)
240  USE dimphy
241  USE yomcst_mod_h
242IMPLICIT NONE
243
244  ! Auteur(s): Z.X. Li (LMD/ENS)
245
246  ! Objet: calculer le cosinus de l'angle zenithal du soleil en
247  ! connaissant la declinaison du soleil, la latitude et la
248  ! longitude du point sur la terre, et le temps universel
249
250  ! Arguments d'entree:
251  ! longi  : declinaison du soleil (en degres)
252  ! gmtime : temps universel en second qui varie entre 0 et 86400
253  ! lat    : latitude en degres
254  ! long   : longitude en degres
255  ! Arguments de sortie:
256  ! pmu0   : cosinus de l'angle zenithal
257
258  ! ====================================================================
259
260  ! ====================================================================
261  REAL longi, gmtime
262  REAL lat(klon), long(klon), pmu0(klon), fract(klon)
263  ! =====================================================================
264  INTEGER n
265  REAL zpi, zpir, omega, zgmtime
266  REAL incl, lat_sun, lon_sun
267  ! ----------------------------------------------------------------------
268  zpi = 4.0*atan(1.0)
269  zpir = zpi/180.0
270  zgmtime = gmtime*86400.
271
272  incl = r_incl*zpir
273
274  lon_sun = longi*zpir
275  lat_sun = asin(sin(lon_sun)*sin(incl))
276
277  ! --initialisation a la nuit
278
279  DO n = 1, klon
280    pmu0(n) = 0.
281    fract(n) = 0.0
282  END DO
283
284  ! 1 degre en longitude = 240 secondes en temps
285
286  DO n = 1, klon
287    omega = zgmtime + long(n)*86400.0/360.0
288    omega = omega/86400.0*2.0*zpi
289    omega = mod(omega+2.0*zpi, 2.0*zpi)
290    omega = omega - zpi
291    pmu0(n) = sin(lat(n)*zpir)*sin(lat_sun) + cos(lat(n)*zpir)*cos(lat_sun)* &
292      cos(omega)
293    pmu0(n) = max(pmu0(n), 0.0)
294    IF (pmu0(n)>1.E-6) fract(n) = 1.0
295  END DO
296
297  RETURN
298END SUBROUTINE zenith
Note: See TracBrowser for help on using the repository browser.