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

Contents of /trunk/libf/phylmd/ini_histday.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (show annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years ago) by guez
File size: 1771 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 module ini_histday_m
2
3 implicit none
4
5 contains
6
7 subroutine ini_histday(dtime, ok_journe, nid_day, nq)
8
9 ! From phylmd/ini_histday.h, v 1.3 2005/05/25 13:10:09
10
11 use dimens_m, only: iim, jjm, llm
12 use temps, only: itau_phy, day_ref, annee_ref
13 USE calendar, only: ymds2ju
14 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 use phyetat0_m, only: rlon, rlat
19 use clesphys, only: ecrit_day
20 use grid_change, only: gr_phy_write_2d
21 use comvert, only: presnivs
22
23 REAL, intent(in):: dtime ! pas temporel de la physique (s)
24 logical, intent(in):: ok_journe
25 integer, intent(out):: nid_day
26 INTEGER, intent(in):: nq ! nombre de traceurs (y compris vapeur d'eau)
27
28 ! Variables local to the procedure:
29 REAL zx_lat(iim, jjm + 1)
30 integer nhori, nvert
31 real zjulian
32
33 !--------------------------------
34
35 IF (ok_journe) THEN
36 CALL ymds2ju(annee_ref, 1, day_ref, 0., zjulian)
37 zx_lat = gr_phy_write_2d(rlat)
38 CALL histbeg_totreg("histday", rlon(2: iim+1), zx_lat(1, :), 1, iim, &
39 1, jjm + 1, itau_phy, zjulian, dtime, nhori, nid_day)
40 CALL histvert(nid_day, "presnivs", "Vertical levels", "mb", &
41 llm, presnivs/100., nvert)
42 if (nq <= 4) then
43 call histdef(nid_day, "Sigma_O3_Royer", &
44 "column-density of ozone, in a cell, from Royer", "DU", &
45 xsize=iim, ysize=jjm+1, horiid=nhori, pzsize=llm, par_oriz=1, &
46 par_szz=llm, pzid=nvert, popp="ave(X)", pfreq_opp=dtime, &
47 pfreq_wrt=real(ecrit_day))
48 end if
49 CALL histend(nid_day)
50 ENDIF
51
52 end subroutine ini_histday
53
54 end module ini_histday_m

  ViewVC Help
Powered by ViewVC 1.1.21