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

Annotation of /trunk/libf/bibio/writehist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 61 - (hide annotations)
Fri Apr 20 14:58:43 2012 UTC (12 years, 1 month ago) by guez
File size: 2648 byte(s)
No more included file in LMDZE, not even "netcdf.inc".

Created a variable containing the list of common source files in
GNUmakefile. So we now also see clearly files that are specific to
each program.

Split module "histcom". Assembled resulting files in directory
"Histcom".

Removed aliasing in calls to "laplacien".

1 guez 56 !
2     ! $Id: writehist.F 1403 2010-07-01 09:02:53Z fairhead $
3     !
4     subroutine writehist(time,vcov,ucov,teta,phi,q,masse,ps,phis)
5    
6     use dimens_m, only: nqmx, llm, jjm
7     USE iniadvtrac_m, ONLY: ttext
8     use com_io_dyn, only: histid,histvid,histuid
9     use paramet_m, only: ip1jm, ip1jmp1, iip1, jjp1
10     use temps, only: itau_dyn
11     use histwrite_m, only: histwrite
12 guez 61 use histsync_m, only: histsync
13     use covnat_m, only: covnat
14 guez 56
15     implicit none
16    
17     C
18     C Ecriture du fichier histoire au format IOIPSL
19     C
20     C Appels succesifs des routines: histwrite
21     C
22     C Entree:
23     C time: temps de l'ecriture
24     C vcov: vents v covariants
25     C ucov: vents u covariants
26     C teta: temperature potentielle
27     C phi : geopotentiel instantane
28     C q : traceurs
29     C masse: masse
30     C ps :pression au sol
31     C phis : geopotentiel au sol
32     C
33     C
34     C L. Fairhead, LMD, 03/99
35     C
36     C =====================================================================
37     C
38     C Declarations
39    
40     C
41     C Arguments
42     C
43    
44     REAL vcov(ip1jm,llm),ucov(ip1jmp1,llm)
45     REAL teta(ip1jmp1,llm),phi(ip1jmp1,llm)
46     REAL ps(ip1jmp1),masse(ip1jmp1,llm)
47     REAL phis(ip1jmp1)
48     REAL q(ip1jmp1,llm,nqmx)
49     integer time
50    
51    
52     ! This routine needs IOIPSL to work
53     C Variables locales
54     C
55     integer iq, ii, ll
56     integer ndexu(ip1jmp1*llm),ndexv(ip1jm*llm),ndex2d(ip1jmp1)
57     logical ok_sync
58     integer itau_w
59     REAL vnat(ip1jm,llm),unat(ip1jmp1,llm)
60    
61     C
62     C Initialisations
63     C
64     ndexu = 0
65     ndexv = 0
66     ndex2d = 0
67     ok_sync =.TRUE.
68     itau_w = itau_dyn + time
69     ! Passage aux composantes naturelles du vent
70     call covnat(llm, ucov, vcov, unat, vnat)
71     C
72     C Appels a histwrite pour l'ecriture des variables a sauvegarder
73     C
74     C Vents U
75     C
76     call histwrite(histuid, 'u', itau_w, unat)
77     C
78     C Vents V
79     C
80     call histwrite(histvid, 'v', itau_w, vnat)
81    
82     C
83     C Temperature potentielle
84     C
85     call histwrite(histid, 'teta', itau_w, teta)
86     C
87     C Geopotentiel
88     C
89     call histwrite(histid, 'phi', itau_w, phi)
90     C
91     C Traceurs
92     C
93     ! DO iq=1,nqmx
94     ! call histwrite(histid, ttext(iq), itau_w, q(:,:,iq),
95     ! . iip1*jjp1*llm, ndexu)
96     ! enddo
97     !C
98     C Masse
99     C
100     call histwrite(histid,'masse',itau_w, masse)
101     C
102     C Pression au sol
103     C
104     call histwrite(histid, 'ps', itau_w, ps)
105     C
106     C Geopotentiel au sol
107     C
108     ! call histwrite(histid, 'phis', itau_w, phis, iip1*jjp1, ndex2d)
109     C
110     C Fin
111     C
112     if (ok_sync) then
113     call histsync(histid)
114     call histsync(histvid)
115     call histsync(histuid)
116     endif
117     return
118     end

  ViewVC Help
Powered by ViewVC 1.1.21