/[lmdze]/trunk/libf/bibio/writedynav.f90
ViewVC logotype

Annotation of /trunk/libf/bibio/writedynav.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 44 - (hide annotations)
Wed Apr 13 12:29:18 2011 UTC (13 years, 1 month ago) by guez
Original Path: trunk/libf/bibio/writedynav.f
File size: 2998 byte(s)
Removed argument "pdteta" of "calfis", because it was not used.

Created module "conf_guide_m", containing procedure
"conf_guide". Moved module variables from "guide_m" to "conf_guide_m".

In module "getparam", removed "ini_getparam" and "fin_getparam" from
generic interface "getpar".

Created module variables in "tau2alpha_m" to replace common "comdxdy".

1 guez 3 !
2     ! $Header: /home/cvsroot/LMDZ4/libf/bibio/writedynav.F,v 1.1.1.1 2004/05/19 12:53:05 lmdzadmin Exp $
3     !
4     subroutine writedynav( histid, nq, time, vcov,
5     , ucov,teta,ppk,phi,q,masse,ps,phis)
6    
7     C
8     C Ecriture du fichier histoire au format IOIPSL
9     C
10     C Appels succesifs des routines: histwrite
11     C
12     C Entree:
13     C histid: ID du fichier histoire
14     C time: temps de l'ecriture
15     C vcov: vents v covariants
16     C ucov: vents u covariants
17     C teta: temperature potentielle
18     C phi : geopotentiel instantane
19     C q : traceurs
20     C masse: masse
21     C ps :pression au sol
22     C phis : geopotentiel au sol
23     C
24     C
25     C Sortie:
26     C fileid: ID du fichier netcdf cree
27     C
28     C L. Fairhead, LMD, 03/99
29     C
30     C =====================================================================
31     C
32     C Declarations
33 guez 30 USE histwrite_m
34     use histcom
35 guez 3 use dimens_m
36     use paramet_m
37     use comconst
38     use comvert
39     use logic
40     use comgeom
41     use serre
42     use temps
43     use ener
44 guez 18 use iniadvtrac_m
45 guez 3 implicit none
46    
47    
48     C
49     C Arguments
50     C
51    
52     INTEGER histid, nq
53     REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm)
54 guez 44 REAL, intent(in):: teta(ip1jmp1*llm)
55     real phi(ip1jmp1,llm),ppk(ip1jmp1*llm)
56 guez 3 REAL ps(ip1jmp1),masse(ip1jmp1,llm)
57     REAL phis(ip1jmp1)
58     REAL q(ip1jmp1,llm,nq)
59     integer, intent(in):: time
60    
61    
62     C Variables locales
63     C
64     integer ndex2d(iip1*jjp1),ndex3d(iip1*jjp1*llm),iq, ii, ll
65     real us(ip1jmp1*llm), vs(ip1jmp1*llm)
66     real tm(ip1jmp1*llm)
67     REAL vnat(ip1jm,llm),unat(ip1jmp1,llm)
68     logical ok_sync
69     integer itau_w
70     C
71     C Initialisations
72     C
73     ndex3d = 0
74     ndex2d = 0
75     ok_sync = .TRUE.
76     us = 999.999
77     vs = 999.999
78     tm = 999.999
79     vnat = 999.999
80     unat = 999.999
81     itau_w = itau_dyn + time
82    
83     C Passage aux composantes naturelles du vent
84     call covnat(llm, ucov, vcov, unat, vnat)
85    
86     C
87     C Appels a histwrite pour l'ecriture des variables a sauvegarder
88     C
89     C Vents U scalaire
90     C
91     call gr_u_scal(llm, unat, us)
92 guez 15 call histwrite(histid, 'u', itau_w, us)
93 guez 3 C
94     C Vents V scalaire
95     C
96     call gr_v_scal(llm, vnat, vs)
97 guez 15 call histwrite(histid, 'v', itau_w, vs)
98 guez 3 C
99     C Temperature potentielle moyennee
100     C
101 guez 15 call histwrite(histid, 'theta', itau_w, teta)
102 guez 3 C
103     C Temperature moyennee
104     C
105     do ii = 1, ijp1llm
106     tm(ii) = teta(ii) * ppk(ii)/cpp
107     enddo
108 guez 15 call histwrite(histid, 'temp', itau_w, tm)
109 guez 3 C
110     C Geopotentiel
111     C
112 guez 15 call histwrite(histid, 'phi', itau_w, phi)
113 guez 3 C
114     C Traceurs
115     C
116     DO iq=1,nq
117 guez 15 call histwrite(histid, ttext(iq), itau_w, q(:,:,iq))
118 guez 3 enddo
119     C
120     C Masse
121     C
122 guez 15 call histwrite(histid, 'masse', itau_w, masse)
123 guez 3 C
124     C Pression au sol
125     C
126 guez 15 call histwrite(histid, 'ps', itau_w, ps)
127 guez 3 C
128     C Geopotentiel au sol
129     C
130 guez 15 call histwrite(histid, 'phis', itau_w, phis)
131 guez 3 C
132     C Fin
133     C
134     if (ok_sync) call histsync(histid)
135     return
136     end

  ViewVC Help
Powered by ViewVC 1.1.21