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

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

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

revision 30 by guez, Thu Apr 1 09:07:28 2010 UTC revision 52 by guez, Fri Sep 23 12:28:01 2011 UTC
# Line 1  Line 1 
1  module inithist_m  module inithist_m
2    
   ! This module is clean: no C preprocessor directive, no include line  
   
3    implicit none    implicit none
4    
5  contains  contains
6    
7    subroutine inithist(day0, anne0, tstep, nq, fileid, filevid, t_ops, &    subroutine inithist(day0, anne0, tstep, nq, fileid, filevid, t_ops, t_wrt)
        t_wrt)  
   
     ! From inithist.F,v 1.1.1.1 2004/05/19 12:53:05  
8    
9      !   Routine d'initialisation des ecritures des fichiers histoires LMDZ      ! From inithist.F, version 1.1.1.1 2004/05/19 12:53:05
     !   au format IOIPSL  
10    
11      !   Appels succesifs des routines: histbeg      ! Routine d'initialisation des écritures des fichiers histoires LMDZ
12      !                                  histhori      ! au format IOIPSL
13      !                                  histver      ! Appels successifs des routines : histbeg, histhori, histver,
14      !                                  histdef      ! histdef, histend
15      !                                  histend  
16        ! Entrées :
17      !   Entree:      ! day0, anne0: date de référence
18      !      day0,anne0: date de reference      ! tstep : durée du pas de temps en secondes
19      !      tstep: duree du pas de temps en seconde      ! t_ops : fréquence de l'opération pour IOIPSL
20      !      t_ops: frequence de l'operation pour IOIPSL      ! t_wrt : fréquence d'écriture sur le fichier
21      !      t_wrt: frequence d'ecriture sur le fichier      ! nq : nombre de traceurs
22      !      nq: nombre de traceurs  
23        ! Sorties :
24      !   Sortie:      ! fileid : ID du fichier Netcdf créé
25      !      fileid: ID du fichier netcdf cree      ! filevid : ID du fichier Netcdf pour la grille v
     !      filevid:ID du fichier netcdf pour la grille v  
26    
27      !   L. Fairhead, LMD, 03/99      ! L. Fairhead, LMD, 03/99
28    
29      USE calendar      USE calendar
30      use histcom      use histcom
# Line 45  contains Line 38  contains
38      use temps      use temps
39      use ener      use ener
40      use iniadvtrac_m      use iniadvtrac_m
41        use nr_util, only: pi
42    
43      !   Arguments      ! Arguments
44      integer day0, anne0      integer day0, anne0
45      real, intent(in):: tstep, t_ops, t_wrt      real, intent(in):: tstep, t_ops, t_wrt
46      integer fileid, filevid      integer fileid, filevid
47      integer nq      integer nq
48    
49      !   Variables locales      ! Variables locales
50      real zjulian      real zjulian
51      integer iq      integer iq
52      real rlong(iip1,jjp1), rlat(iip1,jjp1)      real rlong(iip1, jjp1), rlat(iip1, jjp1)
53      integer uhoriid, vhoriid, thoriid, zvertiid      integer uhoriid, vhoriid, thoriid, zvertiid
54      integer ii,jj      integer ii, jj
55      integer zan, dayref      integer zan, dayref
56    
57      !-----------------------------------------------------------------------      !-----------------------------------------------------------------------
58    
59      !  Initialisations      ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
   
     pi = 4. * atan (1.)  
   
     !  Appel a histbeg: creation du fichier netcdf et initialisations diverses  
