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,L_NSPECTI,L_NSPECTV use watercommon_h, only : mx_eau_sol, To, RLVTT, mH2O use gases_h 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 PROBLEMS WITH ALLOCATED ARRAYS !!! (temporary solution in callcorrk: do not deallocate if CLFvarying ...) clearsky=.true. call callcorrk(ngrid,nlayer,pq,nq,qsurf, & albedo,emis,mu0,pplev,pplay,pt, & tsurf,fract,dist_star,aerosol,cpp3D,muvar, & zdtlw1,zdtsw1,fluxsurf_lw1,fluxsurf_sw1,fluxtop_lw1, & fluxabs_sw1,fluxtop_dn,OLR_nu1,OSR_nu1, & reffrad,tau_col1,cloudfrac,totcloudfrac, & clearsky,firstcall,lastcall) clearsky = .false. !! just in case. ! Sum the fluxes and heating rates from cloudy/clear cases do ig=1,ngrid tf=totcloudfrac(ig) ntf=1.-tf fluxsurf_lw(ig) = ntf*fluxsurf_lw1(ig) + tf*fluxsurf_lw(ig) fluxsurf_sw(ig) = ntf*fluxsurf_sw1(ig) + tf*fluxsurf_sw(ig) fluxtop_lw(ig) = ntf*fluxtop_lw1(ig) + tf*fluxtop_lw(ig) fluxabs_sw(ig) = ntf*fluxabs_sw1(ig) + tf*fluxabs_sw(ig) tau_col(ig) = ntf*tau_col1(ig) + tf*tau_col(ig) do l=1,nlayer zdtlw(ig,l) = ntf*zdtlw1(ig,l) + tf*zdtlw(ig,l) zdtsw(ig,l) = ntf*zdtsw1(ig,l) + tf*zdtsw(ig,l) enddo do nw=1,L_NSPECTV OSR_nu(ig,nw) = ntf*OSR_nu1(ig,nw) + tf*OSR_nu(ig,nw) enddo do nw=1,L_NSPECTI OLR_nu(ig,nw) = ntf*OLR_nu1(ig,nw) + tf*OLR_nu(ig,nw) enddo enddo endif !CLFvarying ! Radiative flux from the sky absorbed by the surface (W.m-2) GSR=0.0 do ig=1,ngrid fluxrad_sky(ig)=emis(ig)*fluxsurf_lw(ig) & +fluxsurf_sw(ig)*(1.-albedo(ig)) if(noradsurf)then ! no lower surface; SW flux just disappears GSR = GSR + fluxsurf_sw(ig)*area(ig) fluxrad_sky(ig)=emis(ig)*fluxsurf_lw(ig) endif enddo if(noradsurf)then print*,'SW lost in deep atmosphere = ',GSR/totarea,' W m^-2' endif ! Net atmospheric radiative heating rate (K.s-1) do l=1,nlayer do ig=1,ngrid dtrad(ig,l)=zdtsw(ig,l)+zdtlw(ig,l) enddo enddo elseif(newtonian)then ! b) Call Newtonian cooling scheme ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ call newtrelax(mu0,sinlat,zpopsk,pt,pplay,pplev,dtrad,firstcall) do ig=1,ngrid zdtsurf(ig) = +(pt(ig,1)-tsurf(ig))/ptimestep enddo ! e.g. surface becomes proxy for 1st atmospheric layer ? else ! c) Atmosphere has no radiative effect ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ do ig=1,ngrid fluxtop_dn(ig) = fract(ig)*mu0(ig)*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(ig) = fluxtop_dn(ig) fluxrad_sky(ig) = fluxtop_dn(ig)*(1.-albedo(ig)) fluxtop_lw(ig) = emis(ig)*stephan*tsurf(ig)**4 enddo ! radiation skips the atmosphere entirely do l=1,nlayer do ig=1,ngrid dtrad(ig,l)=0.0 enddo enddo ! hence no atmospheric radiative heating endif ! if corrk endif ! of if(mod(icount-1,iradia).eq.0) ! Transformation of the radiative tendencies ! ------------------------------------------ do ig=1,ngrid zplanck(ig)=tsurf(ig)*tsurf(ig) zplanck(ig)=emis(ig)*stephan*zplanck(ig)*zplanck(ig) fluxrad(ig)=fluxrad_sky(ig)-zplanck(ig) enddo do l=1,nlayer do ig=1,ngrid pdt(ig,l)=pdt(ig,l)+dtrad(ig,l) enddo enddo !------------------------- ! test energy conservation if(enertest)then dEtotSW = 0.0 dEtotLW = 0.0 dEtotsSW = 0.0 dEtotsLW = 0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dEtotSW = dEtotSW + cpp3D(ig,l)*masse*zdtsw(ig,l)*area(ig) dEtotLW = dEtotLW + cpp3D(ig,l)*masse*zdtlw(ig,l)*area(ig) enddo dEtotsSW = dEtotsSW + fluxsurf_sw(ig)*(1.-albedo(ig))*area(ig) dEtotsLW = dEtotsLW + emis(ig)*fluxsurf_lw(ig)*area(ig)-zplanck(ig)*area(ig) enddo dEtotSW = dEtotSW/totarea dEtotLW = dEtotLW/totarea dEtotsSW = dEtotsSW/totarea dEtotsLW = dEtotsLW/totarea print*,'---------------------------------------------------------------' print*,'In corrk SW atmospheric energy change =',dEtotSW,' W m-2' print*,'In corrk LW atmospheric energy change =',dEtotLW,' W m-2' print*,'atmospheric energy change (SW+LW) =',dEtotLW+dEtotSW,' W m-2' print*,'In corrk SW surface energy change =',dEtotsSW,' W m-2' print*,'In corrk LW surface energy change =',dEtotsLW,' W m-2' print*,'surface energy change (SW+LW) =',dEtotsLW+dEtotsSW,' W m-2' endif !------------------------- endif ! of if (callrad) !----------------------------------------------------------------------- ! 4. Vertical diffusion (turbulent mixing): ! ----------------------------------------- if (calldifv) then do ig=1,ngrid zflubid(ig)=fluxrad(ig)+fluxgrd(ig) enddo zdum1(:,:)=0.0 zdum2(:,:)=0.0 do l=1,nlayer do ig=1,ngrid zdh(ig,l)=pdt(ig,l)/zpopsk(ig,l) enddo enddo call vdifc(ngrid,nlayer,nq,rnat,zpopsk, & ptimestep,capcal,lwrite, & pplay,pplev,zzlay,zzlev,z0, & pu,pv,zh,pq,tsurf,emis,qsurf, & zdum1,zdum2,zdh,pdq,zflubid, & zdudif,zdvdif,zdhdif,zdtsdif,q2, & zdqdif,zdqsdif,lastcall) do l=1,nlayer do ig=1,ngrid pdv(ig,l)=pdv(ig,l)+zdvdif(ig,l) pdu(ig,l)=pdu(ig,l)+zdudif(ig,l) pdt(ig,l)=pdt(ig,l)+zdhdif(ig,l)*zpopsk(ig,l) zdtdif(ig,l)=zdhdif(ig,l)*zpopsk(ig,l) ! for diagnostic only enddo enddo do ig=1,ngrid zdtsurf(ig)=zdtsurf(ig)+zdtsdif(ig) enddo if (tracer) then do iq=1, nq do l=1,nlayer do ig=1,ngrid pdq(ig,l,iq)=pdq(ig,l,iq)+ zdqdif(ig,l,iq) enddo enddo enddo do iq=1, nq do ig=1,ngrid dqsurf(ig,iq)=dqsurf(ig,iq) + zdqsdif(ig,iq) enddo enddo end if ! of if (tracer) !------------------------- ! test energy conservation if(enertest)then dEtot=0.0 dEtots=0.0 vdifcdE(:)=0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dEtot = dEtot + cpp3D(ig,l)*masse*zdtdif(ig,l)*area(ig) vabs = sqrt(pdu(ig,l)**2 + pdv(ig,l)**2) dvabs = sqrt(zdudif(ig,l)**2 + zdvdif(ig,l)**2) dEtot = dEtot + masse*vabs*dvabs*area(ig) vdifcdE(ig) = vdifcdE(ig) + masse*vabs*dvabs enddo dEtot = dEtot - zflubid(ig)*area(ig) ! subtract flux from ground dEtots = dEtots + capcal(ig)*zdtsdif(ig)*area(ig) vdifcdE(ig) = vdifcdE(ig) + capcal(ig)*zdtsdif(ig) enddo dEtot = dEtot/totarea dEtots = dEtots/totarea print*,'In difv atmospheric energy change =',dEtot,' W m-2' print*,'In difv surface energy change =',dEtots,' W m-2' !print*,'Note we ignore the wind change...' ! and latent heat for that matter... endif !------------------------- !------------------------- ! test water conservation if(watertest.and.water)then dWtot=0.0 dWtots=0.0 nconsMAX=0.0 do ig = 1, ngrid vdifcncons(ig)=0.0 do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g iq = igcm_h2o_vap dWtot = dWtot + masse*zdqdif(ig,l,iq)*ptimestep*area(ig) vdifcncons(ig)=vdifcncons(ig) + masse*zdqdif(ig,l,iq) iq = igcm_h2o_ice dWtot = dWtot + masse*zdqdif(ig,l,iq)*ptimestep*area(ig) vdifcncons(ig)=vdifcncons(ig) + masse*zdqdif(ig,l,iq) enddo iq = igcm_h2o_vap dWtots = dWtots + zdqsdif(ig,iq)*ptimestep*area(ig) vdifcncons(ig)=vdifcncons(ig)+zdqsdif(ig,iq) iq = igcm_h2o_ice dWtots = dWtots + zdqsdif(ig,iq)*ptimestep*area(ig) vdifcncons(ig)=vdifcncons(ig)+zdqsdif(ig,iq) if(vdifcncons(ig).gt.nconsMAX)then nconsMAX=vdifcncons(ig) endif enddo dWtot = dWtot/totarea dWtots = dWtots/totarea print*,'---------------------------------------------------------------' print*,'In difv atmospheric water change =',dWtot,' kg m-2' print*,'In difv surface water change =',dWtots,' kg m-2' print*,'In difv non-cons factor =',dWtot+dWtots,' kg m-2' print*,'In difv MAX non-cons factor =',nconsMAX,' kg m-2 s-1' endif !------------------------- else if(.not.newtonian)then do ig=1,ngrid zdtsurf(ig) = zdtsurf(ig) + (fluxrad(ig) + fluxgrd(ig))/capcal(ig) enddo endif endif ! of if (calldifv) !----------------------------------------------------------------------- ! 5. Dry convective adjustment: ! ----------------------------- if(calladj) then do l=1,nlayer do ig=1,ngrid if(nonideal)then zdh(ig,l) = pdt(ig,l)/zpopskNI(ig,l) else zdh(ig,l) = pdt(ig,l)/zpopsk(ig,l) endif enddo enddo zduadj(:,:)=0.0 zdvadj(:,:)=0.0 zdhadj(:,:)=0.0 if(nonideal)then print*,'Nonideal is not used at the moment' call abort call convadj(ngrid,nlayer,nq,ptimestep, & pplay,pplev,zpopskNI, & pu,pv,zh,pq, & pdu,pdv,zdh,pdq, & zduadj,zdvadj,zdhadj, & zdqadj) else call convadj(ngrid,nlayer,nq,ptimestep, & pplay,pplev,zpopsk, & pu,pv,zh,pq, & pdu,pdv,zdh,pdq, & zduadj,zdvadj,zdhadj, & zdqadj) endif do l=1,nlayer do ig=1,ngrid pdu(ig,l) = pdu(ig,l) + zduadj(ig,l) pdv(ig,l) = pdv(ig,l) + zdvadj(ig,l) if(nonideal)then pdt(ig,l) = pdt(ig,l) + zdhadj(ig,l)*zpopskNI(ig,l) zdtadj(ig,l) = zdhadj(ig,l)*zpopskNI(ig,l) ! for diagnostic only else pdt(ig,l) = pdt(ig,l) + zdhadj(ig,l)*zpopsk(ig,l) zdtadj(ig,l) = zdhadj(ig,l)*zpopsk(ig,l) ! for diagnostic only endif enddo enddo if(tracer) then do iq=1, nq do l=1,nlayer do ig=1,ngrid pdq(ig,l,iq) = pdq(ig,l,iq) + zdqadj(ig,l,iq) enddo enddo enddo end if !------------------------- ! test energy conservation if(enertest)then dEtot=0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dEtot = dEtot + cpp3D(ig,l)*masse*zdtadj(ig,l)*area(ig) enddo enddo dEtot=dEtot/totarea print*,'In convadj atmospheric energy change =',dEtot,' W m-2' endif !------------------------- !------------------------- ! test water conservation if(watertest)then dWtot=0.0 do ig = 1, ngrid cadjncons(ig)=0.0 do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g iq = igcm_h2o_vap dWtot = dWtot + masse*zdqadj(ig,l,iq)*ptimestep*area(ig) cadjncons(ig)=cadjncons(ig) + masse*zdqadj(ig,l,iq) iq = igcm_h2o_ice dWtot = dWtot + masse*zdqadj(ig,l,iq)*ptimestep*area(ig) cadjncons(ig)=cadjncons(ig) + masse*zdqadj(ig,l,iq) enddo enddo dWtot=dWtot/totarea print*,'In convadj atmospheric water change =',dWtot,' kg m-2' endif !------------------------- endif ! of if(calladj) !----------------------------------------------------------------------- ! 6. Carbon dioxide condensation-sublimation: ! ------------------------------------------- if (co2cond) then if (.not.tracer) then print*,'We need a CO2 ice tracer to condense CO2' call abort endif call condense_cloud(ngrid,nlayer,nq,ptimestep, & capcal,pplay,pplev,tsurf,pt, & pphi,pdt,pdu,pdv,zdtsurf,pu,pv,pq,pdq, & qsurf(1,igcm_co2_ice),albedo,emis, & zdtc,zdtsurfc,pdpsrf,zduc,zdvc, & zdqc,reffrad,cpp3D) do l=1,nlayer do ig=1,ngrid pdt(ig,l)=pdt(ig,l)+zdtc(ig,l) pdv(ig,l)=pdv(ig,l)+zdvc(ig,l) pdu(ig,l)=pdu(ig,l)+zduc(ig,l) enddo enddo do ig=1,ngrid zdtsurf(ig) = zdtsurf(ig) + zdtsurfc(ig) enddo do iq=1,nq ! should use new notation here ! do l=1,nlayer do ig=1,ngrid pdq(ig,l,iq)=pdq(ig,l,iq)+ zdqc(ig,l,iq) enddo enddo enddo ! Note: we do not add surface co2ice tendency ! because qsurf(:,igcm_co2_ice) is updated in condens_co2cloud !------------------------- ! test energy conservation if(enertest)then dEtot=0.0 dEtots=0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dEtot = dEtot + cpp3D(ig,l)*masse*zdtc(ig,l)*area(ig) enddo dEtots = dEtots + capcal(ig)*zdtsurfc(ig)*area(ig) enddo dEtot=dEtot/totarea dEtots=dEtots/totarea print*,'In co2cloud atmospheric energy change =',dEtot,' W m-2' print*,'In co2cloud surface energy change =',dEtots,' W m-2' endif !------------------------- endif ! of if (co2cond) !----------------------------------------------------------------------- ! 7. Specific parameterizations for tracers ! ----------------------------------------- if (tracer) then ! 7a. Water and ice ! --------------- if (water) then ! ---------------------------------------- ! Water ice condensation in the atmosphere ! ---------------------------------------- if(watercond)then if(RLVTT.gt.1.e-8)then ! Re-evaporate cloud water/ice call evap(ptimestep,pt,pq,pdq,pdt,dqevap,dtevap,qevap,tevap) DO l = 1, nlayer DO ig = 1, ngrid pdq(ig,l,igcm_h2o_vap) = pdq(ig,l,igcm_h2o_vap)+dqevap(ig,l) pdq(ig,l,igcm_h2o_ice) = pdq(ig,l,igcm_h2o_ice)-dqevap(ig,l) pdt(ig,l) = pdt(ig,l)+dtevap(ig,l) enddo enddo call moistadj(pt,qevap,pplev,pplay,dtmoist,dqmoist,ptimestep,rneb_man) do l=1,nlayer do ig=1,ngrid pdq(ig,l,igcm_h2o_vap) = pdq(ig,l,igcm_h2o_vap)+dqmoist(ig,l,igcm_h2o_vap) pdq(ig,l,igcm_h2o_ice) = pdq(ig,l,igcm_h2o_ice)+dqmoist(ig,l,igcm_h2o_ice) pdt(ig,l) = pdt(ig,l)+dtmoist(ig,l) enddo enddo !------------------------- ! test energy conservation if(enertest)then dEtot=0.0 madjdE(:)=0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dEtot = dEtot + cpp3D(ig,l)*masse*(dtmoist(ig,l)+dtevap(ig,l))*area(ig) madjdE(ig) = madjdE(ig) + cpp3D(ig,l)*masse*(dtmoist(ig,l)+dtevap(ig,l)) enddo enddo dEtot=dEtot/totarea print*,'In moistadj atmospheric energy change =',dEtot,' W m-2' endif !------------------------- !------------------------- ! test water conservation if(watertest)then dWtot=0.0 do ig = 1, ngrid do iq = 1 , nq do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dWtot = dWtot + masse*dqmoist(ig,l,igcm_h2o_vap)*area(ig)*ptimestep dWtot = dWtot + masse*dqmoist(ig,l,igcm_h2o_ice)*area(ig)*ptimestep enddo enddo enddo dWtot=dWtot/totarea print*,'In moistadj atmospheric water change =',dWtot,' kg m-2' endif !------------------------- endif ! Re-evaporate cloud water/ice call evap(ptimestep,pt,pq,pdq,pdt,dqevap,dtevap,qevap,tevap) do l = 1, nlayer do ig = 1, ngrid pdq(ig,l,igcm_h2o_vap) = pdq(ig,l,igcm_h2o_vap)+dqevap(ig,l) pdq(ig,l,igcm_h2o_ice) = pdq(ig,l,igcm_h2o_ice)-dqevap(ig,l) pdt(ig,l) = pdt(ig,l)+dtevap(ig,l) enddo enddo ! note: we use qevap but not tevap in largescale/moistadj ! otherwise is a big mess call largescale(ptimestep,pplev,pplay,pt,qevap, & ! a bug was here! pdt,dtlscale,dqvaplscale,dqcldlscale,rneb_lsc,reffH2O) do l=1,nlayer do ig=1,ngrid pdq(ig,l,igcm_h2o_vap) = pdq(ig,l,igcm_h2o_vap)+dqvaplscale(ig,l) pdq(ig,l,igcm_h2o_ice) = pdq(ig,l,igcm_h2o_ice)+dqcldlscale(ig,l) pdt(ig,l) = pdt(ig,l)+dtlscale(ig,l) if(.not.aerofixed)then reffrad(ig,l,2)=reffH2O(ig,l) endif enddo enddo !------------------------- ! test energy conservation if(enertest)then dEtot=0.0 lscaledE(:)=0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dEtot = dEtot + cpp3D(ig,l)*masse*(dtlscale(ig,l)+dtevap(ig,l))*area(ig) lscaledE(ig) = lscaledE(ig) + cpp3D(ig,l)*masse*(dtlscale(ig,l)+dtevap(ig,l)) enddo enddo dEtot=dEtot/totarea print*,'In largescale atmospheric energy change =',dEtot,' W m-2' endif !------------------------- !------------------------- ! test water conservation if(watertest)then dWtot=0.0 do ig = 1, ngrid do iq = 1 , nq do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dWtot = dWtot + masse*dqvaplscale(ig,l)*area(ig)*ptimestep dWtot = dWtot + masse*dqcldlscale(ig,l)*area(ig)*ptimestep enddo enddo enddo dWtot=dWtot/totarea print*,'In largescale atmospheric water change =',dWtot,' kg m-2' endif !------------------------- ! compute cloud fraction do l = 1, nlayer do ig = 1,ngrid cloudfrac(ig,l)=MAX(rneb_lsc(ig,l),rneb_man(ig,l)) enddo enddo ! compute total cloud fraction in column call totalcloudfrac(cloudfrac,totcloudfrac) endif ! of if (watercondense) ! -------------------------------- ! Water ice / liquid precipitation ! -------------------------------- if(waterrain)then zdqrain(:,:,:) = 0.0 zdqsrain(:) = 0.0 zdqssnow(:) = 0.0 call rain(ptimestep,pplev,pplay,pt,pdt,pq,pdq, & zdtrain,zdqrain,zdqsrain,zdqssnow,cloudfrac) do l=1,nlayer do ig=1,ngrid pdq(ig,l,igcm_h2o_vap) = pdq(ig,l,igcm_h2o_vap)+zdqrain(ig,l,igcm_h2o_vap) pdq(ig,l,igcm_h2o_ice) = pdq(ig,l,igcm_h2o_ice)+zdqrain(ig,l,igcm_h2o_ice) pdt(ig,l) = pdt(ig,l)+zdtrain(ig,l) enddo enddo do ig=1,ngrid dqsurf(ig,igcm_h2o_vap) = dqsurf(ig,igcm_h2o_vap)+zdqsrain(ig) ! a bug was here dqsurf(ig,igcm_h2o_ice) = dqsurf(ig,igcm_h2o_ice)+zdqssnow(ig) rainout(ig) = zdqsrain(ig)+zdqssnow(ig) ! diagnostic enddo !------------------------- ! test energy conservation if(enertest)then dEtot=0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dEtot = dEtot + cpp3D(ig,l)*masse*zdtrain(ig,l)*area(ig) enddo enddo dEtot=dEtot/totarea print*,'In rain atmospheric energy change =',dEtot,' W m-2' endif !------------------------- !------------------------- ! test energy conservation if(enertest)then dEtot=0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dEtot = dEtot + cpp3D(ig,l)*masse*zdtrain(ig,l)*area(ig) enddo enddo dEtot=dEtot/totarea print*,'In rain atmospheric T energy change =',dEtot,' W m-2' dEtot=0.0 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g dItot = dItot + masse*zdqrain(ig,l,igcm_h2o_ice)*area(ig) dVtot = dVtot + masse*zdqrain(ig,l,igcm_h2o_vap)*area(ig) enddo dItot = dItot + zdqssnow(ig)*area(ig) dVtot = dVtot + zdqsrain(ig)*area(ig) enddo dEtot=(dItot*RLVTT/cpp + dVtot*RLVTT/cpp)/totarea print*,'In rain dItot =',dItot*RLVTT/(cpp*totarea),' W m-2' print*,'In rain dVtot =',dVtot*RLVTT/(cpp*totarea),' W m-2' print*,'In rain atmospheric L energy change =',dEtot,' W m-2' endif !------------------------- !------------------------- ! test water conservation if(watertest)then dWtot=0.0 dWtots=0.0 do ig = 1, ngrid !do iq = 1 , nq do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g ! equiv to l2c in rain dWtot = dWtot + masse*zdqrain(ig,l,igcm_h2o_vap)*area(ig)*ptimestep dWtot = dWtot + masse*zdqrain(ig,l,igcm_h2o_ice)*area(ig)*ptimestep enddo !enddo dWtots = dWtots + (zdqsrain(ig)+zdqssnow(ig))*area(ig)*ptimestep enddo dWtot=dWtot/totarea dWtots=dWtots/totarea print*,'In rain atmospheric water change =',dWtot,' kg m-2' print*,'In rain surface water change =',dWtots,' kg m-2' print*,'In rain non-cons factor =',dWtot+dWtots,' kg m-2' endif !------------------------- end if ! of if (waterrain) end if ! of if (water) ! 7c. Aerosol particles ! ------------------- ! ------------- ! Sedimentation ! ------------- if (sedimentation) then zdqsed(:,:,:) = 0.0 zdqssed(:,:) = 0.0 !------------------------- ! find qtot if(watertest)then dWtot=0.0 dWtots=0.0 iq=3 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g ! equiv to l2c in rain dWtot = dWtot + masse*pq(ig,l,iq)*area(ig)*ptimestep dWtots = dWtots + masse*pdq(ig,l,iq)*area(ig)*ptimestep enddo enddo dWtot=dWtot/totarea dWtots=dWtots/totarea print*,'Before sedim pq =',dWtot,' kg m-2' print*,'Before sedim pdq =',dWtots,' kg m-2' endif !------------------------- call callsedim(ngrid,nlayer,ptimestep, & pplev,zzlev,pt,pq,pdq,zdqsed,zdqssed,nq,reffH2O) !------------------------- ! find qtot if(watertest)then dWtot=0.0 dWtots=0.0 iq=3 do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g ! equiv to l2c in rain dWtot = dWtot + masse*pq(ig,l,iq)*area(ig)*ptimestep dWtots = dWtots + masse*pdq(ig,l,iq)*area(ig)*ptimestep enddo enddo dWtot=dWtot/totarea dWtots=dWtots/totarea print*,'After sedim pq =',dWtot,' kg m-2' print*,'After sedim pdq =',dWtots,' kg m-2' endif !------------------------- do iq=1,nq ! for now, we only allow H2O ice to sediment ! and as in rain.F90, whether it falls as rain or snow depends ! only on the surface temperature do ig=1,ngrid do l=1,nlayer pdq(ig,l,iq) = pdq(ig,l,iq) + zdqsed(ig,l,iq) enddo dqsurf(ig,iq) = dqsurf(ig,iq) + zdqssed(ig,iq) enddo enddo !------------------------- ! test water conservation if(watertest)then dWtot=0.0 dWtots=0.0 do iq=1,nq do ig = 1, ngrid do l = 1, nlayer masse = (pplev(ig,l) - pplev(ig,l+1))/g ! equiv to l2c in rain dWtot = dWtot + masse*zdqsed(ig,l,iq)*area(ig)*ptimestep enddo dWtots = dWtots + zdqssed(ig,iq)*area(ig)*ptimestep enddo enddo dWtot=dWtot/totarea dWtots=dWtots/totarea print*,'In sedim atmospheric ice change =',dWtot,' kg m-2' print*,'In sedim surface ice change =',dWtots,' kg m-2' print*,'In sedim non-cons factor =',dWtot+dWtots,' kg m-2' endif !------------------------- end if ! of if (sedimentation) ! 7d. Updates ! --------- ! --------------------------------- ! Updating tracer budget on surface ! --------------------------------- if(hydrology)then call hydrol(ptimestep,rnat,tsurf,qsurf,dqsurf,dqs_hyd, & capcal,albedo0,albedo,mu0,zdtsurf,zdtsurf_hyd,hice) ! note: for now, also changes albedo in the subroutine do ig=1,ngrid zdtsurf(ig) = zdtsurf(ig) + zdtsurf_hyd(ig) do iq=1,nq qsurf(ig,iq) = qsurf(ig,iq)+ptimestep*dqs_hyd(ig,iq) enddo enddo ! when hydrology is used, other dqsurf tendencies are all added to dqs_hyd inside !------------------------- ! test energy conservation if(enertest)then dEtot=0.0 do ig = 1, ngrid dEtots = dEtots + capcal(ig)*zdtsurf_hyd(ig)*area(ig) enddo dEtot=dEtot/totarea print*,'In hydrol atmospheric energy change =',dEtot,' W m-2' endif !------------------------- !------------------------- ! test water conservation if(watertest)then dWtots=0.0 do ig = 1, ngrid dWtots = dWtots + dqs_hyd(ig,igcm_h2o_ice)*area(ig)*ptimestep enddo dWtots=dWtots/totarea print*,'In hydrol surface ice change =',dWtots,' kg m-2' dWtots=0.0 do ig = 1, ngrid dWtots = dWtots + dqs_hyd(ig,igcm_h2o_vap)*area(ig)*ptimestep enddo dWtots=dWtots/totarea print*,'In hydrol surface water change =',dWtots,' kg m-2' print*,'---------------------------------------------------------------' endif !------------------------- ELSE ! of if (hydrology) do iq=1,nq do ig=1,ngrid qsurf(ig,iq)=qsurf(ig,iq)+ptimestep*dqsurf(ig,iq) enddo enddo END IF ! of if (hydrology) ! Add qsurf to qsurf_hist, which is what we save in ! diagfi.nc etc. At the same time, we set the water ! content of ocean gridpoints back to zero, in order ! to avoid rounding errors in vdifc, rain do ig = 1, ngrid do iq = 1, nq if(iq.eq.igcm_h2o_vap .and. rnat(ig).eq.0)then ! if liquid water and terrain = ocean qsurf_hist(ig,iq) = qsurf(ig,iq) !qsurf(ig,iq) = qcol(ig,iq) ! the value of qsurf we choose here makes NO DIFFERENCE TO ANYTHING AT ALL else qsurf_hist(ig,iq) = qsurf(ig,iq) endif enddo enddo if(ice_update)then do ig = 1, ngrid ice_min(ig)=min(ice_min(ig),qsurf(ig,igcm_h2o_ice)) enddo endif endif ! of if (tracer) !----------------------------------------------------------------------- ! 9. Surface and sub-surface soil temperature !----------------------------------------------------------------------- ! Increment surface temperature do ig=1,ngrid tsurf(ig)=tsurf(ig)+ptimestep*zdtsurf(ig) enddo ! Compute soil temperatures and subsurface heat flux if (callsoil) then call soil(ngrid,nsoilmx,.false.,inertiedat, & ptimestep,tsurf,tsoil,capcal,fluxgrd) endif !------------------------- ! test energy conservation if(enertest)then dEtots=0.0 do ig = 1, ngrid dEtots = dEtots + capcal(ig)*zdtsurf(ig)*area(ig) enddo dEtots=dEtots/totarea print*,'Surface energy change=',dEtots,' W m-2' endif !------------------------- !----------------------------------------------------------------------- ! 10. Perform diagnostics and write output files !----------------------------------------------------------------------- ! ------------------------------- ! Dynamical fields incrementation ! ------------------------------- ! For output only: the actual model integration is performed in the dynamics ! temperature, zonal and meridional wind do l=1,nlayer do ig=1,ngrid zt(ig,l) = pt(ig,l) + pdt(ig,l)*ptimestep zu(ig,l) = pu(ig,l) + pdu(ig,l)*ptimestep zv(ig,l) = pv(ig,l) + pdv(ig,l)*ptimestep ! diagnostic zdtdyn(ig,l) = ztprevious(ig,l)-pt(ig,l) ztprevious(ig,l) = zt(ig,l) enddo enddo if(firstcall)then zdtdyn(:,:)=0.0 endif ! dynamical heating diagnostic fluxdyn(:)=0. do ig=1,ngrid do l=1,nlayer fluxdyn(ig)=fluxdyn(ig) - (zdtdyn(ig,l)/ptimestep) & *(pplev(ig,l)-pplev(ig,l+1))*cpp3D(ig,l)/g enddo enddo ! tracers do iq=1, nq do l=1,nlayer do ig=1,ngrid zq(ig,l,iq) = pq(ig,l,iq) + pdq(ig,l,iq)*ptimestep enddo enddo enddo ! surface pressure do ig=1,ngrid ps(ig) = pplev(ig,1) + pdpsrf(ig)*ptimestep enddo ! pressure do l=1,nlayer do ig=1,ngrid zplev(ig,l) = pplev(ig,l)/pplev(ig,1)*ps(ig) zplay(ig,l) = pplay(ig,l)/pplev(ig,1)*ps(ig) enddo enddo ! --------------------------------------------------------- ! Surface and soil temperature information ! --------------------------------------------------------- Ts1 = 0.0 Ts2 = 99999.9 Ts3 = 0.0 TsS = 0.0 ! mean temperature at bottom soil layer do ig=1,ngrid Ts1 = Ts1 + area(ig)*tsurf(ig) Ts2 = min(Ts2,tsurf(ig)) Ts3 = max(Ts3,tsurf(ig)) TsS = TsS + area(ig)*tsoil(ig,nsoilmx) end do Ts1=Ts1/totarea TsS=TsS/totarea if(callsoil)then print*,' ave[Tsurf] min[Tsurf] max[Tsurf] ave[Tdeep]' print*,Ts1,Ts2,Ts3,TsS else print*,' ave[Tsurf] min[Tsurf] max[Tsurf]' print*,Ts1,Ts2,Ts3 endif ! --------------------------------------------------------- ! Check the energy balance of the simulation during the run ! --------------------------------------------------------- if(corrk)then ISR = 0.0 ASR = 0.0 OLR = 0.0 GND = 0.0 DYN = 0.0 do ig=1,ngrid ISR = ISR + area(ig)*fluxtop_dn(ig) ASR = ASR + area(ig)*fluxabs_sw(ig) OLR = OLR + area(ig)*fluxtop_lw(ig) GND = GND + area(ig)*fluxgrd(ig) if(.not.callsoil) GND=GND+ area(ig)*fluxrad(ig) DYN = DYN + area(ig)*fluxdyn(ig) 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(ngridmx.eq.1)then DYN=0.0 endif print*,' ISR ASR OLR GND DYN [W m^-2]' print*, ISR/totarea,ASR/totarea,OLR/totarea,GND/totarea,DYN/totarea if(enertest)then print*,'SW energy balance SW++ - ASR = ',dEtotSW+dEtotsSW-ASR/totarea,' W m-2' print*,'LW energy balance LW++ + ***ASR*** = ',dEtotLW+dEtotsLW+ASR/totarea,' W m-2' endif if(meanOLR)then if((ngridmx.gt.1) .or. (mod(icount-1,nint(ecritphy)).eq.0))then ! to record global radiative balance open(92,file="rad_bal.out",form='formatted',access='append') write(92,*) zday,ISR/totarea,ASR/totarea,OLR/totarea close(92) open(93,file="tem_bal.out",form='formatted',access='append') write(93,*) zday,Ts1,Ts2,Ts3,TsS close(93) endif endif endif ! ------------------------------------------------------------------ ! Diagnostic to test radiative-convective timescales in code ! ------------------------------------------------------------------ if(testradtimes)then open(38,file="tau_phys.out",form='formatted',access='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, exiting physics on first call' call abort endif ! --------------------------------------------------------- ! Compute column amounts (kg m-2) if tracers are enabled ! --------------------------------------------------------- if(tracer)then qcol(:,:)=0.0 do iq=1,nq do ig=1,ngrid do l=1,nlayer qcol(ig,iq) = qcol(ig,iq) + zq(ig,l,iq) * & (pplev(ig,l) - pplev(ig,l+1)) / g enddo enddo enddo ! not generalised for arbitrary aerosols yet!!! reffcol(:,:)=0.0 do ig=1,ngrid do l=1,nlayer if(co2cond)then reffcol(ig,1) = reffcol(ig,1) + zq(ig,l,igcm_co2_ice) * & reffrad(ig,l,1) * & (pplev(ig,l) - pplev(ig,l+1)) / g endif if(water)then reffcol(ig,2) = reffcol(ig,2) + zq(ig,l,igcm_h2o_ice) * & reffrad(ig,l,2) * & (pplev(ig,l) - pplev(ig,l+1)) / g endif enddo enddo endif ! --------------------------------------------------------- ! Test for water conservation if water is enabled ! --------------------------------------------------------- if(water)then icesrf = 0.0 liqsrf = 0.0 icecol = 0.0 vapcol = 0.0 h2otot = 0.0 do ig=1,ngrid icesrf = icesrf + area(ig)*qsurf_hist(ig,igcm_h2o_ice) liqsrf = liqsrf + area(ig)*qsurf_hist(ig,igcm_h2o_vap) icecol = icecol + area(ig)*qcol(ig,igcm_h2o_ice) vapcol = vapcol + area(ig)*qcol(ig,igcm_h2o_vap) h2otot = h2otot + area(ig)* & (qcol(ig,igcm_h2o_ice)+qcol(ig,igcm_h2o_vap) & +qsurf_hist(ig,igcm_h2o_ice)+qsurf_hist(ig,igcm_h2o_vap)) end do print*,' Total water amount [kg m^-2]: ',h2otot/totarea print*,' Surface ice Surface liq. Atmos. con. Atmos. vap. [kg m^-2] ' print*, icesrf/totarea,liqsrf/totarea,icecol/totarea,vapcol/totarea if(meanOLR)then if((ngridmx.gt.1) .or. (mod(icount-1,nint(ecritphy)).eq.0))then ! to record global water balance open(98,file="h2o_bal.out",form='formatted',access='append') write(98,*) zday,icesrf/totarea,liqsrf/totarea,icecol/totarea,vapcol/totarea close(98) endif endif endif ! --------------------------------------------------------- ! Calculate RH for diagnostic if water is enabled ! --------------------------------------------------------- if(water)then do l = 1, nlayer do ig = 1, ngrid call watersat(pt(ig,l),pplay(ig,l),qsat(ig,l)) RH(ig,l) = zq(ig,l,igcm_h2o_vap) / qsat(ig,l) enddo enddo ! compute maximum possible H2O column amount (100% saturation) do ig=1,ngrid H2Omaxcol(ig)=0.0 do l=1,nlayer H2Omaxcol(ig) = H2Omaxcol(ig) + qsat(ig,l) * & (pplev(ig,l) - pplev(ig,l+1))/g enddo enddo endif print*,'' 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 if(ice_update)then do ig = 1, ngrid delta_ice = (qsurf(ig,igcm_h2o_ice)-ice_initial(ig)) ! add multiple years of evolution qsurf_hist(ig,igcm_h2o_ice) = & !qsurf_hist(ig,igcm_h2o_ice) + delta_ice*100.0 qsurf_hist(ig,igcm_h2o_ice) + delta_ice*icetstep ! if ice has gone -ve, set to zero if(qsurf_hist(ig,igcm_h2o_ice).lt.0.0)then qsurf_hist(ig,igcm_h2o_ice) = 0.0 !qsurf_hist(ig,igcm_h2o_vap) = 0.0 endif ! if ice is seasonal, set to zero (NEW) if(ice_min(ig).lt.0.01)then qsurf_hist(ig,igcm_h2o_ice) = 0.0 !qsurf_hist(ig,igcm_h2o_vap) = 0.0 endif enddo ! enforce ice conservation ice_tot=0.0 do ig = 1, ngrid ice_tot = ice_tot + qsurf_hist(ig,igcm_h2o_ice)*area(ig) enddo do ig = 1, ngrid qsurf_hist(ig,igcm_h2o_ice) = qsurf_hist(ig,igcm_h2o_ice)*(icesrf/ice_tot) 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_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,"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,"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(ngridmx,trim(noms(iq))//'_surf',trim(noms(iq))//'_surf', & 'kg m^-2',2,qsurf(1,iq) ) call wstats(ngridmx,trim(noms(iq))//'_col',trim(noms(iq))//'_col', & 'kg m^-2',2,qcol(1,iq) ) call wstats(ngridmx,trim(noms(iq))//'_reff', & trim(noms(iq))//'_reff', & 'm',3,reffrad(1,1,iq)) end do 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) 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) ! 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) ! Output aerosols if (igcm_co2_ice.ne.0) call writediagfi(ngridmx,'CO2ice_reff','CO2ice_reff','m',3,reffrad(1,1,1)) if (igcm_h2o_ice.ne.0) call writediagfi(ngridmx,'H2Oice_reff','H2Oice_reff','m',3,reffrad(1,1,2)) if (igcm_co2_ice.ne.0) call writediagfi(ngridmx,'CO2ice_reffcol','CO2ice_reffcol','um kg m^-2',2,reffcol(1,1)) if (igcm_h2o_ice.ne.0) 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 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) endif icount=icount+1 !!! DEALLOCATE STUFF if (lastcall) then IF ( ALLOCATED( gnom ) ) DEALLOCATE( gnom ) !! this was allocated in su_gases.F90 IF ( ALLOCATED( gfrac ) ) DEALLOCATE( gfrac ) !! this was allocated in su_gases.F90 endif return end subroutine physiq