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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 134 - (show annotations)
Wed Apr 29 15:47:56 2015 UTC (9 years ago) by guez
File size: 3356 byte(s)
Sources inside, compilation outside.
1 module writedynav_m
2
3 implicit none
4
5 contains
6
7 subroutine writedynav(vcov, ucov, teta, pk, phi, q, masse, ps, phis, 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 ! Appels successifs des routines histwrite
14
15 USE comconst, ONLY: cpp
16 use covnat_m, only: covnat
17 USE dimens_m, ONLY: iim, jjm, llm, nqmx
18 USE histsync_m, ONLY: histsync
19 USE histwrite_m, ONLY: histwrite
20 USE iniadvtrac_m, ONLY: ttext
21 use initdynav_m, only: histaveid
22 use nr_util, only: assert
23 USE paramet_m, ONLY: iip1, ip1jm, ip1jmp1, jjp1
24 USE temps, ONLY: itau_dyn
25
26 ! Vents covariants :
27 REAL, intent(in):: vcov(:, :, :) ! (iim + 1, jjm, llm)
28 REAL, intent(in):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm)
29
30 REAL, intent(in):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
31 ! temperature potentielle
32
33 real, intent(in):: pk(:, :, :) ! (iim + 1, jjm + 1, llm)
34
35 real, intent(in):: phi(:, :, :) ! (iim + 1, jjm + 1, llm)
36 ! geopotentiel instantane
37
38 REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) traceurs
39 real, intent(in):: masse(:, :, :) ! (iim + 1, jjm + 1, llm)
40 REAL, intent(in):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol
41 REAL, intent(in):: phis(:, :) ! (iim + 1, jjm + 1) geopotentiel au sol
42 integer, intent(in):: time ! temps de l'ecriture
43
44 ! Variables locales
45 integer iq
46 real us(ip1jmp1*llm), vs(ip1jmp1*llm)
47 REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
48 integer itau_w
49
50 !---------------------------------------------------------------
51
52 call assert((/size(vcov, 1), size(ucov, 1), size(teta, 1), size(phi, 1), &
53 size(pk, 1), size(ps, 1), size(masse, 1), size(phis, 1), &
54 size(q, 1)/) == iim + 1, "writedynav iim")
55 call assert((/size(vcov, 2) + 1, size(ucov, 2), size(teta, 2), &
56 size(phi, 2), size(pk, 2), size(ps, 2), size(masse, 2), &
57 size(phis, 2), size(q, 2)/) == jjm + 1, "writedynav jjm")
58 call assert((/size(vcov, 3), size(ucov, 3), size(teta, 3), size(phi, 3), &
59 size(pk, 3), size(masse, 3), size(q, 3)/) == llm, "writedynav llm")
60 call assert(size(q, 4) == nqmx, "writedynav nqmx")
61
62 ! Initialisations
63 us = 999.999
64 vs = 999.999
65 vnat = 999.999
66 unat = 999.999
67 itau_w = itau_dyn + time
68
69 ! Passage aux composantes naturelles du vent
70 call covnat(llm, ucov, vcov, unat, vnat)
71
72 ! Appels a histwrite pour l'ecriture des variables a sauvegarder
73
74 ! Vents U scalaire
75 call gr_u_scal(llm, unat, us)
76 call histwrite(histaveid, 'u', itau_w, us)
77
78 ! Vents V scalaire
79 call gr_v_scal(llm, vnat, vs)
80 call histwrite(histaveid, 'v', itau_w, vs)
81
82 ! Temperature potentielle moyennee
83 call histwrite(histaveid, 'theta', itau_w, teta)
84
85 ! Temperature moyennee
86 call histwrite(histaveid, 'temp', itau_w, teta * pk / cpp)
87
88 call histwrite(histaveid, 'phi', itau_w, phi)
89
90 ! Traceurs
91 DO iq = 1, size(q, 4)
92 call histwrite(histaveid, ttext(iq), itau_w, q(:, :, :, iq))
93 enddo
94
95 call histwrite(histaveid, 'masse', itau_w, masse)
96 call histwrite(histaveid, 'ps', itau_w, ps)
97 call histwrite(histaveid, 'phis', itau_w, phis)
98
99 call histsync(histaveid)
100
101 end subroutine writedynav
102
103 end module writedynav_m

  ViewVC Help
Powered by ViewVC 1.1.21