/[lmdze]/trunk/Sources/dyn3d/writedynav.f
ViewVC logotype

Diff of /trunk/Sources/dyn3d/writedynav.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

trunk/libf/bibio/writedynav.f90 revision 45 by guez, Wed Apr 27 13:00:12 2011 UTC trunk/bibio/writedynav.f90 revision 76 by guez, Fri Nov 15 18:45:49 2013 UTC
# Line 1  Line 1 
1  subroutine writedynav(histid, nq, time, vcov, ucov, teta, ppk, phi, q, masse, &  module writedynav_m
      ps, phis)  
2    
3    ! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05    implicit none
   ! Ecriture du fichier histoire au format IOIPSL  
4    
5    ! Appels successifs des routines histwrite  contains
6    
7    ! Entree:    subroutine writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, time)
   ! 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  
8    
9    implicit none      ! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05
10        ! Écriture du fichier histoire au format IOIPSL
11        ! L. Fairhead, LMD, 03/99
12    
13        ! Appels successifs des routines histwrite
14    
15        USE comconst, ONLY: cpp
16        use covnat_m, only: covnat
17        USE dimens_m, ONLY: iim, jjm, llm, nqmx
18        USE histsync_m, ONLY: histsync
19        USE histwrite_m, ONLY: histwrite
20        USE iniadvtrac_m, ONLY: ttext
21        use initdynav_m, only: histaveid
22        use nr_util, only: assert
23        USE paramet_m, ONLY: iip1, ip1jm, ip1jmp1, jjp1
24        USE temps, ONLY: itau_dyn
25    
26        ! Vents covariants :
27        REAL, intent(in):: vcov(:, :, :) ! (iim + 1, jjm, llm)
28        REAL, intent(in):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm)
29    
30        REAL, intent(in):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
31        ! temperature potentielle
32    
33        real, intent(in):: pk(:, :, :) ! (iim + 1, jjm + 1, llm)
34    
35        real, intent(in):: phi(:, :, :) ! (iim + 1, jjm + 1, llm)
36        ! geopotentiel instantane
37    
38        REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) traceurs
39        real, intent(in):: masse(:, :, :) ! (iim + 1, jjm + 1, llm)
40        REAL, intent(in):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol
41        REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) geopotentiel au sol
42        integer, intent(in):: time ! temps de l'ecriture
43    
44        ! Variables locales
45        integer iq
46        real us(ip1jmp1*llm), vs(ip1jmp1*llm)
47        REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
48        integer itau_w
49    
50        !---------------------------------------------------------------
51    
52        call assert((/size(vcov, 1), size(ucov, 1), size(teta, 1), size(phi, 1), &
53             size(pk, 1), size(ps, 1), size(masse, 1), size(phis, 1), &
54             size(q, 1)/) == iim + 1, "writedynav iim")
55        call assert((/size(vcov, 2) + 1, size(ucov, 2), size(teta, 2), &
56             size(phi, 2), size(pk, 2), size(ps, 2), size(masse, 2), &
57             size(phis, 2), size(q, 2)/) == jjm + 1, "writedynav jjm")
58        call assert((/size(vcov, 3), size(ucov, 3), size(teta, 3), size(phi, 3), &
59             size(pk, 3), size(masse, 3), size(q, 3)/) == llm, "writedynav llm")
60        call assert(size(q, 4) == nqmx, "writedynav nqmx")
61    
62        ! Initialisations
63        us = 999.999
64        vs = 999.999
65        vnat = 999.999
66        unat = 999.999
67        itau_w = itau_dyn + time
68    
69        ! Passage aux composantes naturelles du vent
70        call covnat(llm, ucov, vcov, unat, vnat)
71    
72        ! Appels a histwrite pour l'ecriture des variables a sauvegarder
73    
74        ! Vents U scalaire
75        call gr_u_scal(llm, unat, us)
76        call histwrite(histaveid, 'u', itau_w, us)
77    
78        ! Vents V scalaire
79        call gr_v_scal(llm, vnat, vs)
80        call histwrite(histaveid, 'v', itau_w, vs)
81    
82        ! Temperature potentielle moyennee
83        call histwrite(histaveid, 'theta', itau_w, teta)
84    
85        ! Temperature moyennee
86        call histwrite(histaveid, 'temp', itau_w, teta * pk / cpp)
87    
88        call histwrite(histaveid, 'phi', itau_w, phi)
89    
90        ! Traceurs
91        DO iq = 1, size(q, 4)
92           call histwrite(histaveid, ttext(iq), itau_w, q(:, :, :, iq))
93        enddo
94    
95    INTEGER histid, nq      call histwrite(histaveid, 'masse', itau_w, masse)
96    REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm)      call histwrite(histaveid, 'ps', itau_w, ps)
97    REAL, intent(in):: teta(ip1jmp1*llm) ! temperature potentielle      call histwrite(histaveid, 'phis', itau_w, phis)
   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)  
98    
99    ! Geopotentiel au sol      call histsync(histaveid)
   call histwrite(histid, 'phis', itau_w, phis)  
100    
101    if (ok_sync) call histsync(histid)    end subroutine writedynav
102    
103  end subroutine writedynav  end module writedynav_m

Legend:
Removed from v.45  
changed lines
  Added in v.76

  ViewVC Help
Powered by ViewVC 1.1.21