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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 39 - (hide annotations)
Tue Jan 25 15:11:05 2011 UTC (13 years, 3 months ago) by guez
File size: 4708 byte(s)
"pi" comes from "nr_util". Removed subroutine "initialize" in module
"comconst".

Copied the content of "fxy_sin.h" into "fxysinus", instead of getting
it from an "include" line. Removed file "fxy_sin.h".

"ps" has rank 2 in "gcm" and "dynetat0".

Assumed-shape for argument "q" of "integrd".

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

  ViewVC Help
Powered by ViewVC 1.1.21