subroutine physiq(ngrid,nlayer,nq, & firstcall,lastcall, & pday,ptime,ptimestep, & pplev,pplay,pphi, & pu,pv,pt,pq, & pw, & pdu,pdv,pdt,pdq,pdpsrf,tracerdyn) use radinc_h, only : naerkind use watercommon_h, only : mx_eau_sol, To, RLVTT 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: ! ! 1. Initialization: ! 1.1 Firstcall initializations ! 1.2 Initialization for every call to physiq ! 1.2.5 Compute mean mass and cp, R and thermal conduction coeff. ! 2. Compute radiative transfer tendencies ! (longwave and shortwave). ! 4. Vertical diffusion (turbulent mixing): ! 5. Convective adjustment ! 6. Condensation and sublimation of gases (currently just CO2). ! 7. TRACERS ! 7a. water and water ice ! 7c. other schemes for tracer transport (lifting, sedimentation) ! 7d. updates (pressure variations, surface budget) ! 9. Surface and sub-surface temperature calculations ! 10. Write outputs : ! - "startfi", "histfi" if it's time ! - Saving statistics if "callstats = .true." ! - Output any needed variables in "diagfi" ! 10. Diagnostics: mass conservation of tracers, radiative energy balance etc. ! ! arguments ! --------- ! ! input ! ----- ! ecri period (in dynamical timestep) to write output ! ngrid Size of the horizontal grid. ! All internal loops are performed on that 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 (0 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 if(ice_update)then do ig = 1, ngrid ! if ice amount is declining, set it to zero if(qsurf_hist(ig,igcm_h2o_ice).lt.ice_initial(ig))then qsurf_hist(ig,igcm_h2o_ice)=0.0 qsurf_hist(ig,igcm_h2o_vap)=0.0 endif ! if ice amount is increasing and was never less than 0.5 kg/m2, ! set value in gridbox to 1x10^3 kg/m2 !if(qsurf_hist(ig,igcm_h2o_ice).gt.0.0 .and. ice_min(ig).gt.0.5)then if(qsurf_hist(ig,igcm_h2o_ice).gt.ice_initial(ig))then if(ice_min(ig).gt.0.5)then ! ignore seasonal ice qsurf_hist(ig,igcm_h2o_ice)=1.e3 qsurf_hist(ig,igcm_h2o_vap)=0.0 endif endif enddo endif write(*,*)'PHYSIQ: for physdem ztime_fin =',ztime_fin call physdem1("restartfi.nc",long,lati,nsoilmx,nq, & ! ptimestep,pday,ztime_fin,tsurf,tsoil,emis,q2,qsurf, & ptimestep,pday,ztime_fin,tsurf,tsoil,emis,q2,qsurf_hist, & area,albedodat,inertiedat,zmea,zstd,zsig,zgam,zthe, & cloudfrac,totcloudfrac,hice) endif ! ----------------------------------------------------------------- ! Saving statistics : ! ----------------------------------------------------------------- ! ("stats" stores and accumulates 8 key variables in file "stats.nc" ! which can later be used to make the statistic files of the run: ! "stats") only possible in 3D runs ! if (callstats) then 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,"fluxsurf_sw", & ! "Solar radiative flux to surface","W.m-2",2, & ! fluxsurf_sw_tot) call wstats(ngrid,"fluxtop_lw", & "Thermal IR radiative flux to space","W.m-2",2, & fluxtop_lw) ! call wstats(ngrid,"fluxtop_sw", & ! "Solar radiative flux to space","W.m-2",2, & ! fluxtop_sw_tot) 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 if (water) then vmr=zq(1:ngridmx,1:nlayermx,igcm_h2o_vap)*mugaz/mmol(igcm_h2o_vap) call wstats(ngrid,"vmr_h2ovapor", & "H2O vapour volume mixing ratio","mol/mol", & 3,vmr) endif ! of if (water) endif !tracer if(lastcall) then write (*,*) "Writing stats..." call mkstats(ierr) endif endif !if callstats ! ---------------------------------------------------------------------- ! output in netcdf file "DIAGFI", containing any variable for diagnostic ! (output with period "ecritphy", set in "run.def") ! ---------------------------------------------------------------------- ! 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,"tsurf","Surface temperature","K",2,tsurf) call writediagfi(ngrid,"ps","Surface pressure","Pa",2,ps) if(1.eq.2)then call writediagfi(ngrid,"temp","temperature","K",3,zt) 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) endif ! Total energy balance diagnostics if(callrad.and.(.not.newtonian))then call writediagfi(ngrid,'ALB','Surface albedo',' ',2,albedo) 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,"GND","heat flux from ground","W m-2",2,fluxgrd) call writediagfi(ngrid,"DYN","dynamical heat input","W m-2",2,fluxdyn) !call writediagfi(ngrid,"lscaledE","heat from largescale","W m-2",2,lscaledE) !call writediagfi(ngrid,"madjdE","heat from moistadj","W m-2",2,madjdE) !call writediagfi(ngrid,"vdifcdE","heat from vdifc surface","W m-2",2,vdifcdE) endif ! Temporary inclusions for heating diagnostics ! call writediagfi(ngrid,"zdtdyn","Dyn. heating","T s-1",3,zdtdyn) ! 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) ! debugging !call writediagfi(ngrid,"vdifNC","H2O loss vdifc","kg m-2 s-1",2,vdifcncons) !call writediagfi(ngrid,"cadjNC","H2O loss convadj","kg m-2 s-1",2,cadjncons) call writediagfi(ngrid,'rnat','Terrain type',' ',2,real(rnat)) !call writediagfi(ngrid,'pphi','Geopotential',' ',3,pphi) !print*,vdifcncons !call abort ! Output aerosols ! call writediagfi(ngridmx,'CO2ice_reff','CO2ice_reff','m',3,reffrad(1,1,1)) ! call writediagfi(ngridmx,'H2Oice_reff','H2Oice_reff','m',3,reffrad(1,1,2)) ! call writediagfi(ngridmx,'CO2ice_reffcol','CO2ice_reffcol','um kg m^-2',2,reffcol(1,1)) ! call writediagfi(ngridmx,'H2Oice_reffcol','H2Oice_reffcol','um kg m^-2',2,reffcol(1,2)) ! 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(ngridmx,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf', & ! 'kg m^-2',2,qsurf(1,iq) ) call writediagfi(ngridmx,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf', & 'kg m^-2',2,qsurf_hist(1,iq) ) call writediagfi(ngridmx,trim(noms(iq))//'_col',trim(noms(iq))//'_col', & 'kg m^-2',2,qcol(1,iq) ) if(water)then call writediagfi(ngridmx,"H2Omaxcol","max. poss. H2O column","kg m^-2",2,H2Omaxcol) endif if(watercond)then !call writediagfi(ngrid,"CLF","H2O cloud fraction"," ",3,cloudfrac) call writediagfi(ngrid,"CLFt","H2O column cloud fraction"," ",2,totcloudfrac) endif if(waterrain)then call writediagfi(ngridmx,"rain","rainfall","kg m-2 s-1",2,zdqsrain) call writediagfi(ngridmx,"snow","snowfall","kg m-2 s-1",2,zdqssnow) endif if(hydrology)then call writediagfi(ngridmx,"hice","oceanic ice height","m",2,hice) endif if(ice_update)then call writediagfi(ngridmx,"ice_min","min annual ice","m",2,ice_min) call writediagfi(ngridmx,"ice_ini","initial annual ice","m",2,ice_initial) endif do ig=1,ngrid if(tau_col(ig).gt.1.e3)then print*,'WARNING: tau_col=',tau_col(ig) print*,'at ig=',ig,'in PHYSIQ' !call abort endif end do call writediagfi(ngridmx,"tau_col","Total aerosol optical depth","[]",2,tau_col) enddo endif ! ---------------------------------------------------------- ! Output in netcdf file "diagsoil.nc" for subterranean ! variables (output every "ecritphy", as for writediagfi) ! ---------------------------------------------------------- ! Write soil temperature !call writediagsoil(ngrid,"soiltemp","Soil temperature","K",3,tsoil) else ! if(ngrid.eq.1) ! ---------------------------------------------------------------------- ! Output in grads file "g1d" (ONLY when using testphys1d) ! ---------------------------------------------------------------------- if(countG1D.eq.saveG1D)then call writeg1d(ngrid,1,tsurf,'tsurf','K') call writeg1d(ngrid,1,ps,'ps','Pa') call writeg1d(ngrid,nlayer,zt,'T','K') call writeg1d(ngrid,nlayer,pq(1,1,1),'q','kg/kg') call writeg1d(ngrid,nlayer,aerosol,'aerosol','SI') call writeg1d(ngrid,nlayer,reffrad,'reffrad','SI') call writeg1d(ngrid,nlayer,zdtlw,'dtlw','SI') call writeg1d(ngrid,nlayer,zdtsw,'dtsw','SI') call writeg1d(ngrid,nlayer,zdtdyn,'dtdyn','SI') ! radiation balance call writeg1d(ngrid,1,ISR/totarea,'ISR','W m-2') call writeg1d(ngrid,1,ASR/totarea,'ASR','W m-2') call writeg1d(ngrid,1,OLR/totarea,'OLR','W m-2') if(tracer) then do iq=1,nq call writeg1d(ngrid,1,qsurf(1,iq),trim(noms(iq))//'_s','kg m^-2') call writeg1d(ngrid,1,qcol(1,iq),trim(noms(iq))//'_c','kg m^-2') call writeg1d(ngrid,nlayer,zq(1,1,iq),noms(iq),'kg/kg') end do if(waterrain)then CALL writeg1d(ngrid,1,zdqsrain,'rainfall','kg m-2 s-1') CALL writeg1d(ngrid,1,zdqssnow,'snowfall','kg m-2 s-1') endif end if countG1D=1 else countG1D=countG1D+1 endif ! if time to save endif ! if(ngrid.ne.1) icount=icount+1 return end subroutine physiq