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

Annotation of /trunk/dyn3d/writehist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 261 - (hide annotations)
Wed Mar 7 13:33:15 2018 UTC (6 years, 2 months ago) by guez
File size: 2571 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 guez 69 module writehist_m
2 guez 56
3 guez 69 implicit none
4 guez 56
5 guez 69 contains
6 guez 56
7 guez 261 subroutine writehist(vcov, ucov, teta, pk, phi, q, masse, ps, time)
8 guez 56
9 guez 261 ! From writehist.F, revision 1403, 2010-07-01 09:02:53
10     ! Écriture du fichier histoire au format IOIPSL
11 guez 69 ! L. Fairhead, LMD, 03/99
12 guez 56
13 guez 261 USE comconst, ONLY: cpp
14 guez 212 use covnat_m, only: covnat
15 guez 178 use dimens_m, only: llm
16 guez 212 use histsync_m, only: histsync
17     use histwrite_m, only: histwrite
18 guez 261 use iniadvtrac_m, only: ttext
19 guez 212 use inithist_m, only: histid, histvid, histuid
20 guez 261 use nr_util, only: assert
21 guez 178 use paramet_m, only: ip1jm, ip1jmp1
22 guez 69 use temps, only: itau_dyn
23 guez 56
24 guez 261 ! 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     real, intent(in):: phi(:, :, :) ! (iim + 1, jjm + 1, llm) ! geopotential
33     REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) traceurs
34     real, intent(in):: masse(:, :, :) ! (iim + 1, jjm + 1, llm)
35     REAL, intent(in):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol
36 guez 212 integer, intent(in):: time ! temps de l'ecriture
37 guez 56
38 guez 212 ! Local:
39 guez 261 integer iq, itau_w
40 guez 69 REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
41    
42     !---------------------------------------------------------------------
43    
44 guez 261 call assert([size(vcov, 1), size(ucov, 1), size(teta, 1), size(phi, 1), &
45     size(pk, 1), size(ps, 1), size(masse, 1)] == size(q, 1), &
46     "writedynav iim")
47     call assert([size(vcov, 2) + 1, size(ucov, 2), size(teta, 2), &
48     size(phi, 2), size(pk, 2), size(ps, 2), size(masse, 2)] &
49     == size(q, 2), "writedynav jjm")
50     call assert([size(vcov, 3), size(ucov, 3), size(teta, 3), size(phi, 3), &
51     size(pk, 3), size(masse, 3), size(q, 3)] == llm, "writedynav llm")
52 guez 69
53     itau_w = itau_dyn + time
54     call covnat(llm, ucov, vcov, unat, vnat)
55    
56     call histwrite(histuid, 'u', itau_w, unat)
57     call histwrite(histvid, 'v', itau_w, vnat)
58 guez 261 call histwrite(histid, 'theta', itau_w, teta)
59     call histwrite(histid, 'temp', itau_w, teta * pk / cpp)
60 guez 69 call histwrite(histid, 'phi', itau_w, phi)
61    
62 guez 261 DO iq = 1, size(q, 4)
63     call histwrite(histid, ttext(iq), itau_w, q(:, :, :, iq))
64     enddo
65 guez 69
66     call histwrite(histid, 'masse', itau_w, masse)
67     call histwrite(histid, 'ps', itau_w, ps)
68    
69 guez 261 call histsync(histid)
70     call histsync(histvid)
71     call histsync(histuid)
72 guez 69
73     end subroutine writehist
74    
75     end module writehist_m

  ViewVC Help
Powered by ViewVC 1.1.21