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

Legend:
Removed from v.15  
changed lines
  Added in v.57

  ViewVC Help
Powered by ViewVC 1.1.21