1 |
module writedynav_m |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
subroutine writedynav(vcov, ucov, teta, pk, 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 comconst, ONLY: cpp |
16 |
use covnat_m, only: covnat |
17 |
USE dimens_m, ONLY: iim, jjm, llm, nqmx |
18 |
USE histsync_m, ONLY: histsync |
19 |
USE histwrite_m, ONLY: histwrite |
20 |
USE iniadvtrac_m, ONLY: ttext |
21 |
use initdynav_m, only: histaveid |
22 |
use nr_util, only: assert |
23 |
USE paramet_m, ONLY: iip1, ip1jm, ip1jmp1, jjp1 |
24 |
USE temps, ONLY: itau_dyn |
25 |
|
26 |
! Vents covariants : |
27 |
REAL, intent(in):: vcov(:, :, :) ! (iim + 1, jjm, llm) |
28 |
REAL, intent(in):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm) |
29 |
|
30 |
REAL, intent(in):: teta(:, :, :) ! (iim + 1, jjm + 1, llm) |
31 |
! temperature potentielle |
32 |
|
33 |
real, intent(in):: phi(:, :, :) ! (iim + 1, jjm + 1, llm) |
34 |
! geopotentiel instantane |
35 |
|
36 |
real, intent(in):: pk(:, :, :) ! (iim + 1, jjm + 1, llm) |
37 |
REAL, intent(in):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol |
38 |
real, intent(in):: masse(:, :, :) ! (iim + 1, jjm + 1, llm) |
39 |
REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) geopotentiel au sol |
40 |
REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) traceurs |
41 |
integer, intent(in):: time ! temps de l'ecriture |
42 |
|
43 |
! Variables locales |
44 |
integer iq |
45 |
real us(ip1jmp1*llm), vs(ip1jmp1*llm) |
46 |
REAL vnat(ip1jm, llm), unat(ip1jmp1, llm) |
47 |
integer itau_w |
48 |
|
49 |
!--------------------------------------------------------------- |
50 |
|
51 |
call assert((/size(vcov, 1), size(ucov, 1), size(teta, 1), size(phi, 1), & |
52 |
size(pk, 1), size(ps, 1), size(masse, 1), size(phis, 1), & |
53 |
size(q, 1)/) == iim + 1, "writedynav iim") |
54 |
call assert((/size(vcov, 2) + 1, size(ucov, 2), size(teta, 2), & |
55 |
size(phi, 2), size(pk, 2), size(ps, 2), size(masse, 2), & |
56 |
size(phis, 2), size(q, 2)/) == jjm + 1, "writedynav jjm") |
57 |
call assert((/size(vcov, 3), size(ucov, 3), size(teta, 3), size(phi, 3), & |
58 |
size(pk, 3), size(masse, 3), size(q, 3)/) == llm, "writedynav llm") |
59 |
call assert(size(q, 4) == nqmx, "writedynav nqmx") |
60 |
|
61 |
! Initialisations |
62 |
us = 999.999 |
63 |
vs = 999.999 |
64 |
vnat = 999.999 |
65 |
unat = 999.999 |
66 |
itau_w = itau_dyn + time |
67 |
|
68 |
! Passage aux composantes naturelles du vent |
69 |
call covnat(llm, ucov, vcov, unat, vnat) |
70 |
|
71 |
! Appels a histwrite pour l'ecriture des variables a sauvegarder |
72 |
|
73 |
! Vents U scalaire |
74 |
call gr_u_scal(llm, unat, us) |
75 |
call histwrite(histaveid, 'u', itau_w, us) |
76 |
|
77 |
! Vents V scalaire |
78 |
call gr_v_scal(llm, vnat, vs) |
79 |
call histwrite(histaveid, 'v', itau_w, vs) |
80 |
|
81 |
! Temperature potentielle moyennee |
82 |
call histwrite(histaveid, 'theta', itau_w, teta) |
83 |
|
84 |
! Temperature moyennee |
85 |
call histwrite(histaveid, 'temp', itau_w, teta * pk / cpp) |
86 |
|
87 |
call histwrite(histaveid, 'phi', itau_w, phi) |
88 |
|
89 |
! Traceurs |
90 |
DO iq = 1, size(q, 4) |
91 |
call histwrite(histaveid, ttext(iq), itau_w, q(:, :, :, iq)) |
92 |
enddo |
93 |
|
94 |
call histwrite(histaveid, 'masse', itau_w, masse) |
95 |
call histwrite(histaveid, 'ps', itau_w, ps) |
96 |
call histwrite(histaveid, 'phis', itau_w, phis) |
97 |
|
98 |
call histsync(histaveid) |
99 |
|
100 |
end subroutine writedynav |
101 |
|
102 |
end module writedynav_m |