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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 62 - (hide annotations)
Thu Jul 26 14:37:37 2012 UTC (11 years, 10 months ago) by guez
Original Path: trunk/libf/bibio/writedynav.f90
File size: 2827 byte(s)
Changed handling of compiler in compilation system.

Removed the prefix letters "y", "p", "t" or "z" in some names of variables.

Replaced calls to NetCDF by calls to NetCDF95.

Extracted "ioget_calendar" procedures from "calendar.f90" into a
separate file.

Extracted to a separate file, "mathop2.f90", procedures that were not
part of the generic interface "mathop" in "mathop.f90".

Removed computation of "dq" in "bilan_dyn", which was not used.

In "iniadvtrac", removed schemes 20 Slopes and 30 Prather. Was not
compatible with declarations of array sizes.

In "clcdrag", "ustarhb", "vdif_kcay", "yamada4" and "coefkz", changed
the size of some arrays from "klon" to "knon".

Removed possible call to "conema3" in "physiq".

Removed unused argument "cd" in "yamada".

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

  ViewVC Help
Powered by ViewVC 1.1.21