Changeset 4100 for trunk/LMDZ.VENUS/libf


Ignore:
Timestamp:
Mar 4, 2026, 2:16:03 PM (4 weeks ago)
Author:
emillour
Message:

Venus PCM:
Code tidying: remove common diffusion.h, make these variables module
variables in moldiff_MPF.F90 (turn it into a module in the process).
Turn moldiffcoeff_red.F into a module.
EM

Location:
trunk/LMDZ.VENUS/libf/phyvenus
Files:
1 deleted
3 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.VENUS/libf/phyvenus/moldiff_MPF.F90

    r3956 r4100  
     1MODULE moldiff_MPF_mod
     2
     3  real*8,parameter :: Pdiff=1.e-1 ! pressure (Pa) below which diffusion is computed
     4  real*8,parameter :: tdiffmin=5d0
     5  real*8,parameter :: dzres=0.2d0 ! grid resolution (km) for diffusion
     6
     7CONTAINS
     8
    19subroutine moldiff_MPF(ngrid,nlayer,nq,pplay,pplev,pt,pq,&
    210                       ptimestep,pdteuv,pdtconduc,pdqdiff)
    311
    4 USE chemparam_mod
    5 USE infotrac_phy
    6 USE dimphy   
    7 USE comcstfi_mod
     12use chemparam_mod, only: M_tr
     13use infotrac_phy, only: tname
     14use comcstfi_mod, only: g
     15use moldiffcoeff_red_mod, only: moldiffcoeff_red
    816
    917implicit none
    10 
    11 #include "diffusion.h"
    1218
    1319! June 2023 JYC New method based on the modified pass flow (Parshev et al. 1987)
     
    1925      integer,intent(in) :: nlayer ! number of atmospheric layers
    2026      integer,intent(in) :: nq ! number of advected tracers
    21       real ptimestep
    22       real pplay(ngrid,nlayer)
    23       real pplev(ngrid,nlayer+1)
    24       real pq(ngrid,nlayer,nq)
     27      real,intent(in) :: ptimestep ! physics time step (s)
     28      real,intent(in) :: pplay(ngrid,nlayer) ! mid-layer pressure (Pa)
     29      real,intent(in) :: pplev(ngrid,nlayer+1) ! pressure (Pa) at layer boundaries
     30      real,intent(in) :: pq(ngrid,nlayer,nq) ! input tracer mmr (kg/kg_air)
    2531!      real pdq(ngrid,nlayer,nq)
    26       real pt(ngrid,nlayer)
     32      real,intent(in) :: pt(ngrid,nlayer) ! input temperature (K)
    2733!      real pdt(ngrid,nlayer)
    28       real pdteuv(ngrid,nlayer)
    29       real pdtconduc(ngrid,nlayer)
    30       real pdqdiff(ngrid,nlayer,nq)
     34      real,intent(in) :: pdteuv(ngrid,nlayer) ! tendency on temperature (K/s) dur to EUV
     35      real,intent(in) :: pdtconduc(ngrid,nlayer) ! tendency on temperature (K/s) due to conduction
     36      real,intent(out) :: pdqdiff(ngrid,nlayer,nq) ! tendency on tracers (kg/kg_air/s)
    3137!
    3238! Local
     
    698704
    699705
    700       return
    701       end
     706end subroutine moldiff_MPF
    702707
    703708!    ********************************************************************
     
    804809      enddo
    805810
    806       end
     811      end subroutine tridagbloc
    807812
    808813      subroutine tridag(a,b,c,r,u,n)
     
    82783212    continue
    828833      return
    829       end
     834      end subroutine tridag
    830835
    831836!    ********************************************************************
     
    86687114    CONTINUE
    867872      RETURN
    868       END
     873      END SUBROUTINE LUBKSB
    869874
    870875!    ********************************************************************
     
    952957      ierr =0
    953958      RETURN
    954       END
     959      END SUBROUTINE LUDCMP
    955960
    956961        SUBROUTINE TMNEW(T1,DT1,DT2,DT3,T2,dtime,nl,ig)
     
    972977
    973978        ENDDO
    974         END
     979        END SUBROUTINE TMNEW
    975980
    976981        SUBROUTINE QMNEW(Q1,DQ,Q2,dtime,nl,nq,gc,ig)
     
    992997        ENDDO
    993998        ENDDO
    994         END
     999        END SUBROUTINE QMNEW
    9951000
    9961001        SUBROUTINE HSCALE(p,hp,nl)
     
    10081013        hp(l)=-log(P(l+1)/P(l-1))
    10091014        ENDDO
    1010         END
     1015        END SUBROUTINE HSCALE
    10111016
    10121017        SUBROUTINE MMOY(massemoy,mol_tr,qq,gc,nl,nq)
     
    10261031        enddo
    10271032
    1028         END
     1033        END SUBROUTINE MMOY
    10291034
    10301035        SUBROUTINE DMMOY(M,H,DM,nl)
     
    10401045        enddo
    10411046
    1042         END
     1047        END SUBROUTINE DMMOY
    10431048
    10441049        SUBROUTINE ZVERT(P,T,M,Z,nl,ig)
     
    10731078        enddo
    10741079
    1075         END
     1080        END SUBROUTINE ZVERT
    10761081
    10771082        SUBROUTINE RHOTOT(P,T,M,qq,rhoN,rhoK,nl,nq)
     
    10941099        enddo
    10951100
    1096         END
     1101        END SUBROUTINE RHOTOT
    10971102
    10981103        SUBROUTINE UPPER_RESOL(P,T,Z,M,RR,Rk,                            &
     
    11871192        Praf(nl)=Nraf(nl)*kbolt*Traf(nl)
    11881193        Mraf(nl)=1D0/sum(Qraf(nl,:)/dble(mol_tr(:)))
    1189         END
     1194        END SUBROUTINE UPPER_RESOL
    11901195
    11911196        SUBROUTINE CORRMASS(qq,qint,FacMass,nl,nq)
     
    12001205        enddo
    12011206
    1202         END     
     1207        END SUBROUTINE CORRMASS
    12031208
    12041209
     
    12271232        D(l)=(1d0-Nk(l,nn)/N(l))/interm
    12281233        enddo
    1229         END
     1234        END SUBROUTINE DCOEFF
    12301235 
    12311236        SUBROUTINE HSCALEREAL(nn,Nk,Dz,H,nl,nq)
     
    12561261!       enddo
    12571262
    1258         END
     1263        END SUBROUTINE HSCALEREAL
    12591264       
    12601265        SUBROUTINE VELVERT(nn,T,H,D,Dz,masse,W,nl)
     
    12891294!       enddo
    12901295
    1291         END
     1296        END SUBROUTINE VELVERT
    12921297
    12931298        SUBROUTINE TIMEDIFF(nn,H,W,TIME,nl)
     
    13031308        ENDDO
    13041309
    1305         END
     1310        END SUBROUTINE TIMEDIFF
    13061311
    13071312
     
    13371342        loss(nl)=1D0/dtime       
    13381343
    1339         END
     1344        END SUBROUTINE DIFFPARAM
    13401345
    13411346        SUBROUTINE SEQUENCY(alpha,beta,delta,ksi,eps,zeta,Dad,Kad,rhoad,Loss,Prod,        &
     
    13561361        ENDDO
    13571362
    1358         END
     1363        END SUBROUTINE SEQUENCY
    13591364
    13601365        SUBROUTINE MATCOEFF(alpha,beta,gama,delta,eps,Dad,rhoad,        &
     
    13961401
    13971402
    1398         END
     1403        END SUBROUTINE MATCOEFF
    13991404
    14001405        SUBROUTINE Checkmass(X,Y,nl,nn)
     
    14111416        print*,'no conservation for mass',Xtot,Ytot,nn
    14121417        endif
    1413         END
     1418        END SUBROUTINE Checkmass
    14141419
    14151420        SUBROUTINE Checkmass2(qold,qnew,P,il,nl,nn,nq)
     
    14351440        ENDIF
    14361441
    1437         END
     1442        END SUBROUTINE Checkmass2
    14381443
    14391444        SUBROUTINE GCMGRID_P(Z,P,Q,T,Nk,Rk,qq,qnew,tt,tnew,             &
     
    15391544        enddo
    15401545
    1541         END
     1546        END SUBROUTINE GCMGRID_P
    15421547
    15431548        SUBROUTINE GCMGRID_P2(Z,P,Q,T,Nk,Rk,qq,qnew,tt,tnew             &
     
    16461651!       write(*,*), ' -- Sortie moldiff_red -- '
    16471652
    1648         END
     1653        END SUBROUTINE GCMGRID_P2
    16491654       
     1655END MODULE moldiff_MPF_mod
  • trunk/LMDZ.VENUS/libf/phyvenus/moldiffcoeff_red.F

    r2523 r4100  
     1      MODULE moldiffcoeff_red_mod
     2     
     3      IMPLICIT NONE
     4     
     5      CONTAINS
     6
    17      subroutine moldiffcoeff_red(dij,indic,gcmind,ncompdiff2)
    28
    3        USE chemparam_mod
    4        USE infotrac_phy
    5        USE dimphy   
     9       use chemparam_mod, only: M_tr
     10       use infotrac_phy, only: nqtot, tname
    611
    712       IMPLICIT NONE
     
    1520c
    1621c=======================================================================
    17 #include "diffusion.h"
    1822
    1923c-----------------------------------------------------------------------
     
    304308
    305309
    306       return   
    307       end
     310      end subroutine moldiffcoeff_red
    308311     
     312      END MODULE moldiffcoeff_red_mod
  • trunk/LMDZ.VENUS/libf/phyvenus/physiq_mod.F

    r4092 r4100  
    7979      use compo_hedin83_mod2, only: compo_hedin83_init2
    8080      use compo_hedin83_mod2, only: compo_hedin83_mod
     81      use moldiff_mpf_mod, only: moldiff_mpf
    8182      use nirdata_mod, only: nir_leedat
    8283      use radlwsw_newtoncool_mod, only: radlwsw_newtoncool
     
    239240c
    240241      REAL Fsedim(klon,klev+1)  ! Flux de sedimentation (kg.m-2)
    241       EXTERNAL moldiff_MPF
    242242
    243243c
Note: See TracChangeset for help on using the changeset viewer.