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

Last change on this file since 706 was 38, checked in by emillour, 14 years ago

Ajout du modè Martien (mon LMDZ.MARS.BETA, du 28/01/2011) dans le rértoire mars, pour pouvoir suivre plus facilement les modifs.
EM

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