1 |
subroutine writedynav(histid, nq, time, vcov, ucov, teta, ppk, phi, q, masse, & |
2 |
ps, phis) |
3 |
|
4 |
! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05 |
5 |
! Ecriture du fichier histoire au format IOIPSL |
6 |
|
7 |
! Appels successifs des routines histwrite |
8 |
|
9 |
! Entree: |
10 |
! histid: ID du fichier histoire |
11 |
! time: temps de l'ecriture |
12 |
! vcov: vents v covariants |
13 |
! ucov: vents u covariants |
14 |
! phi : geopotentiel instantane |
15 |
! q : traceurs |
16 |
! ps :pression au sol |
17 |
! phis : geopotentiel au sol |
18 |
|
19 |
! L. Fairhead, LMD, 03/99 |
20 |
|
21 |
USE histwrite_m, ONLY : histwrite |
22 |
USE histcom, ONLY : histsync |
23 |
USE dimens_m, ONLY : llm |
24 |
USE paramet_m, ONLY : iip1, ijp1llm, ip1jm, ip1jmp1, jjp1 |
25 |
USE comconst, ONLY : cpp |
26 |
USE temps, ONLY : itau_dyn |
27 |
USE iniadvtrac_m, ONLY : ttext |
28 |
|
29 |
implicit none |
30 |
|
31 |
INTEGER histid, nq |
32 |
REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm) |
33 |
REAL, intent(in):: teta(ip1jmp1*llm) ! temperature potentielle |
34 |
real phi(ip1jmp1, llm), ppk(ip1jmp1*llm) |
35 |
REAL ps(ip1jmp1) |
36 |
real, intent(in):: masse(ip1jmp1, llm) |
37 |
REAL phis(ip1jmp1) |
38 |
REAL q(ip1jmp1, llm, nq) |
39 |
integer, intent(in):: time |
40 |
|
41 |
! Variables locales |
42 |
integer ndex2d(iip1*jjp1), ndex3d(iip1*jjp1*llm), iq, ii, ll |
43 |
real us(ip1jmp1*llm), vs(ip1jmp1*llm) |
44 |
real tm(ip1jmp1*llm) |
45 |
REAL vnat(ip1jm, llm), unat(ip1jmp1, llm) |
46 |
logical ok_sync |
47 |
integer itau_w |
48 |
|
49 |
!--------------------------------------------------------------- |
50 |
|
51 |
! Initialisations |
52 |
ndex3d = 0 |
53 |
ndex2d = 0 |
54 |
ok_sync = .TRUE. |
55 |
us = 999.999 |
56 |
vs = 999.999 |
57 |
tm = 999.999 |
58 |
vnat = 999.999 |
59 |
unat = 999.999 |
60 |
itau_w = itau_dyn + time |
61 |
|
62 |
! Passage aux composantes naturelles du vent |
63 |
call covnat(llm, ucov, vcov, unat, vnat) |
64 |
|
65 |
! Appels a histwrite pour l'ecriture des variables a sauvegarder |
66 |
|
67 |
! Vents U scalaire |
68 |
call gr_u_scal(llm, unat, us) |
69 |
call histwrite(histid, 'u', itau_w, us) |
70 |
|
71 |
! Vents V scalaire |
72 |
call gr_v_scal(llm, vnat, vs) |
73 |
call histwrite(histid, 'v', itau_w, vs) |
74 |
|
75 |
! Temperature potentielle moyennee |
76 |
call histwrite(histid, 'theta', itau_w, teta) |
77 |
|
78 |
! Temperature moyennee |
79 |
do ii = 1, ijp1llm |
80 |
tm(ii) = teta(ii) * ppk(ii)/cpp |
81 |
enddo |
82 |
call histwrite(histid, 'temp', itau_w, tm) |
83 |
|
84 |
! Geopotentiel |
85 |
call histwrite(histid, 'phi', itau_w, phi) |
86 |
|
87 |
! Traceurs |
88 |
DO iq=1, nq |
89 |
call histwrite(histid, ttext(iq), itau_w, q(:, :, iq)) |
90 |
enddo |
91 |
|
92 |
! Masse |
93 |
call histwrite(histid, 'masse', itau_w, masse) |
94 |
|
95 |
! Pression au sol |
96 |
call histwrite(histid, 'ps', itau_w, ps) |
97 |
|
98 |
! Geopotentiel au sol |
99 |
call histwrite(histid, 'phis', itau_w, phis) |
100 |
|
101 |
if (ok_sync) call histsync(histid) |
102 |
|
103 |
end subroutine writedynav |