| 1 | MODULE recomp_orb_param_mod |
|---|
| 2 | !======================================================================= |
|---|
| 3 | ! |
|---|
| 4 | ! Purpose: Recompute orbit parameters based on values by Laskar et al., |
|---|
| 5 | ! |
|---|
| 6 | ! Author: RV, JBC |
|---|
| 7 | !======================================================================= |
|---|
| 8 | IMPLICIT NONE |
|---|
| 9 | |
|---|
| 10 | CONTAINS |
|---|
| 11 | |
|---|
| 12 | SUBROUTINE recomp_orb_param(i_myear,year_iter) |
|---|
| 13 | |
|---|
| 14 | use time_evol_mod, only: year_bp_ini, var_obl, var_ecc, var_lsp, convert_years |
|---|
| 15 | #ifndef CPP_STD |
|---|
| 16 | use comconst_mod, only: pi |
|---|
| 17 | use planete_h, only: e_elips, obliquit, lsperi, periheli, aphelie, p_elips, peri_day, year_day |
|---|
| 18 | #else |
|---|
| 19 | use planete_mod, only: e_elips, obliquit, lsperi, periastr, apoastr, p_elips, peri_day, year_day |
|---|
| 20 | use comcstfi_mod, only: pi |
|---|
| 21 | #endif |
|---|
| 22 | |
|---|
| 23 | use lask_param_mod, only: yearlask, obllask, ecclask, lsplask, end_lask_param_mod, last_ilask |
|---|
| 24 | |
|---|
| 25 | IMPLICIT NONE |
|---|
| 26 | |
|---|
| 27 | !-------------------------------------------------------- |
|---|
| 28 | ! Input Variables |
|---|
| 29 | !-------------------------------------------------------- |
|---|
| 30 | integer, intent(in) :: i_myear ! Number of simulated Martian years |
|---|
| 31 | integer, intent(in) :: year_iter ! Number of iterations of the PEM |
|---|
| 32 | |
|---|
| 33 | !-------------------------------------------------------- |
|---|
| 34 | ! Output Variables |
|---|
| 35 | !-------------------------------------------------------- |
|---|
| 36 | |
|---|
| 37 | !-------------------------------------------------------- |
|---|
| 38 | ! Local variables |
|---|
| 39 | !-------------------------------------------------------- |
|---|
| 40 | integer :: Year ! Year of the simulation |
|---|
| 41 | real :: Lsp ! Ls perihelion |
|---|
| 42 | integer :: ilask ! Loop variable |
|---|
| 43 | real :: halfaxe ! Million of km |
|---|
| 44 | real :: unitastr |
|---|
| 45 | real :: xa, xb, ya, yb ! Variables for interpolation |
|---|
| 46 | logical :: found_year ! Flag variable |
|---|
| 47 | |
|---|
| 48 | ! ********************************************************************** |
|---|
| 49 | ! 0. Initializations |
|---|
| 50 | ! ********************************************************************** |
|---|
| 51 | #ifdef CPP_STD |
|---|
| 52 | real aphelie |
|---|
| 53 | real periheli |
|---|
| 54 | |
|---|
| 55 | aphelie=apoastr |
|---|
| 56 | periheli=periastr |
|---|
| 57 | #endif |
|---|
| 58 | |
|---|
| 59 | Year = year_bp_ini + i_myear ! We compute the current year |
|---|
| 60 | Lsp = lsperi*180./pi ! We convert in degree |
|---|
| 61 | |
|---|
| 62 | write(*,*) 'recomp_orb_param, Old year =', Year - year_iter |
|---|
| 63 | write(*,*) 'recomp_orb_param, Old obl. =', obliquit |
|---|
| 64 | write(*,*) 'recomp_orb_param, Old ecc. =', e_elips |
|---|
| 65 | write(*,*) 'recomp_orb_param, Old Lsp =', Lsp |
|---|
| 66 | |
|---|
| 67 | found_year = .false. |
|---|
| 68 | do ilask = last_ilask,2,-1 |
|---|
| 69 | xa = yearlask(ilask) |
|---|
| 70 | xb = yearlask(ilask - 1) |
|---|
| 71 | if (xa <= Year .and. Year < xb) then |
|---|
| 72 | ! Obliquity |
|---|
| 73 | if (var_obl) then |
|---|
| 74 | ya = obllask(ilask) |
|---|
| 75 | yb = obllask(ilask - 1) |
|---|
| 76 | obliquit = (Year - xa)*(yb - ya)/(xb - xa) + ya |
|---|
| 77 | endif |
|---|
| 78 | |
|---|
| 79 | ! Eccentricity |
|---|
| 80 | if (var_ecc) then |
|---|
| 81 | ya = ecclask(ilask) |
|---|
| 82 | yb = ecclask(ilask - 1) |
|---|
| 83 | e_elips = (Year - xa)*(yb - ya)/(xb - xa) + ya |
|---|
| 84 | endif |
|---|
| 85 | |
|---|
| 86 | ! Lsp |
|---|
| 87 | if (var_lsp) then |
|---|
| 88 | ya = lsplask(ilask) |
|---|
| 89 | yb = lsplask(ilask - 1) |
|---|
| 90 | if (abs(yb - ya) > 300.) then ! If modulo is "crossed" through the interval |
|---|
| 91 | if (ya < yb) then ! Lsp should be decreasing |
|---|
| 92 | ya = ya + 360. |
|---|
| 93 | else ! Lsp should be increasing |
|---|
| 94 | yb = yb + 360. |
|---|
| 95 | endif |
|---|
| 96 | endif |
|---|
| 97 | Lsp = modulo((Year - xa)*(yb - ya)/(xb - xa) + ya,360.) |
|---|
| 98 | endif |
|---|
| 99 | found_year = .true. |
|---|
| 100 | exit ! The loop is left as soon as the right interval is found |
|---|
| 101 | endif |
|---|
| 102 | enddo |
|---|
| 103 | |
|---|
| 104 | if (.not. found_year) error stop 'The new year could not be found in the file "obl_ecc_lsp.asc".' |
|---|
| 105 | |
|---|
| 106 | halfaxe = 227.94 |
|---|
| 107 | Lsp = Lsp*pi/180. |
|---|
| 108 | periheli = halfaxe*(1 - e_elips) |
|---|
| 109 | aphelie = halfaxe*(1 + e_elips) |
|---|
| 110 | unitastr = 149.597927 |
|---|
| 111 | p_elips = 0.5*(periheli + aphelie)*(1. - e_elips*e_elips)/unitastr |
|---|
| 112 | |
|---|
| 113 | #ifndef CPP_STD |
|---|
| 114 | call call_dayperi(Lsp,e_elips,peri_day,year_day) |
|---|
| 115 | #endif |
|---|
| 116 | |
|---|
| 117 | write(*,*) 'recomp_orb_param, New year =', Year |
|---|
| 118 | write(*,*) 'recomp_orb_param, New obl. =', obliquit |
|---|
| 119 | write(*,*) 'recomp_orb_param, New ecc. =', e_elips |
|---|
| 120 | write(*,*) 'recomp_orb_param, New Lsp =', Lsp |
|---|
| 121 | |
|---|
| 122 | END SUBROUTINE recomp_orb_param |
|---|
| 123 | |
|---|
| 124 | !******************************************************************************** |
|---|
| 125 | |
|---|
| 126 | END MODULE recomp_orb_param_mod |
|---|
| 127 | |
|---|