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