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

Annotation of /trunk/Sources/bibio/writedynav.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
Original Path: trunk/libf/bibio/writedynav.f90
File size: 2848 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 module writedynav_m
2 guez 3
3 guez 56 implicit none
4 guez 3
5 guez 56 contains
6 guez 3
7 guez 61 subroutine writedynav(vcov, ucov, teta, ppk, phi, q, masse, ps, phis, time)
8 guez 3
9 guez 56 ! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05
10     ! Ecriture du fichier histoire au format IOIPSL
11 guez 3
12 guez 56 ! Appels successifs des routines histwrite
13 guez 3
14 guez 56 ! Entree:
15     ! vcov: vents v covariants
16     ! ucov: vents u covariants
17     ! phi : geopotentiel instantane
18     ! q : traceurs
19     ! ps :pression au sol
20     ! phis : geopotentiel au sol
21 guez 3
22 guez 56 ! L. Fairhead, LMD, 03/99
23 guez 3
24 guez 61 use covnat_m, only: covnat
25 guez 56 USE histwrite_m, ONLY: histwrite
26 guez 61 USE histsync_m, ONLY: histsync
27 guez 56 USE dimens_m, ONLY: llm
28     USE paramet_m, ONLY: iip1, ijp1llm, ip1jm, ip1jmp1, jjp1
29     USE comconst, ONLY: cpp
30     USE temps, ONLY: itau_dyn
31     USE iniadvtrac_m, ONLY: ttext
32     use initdynav_m, only: histaveid
33 guez 45
34 guez 56 REAL vcov(ip1jm, llm), ucov(ip1jmp1, llm)
35     REAL, intent(in):: teta(ip1jmp1*llm) ! temperature potentielle
36     real phi(ip1jmp1, llm), ppk(ip1jmp1*llm)
37     REAL ps(ip1jmp1)
38     real, intent(in):: masse(ip1jmp1, llm)
39     REAL phis(ip1jmp1)
40 guez 61 REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
41     integer, intent(in):: time ! temps de l'ecriture
42 guez 45
43 guez 56 ! Variables locales
44     integer ndex2d(iip1*jjp1), ndex3d(iip1*jjp1*llm), iq, ii, ll
45     real us(ip1jmp1*llm), vs(ip1jmp1*llm)
46     real tm(ip1jmp1*llm)
47     REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
48     logical ok_sync
49     integer itau_w
50 guez 45
51 guez 56 !---------------------------------------------------------------
52 guez 45
53 guez 56 ! Initialisations
54     ndex3d = 0
55     ndex2d = 0
56     ok_sync = .TRUE.
57     us = 999.999
58     vs = 999.999
59     tm = 999.999
60     vnat = 999.999
61     unat = 999.999
62     itau_w = itau_dyn + time
63 guez 45
64 guez 56 ! Passage aux composantes naturelles du vent
65     call covnat(llm, ucov, vcov, unat, vnat)
66 guez 45
67 guez 56 ! Appels a histwrite pour l'ecriture des variables a sauvegarder
68 guez 45
69 guez 56 ! Vents U scalaire
70     call gr_u_scal(llm, unat, us)
71     call histwrite(histaveid, 'u', itau_w, us)
72 guez 45
73 guez 56 ! Vents V scalaire
74     call gr_v_scal(llm, vnat, vs)
75     call histwrite(histaveid, 'v', itau_w, vs)
76 guez 45
77 guez 56 ! Temperature potentielle moyennee
78     call histwrite(histaveid, 'theta', itau_w, teta)
79 guez 45
80 guez 56 ! Temperature moyennee
81     do ii = 1, ijp1llm
82     tm(ii) = teta(ii) * ppk(ii)/cpp
83     enddo
84     call histwrite(histaveid, 'temp', itau_w, tm)
85 guez 45
86 guez 56 ! Geopotentiel
87     call histwrite(histaveid, 'phi', itau_w, phi)
88 guez 45
89 guez 56 ! Traceurs
90 guez 61 DO iq = 1, size(q, 4)
91     call histwrite(histaveid, ttext(iq), itau_w, q(:, :, :, iq))
92 guez 56 enddo
93 guez 45
94 guez 56 ! Masse
95     call histwrite(histaveid, 'masse', itau_w, masse)
96 guez 45
97 guez 56 ! Pression au sol
98     call histwrite(histaveid, 'ps', itau_w, ps)
99 guez 45
100 guez 56 ! Geopotentiel au sol
101     call histwrite(histaveid, 'phis', itau_w, phis)
102    
103     if (ok_sync) call histsync(histaveid)
104    
105     end subroutine writedynav
106    
107     end module writedynav_m

  ViewVC Help
Powered by ViewVC 1.1.21