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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 212 - (hide annotations)
Thu Jan 12 12:31:31 2017 UTC (7 years, 4 months ago) by guez
File size: 4094 byte(s)
Moved variables from module com_io_dyn to module inithist_m, where
they are defined.

Split grid_atob.f into grille_m.f and dist_sphe.f. Extracted ASCCI art
to documentation. In grille_m, use automatic arrays instead of maximum
size. In grille_m, instead of printing data for every problematic
point, print a single diagnostic message.

Removed variables top_height, overlap, lev_histhf, lev_histday,
lev_histmth, type_run, ok_isccp, ok_regdyn, lonmin_ins, lonmax_ins,
latmin_ins, latmax_ins of module clesphys, not used.

Removed variable itap of module histwrite_phy_m, not used. There is a
variable itap in module time_phylmdz.

Added output of tro3.

In physiq, no need to compute wo at every time-step, since we only use
it in radlwsw.

1 guez 3 module inithist_m
2    
3     implicit none
4    
5 guez 212 integer histid, histvid, histuid
6    
7 guez 3 contains
8    
9 guez 129 subroutine inithist(tstep, nq, 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 61 ! L. Fairhead, LMD, 03/99
13 guez 3
14 guez 61 ! Routine d'initialisation des écritures des fichiers histoires
15     ! LMDZ au format IOIPSL.
16 guez 3
17 guez 92 USE dimens_m, ONLY: jjm, llm
18     USE disvert_m, ONLY: presnivs
19 guez 139 use dynetat0_m, only: day_ref, annee_ref, rlatu, rlatv, rlonu, rlonv
20 guez 61 USE histbeg_totreg_m, ONLY : histbeg_totreg
21     USE histdef_m, ONLY : histdef
22     USE histend_m, ONLY : histend
23     USE histvert_m, ONLY : histvert
24 guez 92 USE iniadvtrac_m, ONLY: ttext
25     USE nr_util, ONLY: pi
26 guez 61 USE paramet_m, ONLY: iip1, jjp1
27     USE temps, ONLY: itau_dyn
28 guez 92 USE ymds2ju_m, ONLY: ymds2ju
29 guez 3
30 guez 61 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
42     !-----------------------------------------------------------------------
43    
44 guez 129 CALL ymds2ju(annee_ref, 1, day_ref, 0.0, zjulian)
45 guez 3
46     do jj = 1, jjp1
47     do ii = 1, iip1
48 guez 39 rlong(ii, jj) = rlonu(ii) * 180. / pi
49     rlat(ii, jj) = rlatu(jj) * 180. / pi
50 guez 3 enddo
51     enddo
52    
53 guez 56 call histbeg_totreg("dyn_histu.nc", rlong(:,1), rlat(1,:), 1, iip1, 1, &
54     jjp1, itau_dyn, zjulian, tstep, uhoriid, histuid)
55 guez 3
56     do jj = 1, jjm
57     do ii = 1, iip1
58 guez 39 rlong(ii, jj) = rlonv(ii) * 180. / pi
59     rlat(ii, jj) = rlatv(jj) * 180. / pi
60 guez 3 enddo
61     enddo
62    
63 guez 56 ! Creation du fichier histoire pour la grille en V (oblige pour l'instant,
64     ! IOIPSL ne permet pas de grilles avec des nombres de point differents dans
65     ! un meme fichier)
66    
67 guez 61 call histbeg_totreg('dyn_histv.nc', rlong(:, 1), rlat(1, :jjm), 1, iip1, &
68     1, jjm, itau_dyn, zjulian, tstep, vhoriid, histvid)
69    
70 guez 3 do jj = 1, jjp1
71     do ii = 1, iip1
72 guez 39 rlong(ii, jj) = rlonv(ii) * 180. / pi
73     rlat(ii, jj) = rlatu(jj) * 180. / pi
74 guez 3 enddo
75     enddo
76    
77 guez 56 call histbeg_totreg("dyn_hist.nc", rlong(:, 1), rlat(1, :), 1, iip1, 1, &
78     jjp1, itau_dyn, zjulian, tstep, thoriid, histid)
79 guez 61
80 guez 39 ! Appel a histvert pour la grille verticale
81 guez 61
82 guez 67 call histvert(histid, 'presnivs', 'Niveaux pression','mb', presnivs/100., &
83     zvertiid,'down')
84     call histvert(histvid, 'presnivs', 'Niveaux pression','mb', &
85 guez 61 presnivs/100., zvertiid,'down')
86 guez 67 call histvert(histuid, 'presnivs', 'Niveaux pression','mb', &
87 guez 61 presnivs/100., zvertiid,'down')
88    
89 guez 39 ! Appels a histdef pour la definition des variables a sauvegarder
90 guez 3
91 guez 61 call histdef(histuid, 'u', 'vent u', 'm/s', iip1, jjp1, uhoriid, llm, 1, &
92     llm, zvertiid, 'inst(X)', t_ops, t_wrt)
93     call histdef(histvid, 'v', 'vent v', 'm/s', iip1, jjm, vhoriid, llm, 1, &
94     llm, zvertiid, 'inst(X)', t_ops, t_wrt)
95     call histdef(histid, 'teta', 'temperature potentielle', '-', iip1, jjp1, &
96     thoriid, llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
97     call histdef(histid, 'phi', 'geopotentiel', '-', iip1, jjp1, thoriid, &
98     llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
99    
100 guez 39 ! Traceurs
101     DO iq=1, nq
102 guez 61 call histdef(histid, ttext(iq), ttext(iq), '-', iip1, jjp1, thoriid, &
103     llm, 1, llm, zvertiid, 'inst(X)', t_ops, t_wrt)
104 guez 3 enddo
105 guez 61
106     call histdef(histid, 'masse', 'masse', 'kg', iip1, jjp1, thoriid, llm, 1, &
107     llm, zvertiid, 'inst(X)', t_ops, t_wrt)
108     call histdef(histid, 'ps', 'pression naturelle au sol', 'Pa', iip1, jjp1, &
109     thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
110     call histdef(histid, 'phis', 'geopotentiel au sol', '-', iip1, jjp1, &
111     thoriid, 1, 1, 1, -99, 'inst(X)', t_ops, t_wrt)
112    
113 guez 56 call histend(histid)
114     call histend(histuid)
115     call histend(histvid)
116 guez 3
117     end subroutine inithist
118    
119     end module inithist_m

  ViewVC Help
Powered by ViewVC 1.1.21