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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 30 - (hide annotations)
Thu Apr 1 09:07:28 2010 UTC (14 years, 2 months ago) by guez
Original Path: trunk/libf/bibio/writedynav.f
File size: 2973 byte(s)
Imported Source files of the external library "IOIPSL_Lionel" into
"libf/IOIPSL".

Split "cray.f90" into "scopy.f90" and "ssum.f90".

Rewrote "leapfrog" in order to have a clearer algorithmic structure.

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     REAL teta(ip1jmp1*llm),phi(ip1jmp1,llm),ppk(ip1jmp1*llm)
55     REAL ps(ip1jmp1),masse(ip1jmp1,llm)
56     REAL phis(ip1jmp1)
57     REAL q(ip1jmp1,llm,nq)
58     integer, intent(in):: time
59    
60    
61     C Variables locales
62     C
63     integer ndex2d(iip1*jjp1),ndex3d(iip1*jjp1*llm),iq, ii, ll
64     real us(ip1jmp1*llm), vs(ip1jmp1*llm)
65     real tm(ip1jmp1*llm)
66     REAL vnat(ip1jm,llm),unat(ip1jmp1,llm)
67     logical ok_sync
68     integer itau_w
69     C
70     C Initialisations
71     C
72     ndex3d = 0
73     ndex2d = 0
74     ok_sync = .TRUE.
75     us = 999.999
76     vs = 999.999
77     tm = 999.999
78     vnat = 999.999
79     unat = 999.999
80     itau_w = itau_dyn + time
81    
82     C Passage aux composantes naturelles du vent
83     call covnat(llm, ucov, vcov, unat, vnat)
84    
85     C
86     C Appels a histwrite pour l'ecriture des variables a sauvegarder
87     C
88     C Vents U scalaire
89     C
90     call gr_u_scal(llm, unat, us)
91 guez 15 call histwrite(histid, 'u', itau_w, us)
92 guez 3 C
93     C Vents V scalaire
94     C
95     call gr_v_scal(llm, vnat, vs)
96 guez 15 call histwrite(histid, 'v', itau_w, vs)
97 guez 3 C
98     C Temperature potentielle moyennee
99     C
100 guez 15 call histwrite(histid, 'theta', itau_w, teta)
101 guez 3 C
102     C Temperature moyennee
103     C
104     do ii = 1, ijp1llm
105     tm(ii) = teta(ii) * ppk(ii)/cpp
106     enddo
107 guez 15 call histwrite(histid, 'temp', itau_w, tm)
108 guez 3 C
109     C Geopotentiel
110     C
111 guez 15 call histwrite(histid, 'phi', itau_w, phi)
112 guez 3 C
113     C Traceurs
114     C
115     DO iq=1,nq
116 guez 15 call histwrite(histid, ttext(iq), itau_w, q(:,:,iq))
117 guez 3 enddo
118     C
119     C Masse
120     C
121 guez 15 call histwrite(histid, 'masse', itau_w, masse)
122 guez 3 C
123     C Pression au sol
124     C
125 guez 15 call histwrite(histid, 'ps', itau_w, ps)
126 guez 3 C
127     C Geopotentiel au sol
128     C
129 guez 15 call histwrite(histid, 'phis', itau_w, phis)
130 guez 3 C
131     C Fin
132     C
133     if (ok_sync) call histsync(histid)
134     return
135     end

  ViewVC Help
Powered by ViewVC 1.1.21