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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 34 - (hide annotations)
Wed Jun 2 11:01:12 2010 UTC (13 years, 11 months ago) by guez
File size: 1415 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_histhf_m
2 guez 3
3 guez 34 implicit none
4 guez 3
5     contains
6    
7 guez 20 subroutine ini_histhf(dtime, nid_hf, nid_hf3d)
8 guez 3
9     ! From phylmd/ini_histhf.h, version 1.3 2005/05/25 13:10:09
10    
11     use dimens_m, only: iim, jjm, llm
12     use temps, only: day_ref, annee_ref, itau_phy
13     use dimphy, only: klon
14 guez 30 USE calendar, only: ymds2ju
15     use histcom, only: histbeg_totreg, histvert, histend
16 guez 3 use phyetat0_m, only: rlon, rlat
17 guez 20 use comvert, only: presnivs
18 guez 34 use ini_histhf3d_m, only: ini_histhf3d
19 guez 3
20     REAL, intent(in):: dtime ! pas temporel de la physique (s)
21     integer, intent(out):: nid_hf, nid_hf3d
22    
23     REAL zx_lon(iim, jjm + 1), zx_lat(iim, jjm + 1)
24     integer idayref
25     real zjulian
26     integer i, nhori, nvert
27    
28     !-----------------------------------------------
29    
30     idayref = day_ref
31     CALL ymds2ju(annee_ref, 1, idayref, 0.0, zjulian)
32    
33     CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlon, zx_lon)
34     DO i = 1, iim
35     zx_lon(i, 1) = rlon(i+1)
36     zx_lon(i, (jjm + 1)) = rlon(i+1)
37     ENDDO
38    
39     CALL gr_fi_ecrit(1, klon, iim, (jjm + 1), rlat, zx_lat)
40 guez 15 CALL histbeg_totreg("histhf", zx_lon(:, 1), zx_lat(1, :), 1, iim, 1, &
41     (jjm + 1), itau_phy, zjulian, dtime, nhori, nid_hf)
42 guez 3
43     CALL histvert(nid_hf, "presnivs", "Vertical levels", "mb", &
44     llm, presnivs/100., nvert)
45    
46 guez 20 call ini_histhf3d(dtime, nid_hf3d)
47 guez 3 CALL histend(nid_hf)
48    
49     end subroutine ini_histhf
50    
51 guez 34 end module ini_histhf_m

  ViewVC Help
Powered by ViewVC 1.1.21