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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 138 - (show annotations)
Fri May 22 23:13:19 2015 UTC (8 years, 11 months ago) by guez
File size: 2134 byte(s)
Moved variable nb_files from module histcom_var to module
histbeg_totreg_m.

Removed unused argument q of writehist.

No history file is created in program ce0l so there is no need to call
histclo in etat0.

In phyredem, access variables rlat and rlon directly from module
phyetat0_m instead of having them as arguments. This is clearer for
the program gcm. There are bad side effects for the program ce0l: we
have to modify the module variables rlat and rlon in procedure etat0,
and we need the additional file phyetat0.f to compile ce0l.

1 module writehist_m
2
3 implicit none
4
5 contains
6
7 subroutine writehist(time, vcov, ucov, teta, phi, masse, ps)
8
9 ! From writehist.F 1403 2010-07-01 09:02:53Z
10 ! Écriture du fichier histoire au format IOIPSL
11 ! Appels successifs des routines histwrite
12 ! L. Fairhead, LMD, 03/99
13
14 use dimens_m, only: nqmx, llm, jjm
15 use com_io_dyn, only: histid, histvid, histuid
16 use paramet_m, only: ip1jm, ip1jmp1, iip1, jjp1
17 use temps, only: itau_dyn
18 use histwrite_m, only: histwrite
19 use histsync_m, only: histsync
20 use covnat_m, only: covnat
21
22 ! Entree:
23 ! time: temps de l'ecriture
24 ! vcov: vents v covariants
25 ! ucov: vents u covariants
26 ! teta: temperature potentielle
27 ! phi : geopotentiel instantane
28 ! masse: masse
29 ! ps :pression au sol
30
31 ! Arguments
32
33 REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm)
34 REAL teta(ip1jmp1, llm), phi(ip1jmp1, llm)
35 REAL, intent(in):: ps(ip1jmp1), masse(ip1jmp1, llm)
36 integer time
37
38 ! This routine needs IOIPSL to work
39 ! Variables locales
40
41 integer ndexu(ip1jmp1*llm), ndexv(ip1jm*llm), ndex2d(ip1jmp1)
42 logical ok_sync
43 integer itau_w
44 REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
45
46 !---------------------------------------------------------------------
47
48 ! Initialisations
49
50 ndexu = 0
51 ndexv = 0
52 ndex2d = 0
53 ok_sync =.TRUE.
54 itau_w = itau_dyn + time
55 ! Passage aux composantes naturelles du vent
56 call covnat(llm, ucov, vcov, unat, vnat)
57
58 ! Appels a histwrite pour l'ecriture des variables a sauvegarder
59
60 ! Vents U
61
62 call histwrite(histuid, 'u', itau_w, unat)
63
64 ! Vents V
65
66 call histwrite(histvid, 'v', itau_w, vnat)
67
68 ! Temperature potentielle
69
70 call histwrite(histid, 'teta', itau_w, teta)
71
72 ! Geopotentiel
73
74 call histwrite(histid, 'phi', itau_w, phi)
75
76 ! Masse
77
78 call histwrite(histid, 'masse', itau_w, masse)
79
80 ! Pression au sol
81 call histwrite(histid, 'ps', itau_w, ps)
82
83 if (ok_sync) then
84 call histsync(histid)
85 call histsync(histvid)
86 call histsync(histuid)
87 endif
88
89 end subroutine writehist
90
91 end module writehist_m

  ViewVC Help
Powered by ViewVC 1.1.21