- Timestamp:
- Jan 21, 2014, 4:07:42 PM (11 years ago)
- Location:
- trunk/LMDZ.GENERIC/libf/phystd
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/callcorrk.F90
r1147 r1161 718 718 719 719 if(fract(ig) .ge. 1.0e-4) then ! only during daylight! 720 720 ! if(sum(eclipse) .eq. 0.) then 721 ! write(*,*) 'but eclipse is still 0 in physiq !' 722 ! endif 721 723 if((ngrid.eq.1).and.(global1d))then 722 724 do nw=1,L_NSPECTV 723 stel_fract(nw)= stel(nw) 725 stel_fract(nw)= stel(nw)* 0.25 / acosz 724 726 ! globally averaged = divide by 4 725 727 ! but we correct for solar zenith angle … … 727 729 else 728 730 do nw=1,L_NSPECTV 729 stel_fract(nw)= stel(nw) * fract(ig) *(1.-eclipse(ig))731 stel_fract(nw)= stel(nw) * fract(ig) 730 732 end do 731 endif 732 733 endif 733 734 call optcv(dtauv,tauv,taucumv,plevrad, & 734 735 qxvaer,qsvaer,gvaer,wbarv,cosbv,tauray,tauaero, & -
trunk/LMDZ.GENERIC/libf/phystd/callkeys.h
r1151 r1161 5 5 COMMON/callkeys/callrad,corrk,calldifv,UseTurbDiff,calladj & 6 6 & , co2cond,callsoil & 7 & , season,diurnal,tlocked, iradia,lwrite&7 & , season,diurnal,tlocked,rings_shadow,iradia,lwrite & 8 8 & , iaervar,iddist,topdustref,callstats,calleofdump & 9 9 & , enertest & … … 39 39 & , cloudlvl & 40 40 & , pceil & 41 & , strictboundcorrk & 42 & , rings_shadow 41 & , strictboundcorrk 43 42 44 43 logical callrad,corrk,calldifv,UseTurbDiff & 45 44 & , calladj,co2cond,callsoil & 46 & , season,diurnal,tlocked, lwrite &45 & , season,diurnal,tlocked,rings_shadow,lwrite & 47 46 & , callstats,calleofdump & 48 47 & , callgasvis,continuum,H2Ocont_simple,graybody & 49 & , strictboundcorrk & 50 & , rings_shadow 48 & , strictboundcorrk 51 49 52 50 logical enertest -
trunk/LMDZ.GENERIC/libf/phystd/physiq.F90
r1157 r1161 199 199 character*80 fichier 200 200 integer l,ig,ierr,iq,iaer 201 201 202 202 !!! this is saved for diagnostic 203 203 real,dimension(:),allocatable,save :: fluxsurf_lw ! incident LW (IR) surface flux (W.m-2) … … 411 411 logical,save :: ice_update 412 412 413 ! included by MS to compute the daily average of rings shadowing 414 integer, parameter :: nb_hours = 1536 ! set how many times per day are used 415 real :: pas 416 integer :: m 417 ! temporary variables computed at each step of this average 418 real :: ptime_day ! Universal time in sol fraction 419 real :: tmp_zls ! solar longitude 420 real, dimension(:), allocatable :: tmp_fract ! day fraction of the time interval 421 real, dimension(:), allocatable :: tmp_mu0 ! equivalent solar angle 413 422 !======================================================================= 414 423 … … 455 464 ALLOCATE(tau_col(ngrid)) 456 465 457 !! this is defined in comsaison_h 458 ALLOCATE(mu0(ngrid)) 459 ALLOCATE(fract(ngrid)) 460 466 !! this is defined in comsaison_h 467 ALLOCATE(mu0(ngrid)) 468 ALLOCATE(fract(ngrid)) 469 470 461 471 !! this is defined in radcommon_h 462 472 ALLOCATE(eclipse(ngrid)) … … 681 691 ztim1=SIN(declin) 682 692 ! ztim2=COS(declin)*COS(2.*pi*(zday/year_day) - zls*nres) 683 ! ztim3=-COS(declin)*SIN(2. *pi*(zday/year_day) - zls*nres)693 ! ztim3=-COS(declin)*SIN(2. zday/year_day) - zls*nres) 684 694 ! JL12 corrects tidally resonant cases. nres=omega_rot/omega_orb 685 695 ztim2=COS(declin)*COS(2.*pi*(zday/year_day)*nres - zls) … … 696 706 call stelang(ngrid,sinlon,coslon,sinlat,coslat, & 697 707 ztim1,ztim2,ztim3,mu0,fract) 698 699 else 708 else if(diurnal .eq. .false.) then 700 709 701 710 call mucorr(ngrid,declin,lati,mu0,fract,10000.,rad) 702 703 endif 704 711 ! write(*,*) 'mu0 =', mu0 712 ! write(*,*) 'hello, l 738' 713 ! WARNING: this function appears not to work in 1D 714 715 endif 716 705 717 !! Eclipse incoming sunlight (e.g. Saturn ring shadowing) 706 if (rings_shadow) then 707 if (diurnal) then 708 call rings(ngrid, declin, ptime, rad, 0., eclipse) 709 call writediagfi(ngrid,"shad","rings"," ", 2, 1.-eclipse) 710 write(*,*) "Ring shadow activated." 711 else 712 write(*,*) "Ring shadow + diurnal=F not supported yet." 713 endif 714 else 715 eclipse(:) = 0.0 716 endif 718 719 if(rings_shadow) then 720 write(*,*) 'Rings shadow activated' 721 if(diurnal .eq. .false.) then ! we need to compute the daily average insolation 722 pas = 1./nb_hours 723 ptime_day = 0. 724 fract(:) = 0. 725 ALLOCATE(tmp_fract(ngrid)) 726 ALLOCATE(tmp_mu0(ngrid)) 727 tmp_fract(:) = 0. 728 eclipse(:) = 0. 729 tmp_mu0(:) = 0. 730 731 do m=1, nb_hours 732 ptime_day = m*pas 733 call stellarlong(pday+ptime_day,tmp_zls) 734 call orbite(tmp_zls,dist_star,declin) 735 ztim1=SIN(declin) 736 ztim2=COS(declin)*COS(2.*pi*(pday+ptime_day-.5)) 737 ztim3=-COS(declin)*SIN(2.*pi*(pday+ptime_day-.5)) 738 739 call stelang(ngrid,sinlon,coslon,sinlat,coslat, & 740 ztim1,ztim2,ztim3,tmp_mu0,tmp_fract) 741 742 call rings(ngrid, declin, ptime_day, rad, 0., eclipse) 743 744 fract(:) = fract(:) + (1.-eclipse(:))*tmp_fract(:) !! fract prend en compte l'ombre des anneaux et l'alternance jour/nuit 745 enddo 746 747 DEALLOCATE(tmp_fract) 748 DEALLOCATE(tmp_mu0) 749 750 fract(:) = fract(:)/nb_hours 751 752 else 753 call rings(ngrid, declin, ptime, rad, 0., eclipse) 754 fract(:) = fract(:) * (1.-eclipse) 755 endif 756 else 757 eclipse(:) = 0.0 758 endif 759 717 760 718 761 if (corrk) then … … 1431 1474 if(corrk)then 1432 1475 1476 if(ISNAN(sum(fluxtop_dn))) then 1477 write(*,*) 'err fluxtop_dn' 1478 else if(ISNAN(sum(area))) then 1479 write(*,*) 'err area' 1480 else if(ISNAN(totarea)) then 1481 write(*,*) 'area is nan' 1482 endif 1483 1433 1484 ISR = SUM(area(:)*fluxtop_dn(:))/totarea 1434 1485 ASR = SUM(area(:)*fluxabs_sw(:))/totarea … … 1717 1768 call writediagfi(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw) 1718 1769 call writediagfi(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw) 1770 call writediagfi(ngrid,"shad","rings"," ", 2, fract) 1719 1771 ! call writediagfi(ngrid,"ASRcs","absorbed stellar rad (cs).","W m-2",2,fluxabs_sw1) 1720 1772 ! call writediagfi(ngrid,"OLRcs","outgoing longwave rad (cs).","W m-2",2,fluxtop_lw1) … … 1868 1920 !! this is defined in comsaison_h 1869 1921 IF ( ALLOCATED(mu0)) DEALLOCATE(mu0) 1922 1870 1923 IF ( ALLOCATED(fract)) DEALLOCATE(fract) 1924 1925 1926 !! this is defined in radcommon_h 1927 IF ( ALLOCATED(eclipse)) DEALLOCATE(eclipse) 1871 1928 1872 1929 !! this is defined in comsoil_h -
trunk/LMDZ.GENERIC/libf/phystd/rings.F90
r1133 r1161 40 40 REAL, DIMENSION(:), ALLOCATABLE:: x, y, z ! cartesian coordinates of the points on the planet 41 41 REAL :: xs, ys, zs ! cartesian coordinates of the points of the subsolar point 42 REAL, DIMENSION(:), ALLOCATABLE :: k ! parameter (?)42 REAL, DIMENSION(:), ALLOCATABLE :: k 43 43 REAL, DIMENSION(:), ALLOCATABLE :: N ! parameter to compute cartesian coordinates on a ellipsoidal planet 44 44 REAL, DIMENSION(:), ALLOCATABLE :: r ! distance at which the incident ray of sun crosses the equatorial plane … … 48 48 ! equinox --> no shadow (AS: why is this needed?) 49 49 if(declin .eq. 0.) then 50 eclipse = 0.50 eclipse(:) = 0. 51 51 return 52 52 endif … … 63 63 ALLOCATE(N(ngrid)) 64 64 ALLOCATE(r(ngrid)) 65 eclipse = 2000.65 eclipse(:) = 2000. 66 66 67 67 ! Size of the rings … … 85 85 86 86 ! Convert to cartesian coordinates 87 N = rad/sqrt(1-(e**2)*sinlat**2)88 x = N*coslat*coslon89 y = N*coslat*sinlon90 z = N*(1-e**2)*sinlat87 N(:) = rad / sqrt(1-(e**2)*sinlat(:)**2) 88 x(:) = N(:)*coslat(:)*coslon(:) 89 y(:) = N(:)*coslat(:)*sinlon(:) 90 z(:) = N(:)*(1-e**2)*sinlat(:) 91 91 92 92 ! 2) LOCATION OF THE SUBSOLAR POINT … … 95 95 ! SG: the minus sign is important! ... otherwise subsolar point adopts a reverse rotation 96 96 phi_S = -(ptime - 0.5)*2.*pi 97 write(*,*) 'subsol point coords : ', declin*180./pi, phi_S*180./pi97 ! write(*,*) 'subsol point coords : ', declin*180./pi, phi_S*180./pi 98 98 99 99 ! subsolar latitude is declin (declination of the sun) … … 106 106 ! 3) WHERE DOES THE INCIDENT RAY OF SUN CROSS THE EQUATORIAL PLAN ? 107 107 108 k = -z/zs109 r = (k*xs + x)**2 + (k*ys + y)**2110 r = sqrt(r)108 k(:) = -z(:)/zs 109 r(:) = (k(:)*xs + x(:))**2 + (k(:)*ys + y(:))**2 110 r(:) = sqrt(r(:)) 111 111 112 112 ! 4) SO WHERE ARE THE SHADOW OF THESE RINGS ? 113 113 114 114 ! Summer hemisphere is not under the shadow of the rings 115 where(lati *declin .gt. 0.)116 eclipse = 1000.115 where(lati(:)*declin .gt. 0.) 116 eclipse(:) = 1000. 117 117 end where 118 118 119 119 ! No shadow of the rings by night 120 where(x *xs + y*ys + z*zs .lt. 0.)121 eclipse = 1000.120 where(x(:)*xs + y(:)*ys + z(:)*zs .lt. 0.) 121 eclipse(:) = 1000. 122 122 end where 123 123 124 124 ! if the incident rays of sun cross a ring, there is a shadow 125 125 do i=1, nb_A 126 where(r .ge. A_Rint(i) .and. r .le. A_Rext(i) .and. eclipse.ne. 1000.)127 eclipse = 1. - exp(-tau_A(i)/cos(declin))126 where(r(:) .ge. A_Rint(i) .and. r(:) .le. A_Rext(i) .and. eclipse(:) .ne. 1000.) 127 eclipse(:) = 1. - exp(-tau_A(i)/cos(declin)) 128 128 end where 129 129 end do 130 130 131 131 do i=1, nb_B 132 where(r .ge. B_Rint(i) .and. r .le. B_Rext(i) .and. eclipse.ne. 1000.)133 eclipse = 1. - exp(-tau_B(i)/cos(declin))132 where(r(:) .ge. B_Rint(i) .and. r(:) .le. B_Rext(i) .and. eclipse(:) .ne. 1000.) 133 eclipse(:) = 1. - exp(-tau_B(i)/cos(declin)) 134 134 end where 135 135 enddo 136 136 137 137 do i=1, nb_C 138 where(r .ge. C_Rint(i) .and. r .le. C_Rext(i) .and. eclipse.ne. 1000.)139 eclipse = 1. - exp(-tau_C(i)/cos(declin))138 where(r(:) .ge. C_Rint(i) .and. r(:) .le. C_Rext(i) .and. eclipse(:) .ne. 1000.) 139 eclipse(:) = 1. - exp(-tau_C(i)/cos(declin)) 140 140 end where 141 141 enddo 142 142 143 143 ! At the other places and the excluded ones, eclipse is 0. 144 where(eclipse .eq. 2000. .or. eclipse.eq. 1000.)145 eclipse = 0.144 where(eclipse(:) .eq. 2000. .or. eclipse(:) .eq. 1000.) 145 eclipse(:) = 0. 146 146 end where 147 147
Note: See TracChangeset
for help on using the changeset viewer.