Changeset 1329 for trunk/LMDZ.GENERIC
- Timestamp:
- Aug 21, 2014, 6:20:57 PM (10 years ago)
- Location:
- trunk/LMDZ.GENERIC
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/README
r1326 r1329 1046 1046 - fixed variable type declaration bug in ave_stelspec 1047 1047 - 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 452 452 ! temporary variables computed at each step of this average 453 453 real :: ptime_day ! Universal time in sol fraction 454 real :: tmp_zls ! solar longitude454 real :: tmp_zls,tmp_dist_star, tmp_declin, tmp_right_ascen ! tmp solar longitude, stellar dist, declin and RA 455 455 real, dimension(:), allocatable :: tmp_fract ! day fraction of the time interval 456 456 real, dimension(:), allocatable :: tmp_mu0 ! equivalent solar angle … … 762 762 zday=pday+ptime ! compute time, in sols (and fraction thereof) 763 763 764 ! Compute Stellar Longitude (Ls) 764 ! Compute Stellar Longitude (Ls), and orbital parameters 765 765 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 766 766 if (season) then … … 770 770 end if 771 771 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 772 780 773 781 … … 845 853 ! Compute local stellar zenith angles 846 854 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ 847 call orbite(zls,dist_star,declin,right_ascen)848 855 849 856 if (tlocked) then 850 857 ! 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)852 858 ztim1=SIN(declin) 853 859 ztim2=COS(declin)*COS(zlss) 854 860 ztim3=COS(declin)*SIN(zlss) 855 861 856 857 862 call stelang(ngrid,sinlon,coslon,sinlat,coslat, & 858 863 ztim1,ztim2,ztim3,mu0,fract, flatten) 859 864 860 865 elseif (diurnal) then 861 zlss=-2.*pi*(zday-.5)862 866 ztim1=SIN(declin) 863 867 ztim2=COS(declin)*COS(2.*pi*(zday-.5)) … … 890 894 ptime_day = m*pas 891 895 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 896 901 897 902 call stelang(ngrid,sinlon,coslon,sinlat,coslat, & 898 903 ztim1,ztim2,ztim3,tmp_mu0,tmp_fract, flatten) 899 904 900 call rings(ngrid, declin, ptime_day, rad, flatten, eclipse)905 call rings(ngrid, tmp_declin, ptime_day, rad, flatten, eclipse) 901 906 902 907 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.