/[lmdze]/trunk/Sources/bibio/writedynav.f
ViewVC logotype

Annotation of /trunk/Sources/bibio/writedynav.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 56 - (hide annotations)
Tue Jan 10 19:02:02 2012 UTC (12 years, 4 months ago) by guez
Original Path: trunk/libf/bibio/writedynav.f90
File size: 2799 byte(s)
Imported "writehist.f" from LMDZ.

Moved module variable "histaveid" from "com_io_dyn" to "initdynav_m".

In "inithist", access directly module variables from "com_io_dyn"
instead of going through the arguments. Copying from LMDZ, write "u"
and scalar variables to separate files. Create a new variable for the
new file in "com_io_dyn". Copying from LMDZ, change the vertical axes
of the three files.

Removed some useless initializations in "dissip".

In "bilan_dyn", removed useless variable "time". Avoiding the
approximate test on "dt_cum" being a multiple of "dt_app", just
compute "ncum" from known usage of "bilan_dyn" and compute "dt_cum"
from "ncum". Change "periodav" from real to integer in
"conf_gcm_m". Since "day_step" is required to be a multiple of
"iperiod", so is "ncum".

1 guez 56 module writedynav_m
2 guez 3
3 guez 56 implicit none
4 guez 3
5 guez 56 contains
6 guez 3
7 guez 56 subroutine writedynav(nq, time, vcov, ucov, teta, ppk, phi, q, masse, ps, &
8     phis)
9 guez 3
10 guez 56 ! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05
11     ! Ecriture du fichier histoire au format IOIPSL
12 guez 3
13 guez 56 ! Appels successifs des routines histwrite
14 guez 3
15 guez 56 ! Entree:
16     ! time: temps de l'ecriture
17     ! vcov: vents v covariants
18     ! ucov: vents u covariants
19     ! phi : geopotentiel instantane
20     ! q : traceurs
21     ! ps :pression au sol
22     ! phis : geopotentiel au sol
23 guez 3
24 guez 56 ! L. Fairhead, LMD, 03/99
25 guez 3
26 guez 56 USE histwrite_m, ONLY: histwrite
27     USE histcom, ONLY: histsync
28     USE dimens_m, ONLY: llm
29     USE paramet_m, ONLY: iip1, ijp1llm, ip1jm, ip1jmp1, jjp1
30     USE comconst, ONLY: cpp
31     USE temps, ONLY: itau_dyn
32     USE iniadvtrac_m, ONLY: ttext
33     use initdynav_m, only: histaveid
34 guez 45
35 guez 56 INTEGER nq
36     REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm)
37     REAL, intent(in):: teta(ip1jmp1*llm) ! temperature potentielle
38     real phi(ip1jmp1, llm), ppk(ip1jmp1*llm)
39     REAL ps(ip1jmp1)
40     real, intent(in):: masse(ip1jmp1, llm)
41     REAL phis(ip1jmp1)
42     REAL q(ip1jmp1, llm, nq)
43     integer, intent(in):: time
44 guez 45
45 guez 56 ! Variables locales
46     integer ndex2d(iip1*jjp1), ndex3d(iip1*jjp1*llm), iq, ii, ll
47     real us(ip1jmp1*llm), vs(ip1jmp1*llm)
48     real tm(ip1jmp1*llm)
49     REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
50     logical ok_sync
51     integer itau_w
52 guez 45
53 guez 56 !---------------------------------------------------------------
54 guez 45
55 guez 56 ! Initialisations
56     ndex3d = 0
57     ndex2d = 0
58     ok_sync = .TRUE.
59     us = 999.999
60     vs = 999.999
61     tm = 999.999
62     vnat = 999.999
63     unat = 999.999
64     itau_w = itau_dyn + time
65 guez 45
66 guez 56 ! Passage aux composantes naturelles du vent
67     call covnat(llm, ucov, vcov, unat, vnat)
68 guez 45
69 guez 56 ! Appels a histwrite pour l'ecriture des variables a sauvegarder
70 guez 45
71 guez 56 ! Vents U scalaire
72     call gr_u_scal(llm, unat, us)
73     call histwrite(histaveid, 'u', itau_w, us)
74 guez 45
75 guez 56 ! Vents V scalaire
76     call gr_v_scal(llm, vnat, vs)
77     call histwrite(histaveid, 'v', itau_w, vs)
78 guez 45
79 guez 56 ! Temperature potentielle moyennee
80     call histwrite(histaveid, 'theta', itau_w, teta)
81 guez 45
82 guez 56 ! Temperature moyennee
83     do ii = 1, ijp1llm
84     tm(ii) = teta(ii) * ppk(ii)/cpp
85     enddo
86     call histwrite(histaveid, 'temp', itau_w, tm)
87 guez 45
88 guez 56 ! Geopotentiel
89     call histwrite(histaveid, 'phi', itau_w, phi)
90 guez 45
91 guez 56 ! Traceurs
92     DO iq=1, nq
93     call histwrite(histaveid, ttext(iq), itau_w, q(:, :, iq))
94     enddo
95 guez 45
96 guez 56 ! Masse
97     call histwrite(histaveid, 'masse', itau_w, masse)
98 guez 45
99 guez 56 ! Pression au sol
100     call histwrite(histaveid, 'ps', itau_w, ps)
101 guez 45
102 guez 56 ! Geopotentiel au sol
103     call histwrite(histaveid, 'phis', itau_w, phis)
104    
105     if (ok_sync) call histsync(histaveid)
106    
107     end subroutine writedynav
108    
109     end module writedynav_m

  ViewVC Help
Powered by ViewVC 1.1.21