source: LMDZ5/branches/IPSLCM5A2.1/libf/dyn3d_common/initdynav.F90 @ 3183

Last change on this file since 3183 was 2603, checked in by Ehouarn Millour, 8 years ago

Cleanup in the dynamics: turn logic.h into module logic_mod.F90
EM

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