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

Contents of /trunk/libf/bibio/writedynav.f

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21