/[lmdze]/trunk/Sources/dyn3d/writehist.f
ViewVC logotype

Annotation of /trunk/Sources/dyn3d/writehist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 69 - (hide annotations)
Mon Feb 18 16:33:12 2013 UTC (11 years, 3 months ago) by guez
Original Path: trunk/libf/bibio/writehist.f90
File size: 2377 byte(s)
Deleted files cvparam3.f90 and nuagecom.f90. Moved variables from
module cvparam3 to module cv3_param_m. Moved variables rad_chau1 and
rad_chau2 from module nuagecom to module conf_phys_m.

Read clesphys2_nml from conf_phys instead of gcm.

Removed argument iflag_con from several procedures. Access module
variable instead.

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

  ViewVC Help
Powered by ViewVC 1.1.21