source: trunk/LMDZ.MARS/libf/phymars/tabfi.F @ 672

Last change on this file since 672 was 672, checked in by tnavarro, 13 years ago

change perihelion date in solar longitude for any orbit, not only present day one

File size: 22.5 KB
Line 
1c=======================================================================
2      SUBROUTINE tabfi(nid,Lmodif,tab0,day_ini,lmax,p_rad,
3     .                 p_omeg,p_g,p_mugaz,p_daysec,time)
4c=======================================================================
5c
6c   C. Hourdin 15/11/96
7c
8c   Object:        Lecture du tab_cntrl physique dans un fichier
9c   ------            et initialisation des constantes physiques
10c
11c   Arguments:
12c   ----------
13c
14c     Inputs:
15c     ------
16c
17c      - nid:    unitne logique du fichier ou on va lire le tab_cntrl   
18c                      (ouvert dans le programme appellant)
19c
20c                 si nid=0:
21c                       pas de lecture du tab_cntrl mais
22c                       Valeurs par default des constantes physiques
23c       
24c      - tab0:    Offset de tab_cntrl a partir duquel sont ranges
25c                  les parametres physiques (50 pour start_archive)
26c
27c      - Lmodif:  si on souhaite modifier les constantes  Lmodif = 1 = TRUE
28c
29c
30c     Outputs:
31c     --------
32c
33c      - day_ini: tab_cntrl(tab0+3) (Dans les cas ou l'on souhaite
34c                              comparer avec le day_ini dynamique)
35c
36c      - lmax:    tab_cntrl(tab0+2) (pour test avec nlayermx)
37c
38c      - p_rad
39c      - p_omeg   !
40c      - p_g      ! Constantes physiques ayant des
41c      - p_mugaz  ! homonymes dynamiques
42c      - p_daysec !
43c
44c=======================================================================
45 
46      implicit none
47 
48#include "dimensions.h"
49#include "dimphys.h"
50#include "comcstfi.h"
51#include "comgeomfi.h"
52#include "planete.h"
53#include "surfdat.h"
54#include "comsoil.h"
55#include "netcdf.inc"
56#include "dimradmars.h"
57#include "yomaer.h"
58
59c-----------------------------------------------------------------------
60c   Declarations
61c-----------------------------------------------------------------------
62
63c Arguments
64c ---------
65      INTEGER nid,nvarid,tab0
66      INTEGER*4 day_ini
67      INTEGER Lmodif
68      INTEGER lmax
69      REAL p_rad,p_omeg,p_g,p_mugaz,p_daysec,time,peri_ls
70
71c Variables
72c ---------
73      INTEGER length
74      parameter (length = 100)
75      REAL tab_cntrl(length) ! array in which are stored the run's parameters
76      INTEGER  ierr
77      INTEGER size
78      CHARACTER modif*20
79
80c Fonctions DRS et autres
81c -----------------------
82      INTEGER setname, cluvdb, getdat
83
84c-----------------------------------------------------------------------
85c  Initialization of various physical constants to defaut values (nid = 0 case)
86c-----------------------------------------------------------------------
87      IF (nid.eq.0) then
88 
89c Reference pressure
90c-------------------------------------
91c     pressrf = 670.            ! Pression de reference (Pa) ~650
92
93c Infos about Mars for the dynamics and physics
94c----------------------------------------------------------
95      rad=3397200.          ! radius of Mars (m)  ~3397200 m
96      daysec=88775.         ! length of a sol (s)  ~88775 s
97      omeg=4.*asin(1.)/(daysec)       ! rotation rate  (rad.s-1)
98      g=3.72                ! gravity (m.s-2) ~3.72
99      mugaz=43.49           ! Molar mass of the atmosphere (g.mol-1) ~43.49
100      rcp=.256793         ! = r/cp  ~0.256793
101
102c Informations about Mars, only for physics
103c-----------------------------------------------------
104      year_day = 669.       !Modif FH: length of year (sols) ~668.6
105      periheli = 206.66         ! min. Sun-Mars distance (Mkm) ~206.66
106      aphelie = 249.22          ! max. Sun-Mars distance (Mkm) ~249.22
107      peri_day =  485.    ! date of perihelion (sols since N. spring)
108      obliquit = 25.19         ! Obliquity of the planet (deg) ~25.19
109
110c Boundary layer and turbulence
111c----------------------------
112      z0_default =  1.e-2       ! surface roughness (m) ~0.01
113      emin_turb = 1.e-6         ! minimal energy ~1.e-8
114      lmixmin = 30              ! mixing length ~100
115
116c Optical properties of polar caps and ground emissivity
117c-----------------------------------------------------
118      emissiv=.95               ! Emissivity of martian soil ~.95
119      emisice(1)=0.95           ! Emissivity of northern cap
120      emisice(2)=0.95           ! Emissivity of southern cap
121      albedice(1)=0.65          ! Albedo of northern cap
122      albedice(2)=0.65          ! Albedo of southern cap
123      iceradius(1) = 100.e-6    ! mean scat radius of CO2 snow (north)
124      iceradius(2) = 100.e-6    ! mean scat radius of CO2 snow (south)
125      dtemisice(1) = 0.4   ! time scale for snow metamorphism (north)
126      dtemisice(2) = 0.4   ! time scale for snow metamorphism (south)
127
128c dust aerosol properties
129c---------------------------------
130      tauvis= 0.2          ! mean visible optical depth
131
132c  Ancien code radiatif (non utilise avec le code d'apres 03/96)
133c---------------------------------------------------------------
134c     tauir= 0.  ! .2  ratio (mean IR opt.depth)/Visible
135c     scatalb=0. ! .86 scaterring albedo visible (~.86)
136c     asfact=0.  ! .79 assymetrie factor visible   (~.79)
137c     day0 = 0   ! = 0 en general !!!
138
139c soil properties
140      volcapa = 1.e6 ! soil volumetric heat capacity (in comsoil.h)
141      ELSE
142c-----------------------------------------------------------------------
143c  Initialization of physical constants by reading array tab_cntrl(:)
144c               which contains these parameters (nid != 0 case)
145c-----------------------------------------------------------------------
146c Read 'controle' array
147c
148      ierr = NF_INQ_VARID (nid, "controle", nvarid)
149      IF (ierr .NE. NF_NOERR) THEN
150         PRINT*, "Tabfi: Could not find <controle> data"
151         CALL abort
152      ENDIF
153#ifdef NC_DOUBLE
154      ierr = NF_GET_VAR_DOUBLE(nid, nvarid, tab_cntrl)
155#else
156      ierr = NF_GET_VAR_REAL(nid, nvarid, tab_cntrl)
157#endif
158      IF (ierr .NE. NF_NOERR) THEN
159         PRINT*, "Tabfi: Failed reading <controle> array"
160         CALL abort
161      ENDIF
162
163       print*,'tabfi: tab_cntrl',tab_cntrl
164c
165c  Initialization of some physical constants
166c informations on physics grid
167      if(ngridmx.ne.tab_cntrl(tab0+1)) then
168         print*,'tabfi: WARNING !!! tab_cntrl(tab0+1).ne.ngridmx'
169         print*,tab_cntrl(tab0+1),ngridmx
170      endif
171      lmax = nint(tab_cntrl(tab0+2))
172      day_ini = tab_cntrl(tab0+3)
173      time = tab_cntrl(tab0+4)
174      write (*,*) 'IN tabfi day_ini=',day_ini
175c Informations about planet Mars for dynamics and physics
176      rad = tab_cntrl(tab0+5)
177      omeg = tab_cntrl(tab0+6)
178      g = tab_cntrl(tab0+7)
179      mugaz = tab_cntrl(tab0+8)
180      rcp = tab_cntrl(tab0+9)
181      daysec = tab_cntrl(tab0+10)
182      dtphys = tab_cntrl(tab0+11)
183c Informations about planet Mars for the physics only
184      year_day = tab_cntrl(tab0+14)
185      periheli = tab_cntrl(tab0+15)
186      aphelie = tab_cntrl(tab0+16)
187      peri_day = tab_cntrl(tab0+17)
188      obliquit = tab_cntrl(tab0+18)
189c boundary layer and turbeulence
190      z0_default = tab_cntrl(tab0+19)
191      lmixmin = tab_cntrl(tab0+20)
192      emin_turb = tab_cntrl(tab0+21)
193c optical properties of polar caps and ground emissivity
194      albedice(1)= tab_cntrl(tab0+22)
195      albedice(2)= tab_cntrl(tab0+23)
196      emisice(1) = tab_cntrl(tab0+24)
197      emisice(2) = tab_cntrl(tab0+25)
198      emissiv    = tab_cntrl(tab0+26)
199      tauvis     = tab_cntrl(tab0+27)  ! dust opt. depth vis.
200      iceradius(1)= tab_cntrl(tab0+31) ! mean scat radius of CO2 snow (north)
201      iceradius(2)= tab_cntrl(tab0+32) ! mean scat radius of CO2 snow (south)
202      dtemisice(1)= tab_cntrl(tab0+33) !time scale for snow metamorphism (north)
203      dtemisice(2)= tab_cntrl(tab0+34) !time scale for snow metamorphism (south)
204c soil properties
205      volcapa = tab_cntrl(tab0+35) ! volumetric heat capacity
206c-----------------------------------------------------------------------
207c       Save some constants for later use (as routine arguments)
208c-----------------------------------------------------------------------
209      p_omeg = omeg
210      p_g = g
211      p_mugaz = mugaz
212      p_daysec = daysec
213      p_rad=rad
214
215      ENDIF    ! end of (nid = 0)
216
217c-----------------------------------------------------------------------
218c       Write physical constants to output before modifying them
219c-----------------------------------------------------------------------
220 
221   6  FORMAT(a20,e15.6,e15.6)
222   5  FORMAT(a20,f12.2,f12.2)
223 
224      write(*,*) '*****************************************************'
225      write(*,*) 'Reading tab_cntrl when calling tabfi before changes'
226      write(*,*) '*****************************************************'
227      write(*,5) '(1)      = ngridmx?',tab_cntrl(tab0+1),real(ngridmx)
228      write(*,5) '(2)            lmax',tab_cntrl(tab0+2),real(lmax)
229      write(*,5) '(3)         day_ini',tab_cntrl(tab0+3),real(day_ini)
230      write(*,5) '(5)             rad',tab_cntrl(tab0+5),rad
231      write(*,5) '(10)         daysec',tab_cntrl(tab0+10),daysec
232      write(*,6) '(6)            omeg',tab_cntrl(tab0+6),omeg
233      write(*,5) '(7)               g',tab_cntrl(tab0+7),g
234      write(*,5) '(8)           mugaz',tab_cntrl(tab0+8),mugaz
235      write(*,5) '(9)             rcp',tab_cntrl(tab0+9),rcp
236      write(*,6) '(11)        dtphys?',tab_cntrl(tab0+11),dtphys
237
238      write(*,5) '(14)       year_day',tab_cntrl(tab0+14),year_day
239      write(*,5) '(15)       periheli',tab_cntrl(tab0+15),periheli
240      write(*,5) '(16)        aphelie',tab_cntrl(tab0+16),aphelie
241      write(*,5) '(17)       peri_day',tab_cntrl(tab0+17),peri_day
242      write(*,5) '(18)       obliquit',tab_cntrl(tab0+18),obliquit
243
244      write(*,6) '(19)     z0_default',tab_cntrl(tab0+19),z0_default
245      write(*,6) '(21)      emin_turb',tab_cntrl(tab0+21),emin_turb
246      write(*,5) '(20)        lmixmin',tab_cntrl(tab0+20),lmixmin
247
248      write(*,5) '(26)        emissiv',tab_cntrl(tab0+26),emissiv
249      write(*,5) '(24)     emisice(1)',tab_cntrl(tab0+24),emisice(1)
250      write(*,5) '(25)     emisice(2)',tab_cntrl(tab0+25),emisice(2)
251      write(*,5) '(22)    albedice(1)',tab_cntrl(tab0+22),albedice(1)
252      write(*,5) '(23)    albedice(2)',tab_cntrl(tab0+23),albedice(2)
253      write(*,6) '(31)   iceradius(1)',tab_cntrl(tab0+31),iceradius(1)
254      write(*,6) '(32)   iceradius(2)',tab_cntrl(tab0+32),iceradius(2)
255      write(*,5) '(33)   dtemisice(1)',tab_cntrl(tab0+33),dtemisice(1)
256      write(*,5) '(34)   dtemisice(2)',tab_cntrl(tab0+34),dtemisice(2)
257
258      write(*,5) '(27)         tauvis',tab_cntrl(tab0+27),tauvis
259
260      write(*,5) '(35)        volcapa',tab_cntrl(tab0+35),volcapa
261
262      write(*,*)
263      write(*,*) 'Lmodif in tabfi!!!!!!!',Lmodif
264
265c-----------------------------------------------------------------------
266c        Modifications...
267c-----------------------------------------------------------------------
268
269      IF(Lmodif.eq.1) then
270
271      write(*,*)
272      write(*,*) 'Change values in tab_cntrl ? :'
273      write(*,*) '~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
274      write(*,*) '(Current values given above)'
275      write(*,*)
276      write(*,*) '(3)          day_ini : Initial day (=0 at Ls=0)'
277      write(*,*) '(19)              z0 : default surface roughness (m)'
278      write(*,*) '(21)       emin_turb :  minimal energy (PBL)'
279      write(*,*) '(20)         lmixmin : mixing length (PBL)'
280      write(*,*) '(26)         emissiv : ground emissivity'
281      write(*,*) '(24 et 25)   emisice : CO2 ice max emissivity '
282      write(*,*) '(22 et 23)  albedice : CO2 ice cap albedos'
283      write(*,*) '(31 et 32) iceradius : mean scat radius of CO2 snow'
284      write(*,*) '(33 et 34) dtemisice : time scale for snow',
285     &           ' metamorphism'
286      write(*,*) '(27)        tauvis : mean dust vis. reference ',
287     &           'opacity'
288      write(*,*) '(35)         volcapa : soil volumetric heat capacity'
289      write(*,*) '(18)        obliquit : planet obliquity (deg)'
290      write(*,*) '(17)      peri_day : perihelion date (sol since Ls=0)'
291      write(*,*) '(  )      peri_ls : perihelion date (Ls since Ls=0)'
292      write(*,*) '(15)      periheli : min. sun-mars dist (Mkm)'
293      write(*,*) '(16)      aphelie  : max. sun-mars dist (Mkm)'
294      write(*,*)
295 
296 
297      do ! neverending loop
298        write(*,*)
299        write(*,*)
300        write(*,*) 'Changes to perform ?'
301        write(*,*) '   (enter keyword or return )'
302        write(*,*)
303        read(*,fmt='(a20)') modif
304        if (modif(1:1) .eq. ' ') goto 999
305 
306        write(*,*)
307        write(*,*) trim(modif) , ' : '
308
309        if (trim(modif) .eq. 'day_ini') then
310          write(*,*) 'current value:',day_ini
311          write(*,*) 'enter new value:'
312 101      read(*,*,iostat=ierr) day_ini
313          if(ierr.ne.0) goto 101
314          write(*,*) ' '
315          write(*,*) 'day_ini (new value):',day_ini
316
317        else if (trim(modif) .eq. 'z0') then
318          write(*,*) 'current value (m):',z0_default
319          write(*,*) 'enter new value (m):'
320 102      read(*,*,iostat=ierr) z0_default
321          if(ierr.ne.0) goto 102
322          write(*,*) ' '
323          write(*,*) ' z0 (new value):',z0_default
324
325        else if (trim(modif) .eq. 'emin_turb') then
326          write(*,*) 'current value:',emin_turb
327          write(*,*) 'enter new value:'
328 103      read(*,*,iostat=ierr) emin_turb
329          if(ierr.ne.0) goto 103
330          write(*,*) ' '
331          write(*,*) ' emin_turb (new value):',emin_turb
332
333        else if (trim(modif) .eq. 'lmixmin') then
334          write(*,*) 'current value:',lmixmin
335          write(*,*) 'enter new value:'
336 104      read(*,*,iostat=ierr) lmixmin
337          if(ierr.ne.0) goto 104
338          write(*,*) ' '
339          write(*,*) ' lmixmin (new value):',lmixmin
340
341        else if (trim(modif) .eq. 'emissiv') then
342          write(*,*) 'current value:',emissiv
343          write(*,*) 'enter new value:'
344 105      read(*,*,iostat=ierr) emissiv
345          if(ierr.ne.0) goto 105
346          write(*,*) ' '
347          write(*,*) ' emissiv (new value):',emissiv
348
349        else if (trim(modif) .eq. 'emisice') then
350          write(*,*) 'current value emisice(1) North:',emisice(1)
351          write(*,*) 'enter new value:'
352 106      read(*,*,iostat=ierr) emisice(1)
353          if(ierr.ne.0) goto 106
354          write(*,*)
355          write(*,*) ' emisice(1) (new value):',emisice(1)
356          write(*,*)
357
358          write(*,*) 'current value emisice(2) South:',emisice(2)
359          write(*,*) 'enter new value:'
360 107      read(*,*,iostat=ierr) emisice(2)
361          if(ierr.ne.0) goto 107
362          write(*,*)
363          write(*,*) ' emisice(2) (new value):',emisice(2)
364
365        else if (trim(modif) .eq. 'albedice') then
366          write(*,*) 'current value albedice(1) North:',albedice(1)
367          write(*,*) 'enter new value:'
368 108      read(*,*,iostat=ierr) albedice(1)
369          if(ierr.ne.0) goto 108
370          write(*,*)
371          write(*,*) ' albedice(1) (new value):',albedice(1)
372          write(*,*)
373
374          write(*,*) 'current value albedice(2) South:',albedice(2)
375          write(*,*) 'enter new value:'
376 109      read(*,*,iostat=ierr) albedice(2)
377          if(ierr.ne.0) goto 109
378          write(*,*)
379          write(*,*) ' albedice(2) (new value):',albedice(2)
380
381        else if (trim(modif) .eq. 'iceradius') then
382          write(*,*) 'current value iceradius(1) North:',iceradius(1)
383          write(*,*) 'enter new value:'
384 110      read(*,*,iostat=ierr) iceradius(1)
385          if(ierr.ne.0) goto 110
386          write(*,*)
387          write(*,*) ' iceradius(1) (new value):',iceradius(1)
388          write(*,*)
389
390          write(*,*) 'current value iceradius(2) South:',iceradius(2)
391          write(*,*) 'enter new value:'
392 111      read(*,*,iostat=ierr) iceradius(2)
393          if(ierr.ne.0) goto 111
394          write(*,*)
395          write(*,*) ' iceradius(2) (new value):',iceradius(2)
396
397        else if (trim(modif) .eq. 'dtemisice') then
398          write(*,*) 'current value dtemisice(1) North:',dtemisice(1)
399          write(*,*) 'enter new value:'
400 112      read(*,*,iostat=ierr) dtemisice(1)
401          if(ierr.ne.0) goto 112
402          write(*,*)
403          write(*,*) ' dtemisice(1) (new value):',dtemisice(1)
404          write(*,*)
405
406          write(*,*) 'current value dtemisice(2) South:',dtemisice(2)
407          write(*,*) 'enter new value:'
408 113      read(*,*,iostat=ierr) dtemisice(2)
409          if(ierr.ne.0) goto 113
410          write(*,*)
411          write(*,*) ' dtemisice(2) (new value):',dtemisice(2)
412
413        else if (trim(modif) .eq. 'tauvis') then
414          write(*,*) 'current value:',tauvis
415          write(*,*) 'enter new value:'
416 114      read(*,*,iostat=ierr) tauvis
417          if(ierr.ne.0) goto 114
418          write(*,*)
419          write(*,*) ' tauvis (new value):',tauvis
420
421        else if (trim(modif) .eq. 'obliquit') then
422          write(*,*) 'current value:',obliquit
423          write(*,*) 'obliquit should be 25.19 on current Mars'
424          write(*,*) 'enter new value:'
425 115      read(*,*,iostat=ierr) obliquit
426          if(ierr.ne.0) goto 115
427          write(*,*)
428          write(*,*) ' obliquit (new value):',obliquit
429
430        else if (trim(modif) .eq. 'peri_day') then
431          write(*,*) 'current value:',peri_day
432          write(*,*) 'peri_day should be 485 sols on current Mars'
433          write(*,*) 'enter new value:'
434 116      read(*,*,iostat=ierr) peri_day
435          if(ierr.ne.0) goto 116
436          write(*,*)
437          write(*,*) ' peri_day (new value):',peri_day
438         
439        else if (trim(modif) .eq. 'peri_ls') then
440          write(*,*) 'peri_ls value is not stored in start files,'
441          write(*,*) 'but it should be 251 degrees on current Mars'
442          write(*,*) '(peri_day should be 485 sols on current Mars)'
443          write(*,*) 'enter new value:'
444 1160     read(*,*,iostat=ierr) peri_ls
445          if(ierr.ne.0) goto 1160
446          write(*,*)
447          write(*,*) 'peri_ls asked:',peri_ls
448          write(*,*) 'for aphelion =',aphelie
449          write(*,*) 'perihelion =',periheli
450          write(*,*) 'and',year_day,'sols/year'
451          call lsp2solp(peri_ls,peri_day,aphelie,periheli,year_day)
452          write(*,*) 'peri_day (new value):',peri_day
453
454
455        else if (trim(modif) .eq. 'periheli') then
456          write(*,*) 'current value:',periheli
457          write(*,*) 'perihelion should be 206.66 on current Mars'
458          write(*,*) 'enter new value:'
459 117      read(*,*,iostat=ierr) periheli
460          if(ierr.ne.0) goto 117
461          write(*,*)
462          write(*,*) ' periheli (new value):',periheli
463 
464        else if (trim(modif) .eq. 'aphelie') then
465          write(*,*) 'current value:',aphelie
466          write(*,*) 'aphelion should be 249.22 on current Mars'
467          write(*,*) 'enter new value:'
468 118      read(*,*,iostat=ierr) aphelie
469          if(ierr.ne.0) goto 118
470          write(*,*)
471          write(*,*) ' aphelie (new value):',aphelie
472 
473        else if (trim(modif) .eq. 'volcapa') then
474          write(*,*) 'current value:',volcapa
475          write(*,*) 'enter new value:'
476 119      read(*,*,iostat=ierr) volcapa
477          if(ierr.ne.0) goto 119
478          write(*,*)
479          write(*,*) ' volcapa (new value):',volcapa
480 
481        endif
482      enddo ! of do ! neverending loop
483
484 999  continue
485
486c-----------------------------------------------------------------------
487c       Write values of physical constants after modifications
488c-----------------------------------------------------------------------
489 
490      write(*,*) '*****************************************************'
491      write(*,*) 'Reading tab_cntrl when calling tabfi AFTER changes'
492      write(*,*) '*****************************************************'
493      write(*,5) '(1)      = ngridmx?',tab_cntrl(tab0+1),real(ngridmx)
494      write(*,5) '(2)            lmax',tab_cntrl(tab0+2),real(lmax)
495      write(*,5) '(3)         day_ini',tab_cntrl(tab0+3),real(day_ini)
496      write(*,5) '(5)             rad',tab_cntrl(tab0+5),rad
497      write(*,5) '(10)         daysec',tab_cntrl(tab0+10),daysec
498      write(*,6) '(6)            omeg',tab_cntrl(tab0+6),omeg
499      write(*,5) '(7)               g',tab_cntrl(tab0+7),g
500      write(*,5) '(8)           mugaz',tab_cntrl(tab0+8),mugaz
501      write(*,5) '(9)             rcp',tab_cntrl(tab0+9),rcp
502      write(*,6) '(11)        dtphys?',tab_cntrl(tab0+11),dtphys
503 
504      write(*,5) '(14)       year_day',tab_cntrl(tab0+14),year_day
505      write(*,5) '(15)       periheli',tab_cntrl(tab0+15),periheli
506      write(*,5) '(16)        aphelie',tab_cntrl(tab0+16),aphelie
507      write(*,5) '(17)       peri_day',tab_cntrl(tab0+17),peri_day
508      write(*,5) '(18)       obliquit',tab_cntrl(tab0+18),obliquit
509 
510      write(*,6) '(19)     z0_default',tab_cntrl(tab0+19),z0_default
511      write(*,6) '(21)      emin_turb',tab_cntrl(tab0+21),emin_turb
512      write(*,5) '(20)        lmixmin',tab_cntrl(tab0+20),lmixmin
513 
514      write(*,5) '(26)        emissiv',tab_cntrl(tab0+26),emissiv
515      write(*,5) '(24)     emisice(1)',tab_cntrl(tab0+24),emisice(1)
516      write(*,5) '(25)     emisice(2)',tab_cntrl(tab0+25),emisice(2)
517      write(*,5) '(22)    albedice(1)',tab_cntrl(tab0+22),albedice(1)
518      write(*,5) '(23)    albedice(2)',tab_cntrl(tab0+23),albedice(2)
519      write(*,6) '(31)   iceradius(1)',tab_cntrl(tab0+31),iceradius(1)
520      write(*,6) '(32)   iceradius(2)',tab_cntrl(tab0+32),iceradius(2)
521      write(*,5) '(33)   dtemisice(1)',tab_cntrl(tab0+33),dtemisice(1)
522      write(*,5) '(34)   dtemisice(2)',tab_cntrl(tab0+34),dtemisice(2)
523 
524      write(*,5) '(27)         tauvis',tab_cntrl(tab0+27),tauvis
525
526      write(*,5) '(35)        volcapa',tab_cntrl(tab0+35),volcapa
527
528      write(*,*) 
529      write(*,*)
530
531      ENDIF                     !       of if (Lmodif == 1)
532
533c-----------------------------------------------------------------------
534c       Case when using a start file from before March 1996 (without iceradius...
535c-----------------------------------------------------------------------
536      if (iceradius(1).eq.0) then
537         iceradius(1) = 100.e-6
538         iceradius(2) = 100.e-6
539         dtemisice(1) = 0.4
540         dtemisice(2) = 0.4
541         write (*,*) ' tabfi: WARNING : old initialisation file'
542         write (*,*) 'iceradius set to',iceradius(1),iceradius(2) 
543         write (*,*) 'dtemisice set to',dtemisice(1),dtemisice(2) 
544       end if
545
546c-----------------------------------------------------------------------
547      end
548
549
550
551     
552!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
553! gives sol at perihelion for ls at perihelion (for precession cycles)
554      subroutine lsp2solp(lsp,solp,aphelie,periheli,year_day)
555
556      implicit none
557!  Arguments:
558      real lsp     ! Input: ls at perihelion
559      real solp    ! Output: sol at perihelion
560      real aphelie,periheli,year_day ! Input: parameters
561 
562!  Local:
563      double precision zx0 ! eccentric anomaly at Ls=0
564      double precision e_elips
565      double precision pi,degrad
566     
567      parameter (pi=3.14159265358979d0)
568      parameter (degrad=57.2957795130823d0)
569
570      e_elips=(aphelie-periheli)/(aphelie+periheli)     
571      zx0 = -2.0*datan(dtan(0.5*lsp/degrad)
572     .          *dsqrt((1.-e_elips)/(1.+e_elips)))
573      if (zx0 .le. 0.) zx0 = zx0 + 2.*pi
574     
575      solp  = year_day*(1.-(zx0-e_elips*dsin(zx0))/(2.*pi))
576
577
578      end subroutine lsp2solp
579
580
581
Note: See TracBrowser for help on using the repository browser.