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

Annotation of /trunk/phylmd/ini_histrac.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 91 - (hide annotations)
Wed Mar 26 17:18:58 2014 UTC (10 years, 2 months ago) by guez
File size: 3329 byte(s)
Removed unused variables lock_startdate and time_stamp of module
calendar.

Noticed that physiq does not change the surface pressure. So removed
arguments ps and dpfi of subroutine addfi. dpfi was always 0. The
computation of ps in addfi included some averaging at the poles. In
principle, this does not change ps but in practice it does because of
finite numerical precision. So the results of the simulation are
changed. Removed arguments ps and dpfi of calfis. Removed argument
d_ps of physiq.

du at the poles is not computed by dudv1, so declare only the
corresponding latitudes in dudv1. caldyn passes only a section of the
array dudyn as argument.

Removed variable niadv of module iniadvtrac_m.

Declared arguments of exner_hyb as assumed-shape arrays and made all
other horizontal sizes in exner_hyb dynamic. This allows the external
program test_disvert to use exner_hyb at a single horizontal position.

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 91 USE histbeg_totreg_m, ONLY: histbeg_totreg
15 guez 61 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 guez 91 use iniadvtrac_m, only: tnom, ttext
20 guez 34 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 guez 91 integer it, iq
41 guez 34
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 guez 91 CALL histdef(nid_tra, tnom(iq), ttext(iq), "U/kga", iim, jjm+1, &
67 guez 34 nhori, llm, 1, llm, nvert, "ave(X)", zsto, zout)
68     if (lessivage) THEN
69 guez 91 CALL histdef(nid_tra, "fl"//tnom(iq), "Flux "//ttext(iq), &
70 guez 34 "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 guez 91 "tendance thermique"// ttext(iq), "?", &
77 guez 34 iim, jjm+1, nhori, llm, 1, llm, nvert, &
78     "ave(X)", zsto, zout)
79     CALL histdef(nid_tra, "d_tr_cv_"//tnom(iq), &
80 guez 91 "tendance convection"// ttext(iq), "?", &
81 guez 34 iim, jjm+1, nhori, llm, 1, llm, nvert, &
82     "ave(X)", zsto, zout)
83     CALL histdef(nid_tra, "d_tr_cl_"//tnom(iq), &
84 guez 91 "tendance couche limite"// ttext(iq), "?", &
85 guez 34 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