SFEMaNS  version 4.1 (work in progress)
Reference documentation for SFEMaNS
 All Classes Files Functions Variables Groups Pages
sub_plot.f90
Go to the documentation of this file.
1 !
2 !Authors: Jean-Luc Guermond, Lugi Quartapelle, Copyright 1994
3 !
4 MODULE sub_plot
5 
6 CONTAINS
7 
8  SUBROUTINE plot_vit_2d(jj, rr, uu)
9  !---FORMAT PLTMTV
10 
11  IMPLICIT NONE
12 
13  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
14  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
15  REAL(KIND=8), DIMENSION(:), INTENT(INOUT) :: uu ! rigid translation
16 
17  INTEGER :: m, n1, n2, n3
18 
19  OPEN(unit = 20, file ='vit_u', form ='formatted', status ='unknown')
20 
21  WRITE (20, *) '$ DATA = CONTCURVE'
22  WRITE (20, *) '% contstyle = 2'
23  WRITE (20,*) '% meshplot = False'
24 
25  DO m = 1, SIZE(jj, 2)
26  ! n1 = jj(1, m)
27  ! n2 = jj(2, m)
28  ! n3 = jj(3, m)
29  ! WRITE (20, *) rr(1,n1), rr(2,n1), uu(n1), n1
30  ! WRITE (20, *) rr(1,n2), rr(2,n2), uu(n2), n2
31  ! WRITE (20, *) rr(1,n3), rr(2,n3), uu(n3), n3
32  ! WRITE (20, *)
33 
34  n1 = jj(1, m)
35  n2 = jj(6, m)
36  n3 = jj(5, m)
37  WRITE (20, *) rr(1,n1), rr(2,n1), uu(n1), n1
38  WRITE (20, *) rr(1,n2), rr(2,n2), uu(n2), n2
39  WRITE (20, *) rr(1,n3), rr(2,n3), uu(n3), n3
40  WRITE (20, *)
41 
42 
43  n1 = jj(2, m)
44  n2 = jj(4, m)
45  n3 = jj(6, m)
46  WRITE (20, *) rr(1,n1), rr(2,n1), uu(n1), n1
47  WRITE (20, *) rr(1,n2), rr(2,n2), uu(n2), n2
48  WRITE (20, *) rr(1,n3), rr(2,n3), uu(n3), n3
49  WRITE (20, *)
50 
51  n1 = jj(4, m)
52  n2 = jj(3, m)
53  n3 = jj(5, m)
54  WRITE (20, *) rr(1,n1), rr(2,n1), uu(n1), n1
55  WRITE (20, *) rr(1,n2), rr(2,n2), uu(n2), n2
56  WRITE (20, *) rr(1,n3), rr(2,n3), uu(n3), n3
57  WRITE (20, *)
58 
59  n1 = jj(5, m)
60  n2 = jj(6, m)
61  n3 = jj(4, m)
62  WRITE (20, *) rr(1,n1), rr(2,n1), uu(n1), n1
63  WRITE (20, *) rr(1,n2), rr(2,n2), uu(n2), n2
64  WRITE (20, *) rr(1,n3), rr(2,n3), uu(n3), n3
65  WRITE (20, *)
66  ENDDO
67 
68  END SUBROUTINE plot_vit_2d
69 
70  !----------------------------------------------------------------
71 
72  SUBROUTINE plot_arrow(jj, rr, vv)
73  !---FORMAT PLTMTV
74 
75  IMPLICIT NONE
76 
77  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
78  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr, vv
79 
80  INTEGER :: n
81  REAL(KIND=8) :: z_d, w_d
82 
83  OPEN(unit = 20,file ='velocity',form ='formatted',status ='unknown')
84 
85  WRITE (20, 9000) SIZE(rr,2), SIZE(jj,1), SIZE(jj,2)
86 9000 FORMAT('#',3i10)
87  WRITE (20, *) '$ DATA = VECTOR'
88  WRITE (20, *) '% axisscale = FALSE'
89  WRITE (20, *) '% vscale = 1.'
90 
91  ! WRITE (20, *) '% contstyle = 2'
92  ! WRITE (20, *) '% meshplot = True'
93 
94  SELECT CASE(SIZE(rr,1))
95 
96  CASE(2)
97  z_d = 0
98  w_d = 0
99  DO n = 1, SIZE(rr, 2)
100  WRITE (20,*) rr(1,n), rr(2,n), z_d, vv(1,n), vv(2,n), w_d
101  ENDDO
102 
103  CASE(3)
104  DO n = 1, SIZE(rr, 2)
105  WRITE (20,*) rr(:,n), vv(:,n)
106  ENDDO
107 
108  END SELECT
109 
110  CLOSE(20)
111 
112  END SUBROUTINE plot_arrow
113 
114  !------------------------------------------------------------------------------
115 
116  SUBROUTINE plot_pressure_2d(jj, rr, uu)
117  !---FORMAT PLTMTV
118 
119  IMPLICIT NONE
120 
121  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
122  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
123  REAL(KIND=8), DIMENSION(:), INTENT(INOUT) :: uu ! rigid translation
124 
125  INTEGER :: na, nb, n, m, n1, n2, n3
126  REAL(KIND=8) :: xa, xb, xm, r, p_midpoint
127 
128  OPEN (unit = 20, file = 'pressure_2D.plt', form = 'formatted', status = 'unknown')
129 
130  ! set pressure field so as to have a zero pressure at
131  ! the mid-point of the bottom unit side of the cavity
132 
133  ! find the mid-point
134 
135  xa = 0; xb = 1; xm = 0.5
136  na = 1; nb = 1
137  DO n = 1, SIZE(uu)
138 
139  IF (rr(2,n) < 1.0d-4) THEN
140 
141  r = rr(1,n)
142  IF (r < xm .AND. r > xa) then; xa = r; na = n
143  ENDIF
144  IF (r > xm .AND. r < xb) then; xb = r; nb = n
145  ENDIF
146 
147  ENDIF
148 
149  ENDDO
150 
151  p_midpoint = uu(na) + (uu(nb) - uu(na)) * (xm-xa)/(xb-xa)
152 
153  ! reset pressure field
154 
155  uu = uu - p_midpoint
156 
157  WRITE (20, *) '$ DATA = CONTCURVE'
158  WRITE (20, *) '% contstyle = 1'
159  WRITE (20, *) '% contours = "( -.12 -.06 -.03 )"'
160  WRITE (20, *) '% contours = "( -.015 -.0075 -.00375 )"'
161  WRITE (20, *) '% contours = "( 0.0 )"'
162  WRITE (20, *) '% contours = "( .12 .06 .03 )"'
163  WRITE (20, *) '% contours = "( .015 .0075 .00375 )"'
164 
165  WRITE (20,*) '% meshplot = False'
166 
167  DO m = 1, SIZE(jj, 2)
168  n1 = jj(1, m)
169  n2 = jj(2, m)
170  n3 = jj(3, m)
171  WRITE (20, *) rr(1,n1), rr(2,n1), uu(n1), n1
172  WRITE (20, *) rr(1,n2), rr(2,n2), uu(n2), n2
173  WRITE (20, *) rr(1,n3), rr(2,n3), uu(n3), n3
174  WRITE (20, *)
175  ENDDO
176  CLOSE(20)
177 
178  END SUBROUTINE plot_pressure_2d
179 
180  !------------------------------------------------------------------------------
181 
182  SUBROUTINE plot_loc_rel_var(jj, rr, vo, vv)
183  !---FORMAT PLTMTV
184 
185  ! plot local relative variation of speed
186 
187 
188  IMPLICIT NONE
189 
190  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
191  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr, vo, vv
192 
193  REAL(KIND=8), DIMENSION(SIZE(vv,2)) :: uu
194  REAL(KIND=8) :: vo_n_mod
195  INTEGER :: n, m, n1, n2, n3
196 
197  OPEN (unit=20, file ='speed_var',form ='formatted',status ='unknown')
198 
199 
200  DO n = 1, SIZE(vv,2)
201 
202  vo_n_mod = sqrt(sum(vo(:,n)**2))
203 
204  IF (vo_n_mod < 1.0d-15) THEN
205 
206  uu(n) = 0
207 
208  ELSE
209 
210  uu(n) = sqrt(sum( (vo(:,n) - vv(:,n))**2 ))/vo_n_mod
211 
212  ENDIF
213 
214  ENDDO
215 
216 
217  WRITE (20, *) '$ DATA = CONTCURVE'
218  WRITE (20, *) '% contstyle = 2'
219  WRITE (20, *) '% meshplot = True'
220 
221  DO m = 1, SIZE(jj, 2)
222  n1 = jj(1, m)
223  n2 = jj(2, m)
224  n3 = jj(3, m)
225  WRITE (20, *) rr(1,n1), rr(2,n1), uu(n1), n1
226  WRITE (20, *) rr(1,n2), rr(2,n2), uu(n2), n2
227  WRITE (20, *) rr(1,n3), rr(2,n3), uu(n3), n3
228  WRITE (20, *)
229  ENDDO
230 
231  END SUBROUTINE plot_loc_rel_var
232 
233  !------------------------------------------------------------------------------
234 
235  SUBROUTINE plot_pressure(jj, rr, uu)
236  !---FORMAT PLTMTV
237 
238  IMPLICIT NONE
239 
240  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
241  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
242  REAL(KIND=8), DIMENSION(:), INTENT(INOUT) :: uu ! rigid translation
243  INTEGER :: m, n1, n2, n3
244 
245  OPEN (unit=20, file='pressure', form='formatted', status='unknown')
246 
247  WRITE (20, *) '$ DATA = CONTCURVE'
248  WRITE (20, *) '% contstyle = 2'
249  WRITE (20, *) '% cstep = 20'
250  WRITE (20, *) '% meshplot = False'
251 
252  DO m = 1, SIZE(jj, 2)
253  n1 = jj(1, m)
254  n2 = jj(2, m)
255  n3 = jj(3, m)
256  WRITE (20, *) rr(1,n1), rr(2,n1), uu(n1), n1
257  WRITE (20, *) rr(1,n2), rr(2,n2), uu(n2), n2
258  WRITE (20, *) rr(1,n3), rr(2,n3), uu(n3), n3
259  WRITE (20, *)
260  ENDDO
261 
262  END SUBROUTINE plot_pressure
263 
264  SUBROUTINE plot_const_p1_label(jj, rr, uu, file_name )
265  !---FORMAT PLTMTV
266 
267  IMPLICIT NONE
268 
269  CHARACTER(*) :: file_name
270  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
271  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
272  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: uu
273  INTEGER :: m, n1, n2, n3
274 
275  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
276 
277  WRITE (20, *) '$ DATA = CONTCURVE'
278  WRITE (20, *) '% contstyle = 2'
279  WRITE (20, *) '% nsteps = 50'
280  WRITE (20, *) '% meshplot = true'
281  WRITE (20, *)
282 
283  DO m = 1, SIZE(jj, 2)
284  n1 = jj(1, m)
285  n2 = jj(2, m)
286  n3 = jj(3, m)
287  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(m), n1
288  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(m), n2
289  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(m), n3
290  WRITE (20, *)
291  ENDDO
292 100 FORMAT(3(e12.5,3x),i5)
293  CLOSE(20)
294 
295  END SUBROUTINE plot_const_p1_label
296 
297  SUBROUTINE plot_pressure_p1_label(jj, rr, uu, file_name )
298  !---FORMAT PLTMTV
299 
300  IMPLICIT NONE
301 
302  CHARACTER(*) :: file_name
303  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
304  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
305  REAL(KIND=8), DIMENSION(:), INTENT(INOUT) :: uu ! rigid translation
306  INTEGER :: m, n1, n2, n3
307 
308  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
309 
310  WRITE (20, *) '$ DATA = CONTCURVE'
311  WRITE (20, *) '% contstyle = 2'
312  WRITE (20, *) '% nsteps = 50'
313  WRITE (20, *) '% meshplot = false'
314  WRITE (20, *)
315 
316  DO m = 1, SIZE(jj, 2)
317 
318  n1 = jj(1, m)
319  n2 = jj(2, m)
320  n3 = jj(3, m)
321  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
322  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
323  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
324  WRITE (20, *)
325  ENDDO
326 100 FORMAT(3(e12.5,3x),i5)
327  CLOSE(20)
328 
329  END SUBROUTINE plot_pressure_p1_label
330 
331  SUBROUTINE plot_p1_cont_label(jj, jjs, sides, list, rr, uu, file_name )
332  !---FORMAT PLOTMTV
333 
334  IMPLICIT NONE
335 
336  CHARACTER(*) :: file_name
337  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj, jjs
338  INTEGER, DIMENSION(:), INTENT(IN) :: sides, list
339  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
340  REAL(KIND=8), DIMENSION(:), INTENT(INOUT) :: uu ! rigid translation
341 
342  INTEGER :: m, n1, n2, n3, ms
343 
344  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
345 
346  WRITE (20, *) '$ DATA = CONTCURVE'
347  WRITE (20, *) '% contstyle = 2'
348  WRITE (20, *) '% nsteps = 50'
349  WRITE (20, *) '% meshplot = f'
350  WRITE (20, *)
351 
352  DO ms = 1, SIZE(sides)
353  IF(minval(abs(list-sides(ms))) /= 0 ) cycle
354  n1 = jjs(1,ms)
355  n2 = jjs(2,ms)
356  WRITE(20,110) '@line x1=',rr(1,n1), 'y1=',rr(2,n1), 'z1=',0., &
357  'x2=',rr(1,n2), 'y2=',rr(2,n2), 'z2=',0.
358  END DO
359 110 FORMAT(6(a,x,e12.5,3x))
360 
361  WRITE (20, *)
362  DO m = 1, SIZE(jj, 2)
363  n1 = jj(1, m)
364  n2 = jj(2, m)
365  n3 = jj(3, m)
366  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
367  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
368  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
369  WRITE (20, *)
370  ENDDO
371 100 FORMAT(3(e12.5,3x),i5)
372  CLOSE(20)
373 
374  END SUBROUTINE plot_p1_cont_label
375 
376  SUBROUTINE plot_p1_matiere_label(jj, neigh, i_d, rr, uu, file_name )
377  !---FORMAT PLOTMTV
378 
379  IMPLICIT NONE
380 
381  CHARACTER(*) :: file_name
382  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj, neigh
383  INTEGER, DIMENSION(:), INTENT(IN) :: i_d
384  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
385  REAL(KIND=8), DIMENSION(:), INTENT(INOUT) :: uu ! rigid translation
386 
387  INTEGER :: n, m, n1, n2, n3, nghm
388 
389  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
390 
391  WRITE (20, *) '$ DATA = CONTCURVE'
392  WRITE (20, *) '% contstyle = 2'
393  WRITE (20, *) '% nsteps = 50'
394  WRITE (20, *) '% meshplot = f'
395  WRITE (20, *)
396 
397  DO m = 1, SIZE(neigh,2)
398  DO n = 1, 3
399  nghm = neigh(n,m)
400 
401  IF (nghm == 0) THEN
402  CONTINUE
403  ELSE IF (i_d(nghm) /= i_d(m)) THEN
404  CONTINUE
405  ELSE
406  cycle
407  END IF
408 
409  n1 = jj(modulo(n,3) + 1,m)
410  n2 = jj(modulo(n+1,3) + 1,m)
411  WRITE(20,110) '@line x1=',rr(1,n1), 'y1=',rr(2,n1), 'z1=',0., &
412  'x2=',rr(1,n2), 'y2=',rr(2,n2), 'z2=',0.
413  END DO
414  END DO
415 110 FORMAT(6(a,x,e12.5,3x))
416 
417  WRITE (20, *)
418  DO m = 1, SIZE(jj, 2)
419  n1 = jj(1, m)
420  n2 = jj(2, m)
421  n3 = jj(3, m)
422  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
423  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
424  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
425  WRITE (20, *)
426  ENDDO
427 100 FORMAT(3(e12.5,3x),i5)
428  CLOSE(20)
429 
430  END SUBROUTINE plot_p1_matiere_label
431 
432  SUBROUTINE plot_pressure_label(jj, rr, uu, file_name )
433  !---FORMAT PLTMTV
434 
435  IMPLICIT NONE
436 
437  CHARACTER(*) :: file_name
438  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
439  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
440  REAL(KIND=8), DIMENSION(:), INTENT(INOUT) :: uu ! rigid translation
441 
442  INTEGER :: m, n1, n2, n3
443 
444  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
445 
446  WRITE (20, *) '$ DATA = CONTCURVE'
447  WRITE (20, *) '% contstyle = 2'
448  WRITE (20, *) '% nsteps = 50'
449  WRITE (20, *) '% meshplot = False'
450  WRITE (20, *)
451 
452  DO m = 1, SIZE(jj, 2)
453  n1 = jj(1, m)
454  n2 = jj(2, m)
455  n3 = jj(3, m)
456  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
457  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
458  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
459  WRITE (20, *)
460  ENDDO
461 100 FORMAT(3(e12.5,3x),i5)
462 
463  CLOSE(20)
464 
465  END SUBROUTINE plot_pressure_label
466 
467  SUBROUTINE plot_arrow_label(jj, rr, vv, file_name)
468  !---FORMAT PLTMTV
469 
470  IMPLICIT NONE
471 
472  CHARACTER(*) :: file_name
473  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
474  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr, vv
475 
476  INTEGER :: n
477  REAL(KIND=8) :: z_d, w_d
478 
479  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
480 
481  WRITE (20, 9000) SIZE(rr,2), SIZE(jj,1), SIZE(jj,2)
482 9000 FORMAT('#',3i10)
483  WRITE (20, *) '$ DATA = VECTOR'
484  WRITE (20, *) '% axisscale = FALSE'
485  WRITE (20, *) '% vscale = 1.'
486 
487  ! WRITE (20, *) '% contstyle = 2'
488  ! WRITE (20, *) '% meshplot = True'
489 
490  SELECT CASE(SIZE(rr,1))
491 
492  CASE(2)
493  z_d = 0
494  w_d = 0
495  IF (SIZE(vv,1)/=2) THEN
496  DO n = 1, SIZE(rr, 2)
497  WRITE (20,100) rr(1,n), rr(2,n), z_d, vv(n,1), vv(n,2), w_d
498  ENDDO
499  ELSE
500  DO n = 1, SIZE(rr, 2)
501  WRITE (20,100) rr(1,n), rr(2,n), z_d, vv(1,n), vv(2,n), w_d
502  ENDDO
503  END IF
504  CASE(3)
505  IF (SIZE(vv,1)/=3) THEN
506  DO n = 1, SIZE(rr, 2)
507  WRITE (20,100) rr(:,n), vv(:,n)
508  ENDDO
509  ELSE
510  DO n = 1, SIZE(rr, 2)
511  WRITE (20,100) rr(:,n), vv(:,n)
512  ENDDO
513  END IF
514  END SELECT
515 
516  CLOSE(20)
517 100 FORMAT(6(e12.5,3x))
518 
519  END SUBROUTINE plot_arrow_label
520 
521  SUBROUTINE plot_pressure_p2_label(jj, rr, uu, file_name)
522  !---FORMAT PLTMTV
523 
524  IMPLICIT NONE
525 
526  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
527  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
528  REAL(KIND=8), DIMENSION(:), INTENT(INOUT) :: uu ! rigid translation
529  CHARACTER(*) :: file_name
530 
531  INTEGER :: m, n1, n2, n3
532 
533  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
534 
535  WRITE (20, *) '$ DATA = CONTCURVE'
536  WRITE (20, *) '% contstyle = 2'
537  WRITE (20, *) '% nsteps = 50'
538  WRITE (20,*) '% meshplot = false'
539 
540  DO m = 1, SIZE(jj, 2)
541 
542  n1 = jj(1, m)
543  n2 = jj(6, m)
544  n3 = jj(5, m)
545  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
546  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
547  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
548  WRITE (20, 100)
549 
550 
551  n1 = jj(2, m)
552  n2 = jj(4, m)
553  n3 = jj(6, m)
554  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
555  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
556  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
557  WRITE (20, 100)
558 
559  n1 = jj(4, m)
560  n2 = jj(3, m)
561  n3 = jj(5, m)
562  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
563  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
564  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
565  WRITE (20, 100)
566 
567  n1 = jj(5, m)
568  n2 = jj(6, m)
569  n3 = jj(4, m)
570  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
571  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
572  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
573  WRITE (20, 100)
574  ENDDO
575 100 FORMAT(3(e12.5,3x),i5)
576 
577  CLOSE(20)
578 
579  END SUBROUTINE plot_pressure_p2_label
580 
581  !----------------------------------------------------------------
582 
583  SUBROUTINE plot_ensight_vecteur(u8, vit)
584  !--FORMAT ENSIGHT
585 
586  IMPLICIT NONE
587 
588  CHARACTER(LEN=*), INTENT(IN) :: vit
589  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: u8
590  CHARACTER(LEN=80) :: sketuve
591  INTEGER :: i, j
592 
593  OPEN(unit=33,file=vit,status='unknown',form='UNFORMATTED')
594 
595  sketuve = 'rien'
596  WRITE(33) sketuve
597 
598  WRITE(33) ((REAL(u8(i,j),KIND=4),i=1,size(u8,1)),j=1,size(u8,2))
599 
600  CLOSE(33)
601 
602  END SUBROUTINE plot_ensight_vecteur
603 
604  !----------------------------------------------------------------
605 
606  SUBROUTINE plot_ensight_scalaire(p8, pres)
607  !--FORMAT ENSIGHT
608 
609  IMPLICIT NONE
610 
611  CHARACTER(LEN=*), INTENT(IN) :: pres
612  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: p8
613  CHARACTER(LEN=80) :: sketuve
614  INTEGER :: i
615 
616  OPEN(unit=44,file=pres,status='unknown',form='UNFORMATTED')
617 
618  sketuve = 'rien'
619  WRITE(44) sketuve
620 
621  WRITE(44) (REAL(p8(i),KIND =4),i=1,size(p8))
622 
623  CLOSE(44)
624 
625  END SUBROUTINE plot_ensight_scalaire
626 
627  SUBROUTINE plot_vorticity (jj, rr, zz, i)
628 
629  IMPLICIT NONE
630 
631  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
632  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
633  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: zz
634  INTEGER, INTENT(IN) :: i
635 
636  INTEGER :: m, n1, n2, n3
637 
638  ! OPEN (UNIT = 20, FILE = 'vorticity', FORM = 'formatted', STATUS = 'unknown')
639 
640  WRITE (i, *) '$ DATA = CONTCURVE'
641  WRITE (i, *) '% contstyle = 1'
642  WRITE (i, *) '% nsteps = 40'
643 
644  ! WRITE (i, *) '% contours = "(-6.0 -5.0 -4.0 -3.0 -2.0 -1.0)"'
645  ! WRITE (i, *) '% contours = "(-0.5 0.0 0.5)"'
646  ! WRITE (i, *) '% contours = "( 1.0 2.0 3.0 4.0 5.0 6.0)"'
647 
648  WRITE (i, *) '% meshplot = False'
649 
650  DO m = 1, SIZE(jj, 2)
651  n1 = jj(1, m)
652  n2 = jj(2, m)
653  n3 = jj(3, m)
654  WRITE (i, *) rr(1,n1), rr(2,n1), zz(n1), n1
655  WRITE (i, *) rr(1,n2), rr(2,n2), zz(n2), n2
656  WRITE (i, *) rr(1,n3), rr(2,n3), zz(n3), n3
657  WRITE (i, *)
658  ENDDO
659 
660  CLOSE (unit = i)
661 
662  END SUBROUTINE plot_vorticity
663 
664  SUBROUTINE plot_stream (jj, rr, pp, i)
665 
666  IMPLICIT NONE
667 
668  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
669  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
670  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: pp
671  INTEGER, INTENT(IN) :: i
672 
673  INTEGER :: m, n1, n2, n3
674 
675  ! OPEN (UNIT = 20, FILE = 'stream', FORM = 'formatted', STATUS = 'unknown')
676 
677  WRITE (i, *) '$ DATA = CONTCURVE'
678  WRITE (i, *) '% contstyle = 1'
679  WRITE (i, *) '% nsteps = 30'
680 
681  ! WRITE (i, *) '% contours = "(-1.e-10 -1.e-7 -1.e-5 -1.e-4 )"'
682  ! WRITE (i, *) '% contours = "(-0.01 -0.03 -0.05 -0.07 -0.09 )"'
683  ! WRITE (i, *) '% contours = "(-0.100 -0.110 -0.115 -0.1175 )"'
684  ! WRITE (i, *) '% contours = "(+1.0e-8 +1.0e-7 +1.0e-6 +1.0e-5 )"'
685  ! WRITE (i, *) '% contours = "(+5.0e-5 +1.0e-4 +2.5e-4 +5.0e-4 )"'
686  ! WRITE (i, *) '% contours = "(+1.0e-3 +1.5e-3 +3.0e-3 )"'
687  WRITE (i, *) '% meshplot = False'
688 
689  DO m = 1, SIZE(jj, 2)
690  n1 = jj(1, m)
691  n2 = jj(2, m)
692  n3 = jj(3, m)
693  WRITE (i, *) rr(1,n1), rr(2,n1), pp(n1), n1
694  WRITE (i, *) rr(1,n2), rr(2,n2), pp(n2), n2
695  WRITE (i, *) rr(1,n3), rr(2,n3), pp(n3), n3
696  WRITE (i, *)
697  ENDDO
698 
699  CLOSE (unit = i)
700 
701  END SUBROUTINE plot_stream
702 
703  SUBROUTINE plot_scalar_field(jj, rr, uu, file_name)
704  !---FORMAT PLTMTV
705 
706  IMPLICIT NONE
707 
708  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
709  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
710  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: uu ! rigid translation
711  CHARACTER(*) :: file_name
712 
713  INTEGER :: m, n1, n2, n3, unit_w=47
714 
715 100 FORMAT(3(f15.8,3x),i5)
716  OPEN (unit=unit_w, file=file_name, form='formatted', status='unknown')
717 
718  WRITE (unit_w, '(A)') '$ DATA = CONTCURVE'
719  WRITE (unit_w, '(A)') '% contstyle = 2'
720  WRITE (unit_w, '(A)') '% nsteps = 50'
721  WRITE (unit_w, '(A)') '% meshplot = false'
722  WRITE (unit_w, '(A)')
723 
724  IF (SIZE(jj,1)==3) THEN
725  DO m = 1, SIZE(jj, 2)
726  n1 = jj(1, m)
727  n2 = jj(2, m)
728  n3 = jj(3, m)
729  WRITE (unit_w, 100) rr(1,n1), rr(2,n1), uu(n1), n1
730  WRITE (unit_w, 100) rr(1,n2), rr(2,n2), uu(n2), n2
731  WRITE (unit_w, 100) rr(1,n3), rr(2,n3), uu(n3), n3
732  WRITE (unit_w, 100)
733  ENDDO
734  ELSE IF (SIZE(jj,1)==6) THEN
735  DO m = 1, SIZE(jj, 2)
736  n1 = jj(1, m)
737  n2 = jj(6, m)
738  n3 = jj(5, m)
739  WRITE (unit_w, 100) rr(1,n1), rr(2,n1), uu(n1), n1
740  WRITE (unit_w, 100) rr(1,n2), rr(2,n2), uu(n2), n2
741  WRITE (unit_w, 100) rr(1,n3), rr(2,n3), uu(n3), n3
742  WRITE (unit_w, 100)
743 
744  n1 = jj(2, m)
745  n2 = jj(4, m)
746  n3 = jj(6, m)
747  WRITE (unit_w, 100) rr(1,n1), rr(2,n1), uu(n1), n1
748  WRITE (unit_w, 100) rr(1,n2), rr(2,n2), uu(n2), n2
749  WRITE (unit_w, 100) rr(1,n3), rr(2,n3), uu(n3), n3
750  WRITE (unit_w, 100)
751 
752  n1 = jj(4, m)
753  n2 = jj(3, m)
754  n3 = jj(5, m)
755  WRITE (unit_w, 100) rr(1,n1), rr(2,n1), uu(n1), n1
756  WRITE (unit_w, 100) rr(1,n2), rr(2,n2), uu(n2), n2
757  WRITE (unit_w, 100) rr(1,n3), rr(2,n3), uu(n3), n3
758  WRITE (unit_w, 100)
759 
760  n1 = jj(5, m)
761  n2 = jj(6, m)
762  n3 = jj(4, m)
763  WRITE (unit_w, 100) rr(1,n1), rr(2,n1), uu(n1), n1
764  WRITE (unit_w, 100) rr(1,n2), rr(2,n2), uu(n2), n2
765  WRITE (unit_w, 100) rr(1,n3), rr(2,n3), uu(n3), n3
766  WRITE (unit_w, 100)
767  ENDDO
768 
769  ELSE
770  WRITE(*,*) ' Problem in plot_scalar_field '
771  stop
772  END IF
773 
774  CLOSE(unit=unit_w)
775 
776  END SUBROUTINE plot_scalar_field
777 
778  SUBROUTINE plot_two_scalar_field(jj, rr, uu, jj2, rr2, uu2, file_name)
779  !---FORMAT PLTMTV
780 
781  IMPLICIT NONE
782 
783  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj, jj2
784  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr, rr2
785  REAL(KIND=8), DIMENSION(:) :: uu, uu2
786  CHARACTER(*) :: file_name
787 
788  INTEGER :: m, n1, n2, n3
789 
790  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
791 
792  WRITE (20, *) '$ DATA = CONTCURVE'
793  WRITE (20, *) '% contstyle = 2'
794  WRITE (20, *) '% nsteps = 50'
795  WRITE (20, *) '% meshplot = false'
796  WRITE (20, *)
797 
798  IF (SIZE(jj,1)==3) THEN
799  DO m = 1, SIZE(jj, 2)
800 
801  n1 = jj(1, m)
802  n2 = jj(2, m)
803  n3 = jj(3, m)
804  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
805  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
806  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
807  WRITE (20, *)
808  ENDDO
809  DO m = 1, SIZE(jj2, 2)
810 
811  n1 = jj2(1, m)
812  n2 = jj2(2, m)
813  n3 = jj2(3, m)
814  WRITE (20, 100) rr2(1,n1), rr2(2,n1), uu2(n1), n1
815  WRITE (20, 100) rr2(1,n2), rr2(2,n2), uu2(n2), n2
816  WRITE (20, 100) rr2(1,n3), rr2(2,n3), uu2(n3), n3
817  WRITE (20, *)
818  ENDDO
819 
820 
821  ELSE IF (SIZE(jj,1)==6) THEN
822 
823  DO m = 1, SIZE(jj, 2)
824 
825  n1 = jj(1, m)
826  n2 = jj(6, m)
827  n3 = jj(5, m)
828  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
829  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
830  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
831  WRITE (20, 100)
832 
833 
834  n1 = jj(2, m)
835  n2 = jj(4, m)
836  n3 = jj(6, m)
837  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
838  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
839  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
840  WRITE (20, 100)
841 
842  n1 = jj(4, m)
843  n2 = jj(3, m)
844  n3 = jj(5, m)
845  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
846  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
847  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
848  WRITE (20, 100)
849 
850  n1 = jj(5, m)
851  n2 = jj(6, m)
852  n3 = jj(4, m)
853  WRITE (20, 100) rr(1,n1), rr(2,n1), uu(n1), n1
854  WRITE (20, 100) rr(1,n2), rr(2,n2), uu(n2), n2
855  WRITE (20, 100) rr(1,n3), rr(2,n3), uu(n3), n3
856  WRITE (20, 100)
857  ENDDO
858  DO m = 1, SIZE(jj2, 2)
859 
860  n1 = jj2(1, m)
861  n2 = jj2(6, m)
862  n3 = jj2(5, m)
863  WRITE (20, 100) rr2(1,n1), rr2(2,n1), uu2(n1), n1
864  WRITE (20, 100) rr2(1,n2), rr2(2,n2), uu2(n2), n2
865  WRITE (20, 100) rr2(1,n3), rr2(2,n3), uu2(n3), n3
866  WRITE (20, 100)
867 
868 
869  n1 = jj2(2, m)
870  n2 = jj2(4, m)
871  n3 = jj2(6, m)
872  WRITE (20, 100) rr2(1,n1), rr2(2,n1), uu2(n1), n1
873  WRITE (20, 100) rr2(1,n2), rr2(2,n2), uu2(n2), n2
874  WRITE (20, 100) rr2(1,n3), rr2(2,n3), uu2(n3), n3
875  WRITE (20, 100)
876 
877  n1 = jj2(4, m)
878  n2 = jj2(3, m)
879  n3 = jj2(5, m)
880  WRITE (20, 100) rr2(1,n1), rr2(2,n1), uu2(n1), n1
881  WRITE (20, 100) rr2(1,n2), rr2(2,n2), uu2(n2), n2
882  WRITE (20, 100) rr2(1,n3), rr2(2,n3), uu2(n3), n3
883  WRITE (20, 100)
884 
885  n1 = jj2(5, m)
886  n2 = jj2(6, m)
887  n3 = jj2(4, m)
888  WRITE (20, 100) rr2(1,n1), rr2(2,n1), uu2(n1), n1
889  WRITE (20, 100) rr2(1,n2), rr2(2,n2), uu2(n2), n2
890  WRITE (20, 100) rr2(1,n3), rr2(2,n3), uu2(n3), n3
891  WRITE (20, 100)
892  ENDDO
893 
894 
895  ELSE
896  WRITE(*,*) ' Problem in plot_two_scalar_field '
897  stop
898  END IF
899 
900 100 FORMAT(3(e12.5,3x),i5)
901 
902  CLOSE(20)
903 
904  END SUBROUTINE plot_two_scalar_field
905 
906  SUBROUTINE plot_scalar_field_domain(jj, rr, id, index_dom, file_name)
907  !---FORMAT PLTMTV
908 
909  IMPLICIT NONE
910 
911  INTEGER, DIMENSION(:,:), INTENT(IN) :: jj
912  REAL(KIND=8), DIMENSION(:,:), INTENT(IN) :: rr
913  INTEGER, DIMENSION(:), INTENT(IN) :: id
914  INTEGER, INTENT(IN) :: index_dom
915  CHARACTER(*) :: file_name
916 
917  INTEGER :: m, n1, n2, n3
918 
919  OPEN (unit=20, file=file_name, form='formatted', status='unknown')
920 
921  WRITE (20, *) '$ DATA = CONTCURVE'
922  WRITE (20, *) '% contstyle = 2'
923  WRITE (20, *) '% nsteps = 50'
924  WRITE (20, *) '% meshplot = true'
925  WRITE (20, *)
926 
927  IF (SIZE(jj,1)==3) THEN
928  DO m = 1, SIZE(jj, 2)
929  IF (id(m)/=index_dom) cycle
930  n1 = jj(1, m)
931  n2 = jj(2, m)
932  n3 = jj(3, m)
933  WRITE (20, 100) rr(1,n1), rr(2,n1), id(m), n1
934  WRITE (20, 100) rr(1,n2), rr(2,n2), id(m), n2
935  WRITE (20, 100) rr(1,n3), rr(2,n3), id(m), n3
936  WRITE (20, *)
937  ENDDO
938 
939  ELSE IF (SIZE(jj,1)==6) THEN
940 
941  DO m = 1, SIZE(jj, 2)
942  IF (id(m)/=index_dom) cycle
943  n1 = jj(1, m)
944  n2 = jj(6, m)
945  n3 = jj(5, m)
946  WRITE (20, 100) rr(1,n1), rr(2,n1), id(m), n1
947  WRITE (20, 100) rr(1,n2), rr(2,n2), id(m), n2
948  WRITE (20, 100) rr(1,n3), rr(2,n3), id(m), n3
949  WRITE (20, 100)
950 
951  n1 = jj(2, m)
952  n2 = jj(4, m)
953  n3 = jj(6, m)
954  WRITE (20, 100) rr(1,n1), rr(2,n1), id(m), n1
955  WRITE (20, 100) rr(1,n2), rr(2,n2), id(m), n2
956  WRITE (20, 100) rr(1,n3), rr(2,n3), id(m), n3
957  WRITE (20, 100)
958 
959  n1 = jj(4, m)
960  n2 = jj(3, m)
961  n3 = jj(5, m)
962  WRITE (20, 100) rr(1,n1), rr(2,n1), id(m), n1
963  WRITE (20, 100) rr(1,n2), rr(2,n2), id(m), n2
964  WRITE (20, 100) rr(1,n3), rr(2,n3), id(m), n3
965  WRITE (20, 100)
966 
967  n1 = jj(5, m)
968  n2 = jj(6, m)
969  n3 = jj(4, m)
970  WRITE (20, 100) rr(1,n1), rr(2,n1), id(m), n1
971  WRITE (20, 100) rr(1,n2), rr(2,n2), id(m), n2
972  WRITE (20, 100) rr(1,n3), rr(2,n3), id(m), n3
973  WRITE (20, 100)
974  ENDDO
975 
976  ELSE
977  WRITE(*,*) ' Problem in plot_scalar_field '
978  stop
979  END IF
980 
981 100 FORMAT(2(e12.5,3x),i5,3x,i5)
982 
983  CLOSE(20)
984 
985  END SUBROUTINE plot_scalar_field_domain
986 
987  SUBROUTINE vtk_p1_2d(mesh, champ, unit_file, file_name)
988 
989  USE def_type_mesh
990 
991  IMPLICIT NONE
992 
993  TYPE(mesh_type), TARGET :: mesh
994  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: champ
995  INTEGER, INTENT(IN) :: unit_file
996  CHARACTER(*), INTENT(IN) :: file_name
997  INTEGER :: i
998 
999  write(*,*) 'ecriture header'
1000  OPEN (unit=unit_file,file=file_name,form = 'formatted',&
1001  status = 'unknown')
1002 
1003  write(unit_file,'(A)') '# vtk DataFile Version 3.0'
1004  write(unit_file,'(A)') 'vtk '//file_name//''
1005  write(unit_file,'(A)')'ASCII'
1006  write(unit_file,'(A)')'DATASET UNSTRUCTURED_GRID'
1007  write(unit_file,'(A,I7,A)')'POINTS ', mesh%np, ' float'
1008  write(*,*) 'points ...'
1009  DO i=1, mesh%np
1010  write(unit_file,'(2(e14.7,2x),A)') mesh%rr(1,i), &
1011  mesh%rr(2,i), ' 0.0 '
1012  ENDDO
1013  write(*,*) 'cells ...'
1014  write(unit_file,'(A,I7,I8)') 'CELLS ', mesh%me, mesh%me*4
1015  DO i=1, mesh%me
1016  write(unit_file,'(A,3(I8,1x))') '3 ', mesh%jj(1,i)-1, &
1017  mesh%jj(2,i)-1, mesh%jj(3,i)-1
1018  ENDDO
1019  write(unit_file,'(A,I7)') 'CELL_TYPES ', mesh%me
1020  DO i=1, mesh%me
1021  write(unit_file,'(A)') '5'
1022  ENDDO
1023 
1024  write(*,*) 'data ...'
1025  write(unit_file,'(A,I7)') 'POINT_DATA ',mesh%np
1026  write(unit_file,'(A)') 'SCALARS scalars float 1'
1027  write(unit_file,'(A)') 'LOOKUP_TABLE default'
1028  DO i=1, mesh%np
1029  write(unit_file,'(e14.7,2x)') champ(i)
1030  ENDDO
1031  CLOSE(unit_file)
1032  END SUBROUTINE vtk_p1_2d
1033 
1034  SUBROUTINE vtk_2d(mesh, champ, unit_file, file_name)
1035 
1036  USE def_type_mesh
1037 
1038  IMPLICIT NONE
1039 
1040  TYPE(mesh_type), TARGET :: mesh
1041  REAL(KIND=8), DIMENSION(:), INTENT(IN) :: champ
1042  INTEGER, INTENT(IN) :: unit_file
1043  CHARACTER(*), INTENT(IN) :: file_name
1044  INTEGER :: i, vtk_cell, n, m
1045 
1046 
1047  write(*,*) 'ecriture header'
1048  OPEN (unit=unit_file,file=file_name,form = 'formatted',&
1049  status = 'unknown')
1050 
1051  write(unit_file,'(A)') '# vtk DataFile Version 3.0'
1052  write(unit_file,'(A)') 'vtk '//file_name//''
1053  write(unit_file,'(A)')'ASCII'
1054  write(unit_file,'(A)')'DATASET UNSTRUCTURED_GRID'
1055  write(unit_file,'(A,I7,A)')'POINTS ', mesh%np, ' float'
1056  write(*,*) 'points ...'
1057  DO i=1, mesh%np
1058  write(unit_file,'(2(e14.7,2x),A)') mesh%rr(1,i), &
1059  mesh%rr(2,i), ' 0.0 '
1060  ENDDO
1061  write(*,*) 'cells ...'
1062 
1063  write(unit_file,'(A,I7,I8)') 'CELLS ', mesh%me, mesh%me*(mesh%gauss%n_w+1)
1064  IF (mesh%gauss%n_w==3) THEN
1065  vtk_cell = 5
1066  DO m=1, mesh%me
1067  write(unit_file,'(I2,6(I8,1x))') mesh%gauss%n_w, (mesh%jj(n,m)-1, n=1,mesh%gauss%n_w)
1068  ENDDO
1069  ELSE IF (mesh%gauss%n_w==6) THEN
1070  vtk_cell = 22
1071  DO m=1, mesh%me
1072  write(unit_file,'(I2,6(I8,1x))') mesh%gauss%n_w, &
1073  mesh%jj(1,m)-1, mesh%jj(2,m)-1, mesh%jj(3,m)-1, mesh%jj(6,m)-1, mesh%jj(4,m)-1, mesh%jj(5,m)-1
1074  ENDDO
1075  ELSE
1076  WRITE(*,*) ' BUG in vtk_2d '
1077  END IF
1078  write(unit_file,'(A,I7)') 'CELL_TYPES ', mesh%me
1079  DO i=1, mesh%me
1080  write(unit_file,'(i2)') vtk_cell
1081  ENDDO
1082 
1083  write(*,*) 'data ...'
1084  write(unit_file,'(A,I7)') 'POINT_DATA ',mesh%np
1085  write(unit_file,'(A)') 'SCALARS scalars float 1'
1086  write(unit_file,'(A)') 'LOOKUP_TABLE default'
1087  DO i=1, mesh%np
1088  write(unit_file,'(e14.7,2x)') champ(i)
1089  ENDDO
1090  CLOSE(unit_file)
1091 
1092 
1093  write(*,*) 'ecriture header'
1094  OPEN (unit=unit_file,file=file_name,form = 'formatted',&
1095  status = 'unknown')
1096 
1097 !!$ write(unit_file,'(A)') '# vtk DataFile Version 3.0'
1098 !!$ write(unit_file,'(A)') 'vtk '//file_name//''
1099 !!$ write(unit_file,'(A)')'ASCII'
1100 !!$ write(unit_file,'(A)')'DATASET UNSTRUCTURED_GRID'
1101 !!$ write(unit_file,'(A,I7,A)')'POINTS ', mesh%np, ' float'
1102 !!$ write(*,*) 'points ...'
1103 !!$ DO i = 1, mesh%np
1104 !!$ write(unit_file,'(2(e14.7,2x),A)') mesh%rr(1,i), &
1105 !!$ mesh%rr(2,i), ' 0.0 '
1106 !!$ ENDDO
1107 !!$ write(*,*) 'cells ...'
1108 !!$ write(unit_file,'(A,I7,I8)') 'CELLS ', mesh%me, mesh%me*(mesh%gauss%n_w+1)
1109 !!$ DO m = 1, mesh%me
1110 !!$ write(unit_file,'(A,3(I8,1x))') mesh%gauss%n_w, (mesh%jj(n,m)-1, n=1,mesh%gauss%n_w)
1111 !!$ ENDDO
1112 !!$ IF (mesh%gauss%n_w==3) THEN
1113 !!$ vtK_cell = 5
1114 !!$ ELSE IF (mesh%gauss%n_w==6) THEN
1115 !!$ vtK_cell = 22
1116 !!$ ELSE
1117 !!$ WRITE(*,*) ' BUG in vtk_2d '
1118 !!$ END IF
1119 !!$ write(unit_file,'(A,I7)') 'CELL_TYPES ', mesh%me
1120 !!$ DO m = 1, mesh%me
1121 !!$ write(unit_file,'(I2)') vtk_cell
1122 !!$ ENDDO
1123 !!$
1124 !!$ write(*,*) 'data ...'
1125 !!$ write(unit_file,'(A,I7)') 'POINT_DATA ',mesh%np
1126 !!$ write(unit_file,'(A)') 'SCALARS scalars float 1'
1127 !!$ write(unit_file,'(A)') 'LOOKUP_TABLE default'
1128 !!$ DO i=1, mesh%np
1129 !!$ write(unit_file,'(e14.7,2x)') champ(i)
1130 !!$ ENDDO
1131 !!$ CLOSE(unit_file)
1132  END SUBROUTINE vtk_2d
1133 
1134  SUBROUTINE trace_profile(mesh, v, it, freq_plot, list_mode, nom_champ, num_dom)
1135 
1136  USE chaine_caractere
1137  USE def_type_mesh
1138 
1139  IMPLICIT NONE
1140  TYPE(mesh_type) :: mesh
1141  REAL(KIND=8) , DIMENSION(:,:,:), INTENT(INOUT):: v
1142  INTEGER, INTENT(IN) :: it, freq_plot
1143  INTEGER, DIMENSION(:), INTENT(IN) :: list_mode
1144  CHARACTER(len=3), INTENT(IN) :: nom_champ
1145  INTEGER, OPTIONAL, INTENT(IN) :: num_dom
1146 
1147  INTEGER :: l, lblank, rang_s
1148  CHARACTER(len=3) :: tit, tmode, tit_s
1149  CHARACTER(len=1) :: tcomp
1150  INTEGER :: i, k
1151 
1152  IF (present(num_dom)) THEN
1153  rang_s = num_dom
1154  ELSE
1155  rang_s = 0
1156  END IF
1157 
1158  WRITE(tit_s,'(i3)') rang_s
1159  lblank = eval_blank(3,tit_s)
1160  DO l = 1, lblank - 1
1161  tit_s(l:l) = '0'
1162  END DO
1163 
1164  WRITE(tit,'(i3)') it/freq_plot
1165  lblank = eval_blank(3,tit)
1166  DO l = 1, lblank - 1
1167  tit(l:l) = '0'
1168  END DO
1169 
1170  DO i= 1, SIZE(v,3)
1171 
1172  WRITE(tmode,'(i3)') list_mode(i)
1173  lblank = eval_blank(3,tmode)
1174  DO l = 1, lblank - 1
1175  tmode(l:l) = '0'
1176  END DO
1177 
1178  DO k= 1, size(v,2)
1179  WRITE(tcomp,'(i1)') k
1180  CALL plot_scalar_field(mesh%jj, mesh%rr, &
1181  v(:,k,i) , nom_champ//tcomp//'_m='//tmode//'_'//tit_s//'_'//tit//'.plt')
1182  ENDDO
1183 
1184  END DO
1185 
1186  END SUBROUTINE trace_profile
1187 
1188 END MODULE sub_plot
integer function eval_blank(len_str, string)
subroutine plot_p1_cont_label(jj, jjs, sides, list, rr, uu, file_name)
Definition: sub_plot.f90:331
subroutine plot_pressure(jj, rr, uu)
Definition: sub_plot.f90:235
subroutine plot_scalar_field_domain(jj, rr, id, index_dom, file_name)
Definition: sub_plot.f90:906
subroutine plot_arrow(jj, rr, vv)
Definition: sub_plot.f90:72
subroutine plot_vorticity(jj, rr, zz, i)
Definition: sub_plot.f90:627
subroutine trace_profile(mesh, v, it, freq_plot, list_mode, nom_champ, num_dom)
Definition: sub_plot.f90:1134
subroutine plot_ensight_scalaire(p8, pres)
Definition: sub_plot.f90:606
subroutine plot_p1_matiere_label(jj, neigh, i_d, rr, uu, file_name)
Definition: sub_plot.f90:376
subroutine plot_pressure_p1_label(jj, rr, uu, file_name)
Definition: sub_plot.f90:297
subroutine vtk_p1_2d(mesh, champ, unit_file, file_name)
Definition: sub_plot.f90:987
subroutine plot_arrow_label(jj, rr, vv, file_name)
Definition: sub_plot.f90:467
subroutine plot_scalar_field(jj, rr, uu, file_name)
Definition: sub_plot.f90:703
subroutine plot_stream(jj, rr, pp, i)
Definition: sub_plot.f90:664
subroutine plot_const_p1_label(jj, rr, uu, file_name)
Definition: sub_plot.f90:264
subroutine plot_ensight_vecteur(u8, vit)
Definition: sub_plot.f90:583
subroutine plot_two_scalar_field(jj, rr, uu, jj2, rr2, uu2, file_name)
Definition: sub_plot.f90:778
subroutine plot_loc_rel_var(jj, rr, vo, vv)
Definition: sub_plot.f90:182
subroutine plot_pressure_label(jj, rr, uu, file_name)
Definition: sub_plot.f90:432
subroutine plot_vit_2d(jj, rr, uu)
Definition: sub_plot.f90:8
subroutine vtk_2d(mesh, champ, unit_file, file_name)
Definition: sub_plot.f90:1034
subroutine plot_pressure_p2_label(jj, rr, uu, file_name)
Definition: sub_plot.f90:521
subroutine plot_pressure_2d(jj, rr, uu)
Definition: sub_plot.f90:116
section doc_intro_frame_work_num_app Numerical approximation subsection doc_intro_fram_work_num_app_Fourier_FEM Fourier Finite element representation The SFEMaNS code uses a hybrid Fourier Finite element formulation The Fourier decomposition allows to approximate the problem’s solutions for each Fourier mode modulo nonlinear terms that are made explicit The variables are then approximated on a meridian section of the domain with a finite element method The numerical approximation of a function f $f f is written in the following generic form
Definition: doc_intro.h:193