source: LMDZ5/branches/LMDZ6_rc0/libf/bibio/initdynav.F90 @ 5394

Last change on this file since 5394 was 1910, checked in by Laurent Fairhead, 11 years ago

Merged trunk changes r1860:1909 into testing branch

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 5.3 KB
Line 
1! $Id: initdynav.F90 1910 2013-11-29 08:40:25Z evignon $
2
3subroutine initdynav(day0,anne0,tstep,t_ops,t_wrt)
4
5#ifdef CPP_IOIPSL
6  USE IOIPSL
7#endif
8  USE infotrac, ONLY : nqtot, ttext
9  use com_io_dyn_mod, only : histaveid,histvaveid,histuaveid, &
10       dynhistave_file,dynhistvave_file,dynhistuave_file
11  implicit none
12
13
14  !   Routine d'initialisation des ecritures des fichiers histoires LMDZ
15  !   au format IOIPSL. Initialisation du fichier histoire moyenne.
16
17  !   Appels succesifs des routines: histbeg
18  !                                  histhori
19  !                                  histver
20  !                                  histdef
21  !                                  histend
22
23  !   Entree:
24
25  !      infile: nom du fichier histoire a creer
26  !      day0,anne0: date de reference
27  !      tstep : frequence d'ecriture
28  !      t_ops: frequence de l'operation pour IOIPSL
29  !      t_wrt: frequence d'ecriture sur le fichier
30
31
32  !   L. Fairhead, LMD, 03/99
33
34  include "dimensions.h"
35  include "paramet.h"
36  include "comconst.h"
37  include "comvert.h"
38  include "comgeom.h"
39  include "temps.h"
40  include "ener.h"
41  include "logic.h"
42  include "description.h"
43  include "serre.h"
44  include "iniprint.h"
45
46  !   Arguments
47
48  integer day0, anne0
49  real tstep, t_ops, t_wrt
50
51#ifdef CPP_IOIPSL
52  ! This routine needs IOIPSL to work
53  !   Variables locales
54
55  integer tau0
56  real zjulian
57  integer iq
58  real rlong(iip1,jjp1), rlat(iip1,jjp1)
59  integer uhoriid, vhoriid, thoriid, zvertiid
60  integer ii,jj
61  integer zan, dayref
62
63  !--------------------------------------------------------------------
64
65  !  Initialisations
66
67  pi = 4. * atan (1.)
68
69  !  Appel a histbeg: creation du fichier netcdf et initialisations diverses
70
71
72  zan = anne0
73  dayref = day0
74  CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
75  tau0 = itau_dyn
76
77  do jj = 1, jjp1
78     do ii = 1, iip1
79        rlong(ii,jj) = rlonv(ii) * 180. / pi
80        rlat(ii,jj)  = rlatu(jj) * 180. / pi
81     enddo
82  enddo
83
84  ! Creation de 3 fichiers pour les differentes grilles horizontales
85  ! Restriction de IOIPSL: seulement 2 coordonnees dans le meme fichier
86  ! Grille Scalaire       
87  call histbeg(dynhistave_file, iip1, rlong(:,1), jjp1, rlat(1,:), &
88       1, iip1, 1, jjp1, &
89       tau0, zjulian, tstep, thoriid,histaveid)
90
91  ! Creation du fichier histoire pour les grilles en V et U (oblige
92  ! pour l'instant, IOIPSL ne permet pas de grilles avec des nombres
93  ! de point differents dans  un meme fichier)
94  ! Grille V
95  do jj = 1, jjm
96     do ii = 1, iip1
97        rlong(ii,jj) = rlonv(ii) * 180. / pi
98        rlat(ii,jj) = rlatv(jj) * 180. / pi
99     enddo
100  enddo
101
102  call histbeg(dynhistvave_file, iip1, rlong(:,1), jjm, rlat(1,:), &
103       1, iip1, 1, jjm, &
104       tau0, zjulian, tstep, vhoriid,histvaveid)
105  ! Grille U
106  do jj = 1, jjp1
107     do ii = 1, iip1
108        rlong(ii,jj) = rlonu(ii) * 180. / pi
109        rlat(ii,jj) = rlatu(jj) * 180. / pi
110     enddo
111  enddo
112
113  call histbeg(dynhistuave_file, iip1, rlong(:,1),jjp1, rlat(1,:), &
114       1, iip1, 1, jjp1, &
115       tau0, zjulian, tstep, uhoriid,histuaveid)
116
117  !  Appel a histvert pour la grille verticale
118
119  call histvert(histaveid,'presnivs','Niveaux Pression approximatifs','mb', &
120       llm, presnivs/100., zvertiid,'down')
121  call histvert(histuaveid,'presnivs','Niveaux Pression approximatifs','mb', &
122       llm, presnivs/100., zvertiid,'down')
123  call histvert(histvaveid,'presnivs','Niveaux Pression approximatifs','mb', &
124       llm, presnivs/100., zvertiid,'down')
125
126  !  Appels a histdef pour la definition des variables a sauvegarder
127
128  !  Vents U
129
130  call histdef(histuaveid, 'u', 'vent u moyen ', &
131       'm/s', iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, &
132       32, 'ave(X)', t_ops, t_wrt)
133
134  !  Vents V
135
136  call histdef(histvaveid, 'v', 'vent v moyen', &
137       'm/s', iip1, jjm, vhoriid, llm, 1, llm, zvertiid, &
138       32, 'ave(X)', t_ops, t_wrt)
139
140
141  !  Temperature
142
143  call histdef(histaveid, 'temp', 'temperature moyenne', 'K', &
144       iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
145       32, 'ave(X)', t_ops, t_wrt)
146
147  !  Temperature potentielle
148
149  call histdef(histaveid, 'theta', 'temperature potentielle', 'K', &
150       iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
151       32, 'ave(X)', t_ops, t_wrt)
152
153  !  Geopotentiel
154
155  call histdef(histaveid, 'phi', 'geopotentiel moyen', '-', &
156       iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
157       32, 'ave(X)', t_ops, t_wrt)
158
159  !  Traceurs
160
161  !        DO iq=1,nqtot
162  !          call histdef(histaveid, ttext(iq), ttext(iq), '-', &
163  !                  iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
164  !                  32, 'ave(X)', t_ops, t_wrt)
165  !        enddo
166
167  !  Masse
168
169  call histdef(histaveid, 'masse', 'masse', 'kg', &
170       iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
171       32, 'ave(X)', t_ops, t_wrt)
172
173  !  Pression au sol
174
175  call histdef(histaveid, 'ps', 'pression naturelle au sol', 'Pa', &
176       iip1, jjp1, thoriid, 1, 1, 1, -99, &
177       32, 'ave(X)', t_ops, t_wrt)
178
179  !  Geopotentiel au sol
180
181  !      call histdef(histaveid, 'phis', 'geopotentiel au sol', '-', &
182  !                  iip1, jjp1, thoriid, 1, 1, 1, -99, &
183  !                  32, 'ave(X)', t_ops, t_wrt)
184
185  call histend(histaveid)
186  call histend(histuaveid)
187  call histend(histvaveid)
188#else
189  write(lunout,*)"initdynav: Warning this routine should not be", &
190       " used without ioipsl"
191#endif
192  ! of #ifdef CPP_IOIPSL
193
194end subroutine initdynav
Note: See TracBrowser for help on using the repository browser.