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

Annotation of /trunk/phylmd/ini_histrac.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 76 - (hide annotations)
Fri Nov 15 18:45:49 2013 UTC (10 years, 6 months ago) by guez
Original Path: trunk/phylmd/ini_histrac.f90
File size: 3368 byte(s)
Moved everything out of libf.
1 guez 34 module ini_histrac_m
2    
3     implicit none
4    
5     contains
6    
7 guez 61
8 guez 34 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 dimens_m, only: iim, jjm, llm
13     USE calendar, only: ymds2ju
14 guez 61 USE histbeg_totreg_m, ONLY : histbeg_totreg
15     USE histdef_m, ONLY : histdef
16     USE histend_m, ONLY : histend
17     USE histvert_m, ONLY : histvert
18 guez 34 use temps, only: annee_ref, day_ref, itau_phy
19     use iniadvtrac_m, only: niadv, tnom, ttext
20     use dimphy, only: klon
21     use clesphys, only: ecrit_tra
22     use grid_change, only: gr_phy_write_2d
23     use phyetat0_m, only: rlon, rlat
24 guez 66 use disvert_m, only: presnivs
25 guez 34
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, iiq
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 guez 67 CALL histvert(nid_tra, "presnivs", "Vertical levels", "mb", presnivs, nvert)
49 guez 34
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     iiq=niadv(iq)
67     CALL histdef(nid_tra, tnom(iq), ttext(iiq), "U/kga", iim, jjm+1, &
68     nhori, llm, 1, llm, nvert, "ave(X)", zsto, zout)
69     if (lessivage) THEN
70     CALL histdef(nid_tra, "fl"//tnom(iq), "Flux "//ttext(iiq), &
71     "U/m2/s", iim, jjm+1, nhori, llm, 1, llm, nvert, &
72     "ave(X)", zsto, zout)
73     endif
74    
75     !---Ajout Olivia
76     CALL histdef(nid_tra, "d_tr_th_"//tnom(iq), &
77     "tendance thermique"// ttext(iiq), "?", &
78     iim, jjm+1, nhori, llm, 1, llm, nvert, &
79     "ave(X)", zsto, zout)
80     CALL histdef(nid_tra, "d_tr_cv_"//tnom(iq), &
81     "tendance convection"// ttext(iiq), "?", &
82     iim, jjm+1, nhori, llm, 1, llm, nvert, &
83     "ave(X)", zsto, zout)
84     CALL histdef(nid_tra, "d_tr_cl_"//tnom(iq), &
85     "tendance couche limite"// ttext(iiq), "?", &
86     iim, jjm+1, nhori, llm, 1, llm, nvert, &
87     "ave(X)", zsto, zout)
88     !---fin Olivia
89    
90     ENDDO
91    
92     CALL histdef(nid_tra, "pplay", "", "-", &
93     iim, jjm+1, nhori, llm, 1, llm, nvert, &
94     "inst(X)", zout, zout)
95     CALL histdef(nid_tra, "T", "temperature", "K", iim, jjm+1, nhori, llm, &
96     1, llm, nvert, "inst(X)", zout, zout)
97    
98     CALL histend(nid_tra)
99    
100     end subroutine ini_histrac
101    
102     end module ini_histrac_m

  ViewVC Help
Powered by ViewVC 1.1.21