/[lmdze]/trunk/phylmd/ini_histrac.f
ViewVC logotype

Contents of /trunk/phylmd/ini_histrac.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 92 - (show annotations)
Wed Mar 26 18:16:05 2014 UTC (10 years, 1 month ago) by guez
File size: 3330 byte(s)
Extracted procedures that were in module calendar into separate files.

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

  ViewVC Help
Powered by ViewVC 1.1.21