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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 115 - (hide annotations)
Fri Sep 19 17:36:20 2014 UTC (9 years, 8 months ago) by guez
Original Path: trunk/bibio/writehist.f
File size: 2368 byte(s)
Extracted code from tau2alpha for first call into new procedure
init_tau2alpha. dxdys, dxdyu, dxdyv are now local variables if guide
computed by init_tau2alpha. This allows us to remove terrible argument
type of tau2alpha: we just give to tau2alpha the right dxdy and
rlat.

In module conf_guide_m, changed default values of tau_min_*, because
0.02 is too small for the default daystep = 240, iperiod = 5. Changed
default values of guide_[uv] to false. Moved variable ok_guide from
conf_gcm_m to conf_guide_m, ok_guide is no longer an input parameter,
it is computed from guide_*. Had then to move test on ok_guide and
day_step from conf_gcm_m to conf_guide_m. Added checks on input
nudging parameters in procedure conf_guide. Upgraded variable factt to
module conf_guide_m in order to check nudging parameters. Bug fix:
variable guide_q was not in namelist conf_guide_nml.

Removed unused variables aire_min, aire_max of MODULE guide_m.

Moved the call to conf_guide from guide to gcm. This was needed to
define ok_guide before getting into guide (since ok_guide is no longer
in conf_gcm_m). Moved test on grossismx and grossismy from tau2alpha
to guide. It is clearer now that only tau_max is used for a regular
grid, and we do not have to repeat this test in each call to
tau2alpha. In guide, we only call writefield when alpha is not a
constant.

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 guez 83 REAL, intent(in):: ps(ip1jmp1), masse(ip1jmp1, llm)
38 guez 69 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 ndexu(ip1jmp1*llm), ndexv(ip1jm*llm), ndex2d(ip1jmp1)
45     logical ok_sync
46     integer itau_w
47     REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
48    
49     !---------------------------------------------------------------------
50    
51     ! Initialisations
52    
53     ndexu = 0
54     ndexv = 0
55     ndex2d = 0
56     ok_sync =.TRUE.
57     itau_w = itau_dyn + time
58     ! Passage aux composantes naturelles du vent
59     call covnat(llm, ucov, vcov, unat, vnat)
60    
61     ! Appels a histwrite pour l'ecriture des variables a sauvegarder
62    
63     ! Vents U
64    
65     call histwrite(histuid, 'u', itau_w, unat)
66    
67     ! Vents V
68    
69     call histwrite(histvid, 'v', itau_w, vnat)
70    
71     ! Temperature potentielle
72    
73     call histwrite(histid, 'teta', itau_w, teta)
74    
75     ! Geopotentiel
76    
77     call histwrite(histid, 'phi', itau_w, phi)
78    
79     ! Traceurs
80    
81     ! DO iq=1, nqmx
82     ! call histwrite(histid, ttext(iq), itau_w, q(:, :, iq),
83     ! . iip1*jjp1*llm, ndexu)
84     ! enddo
85     !C
86     ! Masse
87    
88     call histwrite(histid, 'masse', itau_w, masse)
89    
90     ! Pression au sol
91     call histwrite(histid, 'ps', itau_w, ps)
92    
93     if (ok_sync) then
94     call histsync(histid)
95     call histsync(histvid)
96     call histsync(histuid)
97     endif
98    
99     end subroutine writehist
100    
101     end module writehist_m

  ViewVC Help
Powered by ViewVC 1.1.21