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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 15 - (hide annotations)
Fri Aug 1 15:24:12 2008 UTC (15 years, 9 months ago) by guez
File size: 4996 byte(s)
-- Minor modification of input/output:

Added variable "Sigma_O3_Royer" to "histday.nc". "ecrit_day" is not
modified in "physiq". Removed variables "pyu1", "pyv1", "ftsol1",
"ftsol2", "ftsol3", "ftsol4", "psrf1", "psrf2", "psrf3", "psrf4"
"mfu", "mfd", "en_u", "en_d", "de_d", "de_u", "coefh" from
"histrac.nc".

Variable "raz_date" of module "conf_gcm_m" has logical type instead of
integer type.

-- Should not change any result at run time:

Modified calls to "IOIPSL_Lionel" procedures because the interfaces of
these procedures have been simplified.

Changed name of variable in module "start_init_orog_m": "masque" to
"mask".

Created a module containing procedure "phyredem".

Removed arguments "punjours", "pdayref" and "ptimestep" of procedure
"iniphysiq".

Renamed procedure "gr_phy_write" to "gr_phy_write_2d". Created
procedure "gr_phy_write_3d".

Removed procedures "ini_undefstd", "moy_undefSTD", "calcul_STDlev",
"calcul_divers".

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

  ViewVC Help
Powered by ViewVC 1.1.21