module physiq_mod implicit none contains subroutine physiq(ngrid,nlayer,nq, & firstcall,lastcall, & pday,ptime,ptimestep, & pplev,pplay,pphi, & pu,pv,pt,pq, & flxw, & pdu,pdv,pdt,pdq,pdpsrf) !! use write_field_phy, only: Writefield_phy !! use ioipsl_getin_p_mod, only: getin_p use radinc_h, only : L_NSPECTI,L_NSPECTV,naerkind, corrkdir, banddir use generic_cloud_common_h, only : epsi_generic, Psat_generic use thermcell_mod, only: init_thermcell_mod use gases_h, only: gnom, gfrac use radcommon_h, only: sigma, glat, grav, BWNV, WNOI, DWNI, DWNV, WNOV use suaer_corrk_mod, only: suaer_corrk use radii_mod, only: n2_reffrad use aerosol_mod, only: iniaerosol, iaero_n2 use surfdat_h, only: phisfi, zmea, zstd, zsig, zgam, zthe, & dryness use comdiurn_h, only: coslat, sinlat, coslon, sinlon use comsaison_h, only: mu0, fract, dist_star, declin, right_ascen use comsoil_h, only: nsoilmx, layer, mlayer, inertiedat use geometry_mod, only: latitude, longitude, cell_area USE comgeomfi_h, only: totarea, totarea_planet USE tracer_h, only: noms, mmol, radius, rho_q, qext, & alpha_lift, alpha_devil, qextrhor, & igcm_dustbin, & igcm_n2_ice, nesp, is_chim, is_condensable,constants_epsi_generic use time_phylmdz_mod, only: ecritphy, iphysiq, nday use phyetat0_mod, only: phyetat0 use wstats_mod, only: callstats, wstats, mkstats use phyredem, only: physdem0, physdem1 use planetwide_mod, only: planetwide_minval,planetwide_maxval,planetwide_sumval use mod_phys_lmdz_para, only : is_master use planete_mod, only: apoastr, periastr, year_day, peri_day, & obliquit, nres, z0 use comcstfi_mod, only: pi, g, rcp, r, rad, mugaz, cpp use time_phylmdz_mod, only: daysec use callkeys_mod, only: albedo_spectral_mode, calladj, calldifv, & callrad, callsoil, nosurf, & n2cond, corrk, diagdtau, & diurnal, enertest, fat1au, & icetstep, intheat, iradia, kastprof, & lwrite, mass_redistrib, meanOLR, & nearn2cond, noseason_day, & season, sedimentation,generic_condensation, & specOLR, & startphy_file, testradtimes, & tracer, UseTurbDiff, & global1d, szangle use generic_tracer_index_mod, only: generic_tracer_index use check_fields_mod, only: check_physics_fields use conc_mod, only: rnew, cpnew, ini_conc_mod use phys_state_var_mod use callcorrk_mod, only: callcorrk use vdifc_mod, only: vdifc use turbdiff_mod, only: turbdiff use turb_mod, only : q2,sensibFlux,turb_resolved use mass_redistribution_mod, only: mass_redistribution use condensation_generic_mod, only: condensation_generic use datafile_mod, only: datadir #ifndef MESOSCALE use vertical_layers_mod, only: presnivs, pseudoalt use mod_phys_lmdz_omp_data, ONLY: is_omp_master #else use comm_wrf, only : comm_HR_SW, comm_HR_LW, & comm_CLOUDFRAC,comm_TOTCLOUDFRAC, comm_RH, & comm_HR_DYN, & comm_DQICE,comm_DQVAP,comm_ALBEQ, & comm_FLUXTOP_DN,comm_FLUXABS_SW, & comm_FLUXTOP_LW,comm_FLUXSURF_SW, & comm_FLUXSURF_LW,comm_FLXGRD, & comm_DTRAIN,comm_DTLSC, & comm_LATENT_HF #endif #ifdef CPP_XIOS use xios_output_mod, only: initialize_xios_output, & update_xios_timestep, & send_xios_field use wxios, only: wxios_context_init, xios_context_finalize #endif implicit none !================================================================== ! ! Purpose ! ------- ! Central subroutine for all the physics parameterisations in the ! universal model. Originally adapted from the Mars LMDZ model. ! ! The model can be run without or with tracer transport ! depending on the value of "tracer" in file "callphys.def". ! ! ! It includes: ! ! I. Initialization : ! I.1 Firstcall initializations. ! I.2 Initialization for every call to physiq. ! ! II. Compute radiative transfer tendencies (longwave and shortwave) : ! II.a Option 1 : Call correlated-k radiative transfer scheme. ! II.b Option 2 : Call Newtonian cooling scheme. ! II.c Option 3 : Atmosphere has no radiative effect. ! ! III. Vertical diffusion (turbulent mixing) : ! ! IV. Convection : ! IV.a Thermal plume model ! IV.b Dry convective adjusment ! ! V. Condensation and sublimation of gases (currently just N2). ! ! VI. Tracers ! VI.1. Water and water ice. ! VI.2 Photochemistry ! VI.3. Aerosols and particles. ! VI.4. Updates (pressure variations, surface budget). ! VI.5. Slab Ocean. ! VI.6. Surface Tracer Update. ! ! VII. Surface and sub-surface soil temperature. ! ! VIII. Perform diagnostics and write output files. ! ! ! arguments ! --------- ! ! INPUT ! ----- ! ! ngrid Size of the horizontal grid. ! nlayer Number of vertical layers. ! nq Number of advected fields. ! ! firstcall True at the first call. ! lastcall True at the last call. ! ! pday Number of days counted from the North. Spring equinoxe. ! ptime Universal time (00 WHEN DOWNWARDS !!) real omega(ngrid,nlayer) ! omega velocity (Pa/s, >0 when downward) integer l,ig,ierr,iq,nw,isoil,iesp, igcm_generic_vap, igcm_generic_ice logical call_ice_vap_generic ! to call only one time the ice/vap pair of a tracer real zls ! Solar longitude (radians). real zlss ! Sub solar point longitude (radians). real zday ! Date (time since Ls=0, calculated in sols). real zzlay(ngrid,nlayer) ! Altitude at the middle of the atmospheric layers. real zzlev(ngrid,nlayer+1) ! Altitude at the atmospheric layer boundaries. ! VARIABLES for the thermal plume model (AF24: deleted) ! TENDENCIES due to various processes : ! For Surface Temperature : (K/s) real zdtsurf(ngrid) ! Cumulated tendencies. real zdtsurfmr(ngrid) ! Mass_redistribution routine. real zdtsurfc(ngrid) ! Condense_n2 routine. real zdtsdif(ngrid) ! Turbdiff/vdifc routines. ! real zdtsurf_hyd(ngrid) ! Hydrol routine. ! For Atmospheric Temperatures : (K/s) real dtlscale(ngrid,nlayer) ! Largescale routine. real dt_generic_condensation(ngrid,nlayer) ! condensation_generic routine. real zdtc(ngrid,nlayer) ! Condense_n2 routine. real zdtdif(ngrid,nlayer) ! Turbdiff/vdifc routines. real zdtmr(ngrid,nlayer) ! Mass_redistribution routine. real zdtsw1(ngrid,nlayer), zdtlw1(ngrid,nlayer) ! Callcorrk routine. real zdtchim(ngrid,nlayer) ! Calchim routine. ! For Surface Tracers : (kg/m2/s) real dqsurf(ngrid,nq) ! Cumulated tendencies. real zdqsurfc(ngrid) ! Condense_n2 routine. real zdqsdif(ngrid,nq) ! Turbdiff/vdifc routines. real zdqssed(ngrid,nq) ! Callsedim routine. real zdqsurfmr(ngrid,nq) ! Mass_redistribution routine. ! For Tracers : (kg/kg_of_air/s) real zdqc(ngrid,nlayer,nq) ! Condense_n2 routine. real zdqadj(ngrid,nlayer,nq) ! Convadj routine. real zdqdif(ngrid,nlayer,nq) ! Turbdiff/vdifc routines. real zdqevap(ngrid,nlayer) ! Turbdiff routine. real zdqsed(ngrid,nlayer,nq) ! Callsedim routine. real zdqmr(ngrid,nlayer,nq) ! Mass_redistribution routine. real dqvaplscale_generic(ngrid,nlayer,nq) ! condensation_generic routine. real dqcldlscale_generic(ngrid,nlayer,nq) ! condensation_generic routine. REAL,allocatable,save :: zdqchim(:,:,:) ! Calchim_asis routine REAL,allocatable,save :: zdqschim(:,:) ! Calchim_asis routine !$OMP THREADPRIVATE(zdqchim,zdqschim) REAL array_zero1(ngrid) REAL array_zero2(ngrid,nlayer) ! For Winds : (m/s/s) real zdvadj(ngrid,nlayer), zduadj(ngrid,nlayer) ! Convadj routine. real zdumr(ngrid,nlayer), zdvmr(ngrid,nlayer) ! Mass_redistribution routine. real zdvdif(ngrid,nlayer), zdudif(ngrid,nlayer) ! Turbdiff/vdifc routines. real zdhdif(ngrid,nlayer) ! Turbdiff/vdifc routines. real zdhadj(ngrid,nlayer) ! Convadj routine. ! For Pressure and Mass : real zdmassmr(ngrid,nlayer) ! Atmospheric Mass tendency for mass_redistribution (kg_of_air/m2/s). real zdmassmr_col(ngrid) ! Atmospheric Column Mass tendency for mass_redistribution (kg_of_air/m2/s). real zdpsrfmr(ngrid) ! Pressure tendency for mass_redistribution routine (Pa/s). ! Local variables for LOCAL CALCULATIONS: ! --------------------------------------- real zflubid(ngrid) real zplanck(ngrid),zpopsk(ngrid,nlayer) real ztim1,ztim2,ztim3, z1,z2 real ztime_fin real zdh(ngrid,nlayer) real gmplanet real taux(ngrid),tauy(ngrid) ! local variables for DIAGNOSTICS : (diagfi & stat) ! ------------------------------------------------- real ps(ngrid) ! Surface Pressure. real zt(ngrid,nlayer) ! Atmospheric Temperature. real zu(ngrid,nlayer),zv(ngrid,nlayer) ! Zonal and Meridional Winds. real zq(ngrid,nlayer,nq) ! Atmospheric Tracers. real zdtadj(ngrid,nlayer) ! Convadj Diagnostic. real zdtdyn(ngrid,nlayer) ! Dynamical Heating (K/s). real zdudyn(ngrid,nlayer) ! Dynamical Zonal Wind tendency (m.s-2). real reff(ngrid,nlayer) ! Effective dust radius (used if doubleq=T). real vmr(ngrid,nlayer) ! volume mixing ratio real time_phys real ISR,ASR,OLR,GND,DYN,GSR,Ts1,Ts2,Ts3,TsS ! for Diagnostic. real qcol(ngrid,nq) ! Tracer Column Mass (kg/m2). ! included by RW for H2O Manabe scheme real rneb_man(ngrid,nlayer) ! H2O cloud fraction (moistadj). real rneb_lsc(ngrid,nlayer) ! H2O cloud fraction (large scale). ! to test energy conservation (RW+JL) real mass(ngrid,nlayer),massarea(ngrid,nlayer) real dEtot, dEtots, AtmToSurf_TurbFlux real,save :: dEtotSW, dEtotsSW, dEtotLW, dEtotsLW !$OMP THREADPRIVATE(dEtotSW, dEtotsSW, dEtotLW, dEtotsLW) !JL12 conservation test for mean flow kinetic energy has been disabled temporarily real dtmoist_max,dtmoist_min real dItot, dItot_tmp, dVtot, dVtot_tmp real dWtot, dWtot_tmp, dWtots, dWtots_tmp real psat_tmp ! AF24: to remove? real qsat_generic(ngrid,nlayer,nq) ! generic condensable tracers (GCS) specific concentration at saturation (kg/kg_of_air). real RH_generic(ngrid,nlayer,nq) ! generic condensable tracers (GCS) Relative humidity. real rneb_generic(ngrid,nlayer,nq) ! GCS cloud fraction (generic condensation). real psat_tmp_generic real, save :: metallicity ! metallicity of planet --- is not used here, but necessary to call function Psat_generic !$OMP THREADPRIVATE(metallicity) real reffrad_generic_zeros_for_wrf(ngrid,nlayer) ! !!! this is temporary, it is only a list of zeros, it will be replaced when a generic aerosol will be implemented ! For Clear Sky Case. (AF24: deleted) real nconsMAX, vdifcncons(ngrid), cadjncons(ngrid) ! Vdfic water conservation test. By RW real muvar(ngrid,nlayer+1) ! For Runaway Greenhouse 1D study. By RW real,save,allocatable :: reffcol(:,:) !$OMP THREADPRIVATE(reffcol) ! Non-oro GW tendencies REAL d_u_hin(ngrid,nlayer), d_v_hin(ngrid,nlayer) REAL d_t_hin(ngrid,nlayer) ! Diagnostics 2D of gw_nonoro REAL zustrhi(ngrid), zvstrhi(ngrid) real :: tsurf2(ngrid) !! real :: flux_o(ngrid),flux_g(ngrid) real :: flux_g(ngrid) real :: flux_sens_lat(ngrid) real :: qsurfint(ngrid,nq) #ifdef MESOSCALE REAL :: lsf_dt(nlayer) REAL :: lsf_dq(nlayer) #endif ! flags to trigger extra sanity checks logical, save :: check_physics_inputs=.false. logical, save :: check_physics_outputs=.false. !$OPM THREADPRIVATE(check_physics_inputs,check_physics_outputs) ! Misc character*2 :: str2 character(len=10) :: tmp1 character(len=10) :: tmp2 !================================================================================================== ! ----------------- ! I. INITIALISATION ! ----------------- ! -------------------------------- ! I.1 First Call Initialisation. ! -------------------------------- if (firstcall) then call getin_p("check_physics_inputs", check_physics_inputs) call getin_p("check_physics_outputs", check_physics_outputs) ! Allocate saved arrays (except for 1D model, where this has already ! been done) #ifndef MESOSCALE if (ngrid>1) call phys_state_var_init(nq) #endif ! Variables set to 0 ! ~~~~~~~~~~~~~~~~~~ dtrad(:,:) = 0.0 fluxrad(:) = 0.0 tau_col(:) = 0.0 zdtsw(:,:) = 0.0 zdtlw(:,:) = 0.0 ! Initialize tracer names, indexes and properties. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ IF (.NOT.ALLOCATED(noms)) ALLOCATE(noms(nq)) ! (because noms is an argument of physdem1 whether or not tracer is on) if (tracer) then call initracer(ngrid,nq) ! if(photochem) then !AF24: removed endif ! Initialize aerosol indexes. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~ call iniaerosol ! allocate related local arrays ! (need be allocated instead of automatic because of "naerkind") allocate(aerosol(ngrid,nlayer,naerkind)) allocate(reffcol(ngrid,naerkind)) #ifdef CPP_XIOS ! Initialize XIOS context write(*,*) "physiq: call wxios_context_init" CALL wxios_context_init #endif ! Read 'startfi.nc' file. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ #ifndef MESOSCALE call phyetat0(startphy_file, & ngrid,nlayer,"startfi.nc",0,0,nsoilmx,nq, & day_ini,time_phys,tsurf,tsoil,emis,q2,qsurf) #else day_ini = pday #endif #ifndef MESOSCALE if (.not.startphy_file) then ! additionnal "academic" initialization of physics if (is_master) write(*,*) "Physiq: initializing tsurf(:) to pt(:,1) !!" tsurf(:)=pt(:,1) if (is_master) write(*,*) "Physiq: initializing tsoil(:) to pt(:,1) !!" do isoil=1,nsoilmx tsoil(1:ngrid,isoil)=tsurf(1:ngrid) enddo if (is_master) write(*,*) "Physiq: initializing day_ini to pdat !" day_ini=pday endif #endif if (pday.ne.day_ini) then write(*,*) "ERROR in physiq.F90:" write(*,*) "bad synchronization between physics and dynamics" write(*,*) "dynamics day: ",pday write(*,*) "physics day: ",day_ini stop endif write (*,*) 'In physiq day_ini =', day_ini ! Initialize albedo calculation. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ albedo(:,:)=0.0 albedo_bareground(:)=0.0 albedo_snow_SPECTV(:)=0.0 albedo_n2_ice_SPECTV(:)=0.0 call surfini(ngrid,nq,qsurf,albedo,albedo_bareground,albedo_snow_SPECTV,albedo_n2_ice_SPECTV) ! Initialize orbital calculation. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ call iniorbit(apoastr,periastr,year_day,peri_day,obliquit) ! Initialize soil. ! ~~~~~~~~~~~~~~~~ if (callsoil) then call soil(ngrid,nsoilmx,firstcall,lastcall,inertiedat, & ptimestep,tsurf,tsoil,capcal,fluxgrd) else ! else of 'callsoil'. print*,'WARNING! Thermal conduction in the soil turned off' capcal(:)=1.e6 fluxgrd(:)=intheat print*,'Flux from ground = ',intheat,' W m^-2' endif ! end of 'callsoil'. icount=1 ! Initialize surface history variable. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ qsurf_hist(:,:)=qsurf(:,:) !! call WriteField_phy("post_qsurf_hist_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1) ! Initialize variable for dynamical heating and zonal wind tendency diagnostic ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ ztprevious(:,:)=pt(:,:) zuprevious(:,:)=pu(:,:) ! Set temperature just above condensation temperature (for Early Mars) ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ if(nearn2cond) then write(*,*)' WARNING! Starting at Tcond+1K' do l=1, nlayer do ig=1,ngrid pdt(ig,l)= ((-3167.8)/(log(.01*pplay(ig,l))-23.23)+4 & -pt(ig,l)) / ptimestep enddo enddo endif if(meanOLR)then call system('rm -f rad_bal.out') ! to record global radiative balance. call system('rm -f tem_bal.out') ! to record global mean/max/min temperatures. call system('rm -f h2o_bal.out') ! to record global hydrological balance. endif !! Initialize variables for water cycle ! AF24: removed ! Set metallicity for GCS ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ metallicity=0.0 ! default value --- is not used here but necessary to call function Psat_generic call getin_p("metallicity",metallicity) ! --- is not used here but necessary to call function Psat_generic ! Set some parameters for the thermal plume model !AF24: removed ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ #ifndef MESOSCALE if (ngrid.ne.1) then ! Note : no need to create a restart file in 1d. call physdem0("restartfi.nc",longitude,latitude,nsoilmx,ngrid,nlayer,nq, & ptimestep,pday+nday,time_phys,cell_area, & albedo_bareground,inertiedat,zmea,zstd,zsig,zgam,zthe) endif !! call WriteField_phy("post_physdem_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1) #endif if (corrk) then ! We initialise the spectral grid here instead of ! at firstcall of callcorrk so we can output XspecIR, XspecVI ! when using Dynamico print*, "physiq_mod: Correlated-k data base folder:",trim(datadir) call getin_p("corrkdir",corrkdir) print*,"corrkdir = ", corrkdir write (tmp1, '(i4)') L_NSPECTI write (tmp2, '(i4)') L_NSPECTV banddir=trim(trim(adjustl(tmp1))//'x'//trim(adjustl(tmp2))) banddir=trim(trim(adjustl(corrkdir))//'/'//trim(adjustl(banddir))) call setspi !Basic infrared properties. call setspv ! Basic visible properties. call sugas_corrk ! Set up gaseous absorption properties. call suaer_corrk ! Set up aerosol optical properties. endif !! call WriteField_phy("post_corrk_firstcall_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1) ! XIOS outputs #ifdef CPP_XIOS write(*,*) "physiq: call initialize_xios_output" call initialize_xios_output(pday,ptime,ptimestep,daysec, & year_day,presnivs,pseudoalt,WNOI,WNOV) #endif !! call WriteField_phy("post_xios_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1) write(*,*) "physiq: end of firstcall" endif ! end of 'firstcall' !! call WriteField_phy("post_firstcall_qsurf",qsurf(1:ngrid,igcm_h2o_vap),1) !! call writediagfi(ngrid,"firstcall_post_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap)) if (check_physics_inputs) then !check the validity of input fields coming from the dynamics call check_physics_fields("begin physiq:", pt, pu, pv, pplev, pq) endif ! call writediagfi(ngrid,"pre_physical_rnat"," "," ",2,rnat) ! call writediagfi(ngrid,"pre_physical_capcal"," "," ",2,capcal) ! ------------------------------------------------------ ! I.2 Initializations done at every physical timestep: ! ------------------------------------------------------ #ifdef CPP_XIOS ! update XIOS time/calendar call update_xios_timestep #endif ! Initialize various variables ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~ if ( .not.nearn2cond ) then pdt(1:ngrid,1:nlayer) = 0.0 endif zdtsurf(1:ngrid) = 0.0 pdq(1:ngrid,1:nlayer,1:nq) = 0.0 dqsurf(1:ngrid,1:nq)= 0.0 pdu(1:ngrid,1:nlayer) = 0.0 pdv(1:ngrid,1:nlayer) = 0.0 pdpsrf(1:ngrid) = 0.0 zflubid(1:ngrid) = 0.0 flux_sens_lat(1:ngrid) = 0.0 taux(1:ngrid) = 0.0 tauy(1:ngrid) = 0.0 zday=pday+ptime ! Compute time, in sols (and fraction thereof). ! Compute Stellar Longitude (Ls), and orbital parameters. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ if (season) then call stellarlong(zday,zls) else call stellarlong(noseason_day,zls) end if call orbite(zls,dist_star,declin,right_ascen) if (diurnal) then zlss=-2.*pi*(zday-.5) else if(diurnal .eqv. .false.) then zlss=9999. endif glat(:) = g !AF24: removed oblateness ! Compute geopotential between layers. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ zzlay(1:ngrid,1:nlayer)=pphi(1:ngrid,1:nlayer) do l=1,nlayer zzlay(1:ngrid,l)= zzlay(1:ngrid,l)/glat(1:ngrid) enddo zzlev(1:ngrid,1)=0. do l=2,nlayer do ig=1,ngrid z1=(pplay(ig,l-1)+pplev(ig,l))/(pplay(ig,l-1)-pplev(ig,l)) z2=(pplev(ig,l)+pplay(ig,l))/(pplev(ig,l)-pplay(ig,l)) zzlev(ig,l)=(z1*zzlay(ig,l-1)+z2*zzlay(ig,l))/(z1+z2) enddo enddo !Altitude of top interface (nlayer+1), using the thicknesss of the level below the top one. LT22 zzlev(1:ngrid,nlayer+1) = 2*zzlev(1:ngrid,nlayer)-zzlev(1:ngrid,nlayer-1) ! Compute potential temperature ! Note : Potential temperature calculation may not be the same in physiq and dynamic... ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ do l=1,nlayer do ig=1,ngrid zpopsk(ig,l)=(pplay(ig,l)/pplev(ig,1))**rcp zh(ig,l)=pt(ig,l)/zpopsk(ig,l) mass(ig,l) = (pplev(ig,l) - pplev(ig,l+1))/glat(ig) massarea(ig,l)=mass(ig,l)*cell_area(ig) enddo enddo ! Compute vertical velocity (m/s) from vertical mass flux ! w = F / (rho*area) and rho = P/(r*T) ! But first linearly interpolate mass flux to mid-layers do l=1,nlayer-1 pw(1:ngrid,l)=0.5*(flxw(1:ngrid,l)+flxw(1:ngrid,l+1)) enddo pw(1:ngrid,nlayer)=0.5*flxw(1:ngrid,nlayer) ! since flxw(nlayer+1)=0 do l=1,nlayer pw(1:ngrid,l)=(pw(1:ngrid,l)*r*pt(1:ngrid,l)) / & (pplay(1:ngrid,l)*cell_area(1:ngrid)) enddo ! omega in Pa/s do l=1,nlayer-1 omega(1:ngrid,l)=0.5*(flxw(1:ngrid,l)+flxw(1:ngrid,l+1)) enddo omega(1:ngrid,nlayer)=0.5*flxw(1:ngrid,nlayer) ! since flxw(nlayer+1)=0 do l=1,nlayer omega(1:ngrid,l)=g*omega(1:ngrid,l)/cell_area(1:ngrid) enddo !--------------------------------- ! II. Compute radiative tendencies !--------------------------------- ! Compute local stellar zenith angles ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ if (diurnal) then ztim1=SIN(declin) ztim2=COS(declin)*COS(2.*pi*(zday-.5)) ztim3=-COS(declin)*SIN(2.*pi*(zday-.5)) call stelang(ngrid,sinlon,coslon,sinlat,coslat, & ztim1,ztim2,ztim3,mu0,fract) else if(diurnal .eqv. .false.) then call mucorr(ngrid,declin,latitude,mu0,fract,10000.,rad) ! WARNING: this function appears not to work in 1D if ((ngrid.eq.1).and.(global1d)) then ! Fixed zenith angle 'szangle' in 1D simulations w/ globally-averaged sunlight. mu0 = cos(pi*szangle/180.0) endif endif if (callrad) then if( mod(icount-1,iradia).eq.0.or.lastcall) then ! Eclipse incoming sunlight !AF24: removed !! call writediagfi(ngrid,"corrk_pre_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap)) !! call writediagfi(ngrid,"corrk_pre_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap)) if (corrk) then ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ ! II.a Call correlated-k radiative transfer scheme ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ if(kastprof)then print*,'kastprof should not = true here' call abort endif ! if(water) then !AF24: removed if(generic_condensation) then do iq=1,nq call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic) if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer epsi_generic=constants_epsi_generic(iq) muvar(1:ngrid,1:nlayer)=mugaz/(1.e0+(1.e0/epsi_generic-1.e0)*pq(1:ngrid,1:nlayer,igcm_generic_vap)) muvar(1:ngrid,nlayer+1)=mugaz/(1.e0+(1.e0/epsi_generic-1.e0)*pq(1:ngrid,nlayer,igcm_generic_vap)) endif end do ! do iq=1,nq loop on tracers ! take into account generic condensable specie (GCS) effect on mean molecular weight else muvar(1:ngrid,1:nlayer+1)=mugaz endif ! if(ok_slab_ocean) then !AF24: removed ! standard callcorrk ! clearsky=.false. call callcorrk(ngrid,nlayer,pq,nq,qsurf, & albedo,albedo_equivalent,emis,mu0,pplev,pplay,pt, & tsurf,fract,dist_star,aerosol,muvar, & zdtlw,zdtsw,fluxsurf_lw,fluxsurf_sw, & fluxsurfabs_sw,fluxtop_lw, & fluxabs_sw,fluxtop_dn,OLR_nu,OSR_nu,GSR_nu, & int_dtaui,int_dtauv, & tau_col,cloudfrac,totcloudfrac, & .false.,firstcall,lastcall) !GG (feb2021): Option to "artificially" decrease the raditive time scale in !the deep atmosphere press > 0.1 bar. Suggested by MT. !! COEFF_RAD to be "tuned" to facilitate convergence of tendency !coeff_rad=0. ! 0 values, it doesn't accelerate the convergence !coeff_rad=0.5 !coeff_rad=1. !do l=1, nlayer ! do ig=1,ngrid ! if(pplay(ig,l).ge.1.d4) then ! zdtsw(ig,l)=zdtsw(ig,l)*(pplay(ig,l)/1.d4)**coeff_rad ! zdtlw(ig,l)=zdtlw(ig,l)*(pplay(ig,l)/1.d4)**coeff_rad ! endif ! enddo !enddo ! AF24: removed CLFvarying Option ! if(ok_slab_ocean) then ! tsurf(:)=tsurf2(:) ! endif ! Radiative flux from the sky absorbed by the surface (W.m-2). GSR=0.0 fluxrad_sky(1:ngrid)=emis(1:ngrid)*fluxsurf_lw(1:ngrid)+fluxsurfabs_sw(1:ngrid) !if(noradsurf)then ! no lower surface; SW flux just disappears ! GSR = SUM(fluxsurf_sw(1:ngrid)*cell_area(1:ngrid))/totarea ! fluxrad_sky(1:ngrid)=emis(1:ngrid)*fluxsurf_lw(1:ngrid) ! print*,'SW lost in deep atmosphere = ',GSR,' W m^-2' !endif ! Net atmospheric radiative heating rate (K.s-1) dtrad(1:ngrid,1:nlayer)=zdtsw(1:ngrid,1:nlayer)+zdtlw(1:ngrid,1:nlayer) ! Late initialization of the Ice Spectral Albedo. We needed the visible bands to do that ! if (firstcall .and. albedo_spectral_mode) then call spectral_albedo_calc(albedo_snow_SPECTV) endif ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ ! II.b Call Newtonian cooling scheme !AF24: removed ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ else ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ ! II.c Atmosphere has no radiative effect ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ fluxtop_dn(1:ngrid) = fract(1:ngrid)*mu0(1:ngrid)*Fat1AU/dist_star**2 if(ngrid.eq.1)then ! / by 4 globally in 1D case... fluxtop_dn(1) = fract(1)*Fat1AU/dist_star**2/2.0 endif fluxsurf_sw(1:ngrid) = fluxtop_dn(1:ngrid) print*,'------------WARNING---WARNING------------' ! by MT2015. print*,'You are in corrk=false mode, ' print*,'and the surface albedo is taken equal to the first visible spectral value' fluxsurfabs_sw(1:ngrid) = fluxtop_dn(1:ngrid)*(1.-albedo(1:ngrid,1)) fluxrad_sky(1:ngrid) = fluxsurfabs_sw(1:ngrid) fluxtop_lw(1:ngrid) = emis(1:ngrid)*sigma*tsurf(1:ngrid)**4 dtrad(1:ngrid,1:nlayer)=0.0 ! no atmospheric radiative heating endif ! end of corrk endif ! of if(mod(icount-1,iradia).eq.0) ! Transformation of the radiative tendencies ! ------------------------------------------ zplanck(1:ngrid)=tsurf(1:ngrid)*tsurf(1:ngrid) zplanck(1:ngrid)=emis(1:ngrid)*sigma*zplanck(1:ngrid)*zplanck(1:ngrid) fluxrad(1:ngrid)=fluxrad_sky(1:ngrid)-zplanck(1:ngrid) pdt(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)+dtrad(1:ngrid,1:nlayer) ! Test of energy conservation !---------------------------- if(enertest)then call planetwide_sumval(cpp*massarea(:,:)*zdtsw(:,:)/totarea_planet,dEtotSW) call planetwide_sumval(cpp*massarea(:,:)*zdtlw(:,:)/totarea_planet,dEtotLW) !call planetwide_sumval(fluxsurf_sw(:)*(1.-albedo_equivalent(:))*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk call planetwide_sumval(fluxsurfabs_sw(:)*cell_area(:)/totarea_planet,dEtotsSW) !JL13 carefull, albedo can have changed since the last time we called corrk call planetwide_sumval((fluxsurf_lw(:)*emis(:)-zplanck(:))*cell_area(:)/totarea_planet,dEtotsLW) dEzRadsw(:,:)=cpp*mass(:,:)*zdtsw(:,:) dEzRadlw(:,:)=cpp*mass(:,:)*zdtlw(:,:) if (is_master) then print*,'---------------------------------------------------------------' print*,'In corrk SW atmospheric heating =',dEtotSW,' W m-2' print*,'In corrk LW atmospheric heating =',dEtotLW,' W m-2' print*,'atmospheric net rad heating (SW+LW) =',dEtotLW+dEtotSW,' W m-2' print*,'In corrk SW surface heating =',dEtotsSW,' W m-2' print*,'In corrk LW surface heating =',dEtotsLW,' W m-2' print*,'surface net rad heating (SW+LW) =',dEtotsLW+dEtotsSW,' W m-2' endif endif ! end of 'enertest' endif ! of if (callrad) !! call writediagfi(ngrid,"vdifc_pre_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap)) !! call writediagfi(ngrid,"vdifc_pre_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap)) ! -------------------------------------------- ! III. Vertical diffusion (turbulent mixing) : ! -------------------------------------------- if (calldifv) then zflubid(1:ngrid)=fluxrad(1:ngrid)+fluxgrd(1:ngrid) ! JL12 the following if test is temporarily there to allow us to compare the old vdifc with turbdiff. if (UseTurbDiff) then call turbdiff(ngrid,nlayer,nq, & ptimestep,capcal, & pplay,pplev,zzlay,zzlev,z0, & pu,pv,pt,zpopsk,pq,tsurf,emis,qsurf, & pdt,pdq,zflubid, & zdudif,zdvdif,zdtdif,zdtsdif, & sensibFlux,q2,zdqdif,zdqevap,zdqsdif, & taux,tauy) else zdh(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)/zpopsk(1:ngrid,1:nlayer) call vdifc(ngrid,nlayer,nq,zpopsk, & ptimestep,capcal,lwrite, & pplay,pplev,zzlay,zzlev,z0, & pu,pv,zh,pq,tsurf,emis,qsurf, & zdh,pdq,zflubid, & zdudif,zdvdif,zdhdif,zdtsdif, & sensibFlux,q2,zdqdif,zdqsdif) zdtdif(1:ngrid,1:nlayer)=zdhdif(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer) ! for diagnostic only zdqevap(1:ngrid,1:nlayer)=0. end if !end of 'UseTurbDiff' zdtsurf(1:ngrid)=zdtsurf(1:ngrid)+zdtsdif(1:ngrid) !!! this is always done, except for turbulence-resolving simulations if (.not. turb_resolved) then pdv(1:ngrid,1:nlayer)=pdv(1:ngrid,1:nlayer)+zdvdif(1:ngrid,1:nlayer) pdu(1:ngrid,1:nlayer)=pdu(1:ngrid,1:nlayer)+zdudif(1:ngrid,1:nlayer) pdt(1:ngrid,1:nlayer)=pdt(1:ngrid,1:nlayer)+zdtdif(1:ngrid,1:nlayer) endif ! if(ok_slab_ocean)then !AF24: removed ! flux_sens_lat(1:ngrid)=(zdtsdif(1:ngrid)*capcal(1:ngrid)-fluxrad(1:ngrid)) ! endif !! call writediagfi(ngrid,"vdifc_post_zdqsdif"," "," ",2,zdqsdif(1:ngrid,igcm_h2o_vap)) if (tracer) then pdq(1:ngrid,1:nlayer,1:nq)=pdq(1:ngrid,1:nlayer,1:nq)+ zdqdif(1:ngrid,1:nlayer,1:nq) dqsurf(1:ngrid,1:nq)=dqsurf(1:ngrid,1:nq) + zdqsdif(1:ngrid,1:nq) end if ! of if (tracer) !! call writediagfi(ngrid,"vdifc_post_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap)) !! call writediagfi(ngrid,"vdifc_post_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap)) ! test energy conservation !------------------------- if(enertest)then dEzdiff(:,:)=cpp*mass(:,:)*zdtdif(:,:) do ig = 1, ngrid dEdiff(ig)=SUM(dEzdiff (ig,:))+ sensibFlux(ig)! subtract flux to the ground dEzdiff(ig,1)= dEzdiff(ig,1)+ sensibFlux(ig)! subtract flux to the ground enddo call planetwide_sumval(dEdiff(:)*cell_area(:)/totarea_planet,dEtot) dEdiffs(:)=capcal(:)*zdtsdif(:)-zflubid(:)-sensibFlux(:) call planetwide_sumval(dEdiffs(:)*cell_area(:)/totarea_planet,dEtots) call planetwide_sumval(sensibFlux(:)*cell_area(:)/totarea_planet,AtmToSurf_TurbFlux) if (is_master) then if (UseTurbDiff) then print*,'In TurbDiff sensible flux (atm=>surf) =',AtmToSurf_TurbFlux,' W m-2' print*,'In TurbDiff non-cons atm nrj change =',dEtot,' W m-2' print*,'In TurbDiff (correc rad+latent heat) surf nrj change =',dEtots,' W m-2' else print*,'In vdifc sensible flux (atm=>surf) =',AtmToSurf_TurbFlux,' W m-2' print*,'In vdifc non-cons atm nrj change =',dEtot,' W m-2' print*,'In vdifc (correc rad+latent heat) surf nrj change =',dEtots,' W m-2' end if endif ! end of 'is_master' ! JL12 : note that the black body radiative flux emitted by the surface has been updated by the implicit scheme but not given back elsewhere. endif ! end of 'enertest' ! ! Test water conservation. !AF24: removed else ! calldifv ! if(.not.newtonian)then zdtsurf(1:ngrid) = zdtsurf(1:ngrid) + (fluxrad(1:ngrid) + fluxgrd(1:ngrid))/capcal(1:ngrid) endif ! end of 'calldifv' !------------------- ! IV. Convection : !------------------- ! ~~~~~~~~~~~~~~~~~~~~~~~~~~ ! IV.a Thermal plume model : AF24: removed ! ~~~~~~~~~~~~~~~~~~~~~~~~~~ ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ ! IV.b Dry convective adjustment : ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ if(calladj) then zdh(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)/zpopsk(1:ngrid,1:nlayer) zduadj(1:ngrid,1:nlayer)=0.0 zdvadj(1:ngrid,1:nlayer)=0.0 zdhadj(1:ngrid,1:nlayer)=0.0 call convadj(ngrid,nlayer,nq,ptimestep, & pplay,pplev,zpopsk, & pu,pv,zh,pq, & pdu,pdv,zdh,pdq, & zduadj,zdvadj,zdhadj, & zdqadj) pdu(1:ngrid,1:nlayer) = pdu(1:ngrid,1:nlayer) + zduadj(1:ngrid,1:nlayer) pdv(1:ngrid,1:nlayer) = pdv(1:ngrid,1:nlayer) + zdvadj(1:ngrid,1:nlayer) pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer) + zdhadj(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer) zdtadj(1:ngrid,1:nlayer) = zdhadj(1:ngrid,1:nlayer)*zpopsk(1:ngrid,1:nlayer) ! for diagnostic only if(tracer) then pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqadj(1:ngrid,1:nlayer,1:nq) end if ! Test energy conservation if(enertest)then call planetwide_sumval(cpp*massarea(:,:)*zdtadj(:,:)/totarea_planet,dEtot) if (is_master) print*,'In convadj atmospheric energy change =',dEtot,' W m-2' endif ! ! Test water conservation !AF24: removed endif ! end of 'calladj' !---------------------------------------------- ! Non orographic Gravity Waves: AF24: removed !--------------------------------------------- !----------------------------------------------- ! V. Carbon dioxide condensation-sublimation : !----------------------------------------------- if (n2cond) then if (.not.tracer) then print*,'We need a N2 ice tracer to condense N2' call abort endif call condense_n2(ngrid,nlayer,nq,ptimestep, & capcal,pplay,pplev,tsurf,pt, & pdt,zdtsurf,pq,pdq, & qsurf,zdqsurfc,albedo,emis, & albedo_bareground,albedo_n2_ice_SPECTV, & zdtc,zdtsurfc,pdpsrf,zdqc) pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+zdtc(1:ngrid,1:nlayer) zdtsurf(1:ngrid) = zdtsurf(1:ngrid) + zdtsurfc(1:ngrid) pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq)+ zdqc(1:ngrid,1:nlayer,1:nq) dqsurf(1:ngrid,igcm_n2_ice) = dqsurf(1:ngrid,igcm_n2_ice) + zdqsurfc(1:ngrid) !! call writediagfi(ngrid,"condense_n2_post_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap)) !! call writediagfi(ngrid,"condense_n2_post_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap)) ! test energy conservation if(enertest)then call planetwide_sumval(cpp*massarea(:,:)*zdtc(:,:)/totarea_planet,dEtot) call planetwide_sumval(capcal(:)*zdtsurfc(:)*cell_area(:)/totarea_planet,dEtots) if (is_master) then print*,'In n2cloud atmospheric energy change =',dEtot,' W m-2' print*,'In n2cloud surface energy change =',dEtots,' W m-2' endif endif endif ! end of 'n2cond' !--------------------------------------------- ! VI. Specific parameterizations for tracers !--------------------------------------------- if (tracer) then ! --------------------- ! VI.1. Water and ice !AF24: removed ! --------------------- ! ------------------------- ! VI.2. Photochemistry !AF24: removed ! ------------------------- ! ------------------------- ! VI.3. Aerosol particles ! ------------------------- !Generic Condensation if (generic_condensation) then call condensation_generic(ngrid,nlayer,nq,ptimestep,pplev,pplay, & pt,pq,pdt,pdq,dt_generic_condensation, & dqvaplscale_generic,dqcldlscale_generic,rneb_generic) pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer)+dt_generic_condensation(1:ngrid,1:nlayer) pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq)+dqvaplscale_generic(1:ngrid,1:nlayer,1:nq) pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq)+dqcldlscale_generic(1:ngrid,1:nlayer,1:nq) if(enertest)then do ig=1,ngrid genericconddE(ig) = cpp*SUM(mass(:,:)*dt_generic_condensation(:,:)) enddo call planetwide_sumval(cpp*massarea(:,:)*dt_generic_condensation(:,:)/totarea_planet,dEtot) if (is_master) print*,'In generic condensation atmospheric energy change =',dEtot,' W m-2' end if ! if (.not. water) then ! Compute GCS (Generic Condensable Specie) cloud fraction. For now we can not have both water cloud fraction and GCS cloud fraction ! Water is the priority ! If you have set water and generic_condensation, then cloudfrac will be water cloudfrac ! ! If you have set generic_condensation (and not water) and you have set several GCS ! then cloudfrac will be only the cloudfrac of the last generic tracer ! (Because it is rewritten every tracer in the loop) ! ! Maybe one should create a cloudfrac_generic(ngrid,nlayer,nq) with 3 dimensions, the last one for tracers ! Let's loop on tracers cloudfrac(:,:)=0.0 do iq=1,nq call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic) if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer do l = 1, nlayer do ig=1,ngrid cloudfrac(ig,l)=rneb_generic(ig,l,iq) enddo enddo endif end do ! do iq=1,nq loop on tracers ! endif ! .not. water endif !generic_condensation !Generic Rain !AF24: removed ! Sedimentation. if (sedimentation) then zdqsed(1:ngrid,1:nlayer,1:nq) = 0.0 zdqssed(1:ngrid,1:nq) = 0.0 ! if(watertest)then !AF24: removed call callsedim(ngrid,nlayer,ptimestep, & pplev,zzlev,pt,pdt,pq,pdq, & zdqsed,zdqssed,nq) ! if(watertest)then !AF24: removed ! Whether it falls as rain or snow depends only on the surface temperature pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqsed(1:ngrid,1:nlayer,1:nq) dqsurf(1:ngrid,1:nq) = dqsurf(1:ngrid,1:nq) + zdqssed(1:ngrid,1:nq) !! call writediagfi(ngrid,"callsedim_post_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap)) ! ! Test water conservation !AF24: removed end if ! end of 'sedimentation' !! call writediagfi(ngrid,"mass_redist_pre_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap)) !! call writediagfi(ngrid,"mass_redist_pre_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap)) ! --------------- ! VI.4. Updates ! --------------- ! Updating Atmospheric Mass and Tracers budgets. if(mass_redistrib) then zdmassmr(1:ngrid,1:nlayer) = mass(1:ngrid,1:nlayer) * 0 ! ( zdqevap(1:ngrid,1:nlayer) & ! ! + zdqrain(1:ngrid,1:nlayer,igcm_h2o_vap) & ! ! + dqmoist(1:ngrid,1:nlayer,igcm_h2o_vap) & ! + dqvaplscale(1:ngrid,1:nlayer) ) do ig = 1, ngrid zdmassmr_col(ig)=SUM(zdmassmr(ig,1:nlayer)) enddo ! call writediagfi(ngrid,"mass_evap","mass gain"," ",3,zdmassmr) ! call writediagfi(ngrid,"mass_evap_col","mass gain col"," ",2,zdmassmr_col) call writediagfi(ngrid,"mass","mass","kg/m2",3,mass) call mass_redistribution(ngrid,nlayer,nq,ptimestep, & capcal,pplay,pplev,pt,tsurf,pq,qsurf, & pu,pv,pdt,zdtsurf,pdq,pdu,pdv,zdmassmr, & zdtmr,zdtsurfmr,zdpsrfmr,zdumr,zdvmr,zdqmr,zdqsurfmr) pdq(1:ngrid,1:nlayer,1:nq) = pdq(1:ngrid,1:nlayer,1:nq) + zdqmr(1:ngrid,1:nlayer,1:nq) dqsurf(1:ngrid,1:nq) = dqsurf(1:ngrid,1:nq) + zdqsurfmr(1:ngrid,1:nq) pdt(1:ngrid,1:nlayer) = pdt(1:ngrid,1:nlayer) + zdtmr(1:ngrid,1:nlayer) pdu(1:ngrid,1:nlayer) = pdu(1:ngrid,1:nlayer) + zdumr(1:ngrid,1:nlayer) pdv(1:ngrid,1:nlayer) = pdv(1:ngrid,1:nlayer) + zdvmr(1:ngrid,1:nlayer) pdpsrf(1:ngrid) = pdpsrf(1:ngrid) + zdpsrfmr(1:ngrid) zdtsurf(1:ngrid) = zdtsurf(1:ngrid) + zdtsurfmr(1:ngrid) endif ! call writediagfi(ngrid,"mass_redistribution_post_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap)) !! call writediagfi(ngrid,"slab_pre_dqsurf"," "," ",2,dqsurf(1:ngrid,igcm_h2o_vap)) !! call writediagfi(ngrid,"slab_pre_qsurf"," "," ",2,qsurf(1:ngrid,igcm_h2o_vap)) ! ------------------ ! VI.5. Slab Ocean !AF24: removed ! ------------------ ! ----------------------------- ! VI.6. Surface Tracer Update ! ----------------------------- ! AF24: deleted hydrology qsurf(1:ngrid,1:nq) = qsurf(1:ngrid,1:nq) + ptimestep*dqsurf(1:ngrid,1:nq) ! Add qsurf to qsurf_hist, which is what we save in diagfi.nc. At the same time, we set the water ! content of ocean gridpoints back to zero, in order to avoid rounding errors in vdifc, rain. qsurf_hist(:,:) = qsurf(:,:) ! if(ice_update)then ! ice_min(1:ngrid)=min(ice_min(1:ngrid),qsurf(1:ngrid,igcm_h2o_ice)) ! endif endif! end of if 'tracer' !------------------------------------------------ ! VII. Surface and sub-surface soil temperature !------------------------------------------------ ! ! Increment surface temperature ! if(ok_slab_ocean)then !AF24: removed tsurf(1:ngrid)=tsurf(1:ngrid)+ptimestep*zdtsurf(1:ngrid) ! Compute soil temperatures and subsurface heat flux. if (callsoil) then call soil(ngrid,nsoilmx,.false.,lastcall,inertiedat, & ptimestep,tsurf,tsoil,capcal,fluxgrd) endif ! if (ok_slab_ocean) then !AF24: removed ! Test energy conservation if(enertest)then call planetwide_sumval(cell_area(:)*capcal(:)*zdtsurf(:)/totarea_planet,dEtots) if (is_master) print*,'Surface energy change =',dEtots,' W m-2' endif !--------------------------------------------------- ! VIII. Perform diagnostics and write output files !--------------------------------------------------- ! Note : For output only: the actual model integration is performed in the dynamics. ! Temperature, zonal and meridional winds. zt(1:ngrid,1:nlayer) = pt(1:ngrid,1:nlayer) + pdt(1:ngrid,1:nlayer)*ptimestep zu(1:ngrid,1:nlayer) = pu(1:ngrid,1:nlayer) + pdu(1:ngrid,1:nlayer)*ptimestep zv(1:ngrid,1:nlayer) = pv(1:ngrid,1:nlayer) + pdv(1:ngrid,1:nlayer)*ptimestep !! Recast thermal plume vertical velocity array for outputs !! AF24: removed ! Diagnostic. zdtdyn(1:ngrid,1:nlayer) = (pt(1:ngrid,1:nlayer)-ztprevious(1:ngrid,1:nlayer)) / ptimestep ztprevious(1:ngrid,1:nlayer) = zt(1:ngrid,1:nlayer) zdudyn(1:ngrid,1:nlayer) = (pu(1:ngrid,1:nlayer)-zuprevious(1:ngrid,1:nlayer)) / ptimestep zuprevious(1:ngrid,1:nlayer) = zu(1:ngrid,1:nlayer) if(firstcall)then zdtdyn(1:ngrid,1:nlayer)=0.0 zdudyn(1:ngrid,1:nlayer)=0.0 endif ! Dynamical heating diagnostic. do ig=1,ngrid fluxdyn(ig)= SUM(zdtdyn(ig,:) *mass(ig,:))*cpp enddo ! Tracers. zq(1:ngrid,1:nlayer,1:nq) = pq(1:ngrid,1:nlayer,1:nq) + pdq(1:ngrid,1:nlayer,1:nq)*ptimestep ! Surface pressure. ps(1:ngrid) = pplev(1:ngrid,1) + pdpsrf(1:ngrid)*ptimestep ! Surface and soil temperature information call planetwide_sumval(cell_area(:)*tsurf(:)/totarea_planet,Ts1) call planetwide_minval(tsurf(:),Ts2) call planetwide_maxval(tsurf(:),Ts3) if(callsoil)then TsS = SUM(cell_area(:)*tsoil(:,nsoilmx))/totarea ! mean temperature at bottom soil layer if (is_master) then print*,' ave[Tsurf] min[Tsurf] max[Tsurf] ave[Tdeep]' print*,Ts1,Ts2,Ts3,TsS end if else if (is_master) then print*,' ave[Tsurf] min[Tsurf] max[Tsurf]' print*,Ts1,Ts2,Ts3 endif end if ! Check the energy balance of the simulation during the run if(corrk)then call planetwide_sumval(cell_area(:)*fluxtop_dn(:)/totarea_planet,ISR) call planetwide_sumval(cell_area(:)*fluxabs_sw(:)/totarea_planet,ASR) call planetwide_sumval(cell_area(:)*fluxtop_lw(:)/totarea_planet,OLR) call planetwide_sumval(cell_area(:)*fluxgrd(:)/totarea_planet,GND) call planetwide_sumval(cell_area(:)*fluxdyn(:)/totarea_planet,DYN) do ig=1,ngrid if(fluxtop_dn(ig).lt.0.0)then print*,'fluxtop_dn has gone crazy' print*,'fluxtop_dn=',fluxtop_dn(ig) print*,'tau_col=',tau_col(ig) print*,'aerosol=',aerosol(ig,:,:) print*,'temp= ',pt(ig,:) print*,'pplay= ',pplay(ig,:) call abort endif end do if(ngrid.eq.1)then DYN=0.0 endif if (is_master) then print*,' ISR ASR OLR GND DYN [W m^-2]' print*, ISR,ASR,OLR,GND,DYN endif if(enertest .and. is_master)then print*,'SW flux/heating difference SW++ - ASR = ',dEtotSW+dEtotsSW-ASR,' W m-2' print*,'LW flux/heating difference LW++ - OLR = ',dEtotLW+dEtotsLW+OLR,' W m-2' print*,'LW energy balance LW++ + ASR = ',dEtotLW+dEtotsLW+ASR,' W m-2' endif if(meanOLR .and. is_master)then if((ngrid.gt.1) .or. (mod(icount-1,ecritphy).eq.0))then ! to record global radiative balance open(92,file="rad_bal.out",form='formatted',position='append') write(92,*) zday,ISR,ASR,OLR close(92) open(93,file="tem_bal.out",form='formatted',position='append') if(callsoil)then write(93,*) zday,Ts1,Ts2,Ts3,TsS else write(93,*) zday,Ts1,Ts2,Ts3 endif close(93) endif endif endif ! end of 'corrk' ! Diagnostic to test radiative-convective timescales in code. if(testradtimes)then open(38,file="tau_phys.out",form='formatted',position='append') ig=1 do l=1,nlayer write(38,*) -1./pdt(ig,l),pt(ig,l),pplay(ig,l) enddo close(38) print*,'As testradtimes enabled,' print*,'exiting physics on first call' call abort endif ! Compute column amounts (kg m-2) if tracers are enabled. if(tracer)then qcol(1:ngrid,1:nq)=0.0 do iq=1,nq do ig=1,ngrid qcol(ig,iq) = SUM( zq(ig,1:nlayer,iq) * mass(ig,1:nlayer)) enddo enddo ! Generalised for arbitrary aerosols now. By LK reffcol(1:ngrid,1:naerkind)=0.0 if(n2cond.and.(iaero_n2.ne.0))then call n2_reffrad(ngrid,nlayer,nq,zq,reffrad(1,1,iaero_n2)) do ig=1,ngrid reffcol(ig,iaero_n2) = SUM(zq(ig,1:nlayer,igcm_n2_ice)*reffrad(ig,1:nlayer,iaero_n2)*mass(ig,1:nlayer)) enddo endif ! if(water.and.(iaero_h2o.ne.0))then !AF24: removed endif ! end of 'tracer' ! ! Test for water conservation. !AF24: removed ! Calculate RH_generic (Generic Relative Humidity) for diagnostic. if(generic_condensation)then RH_generic(:,:,:)=0.0 do iq=1,nq call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic) if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer do l = 1, nlayer do ig=1,ngrid call Psat_generic(zt(ig,l),pplay(ig,l),metallicity,psat_tmp_generic,qsat_generic(ig,l,iq)) RH_generic(ig,l,iq) = zq(ig,l,igcm_generic_vap) / qsat_generic(ig,l,iq) enddo enddo end if end do ! iq=1,nq endif ! end of 'generic_condensation' if (is_master) print*,'--> Ls =',zls*180./pi !---------------------------------------------------------------------- ! Writing NetCDF file "RESTARTFI" at the end of the run !---------------------------------------------------------------------- ! Note: 'restartfi' is stored just before dynamics are stored ! in 'restart'. Between now and the writting of 'restart', ! there will have been the itau=itau+1 instruction and ! a reset of 'time' (lastacll = .true. when itau+1= itaufin) ! thus we store for time=time+dtvr if(lastcall) then ztime_fin = ptime + ptimestep/(float(iphysiq)*daysec) !! Update surface ice distribution to iterate to steady state if requested !! AF24: removed ! endif #ifndef MESOSCALE if (ngrid.ne.1) then write(*,*)'PHYSIQ: for physdem ztime_fin =',ztime_fin call physdem1("restartfi.nc",nsoilmx,ngrid,nlayer,nq, & ptimestep,ztime_fin, & tsurf,tsoil,emis,q2,qsurf_hist) endif #endif ! if(ok_slab_ocean) then ! call ocean_slab_final!(tslab, seaice) ! end if endif ! end of 'lastcall' ! ----------------------------------------------------------------- ! WSTATS: Saving statistics ! ----------------------------------------------------------------- ! ("stats" stores and accumulates key variables in file "stats.nc" ! which can later be used to make the statistic files of the run: ! if flag "callstats" from callphys.def is .true.) call wstats(ngrid,"ps","Surface pressure","Pa",2,ps) call wstats(ngrid,"tsurf","Surface temperature","K",2,tsurf) call wstats(ngrid,"fluxsurf_lw", & "Thermal IR radiative flux to surface","W.m-2",2, & fluxsurf_lw) call wstats(ngrid,"fluxtop_lw", & "Thermal IR radiative flux to space","W.m-2",2, & fluxtop_lw) ! call wstats(ngrid,"fluxsurf_sw", & ! "Solar radiative flux to surface","W.m-2",2, & ! fluxsurf_sw_tot) ! call wstats(ngrid,"fluxtop_sw", & ! "Solar radiative flux to space","W.m-2",2, & ! fluxtop_sw_tot) call wstats(ngrid,"ISR","incoming stellar rad.","W m-2",2,fluxtop_dn) call wstats(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw) call wstats(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw) !call wstats(ngrid,"ALB","Surface albedo"," ",2,albedo_equivalent) !call wstats(ngrid,"ALB_1st","First Band Surface albedo"," ",2,albedo(:,1)) call wstats(ngrid,"p","Pressure","Pa",3,pplay) call wstats(ngrid,"temp","Atmospheric temperature","K",3,zt) call wstats(ngrid,"u","Zonal (East-West) wind","m.s-1",3,zu) call wstats(ngrid,"v","Meridional (North-South) wind","m.s-1",3,zv) call wstats(ngrid,"w","Vertical (down-up) wind","m.s-1",3,pw) call wstats(ngrid,"q2","Boundary layer eddy kinetic energy","m2.s-2",3,q2) if (tracer) then do iq=1,nq call wstats(ngrid,noms(iq),noms(iq),'kg/kg',3,zq(1,1,iq)) call wstats(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf', & 'kg m^-2',2,qsurf(1,iq) ) call wstats(ngrid,trim(noms(iq))//'_col',trim(noms(iq))//'_col', & 'kg m^-2',2,qcol(1,iq) ) ! call wstats(ngrid,trim(noms(iq))//'_reff', & ! trim(noms(iq))//'_reff', & ! 'm',3,reffrad(1,1,iq)) end do endif ! end of 'tracer' !AF24: deleted slab ocean and water if(lastcall.and.callstats) then write (*,*) "Writing stats..." call mkstats(ierr) endif #ifndef MESOSCALE !----------------------------------------------------------------------------------------------------- ! OUTPUT in netcdf file "DIAGFI.NC", containing any variable for diagnostic ! ! Note 1 : output with period "ecritphy", set in "run.def" ! ! Note 2 : writediagfi can also be called from any other subroutine for any variable, ! but its preferable to keep all the calls in one place ... !----------------------------------------------------------------------------------------------------- call writediagfi(ngrid,"Ls","solar longitude","deg",0,zls*180./pi) call writediagfi(ngrid,"Lss","sub solar longitude","deg",0,zlss*180./pi) call writediagfi(ngrid,"RA","right ascension","deg",0,right_ascen*180./pi) call writediagfi(ngrid,"Declin","solar declination","deg",0,declin*180./pi) call writediagfi(ngrid,"tsurf","Surface temperature","K",2,tsurf) call writediagfi(ngrid,"ps","Surface pressure","Pa",2,ps) call writediagfi(ngrid,"temp","temperature","K",3,zt) call writediagfi(ngrid,"teta","potential temperature","K",3,zh) call writediagfi(ngrid,"u","Zonal wind","m.s-1",3,zu) call writediagfi(ngrid,"v","Meridional wind","m.s-1",3,zv) call writediagfi(ngrid,"w","Vertical wind","m.s-1",3,pw) call writediagfi(ngrid,"p","Pressure","Pa",3,pplay) ! Subsurface temperatures ! call writediagsoil(ngrid,"tsurf","Surface temperature","K",2,tsurf) ! call writediagsoil(ngrid,"temp","temperature","K",3,tsoil) ! ! Oceanic layers !AF24: removed ! ! Thermal plume model !AF24: removed ! GW non-oro outputs !AF24: removed ! Total energy balance diagnostics if(callrad)then !call writediagfi(ngrid,"ALB","Surface albedo"," ",2,albedo_equivalent) !call writediagfi(ngrid,"ALB_1st","First Band Surface albedo"," ",2,albedo(:,1)) call writediagfi(ngrid,"ISR","incoming stellar rad.","W m-2",2,fluxtop_dn) call writediagfi(ngrid,"ASR","absorbed stellar rad.","W m-2",2,fluxabs_sw) call writediagfi(ngrid,"OLR","outgoing longwave rad.","W m-2",2,fluxtop_lw) call writediagfi(ngrid,"shad","rings"," ", 2, fract) ! call writediagfi(ngrid,"ASRcs","absorbed stellar rad (cs).","W m-2",2,fluxabs_sw1) ! call writediagfi(ngrid,"OLRcs","outgoing longwave rad (cs).","W m-2",2,fluxtop_lw1) ! call writediagfi(ngrid,"fluxsurfsw","sw surface flux.","W m-2",2,fluxsurf_sw) ! call writediagfi(ngrid,"fluxsurflw","lw back radiation.","W m-2",2,fluxsurf_lw) ! call writediagfi(ngrid,"fluxsurfswcs","sw surface flux (cs).","W m-2",2,fluxsurf_sw1) ! call writediagfi(ngrid,"fluxsurflwcs","lw back radiation (cs).","W m-2",2,fluxsurf_lw1) ! if(ok_slab_ocean) then ! call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrdocean) ! else call writediagfi(ngrid,"GND","heat flux from ground","W m-2",2,fluxgrd) ! endif call writediagfi(ngrid,"DYN","dynamical heat input","W m-2",2,fluxdyn) endif ! end of 'callrad' if(enertest) then if (calldifv) then call writediagfi(ngrid,"q2","turbulent kinetic energy","J.kg^-1",3,q2) call writediagfi(ngrid,"sensibFlux","sensible heat flux","w.m^-2",2,sensibFlux) ! call writediagfi(ngrid,"dEzdiff","turbulent diffusion heating (-sensible flux)","w.m^-2",3,dEzdiff) ! call writediagfi(ngrid,"dEdiff","integrated turbulent diffusion heating (-sensible flux)","w.m^-2",2,dEdiff) ! call writediagfi(ngrid,"dEdiffs","In TurbDiff (correc rad+latent heat) surf nrj change","w.m^-2",2,dEdiffs) endif if (corrk) then call writediagfi(ngrid,"dEzradsw","radiative heating","w.m^-2",3,dEzradsw) call writediagfi(ngrid,"dEzradlw","radiative heating","w.m^-2",3,dEzradlw) endif ! if(watercond) then !AF24: removed if (generic_condensation) then call writediagfi(ngrid,"genericconddE","heat from generic condensation","W m-2",2,genericconddE) call writediagfi(ngrid,"dt_generic_condensation","heating from generic condensation","K s-1",3,dt_generic_condensation) endif endif ! end of 'enertest' ! Diagnostics of optical thickness ! Warning this is exp(-tau), I let you postproc with -log to have tau itself - JVO 19 if (diagdtau) then do nw=1,L_NSPECTV write(str2,'(i2.2)') nw call writediagfi(ngrid,'dtauv'//str2,'Layer optical thickness attenuation in VI band '//str2,'',1,int_dtauv(:,nlayer:1:-1,nw)) enddo do nw=1,L_NSPECTI write(str2,'(i2.2)') nw call writediagfi(ngrid,'dtaui'//str2,'Layer optical thickness attenuation in IR band '//str2,'',1,int_dtaui(:,nlayer:1:-1,nw)) enddo endif ! Temporary inclusions for heating diagnostics. call writediagfi(ngrid,"zdtsw","SW heating","T s-1",3,zdtsw) call writediagfi(ngrid,"zdtlw","LW heating","T s-1",3,zdtlw) call writediagfi(ngrid,"dtrad","radiative heating","K s-1",3,dtrad) call writediagfi(ngrid,"zdtdyn","Dyn. heating","T s-1",3,zdtdyn) ! For Debugging. !call writediagfi(ngrid,'rnat','Terrain type',' ',2,real(rnat)) !call writediagfi(ngrid,'pphi','Geopotential',' ',3,pphi) ! Output aerosols. if (igcm_n2_ice.ne.0.and.iaero_n2.ne.0) & call writediagfi(ngrid,'N2ice_reff','N2ice_reff','m',3,reffrad(1,1,iaero_n2)) if (igcm_n2_ice.ne.0.and.iaero_n2.ne.0) & call writediagfi(ngrid,'N2ice_reffcol','N2ice_reffcol','um kg m^-2',2,reffcol(1,iaero_n2)) ! if (igcm_h2o_ice.ne.0.and.iaero_h2o.ne.0) & !AF24: removed ! Output tracers. if (tracer) then do iq=1,nq call writediagfi(ngrid,noms(iq),noms(iq),'kg/kg',3,zq(1,1,iq)) call writediagfi(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf', & 'kg m^-2',2,qsurf_hist(1,iq) ) call writediagfi(ngrid,trim(noms(iq))//'_col',trim(noms(iq))//'_col', & 'kg m^-2',2,qcol(1,iq) ) ! call writediagfi(ngrid,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf', & ! 'kg m^-2',2,qsurf(1,iq) ) ! if(watercond.or.CLFvarying)then !AF24: removed if(generic_condensation)then call writediagfi(ngrid,"rneb_generic","GCS cloud fraction (generic condensation)"," ",3,rneb_generic) call writediagfi(ngrid,"CLF","GCS cloud fraction"," ",3,cloudfrac) call writediagfi(ngrid,"RH_generic","GCS relative humidity"," ",3,RH_generic) endif ! if(generic_rain)then !AF24: removed ! if((hydrology).and.(.not.ok_slab_ocean))then !AF24: removed call writediagfi(ngrid,"tau_col","Total aerosol optical depth","[]",2,tau_col) enddo ! end of 'nq' loop endif ! end of 'tracer' ! Output spectrum. if(specOLR.and.corrk)then call writediagspecIR(ngrid,"OLR3D","OLR(lon,lat,band)","W/m^2/cm^-1",3,OLR_nu) call writediagspecVI(ngrid,"OSR3D","OSR(lon,lat,band)","W/m^2/cm^-1",3,OSR_nu) call writediagspecVI(ngrid,"GSR3D","GSR(lon,lat,band)","W/m^2/cm^-1",3,GSR_nu) endif #else comm_HR_SW(1:ngrid,1:nlayer) = zdtsw(1:ngrid,1:nlayer) comm_HR_LW(1:ngrid,1:nlayer) = zdtlw(1:ngrid,1:nlayer) comm_ALBEQ(1:ngrid)=albedo_equivalent(1:ngrid) if (.not.calldifv) comm_LATENT_HF(:)=0.0 ! if ((tracer).and.(water)) then !AF24: removed if ((tracer).and.(generic_condensation)) then ! .and.(.not. water) ! If you have set generic_condensation (and not water) and you have set several GCS ! then the outputs given to WRF will be only the ones for the last generic tracer ! (Because it is rewritten every tracer in the loop) ! WRF can take only one moist tracer do iq=1,nq call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic) if (call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer reffrad_generic_zeros_for_wrf(:,:) = 1. comm_CLOUDFRAC(1:ngrid,1:nlayer) = cloudfrac(1:ngrid,1:nlayer) comm_TOTCLOUDFRAC(1:ngrid) = totcloudfrac(1:ngrid) !?????? comm_SURFRAIN(1:ngrid) = zdqsrain_generic(1:ngrid,iq) comm_DQVAP(1:ngrid,1:nlayer) = pdq(1:ngrid,1:nlayer,igcm_generic_vap) comm_DQICE(1:ngrid,1:nlayer)=pdq(1:ngrid,1:nlayer,igcm_generic_ice) ! comm_H2OICE_REFF(1:ngrid,1:nlayer) = reffrad_generic_zeros_for_wrf(1:ngrid,1:nlayer) ! for now zeros ! !comm_H2OICE_REFF(1:ngrid,1:nlayer) = 0*zdtrain_generic(1:ngrid,1:nlayer) ! for now zeros ! comm_REEVAP(1:ngrid) = reevap_precip_generic(1:ngrid,iq) comm_DTRAIN(1:ngrid,1:nlayer) = zdtrain_generic(1:ngrid,1:nlayer) comm_DTLSC(1:ngrid,1:nlayer) = dt_generic_condensation(1:ngrid,1:nlayer) comm_RH(1:ngrid,1:nlayer) = RH_generic(1:ngrid,1:nlayer,iq) endif end do ! do iq=1,nq loop on tracers else comm_CLOUDFRAC(1:ngrid,1:nlayer)=0. comm_TOTCLOUDFRAC(1:ngrid)=0. comm_SURFRAIN(1:ngrid)=0. comm_DQVAP(1:ngrid,1:nlayer)=0. comm_DQICE(1:ngrid,1:nlayer)=0. ! comm_H2OICE_REFF(1:ngrid,1:nlayer)=0. comm_REEVAP(1:ngrid)=0. comm_DTRAIN(1:ngrid,1:nlayer)=0. comm_DTLSC(1:ngrid,1:nlayer)=0. comm_RH(1:ngrid,1:nlayer)=0. endif ! if water, if generic_condensation, else comm_FLUXTOP_DN(1:ngrid)=fluxtop_dn(1:ngrid) comm_FLUXABS_SW(1:ngrid)=fluxabs_sw(1:ngrid) comm_FLUXTOP_LW(1:ngrid)=fluxtop_lw(1:ngrid) comm_FLUXSURF_SW(1:ngrid)=fluxsurf_sw(1:ngrid) comm_FLUXSURF_LW(1:ngrid)=fluxsurf_lw(1:ngrid) comm_FLXGRD(1:ngrid)=fluxgrd(1:ngrid) sensibFlux(1:ngrid) = zflubid(1:ngrid) - capcal(1:ngrid)*zdtsdif(1:ngrid) !!! ???? comm_HR_DYN(1:ngrid,1:nlayer) = zdtdyn(1:ngrid,1:nlayer) ! if (turb_resolved) then ! open(17,file='lsf.txt',form='formatted',status='old') ! rewind(17) ! DO l=1,nlayer ! read(17,*) lsf_dt(l),lsf_dq(l) ! ENDDO ! close(17) ! do ig=1,ngrid ! if ((tracer).and.(water)) then ! pdq(ig,:,igcm_h2o_vap) = pdq(ig,:,igcm_h2o_vap) + lsf_dq(:) ! endif ! pdt(ig,:) = pdt(ig,:) + lsf_dt(:) ! comm_HR_DYN(ig,:) = lsf_dt(:) ! enddo ! endif #endif ! XIOS outputs #ifdef CPP_XIOS ! Send fields to XIOS: (NB these fields must also be defined as ! in context_lmdz_physics.xml to be correctly used) CALL send_xios_field("ls",zls) CALL send_xios_field("ps",ps) CALL send_xios_field("area",cell_area) CALL send_xios_field("p",pplay) CALL send_xios_field("temperature",zt) CALL send_xios_field("u",zu) CALL send_xios_field("v",zv) CALL send_xios_field("omega",omega) ! IF (calltherm) THEN !AF24: removed ! IF (water) THEN !AF24: removed CALL send_xios_field("ISR",fluxtop_dn) CALL send_xios_field("OLR",fluxtop_lw) CALL send_xios_field("ASR",fluxabs_sw) if (specOLR .and. corrk) then call send_xios_field("OLR3D",OLR_nu) call send_xios_field("IR_Bandwidth",DWNI) call send_xios_field("VI_Bandwidth",DWNV) call send_xios_field("OSR3D",OSR_nu) call send_xios_field("GSR3D",GSR_nu) endif if (lastcall.and.is_omp_master) then write(*,*) "physiq: call xios_context_finalize" call xios_context_finalize endif #endif if (check_physics_outputs) then ! Check the validity of updated fields at the end of the physics step call check_physics_fields("end of physiq:", zt, zu, zv, pplev, zq) endif icount=icount+1 end subroutine physiq end module physiq_mod