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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 52 - (hide annotations)
Fri Sep 23 12:28:01 2011 UTC (12 years, 7 months ago) by guez
File size: 4636 byte(s)
Split "conflx.f" into single-procedure files in directory "Conflx".

Split "cv_routines.f" into single-procedure files in directory
"CV_routines". Made module "cvparam" from included file
"cvparam.h". No included file other than "netcdf.inc" left in LMDZE.

1 guez 3 module inithist_m
2    
3     implicit none
4    
5     contains
6    
7 guez 39 subroutine inithist(day0, anne0, tstep, nq, fileid, filevid, t_ops, t_wrt)
8 guez 3
9 guez 39 ! From inithist.F, version 1.1.1.1 2004/05/19 12:53:05
10 guez 3
11 guez 39 ! Routine d'initialisation des écritures des fichiers histoires LMDZ
12     ! au format IOIPSL
13     ! Appels successifs des routines : histbeg, histhori, histver,
14     ! histdef, histend
15 guez 3
16 guez 39 ! 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 guez 3
23 guez 39 ! Sorties :
24     ! fileid : ID du fichier Netcdf créé
25     ! filevid : ID du fichier Netcdf pour la grille v
26 guez 3
27 guez 39 ! L. Fairhead, LMD, 03/99
28 guez 3
29 guez 30 USE calendar
30     use histcom
31 guez 3 use dimens_m
32     use paramet_m
33     use comconst
34     use comvert
35     use logic
36     use comgeom
37     use serre
38     use temps
39     use ener
40 guez 18 use iniadvtrac_m
41 guez 39 use nr_util, only: pi
42 guez 3
43 guez 39 ! Arguments
44 guez 3 integer day0, anne0
45     real, intent(in):: tstep, t_ops, t_wrt
46     integer fileid, filevid
47     integer nq
48    
49 guez 39 ! Variables locales
50 guez 3 real zjulian
51     integer iq
52 guez 39 real rlong(iip1, jjp1), rlat(iip1, jjp1)
53 guez 3 integer uhoriid, vhoriid, thoriid, zvertiid
54 guez 39 integer ii, jj
55 guez 3 integer zan, dayref
56    
57     !-----------------------------------------------------------------------
58    
59 guez 39 ! Appel a histbeg: creation du fichier netcdf et initialisations diverses
60 guez 3
61     zan = anne0
62     dayref = day0
63     CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
64    
65     do jj = 1, jjp1
66     do ii = 1, iip1
67 guez 39 rlong(ii, jj) = rlonu(ii) * 180. / pi
68     rlat(ii, jj) = rlatu(jj) * 180. / pi
69 guez 3 enddo
70     enddo
71    
72 guez 39 call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), &
73 guez 3 1, iip1, 1, jjp1, &
74 guez 27 itau_dyn, zjulian, tstep, uhoriid, fileid)
75 guez 3 !
76 guez 39 ! 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
78     ! un meme fichier)
79 guez 3
80     do jj = 1, jjm
81     do ii = 1, iip1
82 guez 39 rlong(ii, jj) = rlonv(ii) * 180. / pi
83     rlat(ii, jj) = rlatv(jj) * 180. / pi
84 guez 3 enddo
85     enddo
86    
87 guez 39 call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), &
88 guez 3 1, iip1, 1, jjm, &
89 guez 27 itau_dyn, zjulian, tstep, vhoriid, filevid)
90 guez 3 !
91 guez 39 ! Appel a histhori pour rajouter les autres grilles horizontales
92 guez 3 !
93     do jj = 1, jjp1
94     do ii = 1, iip1
95 guez 39 rlong(ii, jj) = rlonv(ii) * 180. / pi
96     rlat(ii, jj) = rlatu(jj) * 180. / pi
97 guez 3 enddo
98     enddo
99    
100 guez 15 call histhori_regular(fileid, iip1, rlong, jjp1, rlat, 'scalar', &
101 guez 3 'Grille points scalaires', thoriid)
102     !
103 guez 39 ! Appel a histvert pour la grille verticale
104 guez 3 !
105 guez 39 call histvert(fileid, 'sig_s', 'Niveaux sigma', '-', &
106 guez 3 llm, nivsigs, zvertiid)
107     ! Pour le fichier V
108 guez 39 call histvert(filevid, 'sig_s', 'Niveaux sigma', '-', &
109 guez 3 llm, nivsigs, zvertiid)
110     !
111 guez 39 ! Appels a histdef pour la definition des variables a sauvegarder
112 guez 3 !
113 guez 39 ! Vents U
114 guez 3 !
115     call histdef(fileid, 'ucov', 'vents u covariants', 'm/s', &
116     iip1, jjp1, uhoriid, llm, 1, llm, zvertiid, &
117 guez 15 'inst(X)', t_ops, t_wrt)
118 guez 3 !
119 guez 39 ! Vents V
120 guez 3 !
121     call histdef(filevid, 'vcov', 'vents v covariants', 'm/s', &
122     iip1, jjm, vhoriid, llm, 1, llm, zvertiid, &
123 guez 15 'inst(X)', t_ops, t_wrt)
124 guez 3
125     !
126 guez 39 ! Temperature potentielle
127 guez 3 !
128     call histdef(fileid, 'teta', 'temperature potentielle', '-', &
129     iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
130 guez 15 'inst(X)', t_ops, t_wrt)
131 guez 3 !
132 guez 39 ! Geopotentiel
133 guez 3 !
134     call histdef(fileid, 'phi', 'geopotentiel instantane', '-', &
135     iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
136 guez 15 'inst(X)', t_ops, t_wrt)
137 guez 3 !
138 guez 39 ! Traceurs
139 guez 3 !
140 guez 39 DO iq=1, nq
141     call histdef(fileid, ttext(iq), ttext(iq), '-', &
142 guez 3 iip1, jjp1, thoriid, llm, 1, llm, zvertiid, &
143 guez 15 'inst(X)', t_ops, t_wrt)
144 guez 3 enddo
145     !
146 guez 39 ! Masse
147 guez 3 !
148     call histdef(fileid, 'masse', 'masse', 'kg', &
149     iip1, jjp1, thoriid, 1, 1, 1, -99, &
150 guez 15 'inst(X)', t_ops, t_wrt)
151 guez 3 !
152 guez 39 ! Pression au sol
153 guez 3 !
154     call histdef(fileid, 'ps', 'pression naturelle au sol', 'Pa', &
155     iip1, jjp1, thoriid, 1, 1, 1, -99, &
156 guez 15 'inst(X)', t_ops, t_wrt)
157 guez 3 !
158 guez 39 ! Pression au sol
159 guez 3 !
160     call histdef(fileid, 'phis', 'geopotentiel au sol', '-', &
161     iip1, jjp1, thoriid, 1, 1, 1, -99, &
162 guez 15 'inst(X)', t_ops, t_wrt)
163 guez 3 !
164 guez 39 ! Fin
165 guez 3 !
166     call histend(fileid)
167     call histend(filevid)
168    
169     end subroutine inithist
170    
171     end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21