Changeset 2428
- Timestamp:
- Nov 3, 2020, 5:22:38 PM (4 years ago)
- Location:
- trunk/LMDZ.GENERIC
- Files:
-
- 2 edited
- 1 moved
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/README
r2427 r2428 1602 1602 In the process turned vdifc into a module, as well as turbdiff, 1603 1603 for better control, and removed unused arguments. 1604 1605 == 03/11/2020 == EM + YJ 1606 Bug fix on mass_redistribution; argument rnat should be real, not integer. 1607 Turned it into a mass_redistribution_mod module. -
trunk/LMDZ.GENERIC/libf/phystd/mass_redistribution_mod.F90
r2425 r2428 1 MODULE mass_redistribution_mod 2 3 CONTAINS 4 1 5 SUBROUTINE mass_redistribution(ngrid,nlayer,nq,ptimestep, & 2 6 rnat,pcapcal,pplay,pplev,pt,ptsrf,pq,pqs, & … … 5 9 6 10 USE watercommon_h, Only: Tsat_water,RLVTT 7 use surfdat_h8 11 use radcommon_h, only: glat 9 USE tracer_h 12 USE tracer_h, ONLY: igcm_h2o_vap 10 13 USE planete_mod, only: bp 11 14 use comcstfi_mod, only: g … … 55 58 ! Arguments : 56 59 ! --------- 57 INTEGER ngrid, nlayer, nq58 REAL ptimestep59 REAL pcapcal(ngrid)60 INTEGERrnat(ngrid)61 REAL pplay(ngrid,nlayer),pplev(ngrid,nlayer+1)62 REAL pt(ngrid,nlayer),pdt(ngrid,nlayer)63 REAL ptsrf(ngrid),pdtsrf(ngrid)64 REAL pdtmr(ngrid,nlayer)65 REAL pu(ngrid,nlayer) , pv(ngrid,nlayer)66 REAL pdu(ngrid,nlayer) , pdv(ngrid,nlayer)67 REAL pdmassmr(ngrid,nlayer)68 REAL pdumr(ngrid,nlayer) , pdvmr(ngrid,nlayer)69 REAL pq(ngrid,nlayer,nq),pdq(ngrid,nlayer,nq)70 REAL pqs(ngrid,nq)71 REAL pdqmr(ngrid,nlayer,nq),pdqsmr(ngrid,nq)72 REAL pdpsrfmr(ngrid),pdtsrfmr(ngrid)60 INTEGER,INTENT(IN) :: ngrid, nlayer, nq 61 REAL,INTENT(IN) :: ptimestep 62 REAL,INTENT(IN) :: pcapcal(ngrid) 63 REAL,INTENT(IN) :: rnat(ngrid) 64 REAL,INTENT(IN) :: pplay(ngrid,nlayer),pplev(ngrid,nlayer+1) 65 REAL,INTENT(IN) :: pt(ngrid,nlayer),pdt(ngrid,nlayer) 66 REAL,INTENT(IN) :: ptsrf(ngrid),pdtsrf(ngrid) 67 REAL,INTENT(OUT) :: pdtmr(ngrid,nlayer) 68 REAL,INTENT(IN) :: pu(ngrid,nlayer) , pv(ngrid,nlayer) 69 REAL,INTENT(IN) :: pdu(ngrid,nlayer) , pdv(ngrid,nlayer) 70 REAL,INTENT(IN) :: pdmassmr(ngrid,nlayer) 71 REAL,INTENT(OUT) :: pdumr(ngrid,nlayer) , pdvmr(ngrid,nlayer) 72 REAL,INTENT(IN) :: pq(ngrid,nlayer,nq),pdq(ngrid,nlayer,nq) 73 REAL,INTENT(IN) :: pqs(ngrid,nq) 74 REAL,INTENT(OUT) :: pdqmr(ngrid,nlayer,nq),pdqsmr(ngrid,nq) 75 REAL,INTENT(OUT) :: pdpsrfmr(ngrid),pdtsrfmr(ngrid) 73 76 ! 74 77 ! Local variables : … … 143 146 if (ztsrf(ig).gt.Tsat(ig)) then 144 147 zmassboil(ig)=(ptsrf(ig)-Tsat(ig))*pcapcal(ig)/RLVTT/ptimestep 145 if ((zmassboil(ig)*ptimestep.gt.pqs(ig,igcm_h2o_vap)).and.( rnat(ig).eq.1)) then148 if ((zmassboil(ig)*ptimestep.gt.pqs(ig,igcm_h2o_vap)).and.(nint(rnat(ig)).eq.1)) then 146 149 zmassboil(ig)=pqs(ig,igcm_h2o_vap)/ptimestep 147 150 endif … … 407 410 408 411 END SUBROUTINE mass_redistribution 412 413 END MODULE mass_redistribution_mod -
trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90
r2427 r2428 55 55 use turbdiff_mod, only: turbdiff 56 56 use turb_mod, only : q2,sensibFlux,turb_resolved 57 use mass_redistribution_mod, only: mass_redistribution 57 58 #ifndef MESOSCALE 58 59 use vertical_layers_mod, only: presnivs, pseudoalt
Note: See TracChangeset
for help on using the changeset viewer.