source: trunk/LMDZ.MARS/libf/phymars/iniorbit.F @ 1242

Last change on this file since 1242 was 1226, checked in by aslmd, 11 years ago

LMDZ.MARS : Replaced comcstfi and planete includes by modules.

File size: 2.6 KB
Line 
1      SUBROUTINE iniorbit
2     $     (paphelie,pperiheli,pyear_day,pperi_day,pobliq)
3      use planete_h
4      USE comcstfi_h
5      IMPLICIT NONE
6
7c=======================================================================
8c
9c   Auteur:
10c   -------
11c     Frederic Hourdin      22 Fevrier 1991
12c
13c   Objet:
14c   ------
15c    Initialisation du sous programme orbite qui calcule
16c    a une date donnee de l'annee de duree year_day commencant
17c    a l'equinoxe de printemps et dont le perihelie se situe
18c    a la date peri_day, la distance au soleil et la declinaison.
19c
20c   Interface:
21c   ----------
22c   - Doit etre appele avant d'utiliser orbite.
23c   - initialise une partie du common planete.h
24c
25c   Arguments:
26c   ----------
27c
28c   Input:
29c   ------
30c   aphelie       \   aphelie et perihelie de l'orbite
31c   periheli      /   en millions de kilometres.
32c
33c=======================================================================
34
35c-----------------------------------------------------------------------
36c   Declarations:
37c   -------------
38
39c   Arguments:
40c   ----------
41
42      REAL paphelie,pperiheli,pyear_day,pperi_day,pobliq
43
44c   Local:
45c   ------
46
47      REAL zxref,zanom,zz,zx0,zdx
48      INTEGER iter
49
50c-----------------------------------------------------------------------
51
52      pi=2.*asin(1.)
53
54      aphelie =paphelie
55      periheli=pperiheli
56      year_day=pyear_day
57      obliquit=pobliq
58      peri_day=pperi_day
59
60      PRINT*,'Perihelie en Mkm  ',periheli
61      PRINT*,'Aphelise  en Mkm  ',aphelie
62      PRINT*,'obliquite en degres  :',obliquit
63      unitastr=149.597927
64      e_elips=(aphelie-periheli)/(periheli+aphelie)
65      p_elips=0.5*(periheli+aphelie)*(1-e_elips*e_elips)/unitastr
66
67      print*,'e_elips',e_elips
68      print*,'p_elips',p_elips
69
70c-----------------------------------------------------------------------
71c calcul de l'angle polaire et de la distance au soleil :
72c -------------------------------------------------------
73
74c  calcul de l'zanomalie moyenne
75
76      zz=(year_day-pperi_day)/year_day
77      zanom=2.*pi*(zz-nint(zz))
78      zxref=abs(zanom)
79      PRINT*,'zanom  ',zanom
80
81c  resolution de l'equation horaire  zx0 - e * sin (zx0) = zxref
82c  methode de Newton
83
84      zx0=zxref+e_elips*sin(zxref)
85      DO 110 iter=1,100
86         zdx=-(zx0-e_elips*sin(zx0)-zxref)/(1.-e_elips*cos(zx0))
87         if(abs(zdx).le.(1.e-12)) goto 120
88         zx0=zx0+zdx
89110   continue
90120   continue
91      zx0=zx0+zdx
92      if(zanom.lt.0.) zx0=-zx0
93      PRINT*,'zx0   ',zx0
94
95c zteta est la longitude solaire
96
97      timeperi=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
98      PRINT*,'longitude solaire du perihelie timeperi = ',timeperi
99
100      RETURN
101      END
Note: See TracBrowser for help on using the repository browser.