--- trunk/libf/bibio/writedynav.f90 2011/04/27 13:00:12 45 +++ trunk/libf/bibio/writedynav.f90 2012/01/10 19:02:02 56 @@ -1,103 +1,109 @@ -subroutine writedynav(histid, nq, time, vcov, ucov, teta, ppk, phi, q, masse, & - ps, phis) +module writedynav_m - ! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05 - ! Ecriture du fichier histoire au format IOIPSL + implicit none - ! Appels successifs des routines histwrite +contains - ! Entree: - ! histid: ID du fichier histoire - ! time: temps de l'ecriture - ! vcov: vents v covariants - ! ucov: vents u covariants - ! phi : geopotentiel instantane - ! q : traceurs - ! ps :pression au sol - ! phis : geopotentiel au sol - - ! L. Fairhead, LMD, 03/99 - - USE histwrite_m, ONLY : histwrite - USE histcom, ONLY : histsync - USE dimens_m, ONLY : llm - USE paramet_m, ONLY : iip1, ijp1llm, ip1jm, ip1jmp1, jjp1 - USE comconst, ONLY : cpp - USE temps, ONLY : itau_dyn - USE iniadvtrac_m, ONLY : ttext + subroutine writedynav(nq, time, vcov, ucov, teta, ppk, phi, q, masse, ps, & + phis) - implicit none + ! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05 + ! Ecriture du fichier histoire au format IOIPSL + + ! Appels successifs des routines histwrite + + ! Entree: + ! time: temps de l'ecriture + ! vcov: vents v covariants + ! ucov: vents u covariants + ! phi : geopotentiel instantane + ! q : traceurs + ! ps :pression au sol + ! phis : geopotentiel au sol + + ! L. Fairhead, LMD, 03/99 + + USE histwrite_m, ONLY: histwrite + USE histcom, ONLY: histsync + USE dimens_m, ONLY: llm + USE paramet_m, ONLY: iip1, ijp1llm, ip1jm, ip1jmp1, jjp1 + USE comconst, ONLY: cpp + USE temps, ONLY: itau_dyn + USE iniadvtrac_m, ONLY: ttext + use initdynav_m, only: histaveid + + INTEGER nq + REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm) + REAL, intent(in):: teta(ip1jmp1*llm) ! temperature potentielle + real phi(ip1jmp1, llm), ppk(ip1jmp1*llm) + REAL ps(ip1jmp1) + real, intent(in):: masse(ip1jmp1, llm) + REAL phis(ip1jmp1) + REAL q(ip1jmp1, llm, nq) + integer, intent(in):: time + + ! Variables locales + integer ndex2d(iip1*jjp1), ndex3d(iip1*jjp1*llm), iq, ii, ll + real us(ip1jmp1*llm), vs(ip1jmp1*llm) + real tm(ip1jmp1*llm) + REAL vnat(ip1jm, llm), unat(ip1jmp1, llm) + logical ok_sync + integer itau_w + + !--------------------------------------------------------------- + + ! Initialisations + ndex3d = 0 + ndex2d = 0 + ok_sync = .TRUE. + us = 999.999 + vs = 999.999 + tm = 999.999 + vnat = 999.999 + unat = 999.999 + itau_w = itau_dyn + time + + ! Passage aux composantes naturelles du vent + call covnat(llm, ucov, vcov, unat, vnat) + + ! Appels a histwrite pour l'ecriture des variables a sauvegarder + + ! Vents U scalaire + call gr_u_scal(llm, unat, us) + call histwrite(histaveid, 'u', itau_w, us) + + ! Vents V scalaire + call gr_v_scal(llm, vnat, vs) + call histwrite(histaveid, 'v', itau_w, vs) + + ! Temperature potentielle moyennee + call histwrite(histaveid, 'theta', itau_w, teta) + + ! Temperature moyennee + do ii = 1, ijp1llm + tm(ii) = teta(ii) * ppk(ii)/cpp + enddo + call histwrite(histaveid, 'temp', itau_w, tm) + + ! Geopotentiel + call histwrite(histaveid, 'phi', itau_w, phi) + + ! Traceurs + DO iq=1, nq + call histwrite(histaveid, ttext(iq), itau_w, q(:, :, iq)) + enddo + + ! Masse + call histwrite(histaveid, 'masse', itau_w, masse) + + ! Pression au sol + call histwrite(histaveid, 'ps', itau_w, ps) - INTEGER histid, nq - REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm) - REAL, intent(in):: teta(ip1jmp1*llm) ! temperature potentielle - real phi(ip1jmp1, llm), ppk(ip1jmp1*llm) - REAL ps(ip1jmp1) - real, intent(in):: masse(ip1jmp1, llm) - REAL phis(ip1jmp1) - REAL q(ip1jmp1, llm, nq) - integer, intent(in):: time - - ! Variables locales - integer ndex2d(iip1*jjp1), ndex3d(iip1*jjp1*llm), iq, ii, ll - real us(ip1jmp1*llm), vs(ip1jmp1*llm) - real tm(ip1jmp1*llm) - REAL vnat(ip1jm, llm), unat(ip1jmp1, llm) - logical ok_sync - integer itau_w - - !--------------------------------------------------------------- - - ! Initialisations - ndex3d = 0 - ndex2d = 0 - ok_sync = .TRUE. - us = 999.999 - vs = 999.999 - tm = 999.999 - vnat = 999.999 - unat = 999.999 - itau_w = itau_dyn + time - - ! Passage aux composantes naturelles du vent - call covnat(llm, ucov, vcov, unat, vnat) - - ! Appels a histwrite pour l'ecriture des variables a sauvegarder - - ! Vents U scalaire - call gr_u_scal(llm, unat, us) - call histwrite(histid, 'u', itau_w, us) - - ! Vents V scalaire - call gr_v_scal(llm, vnat, vs) - call histwrite(histid, 'v', itau_w, vs) - - ! Temperature potentielle moyennee - call histwrite(histid, 'theta', itau_w, teta) - - ! Temperature moyennee - do ii = 1, ijp1llm - tm(ii) = teta(ii) * ppk(ii)/cpp - enddo - call histwrite(histid, 'temp', itau_w, tm) - - ! Geopotentiel - call histwrite(histid, 'phi', itau_w, phi) - - ! Traceurs - DO iq=1, nq - call histwrite(histid, ttext(iq), itau_w, q(:, :, iq)) - enddo - - ! Masse - call histwrite(histid, 'masse', itau_w, masse) - - ! Pression au sol - call histwrite(histid, 'ps', itau_w, ps) + ! Geopotentiel au sol + call histwrite(histaveid, 'phis', itau_w, phis) - ! Geopotentiel au sol - call histwrite(histid, 'phis', itau_w, phis) + if (ok_sync) call histsync(histaveid) - if (ok_sync) call histsync(histid) + end subroutine writedynav -end subroutine writedynav +end module writedynav_m