Changeset 1329 for trunk/LMDZ.GENERIC


Ignore:
Timestamp:
Aug 21, 2014, 6:20:57 PM (10 years ago)
Author:
jleconte
Message:

21/08/2014 == JL

  • in the ring shadow calculation, temporary variables for the orbbital parameters are now used (we want the declin and RA output at the same time w and wo rings)
Location:
trunk/LMDZ.GENERIC
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/README

    r1326 r1329  
    10461046- fixed variable type declaration bug in ave_stelspec
    10471047- corrected insolation calculation in locked and inclined case (proper calculation of sub solar longitude and declination)
     1048- in the ring shadow calculation, temporary variables for the orbbital parameters are now used (we want the declin and RA output at the same time w and wo rings)
  • trunk/LMDZ.GENERIC/libf/phystd/physiq.F90

    r1327 r1329  
    452452      ! temporary variables computed at each step of this average
    453453      real :: ptime_day ! Universal time in sol fraction
    454       real :: tmp_zls   ! solar longitude
     454      real :: tmp_zls,tmp_dist_star, tmp_declin, tmp_right_ascen   ! tmp solar longitude, stellar dist, declin and RA
    455455      real, dimension(:), allocatable :: tmp_fract ! day fraction of the time interval
    456456      real, dimension(:), allocatable :: tmp_mu0 ! equivalent solar angle
     
    762762      zday=pday+ptime ! compute time, in sols (and fraction thereof)
    763763
    764 !     Compute Stellar Longitude (Ls)
     764!     Compute Stellar Longitude (Ls), and orbital parameters
    765765!     ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    766766      if (season) then
     
    770770      end if
    771771
     772      call orbite(zls,dist_star,declin,right_ascen)
     773      if (tlocked) then
     774              zlss=Mod(-(2.*pi*(zday/year_day)*nres - right_ascen),2.*pi)
     775      elseif (diurnal) then
     776              zlss=-2.*pi*(zday-.5)
     777      else if(diurnal .eqv. .false.) then
     778              zlss=9999.
     779      endif
    772780
    773781
     
    845853!          Compute local stellar zenith angles
    846854!          ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    847            call orbite(zls,dist_star,declin,right_ascen)
    848855
    849856           if (tlocked) then
    850857! JL14 corrects tidally resonant (and inclined) cases. nres=omega_rot/omega_orb
    851               zlss=Mod(-(2.*pi*(zday/year_day)*nres - right_ascen),2.*pi)
    852858              ztim1=SIN(declin)
    853859              ztim2=COS(declin)*COS(zlss)
    854860              ztim3=COS(declin)*SIN(zlss)
    855861
    856 
    857862              call stelang(ngrid,sinlon,coslon,sinlat,coslat,    &
    858863               ztim1,ztim2,ztim3,mu0,fract, flatten)
    859864
    860865           elseif (diurnal) then
    861               zlss=-2.*pi*(zday-.5)
    862866              ztim1=SIN(declin)
    863867              ztim2=COS(declin)*COS(2.*pi*(zday-.5))
     
    890894                        ptime_day = m*pas
    891895                        call stellarlong(pday+ptime_day,tmp_zls)
    892                         call orbite(tmp_zls,dist_star,declin)
    893                         ztim1=SIN(declin)
    894                         ztim2=COS(declin)*COS(2.*pi*(pday+ptime_day-.5))
    895                         ztim3=-COS(declin)*SIN(2.*pi*(pday+ptime_day-.5))
     896                        call orbite(tmp_zls,tmp_dist_star,tmp_declin,tmp_right_ascen)
     897                        ztim1=SIN(tmp_declin)
     898                        ztim2=COS(tmp_declin)*COS(2.*pi*(pday+ptime_day-.5))
     899                        ztim3=-COS(tmp_declin)*SIN(2.*pi*(pday+ptime_day-.5))
     900
    896901                 
    897902                        call stelang(ngrid,sinlon,coslon,sinlat,coslat,    &
    898903                                 ztim1,ztim2,ztim3,tmp_mu0,tmp_fract, flatten)
    899904                 
    900                         call rings(ngrid, declin, ptime_day, rad, flatten, eclipse)
     905                        call rings(ngrid, tmp_declin, ptime_day, rad, flatten, eclipse)
    901906                 
    902907                        fract(:) = fract(:) + (1.-eclipse(:))*tmp_fract(:) !! fract prend en compte l'ombre des anneaux et l'alternance jour/nuit
Note: See TracChangeset for help on using the changeset viewer.