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