Changeset 1174 for trunk/LMDZ.GENERIC/libf
- Timestamp:
- Feb 3, 2014, 9:39:45 AM (11 years ago)
- Location:
- trunk/LMDZ.GENERIC/libf/phystd
- Files:
-
- 5 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/libf/phystd/callkeys.h
r1161 r1174 12 12 & , meanOLR, specOLR & 13 13 & , kastprof, Tstrat & 14 & , nosurf, intheat 14 & , nosurf, intheat, flatten & 15 15 & , newtonian, tau_relax, testradtimes & 16 16 & , check_cpp_match, force_cpp & … … 30 30 & , tracer, mass_redistrib, varactive, varfixed, satval & 31 31 & , sedimentation,water,watercond,waterrain & 32 & , aeroco2,aeroh2o,aeroh2so4,aeroback2lay 32 & , aeroco2,aeroh2o,aeroh2so4,aeroback2lay, aeronh3 & 33 33 & , aerofixco2,aerofixh2o & 34 34 & , hydrology, sourceevol, icetstep, albedosnow & … … 43 43 logical callrad,corrk,calldifv,UseTurbDiff & 44 44 & , calladj,co2cond,callsoil & 45 & , season,diurnal,tlocked,rings_shadow,lwrite 45 & , season,diurnal,tlocked,rings_shadow,lwrite & 46 46 & , callstats,calleofdump & 47 47 & , callgasvis,continuum,H2Ocont_simple,graybody & … … 69 69 logical water,watercond,waterrain 70 70 logical aeroco2,aeroh2o,aeroh2so4,aeroback2lay 71 logical aerofixco2,aerofixh2o 71 logical aerofixco2,aerofixh2o, aeronh3 72 72 logical hydrology 73 73 logical sourceevol … … 107 107 real icetstep 108 108 real intheat 109 real flatten -
trunk/LMDZ.GENERIC/libf/phystd/inifis.F
r1155 r1174 158 158 write(*,*) "rings_shadow = ", rings_shadow 159 159 160 write(*,*) "Flattening of the planet (a-b)/a ?" 161 flatten = 0.0 162 call getin("flatten", flatten) 163 write(*,*) "flatten = ", flatten 164 160 165 ! Test of incompatibility: 161 166 ! if tlocked, then diurnal should be false … … 448 453 call getin("aeroback2lay",aeroback2lay) 449 454 write(*,*)" aeroback2lay = ",aeroback2lay 455 456 write(*,*)"Radiatively active ammonia cloud?" 457 aeronh3=.false. ! default value 458 call getin("aeronh3",aeronh3) 459 write(*,*)" aeronh3 = ",aeronh3 450 460 451 461 write(*,*)"TWOLAY AEROSOL: total optical depth ", … … 665 675 coslon(ig)=cos(plon(ig)) 666 676 ENDDO 667 677 668 678 pi=2.*asin(1.) ! NB: pi is a common in comcstfi.h 669 679 -
trunk/LMDZ.GENERIC/libf/phystd/mucorr.F
r1152 r1174 1 SUBROUTINE mucorr(npts,pdeclin, plat, pmu, pfract,phaut,prad)1 SUBROUTINE mucorr(npts,pdeclin, plat, pmu,pfract,phaut,prad,pflat) 2 2 IMPLICIT NONE 3 3 … … 34 34 INTEGER npts 35 35 REAL plat(npts),pmu(npts), pfract(npts) 36 REAL phaut,prad,pdeclin 36 REAL phaut,prad,pdeclin, pflat 37 37 c 38 38 c Local variables : … … 41 41 REAL pi 42 42 REAL z,cz,sz,tz,phi,cphi,sphi,tphi 43 REAL ap,a,t,b, tp 43 REAL ap,a,t,b, tp, rap 44 44 REAL alph 45 45 46 46 c----------------------------------------------------------------------- 47 48 c----- SG: geometry adapted to a flattened planet (Feb2014) 47 49 48 50 pi = 4. * atan(1.0) … … 50 52 cz = cos (z) 51 53 sz = sin (z) 54 rap = 1./((1.-pflat)**2) 52 55 53 56 DO 20 j = 1, npts … … 59 62 tphi = sphi / cphi 60 63 b = cphi * cz 61 t = - tphi * sz / cz64 t = -rap*tphi * sz / cz 62 65 a = 1.0 - t*t 63 66 ap = a … … 76 79 t = tp 77 80 78 pmu(j) = (sphi*sz*t ) / pi + b*sin(t)/pi81 pmu(j) = (sphi*sz*t*rap) / pi + b*sin(t)/pi 79 82 pfract(j) = t / pi 80 83 IF (ap .lt.0.) then 81 pmu(j) = sphi * sz 84 pmu(j) = sphi * sz*rap 82 85 pfract(j) = 1.0 83 86 ENDIF … … 85 88 pmu(j) = pmu(j) / pfract(j) 86 89 IF (pmu(j).eq.0.) pfract(j) = 0. 90 91 pmu(j)=pmu(j)/SQRT(cphi**2 + (rap**2) * (sphi**2)) 87 92 88 93 20 CONTINUE -
trunk/LMDZ.GENERIC/libf/phystd/physiq.F90
r1162 r1174 697 697 698 698 call stelang(ngrid,sinlon,coslon,sinlat,coslat, & 699 ztim1,ztim2,ztim3,mu0,fract )699 ztim1,ztim2,ztim3,mu0,fract, flatten) 700 700 701 701 elseif (diurnal) then … … 705 705 706 706 call stelang(ngrid,sinlon,coslon,sinlat,coslat, & 707 ztim1,ztim2,ztim3,mu0,fract )707 ztim1,ztim2,ztim3,mu0,fract, flatten) 708 708 else if(diurnal .eq. .false.) then 709 709 710 call mucorr(ngrid,declin,lati,mu0,fract,10000.,rad )710 call mucorr(ngrid,declin,lati,mu0,fract,10000.,rad,flatten) 711 711 ! WARNING: this function appears not to work in 1D 712 712 … … 736 736 737 737 call stelang(ngrid,sinlon,coslon,sinlat,coslat, & 738 ztim1,ztim2,ztim3,tmp_mu0,tmp_fract )738 ztim1,ztim2,ztim3,tmp_mu0,tmp_fract, flatten) 739 739 740 call rings(ngrid, declin, ptime_day, rad, 0., eclipse)740 call rings(ngrid, declin, ptime_day, rad, flatten, eclipse) 741 741 742 742 fract(:) = fract(:) + (1.-eclipse(:))*tmp_fract(:) !! fract prend en compte l'ombre des anneaux et l'alternance jour/nuit … … 1767 1767 call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrd) 1768 1768 call writediagfi(ngrid,"DYN","dynamical heat input","W m-2",2,fluxdyn) 1769 1769 endif 1770 1770 1771 1771 if(enertest) then -
trunk/LMDZ.GENERIC/libf/phystd/stelang.F
r253 r1174 1 1 subroutine stelang(kgrid,psilon,pcolon,psilat,pcolat, 2 & ptim1,ptim2,ptim3,pmu0,pfract)2 & ptim1,ptim2,ptim3,pmu0,pfract, pflat) 3 3 IMPLICIT NONE 4 4 … … 66 66 C 67 67 INTEGER kgrid 68 REAL ptim1,ptim2,ptim3 68 REAL ptim1,ptim2,ptim3, pflat 69 69 REAL psilon(kgrid),pcolon(kgrid),pmu0(kgrid),pfract(kgrid) 70 70 REAL psilat(kgrid), pcolat(kgrid) 71 71 C 72 72 INTEGER jl 73 REAL ztim1,ztim2,ztim3 73 REAL ztim1,ztim2,ztim3, rap 74 74 C------------------------------------------------------------------------ 75 75 C------------------------------------------------------------------------ … … 81 81 C -------------- 82 82 C 83 c----- SG: geometry adapted to a flattened planet (Feb2014) 84 85 rap = 1./((1.-pflat)**2) 86 83 87 100 CONTINUE 84 88 C … … 92 96 C 93 97 DO jl=1,kgrid 94 ztim1=psilat(jl)*ptim1 98 ztim1=psilat(jl)*ptim1*rap 95 99 ztim2=pcolat(jl)*ptim2 96 100 ztim3=pcolat(jl)*ptim3 97 101 pmu0(jl)=ztim1+ztim2*pcolon(jl)+ztim3*psilon(jl) 102 pmu0(jl)=pmu0(jl)/SQRT(pcolat(jl)**2 + (rap**2) * (psilat(jl)**2)) 103 98 104 ENDDO 99 105 C
Note: See TracChangeset
for help on using the changeset viewer.