Changeset 2228 for trunk/LMDZ.MARS/libf/phymars/dyn1d
- Timestamp:
- Jan 28, 2020, 4:35:21 PM (5 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/dyn1d/testphys1d.F
r2214 r2228 122 122 character(len=44) :: txt 123 123 124 c New flag to compute paleo orbital configurations + few variables JN 125 REAL halfaxe, excentric, Lsperi 126 Logical paleomars 127 124 128 c======================================================================= 125 129 … … 151 155 periheli = 206.66 ! minimum sun-mars distance (Mkm) ~206.66 152 156 aphelie = 249.22 ! maximum sun-mars distance (Mkm) ~249.22 157 halfaxe = 227.94 ! demi-grand axe de l'ellipse 153 158 peri_day = 485. ! perihelion date (sols since N. Spring) 154 159 obliquit = 25.2 ! Obliquity (deg) ~25.2 160 excentric = 0.0934 ! Eccentricity (0.0934) 155 161 156 162 c Planetary Boundary Layer and Turbulence parameters … … 464 470 c Orbital parameters 465 471 c ------------------ 466 print *,'Min. distance Sun-Mars (Mkm)?' 467 call getin("periheli",periheli) 468 write(*,*) " periheli = ",periheli 469 470 print *,'Max. distance Sun-Mars (Mkm)?' 471 call getin("aphelie",aphelie) 472 write(*,*) " aphelie = ",aphelie 473 474 print *,'Day of perihelion?' 475 call getin("periday",peri_day) 476 write(*,*) " periday = ",peri_day 477 478 print *,'Obliquity?' 479 call getin("obliquit",obliquit) 480 write(*,*) " obliquit = ",obliquit 472 paleomars=.false. ! Default: no water ice reservoir 473 call getin("paleomars",paleomars) 474 if (paleomars==.true.) then 475 write(*,*) "paleomars=", paleomars 476 write(*,*) "Orbital parameters from callphys.def" 477 write(*,*) "Enter eccentricity & Lsperi" 478 print *, 'Martian eccentricity (0<e<1) ?' 479 call getin('excentric ',excentric) 480 write(*,*)"excentric =",excentric 481 print *, 'Solar longitude of perihelion (0<Ls<360) ?' 482 call getin('Lsperi',Lsperi ) 483 write(*,*)"Lsperi=",Lsperi 484 Lsperi = Lsperi*pi/180.0 ! Put it in rad for peri_day 485 periheli = halfaxe*(1-excentric) 486 aphelie = halfaxe*(1+excentric) 487 call call_dayperi(Lsperi,excentric,peri_day,year_day) 488 write(*,*) "Corresponding orbital params for GCM" 489 write(*,*) " periheli = ",periheli 490 write(*,*) " aphelie = ",aphelie 491 write(*,*) "date of perihelion (sol)",peri_day 492 else 493 write(*,*) "paleomars=", paleomars 494 write(*,*) "Default present-day orbital parameters" 495 write(*,*) "Unless specified otherwise" 496 print *,'Min. distance Sun-Mars (Mkm)?' 497 call getin("periheli",periheli) 498 write(*,*) " periheli = ",periheli 499 500 print *,'Max. distance Sun-Mars (Mkm)?' 501 call getin("aphelie",aphelie) 502 write(*,*) " aphelie = ",aphelie 503 504 print *,'Day of perihelion?' 505 call getin("periday",peri_day) 506 write(*,*) " periday = ",peri_day 507 508 print *,'Obliquity?' 509 call getin("obliquit",obliquit) 510 write(*,*) " obliquit = ",obliquit 511 end if 481 512 482 513 c latitude/longitude
Note: See TracChangeset
for help on using the changeset viewer.