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