/[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 34 - (show annotations)
Wed Jun 2 11:01:12 2010 UTC (13 years, 11 months ago) by guez
File size: 1688 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 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 histcom, only: histbeg_totreg, histvert, histend, histdef
15 use phyetat0_m, only: rlon, rlat
16 use clesphys, only: ecrit_day
17 use grid_change, only: gr_phy_write_2d
18 use comvert, only: presnivs
19
20 REAL, intent(in):: dtime ! pas temporel de la physique (s)
21 logical, intent(in):: ok_journe
22 integer, intent(out):: nid_day
23 INTEGER, intent(in):: nq ! nombre de traceurs (y compris vapeur d'eau)
24
25 ! Variables local to the procedure:
26 REAL zx_lat(iim, jjm + 1)
27 integer nhori, nvert
28 real zjulian
29
30 !--------------------------------
31
32 IF (ok_journe) THEN
33 CALL ymds2ju(annee_ref, 1, day_ref, 0., zjulian)
34 zx_lat = gr_phy_write_2d(rlat)
35 CALL histbeg_totreg("histday", rlon(2: iim+1), zx_lat(1, :), 1, iim, &
36 1, jjm + 1, itau_phy, zjulian, dtime, nhori, nid_day)
37 CALL histvert(nid_day, "presnivs", "Vertical levels", "mb", &
38 llm, presnivs/100., nvert)
39 if (nq <= 4) then
40 call histdef(nid_day, "Sigma_O3_Royer", &
41 "column-density of ozone, in a cell, from Royer", "DU", &
42 pxsize=iim, pysize=jjm+1, phoriid=nhori, pzsize=llm, &
43 par_oriz=1, par_szz=llm, pzid=nvert, popp="ave(X)", &
44 pfreq_opp=dtime, pfreq_wrt=real(ecrit_day))
45 end if
46 CALL histend(nid_day)
47 ENDIF
48
49 end subroutine ini_histday
50
51 end module ini_histday_m

  ViewVC Help
Powered by ViewVC 1.1.21