Changeset 4127 for trunk/LMDZ.GENERIC


Ignore:
Timestamp:
Mar 13, 2026, 10:58:37 AM (4 weeks ago)
Author:
aslmd
Message:

minor fortran coding error in the Rayleigh scattering routine and more major bug in missing albedo initialization in Mesoscale mode

Location:
trunk/LMDZ.GENERIC/libf/phygeneric
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phygeneric/physiq_mod.F90

    r4126 r4127  
    620620
    621621         day_ini = pday
     622         DO nw=1,L_NSPECTV
     623            albedo(:,nw)=albedodat(:)
     624         ENDDO
    622625#endif
    623626
  • trunk/LMDZ.GENERIC/libf/phygeneric/rad_correlatedk_rayleigh_scattering_opacity.F90

    r4081 r4127  
    4848      real, intent(out) :: tauray(L_LEVELS,L_NSPECTV)
    4949      real*8 wl,wn
    50       integer N,Nfine,ifine,igas,k
    51       parameter(Nfine=500.0)
     50      integer N,ifine,igas,k
     51      integer, parameter :: Nfine = 500
    5252      real*8 :: Fk ! King factor for the depolarization
    5353      real*8 :: ng(L_LEVELS) ! real refractive index
     
    131131         endif
    132132         tauvari(igas,:) = 0.
     133         tauconsti(igas,:) = 0.
    133134      enddo
    134135     
     
    149150         bwidth = (10000.0/BWNV(N)) - (10000.0/BWNV(N+1))
    150151         do ifine=1,Nfine
    151             wl=bstart+dble(ifine)*bwidth/Nfine
    152             wn=BWNV(N)+dble(ifine)*(BWNV(N+1)-BWNV(N))/Nfine
    153 
    154             tauvar(:)=0.0
     152            wl=bstart+dble(ifine)*bwidth/dble(Nfine)
     153            wn=BWNV(N)+dble(ifine)*(BWNV(N+1)-BWNV(N))/dble(Nfine)
     154
     155            tauvar(:)=0.0d0
    155156            do igas=1,ngasmx
    156157               if (maxval(mass_frac(igas,:)).ge.1e-2) then
     
    319320                 ! pmid*scalep -> mbar to Pa
    320321                 ! muvar/1000 -> g/mol to kg/mol
    321                
    322322                 tauvar(:)=tauvar(:)+tauconsti(igas,:)*tauvari(igas,:)
    323                  
     323         
    324324               endif !greater than 0.01
    325325
     
    327327
    328328            call rad_blackbody_planck_law_wavelength(dble(wl*1e-6),dble(tstellar),df)
    329             df=df*bwidth/Nfine
     329            df=df*bwidth/dble(Nfine)
    330330            tauwei=tauwei+df
    331331            tausum(:)=tausum(:)+tauvar(:)*df
Note: See TracChangeset for help on using the changeset viewer.