/[lmdze]/trunk/dyn3d/writedynav.f
ViewVC logotype

Contents of /trunk/dyn3d/writedynav.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 261 - (show annotations)
Wed Mar 7 13:33:15 2018 UTC (6 years, 2 months ago) by guez
File size: 3051 byte(s)
Make procedures writehist and writedynav as identical as
possible. Remove output of phis in writedynav: already in histins and
constant. Add output of temp and q in writehist.

Remove unused variable ndm of module dimens_m. Remove unused variables
kftd, mvar, jcfil and jcfllm of module paramet_m.

1 module writedynav_m
2
3 implicit none
4
5 contains
6
7 subroutine writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, time)
8
9 ! From LMDZ4/libf/bibio/writedynav.F, version 1.1.1.1 2004/05/19 12:53:05
10 ! Écriture du fichier histoire au format IOIPSL
11 ! L. Fairhead, LMD, 03/99
12
13 USE comconst, ONLY: cpp
14 use covnat_m, only: covnat
15 use dimens_m, only: llm, nqmx
16 use histsync_m, only: histsync
17 use histwrite_m, only: histwrite
18 use iniadvtrac_m, only: ttext
19 use initdynav_m, only: histaveid
20 use nr_util, only: assert
21 use paramet_m, only: ip1jm, ip1jmp1
22 use temps, only: itau_dyn
23
24 ! Vent covariant :
25 REAL, intent(in):: vcov(:, :, :) ! (iim + 1, jjm, llm)
26 REAL, intent(in):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm)
27
28 REAL, intent(in):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
29 ! temperature potentielle
30
31 real, intent(in):: pk(:, :, :) ! (iim + 1, jjm + 1, llm)
32
33 real, intent(in):: phi(:, :, :) ! (iim + 1, jjm + 1, llm)
34 ! geopotentiel instantann\'e
35
36 REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) traceurs
37 real, intent(in):: masse(:, :, :) ! (iim + 1, jjm + 1, llm)
38 REAL, intent(in):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol
39 integer, intent(in):: time ! temps de l'ecriture
40
41 ! Local:
42 integer iq, itau_w
43 real us(ip1jmp1*llm), vs(ip1jmp1*llm)
44 REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
45
46 !---------------------------------------------------------------------
47
48 call assert([size(vcov, 1), size(ucov, 1), size(teta, 1), size(phi, 1), &
49 size(pk, 1), size(ps, 1), size(masse, 1)] == size(q, 1), &
50 "writedynav iim")
51 call assert([size(vcov, 2) + 1, size(ucov, 2), size(teta, 2), &
52 size(phi, 2), size(pk, 2), size(ps, 2), size(masse, 2)] &
53 == size(q, 2), "writedynav jjm")
54 call assert([size(vcov, 3), size(ucov, 3), size(teta, 3), size(phi, 3), &
55 size(pk, 3), size(masse, 3), size(q, 3)] == llm, "writedynav llm")
56 call assert(size(q, 4) == nqmx, "writedynav nqmx")
57
58 ! Initialisations
59 us = 999.999
60 vs = 999.999
61
62 itau_w = itau_dyn + time
63
64 ! Passage aux composantes naturelles du vent
65 call covnat(llm, ucov, vcov, unat, vnat)
66
67 ! Appels a histwrite pour l'ecriture des variables a sauvegarder
68
69 ! Vents U scalaire
70 call gr_u_scal(llm, unat, us)
71 call histwrite(histaveid, 'u', itau_w, us)
72
73 ! Vents V scalaire
74 call gr_v_scal(llm, vnat, vs)
75 call histwrite(histaveid, 'v', itau_w, vs)
76
77 ! Temperature potentielle moyennee
78 call histwrite(histaveid, 'theta', itau_w, teta)
79
80 ! Temperature moyennee
81 call histwrite(histaveid, 'temp', itau_w, teta * pk / cpp)
82
83 call histwrite(histaveid, 'phi', itau_w, phi)
84
85 ! Traceurs
86 DO iq = 1, size(q, 4)
87 call histwrite(histaveid, ttext(iq), itau_w, q(:, :, :, iq))
88 enddo
89
90 call histwrite(histaveid, 'masse', itau_w, masse)
91 call histwrite(histaveid, 'ps', itau_w, ps)
92
93 call histsync(histaveid)
94
95 end subroutine writedynav
96
97 end module writedynav_m

  ViewVC Help
Powered by ViewVC 1.1.21