Ignore:
Timestamp:
Feb 18, 2026, 10:28:23 AM (6 weeks ago)
Author:
mturbet
Message:

renaming and merging of radiative routines in phygeneric

File:
1 moved

Legend:

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

    r4076 r4077  
    1 MODULE callcorrk_mod
     1MODULE rad_correlatedk_mod
    22
    33IMPLICIT NONE
     
    55CONTAINS
    66
    7       subroutine callcorrk(ngrid,nlayer,pq,nq,qsurf,zls,       &
     7      subroutine rad_correlatedk(ngrid,nlayer,pq,nq,qsurf,zls,       &
    88          albedo,albedo_equivalent,emis,mu0,pplev,pplay,pt,    &
    99          tsurf,fract,dist_star,aerosol,muvar,                 &
     
    2727      use ioipsl_getin_p_mod, only: getin_p
    2828      use gases_h, only: ngasmx
    29       use radii_mod, only : su_aer_radii,co2_reffrad,h2o_reffrad,dust_reffrad,h2so4_reffrad,back2lay_reffrad
    30       use aerosol_mod, only : iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4, &
     29      use aerosol_radius, only : aerosol_radius_init,aerosol_radius_co2,h2o_reffrad, &
     30                    aerosol_radius_dust,aerosol_radius_h2so4,aerosol_radius_back2lay
     31      use aerosol_global_variables , only : iaero_co2,iaero_h2o,iaero_dust,iaero_h2so4, &
    3132                              iaero_back2lay, iaero_aurora,               &
    3233                              iaero_venus1,iaero_venus2,iaero_venus2p,    &
    3334                              iaero_venus3,iaero_venusUV
    34       use aeropacity_mod, only: aeropacity
    35       use aeroptproperties_mod, only: aeroptproperties
     35      use aerosol_opacity_mod, only: aerosol_opacity
     36      use aerosol_optical_properties_mod, only: aerosol_optical_properties
    3637      use tracer_h, only: igcm_h2o_ice, igcm_h2o_vap, igcm_co2_ice
    3738      use tracer_h, only: constants_epsi_generic
     
    4142                              CLFvarying,tplanckmin,tplanckmax,global1d, &
    4243                              generic_condensation, aerovenus, nvarlayer, varspec
    43       use optcv_mod, only: optcv
    44       use optci_mod, only: optci
    45       use sfluxi_mod, only: sfluxi
    46       use sfluxv_mod, only: sfluxv
    47       use recombin_corrk_mod, only: corrk_recombin, call_recombin
     44      use rad_correlatedk_opacities_stellar_mod, &
     45        only: rad_correlatedk_opacities_stellar
     46      use rad_correlatedk_opacities_thermal_mod, &
     47        only: rad_correlatedk_opacities_thermal
     48      use rad_correlatedk_fluxes_thermal_mod, &
     49        only: rad_correlatedk_fluxes_thermal
     50      use rad_correlatedk_fluxes_stellar_mod, &
     51        only: rad_correlatedk_fluxes_stellar
     52      use rad_correlatedk_online_recombination_mod, &
     53        only: corrk_recombin, &
     54        rad_correlatedk_recombination_main
    4855      use pindex_mod, only: pindex
    4956      use generic_cloud_common_h, only: Psat_generic, epsi_generic
     
    111118      REAL,INTENT(OUT) :: OSR_nu(ngrid,L_NSPECTV)        ! Outgoing SW radiation in each band (Normalized to the band width (W/m2/cm-1).
    112119      REAL,INTENT(OUT) :: GSR_nu(ngrid,L_NSPECTV)        ! Surface SW radiation in each band (Normalized to the band width (W/m2/cm-1).
    113       REAL,INTENT(OUT) :: tau_col(ngrid)                 ! Diagnostic from aeropacity.
     120      REAL,INTENT(OUT) :: tau_col(ngrid)                 ! Diagnostic from aerosol_opacity.
    114121      REAL,INTENT(OUT) :: albedo_equivalent(ngrid)       ! Spectrally Integrated Albedo. For Diagnostic. By MT2015
    115122      REAL,INTENT(OUT) :: totcloudfrac(ngrid)            ! Column Fraction of clouds (%).
     
    148155      REAL*8 tlevrad(L_LEVELS),plevrad(L_LEVELS)
    149156
    150       ! Optical values for the optci/cv subroutines
     157      ! Optical values for the rad_correlatedk_opacities_thermal/cv subroutines
    151158      REAL*8 stel(L_NSPECTV),stel_fract(L_NSPECTV)
    152159      ! NB: Arrays below are "save" to avoid reallocating them at every call
     
    200207      character(len=10) :: tmp2
    201208      character(len=100) :: message
    202       character(len=10),parameter :: subname="callcorrk"
     209      character(len=10),parameter :: subname="rad_correlatedk"
    203210
    204211      ! For fixed water vapour profiles.
     
    242249      if(firstcall) then
    243250
    244         ! test on allocated necessary because of CLFvarying (two calls to callcorrk in physiq)
     251        ! test on allocated necessary because of CLFvarying (two calls to rad_correlatedk in physiq)
    245252        if(.not.allocated(QVISsQREF3d)) then
    246253          allocate(QVISsQREF3d(ngrid,nlayer,L_NSPECTV,naerkind))
     
    308315        endif
    309316
    310          !!! ALLOCATED instances are necessary because of CLFvarying (strategy to call callcorrk twice in physiq...)
     317         !!! ALLOCATED instances are necessary because of CLFvarying (strategy to call rad_correlatedk twice in physiq...)
    311318         IF(.not.ALLOCATED(QREFvis3d))THEN
    312319           ALLOCATE(QREFvis3d(ngrid,nlayer,naerkind), stat=ok)
     
    343350#endif
    344351
    345          call su_aer_radii(ngrid,nlayer,reffrad,nueffrad)
     352         call aerosol_radius_init(ngrid,nlayer,reffrad,nueffrad)
    346353         
    347354         
     
    351358
    352359      !this block is now done at firstcall of physiq_mod
    353          ! print*, "callcorrk: Correlated-k data base folder:",trim(datadir)
     360         ! print*, "rad_correlatedk: Correlated-k data base folder:",trim(datadir)
    354361         ! call getin_p("corrkdir",corrkdir)
    355362         ! print*, "corrkdir = ",corrkdir
     
    359366         ! banddir=trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir))
    360367
    361          ! call setspi            ! Basic infrared properties.
    362          ! call setspv            ! Basic visible properties.
    363          ! call sugas_corrk       ! Set up gaseous absorption properties.
    364          ! call suaer_corrk       ! Set up aerosol optical properties.
     368         ! call rad_correlatedk_init_thermal            ! Basic infrared properties.
     369         ! call rad_correlatedk_init_stellar            ! Basic visible properties.
     370         ! call rad_correlatedk_read_opacity_tables        ! Set up gaseous absorption properties.
     371         ! call rad_correlatedk_ini_aerosol       ! Set up aerosol optical properties.
    365372       
    366373
    367          ! now that L_NGAUSS has been initialized (by sugas_corrk)
     374         ! now that L_NGAUSS has been initialized (by rad_correlatedk_read_opacity_tables )
    368375         ! allocate related arrays
    369376         if(.not.allocated(dtaui)) then
     
    446453
    447454         if((igcm_h2o_vap.eq.0) .and. varactive .and. water)then
    448             message='varactive in callcorrk but no h2o_vap tracer.'
     455            message='varactive in rad_correlatedk but no h2o_vap tracer.'
    449456            call abort_physic(subname,message,1)
    450457         endif
     
    485492
    486493         if ((iaer.eq.iaero_co2).and.tracer.and.(igcm_co2_ice.gt.0)) then ! Treat condensed co2 particles.
    487             call co2_reffrad(ngrid,nlayer,nq,pq,reffrad(1,1,iaero_co2))
     494            call aerosol_radius_co2(ngrid,nlayer,nq,pq,reffrad(1,1,iaero_co2))
    488495
    489496            call planetwide_maxval(reffrad(:,:,iaero_co2),maxvalue)
     
    517524         
    518525         if(iaer.eq.iaero_dust)then
    519             call dust_reffrad(ngrid,nlayer,reffrad(1,1,iaero_dust))
     526            call aerosol_radius_dust(ngrid,nlayer,reffrad(1,1,iaero_dust))
    520527            if (is_master) then
    521528               print*,'Dust particle size = ',reffrad(1,1,iaer)/1.e-6,' um'
     
    524531         
    525532         if(iaer.eq.iaero_h2so4)then
    526             call h2so4_reffrad(ngrid,nlayer,reffrad(1,1,iaero_h2so4))
     533            call aerosol_radius_h2so4(ngrid,nlayer,reffrad(1,1,iaero_h2so4))
    527534            if (is_master) then
    528535               print*,'H2SO4 particle size =',reffrad(1,1,iaer)/1.e-6,' um'
     
    531538         
    532539          if(iaer.eq.iaero_back2lay)then
    533             call back2lay_reffrad(ngrid,reffrad(1,1,iaero_back2lay),nlayer,pplev)
    534          endif
    535 
    536          !  For n-layer aerosol size set once for all at firstcall in su_aer_radii
     540            call aerosol_radius_back2lay(ngrid,reffrad(1,1,iaero_back2lay),nlayer,pplev)
     541         endif
     542
     543         !  For n-layer aerosol size set once for all at firstcall in aerosol_radius_init
    537544
    538545!         if(iaer.eq.iaero_aurora)then
     
    549556
    550557      ! Get 3D aerosol optical properties.
    551       call aeroptproperties(ngrid,nlayer,reffrad,nueffrad,         &
     558      call aerosol_optical_properties(ngrid,nlayer,reffrad,nueffrad,         &
    552559           QVISsQREF3d,omegaVIS3d,gVIS3d,                          &
    553560           QIRsQREF3d,omegaIR3d,gIR3d,                             &
     
    555562
    556563      ! Get aerosol optical depths.
    557       call aeropacity(ngrid,nlayer,nq,pplay,pplev,pt,pq,zls,aerosol,      &
     564      call aerosol_opacity(ngrid,nlayer,nq,pplay,pplev,pt,pq,zls,aerosol,      &
    558565           reffrad,nueffrad,QREFvis3d,QREFir3d,                             &
    559566           tau_col,cloudfrac,totcloudfrac,clearsky)               
     
    692699                       (pplev(ig,L_NLAYRAD-k)-pplev(ig,L_NLAYRAD-k+1))
    693700               ! As 'aerosol' is at reference (visible) wavelenght we scale it as
    694                ! it will be multplied by qxi/v in optci/v
     701               ! it will be multplied by qxi/v in rad_correlatedk_opacities_thermal/v
    695702               temp=aerosol(ig,L_NLAYRAD-k,iaer)/QREFvis3d(ig,L_NLAYRAD-k,iaer)
    696703               tauaero(2*k+2,iaer)=max(temp*pweight,0.d0)
     
    918925            qvar(1)=qvar(2)
    919926
    920             write(*,*)trim(subname),' :Warning: reducing qvar in callcorrk.'
     927            write(*,*)trim(subname),' :Warning: reducing qvar in rad_correlatedk.'
    921928            write(*,*)trim(subname),' :Temperature profile no longer consistent ', &
    922929                   'with saturated H2O. qsat=',satval
     
    10181025      ! -- JVO 20 : Also add a sanity test checking that tlevrad is
    10191026      !             within Planck function temperature boundaries,
    1020       !             which would cause gfluxi/sfluxi to crash.
     1027      !             which would cause rad_correlatedk_fluxes_solver_thermal/rad_correlatedk_fluxes_thermal to crash.
    10211028      do k=1,L_LEVELS
    10221029
     
    11271134! Recombine reference corrk tables if needed - Added by JVO, 2020.
    11281135         if (corrk_recombin) then
    1129            call call_recombin(ig,nlayer,pq(ig,:,:),pplay(ig,:),pt(ig,:),qvar(:),tmid(:),pmid(:))
     1136           call rad_correlatedk_recombination_main(ig,nlayer,pq(ig,:,:),pplay(ig,:),pt(ig,:),qvar(:),tmid(:),pmid(:))
    11301137         endif
    11311138! ----------------------------------------------------------------
     
    11491156            endif
    11501157
    1151             call optcv(dtauv,tauv,taucumv,plevrad,                 &
    1152                  qxvaer,qsvaer,gvaer,wbarv,cosbv,tauaero,          &
     1158            call rad_correlatedk_opacities_stellar(dtauv,    &
     1159                 tauv,taucumv,plevrad,                       &
     1160                 qxvaer,qsvaer,gvaer,wbarv,cosbv,tauaero,    &
    11531161                 tmid,pmid,taugsurf,qvar,muvarrad,fracvari)
    11541162
    1155             call sfluxv(dtauv,tauv,taucumv,albv,dwnv,wbarv,cosbv,  &
    1156                  acosz,stel_fract,                                 &
    1157                  nfluxtopv,fluxtopvdn,nfluxoutv_nu,nfluxgndv_nu,   &
     1163            call rad_correlatedk_fluxes_stellar(dtauv,tauv,  &
     1164                 taucumv,albv,dwnv,wbarv,cosbv,              &
     1165                 acosz,stel_fract,                           &
     1166                 nfluxtopv,fluxtopvdn,                       &
     1167                 nfluxoutv_nu,nfluxgndv_nu,                  &
    11581168                 fmnetv,fluxupv,fluxdnv,fzerov,taugsurf)
    11591169
     
    11921202!-----------------------------------------------------------------------
    11931203
    1194          call optci(plevrad,tlevrad,dtaui,taucumi,                  &
    1195               qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid,    &
     1204         call rad_correlatedk_opacities_thermal(plevrad,          &
     1205              tlevrad,dtaui,taucumi,                              &
     1206              qxiaer,qsiaer,giaer,cosbi,wbari,tauaero,tmid,pmid,  &
    11961207              taugsurfi,qvar,muvarrad,fracvari)
    11971208
    1198          call sfluxi(plevrad,tlevrad,dtaui,taucumi,ubari,albi,      &
    1199               wnoi,dwni,cosbi,wbari,nfluxtopi,nfluxtopi_nu,         &
     1209         call rad_correlatedk_fluxes_thermal(plevrad,             &
     1210              tlevrad,dtaui,taucumi,ubari,albi,                   &
     1211              wnoi,dwni,cosbi,wbari,nfluxtopi,nfluxtopi_nu,       &
    12001212              fmneti,fluxupi,fluxdni,fluxupi_nu,fzeroi,taugsurfi)
    12011213
     
    13461358
    13471359
    1348     end subroutine callcorrk
    1349 
    1350 END MODULE callcorrk_mod
     1360    end subroutine rad_correlatedk
     1361
     1362END MODULE rad_correlatedk_mod
Note: See TracChangeset for help on using the changeset viewer.