source: LMDZ6/trunk/libf/dyn3dmem/initfluxsto_p.F90 @ 5255

Last change on this file since 5255 was 5246, checked in by abarral, 29 hours ago

Convert fixed-form to free-form sources .F -> .{f,F}90
(WIP: some .F remain, will be handled in subsequent commits)

  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 7.0 KB
Line 
1!
2! $Id$
3!
4subroutine initfluxsto_p &
5        (infile,tstep,t_ops,t_wrt, &
6        fileid,filevid,filedid)
7
8#ifdef CPP_IOIPSL
9  ! This routine needs IOIPSL
10   USE IOIPSL
11#endif
12   USE parallel_lmdz
13   use Write_field
14   use misc_mod
15   USE comconst_mod, ONLY: pi
16   USE comvert_mod, ONLY: nivsigs
17   USE temps_mod, ONLY: annee_ref, day_ref, itau_dyn
18
19  implicit none
20
21  !
22  !   Routine d'initialisation des ecritures des fichiers histoires LMDZ
23  !   au format IOIPSL
24  !
25  !   Appels succesifs des routines: histbeg
26  !                              histhori
27  !                              histver
28  !                              histdef
29  !                              histend
30  !
31  !   Entree:
32  !
33  !  infile: nom du fichier histoire a creer
34  !  day0,anne0: date de reference
35  !  tstep: duree du pas de temps en seconde
36  !  t_ops: frequence de l'operation pour IOIPSL
37  !  t_wrt: frequence d'ecriture sur le fichier
38  !
39  !   Sortie:
40  !  fileid: ID du fichier netcdf cree
41  !  filevid:ID du fichier netcdf pour la grille v
42  !
43  !   L. Fairhead, LMD, 03/99
44  !
45  ! =====================================================================
46  !
47  !   Declarations
48  include "dimensions.h"
49  include "paramet.h"
50  include "comgeom.h"
51  include "description.h"
52  include "iniprint.h"
53
54  !   Arguments
55  !
56  character(len=*) :: infile
57  real :: tstep, t_ops, t_wrt
58  integer :: fileid, filevid,filedid
59
60#ifdef CPP_IOIPSL
61  ! This routine needs IOIPSL
62  !   Variables locales
63  !
64  real :: nivd(1)
65  integer :: tau0
66  real :: zjulian
67  character(len=3) :: str
68  character(len=10) :: ctrac
69  integer :: iq
70  real :: rlong(iip1,jjp1), rlat(iip1,jjp1),rl(1,1)
71  integer :: uhoriid, vhoriid, thoriid, zvertiid,dhoriid,dvertiid
72  integer :: ii,jj
73  integer :: zan, idayref
74  logical :: ok_sync
75  integer :: jjb,jje,jjn
76
77  ! definition du domaine d'ecriture pour le rebuild
78
79  INTEGER,DIMENSION(2) :: ddid
80  INTEGER,DIMENSION(2) :: dsg
81  INTEGER,DIMENSION(2) :: dsl
82  INTEGER,DIMENSION(2) :: dpf
83  INTEGER,DIMENSION(2) :: dpl
84  INTEGER,DIMENSION(2) :: dhs
85  INTEGER,DIMENSION(2) :: dhe
86
87  INTEGER :: dynu_domain_id
88  INTEGER :: dynv_domain_id
89
90  !
91  !  Initialisations
92  !
93  pi = 4. * atan (1.)
94  str='q  '
95  ctrac = 'traceur   '
96  ok_sync = .true.
97  !
98  !  Appel a histbeg: creation du fichier netcdf et initialisations diverses
99  !
100
101  zan = annee_ref
102  idayref = day_ref
103  CALL ymds2ju(zan, 1, idayref, 0.0, zjulian)
104  tau0 = itau_dyn
105
106    do jj = 1, jjp1
107    do ii = 1, iip1
108      rlong(ii,jj) = rlonu(ii) * 180. / pi
109      rlat(ii,jj) = rlatu(jj) * 180. / pi
110    enddo
111  enddo
112
113  jjb=jj_begin
114  jje=jj_end
115  jjn=jj_nb
116
117  ddid=(/ 1,2 /)
118  dsg=(/ iip1,jjp1 /)
119  dsl=(/ iip1,jjn /)
120  dpf=(/ 1,jjb /)
121  dpl=(/ iip1,jje /)
122  dhs=(/ 0,0 /)
123  dhe=(/ 0,0 /)
124
125  call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
126        'box',dynu_domain_id)
127
128  call histbeg(trim(infile),iip1, rlong(:,1), jjn, rlat(1,jjb:jje), &
129        1, iip1, 1, jjn, tau0, zjulian, tstep, uhoriid, &
130        fileid,dynu_domain_id)
131  !
132  !  Creation du fichier histoire pour la grille en V (oblige pour l'instant,
133  !  IOIPSL ne permet pas de grilles avec des nombres de point differents dans
134  !  un meme fichier)
135
136
137  do jj = 1, jjm
138    do ii = 1, iip1
139      rlong(ii,jj) = rlonv(ii) * 180. / pi
140      rlat(ii,jj) = rlatv(jj) * 180. / pi
141    enddo
142  enddo
143
144  jjb=jj_begin
145  jje=jj_end
146  jjn=jj_nb
147  if (pole_sud) jje=jj_end-1
148  if (pole_sud) jjn=jj_nb-1
149
150  ddid=(/ 1,2 /)
151  dsg=(/ iip1,jjm /)
152  dsl=(/ iip1,jjn /)
153  dpf=(/ 1,jjb /)
154  dpl=(/ iip1,jje /)
155  dhs=(/ 0,0 /)
156  dhe=(/ 0,0 /)
157
158  call flio_dom_set(mpi_size,mpi_rank,ddid,dsg,dsl,dpf,dpl,dhs,dhe, &
159        'box',dynv_domain_id)
160
161  call histbeg('fluxstokev',iip1, rlong(:,1), jjn, rlat(1,jjb:jje), &
162        1, iip1, 1, jjn,tau0, zjulian, tstep, vhoriid, &
163        filevid,dynv_domain_id)
164
165  rl(1,1) = 1.
166
167  if (mpi_rank==0) then
168
169    call histbeg('defstoke.nc', 1, rl, 1, rl, &
170          1, 1, 1, 1, &
171          tau0, zjulian, tstep, dhoriid, filedid)
172
173  endif
174  !
175  !  Appel a histhori pour rajouter les autres grilles horizontales
176  !
177  do jj = 1, jjp1
178    do ii = 1, iip1
179      rlong(ii,jj) = rlonv(ii) * 180. / pi
180      rlat(ii,jj) = rlatu(jj) * 180. / pi
181    enddo
182  enddo
183
184  jjb=jj_begin
185  jje=jj_end
186  jjn=jj_nb
187
188  call histhori(fileid, iip1, rlong(:,jjb:jje),jjn,rlat(:,jjb:jje), &
189        'scalar','Grille points scalaires', thoriid)
190
191  !
192  !  Appel a histvert pour la grille verticale
193  !
194  call histvert(fileid, 'sig_s', 'Niveaux sigma', &
195        'sigma_level', &
196        llm, nivsigs, zvertiid)
197  ! Pour le fichier V
198  call histvert(filevid, 'sig_s', 'Niveaux sigma', &
199        'sigma_level', &
200        llm, nivsigs, zvertiid)
201  ! pour le fichier def
202  if (mpi_rank==0) then
203     nivd(1) = 1
204     call histvert(filedid, 'sig_s', 'Niveaux sigma', &
205           'sigma_level', &
206           1, nivd, dvertiid)
207  endif
208  !
209  !  Appels a histdef pour la definition des variables a sauvegarder
210
211    CALL histdef(fileid, "phis", "Surface geop. height", "-", &
212          iip1,jjn,thoriid, 1,1,1, -99, 32, &
213          "once", t_ops, t_wrt)
214
215     CALL histdef(fileid, "aire", "Grid area", "-", &
216           iip1,jjn,thoriid, 1,1,1, -99, 32, &
217           "once", t_ops, t_wrt)
218
219    if (mpi_rank==0) then
220
221    CALL histdef(filedid, "dtvr", "tps dyn", "s", &
222          1,1,dhoriid, 1,1,1, -99, 32, &
223          "once", t_ops, t_wrt)
224
225     CALL histdef(filedid, "istdyn", "tps stock", "s", &
226           1,1,dhoriid, 1,1,1, -99, 32, &
227           "once", t_ops, t_wrt)
228
229     CALL histdef(filedid, "istphy", "tps stock phy", "s", &
230           1,1,dhoriid, 1,1,1, -99, 32, &
231           "once", t_ops, t_wrt)
232
233    endif
234  !
235  ! Masse
236  !
237  call histdef(fileid, 'masse', 'Masse', 'kg', &
238        iip1, jjn, thoriid, llm, 1, llm, zvertiid, &
239        32, 'inst(X)', t_ops, t_wrt)
240  !
241  !  Pbaru
242  !
243  call histdef(fileid, 'pbaru', 'flx de masse zonal', 'kg m/s', &
244        iip1, jjn, uhoriid, llm, 1, llm, zvertiid, &
245        32, 'inst(X)', t_ops, t_wrt)
246
247  !
248  !  Pbarv
249  !
250  if (pole_sud) jjn=jj_nb-1
251
252  call histdef(filevid, 'pbarv', 'flx de masse mer', 'kg m/s', &
253        iip1, jjn, vhoriid, llm, 1, llm, zvertiid, &
254        32, 'inst(X)', t_ops, t_wrt)
255  !
256  !  w
257  !
258  if (pole_sud) jjn=jj_nb
259  call histdef(fileid, 'w', 'flx de masse vert', 'kg m/s', &
260        iip1, jjn, thoriid, llm, 1, llm, zvertiid, &
261        32, 'inst(X)', t_ops, t_wrt)
262
263  !
264  !  Temperature potentielle
265  !
266  call histdef(fileid, 'teta', 'temperature potentielle', '-', &
267        iip1, jjn, thoriid, llm, 1, llm, zvertiid, &
268        32, 'inst(X)', t_ops, t_wrt)
269  !
270
271  !
272  ! Geopotentiel
273  !
274  call histdef(fileid, 'phi', 'geopotentiel instantane', '-', &
275        iip1, jjn, thoriid, llm, 1, llm, zvertiid, &
276        32, 'inst(X)', t_ops, t_wrt)
277  !
278  !  Fin
279  !
280  call histend(fileid)
281  call histend(filevid)
282  if (mpi_rank==0) call histend(filedid)
283  if (ok_sync) then
284    call histsync(fileid)
285    call histsync(filevid)
286    if (mpi_rank==0) call histsync(filedid)
287  endif
288
289#else
290  write(lunout,*)'initfluxsto_p: Needs IOIPSL to function'
291#endif
292  ! #endif of #ifdef CPP_IOIPSL
293  return
294end subroutine initfluxsto_p
Note: See TracBrowser for help on using the repository browser.