/[lmdze]/trunk/libf/bibio/writehist.f
ViewVC logotype

Contents of /trunk/libf/bibio/writehist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 56 - (show annotations)
Tue Jan 10 19:02:02 2012 UTC (12 years, 4 months ago) by guez
File size: 2612 byte(s)
Imported "writehist.f" from LMDZ.

Moved module variable "histaveid" from "com_io_dyn" to "initdynav_m".

In "inithist", access directly module variables from "com_io_dyn"
instead of going through the arguments. Copying from LMDZ, write "u"
and scalar variables to separate files. Create a new variable for the
new file in "com_io_dyn". Copying from LMDZ, change the vertical axes
of the three files.

Removed some useless initializations in "dissip".

In "bilan_dyn", removed useless variable "time". Avoiding the
approximate test on "dt_cum" being a multiple of "dt_app", just
compute "ncum" from known usage of "bilan_dyn" and compute "dt_cum"
from "ncum". Change "periodav" from real to integer in
"conf_gcm_m". Since "day_step" is required to be a multiple of
"iperiod", so is "ncum".

1 !
2 ! $Id: writehist.F 1403 2010-07-01 09:02:53Z fairhead $
3 !
4 subroutine writehist(time,vcov,ucov,teta,phi,q,masse,ps,phis)
5
6 use dimens_m, only: nqmx, llm, jjm
7 USE iniadvtrac_m, ONLY: ttext
8 use com_io_dyn, only: histid,histvid,histuid
9 use paramet_m, only: ip1jm, ip1jmp1, iip1, jjp1
10 use temps, only: itau_dyn
11 use histwrite_m, only: histwrite
12 use histcom, only: histsync
13
14 implicit none
15
16 C
17 C Ecriture du fichier histoire au format IOIPSL
18 C
19 C Appels succesifs des routines: histwrite
20 C
21 C Entree:
22 C time: temps de l'ecriture
23 C vcov: vents v covariants
24 C ucov: vents u covariants
25 C teta: temperature potentielle
26 C phi : geopotentiel instantane
27 C q : traceurs
28 C masse: masse
29 C ps :pression au sol
30 C phis : geopotentiel au sol
31 C
32 C
33 C L. Fairhead, LMD, 03/99
34 C
35 C =====================================================================
36 C
37 C Declarations
38
39 C
40 C Arguments
41 C
42
43 REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm)
44 REAL teta(ip1jmp1,llm),phi(ip1jmp1,llm)
45 REAL ps(ip1jmp1),masse(ip1jmp1,llm)
46 REAL phis(ip1jmp1)
47 REAL q(ip1jmp1,llm,nqmx)
48 integer time
49
50
51 ! This routine needs IOIPSL to work
52 C Variables locales
53 C
54 integer iq, ii, ll
55 integer ndexu(ip1jmp1*llm),ndexv(ip1jm*llm),ndex2d(ip1jmp1)
56 logical ok_sync
57 integer itau_w
58 REAL vnat(ip1jm,llm),unat(ip1jmp1,llm)
59
60 C
61 C Initialisations
62 C
63 ndexu = 0
64 ndexv = 0
65 ndex2d = 0
66 ok_sync =.TRUE.
67 itau_w = itau_dyn + time
68 ! Passage aux composantes naturelles du vent
69 call covnat(llm, ucov, vcov, unat, vnat)
70 C
71 C Appels a histwrite pour l'ecriture des variables a sauvegarder
72 C
73 C Vents U
74 C
75 call histwrite(histuid, 'u', itau_w, unat)
76 C
77 C Vents V
78 C
79 call histwrite(histvid, 'v', itau_w, vnat)
80
81 C
82 C Temperature potentielle
83 C
84 call histwrite(histid, 'teta', itau_w, teta)
85 C
86 C Geopotentiel
87 C
88 call histwrite(histid, 'phi', itau_w, phi)
89 C
90 C Traceurs
91 C
92 ! DO iq=1,nqmx
93 ! call histwrite(histid, ttext(iq), itau_w, q(:,:,iq),
94 ! . iip1*jjp1*llm, ndexu)
95 ! enddo
96 !C
97 C Masse
98 C
99 call histwrite(histid,'masse',itau_w, masse)
100 C
101 C Pression au sol
102 C
103 call histwrite(histid, 'ps', itau_w, ps)
104 C
105 C Geopotentiel au sol
106 C
107 ! call histwrite(histid, 'phis', itau_w, phis, iip1*jjp1, ndex2d)
108 C
109 C Fin
110 C
111 if (ok_sync) then
112 call histsync(histid)
113 call histsync(histvid)
114 call histsync(histuid)
115 endif
116 return
117 end

  ViewVC Help
Powered by ViewVC 1.1.21