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

Annotation of /trunk/dyn3d/writehist.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 263 - (hide annotations)
Wed Mar 7 14:41:46 2018 UTC (6 years, 2 months ago) by guez
File size: 2505 byte(s)
Allow output of dynamics at frequency greater than one per day. Change
meaning of iecri.

Remove arguments tstep, nq of procedure inithist: use directly
module variables dtvr and nqmx instead.

Change meaning of last dummy argument of writehist. Simpler and
clearer, because time suggested a dimensional quantity.

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 263 subroutine writehist(vcov, ucov, teta, pk, phi, q, masse, ps, itau_w)
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 56
23 guez 261 ! Vent covariant :
24     REAL, intent(in):: vcov(:, :, :) ! (iim + 1, jjm, llm)
25     REAL, intent(in):: ucov(:, :, :) ! (iim + 1, jjm + 1, llm)
26    
27     REAL, intent(in):: teta(:, :, :) ! (iim + 1, jjm + 1, llm)
28     ! temperature potentielle
29    
30     real, intent(in):: pk(:, :, :) ! (iim + 1, jjm + 1, llm)
31     real, intent(in):: phi(:, :, :) ! (iim + 1, jjm + 1, llm) ! geopotential
32     REAL, intent(in):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx) traceurs
33     real, intent(in):: masse(:, :, :) ! (iim + 1, jjm + 1, llm)
34     REAL, intent(in):: ps(:, :) ! (iim + 1, jjm + 1) pression au sol
35 guez 263 integer, intent(in):: itau_w ! temps de l'ecriture
36 guez 56
37 guez 212 ! Local:
38 guez 263 integer iq
39 guez 69 REAL vnat(ip1jm, llm), unat(ip1jmp1, llm)
40    
41     !---------------------------------------------------------------------
42    
43 guez 261 call assert([size(vcov, 1), size(ucov, 1), size(teta, 1), size(phi, 1), &
44     size(pk, 1), size(ps, 1), size(masse, 1)] == size(q, 1), &
45 guez 262 "writehist iim")
46 guez 261 call assert([size(vcov, 2) + 1, size(ucov, 2), size(teta, 2), &
47     size(phi, 2), size(pk, 2), size(ps, 2), size(masse, 2)] &
48 guez 262 == size(q, 2), "writehist jjm")
49 guez 261 call assert([size(vcov, 3), size(ucov, 3), size(teta, 3), size(phi, 3), &
50 guez 262 size(pk, 3), size(masse, 3), size(q, 3)] == llm, "writehist llm")
51 guez 69
52     call covnat(llm, ucov, vcov, unat, vnat)
53    
54     call histwrite(histuid, 'u', itau_w, unat)
55     call histwrite(histvid, 'v', itau_w, vnat)
56 guez 261 call histwrite(histid, 'theta', itau_w, teta)
57     call histwrite(histid, 'temp', itau_w, teta * pk / cpp)
58 guez 69 call histwrite(histid, 'phi', itau_w, phi)
59    
60 guez 261 DO iq = 1, size(q, 4)
61     call histwrite(histid, ttext(iq), itau_w, q(:, :, :, iq))
62     enddo
63 guez 69
64     call histwrite(histid, 'masse', itau_w, masse)
65     call histwrite(histid, 'ps', itau_w, ps)
66    
67 guez 261 call histsync(histid)
68     call histsync(histvid)
69     call histsync(histuid)
70 guez 69
71     end subroutine writehist
72    
73     end module writehist_m

  ViewVC Help
Powered by ViewVC 1.1.21