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