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