60    
61      zan = anne0      zan = anne0
62      dayref = day0      dayref = day0
# Line 74  contains Line 64  contains
64    
65      do jj = 1, jjp1      do jj = 1, jjp1
66         do ii = 1, iip1         do ii = 1, iip1
67            rlong(ii,jj) = rlonu(ii) * 180. / pi            rlong(ii, jj) = rlonu(ii) * 180. / pi
68            rlat(ii,jj) = rlatu(jj) * 180. / pi            rlat(ii, jj) = rlatu(jj) * 180. / pi
69         enddo         enddo
70      enddo      enddo
71    
72      call histbeg_totreg("dyn_hist.nc", rlong(:,1), rlat(1,:), &      call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), &
73           1, iip1, 1, jjp1, &           1, iip1, 1, jjp1, &
74           itau_dyn, zjulian, tstep, uhoriid, fileid)           itau_dyn, zjulian, tstep, uhoriid, fileid)
75      !      !
76      !  Creation du fichier histoire pour la grille en V (oblige pour l'instant,      ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
77      !  IOIPSL ne permet pas de grilles avec des nombres de point differents dans      ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
78      !  un meme fichier)      ! un meme fichier)
79    
80      do jj = 1, jjm      do jj = 1, jjm
81         do ii = 1, iip1         do ii = 1, iip1
82            rlong(ii,jj) = rlonv(ii) * 180. / pi            rlong(ii, jj) = rlonv(ii) * 180. / pi
83            rlat(ii,jj) = rlatv(jj) * 180. / pi            rlat(ii, jj) = rlatv(jj) * 180. / pi
84         enddo         enddo
85      enddo      enddo
86    
87      call histbeg_totreg('dyn_histv.nc', rlong(:,1), rlat(1,:jjm), &      call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), &
88           1, iip1, 1, jjm, &           1, iip1, 1, jjm, &
89           itau_dyn, zjulian, tstep, vhoriid, filevid)           itau_dyn, zjulian, tstep, vhoriid, filevid)
90      !      !
91      !  Appel a histhori pour rajouter les autres grilles horizontales      ! Appel a histhori pour rajouter les autres grilles horizontales
92      !      !
93      do jj = 1, jjp1      do jj = 1, jjp1
94         do ii = 1, iip1         do ii = 1, iip1
95            rlong(ii,jj) = rlonv(ii) * 180. / pi            rlong(ii, jj) = rlonv(ii) * 180. / pi
96            rlat(ii,jj) = rlatu(jj) * 180. / pi            rlat(ii, jj) = rlatu(jj) * 180. / pi
97         enddo         enddo
98      enddo      enddo
99    
100      call histhori_regular(fileid, iip1, rlong, jjp1, rlat, 'scalar', &      call histhori_regular(fileid, iip1, rlong, jjp1, rlat, 'scalar', &
101           'Grille points scalaires', thoriid)           'Grille points scalaires', thoriid)
102      !      !
103      !  Appel a histvert pour la grille verticale      ! Appel a histvert pour la grille verticale
104      !      !
105      call histvert(fileid, 'sig_s', 'Niveaux sigma','-', &      call histvert(fileid, 'sig_s', 'Niveaux sigma', '-', &
106           llm, nivsigs, zvertiid)           llm, nivsigs, zvertiid)
107      ! Pour le fichier V      ! Pour le fichier V
108      call histvert(filevid, 'sig_s', 'Niveaux sigma','-', &      call histvert(filevid, 'sig_s', 'Niveaux sigma', '-', &
109           llm, nivsigs, zvertiid)           llm, nivsigs, zvertiid)
110      !      !
111      !  Appels a histdef pour la definition des variables a sauvegarder      ! Appels a histdef pour la definition des variables a sauvegarder
112      !      !
113      !  Vents U      ! Vents U
114      !      !
115      call histdef(fileid, 'ucov', 'vents u covariants', 'm/s', &      call histdef(fileid, 'ucov', 'vents u covariants', 'm/s', &
116           iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, &           iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, &
117           'inst(X)', t_ops, t_wrt)           'inst(X)', t_ops, t_wrt)
118      !      !
119      !  Vents V      ! Vents V
120      !      !
121      call histdef(filevid, 'vcov', 'vents v covariants', 'm/s', &      call histdef(filevid, 'vcov', 'vents v covariants', 'm/s', &
122           iip1, jjm, vhoriid, llm, 1, llm, zvertiid, &           iip1, jjm, vhoriid, llm, 1, llm, zvertiid, &
123           'inst(X)', t_ops, t_wrt)           'inst(X)', t_ops, t_wrt)
124    
125      !      !
126      !  Temperature potentielle      ! Temperature potentielle
127      !      !
128      call histdef(fileid, 'teta', 'temperature potentielle', '-', &      call histdef(fileid, 'teta', 'temperature potentielle', '-', &
129           iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &           iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
130           'inst(X)', t_ops, t_wrt)           'inst(X)', t_ops, t_wrt)
131      !      !
132      !  Geopotentiel      ! Geopotentiel
133      !      !
134      call histdef(fileid, 'phi', 'geopotentiel instantane', '-', &      call histdef(fileid, 'phi', 'geopotentiel instantane', '-', &
135           iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &           iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
136           'inst(X)', t_ops, t_wrt)           'inst(X)', t_ops, t_wrt)
137      !      !
138      !  Traceurs      ! Traceurs
139      !      !
140      DO iq=1,nq      DO iq=1, nq
141         call histdef(fileid, ttext(iq),  ttext(iq), '-', &         call histdef(fileid, ttext(iq), ttext(iq), '-', &
142              iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &              iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
143              'inst(X)', t_ops, t_wrt)              'inst(X)', t_ops, t_wrt)
144      enddo      enddo
145      !      !
146      !  Masse      ! Masse
147      !      !
148      call histdef(fileid, 'masse', 'masse', 'kg', &      call histdef(fileid, 'masse', 'masse', 'kg', &
149           iip1, jjp1, thoriid, 1, 1, 1, -99, &           iip1, jjp1, thoriid, 1, 1, 1, -99, &
150           'inst(X)', t_ops, t_wrt)           'inst(X)', t_ops, t_wrt)
151      !      !
152      !  Pression au sol      ! Pression au sol
153      !      !
154      call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &      call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &
155           iip1, jjp1, thoriid, 1, 1, 1, -99, &           iip1, jjp1, thoriid, 1, 1, 1, -99, &
156           'inst(X)', t_ops, t_wrt)           'inst(X)', t_ops, t_wrt)
157      !      !
158      !  Pression au sol      ! Pression au sol
159      !      !
160      call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &      call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &
161           iip1, jjp1, thoriid, 1, 1, 1, -99, &           iip1, jjp1, thoriid, 1, 1, 1, -99, &
162           'inst(X)', t_ops, t_wrt)           'inst(X)', t_ops, t_wrt)
163      !      !
164      !  Fin      ! Fin
165      !      !
166      call histend(fileid)      call histend(fileid)
167      call histend(filevid)      call histend(filevid)

Legend:
Removed from v.30  
changed lines
  Added in v.52

  ViewVC Help
Powered by ViewVC 1.1.21