- Timestamp:
- May 6, 2024, 12:03:23 PM (7 months ago)
- Location:
- trunk/LMDZ.GENERIC
- Files:
-
- 3 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.GENERIC/changelog.txt
r3312 r3322 1921 1921 No correction factor with diurnal equal true for photolysis rate 1922 1922 + writediagspecUV change name to flux_surf for clarity 1923 1924 == 06/05/2024 == EM 1925 Cleanup and cosmetics on convective adjustment routine. Make it a 1926 module in the process. -
trunk/LMDZ.GENERIC/libf/phystd/convadj.F
r3236 r3322 1 module convadj_mod 2 3 implicit none 4 5 contains 6 1 7 subroutine convadj(ngrid,nlay,nq,ptimestep, 2 8 & pplay,pplev,ppopsk, … … 5 11 & pduadj,pdvadj,pdhadj,pdqadj) 6 12 7 USE tracer_h13 use tracer_h, only: igcm_h2o_vap 8 14 use comcstfi_mod, only: g 9 use generic_cloud_common_h 15 use generic_cloud_common_h, only: epsi_generic 10 16 use generic_tracer_index_mod, only: generic_tracer_index 11 17 use callkeys_mod, only: tracer,water,generic_condensation, … … 18 24 ! Purpose 19 25 ! ------- 20 ! Calculates dry convective adjustment. If one tracer is CO2, 21 ! we take into account the molecular mass variation (e.g. when 22 ! CO2 condenses) to trigger convection (F. Forget 01/2005) 23 ! 24 ! Authors 25 ! ------- 26 ! Original author unknown. 27 ! Modif. 2005 by F. Forget. 28 ! 26 ! Compute dry convective adjustment. 27 ! See old reference paper: Hourdin et al. JAS 1993 28 ! "Meteorological Variability and the Annual Surface 29 ! Pressure Cycle on Mars" 30 ! https://doi.org/10.1175/1520-0469(1993)050%3C3625:MVATAS%3E2.0.CO;2 29 31 !================================================================== 30 31 ! ------------32 ! Declarations33 ! ------------34 35 32 36 33 ! Arguments 37 34 ! --------- 38 35 39 INTEGER ngrid,nlay 40 REAL ptimestep 41 REAL ph(ngrid,nlay),pdhfi(ngrid,nlay),pdhadj(ngrid,nlay) 42 REAL pplay(ngrid,nlay),pplev(ngrid,nlay+1),ppopsk(ngrid,nlay) 43 REAL pu(ngrid,nlay),pdufi(ngrid,nlay),pduadj(ngrid,nlay) 44 REAL pv(ngrid,nlay),pdvfi(ngrid,nlay),pdvadj(ngrid,nlay) 45 46 ! Tracers 47 integer nq 48 real pq(ngrid,nlay,nq), pdqfi(ngrid,nlay,nq) 49 real pdqadj(ngrid,nlay,nq) 36 INTEGER,INTENT(IN) :: ngrid ! number of atmospheric columns 37 INTEGER,INTENT(IN) :: nlay ! number of atmospheric layers 38 INTEGER,INTENT(IN) :: nq ! number of tracers 39 REAL,INTENT(IN) :: ptimestep ! physics time step (s) 40 REAL,INTENT(IN) :: pplay(ngrid,nlay) ! mi-layer pressure (Pa) 41 REAL,INTENT(IN) :: pplev(ngrid,nlay+1) ! inter-layer pressure (Pa) 42 REAL,INTENT(IN) :: ppopsk(ngrid,nlay) ! Exner (wrt surface pressure) 43 ! fields from dynamics 44 REAL,INTENT(IN) :: ph(ngrid,nlay) ! potential temperature (K) 45 REAL,INTENT(IN) :: pu(ngrid,nlay) ! zonal wind (m/s) 46 REAL,INTENT(IN) :: pv(ngrid,nlay) ! meridioanl wind (m/s) 47 REAL,INTENT(IN) :: pq(ngrid,nlay,nq) ! tracers (../kg_air) 48 ! tendencies due to previous physical processes 49 REAL,INTENT(IN) :: pdufi(ngrid,nlay) ! on zonal wind (m/s/s) 50 REAL,INTENT(IN) :: pdvfi(ngrid,nlay) ! on meridional wind (m/s/s) 51 REAL,INTENT(IN) :: pdhfi(ngrid,nlay)! on potential temperature (/K/s) 52 REAL,INTENT(IN) :: pdqfi(ngrid,nlay,nq) ! on tracers (../kg_air/s) 53 ! tendencies due to convetive adjustement 54 REAL,INTENT(OUT) :: pduadj(ngrid,nlay) ! on zonal wind (m/s/s) 55 REAL,INTENT(OUT) :: pdvadj(ngrid,nlay) ! on meridinal wind (m/s/s) 56 REAL,INTENT(OUT) :: pdhadj(ngrid,nlay) ! on potential temperature (/K/s) 57 REAL,INTENT(OUT) :: pdqadj(ngrid,nlay,nq) ! on traceurs (../kg_air/s) 50 58 51 59 … … 66 74 INTEGER iq 67 75 REAL zq(ngrid,nlay,nq), zq2(ngrid,nlay,nq) 68 REAL zqm(nq) ,zqco2m76 REAL zqm(nq) 69 77 70 78 integer igcm_generic_vap, igcm_generic_ice … … 75 83 ! for conservation test 76 84 real masse,cadjncons 77 78 EXTERNAL SCOPY79 85 80 86 ! -------------- … … 87 93 if(tracer) then 88 94 zq(:,:,:)=pq(:,:,:)+pdqfi(:,:,:)*ptimestep 89 !zq(:,:,:)=090 95 end if 91 96 92 CALL scopy(ngrid*nlay,zh,1,zh2,1) 93 CALL scopy(ngrid*nlay,zu,1,zu2,1) 94 CALL scopy(ngrid*nlay,zv,1,zv2,1) 95 CALL scopy(ngrid*nlay*nq,zq,1,zq2,1) 96 97 zh2(:,:)=zh(:,:) 98 zu2(:,:)=zu(:,:) 99 zv2(:,:)=zv(:,:) 100 zq2(:,:,:)=zq(:,:,:) 97 101 98 102 ! ----------------------------- … … 115 119 endif 116 120 ENDDO 117 CALL scopy(ngrid*nlay,zvh,1,zvh2,1)118 CALL scopy(ngrid*nlay,zvh2,1,zhc,1)121 zvh2(:,:)=zvh(:,:) 122 zhc(:,:)=zvh2(:,:) 119 123 else 120 CALL scopy(ngrid*nlay,zh2,1,zhc,1)124 zhc(:,:)=zh2(:,:) 121 125 endif 122 126 … … 285 289 cadjncons=0.0 286 290 if(water)then 287 do l = 1, nlay291 do l = 1, nlay 288 292 masse = (pplev(i,l) - pplev(i,l+1))/g 289 293 iq = igcm_h2o_vap 290 294 cadjncons = cadjncons + 291 295 & masse*(zq2(i,l,iq)-zq(i,l,iq))/ptimestep 292 end do296 end do 293 297 endif 294 298 … … 353 357 end if 354 358 355 356 ! output 357 ! if (ngrid.eq.1) then 358 ! ig=1 359 ! iq =1 360 ! write(*,*)'**** l, pq(ig,l,iq),zq(ig,l,iq),zq2(ig,l,iq)' 361 ! do l=nlay,1,-1 362 ! write(*,*) l, pq(ig,l,iq),zq(ig,l,iq),zq2(ig,l,iq) 363 ! end do 364 ! end if 365 366 367 return 368 end 359 end subroutine convadj 360 361 end module convadj_mod -
trunk/LMDZ.GENERIC/libf/phystd/physiq_mod.F90
r3299 r3322 80 80 use volcano_mod, only: volcano 81 81 use pindex_mod, only: pindex 82 use convadj_mod, only: convadj 82 83 use vdifc_mod, only: vdifc 83 84 use turbdiff_mod, only: turbdiff
Note: See TracChangeset
for help on using the changeset viewer.