/[lmdze]/trunk/libf/bibio/writedynav.f90
ViewVC logotype

Diff of /trunk/libf/bibio/writedynav.f90

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

revision 66 by guez, Thu Jul 26 14:37:37 2012 UTC revision 67 by guez, Tue Oct 2 15:50:56 2012 UTC
# Line 4  module writedynav_m Line 4  module writedynav_m
4    
5  contains  contains
6    
7    subroutine writedynav(vcov, ucov, teta, ppk, phi, q, masse, ps, phis, time)    subroutine writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, time)
8    
9      ! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05      ! 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      ! Écriture du fichier histoire au format IOIPSL
# Line 12  contains Line 12  contains
12    
13      ! Appels successifs des routines histwrite      ! Appels successifs des routines histwrite
14    
15        USE comconst, ONLY: cpp
16      use covnat_m, only: covnat      use covnat_m, only: covnat
17      USE histwrite_m, ONLY: histwrite      USE dimens_m, ONLY: iim, jjm, llm, nqmx
18      USE histsync_m, ONLY: histsync      USE histsync_m, ONLY: histsync
19      USE dimens_m, ONLY: llm      USE histwrite_m, ONLY: histwrite
     USE paramet_m, ONLY: iip1, ijp1llm, ip1jm, ip1jmp1, jjp1  
     USE comconst, ONLY: cpp  
     USE temps, ONLY: itau_dyn  
20      USE iniadvtrac_m, ONLY: ttext      USE iniadvtrac_m, ONLY: ttext
21      use initdynav_m, only: histaveid      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      REAL, intent(in):: vcov(ip1jm, llm), ucov(ip1jmp1, llm) ! vents covariants      ! Vents covariants :
27      REAL, intent(in):: teta(ip1jmp1*llm) ! temperature potentielle      REAL, intent(in):: vcov(:, :, :) ! (iim + 1, jjm, llm)
28      real, intent(in):: phi(ip1jmp1, llm) ! geopotentiel instantane      REAL, intent(in):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm)
29      real, intent(in):: ppk(ip1jmp1*llm)  
30      REAL, intent(in):: ps(ip1jmp1) ! pression au sol      REAL, intent(in):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
31      real, intent(in):: masse(ip1jmp1, llm)      ! temperature potentielle
32      REAL, intent(in):: phis(ip1jmp1) ! geopotentiel au sol  
33        real, intent(in):: phi(:, :, :) ! (iim + 1, jjm + 1, llm)
34        ! geopotentiel instantane
35    
36        real, intent(in):: pk(:, :, :) ! (iim + 1, jjm + 1, llm)
37        REAL, intent(in):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol
38        real, intent(in):: masse(:, :, :) ! (iim + 1, jjm + 1, llm)
39        REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) geopotentiel au sol
40      REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) traceurs      REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) traceurs
41      integer, intent(in):: time ! temps de l'ecriture      integer, intent(in):: time ! temps de l'ecriture
42    
43      ! Variables locales      ! Variables locales
44      integer ndex2d(iip1*jjp1), ndex3d(iip1*jjp1*llm), iq, ii, ll      integer iq
45      real us(ip1jmp1*llm), vs(ip1jmp1*llm)      real us(ip1jmp1*llm), vs(ip1jmp1*llm)
     real tm(ip1jmp1*llm)  
46      REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)      REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
     logical ok_sync  
47      integer itau_w      integer itau_w
48    
49      !---------------------------------------------------------------      !---------------------------------------------------------------
50    
51        call assert((/size(vcov, 1), size(ucov, 1), size(teta, 1), size(phi, 1), &
52             size(pk, 1), size(ps, 1), size(masse, 1), size(phis, 1), &
53             size(q, 1)/) == iim + 1, "writedynav iim")
54        call assert((/size(vcov, 2) + 1, size(ucov, 2), size(teta, 2), &
55             size(phi, 2), size(pk, 2), size(ps, 2), size(masse, 2), &
56             size(phis, 2), size(q, 2)/) == jjm + 1, "writedynav jjm")
57        call assert((/size(vcov, 3), size(ucov, 3), size(teta, 3), size(phi, 3), &
58             size(pk, 3), size(masse, 3), size(q, 3)/) == llm, "writedynav llm")
59        call assert(size(q, 4) == nqmx, "writedynav nqmx")
60    
61      ! Initialisations      ! Initialisations
     ndex3d = 0  
     ndex2d = 0  
     ok_sync = .TRUE.  
62      us = 999.999      us = 999.999
63      vs = 999.999      vs = 999.999
     tm = 999.999  
64      vnat = 999.999      vnat = 999.999
65      unat = 999.999      unat = 999.999
66      itau_w = itau_dyn + time      itau_w = itau_dyn + time
# Line 70  contains Line 82  contains
82      call histwrite(histaveid, 'theta', itau_w, teta)      call histwrite(histaveid, 'theta', itau_w, teta)
83    
84      ! Temperature moyennee      ! Temperature moyennee
85      do ii = 1, ijp1llm      call histwrite(histaveid, 'temp', itau_w, teta * pk / cpp)
        tm(ii) = teta(ii) * ppk(ii)/cpp  
     enddo  
     call histwrite(histaveid, 'temp', itau_w, tm)  
86    
     ! Geopotentiel  
87      call histwrite(histaveid, 'phi', itau_w, phi)      call histwrite(histaveid, 'phi', itau_w, phi)
88    
89      ! Traceurs      ! Traceurs
# Line 83  contains Line 91  contains
91         call histwrite(histaveid, ttext(iq), itau_w, q(:, :, :, iq))         call histwrite(histaveid, ttext(iq), itau_w, q(:, :, :, iq))
92      enddo      enddo
93    
     ! Masse  
94      call histwrite(histaveid, 'masse', itau_w, masse)      call histwrite(histaveid, 'masse', itau_w, masse)
   
     ! Pression au sol  
95      call histwrite(histaveid, 'ps', itau_w, ps)      call histwrite(histaveid, 'ps', itau_w, ps)
   
     ! Geopotentiel au sol  
96      call histwrite(histaveid, 'phis', itau_w, phis)      call histwrite(histaveid, 'phis', itau_w, phis)
97    
98      if (ok_sync) call histsync(histaveid)      call histsync(histaveid)
99    
100    end subroutine writedynav    end subroutine writedynav
101    

Legend:
Removed from v.66  
changed lines
  Added in v.67

  ViewVC Help
Powered by ViewVC 1.1.21