source: LMDZ4/trunk/libf/phylmd/solarlong.F @ 2078

Last change on this file since 2078 was 1387, checked in by Laurent Fairhead, 15 years ago

Correction of mixed MPI/OpenMP bug


Correction bug mixte MPI/OpenMP

File size: 3.5 KB
RevLine 
[1201]1      SUBROUTINE solarlong(pday,psollong,pdist_sol)
2
3      USE ioipsl
4
5      IMPLICIT NONE
6
7c=======================================================================
8c
9c   Objet:
10c   ------
11c
12c      Calcul de la distance soleil-planete et de la declinaison
13c   en fonction du jour de l'annee.
14c
15c
16c   Methode:
17c   --------
18c
19c      Calcul complet de l'elipse
20c
21c   Interface:
22c   ----------
23c
24c      Uncommon comprenant les parametres orbitaux.
25c
26c   Arguments:
27c   ----------
28c
29c   Input:
30c   ------
31c   pday          jour de l'annee (le jour 0 correspondant a l'equinoxe)
32c   lwrite        clef logique pour sorties de controle
33c
34c   Output:
35c   -------
36c   pdist_sol     distance entre le soleil et la planete
37c                 ( en unite astronomique pour utiliser la constante
38c                  solaire terrestre 1370 Wm-2 )
39c   pdecli        declinaison ( en radians )
40c
41c=======================================================================
42c-----------------------------------------------------------------------
43c   Declarations:
44c   -------------
45
46#include "planete.h"
47#include "YOMCST.h"
48      include 'iniprint.h'
49
50c arguments:
51c ----------
52
53      REAL pday,pdist_sol,pdecli,psollong
54      LOGICAL lwrite
55
56c Local:
57c ------
58
59      REAL zanom,xref,zx0,zdx,zteta,zz,pi
60      INTEGER iter
61      REAL :: pyear_day,pperi_day
62      REAL :: jD_eq, jD_peri
[1387]63      LOGICAL,SAVE :: first=.TRUE.
64c$OMP THREADPRIVATE(first)
[1201]65
66c-----------------------------------------------------------------------
67c calcul de l'angle polaire et de la distance au soleil :
68c -------------------------------------------------------
69
70c   Initialisation eventuelle:
[1387]71      if(first) then
[1201]72        call ioget_calendar(pyear_day)
73        call ymds2ju(2000, 3, 21, 0., jD_eq)
74        call ymds2ju(2001, 1, 4, 0., jD_peri)
75        pperi_day = jD_peri - jD_eq
76        pperi_day = R_peri + 180.
77        write(lunout,*)' Number of days in a year = ',pyear_day
78c         call iniorbit(249.22,206.66,669.,485.,25.2)
79         call iniorbit(152.59,146.61,pyear_day,pperi_day,R_incl)
[1387]80         first=.FALSE.
[1201]81      endif
82
83c  calcul de l'zanomalie moyenne
84
85      zz=(pday-peri_day)/year_day
86      pi=2.*asin(1.)
87      zanom=2.*pi*(zz-nint(zz))
88      xref=abs(zanom)
89
90c  resolution de l'equation horaire  zx0 - e * sin (zx0) = xref
91c  methode de Newton
92
93!      zx0=xref+e_elips*sin(xref)
94      zx0=xref+R_ecc*sin(xref)
95      DO 110 iter=1,10
96!         zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0))
97         zdx=-(zx0-R_ecc*sin(zx0)-xref)/(1.-R_ecc*cos(zx0))
98         if(abs(zdx).le.(1.e-7)) goto 120
99         zx0=zx0+zdx
100110   continue
101120   continue
102      zx0=zx0+zdx
103      if(zanom.lt.0.) zx0=-zx0
104
105c zteta est la longitude solaire
106
107!      zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
108      zteta=2.*atan(sqrt((1.+R_ecc)/(1.-R_ecc))*tan(zx0/2.))
109
110      psollong=zteta-timeperi
111
112      IF(psollong.LT.0.) psollong=psollong+2.*pi
113      IF(psollong.GT.2.*pi) psollong=psollong-2.*pi
114
115      psollong = psollong * 180. / pi
116
117c distance soleil
118
119      pdist_sol = (1-R_ecc*R_ecc)
120     &      /(1+R_ecc*COS(pi/180.*(psollong-(R_peri+180.0))))
121!      pdist_sol = (1-e_elips*e_elips)
122!     &      /(1+e_elips*COS(pi/180.*(psollong-(R_peri+180.0))))
123c-----------------------------------------------------------------------
124c   sorties eventuelles:
125c   ---------------------
126
127c     IF (lwrite) THEN
128c        PRINT*,'jour de l"annee   :',pday
129c        PRINT*,'distance au soleil (en unite astronomique) :',pdist_sol
130c        PRINT*,'declinaison (en degres) :',pdecli*180./pi
131c     ENDIF
132
133      RETURN
134      END
Note: See TracBrowser for help on using the repository browser.