]> SALOME platform Git repositories - modules/paravis.git/blob - src/PV_SWIG/presentations.py
Salome HOME
Merge branch 'V8_3_BR' into gdd/python3_dev
[modules/paravis.git] / src / PV_SWIG / presentations.py
1 # Copyright (C) 2010-2016  CEA/DEN, EDF R&D
2 #
3 # This library is free software; you can redistribute it and/or
4 # modify it under the terms of the GNU Lesser General Public
5 # License as published by the Free Software Foundation; either
6 # version 2.1 of the License, or (at your option) any later version.
7 #
8 # This library is distributed in the hope that it will be useful,
9 # but WITHOUT ANY WARRANTY; without even the implied warranty of
10 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 # Lesser General Public License for more details.
12 #
13 # You should have received a copy of the GNU Lesser General Public
14 # License along with this library; if not, write to the Free Software
15 # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 #
17 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 #
19
20 """
21 This module is intended to provide Python API for building presentations
22 typical for Post-Pro module (Scalar Map, Deformed Shape, Vectors, etc.)
23 """
24
25
26
27 ##from __future__ import print_function
28
29 import os
30 import re
31 import warnings
32 from math import sqrt, sin, cos, radians
33
34 # Do not use pv as a short name.
35 # It is a name of function from numpy and may be redefined implicitly by 'from numpy import *' call.
36 # import pvsimple as pv
37 import pvsimple as pvs
38 #try:
39 #    # TODO(MZN): to be removed (issue with Point Sprite texture)
40 #    #import paravisSM as sm
41 #except:
42 #    import paraview.simple as pvs
43 #    import paraview.servermanager as sm
44
45
46 # Constants
47 EPS = 1E-3
48 FLT_MIN = 1E-37
49 VTK_LARGE_FLOAT = 1E+38
50 GAP_COEFFICIENT = 0.0001
51
52
53 # Globals
54 _current_bar = None
55 _med_field_sep = '@@][@@'
56
57
58 # Enumerations
59 class PrsTypeEnum:
60     """
61     Post-Pro presentation types.
62     """
63     MESH = 0
64     SCALARMAP = 1
65     ISOSURFACES = 2
66     CUTPLANES = 3
67     CUTLINES = 4
68     DEFORMEDSHAPE = 5
69     DEFORMEDSHAPESCALARMAP = 6
70     VECTORS = 7
71     PLOT3D = 8
72     STREAMLINES = 9
73     GAUSSPOINTS = 10
74
75     _type2name = {MESH: 'Mesh',
76                   SCALARMAP: 'Scalar Map',
77                   ISOSURFACES: 'Iso Surfaces',
78                   CUTPLANES: 'Cut Planes',
79                   CUTLINES: 'Cut Lines',
80                   DEFORMEDSHAPE: 'Deformed Shape',
81                   DEFORMEDSHAPESCALARMAP: 'Deformed Shape And Scalar Map',
82                   VECTORS: 'Vectors',
83                   PLOT3D: 'Plot3D',
84                   STREAMLINES: 'Stream Lines',
85                   GAUSSPOINTS: 'Gauss Points'}
86
87     @classmethod
88     def get_name(cls, type):
89         """Return presentaion name by its type."""
90         return cls._type2name[type]
91
92
93 class EntityType:
94     """
95     Entity types.
96     """
97     NODE = 0
98     CELL = 1
99
100     _type2name = {NODE: 'P1',
101                   CELL: 'P0'}
102
103     _name2type = {'P1': NODE,
104                   'P0': CELL}
105
106     _type2pvtype = {NODE: 'POINT_DATA',
107                     CELL: 'CELL_DATA'}
108
109     @classmethod
110     def get_name(cls, type):
111         """Return entity name (used in full group names) by its type."""
112         return cls._type2name[type]
113
114     @classmethod
115     def get_type(cls, name):
116         """Return entity type by its name (used in full group names)."""
117         return cls._name2type[name]
118
119     @classmethod
120     def get_pvtype(cls, type):
121         """Return entity type from ['CELL_DATA', 'POINT_DATA']"""
122         return cls._type2pvtype[type]
123
124
125 class Orientation:
126     """Orientation types.
127
128     Defines a set of plane orientation possibilities:
129       AUTO: plane orientation should be calculated.
130       XY: plane formed by X and Y axis.
131       YZ: plane formed by Y and Z axis.
132       ZX: plane formed by Z and X axis
133
134     """
135     AUTO = 0
136     XY = 1
137     YZ = 2
138     ZX = 3
139
140
141 class GlyphPos:
142     """Glyph positions.
143
144     Set of elements defining the position of the vector head:
145       CENTER: in the center of the vector
146       TAIL: in the tail of the vector
147       HEAD: in the head of the vector
148
149     """
150     CENTER = 0
151     TAIL = 1
152     HEAD = 2
153
154
155 class GaussType:
156     """
157     Gauss Points primitive types.
158     """
159     SPRITE = 0
160     POINT = 1
161     SPHERE = 2
162
163     _type2mode = {SPRITE: 'Texture',
164                   POINT: 'SimplePoint',
165                   SPHERE: 'Sphere (Texture)'}
166
167     @classmethod
168     def get_mode(cls, type):
169         """Return paraview point sprite mode by the primitive type."""
170         return cls._type2mode[type]
171
172
173 # Auxiliary functions
174
175 def get_field_mesh_name(full_field_name):
176     """Return mesh name of the field by its full name."""
177     aList = full_field_name.split('/')
178     if len(aList) >= 2 :
179         field_name = full_field_name.split('/')[1]
180         return field_name
181
182
183 def get_field_entity(full_field_name):
184     """Return entity type of the field by its full name."""
185     aList = full_field_name.split(_med_field_sep)
186     if len(aList) == 2 :
187         entity_name = full_field_name.split(_med_field_sep)[-1]
188         entity = EntityType.get_type(entity_name)
189         return entity
190
191
192 def get_field_short_name(full_field_name):
193     """Return short name of the field by its full name."""
194     aList = full_field_name.split('/')
195     if len(aList) == 4 :
196         short_name_with_type = full_field_name.split('/')[-1]
197         short_name = short_name_with_type.split(_med_field_sep)[0]
198         return short_name
199
200
201 def find_mesh_full_name(proxy, short_mesh_name):
202     """Return full mesh path by short mesh name, if found"""
203     proxy.UpdatePipeline()
204     all_mesh_names = get_mesh_full_names(proxy)
205     for name in all_mesh_names:
206         if short_mesh_name == get_field_short_name(name):
207             return name
208
209
210 def process_prs_for_test(prs, view, picture_name, show_bar=True):
211     """Show presentation and record snapshot image.
212
213     Arguments:
214       prs: the presentation to show
215       view: the render view
216       picture_name: the full name of the graphics file to save
217       show_bar: to show scalar bar or not
218
219     """
220     # Show the presentation only
221     display_only(prs, view)
222
223     # Show scalar bar
224     global _current_bar
225     if show_bar and _current_bar:
226         _current_bar.Visibility = 1
227
228     # Reset the view
229     reset_view(view)
230
231     # Create a directory for screenshot if necessary
232     file_name = re.sub("\s+", "_", picture_name)
233     pic_dir = os.path.dirname(picture_name)
234     if not os.path.exists(pic_dir):
235         os.makedirs(pic_dir)
236
237     # Save picture
238     print("Write image:", file_name)
239     pvs.WriteImage(file_name, view=view, Magnification=1)
240
241
242 def reset_view(view=None):
243     """Reset the view.
244
245     Set predefined (taken from Post-Pro) camera settings.
246     If the view is not passed, the active view is used.
247
248     """
249     if not view:
250         view = pvs.GetRenderView()
251
252     # Camera preferences
253     view.CameraFocalPoint = [0.0, 0.0, 0.0]
254     view.CameraViewUp = [0.0, 0.0, 1.0]
255     view.CameraPosition = [738.946, -738.946, 738.946]
256
257     # Turn on the headligth
258     view.LightSwitch = 1
259     view.LightIntensity = 0.5
260
261     # Use parallel projection
262     view.CameraParallelProjection = 1
263
264     view.ResetCamera()
265     pvs.Render(view=view)
266
267
268 def hide_all(view, to_remove=False):
269     """Hide all representations in the view."""
270     if not view:
271         view = pvs.GetRenderView()
272
273     rep_list = view.Representations
274     for rep in rep_list:
275         if hasattr(rep, 'Visibility') and rep.Visibility != 0:
276             rep.Visibility = 0
277         if to_remove:
278             view.Representations.remove(rep)
279     pvs.Render(view=view)
280
281
282 def display_only(prs, view=None):
283     """Display only the given presentation in the view."""
284     if not view:
285         view = pvs.GetRenderView()
286
287     rep_list = view.Representations
288     for rep in rep_list:
289         if hasattr(rep, 'Visibility'):
290             rep.Visibility = (rep == prs)
291     pvs.Render(view=view)
292
293
294 def set_visible_lines(xy_prs, lines):
295     """Set visible only the given lines for XYChartRepresentation."""
296     sv = xy_prs.GetProperty("SeriesVisibility").GetData()
297     visible = '0'
298
299     for i in range(0, len(sv)):
300         if i % 2 == 0:
301             line_name = sv[i]
302             if line_name in lines:
303                 visible = '1'
304             else:
305                 visible = '0'
306         else:
307             sv[i] = visible
308
309     xy_prs.SeriesVisibility = sv
310
311
312 def check_vector_mode(vector_mode, nb_components):
313     """Check vector mode.
314
315     Check if vector mode is correct for the data array with the
316     given number of components.
317
318     Arguments:
319       vector_mode: 'Magnitude', 'X', 'Y' or 'Z'
320       nb_components: number of component in the data array
321
322     Raises:
323       ValueError: in case of the vector mode is unexistent
324       or nonapplicable.
325
326     """
327     if vector_mode not in ('Magnitude', 'X', 'Y', 'Z'):
328         raise ValueError("Unexistent vector mode: " + vector_mode)
329
330     if ((nb_components == 1 and (vector_mode == 'Y' or vector_mode == 'Z')) or
331         (nb_components == 2 and  vector_mode == 'Z')):
332         raise ValueError("Incorrect vector mode " + vector_mode + " for " +
333                          nb_components + "-component field")
334
335
336 def get_vector_component(vector_mode):
337     """Get vector component as ineger.
338
339     Translate vector component notation from string
340     to integer:
341       'Magnitude': -1
342       'X': 0
343       'Y': 1
344       'Z': 2
345
346     """
347     vcomponent = -1
348
349     if vector_mode == 'X':
350         vcomponent = 0
351     elif vector_mode == 'Y':
352         vcomponent = 1
353     elif vector_mode == 'Z':
354         vcomponent = 2
355
356     return vcomponent
357
358
359 def get_data_range(proxy, entity, field_name, vector_mode='Magnitude',
360                    cut_off=False):
361     """Get data range for the field.
362
363     Arguments:
364       proxy: the pipeline object, containig data array for the field
365       entity: the field entity
366       field_name: the field name
367       vector_mode: the vector mode ('Magnitude', 'X', 'Y' or 'Z')
368
369     Returns:
370       Data range as [min, max]
371
372     """
373     proxy.UpdatePipeline()
374     entity_data_info = None
375     field_data = proxy.GetFieldDataInformation()
376
377     if field_name in list(field_data.keys()):
378         entity_data_info = field_data
379     elif entity == EntityType.CELL:
380         entity_data_info = proxy.GetCellDataInformation()
381     elif entity == EntityType.NODE:
382         entity_data_info = proxy.GetPointDataInformation()
383
384     data_range = []
385
386     if field_name in list(entity_data_info.keys()):
387         vcomp = get_vector_component(vector_mode)
388         data_range = entity_data_info[field_name].GetComponentRange(vcomp)
389     else:
390         pv_entity = EntityType.get_pvtype(entity)
391         warnings.warn("Field " + field_name +
392                       " is unknown for " + pv_entity + "!")
393
394     # Cut off the range
395     if cut_off and (data_range[0] <= data_range[1]):
396         data_range = list(data_range)
397         delta = abs(data_range[1] - data_range[0]) * GAP_COEFFICIENT
398         data_range[0] += delta
399         data_range[1] -= delta
400
401     return data_range
402
403
404 def get_bounds(proxy):
405     """Get bounds of the proxy in 3D."""
406     proxy.UpdatePipeline()
407     dataInfo = proxy.GetDataInformation()
408     bounds_info = dataInfo.GetBounds()
409     return bounds_info
410
411
412 def get_x_range(proxy):
413     """Get X range of the proxy bounds in 3D."""
414     proxy.UpdatePipeline()
415     bounds_info = get_bounds(proxy)
416     return bounds_info[0:2]
417
418
419 def get_y_range(proxy):
420     """Get Y range of the proxy bounds in 3D."""
421     proxy.UpdatePipeline()
422     bounds_info = get_bounds(proxy)
423     return bounds_info[2:4]
424
425
426 def get_z_range(proxy):
427     """Get Z range of the proxy bounds in 3D."""
428     proxy.UpdatePipeline()
429     bounds_info = get_bounds(proxy)
430     return bounds_info[4:6]
431
432
433 def is_planar_input(proxy):
434     """Check if the given input is planar."""
435     proxy.UpdatePipeline()
436     bounds_info = get_bounds(proxy)
437
438     if (abs(bounds_info[0] - bounds_info[1]) <= FLT_MIN or
439         abs(bounds_info[2] - bounds_info[3]) <= FLT_MIN or
440         abs(bounds_info[4] - bounds_info[5]) <= FLT_MIN):
441         return True
442
443     return False
444
445
446 def is_data_on_cells(proxy, field_name):
447     """Check the existence of a field on cells with the given name."""
448     proxy.UpdatePipeline()
449     cell_data_info = proxy.GetCellDataInformation()
450     return (field_name in list(cell_data_info.keys()))
451
452
453 def is_empty(proxy):
454     """Check if the object contains any points or cells.
455
456     Returns:
457       True: if the given proxy doesn't contain any points or cells
458       False: otherwise
459
460     """
461     proxy.UpdatePipeline()
462     data_info = proxy.GetDataInformation()
463
464     nb_cells = data_info.GetNumberOfCells()
465     nb_points = data_info.GetNumberOfPoints()
466
467     return not(nb_cells + nb_points)
468
469
470 def get_orientation(proxy):
471     """Get the optimum cutting plane orientation for Plot 3D."""
472     proxy.UpdatePipeline()
473     orientation = Orientation.XY
474
475     bounds = get_bounds(proxy)
476     delta = [bounds[1] - bounds[0],
477              bounds[3] - bounds[2],
478              bounds[5] - bounds[4]]
479
480     if (delta[0] >= delta[1] and delta[0] >= delta[2]):
481         if (delta[1] >= delta[2]):
482             orientation = Orientation.XY
483         else:
484             orientation = Orientation.ZX
485     elif (delta[1] >= delta[0] and delta[1] >= delta[2]):
486         if (delta[0] >= delta[2]):
487             orientation = Orientation.XY
488         else:
489             orientation = Orientation.YZ
490     elif (delta[2] >= delta[0] and delta[2] >= delta[1]):
491         if (delta[0] >= delta[1]):
492             orientation = Orientation.ZX
493         else:
494             orientation = Orientation.YZ
495
496     return orientation
497
498
499 def dot_product(a, b):
500     """Dot product of two 3-vectors."""
501     dot = a[0] * b[0] + a[1] * b[1] + a[2] * b[2]
502     return dot
503
504
505 def multiply3x3(a, b):
506     """Mutltiply one 3x3 matrix by another."""
507     c = [[0, 0, 0],
508          [0, 0, 0],
509          [0, 0, 0]]
510
511     for i in range(3):
512         c[0][i] = a[0][0] * b[0][i] + a[0][1] * b[1][i] + a[0][2] * b[2][i]
513         c[1][i] = a[1][0] * b[0][i] + a[1][1] * b[1][i] + a[1][2] * b[2][i]
514         c[2][i] = a[2][0] * b[0][i] + a[2][1] * b[1][i] + a[2][2] * b[2][i]
515
516     return c
517
518
519 def get_rx(ang):
520     """Get X rotation matrix by angle."""
521     rx = [[1.0, 0.0,      0.0],
522           [0.0, cos(ang), -sin(ang)],
523           [0.0, sin(ang), cos(ang)]]
524
525     return rx
526
527
528 def get_ry(ang):
529     """Get Y rotation matrix by angle."""
530     ry = [[cos(ang),  0.0, sin(ang)],
531           [0.0,       1.0, 0.0],
532           [-sin(ang), 0.0, cos(ang)]]
533
534     return ry
535
536
537 def get_rz(ang):
538     """Get Z rotation matrix by angle."""
539     rz = [[cos(ang), -sin(ang), 0.0],
540           [sin(ang), cos(ang),  0.0],
541           [0.0,      0.0,       1.0]]
542
543     return rz
544
545
546 def get_normal_by_orientation(orientation, ang1=0, ang2=0):
547     """Get normal for the plane by its orientation."""
548     i_plane = 0
549     rotation = [[], [], []]
550     rx = ry = rz = [[1.0, 0.0, 0.0],
551                     [0.0, 1.0, 0.0],
552                     [0.0, 0.0, 1.0]]
553
554     normal = [0.0, 0.0, 0.0]
555     if orientation == Orientation.XY:
556         if abs(ang1) > EPS:
557             rx = get_rx(ang1)
558         if abs(ang2) > EPS:
559             ry = get_ry(ang2)
560         rotation = multiply3x3(rx, ry)
561         i_plane = 2
562     elif orientation == Orientation.ZX:
563         if abs(ang1) > EPS:
564             rz = get_rz(ang1)
565         if abs(ang2) > EPS:
566             rx = get_rx(ang2)
567         rotation = multiply3x3(rz, rx)
568         i_plane = 1
569     elif orientation == Orientation.YZ:
570         if abs(ang1) > EPS:
571             ry = get_ry(ang1)
572         if abs(ang2) > EPS:
573             rz = get_rz(ang2)
574         rotation = multiply3x3(ry, rz)
575         i_plane = 0
576
577     for i in range(0, 3):
578         normal[i] = rotation[i][i_plane]
579
580     return normal
581
582
583 def get_bound_project(bound_box, dir):
584     """Get bounds projection"""
585     bound_points = [[bound_box[0], bound_box[2], bound_box[4]],
586                     [bound_box[1], bound_box[2], bound_box[4]],
587                     [bound_box[0], bound_box[3], bound_box[4]],
588                     [bound_box[1], bound_box[3], bound_box[4]],
589                     [bound_box[0], bound_box[2], bound_box[5]],
590                     [bound_box[1], bound_box[2], bound_box[5]],
591                     [bound_box[0], bound_box[3], bound_box[5]],
592                     [bound_box[1], bound_box[3], bound_box[5]]]
593
594     bound_prj = [0, 0, 0]
595     bound_prj[0] = dot_product(dir, bound_points[0])
596     bound_prj[1] = bound_prj[0]
597
598     for i in range(1, 8):
599         tmp = dot_product(dir, bound_points[i])
600         if bound_prj[1] < tmp:
601             bound_prj[1] = tmp
602         if bound_prj[0] > tmp:
603             bound_prj[0] = tmp
604
605     bound_prj[2] = bound_prj[1] - bound_prj[0]
606     bound_prj[1] = bound_prj[0] + (1.0 - EPS) * bound_prj[2]
607     bound_prj[0] = bound_prj[0] + EPS * bound_prj[2]
608     bound_prj[2] = bound_prj[1] - bound_prj[0]
609
610     return bound_prj
611
612
613 def get_positions(nb_planes, dir, bounds, displacement):
614     """Compute plane positions."""
615     positions = []
616     bound_prj = get_bound_project(bounds, dir)
617     if nb_planes > 1:
618         step = bound_prj[2] / (nb_planes - 1)
619         abs_displacement = step * displacement
620         start_pos = bound_prj[0] - 0.5 * step + abs_displacement
621         for i in range(nb_planes):
622             pos = start_pos + i * step
623             positions.append(pos)
624     else:
625         pos = bound_prj[0] + bound_prj[2] * displacement
626         positions.append(pos)
627
628     return positions
629
630
631 def get_contours(scalar_range, nb_contours):
632     """Generate contour values."""
633     contours = []
634     for i in range(nb_contours):
635         pos = scalar_range[0] + i * (
636             scalar_range[1] - scalar_range[0]) / (nb_contours - 1)
637         contours.append(pos)
638
639     return contours
640
641
642 def get_nb_components(proxy, entity, field_name):
643     """Return number of components for the field."""
644     proxy.UpdatePipeline()
645     entity_data_info = None
646     field_data = proxy.GetFieldDataInformation()
647
648     if field_name in list(field_data.keys()):
649         entity_data_info = field_data
650     elif entity == EntityType.CELL:
651         select_cells_with_data(proxy, on_cells=[field_name])
652         entity_data_info = proxy.GetCellDataInformation()
653     elif entity == EntityType.NODE:
654         select_cells_with_data(proxy, on_points=[field_name])
655         entity_data_info = proxy.GetPointDataInformation()
656
657     nb_comp = None
658     if field_name in list(entity_data_info.keys()):
659         nb_comp = entity_data_info[field_name].GetNumberOfComponents()
660     else:
661         pv_entity = EntityType.get_pvtype(entity)
662         raise ValueError("Field " + field_name +
663                          " is unknown for " + pv_entity + "!")
664
665     return nb_comp
666
667
668 def get_scale_factor(proxy):
669     """Compute scale factor."""
670     if not proxy:
671         return 0.0
672
673     proxy.UpdatePipeline()
674     data_info = proxy.GetDataInformation()
675
676     nb_cells = data_info.GetNumberOfCells()
677     nb_points = data_info.GetNumberOfPoints()
678     nb_elements = nb_cells if nb_cells > 0  else nb_points
679     bounds = get_bounds(proxy)
680
681     volume = 1
682     vol = dim = 0
683
684     for i in range(0, 6, 2):
685         vol = abs(bounds[i + 1] - bounds[i])
686         if vol > 0:
687             dim += 1
688             volume *= vol
689
690     if nb_elements == 0 or dim < 1 / VTK_LARGE_FLOAT:
691         return 0
692
693     volume /= nb_elements
694
695     return pow(volume, 1 / dim)
696
697
698 def get_default_scale(prs_type, proxy, entity, field_name):
699     """Get default scale factor."""
700     proxy.UpdatePipeline()
701     data_range = get_data_range(proxy, entity, field_name)
702
703     if prs_type == PrsTypeEnum.DEFORMEDSHAPE:
704         EPS = 1.0 / VTK_LARGE_FLOAT
705         if abs(data_range[1]) > EPS:
706             scale_factor = get_scale_factor(proxy)
707             return scale_factor / data_range[1]
708     elif prs_type == PrsTypeEnum.PLOT3D:
709         bounds = get_bounds(proxy)
710         length = sqrt((bounds[1] - bounds[0]) ** 2 +
711                       (bounds[3] - bounds[2]) ** 2 +
712                       (bounds[5] - bounds[4]) ** 2)
713
714         EPS = 0.3
715         if data_range[1] > 0:
716             return length / data_range[1] * EPS
717
718     return 0
719
720
721 def get_calc_magnitude(proxy, array_entity, array_name):
722     """Compute magnitude for the given vector array via Calculator.
723
724     Returns:
725       the calculator object.
726
727     """
728     proxy.UpdatePipeline()
729     calculator = None
730
731     # Transform vector array to scalar array if possible
732     nb_components = get_nb_components(proxy, array_entity, array_name)
733     if (nb_components > 1):
734         calculator = pvs.Calculator(proxy)
735         attribute_mode = "Point Data"
736         if array_entity != EntityType.NODE:
737             attribute_mode = "Cell Data"
738         calculator.AttributeMode = attribute_mode
739         if (nb_components == 2):
740             # Workaroud: calculator unable to compute magnitude
741             # if number of components equal to 2
742             func = "sqrt(" + array_name + "_X^2+" + array_name + "_Y^2)"
743             calculator.Function = func
744         else:
745             calculator.Function = "mag(" + array_name + ")"
746         calculator.ResultArrayName = array_name + "_magnitude"
747         calculator.UpdatePipeline()
748
749     return calculator
750
751
752 def get_add_component_calc(proxy, array_entity, array_name):
753     """Creates 3-component array from 2-component.
754
755     The first two components is from the original array. The 3rd component
756     is zero.
757     If the number of components is not equal to 2 - return original array name.
758
759     Returns:
760       the calculator object.
761
762     """
763     proxy.UpdatePipeline()
764     calculator = None
765
766     nb_components = get_nb_components(proxy, array_entity, array_name)
767     if nb_components == 2:
768         calculator = pvs.Calculator(proxy)
769         attribute_mode = "Point Data"
770         if array_entity != EntityType.NODE:
771             attribute_mode = "Cell Data"
772         calculator.AttributeMode = attribute_mode
773         expression = "iHat * " + array_name + "_X + jHat * " + array_name + "_Y + kHat * 0"
774         calculator.Function = expression
775         calculator.ResultArrayName = array_name + "_3c"
776         calculator.UpdatePipeline()
777
778     return calculator
779
780
781 def select_all_cells(proxy):
782     """Select all cell types.
783
784     Used in creation of mesh/submesh presentation.
785
786     """
787     proxy.UpdatePipeline()
788     extractCT = pvs.ExtractCellType()
789     extractCT.AllGeoTypes = extractCT.GetProperty("GeoTypesInfo")[::2]
790     extractCT.UpdatePipelineInformation()
791
792
793 def select_cells_with_data(proxy, on_points=[], on_cells=[], on_gauss=[]):
794     """Select cell types with data.
795
796     Only cell types with data for the given fields will be selected.
797     If no fields defined (neither on points nor on cells) only cell
798     types with data for even one field (from available) will be selected.
799
800     """
801     if not proxy.GetProperty("FieldsTreeInfo"):
802         return
803
804     proxy.UpdatePipeline()
805     if not hasattr(proxy, 'Entity'):
806         fields_info = proxy.GetProperty("FieldsTreeInfo")[::2]
807         arr_name_with_dis=[elt.split("/")[-1] for elt in fields_info]
808
809         fields = []
810         for name in on_gauss:
811             fields.append(name+_med_field_sep+'GAUSS')
812         for name in on_cells:
813             fields.append(name+_med_field_sep+'P0')
814         for name in on_points:
815             fields.append(name+_med_field_sep+'P1')
816
817         field_list = []
818         for name in fields:
819             if arr_name_with_dis.count(name) > 0:
820                 index = arr_name_with_dis.index(name)
821                 field_list.append(fields_info[index])
822
823         if field_list:
824             proxy.AllArrays = field_list
825             proxy.UpdatePipeline()
826         return len(field_list) != 0
827
828     # TODO: VTN. Looks like this code is out of date.
829
830     #all_cell_types = proxy.CellTypes.Available
831     all_cell_types = proxy.Entity.Available
832     all_arrays = list(proxy.CellArrays.GetData())
833     all_arrays.extend(proxy.PointArrays.GetData())
834
835     if not all_arrays:
836         file_name = proxy.FileName.split(os.sep)[-1]
837         print("Warning: " + file_name + " doesn't contain any data array.")
838
839     # List of cell types to be selected
840     cell_types_on = []
841
842     for cell_type in all_cell_types:
843         #proxy.CellTypes = [cell_type]
844         proxy.Entity = [cell_type]
845         proxy.UpdatePipeline()
846
847         cell_arrays = list(proxy.GetCellDataInformation().keys())
848         point_arrays = list(proxy.GetPointDataInformation().keys())
849
850         if on_points or on_cells:
851             if on_points is None:
852                 on_points = []
853             if on_cells is None:
854                 on_cells = []
855
856             if (all(array in cell_arrays for array in on_cells) and
857                 all(array in point_arrays for array in on_points)):
858                 # Add cell type to the list
859                 cell_types_on.append(cell_type)
860         else:
861             in_arrays = lambda array: ((array in cell_arrays) or
862                                        (array in point_arrays))
863             if any(in_arrays(array) for array in all_arrays):
864                 cell_types_on.append(cell_type)
865
866     # Select cell types
867     #proxy.CellTypes = cell_types_on
868     proxy.Entity = cell_types_on
869     proxy.UpdatePipeline()
870
871 def if_possible(proxy, field_name, entity, prs_type, extrGrps=None):
872     """Check if the presentation creation is possible on the given field."""
873     proxy.UpdatePipeline()
874     result = True
875     if (prs_type == PrsTypeEnum.DEFORMEDSHAPE or
876         prs_type == PrsTypeEnum.DEFORMEDSHAPESCALARMAP or
877         prs_type == PrsTypeEnum.VECTORS or
878         prs_type == PrsTypeEnum.STREAMLINES):
879         nb_comp = get_nb_components(proxy, entity, field_name)
880         result = (nb_comp > 1)
881     elif (prs_type == PrsTypeEnum.GAUSSPOINTS):
882         result = (entity == EntityType.CELL or
883                   field_name in proxy.QuadraturePointArrays.Available)
884     elif (prs_type == PrsTypeEnum.MESH):
885         result = len(get_group_names(extrGrps)) > 0
886
887     return result
888
889
890 def add_scalar_bar(field_name, nb_components,
891                    vector_mode, lookup_table, time_value):
892     """Add scalar bar with predefined properties."""
893     global _current_bar
894
895     # Construct bar title
896     title = "\n".join([field_name, str(time_value)])
897     if nb_components > 1:
898         title = "\n".join([title, vector_mode])
899
900     # Create scalar bar
901     scalar_bar = pvs.CreateScalarBar(Enabled=1)
902     scalar_bar.Orientation = 'Vertical'
903     scalar_bar.Title = title
904     scalar_bar.LookupTable = lookup_table
905
906     # Set default properties same as in Post-Pro
907     scalar_bar.NumberOfLabels = 5
908     scalar_bar.AutomaticLabelFormat = 0
909     scalar_bar.LabelFormat = '%-#6.6g'
910     # Title
911     scalar_bar.TitleFontFamily = 'Arial'
912     scalar_bar.TitleFontSize = 8
913     scalar_bar.TitleBold = 1
914     scalar_bar.TitleItalic = 1
915     scalar_bar.TitleShadow = 1
916     # Labels
917     scalar_bar.LabelFontFamily = 'Arial'
918     scalar_bar.LabelFontSize = 8
919     scalar_bar.LabelBold = 1
920     scalar_bar.LabelItalic = 1
921     scalar_bar.LabelShadow = 1
922
923     # Add the scalar bar to the view
924     pvs.GetRenderView().Representations.append(scalar_bar)
925
926     # Reassign the current bar
927     _current_bar = scalar_bar
928
929     return _current_bar
930
931
932 def get_bar():
933     """Get current scalar bar."""
934     return _current_bar
935
936
937 def get_lookup_table(field_name, nb_components, vector_mode='Magnitude'):
938     """Get lookup table for the given field."""
939     lookup_table = pvs.GetLookupTableForArray(field_name, nb_components)
940
941     if vector_mode == 'Magnitude':
942         lookup_table.VectorMode = vector_mode
943     elif vector_mode == 'X':
944         lookup_table.VectorMode = 'Component'
945         lookup_table.VectorComponent = 0
946     elif vector_mode == 'Y':
947         lookup_table.VectorMode = 'Component'
948         lookup_table.VectorComponent = 1
949     elif vector_mode == 'Z':
950         lookup_table.VectorMode = 'Component'
951         lookup_table.VectorComponent = 2
952     else:
953         raise ValueError("Incorrect vector mode: " + vector_mode)
954
955     lookup_table.Discretize = 0
956     lookup_table.ColorSpace = 'HSV'
957     if hasattr(lookup_table,"LockDataRange"):
958         lookup_table.LockDataRange = 0
959     elif hasattr(lookup_table,"LockScalarRange"):
960         lookup_table.LockScalarRange = 0
961     else:
962         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
963
964     return lookup_table
965
966
967 def get_group_mesh_name(full_group_name):
968     """Return mesh name of the group by its full name."""
969     aList = full_group_name.split('/')
970     if len(aList) >= 2 :
971         group_name = full_group_name.split('/')[1]
972         return group_name
973
974 def get_group_entity(full_group_name):
975     """Return entity type of the group by its full name."""
976     aList = full_group_name.split('/')
977     if len(aList) >= 3 :
978         entity_name = full_group_name.split('/')[2]
979         entity = EntityType.get_type(entity_name)
980         return entity
981
982
983 def get_group_short_name(full_group_name):
984     """Return short name of the group by its full name."""
985     short_name = re.sub('^GRP_', '', full_group_name)
986     return short_name
987
988
989 def get_mesh_full_names(proxy):
990     """Return all mesh names in the given proxy as a set."""
991     proxy.UpdatePipeline()
992     fields = proxy.GetProperty("FieldsTreeInfo")[::2]
993     mesh_full_names = set([item for item in fields if get_field_mesh_name(item) == get_field_short_name(item)])
994     return mesh_full_names
995
996
997 def get_group_names(extrGrps):
998     """Return full names of all groups of the given 'ExtractGroup' filter object.
999     """
1000     group_names = [x for x in list(extrGrps.GetProperty("GroupsFlagsInfo")[::2]) if x[:4]=="GRP_"]
1001     return group_names
1002
1003
1004 def get_time(proxy, timestamp_nb):
1005     """Get time value by timestamp number."""
1006     #proxy.UpdatePipeline()
1007     # Check timestamp number
1008     timestamps = []
1009
1010     if (hasattr(proxy, 'TimestepValues')):
1011         timestamps = proxy.TimestepValues.GetData()
1012     elif (hasattr(proxy.Input, 'TimestepValues')):
1013         timestamps = proxy.Input.TimestepValues.GetData()
1014
1015     length = len(timestamps)
1016     if (timestamp_nb > 0 and (timestamp_nb - 1) not in range(length) ) or (timestamp_nb < 0 and -timestamp_nb > length):
1017         raise ValueError("Timestamp number is out of range: " + str(timestamp_nb))
1018
1019     if not timestamps:
1020         return 0.0
1021
1022     # Return time value
1023     if timestamp_nb > 0:
1024         return timestamps[timestamp_nb - 1]
1025     else:
1026         return timestamps[timestamp_nb]
1027
1028 def create_prs(prs_type, proxy, field_entity, field_name, timestamp_nb):
1029     """Auxiliary function.
1030
1031     Build presentation of the given type on the given field and
1032     timestamp number.
1033     Set the presentation properties like visu.CreatePrsForResult() do.
1034
1035     """
1036     proxy.UpdatePipeline()
1037     prs = None
1038
1039     if prs_type == PrsTypeEnum.SCALARMAP:
1040         prs = ScalarMapOnField(proxy, field_entity, field_name, timestamp_nb)
1041     elif prs_type == PrsTypeEnum.CUTPLANES:
1042         prs = CutPlanesOnField(proxy, field_entity, field_name, timestamp_nb,
1043                                orientation=Orientation.ZX)
1044     elif prs_type == PrsTypeEnum.CUTLINES:
1045         prs = CutLinesOnField(proxy, field_entity, field_name, timestamp_nb,
1046                               orientation1=Orientation.XY,
1047                               orientation2=Orientation.ZX)
1048     elif prs_type == PrsTypeEnum.DEFORMEDSHAPE:
1049         prs = DeformedShapeOnField(proxy, field_entity,
1050                                    field_name, timestamp_nb)
1051     elif prs_type == PrsTypeEnum.DEFORMEDSHAPESCALARMAP:
1052         prs = DeformedShapeAndScalarMapOnField(proxy, field_entity,
1053                                                field_name, timestamp_nb)
1054     elif prs_type == PrsTypeEnum.VECTORS:
1055         prs = VectorsOnField(proxy, field_entity, field_name, timestamp_nb)
1056     elif prs_type == PrsTypeEnum.PLOT3D:
1057         prs = Plot3DOnField(proxy, field_entity, field_name, timestamp_nb)
1058     elif prs_type == PrsTypeEnum.ISOSURFACES:
1059         prs = IsoSurfacesOnField(proxy, field_entity, field_name, timestamp_nb)
1060     elif prs_type == PrsTypeEnum.GAUSSPOINTS:
1061         prs = GaussPointsOnField(proxy, field_entity, field_name, timestamp_nb)
1062     elif prs_type == PrsTypeEnum.STREAMLINES:
1063         prs = StreamLinesOnField(proxy, field_entity, field_name, timestamp_nb)
1064     else:
1065         raise ValueError("Unexistent presentation type.")
1066
1067     return prs
1068
1069
1070 # Functions for building Post-Pro presentations
1071 def ScalarMapOnField(proxy, entity, field_name, timestamp_nb,
1072                      vector_mode='Magnitude'):
1073     """Creates Scalar Map presentation on the given field.
1074
1075     Arguments:
1076       proxy: the pipeline object, containig data
1077       entity: the entity type from PrsTypeEnum
1078       field_name: the field name
1079       timestamp_nb: the number of time step (1, 2, ...)
1080       vector_mode: the mode of transformation of vector values
1081       into scalar values, applicable only if the field contains vector values.
1082       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1083
1084     Returns:
1085       Scalar Map as representation object.
1086
1087     """
1088     proxy.UpdatePipeline()
1089     # We don't need mesh parts with no data on them
1090     if entity == EntityType.NODE:
1091         select_cells_with_data(proxy, on_points=[field_name])
1092     else:
1093         select_cells_with_data(proxy, on_cells=[field_name])
1094
1095     # Check vector mode
1096     nb_components = get_nb_components(proxy, entity, field_name)
1097     check_vector_mode(vector_mode, nb_components)
1098
1099     # Get time value
1100     time_value = get_time(proxy, timestamp_nb)
1101
1102     # Set timestamp
1103     pvs.GetRenderView().ViewTime = time_value
1104     pvs.UpdatePipeline(time_value, proxy)
1105
1106     # Get Scalar Map representation object
1107     scalarmap = pvs.GetRepresentation(proxy)
1108
1109     # Get lookup table
1110     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1111
1112     # Set field range if necessary
1113     data_range = get_data_range(proxy, entity,
1114                                 field_name, vector_mode)
1115     if hasattr(lookup_table,"LockDataRange"):
1116         lookup_table.LockDataRange = 1
1117     elif hasattr(lookup_table,"LockScalarRange"):
1118         lookup_table.LockScalarRange = 1
1119     else:
1120         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1121
1122     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1123     # Set properties
1124     pvs.ColorBy(scalarmap, (EntityType.get_pvtype(entity), field_name))
1125     scalarmap.LookupTable = lookup_table
1126
1127     # Add scalar bar
1128     bar_title = field_name + ", " + str(time_value)
1129     if (nb_components > 1):
1130         bar_title += "\n" + vector_mode
1131     add_scalar_bar(field_name, nb_components, vector_mode,
1132                    lookup_table, time_value)
1133
1134     return scalarmap
1135
1136
1137 def CutPlanesOnField(proxy, entity, field_name, timestamp_nb,
1138                      nb_planes=10, orientation=Orientation.YZ,
1139                      angle1=0, angle2=0,
1140                      displacement=0.5, vector_mode='Magnitude'):
1141     """Creates Cut Planes presentation on the given field.
1142
1143     Arguments:
1144       proxy: the pipeline object, containig data
1145       entity: the entity type from PrsTypeEnum
1146       field_name: the field name
1147       timestamp_nb: the number of time step (1, 2, ...)
1148       nb_planes: number of cutting planes
1149       orientation: cutting planes orientation in 3D space
1150       angle1: rotation of the planes in 3d space around the first axis of the
1151       selected orientation (X axis for XY, Y axis for YZ, Z axis for ZX).
1152       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1153       angle2: rotation of the planes in 3d space around the second axis of the
1154       selected orientation. Acceptable range: [-45, 45].
1155       displacement: the displacement of the planes into one or another side
1156       vector_mode: the mode of transformation of vector values
1157       into scalar values, applicable only if the field contains vector values.
1158       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1159
1160     Returns:
1161       Cut Planes as representation object.
1162
1163     """
1164     proxy.UpdatePipeline()
1165     if entity == EntityType.NODE:
1166         select_cells_with_data(proxy, on_points=[field_name])
1167     else:
1168         select_cells_with_data(proxy, on_cells=[field_name])
1169
1170     # Check vector mode
1171     nb_components = get_nb_components(proxy, entity, field_name)
1172     check_vector_mode(vector_mode, nb_components)
1173
1174     # Get time value
1175     time_value = get_time(proxy, timestamp_nb)
1176
1177     # Set timestamp
1178     pvs.GetRenderView().ViewTime = time_value
1179     pvs.UpdatePipeline(time_value, proxy)
1180
1181     # Create slice filter
1182     slice_filter = pvs.Slice(proxy)
1183     slice_filter.SliceType = "Plane"
1184
1185     # Set cut planes normal
1186     normal = get_normal_by_orientation(orientation,
1187                                        radians(angle1), radians(angle2))
1188     slice_filter.SliceType.Normal = normal
1189
1190     # Set cut planes positions
1191     positions = get_positions(nb_planes, normal,
1192                               get_bounds(proxy), displacement)
1193     slice_filter.SliceOffsetValues = positions
1194
1195     # Get Cut Planes representation object
1196     cut_planes = pvs.GetRepresentation(slice_filter)
1197
1198     # Get lookup table
1199     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1200
1201     # Set field range if necessary
1202     data_range = get_data_range(proxy, entity,
1203                                 field_name, vector_mode)
1204
1205     if hasattr(lookup_table,"LockDataRange"):
1206         lookup_table.LockDataRange = 1
1207     elif hasattr(lookup_table,"LockScalarRange"):
1208         lookup_table.LockScalarRange = 1
1209     else:
1210         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1211
1212     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1213
1214     # Set properties
1215     pvs.ColorBy(cut_planes, (EntityType.get_pvtype(entity), field_name))
1216     cut_planes.LookupTable = lookup_table
1217
1218     # Add scalar bar
1219     add_scalar_bar(field_name, nb_components,
1220                    vector_mode, lookup_table, time_value)
1221
1222     return cut_planes
1223
1224
1225 def CutLinesOnField(proxy, entity, field_name, timestamp_nb,
1226                     nb_lines=10,
1227                     orientation1=Orientation.XY,
1228                     base_angle1=0, base_angle2=0,
1229                     orientation2=Orientation.YZ,
1230                     cut_angle1=0, cut_angle2=0,
1231                     displacement1=0.5, displacement2=0.5,
1232                     generate_curves=False,
1233                     vector_mode='Magnitude'):
1234     """Creates Cut Lines presentation on the given field.
1235
1236     Arguments:
1237       proxy: the pipeline object, containig data
1238       entity: the entity type from PrsTypeEnum
1239       field_name: the field name
1240       timestamp_nb: the number of time step (1, 2, ...)
1241       nb_lines: number of lines
1242       orientation1: base plane orientation in 3D space
1243       base_angle1: rotation of the base plane in 3d space around the first
1244       axis of the orientation1 (X axis for XY, Y axis for YZ, Z axis for ZX).
1245       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1246       base_angle2: rotation of the base plane in 3d space around the second
1247       axis of the orientation1. Acceptable range: [-45, 45].
1248       orientation2: cutting planes orientation in 3D space
1249       cut_angle1: rotation of the cut planes in 3d space around the first
1250       axis of the orientation2. Acceptable range: [-45, 45].
1251       cut_angle2: rotation of the cuting planes in 3d space around the second
1252       axis of the orientation2. Acceptable range: [-45, 45].
1253       displacement1: base plane displacement
1254       displacement2: cutting planes displacement
1255       generate_curves: if true, 'PlotOverLine' filter will be created
1256       for each cut line
1257       vector_mode: the mode of transformation of vector values
1258       into scalar values, applicable only if the field contains vector values.
1259       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1260
1261     Returns:
1262       Cut Lines as representation object if generate_curves == False,
1263       (Cut Lines as representation object, list of 'PlotOverLine') otherwise
1264
1265     """
1266     proxy.UpdatePipeline()
1267     if entity == EntityType.NODE:
1268         select_cells_with_data(proxy, on_points=[field_name])
1269     else:
1270         select_cells_with_data(proxy, on_cells=[field_name])
1271
1272     # Check vector mode
1273     nb_components = get_nb_components(proxy, entity, field_name)
1274     check_vector_mode(vector_mode, nb_components)
1275
1276     # Get time value
1277     time_value = get_time(proxy, timestamp_nb)
1278
1279     # Set timestamp
1280     pvs.GetRenderView().ViewTime = time_value
1281     pvs.UpdatePipeline(time_value, proxy)
1282
1283     # Create base plane
1284     base_plane = pvs.Slice(proxy)
1285     base_plane.SliceType = "Plane"
1286
1287     # Set base plane normal
1288     base_normal = get_normal_by_orientation(orientation1,
1289                                             radians(base_angle1),
1290                                             radians(base_angle2))
1291     base_plane.SliceType.Normal = base_normal
1292
1293     # Set base plane position
1294     base_position = get_positions(1, base_normal,
1295                                   get_bounds(proxy), displacement1)
1296     base_plane.SliceOffsetValues = base_position
1297
1298     # Check base plane
1299     base_plane.UpdatePipeline()
1300     if (base_plane.GetDataInformation().GetNumberOfCells() == 0):
1301         base_plane = proxy
1302
1303     # Create cutting planes
1304     cut_planes = pvs.Slice(base_plane)
1305     cut_planes.SliceType = "Plane"
1306
1307     # Set cutting planes normal and get positions
1308     cut_normal = get_normal_by_orientation(orientation2,
1309                                            radians(cut_angle1),
1310                                            radians(cut_angle2))
1311     cut_planes.SliceType.Normal = cut_normal
1312
1313     # Set cutting planes position
1314     cut_positions = get_positions(nb_lines, cut_normal,
1315                                   get_bounds(base_plane), displacement2)
1316
1317     # Generate curves
1318     curves = []
1319     if generate_curves:
1320         index = 0
1321         for pos in cut_positions:
1322             # Get points for plot over line objects
1323             cut_planes.SliceOffsetValues = pos
1324             cut_planes.UpdatePipeline()
1325             bounds = get_bounds(cut_planes)
1326             point1 = [bounds[0], bounds[2], bounds[4]]
1327             point2 = [bounds[1], bounds[3], bounds[5]]
1328
1329             # Create plot over line filter
1330             pol = pvs.PlotOverLine(cut_planes,
1331                                   Source="High Resolution Line Source")
1332             pvs.RenameSource('Y' + str(index), pol)
1333             pol.Source.Point1 = point1
1334             pol.Source.Point2 = point2
1335             pol.UpdatePipeline()
1336             curves.append(pol)
1337
1338             index += 1
1339
1340     cut_planes.SliceOffsetValues = cut_positions
1341     cut_planes.UpdatePipeline()
1342
1343     # Get Cut Lines representation object
1344     cut_lines = pvs.GetRepresentation(cut_planes)
1345
1346     # Get lookup table
1347     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1348
1349     # Set field range if necessary
1350     data_range = get_data_range(proxy, entity,
1351                                 field_name, vector_mode)
1352     if hasattr(lookup_table,"LockDataRange"):
1353         lookup_table.LockDataRange = 1
1354     elif hasattr(lookup_table,"LockScalarRange"):
1355         lookup_table.LockScalarRange = 1
1356     else:
1357         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1358
1359     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1360
1361     # Set properties
1362     pvs.ColorBy(cut_lines, (EntityType.get_pvtype(entity), field_name))
1363     cut_lines.LookupTable = lookup_table
1364
1365     # Set wireframe represenatation mode
1366     cut_lines.Representation = 'Wireframe'
1367
1368     # Add scalar bar
1369     add_scalar_bar(field_name, nb_components,
1370                    vector_mode, lookup_table, time_value)
1371
1372     result = cut_lines
1373     # If curves were generated return tuple (cut lines, list of curves)
1374     if curves:
1375         result = cut_lines, curves
1376
1377     return result
1378
1379
1380 def CutSegmentOnField(proxy, entity, field_name, timestamp_nb,
1381                       point1, point2, vector_mode='Magnitude'):
1382     """Creates Cut Segment presentation on the given field.
1383
1384     Arguments:
1385       proxy: the pipeline object, containig data
1386       entity: the entity type from PrsTypeEnum
1387       field_name: the field name
1388       timestamp_nb: the number of time step (1, 2, ...)
1389       point1: set the first point of the segment (as [x, y, z])
1390       point1: set the second point of the segment (as [x, y, z])
1391       vector_mode: the mode of transformation of vector values
1392       into scalar values, applicable only if the field contains vector values.
1393       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1394
1395     Returns:
1396       Cut Segment as 3D representation object.
1397
1398     """
1399     proxy.UpdatePipeline()
1400     if entity == EntityType.NODE:
1401         select_cells_with_data(proxy, on_points=[field_name])
1402     else:
1403         select_cells_with_data(proxy, on_cells=[field_name])
1404
1405     # Check vector mode
1406     nb_components = get_nb_components(proxy, entity, field_name)
1407     check_vector_mode(vector_mode, nb_components)
1408
1409     # Get time value
1410     time_value = get_time(proxy, timestamp_nb)
1411
1412     # Set timestamp
1413     pvs.GetRenderView().ViewTime = time_value
1414     pvs.UpdatePipeline(time_value, proxy)
1415
1416     # Create plot over line filter
1417     pol = pvs.PlotOverLine(proxy, Source="High Resolution Line Source")
1418     pol.Source.Point1 = point1
1419     pol.Source.Point2 = point2
1420     pol.UpdatePipeline()
1421
1422     # Get Cut Segment representation object
1423     cut_segment = pvs.GetRepresentation(pol)
1424
1425     # Get lookup table
1426     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1427
1428     # Set field range if necessary
1429     data_range = get_data_range(proxy, entity,
1430                                 field_name, vector_mode)
1431     if hasattr(lookup_table,"LockDataRange"):
1432         lookup_table.LockDataRange = 1
1433     elif hasattr(lookup_table,"LockScalarRange"):
1434         lookup_table.LockScalarRange = 1
1435     else:
1436         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1437
1438     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1439
1440     # Set properties
1441     pvs.ColorBy(cut_segment, (EntityType.get_pvtype(entity), field_name))
1442     cut_segment.LookupTable = lookup_table
1443
1444     # Set wireframe represenatation mode
1445     cut_segment.Representation = 'Wireframe'
1446
1447     # Add scalar bar
1448     add_scalar_bar(field_name, nb_components,
1449                    vector_mode, lookup_table, time_value)
1450
1451     return cut_segment
1452
1453
1454 def VectorsOnField(proxy, entity, field_name, timestamp_nb,
1455                    scale_factor=None,
1456                    glyph_pos=GlyphPos.TAIL, glyph_type='2D Glyph',
1457                    is_colored=False, vector_mode='Magnitude'):
1458     """Creates Vectors presentation on the given field.
1459
1460     Arguments:
1461       proxy: the pipeline object, containig data
1462       entity: the entity type from PrsTypeEnum
1463       field_name: the field name
1464       timestamp_nb: the number of time step (1, 2, ...)
1465       scale_factor: scale factor
1466       glyph_pos: the position of glyphs
1467       glyph_type: the type of glyphs
1468       is_colored: this option allows to color the presentation according to
1469       the corresponding data array values
1470       vector_mode: the mode of transformation of vector values
1471       into scalar values, applicable only if the field contains vector values.
1472       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1473
1474     Returns:
1475       Vectors as representation object.
1476
1477     """
1478     proxy.UpdatePipeline()
1479     if entity == EntityType.NODE:
1480         select_cells_with_data(proxy, on_points=[field_name])
1481     else:
1482         select_cells_with_data(proxy, on_cells=[field_name])
1483
1484     # Check vector mode
1485     nb_components = get_nb_components(proxy, entity, field_name)
1486     check_vector_mode(vector_mode, nb_components)
1487
1488     # Get time value
1489     time_value = get_time(proxy, timestamp_nb)
1490
1491     # Set timestamp
1492     pvs.GetRenderView().ViewTime = time_value
1493     pvs.UpdatePipeline(time_value, proxy)
1494
1495     # Extract only groups with data for the field
1496     source = proxy
1497
1498     # Cell centers
1499     if is_data_on_cells(proxy, field_name):
1500         cell_centers = pvs.CellCenters(source)
1501         cell_centers.VertexCells = 1
1502         source = cell_centers
1503
1504     vector_array = field_name
1505     # If the given vector array has only 2 components, add the third one
1506     if nb_components == 2:
1507         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1508         vector_array = calc.ResultArrayName
1509         source = calc
1510
1511     # Glyph
1512     glyph = pvs.Glyph(source)
1513     glyph.Vectors = vector_array
1514     glyph.ScaleMode = 'vector'
1515     #glyph.MaskPoints = 0
1516
1517     # Set glyph type
1518     glyph.GlyphType = glyph_type
1519     if glyph_type == '2D Glyph':
1520         glyph.GlyphType.GlyphType = 'Arrow'
1521     elif glyph_type == 'Cone':
1522         glyph.GlyphType.Resolution = 7
1523         glyph.GlyphType.Height = 2
1524         glyph.GlyphType.Radius = 0.2
1525
1526     # Set glyph position if possible
1527     if glyph.GlyphType.GetProperty("Center"):
1528         if (glyph_pos == GlyphPos.TAIL):
1529             glyph.GlyphType.Center = [0.5, 0.0, 0.0]
1530         elif (glyph_pos == GlyphPos.HEAD):
1531             glyph.GlyphType.Center = [-0.5, 0.0, 0.0]
1532         elif (glyph_pos == GlyphPos.CENTER):
1533             glyph.GlyphType.Center = [0.0, 0.0, 0.0]
1534
1535     if scale_factor is not None:
1536         glyph.ScaleFactor = scale_factor
1537     else:
1538         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1539                                       proxy, entity, field_name)
1540         glyph.ScaleFactor = def_scale
1541
1542     glyph.UpdatePipeline()
1543
1544     # Get Vectors representation object
1545     vectors = pvs.GetRepresentation(glyph)
1546
1547     # Get lookup table
1548     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1549
1550     # Set field range if necessary
1551     data_range = get_data_range(proxy, entity,
1552                                 field_name, vector_mode)
1553     if hasattr(lookup_table,"LockDataRange"):
1554         lookup_table.LockDataRange = 1
1555     elif hasattr(lookup_table,"LockScalarRange"):
1556         lookup_table.LockScalarRange = 1
1557     else:
1558         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1559
1560     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1561
1562     # Set properties
1563     if (is_colored):
1564         pvs.ColorBy(vectors, (EntityType.get_pvtype(entity), 'GlyphVector'))
1565     else:
1566         pvs.ColorBy(vectors, (EntityType.get_pvtype(entity), None))
1567     vectors.LookupTable = lookup_table
1568
1569     vectors.LineWidth = 1.0
1570
1571     # Set wireframe represenatation mode
1572     vectors.Representation = 'Wireframe'
1573
1574     # Add scalar bar
1575     add_scalar_bar(field_name, nb_components,
1576                    vector_mode, lookup_table, time_value)
1577
1578     return vectors
1579
1580
1581 def DeformedShapeOnField(proxy, entity, field_name,
1582                          timestamp_nb,
1583                          scale_factor=None, is_colored=False,
1584                          vector_mode='Magnitude'):
1585     """Creates Defromed Shape presentation on the given field.
1586
1587     Arguments:
1588       proxy: the pipeline object, containig data
1589       entity: the entity type from PrsTypeEnum
1590       field_name: the field name
1591       timestamp_nb: the number of time step (1, 2, ...)
1592       scale_factor: scale factor of the deformation
1593       is_colored: this option allows to color the presentation according to
1594       the corresponding data array values
1595       vector_mode: the mode of transformation of vector values
1596       into scalar values, applicable only if the field contains vector values.
1597       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1598
1599     Returns:
1600       Defromed Shape as representation object.
1601
1602     """
1603     proxy.UpdatePipeline()
1604     # We don't need mesh parts with no data on them
1605     if entity == EntityType.NODE:
1606         select_cells_with_data(proxy, on_points=[field_name])
1607     else:
1608         select_cells_with_data(proxy, on_cells=[field_name])
1609
1610     # Check vector mode
1611     nb_components = get_nb_components(proxy, entity, field_name)
1612     check_vector_mode(vector_mode, nb_components)
1613
1614     # Get time value
1615     time_value = get_time(proxy, timestamp_nb)
1616
1617     # Set timestamp
1618     pvs.GetRenderView().ViewTime = time_value
1619     pvs.UpdatePipeline(time_value, proxy)
1620
1621     # Do merge
1622     source = pvs.MergeBlocks(proxy)
1623     pvs.UpdatePipeline()
1624
1625     # Cell data to point data
1626     if is_data_on_cells(proxy, field_name):
1627         cell_to_point = pvs.CellDatatoPointData()
1628         cell_to_point.PassCellData = 1
1629         source = cell_to_point
1630
1631     vector_array = field_name
1632     # If the given vector array has only 2 components, add the third one
1633     if nb_components == 2:
1634         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1635         vector_array = calc.ResultArrayName
1636         source = calc
1637
1638     # Warp by vector
1639     warp_vector = pvs.WarpByVector(source)
1640     warp_vector.Vectors = [vector_array]
1641     if scale_factor is not None:
1642         warp_vector.ScaleFactor = scale_factor
1643     else:
1644         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1645                                       proxy, entity, field_name)
1646         warp_vector.ScaleFactor = def_scale
1647
1648     # Get Deformed Shape representation object
1649     defshape = pvs.GetRepresentation(warp_vector)
1650
1651     # Get lookup table
1652     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1653
1654     # Set field range if necessary
1655     data_range = get_data_range(proxy, entity,
1656                                 field_name, vector_mode)
1657     if hasattr(lookup_table,"LockDataRange"):
1658         lookup_table.LockDataRange = 1
1659     elif hasattr(lookup_table,"LockScalarRange"):
1660         lookup_table.LockScalarRange = 1
1661     else:
1662         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1663
1664     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1665
1666     # Set properties
1667     if is_colored:
1668         pvs.ColorBy(defshape, (EntityType.get_pvtype(entity), field_name))
1669     else:
1670         pvs.ColorBy(defshape, (EntityType.get_pvtype(entity), None))
1671     defshape.LookupTable = lookup_table
1672
1673     # Set wireframe represenatation mode
1674     defshape.Representation = 'Wireframe'
1675
1676     # Add scalar bar
1677     add_scalar_bar(field_name, nb_components,
1678                    vector_mode, lookup_table, time_value)
1679
1680     return defshape
1681
1682
1683 def DeformedShapeAndScalarMapOnField(proxy, entity, field_name,
1684                                      timestamp_nb,
1685                                      scale_factor=None,
1686                                      scalar_entity=None,
1687                                      scalar_field_name=None,
1688                                      vector_mode='Magnitude'):
1689     """Creates Defromed Shape And Scalar Map presentation on the given field.
1690
1691     Arguments:
1692       proxy: the pipeline object, containig data
1693       entity: the entity type from PrsTypeEnum
1694       field_name: the field name
1695       timestamp_nb: the number of time step (1, 2, ...)
1696       scale_factor: scale factor of the deformation
1697       scalar_entity: scalar field entity
1698       scalar_field_name: scalar field, i.e. the field for coloring
1699       vector_mode: the mode of transformation of vector values
1700       into scalar values, applicable only if the field contains vector values.
1701       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1702
1703     Returns:
1704       Defromed Shape And Scalar Map as representation object.
1705
1706     """
1707     proxy.UpdatePipeline()
1708     # We don't need mesh parts with no data on them
1709     on_points = []
1710     on_cells = []
1711
1712     if entity == EntityType.NODE:
1713         on_points.append(field_name)
1714     else:
1715         on_cells.append(field_name)
1716
1717     if scalar_entity and scalar_field_name:
1718         if scalar_entity == EntityType.NODE:
1719             on_points.append(scalar_field_name)
1720         else:
1721             on_cells.append(scalar_field_name)
1722
1723     nb_components = get_nb_components(proxy, entity, field_name)
1724
1725     # Select fields
1726     select_cells_with_data(proxy, on_points, on_cells)
1727
1728     # Check vector mode
1729     check_vector_mode(vector_mode, nb_components)
1730
1731     # Get time value
1732     time_value = get_time(proxy, timestamp_nb)
1733
1734     # Set timestamp
1735     pvs.GetRenderView().ViewTime = time_value
1736     pvs.UpdatePipeline(time_value, proxy)
1737
1738     # Set scalar field by default
1739     scalar_field_entity = scalar_entity
1740     scalar_field = scalar_field_name
1741     if (scalar_field_entity is None) or (scalar_field is None):
1742         scalar_field_entity = entity
1743         scalar_field = field_name
1744
1745     # Do merge
1746     source = pvs.MergeBlocks(proxy)
1747     pvs.UpdatePipeline()
1748
1749     # Cell data to point data
1750     if is_data_on_cells(proxy, field_name):
1751         cell_to_point = pvs.CellDatatoPointData(source)
1752         cell_to_point.PassCellData = 1
1753         source = cell_to_point
1754
1755     vector_array = field_name
1756     # If the given vector array has only 2 components, add the third one
1757     if nb_components == 2:
1758         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1759         vector_array = calc.ResultArrayName
1760         source = calc
1761
1762     # Warp by vector
1763     warp_vector = pvs.WarpByVector(source)
1764     warp_vector.Vectors = [vector_array]
1765     if scale_factor is not None:
1766         warp_vector.ScaleFactor = scale_factor
1767     else:
1768         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1769                                       proxy, entity, field_name)
1770         warp_vector.ScaleFactor = def_scale
1771
1772     # Get Defromed Shape And Scalar Map representation object
1773     defshapemap = pvs.GetRepresentation(warp_vector)
1774
1775     # Get lookup table
1776     lookup_table = get_lookup_table(scalar_field, nb_components, vector_mode)
1777
1778     # Set field range if necessary
1779     data_range = get_data_range(proxy, scalar_field_entity,
1780                                 scalar_field, vector_mode)
1781     if hasattr(lookup_table,"LockDataRange"):
1782         lookup_table.LockDataRange = 1
1783     elif hasattr(lookup_table,"LockScalarRange"):
1784         lookup_table.LockScalarRange = 1
1785     else:
1786         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1787
1788     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1789
1790     # Set properties
1791     pvs.ColorBy(defshapemap, (EntityType.get_pvtype(scalar_field_entity), scalar_field))
1792     defshapemap.LookupTable = lookup_table
1793
1794     # Add scalar bar
1795     add_scalar_bar(field_name, nb_components,
1796                    vector_mode, lookup_table, time_value)
1797
1798     return defshapemap
1799
1800
1801 def Plot3DOnField(proxy, entity, field_name, timestamp_nb,
1802                   orientation=Orientation.AUTO,
1803                   angle1=0, angle2=0,
1804                   position=0.5, is_relative=True,
1805                   scale_factor=None,
1806                   is_contour=False, nb_contours=32,
1807                   vector_mode='Magnitude'):
1808     """Creates Plot 3D presentation on the given field.
1809
1810     Arguments:
1811       proxy: the pipeline object, containig data
1812       entity: the entity type from PrsTypeEnum
1813       field_name: the field name
1814       timestamp_nb: the number of time step (1, 2, ...)
1815       orientation: the cut plane plane orientation in 3D space, if
1816       the input is planar - will not be taken into account
1817       angle1: rotation of the cut plane in 3d space around the first axis
1818       of the selected orientation (X axis for XY, Y axis for YZ,
1819       Z axis for ZX).
1820       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1821       angle2: rotation of the cut plane in 3d space around the second axis
1822       of the selected orientation. Acceptable range: [-45, 45].
1823       position: position of the cut plane in the object (ranging from 0 to 1).
1824       The value 0.5 corresponds to cutting by halves.
1825       is_relative: defines if the cut plane position is relative or absolute
1826       scale_factor: deformation scale factor
1827       is_contour: if True - Plot 3D will be represented with a set of contours,
1828       otherwise - Plot 3D will be represented with a smooth surface
1829       nb_contours: number of contours, applied if is_contour is True
1830       vector_mode: the mode of transformation of vector values
1831       into scalar values, applicable only if the field contains vector values.
1832       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1833
1834     Returns:
1835       Plot 3D as representation object.
1836
1837     """
1838     proxy.UpdatePipeline()
1839     # We don't need mesh parts with no data on them
1840     if entity == EntityType.NODE:
1841         select_cells_with_data(proxy, on_points=[field_name])
1842     else:
1843         select_cells_with_data(proxy, on_cells=[field_name])
1844
1845     # Check vector mode
1846     nb_components = get_nb_components(proxy, entity, field_name)
1847     check_vector_mode(vector_mode, nb_components)
1848
1849     # Get time value
1850     time_value = get_time(proxy, timestamp_nb)
1851
1852     # Set timestamp
1853     pvs.GetRenderView().ViewTime = time_value
1854     pvs.UpdatePipeline(time_value, proxy)
1855
1856     # Do merge
1857     merge_blocks = pvs.MergeBlocks(proxy)
1858     merge_blocks.UpdatePipeline()
1859
1860     poly_data = None
1861
1862     # Cutting plane
1863
1864     # Define orientation if necessary (auto mode)
1865     plane_orientation = orientation
1866     if (orientation == Orientation.AUTO):
1867         plane_orientation = get_orientation(proxy)
1868
1869     # Get cutting plane normal
1870     normal = None
1871
1872     if (not is_planar_input(proxy)):
1873         normal = get_normal_by_orientation(plane_orientation,
1874                                            radians(angle1), radians(angle2))
1875
1876         # Create slice filter
1877         slice_filter = pvs.Slice(merge_blocks)
1878         slice_filter.SliceType = "Plane"
1879
1880         # Set cutting plane normal
1881         slice_filter.SliceType.Normal = normal
1882
1883         # Set cutting plane position
1884         if (is_relative):
1885             base_position = get_positions(1, normal,
1886                                           get_bounds(proxy), position)
1887             slice_filter.SliceOffsetValues = base_position
1888         else:
1889             slice_filter.SliceOffsetValues = position
1890
1891         slice_filter.UpdatePipeline()
1892         poly_data = slice_filter
1893     else:
1894         normal = get_normal_by_orientation(plane_orientation, 0, 0)
1895
1896     use_normal = 0
1897     # Geometry filter
1898     if not poly_data or poly_data.GetDataInformation().GetNumberOfCells() == 0:
1899         geometry_filter = pvs.GeometryFilter(merge_blocks)
1900         poly_data = geometry_filter
1901         use_normal = 1  # TODO(MZN): workaround
1902
1903     warp_scalar = None
1904     plot3d = None
1905     source = poly_data
1906
1907     if is_data_on_cells(poly_data, field_name):
1908         # Cell data to point data
1909         cell_to_point = pvs.CellDatatoPointData(poly_data)
1910         cell_to_point.PassCellData = 1
1911         source = cell_to_point
1912
1913     scalars = ['POINTS', field_name]
1914
1915     # Transform vector array to scalar array if necessary
1916     if (nb_components > 1):
1917         calc = get_calc_magnitude(source, EntityType.NODE, field_name)
1918         scalars = ['POINTS', calc.ResultArrayName]
1919         source = calc
1920
1921     # Warp by scalar
1922     warp_scalar = pvs.WarpByScalar(source)
1923     warp_scalar.Scalars = scalars
1924     warp_scalar.Normal = normal
1925     warp_scalar.UseNormal = use_normal
1926     if scale_factor is not None:
1927         warp_scalar.ScaleFactor = scale_factor
1928     else:
1929         def_scale = get_default_scale(PrsTypeEnum.PLOT3D,
1930                                       proxy, entity, field_name)
1931         warp_scalar.ScaleFactor = def_scale
1932
1933     warp_scalar.UpdatePipeline()
1934     source = warp_scalar
1935
1936     if (is_contour):
1937         # Contours
1938         contour = pvs.Contour(warp_scalar)
1939         contour.PointMergeMethod = "Uniform Binning"
1940         contour.ContourBy = ['POINTS', field_name]
1941         scalar_range = get_data_range(proxy, entity,
1942                                       field_name, vector_mode)
1943         contour.Isosurfaces = get_contours(scalar_range, nb_contours)
1944         contour.UpdatePipeline()
1945         source = contour
1946
1947     # Get Plot 3D representation object
1948     plot3d = pvs.GetRepresentation(source)
1949
1950     # Get lookup table
1951     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1952
1953     # Set field range if necessary
1954     data_range = get_data_range(proxy, entity,
1955                                 field_name, vector_mode)
1956     if hasattr(lookup_table,"LockDataRange"):
1957         lookup_table.LockDataRange = 1
1958     elif hasattr(lookup_table,"LockScalarRange"):
1959         lookup_table.LockScalarRange = 1
1960     else:
1961         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1962
1963     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1964
1965     # Set properties
1966     pvs.ColorBy(plot3d, (EntityType.get_pvtype(entity), field_name))
1967     plot3d.LookupTable = lookup_table
1968
1969     # Add scalar bar
1970     add_scalar_bar(field_name, nb_components,
1971                    vector_mode, lookup_table, time_value)
1972
1973     return plot3d
1974
1975
1976 def IsoSurfacesOnField(proxy, entity, field_name, timestamp_nb,
1977                        custom_range=None, nb_surfaces=10,
1978                        is_colored=True, color=None, vector_mode='Magnitude'):
1979     """Creates Iso Surfaces presentation on the given field.
1980
1981     Arguments:
1982       proxy: the pipeline object, containig data
1983       entity: the entity type from PrsTypeEnum
1984       field_name: the field name
1985       timestamp_nb: the number of time step (1, 2, ...)
1986       custom_range: scalar range, if undefined the source range will be applied
1987       nb_surfaces: number of surfaces, which will be generated
1988       is_colored: this option allows to color the presentation according to
1989       the corresponding data array values. If False - the presentation will
1990       be one-coloured.
1991       color: defines the presentation color as [R, G, B] triple. Taken into
1992       account only if is_colored is False.
1993       vector_mode: the mode of transformation of vector values
1994       into scalar values, applicable only if the field contains vector values.
1995       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1996
1997     Returns:
1998       Iso Surfaces as representation object.
1999
2000     """
2001     proxy.UpdatePipeline()
2002     # We don't need mesh parts with no data on them
2003     if entity == EntityType.NODE:
2004         select_cells_with_data(proxy, on_points=[field_name])
2005     else:
2006         select_cells_with_data(proxy, on_cells=[field_name])
2007
2008     # Check vector mode
2009     nb_components = get_nb_components(proxy, entity, field_name)
2010     check_vector_mode(vector_mode, nb_components)
2011
2012     # Get time value
2013     time_value = get_time(proxy, timestamp_nb)
2014
2015     # Set timestamp
2016     pvs.GetRenderView().ViewTime = time_value
2017     pvs.UpdatePipeline(time_value, proxy)
2018
2019     # Do merge
2020     source = pvs.MergeBlocks(proxy)
2021     pvs.UpdatePipeline()
2022
2023     # Transform cell data into point data if necessary
2024     if is_data_on_cells(proxy, field_name):
2025         cell_to_point = pvs.CellDatatoPointData(source)
2026         cell_to_point.PassCellData = 1
2027         source = cell_to_point
2028
2029     contour_by = ['POINTS', field_name]
2030
2031     # Transform vector array to scalar array if necessary
2032     if (nb_components > 1):
2033         calc = get_calc_magnitude(source, EntityType.NODE, field_name)
2034         contour_by = ['POINTS', calc.ResultArrayName]
2035         source = calc
2036
2037     # Contour filter settings
2038     contour = pvs.Contour(source)
2039     contour.ComputeScalars = 1
2040     contour.ContourBy = contour_by
2041
2042     # Specify the range
2043     scalar_range = custom_range
2044     if (scalar_range is None):
2045         scalar_range = get_data_range(proxy, entity,
2046                                       field_name, cut_off=True)
2047
2048     # Get contour values for the range
2049     surfaces = get_contours(scalar_range, nb_surfaces)
2050
2051     # Set contour values
2052     contour.Isosurfaces = surfaces
2053
2054     # Get Iso Surfaces representation object
2055     isosurfaces = pvs.GetRepresentation(contour)
2056
2057     # Get lookup table
2058     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2059
2060     # Set field range if necessary
2061     data_range = get_data_range(proxy, entity,
2062                                 field_name, vector_mode)
2063     if hasattr(lookup_table,"LockDataRange"):
2064         lookup_table.LockDataRange = 1
2065     elif hasattr(lookup_table,"LockScalarRange"):
2066         lookup_table.LockScalarRange = 1
2067     else:
2068         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
2069
2070     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2071
2072     # Set display properties
2073     if (is_colored):
2074         pvs.ColorBy(isosurfaces, (EntityType.get_pvtype(entity), field_name))
2075     else:
2076         pvs.ColorBy(isosurfaces, (EntityType.get_pvtype(entity), None))
2077         if color:
2078             isosurfaces.DiffuseColor = color
2079     isosurfaces.LookupTable = lookup_table
2080
2081     # Add scalar bar
2082     add_scalar_bar(field_name, nb_components,
2083                    vector_mode, lookup_table, time_value)
2084
2085     return isosurfaces
2086
2087
2088 def GaussPointsOnField(proxy, entity, field_name,
2089                        timestamp_nb,
2090                        is_deformed=True, scale_factor=None,
2091                        is_colored=True, color=None,
2092                        primitive=GaussType.SPRITE,
2093                        is_proportional=True,
2094                        max_pixel_size=256,
2095                        multiplier=None, vector_mode='Magnitude'):
2096     """Creates Gauss Points on the given field.
2097
2098     Arguments:
2099
2100     proxy: the pipeline object, containig data
2101     entity: the field entity type from PrsTypeEnum
2102     field_name: the field name
2103     timestamp_nb: the number of time step (1, 2, ...)
2104     is_deformed: defines whether the Gauss Points will be deformed or not
2105     scale_factor -- the scale factor for deformation. Will be taken into
2106     account only if is_deformed is True.
2107     If not passed by user, default scale will be computed.
2108     is_colored -- defines whether the Gauss Points will be multicolored,
2109     using the corresponding data values
2110     color: defines the presentation color as [R, G, B] triple. Taken into
2111     account only if is_colored is False.
2112     primitive: primitive type from GaussType
2113     is_proportional: if True, the size of primitives will depends on
2114     the gauss point value
2115     max_pixel_size: the maximum sizr of the Gauss Points primitive in pixels
2116     multiplier: coefficient between data values and the size of primitives
2117     If not passed by user, default scale will be computed.
2118     vector_mode: the mode of transformation of vector values into
2119     scalar values, applicable only if the field contains vector values.
2120     Possible modes: 'Magnitude' - vector module;
2121     'X', 'Y', 'Z' - vector components.
2122
2123     Returns:
2124       Gauss Points as representation object.
2125
2126     """
2127     proxy.UpdatePipeline()
2128     # We don't need mesh parts with no data on them
2129     on_gauss = select_cells_with_data(proxy, on_gauss=[field_name])
2130     if not on_gauss:
2131         if entity == EntityType.NODE:
2132             select_cells_with_data(proxy, on_points=[field_name])
2133         else:
2134             select_cells_with_data(proxy, on_cells=[field_name])
2135
2136     # Check vector mode
2137     nb_components = get_nb_components(proxy, entity, field_name)
2138     check_vector_mode(vector_mode, nb_components)
2139
2140     # Get time value
2141     time_value = get_time(proxy, timestamp_nb)
2142
2143     # Set timestamp
2144     pvs.GetRenderView().ViewTime = time_value
2145     pvs.UpdatePipeline(time_value, proxy)
2146
2147     source = proxy
2148
2149     # If no quadrature point array is passed, use cell centers
2150     if on_gauss:
2151         generate_qp = pvs.GenerateQuadraturePoints(source)
2152         generate_qp.QuadratureSchemeDef = ['CELLS', 'ELGA@0']
2153         source = generate_qp
2154     else:
2155         # Cell centers
2156         cell_centers = pvs.CellCenters(source)
2157         cell_centers.VertexCells = 1
2158         source = cell_centers
2159
2160     source.UpdatePipeline()
2161
2162     # Check if deformation enabled
2163     if is_deformed and nb_components > 1:
2164         vector_array = field_name
2165         # If the given vector array has only 2 components, add the third one
2166         if nb_components == 2:
2167             calc = get_add_component_calc(source, EntityType.NODE, field_name)
2168             vector_array = calc.ResultArrayName
2169             source = calc
2170
2171         # Warp by vector
2172         warp_vector = pvs.WarpByVector(source)
2173         warp_vector.Vectors = [vector_array]
2174         if scale_factor is not None:
2175             warp_vector.ScaleFactor = scale_factor
2176         else:
2177             def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE, proxy,
2178                                           entity, field_name)
2179             warp_vector.ScaleFactor = def_scale
2180         warp_vector.UpdatePipeline()
2181         source = warp_vector
2182
2183     # Get Gauss Points representation object
2184     gausspnt = pvs.GetRepresentation(source)
2185
2186     # Get lookup table
2187     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2188
2189     # Set field range if necessary
2190     data_range = get_data_range(proxy, entity,
2191                                 field_name, vector_mode)
2192     if hasattr(lookup_table,"LockDataRange"):
2193         lookup_table.LockDataRange = 1
2194     elif hasattr(lookup_table,"LockScalarRange"):
2195         lookup_table.LockScalarRange = 1
2196     else:
2197         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
2198
2199     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2200
2201     # Set display properties
2202     if is_colored:
2203         pvs.ColorBy(gausspnt, (EntityType.get_pvtype(entity), field_name))
2204     else:
2205         pvs.ColorBy(gausspnt, (EntityType.get_pvtype(entity), None))
2206         if color:
2207             gausspnt.DiffuseColor = color
2208
2209     gausspnt.LookupTable = lookup_table
2210
2211     # Add scalar bar
2212     add_scalar_bar(field_name, nb_components,
2213                    vector_mode, lookup_table, time_value)
2214
2215     # Set point sprite representation
2216     gausspnt.Representation = 'Point Sprite'
2217
2218     # Point sprite settings
2219     gausspnt.InterpolateScalarsBeforeMapping = 0
2220     gausspnt.MaxPixelSize = max_pixel_size
2221
2222     # Render mode
2223     gausspnt.RenderMode = GaussType.get_mode(primitive)
2224
2225     #if primitive == GaussType.SPRITE:
2226         # Set texture
2227         # TODO(MZN): replace with pvsimple high-level interface
2228     #    texture = sm.CreateProxy("textures", "SpriteTexture")
2229     #    alphamprop = texture.GetProperty("AlphaMethod")
2230     #    alphamprop.SetElement(0, 2)  # Clamp
2231     #    alphatprop = texture.GetProperty("AlphaThreshold")
2232     #    alphatprop.SetElement(0, 63)
2233     #    maxprop = texture.GetProperty("Maximum")
2234     #    maxprop.SetElement(0, 255)
2235     #    texture.UpdateVTKObjects()
2236
2237     #    gausspnt.Texture = texture
2238         #gausspnt.Texture.AlphaMethod = 'Clamp'
2239         #gausspnt.Texture.AlphaThreshold = 63
2240         #gausspnt.Texture.Maximum= 255
2241
2242     # Proportional radius
2243     gausspnt.RadiusUseScalarRange = 0
2244     gausspnt.RadiusIsProportional = 0
2245
2246     if is_proportional:
2247         mult = multiplier
2248         if mult is None and data_range[1] != 0:
2249             mult = abs(0.1 / data_range[1])
2250
2251         gausspnt.RadiusScalarRange = data_range
2252         gausspnt.RadiusTransferFunctionEnabled = 1
2253         gausspnt.RadiusMode = 'Scalar'
2254         gausspnt.RadiusArray = ['POINTS', field_name]
2255         if nb_components > 1:
2256             v_comp = get_vector_component(vector_mode)
2257             gausspnt.RadiusVectorComponent = v_comp
2258         gausspnt.RadiusTransferFunctionMode = 'Table'
2259         gausspnt.RadiusScalarRange = data_range
2260         gausspnt.RadiusUseScalarRange = 1
2261         if mult is not None:
2262             gausspnt.RadiusIsProportional = 1
2263             gausspnt.RadiusProportionalFactor = mult
2264     else:
2265         gausspnt.RadiusTransferFunctionEnabled = 0
2266         gausspnt.RadiusMode = 'Constant'
2267         gausspnt.RadiusArray = ['POINTS', 'Constant Radius']
2268
2269     return gausspnt
2270
2271 def GaussPointsOnField1(proxy, entity, field_name,
2272                         timestamp_nb,
2273                         is_colored=True, color=None,
2274                         primitive=GaussType.SPHERE,
2275                         is_proportional=True,
2276                         max_pixel_size=256,
2277                         multiplier=None,
2278                         vector_mode='Magnitude'):
2279     """Creates Gauss Points on the given field. Use GaussPoints() Paraview interface.
2280
2281     Arguments:
2282     proxy: the pipeline object, containig data
2283     entity: the field entity type from PrsTypeEnum
2284     field_name: the field name
2285     timestamp_nb: the number of time step (1, 2, ...)
2286     is_colored -- defines whether the Gauss Points will be multicolored,
2287     using the corresponding data values
2288     color: defines the presentation color as [R, G, B] triple. Taken into
2289     account only if is_colored is False.
2290     primitive: primitive type from GaussType
2291     is_proportional: if True, the size of primitives will depends on
2292     the gauss point value
2293     max_pixel_size: the maximum sizr of the Gauss Points primitive in pixels
2294     multiplier: coefficient between data values and the size of primitives
2295     If not passed by user, default scale will be computed.
2296     vector_mode: the mode of transformation of vector values into
2297     scalar values, applicable only if the field contains vector values.
2298     Possible modes: 'Magnitude' - vector module;
2299     'X', 'Y', 'Z' - vector components.
2300
2301     Returns:
2302       Gauss Points as representation object.
2303
2304     """
2305     proxy.UpdatePipeline()
2306     select_cells_with_data(proxy, on_gauss=[field_name])
2307
2308     nb_components = get_nb_components(proxy, entity, field_name)
2309
2310     # Get time value
2311     time_value = get_time(proxy, timestamp_nb)
2312
2313     # Set timestamp
2314     pvs.GetRenderView().ViewTime = time_value
2315     proxy.UpdatePipeline(time=time_value)
2316
2317     # Create Gauss Points object
2318     source = pvs.ELGAfieldToPointSprite(proxy)
2319     source.UpdatePipeline()
2320
2321     # Get Gauss Points representation object
2322     gausspnt = pvs.GetRepresentation(source)
2323
2324     # Get lookup table
2325     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2326
2327     # Set field range if necessary
2328     data_range = get_data_range(proxy, entity,
2329                                 field_name, vector_mode)
2330     if hasattr(lookup_table,"LockDataRange"):
2331         lookup_table.LockDataRange = 1
2332     elif hasattr(lookup_table,"LockScalarRange"):
2333         lookup_table.LockScalarRange = 1
2334     else:
2335         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
2336
2337     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2338
2339     # Set display properties
2340     if is_colored:
2341         pvs.ColorBy(gausspnt, (EntityType.get_pvtype(entity), field_name))
2342     else:
2343         pvs.ColorBy(gausspnt, (EntityType.get_pvtype(entity), None))
2344         if color:
2345             gausspnt.DiffuseColor = color
2346
2347     gausspnt.LookupTable = lookup_table
2348
2349     # Add scalar bar
2350     add_scalar_bar(field_name, nb_components,
2351                    vector_mode, lookup_table, time_value)
2352
2353     # Set point sprite representation
2354     gausspnt.Representation = 'Point Sprite'
2355
2356     # Point sprite settings
2357     gausspnt.InterpolateScalarsBeforeMapping = 0
2358     gausspnt.MaxPixelSize = max_pixel_size
2359
2360     # Render mode
2361     gausspnt.RenderMode = GaussType.get_mode(primitive)
2362
2363     #if primitive == GaussType.SPRITE:
2364         # Set texture
2365         # TODO(MZN): replace with pvsimple high-level interface
2366     #    texture = sm.CreateProxy("textures", "SpriteTexture")
2367     #    alphamprop = texture.GetProperty("AlphaMethod")
2368     #    alphamprop.SetElement(0, 2)  # Clamp
2369     #    alphatprop = texture.GetProperty("AlphaThreshold")
2370     #    alphatprop.SetElement(0, 63)
2371     #    maxprop = texture.GetProperty("Maximum")
2372     #    maxprop.SetElement(0, 255)
2373     #    texture.UpdateVTKObjects()
2374
2375     #    gausspnt.Texture = texture
2376         #gausspnt.Texture.AlphaMethod = 'Clamp'
2377         #gausspnt.Texture.AlphaThreshold = 63
2378         #gausspnt.Texture.Maximum= 255
2379
2380     # Proportional radius
2381     gausspnt.RadiusUseScalarRange = 0
2382     gausspnt.RadiusIsProportional = 0
2383
2384     if is_proportional:
2385         mult = multiplier
2386         if mult is None and data_range[1] != 0:
2387             mult = abs(0.1 / data_range[1])
2388
2389         gausspnt.RadiusScalarRange = data_range
2390         gausspnt.RadiusTransferFunctionEnabled = 1
2391         gausspnt.RadiusMode = 'Scalar'
2392         gausspnt.RadiusArray = ['POINTS', field_name]
2393         if nb_components > 1:
2394             v_comp = get_vector_component(vector_mode)
2395             gausspnt.RadiusVectorComponent = v_comp
2396         gausspnt.RadiusTransferFunctionMode = 'Table'
2397         gausspnt.RadiusScalarRange = data_range
2398         gausspnt.RadiusUseScalarRange = 1
2399         if mult is not None:
2400             gausspnt.RadiusIsProportional = 1
2401             gausspnt.RadiusProportionalFactor = mult
2402     else:
2403         gausspnt.RadiusTransferFunctionEnabled = 0
2404         gausspnt.RadiusMode = 'Constant'
2405         gausspnt.RadiusArray = ['POINTS', 'Constant Radius']
2406
2407     return gausspnt
2408
2409 def StreamLinesOnField(proxy, entity, field_name, timestamp_nb,
2410                        direction='BOTH', is_colored=False, color=None,
2411                        vector_mode='Magnitude'):
2412     """Creates Stream Lines presentation on the given field.
2413
2414     Arguments:
2415       proxy: the pipeline object, containig data
2416       entity: the entity type from PrsTypeEnum
2417       field_name: the field name
2418       timestamp_nb: the number of time step (1, 2, ...)
2419       direction: the stream lines direction ('FORWARD', 'BACKWARD' or 'BOTH')
2420       is_colored: this option allows to color the presentation according to
2421       the corresponding data values. If False - the presentation will
2422       be one-coloured.
2423       color: defines the presentation color as [R, G, B] triple. Taken into
2424       account only if is_colored is False.
2425       vector_mode: the mode of transformation of vector values
2426       into scalar values, applicable only if the field contains vector values.
2427       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
2428
2429     Returns:
2430       Stream Lines as representation object.
2431
2432     """
2433     proxy.UpdatePipeline()
2434     # We don't need mesh parts with no data on them
2435     if entity == EntityType.NODE:
2436         select_cells_with_data(proxy, on_points=[field_name])
2437     else:
2438         select_cells_with_data(proxy, on_cells=[field_name])
2439
2440     # Check vector mode
2441     nb_components = get_nb_components(proxy, entity, field_name)
2442     check_vector_mode(vector_mode, nb_components)
2443
2444     # Get time value
2445     time_value = get_time(proxy, timestamp_nb)
2446
2447     # Set timestamp
2448     pvs.GetRenderView().ViewTime = time_value
2449     pvs.UpdatePipeline(time_value, proxy)
2450
2451     # Do merge
2452     source = pvs.MergeBlocks(proxy)
2453     pvs.UpdatePipeline()
2454
2455     # Cell data to point data
2456     if is_data_on_cells(proxy, field_name):
2457         cell_to_point = pvs.CellDatatoPointData(source)
2458         cell_to_point.PassCellData = 1
2459         pvs.UpdatePipeline()
2460         source = cell_to_point
2461
2462     vector_array = field_name
2463     # If the given vector array has only 2 components, add the third one
2464     if nb_components == 2:
2465         calc = get_add_component_calc(source, EntityType.NODE, field_name)
2466         vector_array = calc.ResultArrayName
2467         pvs.UpdatePipeline()
2468         source = calc
2469
2470     # Stream Tracer
2471     stream = pvs.StreamTracer(source)
2472     stream.SeedType = "Point Source"
2473     stream.Vectors = ['POINTS', vector_array]
2474     stream.IntegrationDirection = direction
2475     stream.IntegratorType = 'Runge-Kutta 2'
2476     stream.SeedType = 'High Resolution Line Source'
2477     stream.UpdatePipeline()
2478
2479     # Get Stream Lines representation object
2480     if is_empty(stream):
2481         return None
2482     streamlines = pvs.GetRepresentation(stream)
2483
2484     # Get lookup table
2485     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2486
2487     # Set field range if necessary
2488     data_range = get_data_range(proxy, entity,
2489                                 field_name, vector_mode)
2490     if hasattr(lookup_table,"LockDataRange"):
2491         lookup_table.LockDataRange = 1
2492     elif hasattr(lookup_table,"LockScalarRange"):
2493         lookup_table.LockScalarRange = 1
2494     else:
2495         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
2496
2497     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2498
2499     # Set properties
2500     if is_colored:
2501         pvs.ColorBy(streamlines, (EntityType.get_pvtype(entity), field_name))
2502     else:
2503         pvs.ColorBy(streamlines, (EntityType.get_pvtype(entity), None))
2504         if color:
2505             streamlines.DiffuseColor = color
2506
2507     streamlines.LookupTable = lookup_table
2508
2509     # Add scalar bar
2510     add_scalar_bar(field_name, nb_components,
2511                    vector_mode, lookup_table, time_value)
2512
2513     return streamlines
2514
2515
2516 def MeshOnEntity(proxy, mesh_name, entity):
2517     """Creates submesh of the entity type for the mesh.
2518
2519     Arguments:
2520       proxy -- the pipeline object, containig data
2521       mesh_name -- the full or short name of mesh field
2522
2523     Returns:
2524       Submesh as representation object of the given source.
2525
2526     """
2527     proxy.UpdatePipeline()
2528     mesh_full_name = None
2529     aList = mesh_name.split('/')
2530     if len(aList) >= 2:
2531         mesh_full_name = mesh_name
2532     else:
2533         mesh_full_name = find_mesh_full_name(proxy, mesh_name)
2534     if not mesh_full_name:
2535         raise RuntimeError("The given mesh name was not found")
2536     # Select only the given mesh
2537     proxy.AllArrays = [mesh_full_name]
2538     proxy.UpdatePipeline()
2539
2540     # Get representation object if the submesh is not empty
2541     prs = None
2542     if (proxy.GetDataInformation().GetNumberOfPoints() or
2543         proxy.GetDataInformation().GetNumberOfCells()):
2544         my_view = pvs.GetRenderView()
2545         prs = pvs.GetRepresentation(proxy, view=my_view)
2546         prs.ColorArrayName = (None, '')
2547
2548     return prs
2549
2550
2551 def MeshOnGroup(proxy, extrGroups, group_name):
2552     """Creates submesh on the group.
2553
2554     Arguments:
2555       proxy -- the pipeline object, containig data
2556       group_name -- the full group name
2557       extrGroups -- all extracted groups object
2558
2559     Returns:
2560       Representation object of the given source with single group
2561       selected.
2562
2563     """
2564     proxy.UpdatePipeline()
2565     # Deselect all groups
2566     extrGroups.AllGroups = []
2567     extrGroups.UpdatePipelineInformation()
2568     # Select only the group with the given name
2569     extrGroups.AllGroups = [group_name]
2570     extrGroups.UpdatePipelineInformation()
2571
2572     # Get representation object if the submesh is not empty
2573     prs = None
2574
2575     # Check if the group was set
2576     if len(extrGroups.AllGroups) == 1 and \
2577        extrGroups.AllGroups[0] == group_name:
2578         # Check if the submesh is not empty
2579         nb_points = proxy.GetDataInformation().GetNumberOfPoints()
2580         nb_cells = proxy.GetDataInformation().GetNumberOfCells()
2581
2582         if nb_points or nb_cells:
2583 #            prs = pvs.GetRepresentation(proxy)
2584             prs = pvs.Show()
2585             prs.ColorArrayName = (None, '')
2586             display_only(prs)
2587
2588     return prs
2589
2590
2591 def CreatePrsForFile(file_name, prs_types,
2592                      picture_dir, picture_ext):
2593     """Build presentations of the given types for the file.
2594
2595     Build presentations for all fields on all timestamps.
2596
2597     Arguments:
2598       file_name: full path to the MED file
2599       prs_types: the list of presentation types to build
2600       picture_dir: the directory path for saving snapshots
2601       picture_ext: graphics files extension (determines file type)
2602
2603     """
2604     # Import MED file
2605     print("Import " + file_name.split(os.sep)[-1] + "...")
2606
2607     try:
2608         proxy = pvs.MEDReader(FileName=file_name)
2609         if proxy is None:
2610             print("FAILED")
2611         else:
2612             #proxy.UpdatePipeline()
2613             print("OK")
2614     except:
2615         print("FAILED")
2616     else:
2617         # Get view
2618         view = pvs.GetRenderView()
2619         time_value = get_time(proxy, 0)
2620         view.ViewTime = time_value
2621         pvs.UpdatePipeline(time=time_value, proxy=proxy)
2622
2623         # Create required presentations for the proxy
2624         CreatePrsForProxy(proxy, view, prs_types,
2625                           picture_dir, picture_ext)
2626
2627 def CreatePrsForProxy(proxy, view, prs_types, picture_dir, picture_ext):
2628     """Build presentations of the given types for all fields of the proxy.
2629
2630     Save snapshots in graphics files (type depends on the given extension).
2631     Stores the files in the given directory.
2632
2633     Arguments:
2634       proxy: the pipeline object, containig data
2635       view: the render view
2636       prs_types: the list of presentation types to build
2637       picture_dir: the directory path for saving snapshots
2638       picture_ext: graphics files extension (determines file type)
2639
2640     """
2641     proxy.UpdatePipeline()
2642     # List of the field names
2643     fields_info = proxy.GetProperty("FieldsTreeInfo")[::2]
2644
2645     # Add path separator to the end of picture path if necessery
2646     if not picture_dir.endswith(os.sep):
2647         picture_dir += os.sep
2648
2649     # Mesh Presentation
2650     if PrsTypeEnum.MESH in prs_types:
2651         # Iterate on meshes
2652         mesh_names = get_mesh_full_names(proxy)
2653         for mesh_name in mesh_names:
2654             # Build mesh field presentation
2655             print("Creating submesh for '" + get_field_short_name(mesh_name) + "' mesh... ")
2656             prs = MeshOnEntity(proxy, mesh_name, None)
2657             if prs is None:
2658                 print("FAILED")
2659                 continue
2660             else:
2661                 print("OK")
2662             # Construct image file name
2663             pic_name = picture_dir + get_field_short_name(mesh_name) + "." + picture_ext
2664
2665             # Show and dump the presentation into a graphics file
2666             process_prs_for_test(prs, view, pic_name, False)
2667
2668             # Create Mesh presentation. Build all groups.
2669             extGrp = pvs.ExtractGroup()
2670             extGrp.UpdatePipelineInformation()
2671             if if_possible(proxy, None, None, PrsTypeEnum.MESH, extGrp):
2672                 for group in get_group_names(extGrp):
2673                     print("Creating submesh on group " + get_group_short_name(group) + "... ")
2674                     prs = MeshOnGroup(proxy, extGrp, group)
2675                     if prs is None:
2676                         print("FAILED")
2677                         continue
2678                     else:
2679                         print("OK")
2680                     # Construct image file name
2681                     pic_name = picture_dir + get_group_short_name(group) + "." + picture_ext
2682
2683                     # Show and dump the presentation into a graphics file
2684                     process_prs_for_test(prs, view, pic_name, False)
2685
2686     # Presentations on fields
2687     for field in fields_info:
2688         field_name = get_field_short_name(field)
2689         # Ignore mesh presentation
2690         if field_name == get_field_mesh_name(field):
2691             continue
2692         field_entity = get_field_entity(field)
2693         # Select only the current field:
2694         # necessary for getting the right timestamps
2695         proxy.AllArrays = [field]
2696         proxy.UpdatePipeline()
2697
2698         # Get timestamps
2699         timestamps = proxy.TimestepValues.GetData()
2700
2701         for prs_type in prs_types:
2702             # Ignore mesh presentation
2703             if prs_type == PrsTypeEnum.MESH:
2704                 continue
2705
2706             # Get name of presentation type
2707             prs_name = PrsTypeEnum.get_name(prs_type)
2708
2709             # Build the presentation if possible
2710             possible = if_possible(proxy, field_name,
2711                                    field_entity, prs_type)
2712             if possible:
2713                 # Presentation type for graphics file name
2714                 f_prs_type = prs_name.replace(' ', '').upper()
2715
2716                 for timestamp_nb in range(1, len(timestamps) + 1):
2717                     time = timestamps[timestamp_nb - 1]
2718                     if (time == 0.0):
2719                         scalar_range = get_data_range(proxy, field_entity,
2720                                                       field_name, cut_off=True)
2721                         # exclude time stamps with null lenght of scalar range
2722                         if (scalar_range[0] == scalar_range[1]):
2723                             continue
2724                     print("Creating " + prs_name + " on " + field_name + ", time = " + str(time) + "... ")
2725                     try:
2726                         prs = create_prs(prs_type, proxy,
2727                                          field_entity, field_name, timestamp_nb)
2728                     except ValueError:
2729                         """ This exception comes from get_nb_components(...) function.
2730                             The reason of exception is an implementation of MEDReader
2731                             activating the first leaf when reading MED file (refer to
2732                             MEDFileFieldRepresentationTree::activateTheFirst() and
2733                             MEDFileFieldRepresentationTree::getTheSingleActivated(...) methods).
2734                         """
2735                         print("ValueError exception is catched")
2736                         continue
2737                     if prs is None:
2738                         print("FAILED")
2739                         continue
2740                     else:
2741                         print("OK")
2742
2743                     # Construct image file name
2744                     pic_name = picture_dir + field_name + "_" + str(time) + "_" + f_prs_type + "." + picture_ext
2745
2746                     # Show and dump the presentation into a graphics file
2747                     process_prs_for_test(prs, view, pic_name)
2748     return
2749
2750
2751 def delete_pv_object(obj):
2752   # There is a bug when repeating CreateRenderView/Delete calls
2753   # Here is a workaround proposed by KW (#10744)
2754   import gc
2755   del obj
2756   gc.collect()