Changeset 669 for trunk/LMDZ.MARS/libf/phymars/tabfi.F
- Timestamp:
- May 24, 2012, 5:02:41 PM (13 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.MARS/libf/phymars/tabfi.F
r552 r669 67 67 INTEGER Lmodif 68 68 INTEGER lmax 69 REAL p_rad,p_omeg,p_g,p_mugaz,p_daysec,time 69 REAL p_rad,p_omeg,p_g,p_mugaz,p_daysec,time,peri_ls 70 70 71 71 c Variables … … 289 289 write(*,*) '(18) obliquit : planet obliquity (deg)' 290 290 write(*,*) '(17) peri_day : perihelion date (sol since Ls=0)' 291 write(*,*) '( ) peri_ls : perihelion date (Ls since Ls=0)' 291 292 write(*,*) '(15) periheli : min. sun-mars dist (Mkm)' 292 293 write(*,*) '(16) aphelie : max. sun-mars dist (Mkm)' … … 429 430 else if (trim(modif) .eq. 'peri_day') then 430 431 write(*,*) 'current value:',peri_day 431 write(*,*) 'peri_day should be 485 on current Mars'432 write(*,*) 'peri_day should be 485 sols on current Mars' 432 433 write(*,*) 'enter new value:' 433 434 116 read(*,*,iostat=ierr) peri_day … … 435 436 write(*,*) 436 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 call lsp2solp(peri_ls,peri_day) 449 write(*,*) ' peri_day (new value):',peri_day 450 437 451 438 452 else if (trim(modif) .eq. 'periheli') then … … 529 543 c----------------------------------------------------------------------- 530 544 end 545 546 547 548 549 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 550 ! gives sol at perihelion for ls at perihelion (for precession cycles) 551 subroutine lsp2solp(lsp,solp) 552 553 implicit none 554 ! Arguments: 555 real lsp ! Input: ls at perihelion 556 real solp ! Output: sol at perihelion 557 558 ! Local: 559 double precision zx0 ! eccentric anomaly at Ls=0 560 double precision year_day 561 double precision e_elips 562 double precision pi,degrad 563 564 parameter (year_day=668.6d0) ! number of sols in a martian year 565 parameter (e_elips=0.0934d0) ! eccentricity of orbit 566 parameter (pi=3.14159265358979d0) 567 parameter (degrad=57.2957795130823d0) 568 569 zx0 = -2.0*datan(dtan(0.5*lsp/degrad) 570 . *dsqrt((1.-e_elips)/(1.+e_elips))) 571 if (zx0 .le. 0.) zx0 = zx0 + 2.*pi 572 573 solp = year_day*(1.-(zx0-e_elips*dsin(zx0))/(2.*pi)) 574 575 576 end subroutine lsp2solp 577 578 579
Note: See TracChangeset
for help on using the changeset viewer.