The TELEMAC-MASCARET system  trunk
utimp_drogues.f
Go to the documentation of this file.
1 ! ************************
2  SUBROUTINE utimp_drogues
3 ! ************************
4 !
5  &( ltl,atl,npoin2,npoin3, xflot,yflot,zflot,tagflo,clsflo,
6  & nflot,nflot_max,floprd,deja, t2dflo,t2dblo, mardat,martim )
7 !
8 !***********************************************************************
9 ! TELEMAC2D V8P1
10 !***********************************************************************
11 !
12 !brief Specific to drogues.
13 !+ Allow multiple type of file outputs for drogues
14 !
15 !history S.E.BOURBAN (HRW)
16 !+ 13/07/2018
17 !+ V8P0
18 !+ Initial implementation to support TecPlot and BlueKenue parcel
19 !+ files.
20 !
21 !history M.S.TURNBULL (HRW)
22 !+ 04/11/2019
23 !+ V8P2
24 ! Corrections
25 !
26 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
27 !| ATL |-->| TIME OF TIME STEP, IN SECONDS
28 !| CLSFLO |<->| CLASS OF FLOATS
29 !| DEJA |<->| INDICATES FIRST CALL OF SUBROUTINE
30 !| MARDAT |-->| DATE (YEAR, MONTH,DAY)
31 !| MARTIM |-->| TIME (HOUR, MINUTE,SECOND)
32 !| FLOPRD |-->| PERIOD OF LISTING OUTPUTS
33 !| LTL |-->| CURRENT TIME STEP
34 !| LUFLO |-->| LOGICAL UNIT FOR THE ASCII OUTPUT FILE
35 !| LUBLO |-->| LOGICAL UNIT FOR THE BINARY OUTPUT FILE
36 !| FMTBLO |-->| FORMAT FOR THE OUTPUT FILE
37 !| NPOIN2 |-->| NUMBER OF POINTS IN 2D MESH
38 !| NPOIN3 |-->| NUMBER OF POINTS
39 !| NFLOT |-->| NUMBER OF FLOATS.
40 !| NFLOT_MAX |-->| MAXIMUM NUMBER OF FLOATS.
41 !| TAGFLO |<->| TAGS OF FLOATS
42 !| XFLOT |<->| ABSCISSAE OF FLOATS
43 !| YFLOT |<->| ORDINATES OF FLOATS
44 !| ZFLOT |<->| ELEVATIONS OF FLOATS
45 !| T2DFLO |-->| ASCII DROGUES FILE
46 !| T2DBLO |-->| BINARY DROGUES FILE
47 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
48 !
49  USE bief
50 !
52  USE declarations_telemac, ONLY: iframe
53  USE interface_parallel, ONLY : p_sum
54  IMPLICIT NONE
55 !
56 !+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
57 !
58  DOUBLE PRECISION, INTENT(IN) :: ATL
59  INTEGER , INTENT(IN) :: FLOPRD,LTL
60  INTEGER , INTENT(IN) :: MARDAT(3),MARTIM(3)
61  INTEGER , INTENT(IN) :: NFLOT,NFLOT_MAX,NPOIN2,NPOIN3
62  DOUBLE PRECISION, INTENT(INOUT) :: XFLOT(nflot_max)
63  DOUBLE PRECISION, INTENT(INOUT) :: YFLOT(nflot_max)
64  DOUBLE PRECISION, INTENT(INOUT) :: ZFLOT(nflot_max)
65  INTEGER , INTENT(INOUT) :: TAGFLO(nflot_max)
66  INTEGER , INTENT(INOUT) :: CLSFLO(nflot_max)
67  TYPE(bief_file) , INTENT(IN) :: T2DFLO,T2DBLO
68  LOGICAL , INTENT(INOUT) :: DEJA
69 !
70 !+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
71 !
72  INTEGER ID, LUFLO,LUBLO, IPROC,IFLOT, NFLOTG
73  INTEGER I1,I2
74  INTEGER NFLOT_START, NFLOTG_START
75  DOUBLE PRECISION V1
76  LOGICAL YESITDOES
77  CHARACTER(LEN=32) TEXTE(3)
78  REAL R1,R2,R3
79  INTEGER CURDAT(3),CURTIM(3),MILSEC
80  DOUBLE PRECISION JAT
81 !
82  CHARACTER(LEN=11) EXTENS
83  EXTERNAL extens
84 !
85 !+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+
86 !
87 !
88 !-----------------------------------------------------------------------
89 !
90 ! INITIALISING HEADER OF THE OUTPUT FILE
91 !
92  IF(.NOT.deja) THEN
93 !
94 ! WRITING / COMBINING IN ONLY ONE OUTPUT FILE
95  IF(ipid.EQ.0) THEN
96 !
97 ! HEADER OF TECPLOT FILE
98  IF( t2dblo%FMT .EQ.'TECPLOT ' .AND. t2dflo%NAME.NE.'' ) THEN
99  luflo = t2dflo%LU
100  texte(1)='X '
101  texte(2)='Y '
102  texte(3)='ELEVATION Z M '
103  WRITE(luflo,100) 'TITLE = "DROGUES FILE"'
104  IF( npoin2.EQ.npoin3 ) THEN
105  WRITE(luflo,100) 'VARIABLES = "LABELS","'//
106  & texte(1)//'","'//texte(2)//
107  & '","COLOUR"'
108  ELSE
109  WRITE(luflo,100) 'VARIABLES = "LABELS","'//
110  & texte(1)//'","'//texte(2)//'","'//texte(3)//
111  & '","COLOUR"'
112  ENDIF
113  ENDIF
114 !
115 ! HEADER OF BLUE KENUE PARCEL FILE IN ASCII FORM
116  IF( t2dblo%FMT .EQ.'BKBINPCL' .AND. t2dblo%NAME.NE.'' ) THEN
117  iframe = 0
118  lublo = t2dblo%LU
119  CLOSE(lublo)
120  OPEN( lublo, file=trim(t2dblo%TELNAME),form='FORMATTED',
121  & access='STREAM', action='WRITE')
122  WRITE(lublo,100) ':FileType pcl BINARY EnSim 1.0'
123  WRITE(lublo,100) ':WrittenBy TELEMAC'
124  WRITE(lublo,100) ':CreationTime 2018/08/21 14:16:28.000'
125  WRITE(lublo,100) ':Name ' // 'DROGUES_FILE'
126  IF( npoin2.EQ.npoin3 ) THEN
127  WRITE(lublo,100) ':AttributeCount 1 TAG'
128  WRITE(lublo,100) ':AttributeCount 2 CLASS'
129  WRITE(lublo,100) ':AttributeUnits 1 #'
130  WRITE(lublo,100) ':AttributeUnits 2 #'
131  ELSE
132  WRITE(lublo,100) ':AttributeCount 1 Z'
133  WRITE(lublo,100) ':AttributeCount 2 TAG'
134  WRITE(lublo,100) ':AttributeCount 3 CLASS'
135  WRITE(lublo,100) ':AttributeUnits 1 M'
136  WRITE(lublo,100) ':AttributeUnits 2 #'
137  WRITE(lublo,100) ':AttributeUnits 3 #'
138  ENDIF
139  WRITE(lublo,100) ':EndHeader'
140  CLOSE(lublo)
141 !
142 ! CORE OF BLUE KENUE PARCEL FILE IN BINARY LITTLE ENDIAN FORM
143  OPEN( lublo, file=trim(t2dblo%TELNAME),form='UNFORMATTED',
144  & action='WRITE',convert='LITTLE_ENDIAN',
145  & access='STREAM', status='OLD',position='APPEND')
146 !
147  ENDIF
148 !
149  ENDIF
150 !
151  deja = .true.
152  ENDIF
153 !
154 !-----------------------------------------------------------------------
155 !
156  nflot_start = nflot_max
157  nflotg_start = nflot_start
158 !
159  IF( (ltl/floprd)*floprd.NE.ltl ) RETURN
160 !
161  curdat = mardat
162  curtim = martim
163  milsec = int( 1000.d0*( atl - int(atl) ) )
164  jat = jultim(curdat(1),curdat(2),curdat(3),
165  & curtim(1),curtim(2),curtim(3),(1.d0*int(atl)) )
166  CALL gregtim(jat,curdat(1),curdat(2),curdat(3),
167  & curtim(1),curtim(2),curtim(3))
168 !
169 !-----------------------------------------------------------------------
170 !
171 ! TECPLOT AND BLUE KENUE PARCEL FILES
172 !
173  IF( ncsize.GT.1 ) THEN
174 !
175 ! WAITING ALL PROCESSORS (SO THAT NFLOT IS UPDATED FOR ALL
176 ! BEFORE CALLING P_SUM)
177 !
178 ! CALL P_SYNC
179 !
180 ! 1) PARALLEL WRITE UP
181 !
182  nflotg = nflot
183  nflotg = p_sum(nflot)
184  IF( nflotg.GT.0 ) THEN
185 !
186 ! ONE TIME STEP OF THE TECPLOT FILE
187  IF( t2dblo%FMT .EQ.'TECPLOT ' .AND. t2dflo%NAME.NE.'' ) THEN
188  luflo = t2dflo%LU
189 !
190 ! 1) EVERY PROCESSOR WRITES ITS OWN DATA IN A FILE
191 !
192  CALL get_free_id(id)
193  IF(nflot.GT.0) THEN
194  OPEN(id,file=extens(ncsize,ipid+1),form='UNFORMATTED',
195  & status='NEW',action='WRITE',access='STREAM')
196  IF( npoin2.EQ.npoin3 ) THEN
197  DO iflot = 1,nflot
198  WRITE(id) tagflo(iflot),REAL(XFLOT(IFLOT)),
199  & REAL(YFLOT(IFLOT)),1
200  ENDDO
201  ELSE
202  DO iflot = 1,nflot
203  WRITE(id) tagflo(iflot),REAL(XFLOT(IFLOT)),
204  & REAL(YFLOT(IFLOT)),REAL(ZFLOT(IFLOT)),1
205  ENDDO
206  ENDIF
207  CLOSE(id)
208  ENDIF
209 !
210 ! 2) WAITING ALL PROCESSORS
211 !
212  CALL p_sync
213 !
214 ! 3) PROCESSOR 0 READS ALL EXISTING FILES AND MERGES
215 ! THEM IN THE FINAL FILE
216 !
217  IF(ipid.EQ.0) THEN
218 !
219  WRITE(luflo,200) 'ZONE DATAPACKING=POINT, T="G_',atl,
220  & ' SECONDS"',', I=',nflotg,', SOLUTIONTIME=',atl
221  DO iproc = 1,ncsize
222  INQUIRE(file=extens(ncsize,iproc),exist=yesitdoes)
223  IF( yesitdoes ) THEN
224  OPEN(id,file=extens(ncsize,iproc),form='UNFORMATTED',
225  & status='OLD',action='READ',access='STREAM')
226  IF( npoin2.EQ.npoin3 ) THEN
227  22 CONTINUE
228  READ(id,err=24,end=24) i1,r1,r2,i2
229  WRITE(luflo,300) i1,r1,r2,i2
230  GO TO 22
231  ELSE
232  23 CONTINUE
233  READ(id,err=24,end=24) i1,r1,r2,r3,i2
234  WRITE(luflo,301) i1,r1,r2,r3,i2
235  GO TO 23
236  ENDIF
237  24 CONTINUE
238  CLOSE(id,status='DELETE')
239  ENDIF
240  ENDDO
241  ENDIF
242 !
243  ENDIF
244 !
245 ! HEADER OF BLUE KENUE PARCEL FILE IN ASCII FORM
246  IF( t2dblo%FMT .EQ.'BKBINPCL' .AND. t2dblo%NAME.NE.'' ) THEN
247  lublo = t2dblo%LU
248 !
249 ! 1) EVERY PROCESSOR WRITES ITS OWN DATA IN A FILE
250 !
251  CALL get_free_id(id)
252  IF(nflot.GT.0) THEN
253  OPEN(id,file=extens(ncsize,ipid+1),form='UNFORMATTED',
254  & status='NEW',action='WRITE',access='STREAM')
255  WRITE(id) nflot
256  WRITE(id) ( tagflo(iflot), iflot = 1,nflot )
257  WRITE(id) ( clsflo(iflot), iflot = 1,nflot )
258  WRITE(id) ( xflot(iflot), iflot = 1,nflot )
259  WRITE(id) ( yflot(iflot), iflot = 1,nflot )
260  IF( npoin2.NE.npoin3 ) THEN
261  WRITE(id) ( zflot(iflot), iflot = 1,nflot )
262  ENDIF
263  CLOSE(id)
264  ENDIF
265 !
266 ! 2) WAITING ALL PROCESSORS
267 !
268  CALL p_sync
269 !
270 ! 3) PROCESSOR 0 READS ALL EXISTING FILES AND MERGES
271 ! THEM IN ORDER OF THEIR TAGS TO THE FINAL FILE
272 !
273 ! /!\ IT IS IMPORTANT TO NOTE THAT XFLOT,YFLOT, ETC.
274 ! FOR IPID=0 REMAIN UNCHANGED, EVEN IF VALUES
275 ! FROM OTHER PROCESSORS GET COPIED OVER.
276 !
277  IF(ipid.EQ.0) THEN
278 !
279  iframe = iframe + 1
280  WRITE(lublo) iframe,iframe,
281  & curdat(1),curdat(2),curdat(3),
282  & curtim(1),curtim(2),curtim(3),milsec
283  WRITE(lublo) nflotg_start
284  i2 = 0
285  DO iproc = 1,ncsize
286  INQUIRE(file=extens(ncsize,iproc),exist=yesitdoes)
287  IF( yesitdoes ) THEN
288  OPEN(id,file=extens(ncsize,iproc),form='UNFORMATTED',
289  & status='OLD',action='READ',access='STREAM')
290  READ(id) i1
291  READ(id) ( tagflo(i2+iflot), iflot = 1,i1 )
292  READ(id) ( clsflo(i2+iflot), iflot = 1,i1 )
293  READ(id) ( xflot(i2+iflot), iflot = 1,i1 )
294  READ(id) ( yflot(i2+iflot), iflot = 1,i1 )
295  IF( npoin2.NE.npoin3 ) THEN
296  READ(id) ( zflot(i2+iflot), iflot = 1,i1 )
297  ENDIF
298  CLOSE(id,status='DELETE')
299  i2 = i2 + i1
300  ENDIF
301  ENDDO
302 !
303  DO iflot = nflotg+1,nflotg_start
304  v1 = -1.d0
305  xflot(iflot) = v1
306  yflot(iflot) = v1
307  tagflo(iflot) = iflot
308  clsflo(iflot) = 1
309  END DO
310 !
311  WRITE(lublo)
312  & ( REAL(XFLOT(IFLOT)), IFLOT = 1,nflotg_start )
313  WRITE(lublo)
314  & ( REAL(YFLOT(IFLOT)), IFLOT = 1,nflotg_start )
315  IF( npoin2.NE.npoin3 ) THEN
316  WRITE(lublo)
317  & ( REAL(ZFLOT(IFLOT)), IFLOT = 1,nflotg_start )
318  ENDIF
319  WRITE(lublo)
320  & ( REAL(TAGFLO(IFLOT)), IFLOT = 1,nflotg_start )
321  WRITE(lublo)
322  & ( REAL(CLSFLO(IFLOT)), IFLOT = 1,nflotg_start )
323  WRITE(lublo) nflotg_start
324 !
325  ENDIF
326 !
327  ENDIF
328 !
329  ENDIF
330 !
331  ELSE
332 !
333 ! SCALAR VERSION
334 !
335  IF( nflot.GT.0 ) THEN
336 !
337 ! ONE TIME STEP OF THE TECPLOT FILE
338  IF( t2dblo%FMT .EQ.'TECPLOT ' .AND. t2dflo%NAME.NE.'' ) THEN
339  luflo = t2dflo%LU
340 !
341  WRITE(luflo,200) 'ZONE DATAPACKING=POINT, T="G_',atl,
342  & ' SECONDS"',', I=',nflot,', SOLUTIONTIME=',atl
343  IF( npoin2.EQ.npoin3 ) THEN
344  DO iflot = 1,nflot
345  WRITE(luflo,300) tagflo(iflot),
346  & xflot(iflot),yflot(iflot),1
347  ENDDO
348  ELSE
349  DO iflot = 1,nflot
350  WRITE(luflo,301) tagflo(iflot),
351  & xflot(iflot),yflot(iflot),zflot(iflot),1
352  ENDDO
353  ENDIF
354 !
355  ENDIF
356 !
357 ! ONE TIME STEP OF THE BLUE KENUE PARCEL FILE IN ASCII FORM
358  IF( t2dblo%FMT .EQ.'BKBINPCL' .AND. t2dblo%NAME.NE.'' ) THEN
359  lublo = t2dblo%LU
360 !
361  iframe = iframe + 1
362  WRITE(lublo) iframe,iframe,
363  & curdat(1),curdat(2),curdat(3),
364  & curtim(1),curtim(2),curtim(3),milsec
365  WRITE(lublo) nflot_start
366 !
367  DO iflot = nflot+1,nflot_start
368  v1 = -1.d0
369  xflot(iflot) = v1
370  yflot(iflot) = v1
371  END DO
372 !
373  WRITE(lublo)
374  & ( REAL(XFLOT(IFLOT)), IFLOT = 1,nflot_start )
375  WRITE(lublo)
376  & ( REAL(YFLOT(IFLOT)), IFLOT = 1,nflot_start )
377  IF( npoin2.NE.npoin3 ) THEN
378  WRITE(lublo)
379  & ( REAL(ZFLOT(IFLOT)), IFLOT = 1,nflot_start )
380  ENDIF
381  WRITE(lublo)
382  & ( REAL(TAGFLO(IFLOT)), IFLOT = 1,nflot_start )
383  WRITE(lublo)
384  & ( REAL(CLSFLO(IFLOT)), IFLOT = 1,nflot_start )
385  WRITE(lublo) nflot_start
386 !
387  ENDIF
388 !
389  ENDIF
390 !
391  ENDIF
392 !
393 !-----------------------------------------------------------------------
394 !
395  100 FORMAT(a)
396  200 FORMAT(a,f12.4,a,a,i4,a,f12.4)
397  300 FORMAT(i6,',',f16.8,',',f16.8,',',i2)
398  301 FORMAT(i6,',',f16.8,',',f16.8,',',f16.8,',',i2)
399 !
400 !-----------------------------------------------------------------------
401 !
402  RETURN
403  END SUBROUTINE utimp_drogues
subroutine utimp_drogues(LTL, ATL, NPOIN2, NPOIN3, XFLOT, YFLOT, ZFLOT, TAGFLO, CLSFLO, NFLOT, NFLOT_MAX, FLOPRD, DEJA, T2DFLO, T2DBLO, MARDAT, MARTIM)
Definition: utimp_drogues.f:8
subroutine gregtim(JULTIM, YEAR, MONTH, DAY, HOUR, MINU, SEC)
Definition: gregtim.f:7
double precision function jultim(YEAR, MONTH, DAY, HOUR, MINU, SEC, AT)
Definition: jultim.f:7
subroutine p_sync
Definition: p_sync.F:4
Definition: bief.f:3