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

Annotation of /trunk/Sources/bibio/inithist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 92 - (hide annotations)
Wed Mar 26 18:16:05 2014 UTC (10 years, 1 month ago) by guez
Original Path: trunk/bibio/inithist.f
File size: 4210 byte(s)
Extracted procedures that were in module calendar into separate files.

1 guez 3 module inithist_m
2    
3     implicit none
4    
5     contains
6    
7 guez 56 subroutine inithist(day0, anne0, tstep, nq, 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 61 ! L. Fairhead, LMD, 03/99
11 guez 3
12 guez 61 ! Routine d'initialisation des écritures des fichiers histoires
13     ! LMDZ au format IOIPSL.
14 guez 3
15 guez 92 USE comgeom, ONLY: rlatu, rlatv, rlonu, rlonv
16 guez 61 USE com_io_dyn, ONLY: histid, histuid, histvid
17 guez 92 USE dimens_m, ONLY: jjm, llm
18     USE disvert_m, ONLY: presnivs
19 guez 61 USE histbeg_totreg_m, ONLY : histbeg_totreg
20     USE histdef_m, ONLY : histdef
21     USE histend_m, ONLY : histend
22     USE histvert_m, ONLY : histvert
23 guez 92 USE iniadvtrac_m, ONLY: ttext
24     USE nr_util, ONLY: pi
25 guez 61 USE paramet_m, ONLY: iip1, jjp1
26     USE temps, ONLY: itau_dyn
27 guez 92 USE ymds2ju_m, ONLY: ymds2ju
28 guez 3
29 guez 61 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 guez 3
35 guez 61 ! Variables locales:
36 guez 3 real zjulian
37     integer iq
38 guez 39 real rlong(iip1, jjp1), rlat(iip1, jjp1)
39 guez 3 integer uhoriid, vhoriid, thoriid, zvertiid
40 guez 39 integer ii, jj
41 guez 3 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 guez 39 rlong(ii, jj) = rlonu(ii) * 180. / pi
52     rlat(ii, jj) = rlatu(jj) * 180. / pi
53 guez 3 enddo
54     enddo
55    
56 guez 56 call histbeg_totreg("dyn_histu.nc", rlong(:,1), rlat(1,:), 1, iip1, 1, &
57     jjp1, itau_dyn, zjulian, tstep, uhoriid, histuid)
58 guez 3
59     do jj = 1, jjm
60     do ii = 1, iip1
61 guez 39 rlong(ii, jj) = rlonv(ii) * 180. / pi
62     rlat(ii, jj) = rlatv(jj) * 180. / pi
63 guez 3 enddo
64     enddo
65    
66 guez 56 ! 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 guez 61 call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), 1, iip1, &
71     1, jjm, itau_dyn, zjulian, tstep, vhoriid, histvid)
72    
73 guez 3 do jj = 1, jjp1
74     do ii = 1, iip1
75 guez 39 rlong(ii, jj) = rlonv(ii) * 180. / pi
76     rlat(ii, jj) = rlatu(jj) * 180. / pi
77 guez 3 enddo
78     enddo
79    
80 guez 56 call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), 1, iip1, 1, &
81     jjp1, itau_dyn, zjulian, tstep, thoriid, histid)
82 guez 61
83 guez 39 ! Appel a histvert pour la grille verticale
84 guez 61
85 guez 67 call histvert(histid, 'presnivs', 'Niveaux pression','mb', presnivs/100., &
86     zvertiid,'down')
87     call histvert(histvid, 'presnivs', 'Niveaux pression','mb', &
88 guez 61 presnivs/100., zvertiid,'down')
89 guez 67 call histvert(histuid, 'presnivs', 'Niveaux pression','mb', &
90 guez 61 presnivs/100., zvertiid,'down')
91    
92 guez 39 ! Appels a histdef pour la definition des variables a sauvegarder
93 guez 3
94 guez 61 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 guez 39 ! Traceurs
104     DO iq=1, nq
105 guez 61 call histdef(histid, ttext(iq), ttext(iq), '-', iip1, jjp1, thoriid, &
106     llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
107 guez 3 enddo
108 guez 61
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 guez 56 call histend(histid)
117     call histend(histuid)
118     call histend(histvid)
119 guez 3
120     end subroutine inithist
121    
122     end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21