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

Last change on this file since 2282 was 2281, checked in by adelavois, 5 years ago

Mars GCM:
Martian physics is now able to start without startfi.nc
Major update for phyetat0_mod and physiq_mod based on what have been done for the Generic physics

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