1 |
guez |
34 |
module ini_histrac_m |
2 |
|
|
|
3 |
|
|
implicit none |
4 |
|
|
|
5 |
|
|
contains |
6 |
|
|
|
7 |
|
|
subroutine ini_histrac(nid_tra, pdtphys, nq_phys, lessivage) |
8 |
|
|
|
9 |
|
|
! From phylmd/ini_histrac.h, version 1.10 2006/02/21 08:08:30 |
10 |
|
|
|
11 |
|
|
use dimens_m, only: iim, jjm, llm |
12 |
|
|
USE calendar, only: ymds2ju |
13 |
|
|
use histcom, only: histbeg_totreg, histvert, histend, histdef |
14 |
|
|
use temps, only: annee_ref, day_ref, itau_phy |
15 |
|
|
use iniadvtrac_m, only: niadv, tnom, ttext |
16 |
|
|
use dimphy, only: klon |
17 |
|
|
use clesphys, only: ecrit_tra |
18 |
|
|
use grid_change, only: gr_phy_write_2d |
19 |
|
|
use phyetat0_m, only: rlon, rlat |
20 |
|
|
use comvert, only: presnivs |
21 |
|
|
|
22 |
|
|
INTEGER, intent(out):: nid_tra |
23 |
|
|
real, intent(in):: pdtphys ! pas d'integration pour la physique (s) |
24 |
|
|
|
25 |
|
|
integer, intent(in):: nq_phys |
26 |
|
|
! (nombre de traceurs auxquels on applique la physique) |
27 |
|
|
|
28 |
|
|
logical, intent(in):: lessivage |
29 |
|
|
|
30 |
|
|
! Variables local to the procedure: |
31 |
|
|
|
32 |
|
|
REAL zjulian |
33 |
|
|
REAL zx_lat(iim, jjm+1) |
34 |
|
|
INTEGER nhori, nvert |
35 |
|
|
REAL zsto, zout |
36 |
|
|
integer it, iq, iiq |
37 |
|
|
|
38 |
|
|
!--------------------------------------------------------- |
39 |
|
|
|
40 |
|
|
CALL ymds2ju(annee_ref, month=1, day=day_ref, sec=0.0, julian=zjulian) |
41 |
|
|
zx_lat(:, :) = gr_phy_write_2d(rlat) |
42 |
|
|
CALL histbeg_totreg("histrac", rlon(2:iim+1), zx_lat(1, :), & |
43 |
|
|
1, iim, 1, jjm+1, itau_phy, zjulian, pdtphys, nhori, nid_tra) |
44 |
|
|
CALL histvert(nid_tra, "presnivs", "Vertical levels", "mb", llm, & |
45 |
|
|
presnivs, nvert) |
46 |
|
|
|
47 |
|
|
zsto = pdtphys |
48 |
|
|
zout = pdtphys * REAL(ecrit_tra) |
49 |
|
|
|
50 |
|
|
CALL histdef(nid_tra, "phis", "Surface geop. height", "-", & |
51 |
|
|
iim, jjm+1, nhori, 1, 1, 1, -99, & |
52 |
|
|
"once", zsto, zout) |
53 |
|
|
CALL histdef(nid_tra, "aire", "Grid area", "-", & |
54 |
|
|
iim, jjm+1, nhori, 1, 1, 1, -99, & |
55 |
|
|
"once", zsto, zout) |
56 |
|
|
CALL histdef(nid_tra, "zmasse", "column density of air in cell", & |
57 |
|
|
"kg m-2", iim, jjm + 1, nhori, llm, 1, llm, nvert, "ave(X)", & |
58 |
|
|
zsto, zout) |
59 |
|
|
|
60 |
|
|
DO it = 1, nq_phys |
61 |
|
|
! champ 2D |
62 |
|
|
iq=it+2 |
63 |
|
|
iiq=niadv(iq) |
64 |
|
|
CALL histdef(nid_tra, tnom(iq), ttext(iiq), "U/kga", iim, jjm+1, & |
65 |
|
|
nhori, llm, 1, llm, nvert, "ave(X)", zsto, zout) |
66 |
|
|
if (lessivage) THEN |
67 |
|
|
CALL histdef(nid_tra, "fl"//tnom(iq), "Flux "//ttext(iiq), & |
68 |
|
|
"U/m2/s", iim, jjm+1, nhori, llm, 1, llm, nvert, & |
69 |
|
|
"ave(X)", zsto, zout) |
70 |
|
|
endif |
71 |
|
|
|
72 |
|
|
!---Ajout Olivia |
73 |
|
|
CALL histdef(nid_tra, "d_tr_th_"//tnom(iq), & |
74 |
|
|
"tendance thermique"// ttext(iiq), "?", & |
75 |
|
|
iim, jjm+1, nhori, llm, 1, llm, nvert, & |
76 |
|
|
"ave(X)", zsto, zout) |
77 |
|
|
CALL histdef(nid_tra, "d_tr_cv_"//tnom(iq), & |
78 |
|
|
"tendance convection"// ttext(iiq), "?", & |
79 |
|
|
iim, jjm+1, nhori, llm, 1, llm, nvert, & |
80 |
|
|
"ave(X)", zsto, zout) |
81 |
|
|
CALL histdef(nid_tra, "d_tr_cl_"//tnom(iq), & |
82 |
|
|
"tendance couche limite"// ttext(iiq), "?", & |
83 |
|
|
iim, jjm+1, nhori, llm, 1, llm, nvert, & |
84 |
|
|
"ave(X)", zsto, zout) |
85 |
|
|
!---fin Olivia |
86 |
|
|
|
87 |
|
|
ENDDO |
88 |
|
|
|
89 |
|
|
CALL histdef(nid_tra, "pplay", "", "-", & |
90 |
|
|
iim, jjm+1, nhori, llm, 1, llm, nvert, & |
91 |
|
|
"inst(X)", zout, zout) |
92 |
|
|
CALL histdef(nid_tra, "T", "temperature", "K", iim, jjm+1, nhori, llm, & |
93 |
|
|
1, llm, nvert, "inst(X)", zout, zout) |
94 |
|
|
|
95 |
|
|
CALL histend(nid_tra) |
96 |
|
|
|
97 |
|
|
end subroutine ini_histrac |
98 |
|
|
|
99 |
|
|
end module ini_histrac_m |