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

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

Parent Directory Parent Directory | Revision Log Revision Log


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

  ViewVC Help
Powered by ViewVC 1.1.21