! ! $Id: leapfrog.F 1446 2010-10-22 09:27:25Z emillour $ ! c c SUBROUTINE leapfrog_nogcm_pluto(ucov,vcov,teta,ps,masse,phis,q, & time_0) cIM : pour sortir les param. du modele dans un fis. netcdf 110106 #ifdef CPP_IOIPSL use IOIPSL #endif USE infotrac, ONLY: nqtot,ok_iso_verif,tname USE write_field, ONLY: writefield USE control_mod, ONLY: planet_type,nday,day_step,iperiod,iphysiq, & less1day,fractday,ndynstep,iconser, & dissip_period,offline,ip_ebil_dyn, & ok_dynzon,periodav,ok_dyn_ave,iecri, & ok_dyn_ins,output_grads_dyn use exner_hyb_m, only: exner_hyb use exner_milieu_m, only: exner_milieu use cpdet_mod, only: cpdet,tpot2t,t2tpot use sponge_mod, only: callsponge,mode_sponge,sponge use comuforc_h USE comvert_mod, ONLY: ap,bp,pressure_exner,presnivs, & aps,bps,presnivs,pseudoalt,preff,scaleheight USE comconst_mod, ONLY: daysec,dtvr,dtphys,dtdiss, . cpp,ihf,iflag_top_bound,pi,kappa,r USE logic_mod, ONLY: iflag_phys,ok_guide,forward,leapf,apphys, . statcl,conser,purmats,tidal,ok_strato USE temps_mod, ONLY: jD_ref,jH_ref,itaufin,day_ini,day_ref, . start_time,dt USE callkeys_mod USE tracer_h, ONLY: igcm_n2,igcm_ch4_gas,igcm_co_gas, . igcm_prec_haze,igcm_haze IMPLICIT NONE c ...... Version du 30/11/24 .......... c======================================================================= c c Auteur: T. Bertrand, F. Forget c ------- c c Objet: c ------ c c GCM LMD sans dynamique, pour modele saisonnier de surface c c======================================================================= c c----------------------------------------------------------------------- c Declarations: c ------------- #include "dimensions.h" #include "paramet.h" #include "comdissnew.h" #include "comgeom.h" !#include "com_io_dyn.h" #include "iniprint.h" #include "academic.h" REAL,INTENT(IN) :: time_0 ! not used c dynamical variables: REAL,INTENT(INOUT) :: ucov(ip1jmp1,llm) ! zonal covariant wind REAL,INTENT(INOUT) :: vcov(ip1jm,llm) ! meridional covariant wind REAL,INTENT(INOUT) :: teta(ip1jmp1,llm) ! potential temperature REAL,INTENT(INOUT) :: ps(ip1jmp1) ! surface pressure (Pa) REAL,INTENT(INOUT) :: masse(ip1jmp1,llm) ! air mass REAL,INTENT(INOUT) :: phis(ip1jmp1) ! geopotentiat at the surface REAL,INTENT(INOUT) :: q(ip1jmp1,llm,nqtot) ! advected tracers REAL p (ip1jmp1,llmp1 ) ! interlayer pressure REAL pks(ip1jmp1) ! exner at the surface REAL pk(ip1jmp1,llm) ! exner at mid-layer REAL pkf(ip1jmp1,llm) ! filtered exner at mid-layer REAL phi(ip1jmp1,llm) ! geopotential REAL w(ip1jmp1,llm) ! vertical velocity REAL kpd(ip1jmp1) ! TB15 = exp (-z/H) ! ADAPTATION GCM POUR CP(T) REAL temp(ip1jmp1,llm) ! temperature REAL tsurpk(ip1jmp1,llm) ! cpp*T/pk ! nogcm REAL tau_ps real globaverage2d c variables dynamiques intermediaire pour le transport REAL pbaru(ip1jmp1,llm),pbarv(ip1jm,llm) !flux de masse c variables dynamiques au pas -1 REAL vcovm1(ip1jm,llm),ucovm1(ip1jmp1,llm) REAL tetam1(ip1jmp1,llm),psm1(ip1jmp1) REAL massem1(ip1jmp1,llm) c tendances dynamiques en */s REAL dv(ip1jm,llm),du(ip1jmp1,llm) REAL dteta(ip1jmp1,llm),dq(ip1jmp1,llm,nqtot) c tendances physiques */s REAL dvfi(ip1jm,llm),dufi(ip1jmp1,llm) REAL dtetafi(ip1jmp1,llm),dqfi(ip1jmp1,llm,nqtot),dpfi(ip1jmp1) c variables pour le fichier histoire ! REAL dtav ! intervalle de temps elementaire REAL tppn(iim),tpps(iim),tpn,tps c INTEGER itau,itaufinp1,iav ! INTEGER iday ! jour julien REAL time REAL SSUM cym LOGICAL lafin LOGICAL :: lafin=.false. INTEGER ij,iq,l REAL rdaym_ini ! jD_cur: jour julien courant ! jH_cur: heure julienne courante REAL :: jD_cur, jH_cur ! INTEGER :: an, mois, jour ! REAL :: secondes LOGICAL first,callinigrads cIM : pour sortir les param. du modele dans un fis. netcdf 110106 save first data first/.true./ logical ok_sync parameter (ok_sync = .true.) logical physics data callinigrads/.true./ character*10 string10 ! REAL alpha(ip1jmp1,llm),beta(ip1jmp1,llm) REAL :: flxw(ip1jmp1,llm) ! flux de masse verticale REAL psmean ! pression moyenne REAL pqn2mean ! moyenne globale ps*qco REAL pqch4mean ! moyenne globale ps*qch4 REAL pqcomean ! moyenne globale ps*qco REAL pqhazemean REAL pqprechazemean REAL qmean_ch4 ! mass mean mixing ratio vap ch4 REAL qmean_co ! mass mean mixing ratio vap co REAL qmean_haze REAL qmean_prechaze REAL pqch4(ip1jmp1) ! average methane mass index : ps*qch4 REAL pqco(ip1jmp1) ! average CO mass index : ps*q_co REAL pqhaze(ip1jmp1) REAL pqprechaze(ip1jmp1) REAL p0 ! pression de reference REAL p00d ! globalaverage(kpd) REAL qmean_n2,qmean_n2_vert ! mass mean mixing ratio vap n2 REAL pqn2(ip1jmp1) ! average n2 mass index : ps*q_n2 REAL oldps(ip1jmp1) ! saving old pressure ps to calculate qch4 c+jld variables test conservation energie REAL ecin(ip1jmp1,llm),ecin0(ip1jmp1,llm) C Tendance de la temp. potentiel d (theta)/ d t due a la C tansformation d'energie cinetique en energie thermique C cree par la dissipation REAL dtetaecdt(ip1jmp1,llm) REAL vcont(ip1jm,llm),ucont(ip1jmp1,llm) REAL vnat(ip1jm,llm),unat(ip1jmp1,llm) CHARACTER*15 ztit character(len=*),parameter :: modname="leapfrog" character*80 abort_message logical dissip_conservative save dissip_conservative data dissip_conservative/.true./ INTEGER testita PARAMETER (testita = 9) logical , parameter :: flag_verif = .false. ! for CP(T) real :: dtec real :: ztetaec(ip1jmp1,llm) ! TEMP : diagnostic mass real :: n2mass(iip1,jjp1) real :: n2ice_ij(iip1,jjp1) integer :: i,j,ig integer, parameter :: ngrid = 2+(jjm-1)*iim ! local mass !real :: mq(ip1jmp1,llm) if (nday>=0) then itaufin = nday*day_step else ! to run a given (-nday) number of dynamical steps itaufin = -nday endif if (less1day) then c MODIF VENUS: to run less than one day: itaufin = int(fractday*day_step) endif if (ndynstep.gt.0) then ! running a given number of dynamical steps itaufin=ndynstep endif itaufinp1 = itaufin +1 c INITIALISATIONS dufi(:,:) =0. dvfi(:,:) =0. dtetafi(:,:)=0. dqfi(:,:,:) =0. dpfi(:) =0. itau = 0 physics=.true. if (iflag_phys==0.or.iflag_phys==2) physics=.false. ! ED18 nogcm ! firstcall_globaverage2d = .true. c iday = day_ini+itau/day_step c time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0 c IF(time.GT.1.) THEN c time = time-1. c iday = iday+1 c ENDIF c----------------------------------------------------------------------- c On initialise la pression et la fonction d'Exner : c -------------------------------------------------- dq(:,:,:)=0. CALL pression ( ip1jmp1, ap, bp, ps, p ) if (pressure_exner) then CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf ) else CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf ) endif c------------------ c TEST PK MONOTONE c------------------ write(*,*) "Test PK" do ij=1,ip1jmp1 do l=2,llm ! PRINT*,'pk(ij,l) = ',pk(ij,l) ! PRINT*,'pk(ij,l-1) = ',pk(ij,l-1) if(pk(ij,l).gt.pk(ij,l-1)) then c write(*,*) ij,l,pk(ij,l) abort_message = 'PK non strictement decroissante' call abort_gcm(modname,abort_message,1) c write(*,*) "ATTENTION, Test PK deconnecte..." endif enddo enddo write(*,*) "Fin Test PK" c stop c------------------ c Preparing mixing of pressure and tracers in nogcm c kpd=exp(-z/H) Needed to define a reference pressure : c P0=pmean/globalaverage(kpd) c thus P(i) = P0*exp(-z(i)/H) = P0*kpd(i) c It is checked that globalaverage2d(Pi)=pmean DO ij=1,ip1jmp1 kpd(ij) = exp(-phis(ij)/(r*38.)) ENDDO p00d=globaverage2d(kpd) ! mean pres at ref level tau_ps = tau_n2 ! constante de rappel for pressure (s) PRINT*,'igcm_n2 = ',igcm_n2 do iq=1,nqtot if (tname(iq)=="n2") then igcm_n2=iq exit endif enddo c----------------------------------------------------------------------- c Debut de l'integration temporelle: c ---------------------------------- 1 CONTINUE ! Matsuno Forward step begins here c date: (NB: date remains unchanged for Backward step) c ----- jD_cur = jD_ref + day_ini - day_ref + & (itau+1)/day_step jH_cur = jH_ref + start_time + & mod(itau+1,day_step)/float(day_step) jD_cur = jD_cur + int(jH_cur) jH_cur = jH_cur - int(jH_cur) ! Save fields obtained at previous time step as '...m1' CALL SCOPY( ijmllm ,vcov , 1, vcovm1 , 1 ) CALL SCOPY( ijp1llm,ucov , 1, ucovm1 , 1 ) CALL SCOPY( ijp1llm,teta , 1, tetam1 , 1 ) CALL SCOPY( ijp1llm,masse, 1, massem1, 1 ) CALL SCOPY( ip1jmp1, ps , 1, psm1 , 1 ) forward = .TRUE. leapf = .FALSE. dt = dtvr 2 CONTINUE ! Matsuno backward or leapfrog step begins here c gestion des appels de la physique et des dissipations: c ------------------------------------------------------ apphys = .FALSE. statcl = .FALSE. conser = .FALSE. IF( purmats ) THEN ! Purely Matsuno time stepping IF( MOD(itau,iconser) .EQ.0.AND. forward ) conser = .TRUE. IF( MOD(itau,iphysiq ).EQ.0.AND..NOT.forward s .and. physics ) apphys = .TRUE. ELSE ! Leapfrog/Matsuno time stepping IF( MOD(itau ,iconser) .EQ. 0 ) conser = .TRUE. IF( MOD(itau+1,iphysiq).EQ.0.AND.physics) apphys=.TRUE. END IF c----------------------------------------------------------------------- c calcul des tendances dynamiques: c -------------------------------- dv(:,:) = 0. du(:,:) = 0. dteta(:,:) = 0. dq(:,:,:) = 0. c----------------------------------------------------------------------- c calcul des tendances physiques: c ------------------------------- c IF( purmats ) THEN IF( itau.EQ.itaufin.AND..NOT.forward ) lafin = .TRUE. ELSE IF( itau+1. EQ. itaufin ) lafin = .TRUE. ENDIF IF( apphys ) THEN CALL pression ( ip1jmp1, ap, bp, ps, p ) if (pressure_exner) then CALL exner_hyb( ip1jmp1, ps, p,pks, pk, pkf ) else CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf ) endif ! Compute geopotential (physics might need it) CALL geopot ( ip1jmp1, teta , pk , pks, phis , phi ) jD_cur = jD_ref + day_ini - day_ref + & (itau+1)/day_step ! AS: we make jD_cur to be pday jD_cur = int(day_ini + itau/day_step) jH_cur = jH_ref + start_time + & mod(itau+1,day_step)/float(day_step) jH_cur = jH_ref + start_time + & mod(itau,day_step)/float(day_step) jD_cur = jD_cur + int(jH_cur) jH_cur = jH_cur - int(jH_cur) c Interface avec les routines de phylmd (phymars ... ) c ----------------------------------------------------- CALL calfis( lafin , jD_cur, jH_cur, $ ucov,vcov,teta,q,masse,ps,p,pk,phis,phi , $ du,dv,dteta,dq, $ flxw, $ dufi,dvfi,dtetafi,dqfi,dpfi ) c ajout des tendances physiques c ------------------------------ c CALL addfi( dtphys, leapf, forward , c $ ucov, vcov, teta , q ,ps , c $ dufi, dvfi, dtetafi , dqfi ,dpfi ) c CALL pression (ip1jmp1,ap,bp,ps,p) c CALL massdair(p,masse) ! 1) Saving pressure and mass of N2 tracer in each mess (kg) before Horizontal mixing of pressure DO ij=1,ip1jmp1 oldps(ij)=ps(ij) ENDDO !DO l=1, llm ! DO ij=1,ip1jmp1 ! mq(ij,l) = masse(ij,l)*q(ij,l,igcm_n2) ! ENDDO !ENDDO ! 2) Increment q with physical tendancy (mixing ratio kg/kg) DO ij=1,ip1jmp1 if (methane) then q(ij,1,igcm_ch4_gas)=q(ij,1,igcm_ch4_gas)+ & dqfi(ij,1,igcm_ch4_gas)*dtphys endif if (carbox) then q(ij,1,igcm_co_gas)=q(ij,1,igcm_co_gas)+ & dqfi(ij,1,igcm_co_gas)*dtphys endif if (fasthaze) then q(ij,1,igcm_haze)=q(ij,1,igcm_haze)+ & dqfi(ij,1,igcm_haze)*dtphys q(ij,1,igcm_prec_haze)=q(ij,1,igcm_prec_haze)+ & dqfi(ij,1,igcm_prec_haze)*dtphys endif ENDDO ! 3) Update pressure DO ij=1,ip1jmp1 ps(ij) = ps(ij) + dpfi(ij)*dtphys ENDDO c ! 4) Horizontal mixing of pressure ! Rappel newtonien vers psmean psmean= globaverage2d(ps) ! mean pressure p0=psmean/p00d DO ij=1,ip1jmp1 ps(ij)=ps(ij) +(p0*kpd(ij) & -ps(ij))*(1-exp(-dtphys/tau_ps)) ENDDO write(72,*) 'psmean = ',psmean ! mean pressure redistributed ! 5) Security check: test pressure negative DO ij=1,ip1jmp1 IF (ps(ij).le.0.) then !write(*,*) 'Negative pressure :' !write(*,*) 'ps = ',ps(ij),' ig = ',ij !write(*,*) 'pmean = ',p0*kpd(ij) !write(*,*) 'tau_ps = ',tau_ps !STOP ps(ij)=1.e-6*kpd(ij)/p00d ENDIF ENDDO ! 6) Horizontal redistribution of other tracers ! intermediaire de calcul si methane or CO if ((methane).or.(carbox)) then psmean= globaverage2d(ps) endif ! 6a) Update mmr q for new pressure (both q and dqfi were calculated with old pressure in physiq) DO ij=1,ip1jmp1 if (methane) then q(ij,1,igcm_ch4_gas)=q(ij,1,igcm_ch4_gas)*oldps(ij)/ & ps(ij) ! average methane mass index : ps * qch4 => 'pressure of ch4' pqch4(ij)= ps(ij) * q(ij,1,igcm_ch4_gas) endif if (carbox) then q(ij,1,igcm_co_gas)=q(ij,1,igcm_co_gas)*oldps(ij)/ & ps(ij) pqco(ij)= ps(ij) * q(ij,1,igcm_co_gas) endif if (fasthaze) then q(ij,1,igcm_haze)=q(ij,1,igcm_haze)*oldps(ij)/ & ps(ij) pqhaze(ij)= ps(ij) * q(ij,1,igcm_haze) q(ij,1,igcm_prec_haze)=q(ij,1,igcm_prec_haze)*oldps(ij)/ & ps(ij) pqprechaze(ij)= ps(ij) * q(ij,1,igcm_prec_haze) endif ENDDO ! 6b) Rappel newtonien vers q_mean if (methane) then pqch4mean=globaverage2d(pqch4) qmean_ch4= pqch4mean / psmean endif if (carbox) then pqcomean=globaverage2d(pqco) qmean_co= pqcomean / psmean endif if (fasthaze) then pqprechazemean=globaverage2d(pqprechaze) qmean_prechaze= pqprechazemean / psmean !pqhazemean=globaverage2d(pqhaze) !qmean_haze= pqhazemean / psmean endif DO ij=1,ip1jmp1 if (methane) then q(ij,1,igcm_ch4_gas)=q(ij,1,igcm_ch4_gas)+ & (qmean_ch4-q(ij,1,igcm_ch4_gas))* & (1.-exp(-dtphys/tau_ch4)) endif if (carbox) then q(ij,1,igcm_co_gas)=q(ij,1,igcm_co_gas)+ & (qmean_co-q(ij,1,igcm_co_gas))* & (1.-exp(-dtphys/tau_co)) endif ! CO if (fasthaze) then ! q(ij,1,igcm_haze)=q(ij,1,igcm_haze)+ & ! (qmean_haze-q(ij,1,igcm_haze))* & ! (1.-exp(-dtphys/tau_haze)) q(ij,1,igcm_prec_haze)=q(ij,1,igcm_prec_haze)+ & (qmean_prechaze-q(ij,1,igcm_prec_haze))* & (1.-exp(-dtphys/tau_prechaze)) endif ! fasthaze ENDDO ! 7) Update pressure p and pk CALL pression (ip1jmp1,ap,bp,ps,p) CALL massdair(p,masse) if (pressure_exner) then CALL exner_hyb( ip1jmp1, ps, p, pks, pk, pkf ) else CALL exner_milieu( ip1jmp1, ps, p, pks, pk, pkf ) endif ENDIF ! of IF( apphys ) c ******************************************************************** c ******************************************************************** c .... fin de l'integration dynamique et physique pour le pas itau .. c ******************************************************************** c ******************************************************************** IF ( .NOT.purmats ) THEN c ........................................................ c .............. schema matsuno + leapfrog .............. c ........................................................ IF(forward. OR. leapf) THEN itau= itau + 1 c iday= day_ini+itau/day_step c time= REAL(itau-(iday-day_ini)*day_step)/day_step+time_0 c IF(time.GT.1.) THEN c time = time-1. c iday = iday+1 c ENDIF ENDIF IF( itau. EQ. itaufinp1 ) then if (flag_verif) then write(79,*) 'ucov',ucov write(80,*) 'vcov',vcov write(81,*) 'teta',teta write(82,*) 'ps',ps write(83,*) 'q',q WRITE(85,*) 'q1 = ',q(:,:,1) WRITE(86,*) 'q3 = ',q(:,:,3) endif abort_message = 'Simulation finished' call abort_gcm(modname,abort_message,0) ENDIF c----------------------------------------------------------------------- c ecriture du fichier histoire moyenne: c ------------------------------------- IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN IF(itau.EQ.itaufin) THEN iav=1 ELSE iav=0 ENDIF ! ! Ehouarn: re-compute geopotential for outputs CALL geopot(ip1jmp1,teta,pk,pks,phis,phi) IF (ok_dyn_ave) THEN #ifdef CPP_IOIPSL CALL writedynav(itau,vcov, & ucov,teta,pk,phi,q,masse,ps,phis) #endif ENDIF ENDIF ! of IF((MOD(itau,iperiod).EQ.0).OR.(itau.EQ.itaufin)) if (ok_iso_verif) then call check_isotopes_seq(q,ip1jmp1,'leapfrog 1584') endif !if (ok_iso_verif) then c----------------------------------------------------------------------- c ecriture de la bande histoire: c ------------------------------ IF( MOD(itau,iecri).EQ.0) THEN ! Ehouarn: output only during LF or Backward Matsuno if (leapf.or.(.not.leapf.and.(.not.forward))) then ! ADAPTATION GCM POUR CP(T) call tpot2t(ijp1llm,teta,temp,pk) tsurpk = cpp*temp/pk CALL geopot(ip1jmp1,tsurpk,pk,pks,phis,phi) unat=0. do l=1,llm unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm) vnat(:,l)=vcov(:,l)/cv(:) enddo #ifdef CPP_IOIPSL if (ok_dyn_ins) then ! write(lunout,*) "leapfrog: call writehist, itau=",itau CALL writehist(itau,vcov,ucov,teta,phi,q,masse,ps,phis) ! call WriteField('ucov',reshape(ucov,(/iip1,jmp1,llm/))) ! call WriteField('vcov',reshape(vcov,(/iip1,jjm,llm/))) ! call WriteField('teta',reshape(teta,(/iip1,jmp1,llm/))) ! call WriteField('ps',reshape(ps,(/iip1,jmp1/))) ! call WriteField('masse',reshape(masse,(/iip1,jmp1,llm/))) endif ! of if (ok_dyn_ins) #endif ! For some Grads outputs of fields if (output_grads_dyn) then #include "write_grads_dyn.h" endif endif ! of if (leapf.or.(.not.leapf.and.(.not.forward))) ENDIF ! of IF(MOD(itau,iecri).EQ.0) IF(itau.EQ.itaufin) THEN CALL dynredem1("restart.nc",start_time, & vcov,ucov,teta,q,masse,ps) CLOSE(99) !!! Ehouarn: Why not stop here and now? ENDIF ! of IF (itau.EQ.itaufin) c----------------------------------------------------------------------- c gestion de l'integration temporelle: c ------------------------------------ IF( MOD(itau,iperiod).EQ.0 ) THEN GO TO 1 ELSE IF ( MOD(itau-1,iperiod). EQ. 0 ) THEN IF( forward ) THEN c fin du pas forward et debut du pas backward forward = .FALSE. leapf = .FALSE. GO TO 2 ELSE c fin du pas backward et debut du premier pas leapfrog leapf = .TRUE. dt = 2.*dtvr GO TO 2 END IF ! of IF (forward) ELSE c ...... pas leapfrog ..... leapf = .TRUE. dt = 2.*dtvr GO TO 2 END IF ! of IF (MOD(itau,iperiod).EQ.0) ! ELSEIF (MOD(itau-1,iperiod).EQ.0) ELSE ! of IF (.not.purmats) if (ok_iso_verif) then call check_isotopes_seq(q,ip1jmp1,'leapfrog 1664') endif !if (ok_iso_verif) then c ........................................................ c .............. schema matsuno ............... c ........................................................ IF( forward ) THEN itau = itau + 1 c iday = day_ini+itau/day_step c time = REAL(itau-(iday-day_ini)*day_step)/day_step+time_0 c c IF(time.GT.1.) THEN c time = time-1. c iday = iday+1 c ENDIF forward = .FALSE. IF( itau. EQ. itaufinp1 ) then abort_message = 'Simulation finished' call abort_gcm(modname,abort_message,0) ENDIF GO TO 2 ELSE ! of IF(forward) i.e. backward step if (ok_iso_verif) then call check_isotopes_seq(q,ip1jmp1,'leapfrog 1698') endif !if (ok_iso_verif) then IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) THEN IF(itau.EQ.itaufin) THEN iav=1 ELSE iav=0 ENDIF ! ! Ehouarn: re-compute geopotential for outputs CALL geopot(ip1jmp1,teta,pk,pks,phis,phi) IF (ok_dyn_ave) THEN #ifdef CPP_IOIPSL CALL writedynav(itau,vcov, & ucov,teta,pk,phi,q,masse,ps,phis) #endif ENDIF ENDIF ! of IF(MOD(itau,iperiod).EQ.0 .OR. itau.EQ.itaufin) IF(MOD(itau,iecri ).EQ.0) THEN c IF(MOD(itau,iecri*day_step).EQ.0) THEN ! ADAPTATION GCM POUR CP(T) call tpot2t(ijp1llm,teta,temp,pk) tsurpk = cpp*temp/pk CALL geopot(ip1jmp1,tsurpk,pk,pks,phis,phi) unat=0. do l=1,llm unat(iip2:ip1jm,l)=ucov(iip2:ip1jm,l)/cu(iip2:ip1jm) vnat(:,l)=vcov(:,l)/cv(:) enddo #ifdef CPP_IOIPSL if (ok_dyn_ins) then ! write(lunout,*) "leapfrog: call writehist (b)", ! & itau,iecri CALL writehist(itau,vcov,ucov,teta,phi,q,masse,ps,phis) endif ! of if (ok_dyn_ins) #endif ! For some Grads outputs if (output_grads_dyn) then #include "write_grads_dyn.h" endif ENDIF ! of IF(MOD(itau,iecri ).EQ.0) IF(itau.EQ.itaufin) THEN CALL dynredem1("restart.nc",start_time, & vcov,ucov,teta,q,masse,ps) ENDIF ! of IF(itau.EQ.itaufin) forward = .TRUE. GO TO 1 ENDIF ! of IF (forward) END IF ! of IF(.not.purmats) STOP END ! ************************************************************ ! globalaverage VERSION PLuto ! ************************************************************ FUNCTION globaverage2d(var) ! ************************************************************ IMPLICIT NONE #include "dimensions.h" #include "paramet.h" #include "comgeom2.h" !ip1jmp1 called in paramet.h = iip1 x jjp1 REAL var(iip1,jjp1), globaverage2d integer i,j REAL airetot save airetot logical firstcall data firstcall/.true./ save firstcall if (firstcall) then firstcall=.false. airetot =0. DO j=2,jjp1-1 ! lat DO i = 1, iim ! lon airetot = airetot + aire(i,j) END DO END DO DO i=1,iim airetot=airetot + aire(i,1)+aire(i,jjp1) ENDDO end if globaverage2d = 0. DO j=2,jjp1-1 DO i = 1, iim globaverage2d = globaverage2d + var(i,j)*aire(i,j) END DO END DO DO i=1,iim globaverage2d = globaverage2d + var(i,1)*aire(i,1) globaverage2d = globaverage2d + var(i,jjp1)*aire(i,jjp1) ENDDO globaverage2d = globaverage2d/airetot return end