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

Contents of /trunk/dyn3d/writehist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 262 - (show annotations)
Wed Mar 7 13:46:18 2018 UTC (6 years, 2 months ago) by guez
File size: 2568 byte(s)
Remove output file dynhist_ave.nc. As in physics, keep only one output
file for instantaneous output, pending replacement by a more powerful
solution like XIOS.

1 module writehist_m
2
3 implicit none
4
5 contains
6
7 subroutine writehist(vcov, ucov, teta, pk, phi, q, masse, ps, time)
8
9 ! From writehist.F, revision 1403, 2010-07-01 09:02:53
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
16 use histsync_m, only: histsync
17 use histwrite_m, only: histwrite
18 use iniadvtrac_m, only: ttext
19 use inithist_m, only: histid, histvid, histuid
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 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 integer, intent(in):: time ! temps de l'ecriture
37
38 ! Local:
39 integer iq, itau_w
40 REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
41
42 !---------------------------------------------------------------------
43
44 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 "writehist 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), "writehist 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, "writehist llm")
52
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 call histwrite(histid, 'theta', itau_w, teta)
59 call histwrite(histid, 'temp', itau_w, teta * pk / cpp)
60 call histwrite(histid, 'phi', itau_w, phi)
61
62 DO iq = 1, size(q, 4)
63 call histwrite(histid, ttext(iq), itau_w, q(:, :, :, iq))
64 enddo
65
66 call histwrite(histid, 'masse', itau_w, masse)
67 call histwrite(histid, 'ps', itau_w, ps)
68
69 call histsync(histid)
70 call histsync(histvid)
71 call histsync(histuid)
72
73 end subroutine writehist
74
75 end module writehist_m

  ViewVC Help
Powered by ViewVC 1.1.21