/[lmdze]/trunk/libf/phylmd/ini_histhf3d.f90
ViewVC logotype

Annotation of /trunk/libf/phylmd/ini_histhf3d.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 34 - (hide annotations)
Wed Jun 2 11:01:12 2010 UTC (14 years ago) by guez
File size: 2289 byte(s)
Split "ini_hist.f90" into single-procedure files.

In "calfis" and "physiq", removed dummy argument "nq" since "nq" must
be equal to "nqmx".

In "calfis", renamed dummy argument "pq" to "q", same name as actual
argument in "leapfrog". Renamed local variable "zqfi" to "qx", same
name as dummy argument in "physiq".

Removed arguments "itop_con" and "ibas_con" of "phytrac", which were
not used.

1 guez 34 module ini_histhf3d_m
2    
3     implicit none
4    
5     contains
6    
7     subroutine ini_histhf3d(dtime, nid_hf3d)
8    
9     ! From phylmd/ini_histhf3d.h, v 1.2 2005/05/25 13:10:09
10    
11     ! sorties hf 3d
12    
13     use dimens_m, only: iim, jjm, llm
14     use dimphy, only: klon, nbtr
15     use temps, only: itau_phy, day_ref, annee_ref
16     use clesphys, only: ecrit_hf
17     use phyetat0_m, only: rlon, rlat
18     USE calendar, only: ymds2ju
19     use histcom, only: histbeg_totreg, histvert, histend, histdef
20     use comvert, only: presnivs
21    
22     REAL, intent(in):: dtime ! pas temporel de la physique (s)
23     integer, intent(out):: nid_hf3d
24    
25     real zstohf, zout
26     REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1)
27     real zjulian
28     integer i, nhori, nvert, idayref
29    
30     !------------------------------------------
31    
32     zstohf = dtime * REAL(ecrit_hf)
33     zout = dtime * REAL(ecrit_hf)
34    
35     idayref = day_ref
36     CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
37    
38     CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon)
39     DO i = 1, iim
40     zx_lon(i, 1) = rlon(i+1)
41     zx_lon(i, (jjm + 1)) = rlon(i+1)
42     ENDDO
43    
44     CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat)
45     CALL histbeg_totreg("histhf3d", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, &
46     (jjm + 1), itau_phy, zjulian, dtime, nhori, nid_hf3d)
47    
48     CALL histvert(nid_hf3d, "presnivs", "Vertical levels", "mb", &
49     llm, presnivs/100., nvert)
50    
51     ! Champs 3D:
52    
53     CALL histdef(nid_hf3d, "temp", "Air temperature", "K", &
54     iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
55     "ave(X)", zstohf, zout)
56    
57     CALL histdef(nid_hf3d, "ovap", "Specific humidity", "kg/kg", &
58     iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
59     "ave(X)", zstohf, zout)
60    
61     CALL histdef(nid_hf3d, "vitu", "Zonal wind", "m/s", &
62     iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
63     "ave(X)", zstohf, zout)
64    
65     CALL histdef(nid_hf3d, "vitv", "Meridional wind", "m/s", &
66     iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
67     "ave(X)", zstohf, zout)
68    
69     if (nbtr >= 3) then
70     CALL histdef(nid_hf3d, "O3", "Ozone mass fraction", "?", iim, &
71     (jjm + 1), nhori, llm, 1, llm, nvert, "ave(X)", zstohf, &
72     zout)
73     end if
74    
75     CALL histend(nid_hf3d)
76    
77     end subroutine ini_histhf3d
78    
79     end module ini_histhf3d_m

  ViewVC Help
Powered by ViewVC 1.1.21