/[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 61 - (hide annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years, 1 month ago) by guez
File size: 2375 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

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 guez 61 USE histbeg_totreg_m, ONLY : histbeg_totreg
20     USE histdef_m, ONLY : histdef
21     USE histend_m, ONLY : histend
22     USE histvert_m, ONLY : histvert
23 guez 34 use comvert, only: presnivs
24    
25     REAL, intent(in):: dtime ! pas temporel de la physique (s)
26     integer, intent(out):: nid_hf3d
27    
28     real zstohf, zout
29     REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1)
30     real zjulian
31     integer i, nhori, nvert, idayref
32    
33     !------------------------------------------
34    
35     zstohf = dtime * REAL(ecrit_hf)
36     zout = dtime * REAL(ecrit_hf)
37    
38     idayref = day_ref
39     CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
40    
41     CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon)
42     DO i = 1, iim
43     zx_lon(i, 1) = rlon(i+1)
44     zx_lon(i, (jjm + 1)) = rlon(i+1)
45     ENDDO
46    
47     CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat)
48     CALL histbeg_totreg("histhf3d", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, &
49     (jjm + 1), itau_phy, zjulian, dtime, nhori, nid_hf3d)
50    
51     CALL histvert(nid_hf3d, "presnivs", "Vertical levels", "mb", &
52     llm, presnivs/100., nvert)
53    
54     ! Champs 3D:
55    
56     CALL histdef(nid_hf3d, "temp", "Air temperature", "K", &
57     iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
58     "ave(X)", zstohf, zout)
59    
60     CALL histdef(nid_hf3d, "ovap", "Specific humidity", "kg/kg", &
61     iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
62     "ave(X)", zstohf, zout)
63    
64     CALL histdef(nid_hf3d, "vitu", "Zonal wind", "m/s", &
65     iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
66     "ave(X)", zstohf, zout)
67    
68     CALL histdef(nid_hf3d, "vitv", "Meridional wind", "m/s", &
69     iim, (jjm + 1), nhori, llm, 1, llm, nvert, &
70     "ave(X)", zstohf, zout)
71    
72     if (nbtr >= 3) then
73     CALL histdef(nid_hf3d, "O3", "Ozone mass fraction", "?", iim, &
74     (jjm + 1), nhori, llm, 1, llm, nvert, "ave(X)", zstohf, &
75     zout)
76     end if
77    
78     CALL histend(nid_hf3d)
79    
80     end subroutine ini_histhf3d
81    
82     end module ini_histhf3d_m

  ViewVC Help
Powered by ViewVC 1.1.21