/[lmdze]/trunk/bibio/inithist.f
ViewVC logotype

Contents of /trunk/bibio/inithist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 82 - (show annotations)
Wed Mar 5 14:57:53 2014 UTC (10 years, 2 months ago) by guez
File size: 4209 byte(s)
Changed all ".f90" suffixes to ".f".
1 module inithist_m
2
3 implicit none
4
5 contains
6
7 subroutine inithist(day0, anne0, tstep, nq, t_ops, t_wrt)
8
9 ! From inithist.F, version 1.1.1.1 2004/05/19 12:53:05
10 ! L. Fairhead, LMD, 03/99
11
12 ! Routine d'initialisation des écritures des fichiers histoires
13 ! LMDZ au format IOIPSL.
14
15 USE calendar, ONLY: ymds2ju
16 USE com_io_dyn, ONLY: histid, histuid, histvid
17 USE histbeg_totreg_m, ONLY : histbeg_totreg
18 USE histdef_m, ONLY : histdef
19 USE histend_m, ONLY : histend
20 USE histvert_m, ONLY : histvert
21 USE dimens_m, ONLY: jjm, llm
22 USE paramet_m, ONLY: iip1, jjp1
23 USE disvert_m, ONLY: presnivs
24 USE comgeom, ONLY: rlatu, rlatv, rlonu, rlonv
25 USE temps, ONLY: itau_dyn
26 USE iniadvtrac_m, ONLY: ttext
27 USE nr_util, ONLY: pi
28
29 integer, intent(in):: day0, anne0 ! date de référence
30 real, intent(in):: tstep ! durée du pas de temps en secondes
31 integer, intent(in):: nq ! nombre de traceurs
32 real, intent(in):: t_ops ! fréquence de l'opération pour IOIPSL
33 real, intent(in):: t_wrt ! fréquence d'écriture sur le fichier
34
35 ! Variables locales:
36 real zjulian
37 integer iq
38 real rlong(iip1, jjp1), rlat(iip1, jjp1)
39 integer uhoriid, vhoriid, thoriid, zvertiid
40 integer ii, jj
41 integer zan, dayref
42
43 !-----------------------------------------------------------------------
44
45 zan = anne0
46 dayref = day0
47 CALL ymds2ju(zan, 1, dayref, 0.0, zjulian)
48
49 do jj = 1, jjp1
50 do ii = 1, iip1
51 rlong(ii, jj) = rlonu(ii) * 180. / pi
52 rlat(ii, jj) = rlatu(jj) * 180. / pi
53 enddo
54 enddo
55
56 call histbeg_totreg("dyn_histu.nc", rlong(:,1), rlat(1,:), 1, iip1, 1, &
57 jjp1, itau_dyn, zjulian, tstep, uhoriid, histuid)
58
59 do jj = 1, jjm
60 do ii = 1, iip1
61 rlong(ii, jj) = rlonv(ii) * 180. / pi
62 rlat(ii, jj) = rlatv(jj) * 180. / pi
63 enddo
64 enddo
65
66 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
67 ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
68 ! un meme fichier)
69
70 call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), 1, iip1, &
71 1, jjm, itau_dyn, zjulian, tstep, vhoriid, histvid)
72
73 do jj = 1, jjp1
74 do ii = 1, iip1
75 rlong(ii, jj) = rlonv(ii) * 180. / pi
76 rlat(ii, jj) = rlatu(jj) * 180. / pi
77 enddo
78 enddo
79
80 call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), 1, iip1, 1, &
81 jjp1, itau_dyn, zjulian, tstep, thoriid, histid)
82
83 ! Appel a histvert pour la grille verticale
84
85 call histvert(histid, 'presnivs', 'Niveaux pression','mb', presnivs/100., &
86 zvertiid,'down')
87 call histvert(histvid, 'presnivs', 'Niveaux pression','mb', &
88 presnivs/100., zvertiid,'down')
89 call histvert(histuid, 'presnivs', 'Niveaux pression','mb', &
90 presnivs/100., zvertiid,'down')
91
92 ! Appels a histdef pour la definition des variables a sauvegarder
93
94 call histdef(histuid, 'u', 'vent u', 'm/s', iip1, jjp1, uhoriid, llm, 1, &
95 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
96 call histdef(histvid, 'v', 'vent v', 'm/s', iip1, jjm, vhoriid, llm, 1, &
97 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
98 call histdef(histid, 'teta', 'temperature potentielle', '-', iip1, jjp1, &
99 thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
100 call histdef(histid, 'phi', 'geopotentiel', '-', iip1, jjp1, thoriid, &
101 llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
102
103 ! Traceurs
104 DO iq=1, nq
105 call histdef(histid, ttext(iq), ttext(iq), '-', iip1, jjp1, thoriid, &
106 llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
107 enddo
108
109 call histdef(histid, 'masse', 'masse', 'kg', iip1, jjp1, thoriid, llm, 1, &
110 llm, zvertiid, 'inst(X)', t_ops, t_wrt)
111 call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', iip1, jjp1, &
112 thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
113 call histdef(histid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, &
114 thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
115
116 call histend(histid)
117 call histend(histuid)
118 call histend(histvid)
119
120 end subroutine inithist
121
122 end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21