]> SALOME platform Git repositories - modules/paravis.git/blob - src/PV_SWIG/VTKWrapping/presentations.py
Salome HOME
Merge branch 'akl/tests_update'
[modules/paravis.git] / src / PV_SWIG / VTKWrapping / presentations.py
1 # Copyright (C) 2010-2014  CEA/DEN, EDF R&D
2 #
3 # This library is free software; you can redistribute it and/or
4 # modify it under the terms of the GNU Lesser General Public
5 # License as published by the Free Software Foundation; either
6 # version 2.1 of the License, or (at your option) any later version.
7 #
8 # This library is distributed in the hope that it will be useful,
9 # but WITHOUT ANY WARRANTY; without even the implied warranty of
10 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 # Lesser General Public License for more details.
12 #
13 # You should have received a copy of the GNU Lesser General Public
14 # License along with this library; if not, write to the Free Software
15 # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 #
17 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 #
19
20 """
21 This module is intended to provide Python API for building presentations
22 typical for Post-Pro module (Scalar Map, Deformed Shape, Vectors, etc.)
23 """
24
25
26 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=[], on_cells=[], on_gauss=[]):
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         separator = proxy.GetProperty("Separator").GetData()
776
777         fields_info = proxy.GetProperty("FieldsTreeInfo")[::2]
778         arr_name_with_dis=[elt.split("/")[-1] for elt in fields_info]
779
780         proxy.AllArrays = []
781         
782         fields = []
783         for name in on_gauss:
784             fields.append(name+separator+'GAUSS')
785         for name in on_cells:
786             fields.append(name+separator+'P0')
787         for name in on_points:
788             fields.append(name+separator+'P1')
789
790         field_list = []
791         for name in fields:
792             if arr_name_with_dis.count(name) > 0:
793                 index = arr_name_with_dis.index(name)
794                 field_list.append(fields_info[index])
795                 
796         proxy.AllArrays = field_list
797         proxy.UpdatePipeline()
798         return len(field_list) != 0
799
800     # TODO: VTN. Looks like this code is out of date.
801     
802     #all_cell_types = proxy.CellTypes.Available
803     all_cell_types = proxy.Entity.Available
804     all_arrays = list(proxy.CellArrays.GetData())
805     all_arrays.extend(proxy.PointArrays.GetData())
806
807     if not all_arrays:
808         file_name = proxy.FileName.split(os.sep)[-1]
809         print "Warning: " + file_name + " doesn't contain any data array."
810
811     # List of cell types to be selected
812     cell_types_on = []
813
814     for cell_type in all_cell_types:
815         #proxy.CellTypes = [cell_type]
816         proxy.Entity = [cell_type]
817         proxy.UpdatePipeline()
818
819         cell_arrays = proxy.GetCellDataInformation().keys()
820         point_arrays = proxy.GetPointDataInformation().keys()
821
822         if on_points or on_cells:
823             if on_points is None:
824                 on_points = []
825             if on_cells is None:
826                 on_cells = []
827
828             if (all(array in cell_arrays for array in on_cells) and
829                 all(array in point_arrays for array in on_points)):
830                 # Add cell type to the list
831                 cell_types_on.append(cell_type)
832         else:
833             in_arrays = lambda array: ((array in cell_arrays) or
834                                        (array in point_arrays))
835             if any(in_arrays(array) for array in all_arrays):
836                 cell_types_on.append(cell_type)
837
838     # Select cell types
839     #proxy.CellTypes = cell_types_on
840     proxy.Entity = cell_types_on
841     proxy.UpdatePipeline()
842
843 def if_possible(proxy, field_name, entity, prs_type):
844     """Check if the presentation creation is possible on the given field."""
845     result = True
846     if (prs_type == PrsTypeEnum.DEFORMEDSHAPE or
847         prs_type == PrsTypeEnum.DEFORMEDSHAPESCALARMAP or
848         prs_type == PrsTypeEnum.VECTORS or
849         prs_type == PrsTypeEnum.STREAMLINES):
850         nb_comp = get_nb_components(proxy, entity, field_name)
851         result = (nb_comp > 1)
852     elif (prs_type == PrsTypeEnum.GAUSSPOINTS):
853         result = (entity == EntityType.CELL or
854                   field_name in proxy.QuadraturePointArrays.Available)
855     elif (prs_type == PrsTypeEnum.MESH):
856         result = len(get_group_names(proxy, field_name, entity)) > 0
857
858     return result
859
860
861 def add_scalar_bar(field_name, nb_components,
862                    vector_mode, lookup_table, time_value):
863     """Add scalar bar with predefined properties."""
864     global _current_bar
865
866     # Construct bar title
867     title = "\n".join([field_name, str(time_value)])
868     if nb_components > 1:
869         title = "\n".join([title, vector_mode])
870
871     # Create scalar bar
872     scalar_bar = pvs.CreateScalarBar(Enabled=1)
873     scalar_bar.Orientation = 'Vertical'
874     scalar_bar.Title = title
875     scalar_bar.LookupTable = lookup_table
876
877     # Set default properties same as in Post-Pro
878     scalar_bar.NumberOfLabels = 5
879     scalar_bar.AutomaticLabelFormat = 0
880     scalar_bar.LabelFormat = '%-#6.6g'
881     # Title
882     scalar_bar.TitleFontFamily = 'Arial'
883     scalar_bar.TitleFontSize = 8
884     scalar_bar.TitleBold = 1
885     scalar_bar.TitleItalic = 1
886     scalar_bar.TitleShadow = 1
887     # Labels
888     scalar_bar.LabelFontFamily = 'Arial'
889     scalar_bar.LabelFontSize = 8
890     scalar_bar.LabelBold = 1
891     scalar_bar.LabelItalic = 1
892     scalar_bar.LabelShadow = 1
893
894     # Add the scalar bar to the view
895     pvs.GetRenderView().Representations.append(scalar_bar)
896
897     # Reassign the current bar
898     _current_bar = scalar_bar
899
900     return scalar_bar
901
902
903 def get_bar():
904     """Get current scalar bar."""
905     global _current_bar
906
907     return _current_bar
908
909
910 def get_lookup_table(field_name, nb_components, vector_mode='Magnitude'):
911     """Get lookup table for the given field."""
912     lookup_table = pvs.GetLookupTableForArray(field_name, nb_components)
913
914     if vector_mode == 'Magnitude':
915         lookup_table.VectorMode = vector_mode
916     elif vector_mode == 'X':
917         lookup_table.VectorMode = 'Component'
918         lookup_table.VectorComponent = 0
919     elif vector_mode == 'Y':
920         lookup_table.VectorMode = 'Component'
921         lookup_table.VectorComponent = 1
922     elif vector_mode == 'Z':
923         lookup_table.VectorMode = 'Component'
924         lookup_table.VectorComponent = 2
925     else:
926         raise ValueError("Incorrect vector mode: " + vector_mode)
927
928     lookup_table.Discretize = 0
929     lookup_table.ColorSpace = 'HSV'
930     lookup_table.LockScalarRange = 0
931
932     return lookup_table
933
934
935 def get_group_mesh_name(full_group_name):
936     """Return mesh name of the group by its full name."""
937     aList = full_group_name.split('/')
938     if len(aList) >= 2 :
939         group_name = full_group_name.split('/')[1]
940         return group_name
941
942
943 def get_group_entity(full_group_name):
944     """Return entity type of the group by its full name."""
945     aList = full_group_name.split('/')
946     if len(aList) >= 3 :
947         entity_name = full_group_name.split('/')[2]
948         entity = EntityType.get_type(entity_name)
949         return entity
950
951
952 def get_group_short_name(full_group_name):
953     """Return short name of the group by its full name."""
954     aList = full_group_name.split('/')
955     if len(aList) >= 4 :
956         short_name = full_group_name.split('/')[3]
957         return short_name
958
959
960 def get_mesh_names(proxy):
961     """Return all mesh names in the given proxy as a set."""
962     fields = proxy.GetProperty("FieldsTreeInfo")[::2]
963     mesh_names = set([get_field_mesh_name(item) for item in fields])
964     return mesh_names
965
966
967 def get_group_names(proxy, mesh_name, entity, wo_nogroups=False):
968     """Return full names of all groups of the given entity type
969     from the mesh with the given name as a list.
970     """
971     groups = proxy.Groups.Available
972
973     condition = lambda item: (get_group_mesh_name(item) == mesh_name and
974                               get_group_entity(item) == entity)
975     group_names = [item for item in groups if condition(item)]
976
977     if wo_nogroups:
978         # Remove "No_Group" group
979         not_no_group = lambda item: get_group_short_name(item) != "No_Group"
980         group_names = filter(not_no_group, group_names)
981
982     return group_names
983
984
985 def get_time(proxy, timestamp_nb):
986     """Get time value by timestamp number."""
987     # Check timestamp number
988     timestamps = []
989     
990     if (hasattr(proxy, 'TimestepValues')):
991         timestamps = proxy.TimestepValues.GetData()
992     elif (hasattr(proxy.Input, 'TimestepValues')):
993         timestamps = proxy.Input.TimestepValues.GetData()
994
995     length = len(timestamps)
996     if (timestamp_nb > 0 and (timestamp_nb - 1) not in xrange(length) ) or (timestamp_nb < 0 and -timestamp_nb > length):
997         raise ValueError("Timestamp number is out of range: " + str(timestamp_nb))
998
999     # Return time value
1000     if timestamp_nb > 0:
1001         return timestamps[timestamp_nb - 1]
1002     else:
1003         return timestamps[timestamp_nb]
1004
1005 def create_prs(prs_type, proxy, field_entity, field_name, timestamp_nb):
1006     """Auxiliary function.
1007
1008     Build presentation of the given type on the given field and
1009     timestamp number.
1010     Set the presentation properties like visu.CreatePrsForResult() do.
1011
1012     """
1013     prs = None
1014
1015     if prs_type == PrsTypeEnum.SCALARMAP:
1016         prs = ScalarMapOnField(proxy, field_entity, field_name, timestamp_nb)
1017     elif prs_type == PrsTypeEnum.CUTPLANES:
1018         prs = CutPlanesOnField(proxy, field_entity, field_name, timestamp_nb,
1019                                orientation=Orientation.ZX)
1020     elif prs_type == PrsTypeEnum.CUTLINES:
1021         prs = CutLinesOnField(proxy, field_entity, field_name, timestamp_nb,
1022                               orientation1=Orientation.XY,
1023                               orientation2=Orientation.ZX)
1024     elif prs_type == PrsTypeEnum.DEFORMEDSHAPE:
1025         prs = DeformedShapeOnField(proxy, field_entity,
1026                                    field_name, timestamp_nb)
1027     elif prs_type == PrsTypeEnum.DEFORMEDSHAPESCALARMAP:
1028         prs = DeformedShapeAndScalarMapOnField(proxy, field_entity,
1029                                                field_name, timestamp_nb)
1030     elif prs_type == PrsTypeEnum.VECTORS:
1031         prs = VectorsOnField(proxy, field_entity, field_name, timestamp_nb)
1032     elif prs_type == PrsTypeEnum.PLOT3D:
1033         prs = Plot3DOnField(proxy, field_entity, field_name, timestamp_nb)
1034     elif prs_type == PrsTypeEnum.ISOSURFACES:
1035         prs = IsoSurfacesOnField(proxy, field_entity, field_name, timestamp_nb)
1036     elif prs_type == PrsTypeEnum.GAUSSPOINTS:
1037         prs = GaussPointsOnField(proxy, field_entity, field_name, timestamp_nb)
1038     elif prs_type == PrsTypeEnum.STREAMLINES:
1039         prs = StreamLinesOnField(proxy, field_entity, field_name, timestamp_nb)
1040     else:
1041         raise ValueError("Unexistent presentation type.")
1042
1043     return prs
1044
1045
1046 # Functions for building Post-Pro presentations
1047 def ScalarMapOnField(proxy, entity, field_name, timestamp_nb,
1048                      vector_mode='Magnitude'):
1049     """Creates Scalar Map presentation on the given field.
1050
1051     Arguments:
1052       proxy: the pipeline object, containig data
1053       entity: the entity type from PrsTypeEnum
1054       field_name: the field name
1055       timestamp_nb: the number of time step (1, 2, ...)
1056       vector_mode: the mode of transformation of vector values
1057       into scalar values, applicable only if the field contains vector values.
1058       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1059
1060     Returns:
1061       Scalar Map as representation object.
1062
1063     """
1064     # We don't need mesh parts with no data on them
1065     if entity == EntityType.NODE:
1066         select_cells_with_data(proxy, on_points=[field_name])
1067     else:
1068         select_cells_with_data(proxy, on_cells=[field_name])
1069
1070     # Check vector mode
1071     nb_components = get_nb_components(proxy, entity, field_name)
1072     check_vector_mode(vector_mode, nb_components)
1073
1074     # Get time value
1075     time_value = get_time(proxy, timestamp_nb)
1076
1077     # Set timestamp
1078     pvs.GetRenderView().ViewTime = time_value
1079     pvs.UpdatePipeline(time_value, proxy)
1080
1081     # Get Scalar Map representation object
1082     scalarmap = pvs.GetRepresentation(proxy)
1083
1084     # Get lookup table
1085     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1086
1087     # Set field range if necessary
1088     data_range = get_data_range(proxy, entity,
1089                                 field_name, vector_mode)
1090     lookup_table.LockScalarRange = 1
1091     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1092     # Set properties
1093     scalarmap.ColorAttributeType = EntityType.get_pvtype(entity)
1094     scalarmap.ColorArrayName = field_name
1095     scalarmap.LookupTable = lookup_table
1096
1097     # Add scalar bar
1098     bar_title = field_name + ", " + str(time_value)
1099     if (nb_components > 1):
1100         bar_title += "\n" + vector_mode
1101     add_scalar_bar(field_name, nb_components, vector_mode,
1102                    lookup_table, time_value)
1103
1104     return scalarmap
1105
1106
1107 def CutPlanesOnField(proxy, entity, field_name, timestamp_nb,
1108                      nb_planes=10, orientation=Orientation.YZ,
1109                      angle1=0, angle2=0,
1110                      displacement=0.5, vector_mode='Magnitude'):
1111     """Creates Cut Planes presentation on the given field.
1112
1113     Arguments:
1114       proxy: the pipeline object, containig data
1115       entity: the entity type from PrsTypeEnum
1116       field_name: the field name
1117       timestamp_nb: the number of time step (1, 2, ...)
1118       nb_planes: number of cutting planes
1119       orientation: cutting planes orientation in 3D space
1120       angle1: rotation of the planes in 3d space around the first axis of the
1121       selected orientation (X axis for XY, Y axis for YZ, Z axis for ZX).
1122       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1123       angle2: rotation of the planes in 3d space around the second axis of the
1124       selected orientation. Acceptable range: [-45, 45].
1125       displacement: the displacement of the planes into one or another side
1126       vector_mode: the mode of transformation of vector values
1127       into scalar values, applicable only if the field contains vector values.
1128       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1129
1130     Returns:
1131       Cut Planes as representation object.
1132
1133     """
1134     if entity == EntityType.NODE:
1135         select_cells_with_data(proxy, on_points=[field_name])
1136     else:
1137         select_cells_with_data(proxy, on_cells=[field_name])
1138
1139     # Check vector mode
1140     nb_components = get_nb_components(proxy, entity, field_name)
1141     check_vector_mode(vector_mode, nb_components)
1142
1143     # Get time value
1144     time_value = get_time(proxy, timestamp_nb)
1145
1146     # Set timestamp
1147     pvs.GetRenderView().ViewTime = time_value
1148     pvs.UpdatePipeline(time_value, proxy)
1149
1150     # Create slice filter
1151     slice_filter = pvs.Slice(proxy)
1152     slice_filter.SliceType = "Plane"
1153
1154     # Set cut planes normal
1155     normal = get_normal_by_orientation(orientation,
1156                                        radians(angle1), radians(angle2))
1157     slice_filter.SliceType.Normal = normal
1158
1159     # Set cut planes positions
1160     positions = get_positions(nb_planes, normal,
1161                               get_bounds(proxy), displacement)
1162     slice_filter.SliceOffsetValues = positions
1163
1164     # Get Cut Planes representation object
1165     cut_planes = pvs.GetRepresentation(slice_filter)
1166
1167     # Get lookup table
1168     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1169
1170     # Set field range if necessary
1171     data_range = get_data_range(proxy, entity,
1172                                 field_name, vector_mode)
1173     lookup_table.LockScalarRange = 1
1174     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1175
1176     # Set properties
1177     cut_planes.ColorAttributeType = EntityType.get_pvtype(entity)
1178     cut_planes.ColorArrayName = field_name
1179     cut_planes.LookupTable = lookup_table
1180
1181     # Add scalar bar
1182     add_scalar_bar(field_name, nb_components,
1183                    vector_mode, lookup_table, time_value)
1184
1185     return cut_planes
1186
1187
1188 def CutLinesOnField(proxy, entity, field_name, timestamp_nb,
1189                     nb_lines=10,
1190                     orientation1=Orientation.XY,
1191                     base_angle1=0, base_angle2=0,
1192                     orientation2=Orientation.YZ,
1193                     cut_angle1=0, cut_angle2=0,
1194                     displacement1=0.5, displacement2=0.5,
1195                     generate_curves=False,
1196                     vector_mode='Magnitude'):
1197     """Creates Cut Lines presentation on the given field.
1198
1199     Arguments:
1200       proxy: the pipeline object, containig data
1201       entity: the entity type from PrsTypeEnum
1202       field_name: the field name
1203       timestamp_nb: the number of time step (1, 2, ...)
1204       nb_lines: number of lines
1205       orientation1: base plane orientation in 3D space
1206       base_angle1: rotation of the base plane in 3d space around the first
1207       axis of the orientation1 (X axis for XY, Y axis for YZ, Z axis for ZX).
1208       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1209       base_angle2: rotation of the base plane in 3d space around the second
1210       axis of the orientation1. Acceptable range: [-45, 45].
1211       orientation2: cutting planes orientation in 3D space
1212       cut_angle1: rotation of the cut planes in 3d space around the first
1213       axis of the orientation2. Acceptable range: [-45, 45].
1214       cut_angle2: rotation of the cuting planes in 3d space around the second
1215       axis of the orientation2. Acceptable range: [-45, 45].
1216       displacement1: base plane displacement
1217       displacement2: cutting planes displacement
1218       generate_curves: if true, 'PlotOverLine' filter will be created
1219       for each cut line
1220       vector_mode: the mode of transformation of vector values
1221       into scalar values, applicable only if the field contains vector values.
1222       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1223
1224     Returns:
1225       Cut Lines as representation object if generate_curves == False,
1226       (Cut Lines as representation object, list of 'PlotOverLine') otherwise
1227
1228     """
1229     if entity == EntityType.NODE:
1230         select_cells_with_data(proxy, on_points=[field_name])
1231     else:
1232         select_cells_with_data(proxy, on_cells=[field_name])
1233
1234     # Check vector mode
1235     nb_components = get_nb_components(proxy, entity, field_name)
1236     check_vector_mode(vector_mode, nb_components)
1237
1238     # Get time value
1239     time_value = get_time(proxy, timestamp_nb)
1240
1241     # Set timestamp
1242     pvs.GetRenderView().ViewTime = time_value
1243     pvs.UpdatePipeline(time_value, proxy)
1244
1245     # Create base plane
1246     base_plane = pvs.Slice(proxy)
1247     base_plane.SliceType = "Plane"
1248
1249     # Set base plane normal
1250     base_normal = get_normal_by_orientation(orientation1,
1251                                             radians(base_angle1),
1252                                             radians(base_angle2))
1253     base_plane.SliceType.Normal = base_normal
1254
1255     # Set base plane position
1256     base_position = get_positions(1, base_normal,
1257                                   get_bounds(proxy), displacement1)
1258     base_plane.SliceOffsetValues = base_position
1259
1260     # Check base plane
1261     base_plane.UpdatePipeline()
1262     if (base_plane.GetDataInformation().GetNumberOfCells() == 0):
1263         base_plane = proxy
1264
1265     # Create cutting planes
1266     cut_planes = pvs.Slice(base_plane)
1267     cut_planes.SliceType = "Plane"
1268
1269     # Set cutting planes normal and get positions
1270     cut_normal = get_normal_by_orientation(orientation2,
1271                                            radians(cut_angle1),
1272                                            radians(cut_angle2))
1273     cut_planes.SliceType.Normal = cut_normal
1274
1275     # Set cutting planes position
1276     cut_positions = get_positions(nb_lines, cut_normal,
1277                                   get_bounds(base_plane), displacement2)
1278
1279     # Generate curves
1280     curves = []
1281     if generate_curves:
1282         index = 0
1283         for pos in cut_positions:
1284             # Get points for plot over line objects
1285             cut_planes.SliceOffsetValues = pos
1286             cut_planes.UpdatePipeline()
1287             bounds = get_bounds(cut_planes)
1288             point1 = [bounds[0], bounds[2], bounds[4]]
1289             point2 = [bounds[1], bounds[3], bounds[5]]
1290
1291             # Create plot over line filter
1292             pol = pvs.PlotOverLine(cut_planes,
1293                                   Source="High Resolution Line Source")
1294             pvs.RenameSource('Y' + str(index), pol)
1295             pol.Source.Point1 = point1
1296             pol.Source.Point2 = point2
1297             pol.UpdatePipeline()
1298             curves.append(pol)
1299
1300             index += 1
1301
1302     cut_planes.SliceOffsetValues = cut_positions
1303     cut_planes.UpdatePipeline()
1304
1305     # Get Cut Lines representation object
1306     cut_lines = pvs.GetRepresentation(cut_planes)
1307
1308     # Get lookup table
1309     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1310
1311     # Set field range if necessary
1312     data_range = get_data_range(proxy, entity,
1313                                 field_name, vector_mode)
1314     lookup_table.LockScalarRange = 1
1315     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1316
1317     # Set properties
1318     cut_lines.ColorAttributeType = EntityType.get_pvtype(entity)
1319     cut_lines.ColorArrayName = field_name
1320     cut_lines.LookupTable = lookup_table
1321
1322     # Set wireframe represenatation mode
1323     cut_lines.Representation = 'Wireframe'
1324
1325     # Add scalar bar
1326     add_scalar_bar(field_name, nb_components,
1327                    vector_mode, lookup_table, time_value)
1328
1329     result = cut_lines
1330     # If curves were generated return tuple (cut lines, list of curves)
1331     if curves:
1332         result = cut_lines, curves
1333
1334     return result
1335
1336
1337 def CutSegmentOnField(proxy, entity, field_name, timestamp_nb,
1338                       point1, point2, vector_mode='Magnitude'):
1339     """Creates Cut Segment presentation on the given field.
1340
1341     Arguments:
1342       proxy: the pipeline object, containig data
1343       entity: the entity type from PrsTypeEnum
1344       field_name: the field name
1345       timestamp_nb: the number of time step (1, 2, ...)
1346       point1: set the first point of the segment (as [x, y, z])
1347       point1: set the second point of the segment (as [x, y, z])
1348       vector_mode: the mode of transformation of vector values
1349       into scalar values, applicable only if the field contains vector values.
1350       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1351
1352     Returns:
1353       Cut Segment as 3D representation object.
1354
1355     """
1356     if entity == EntityType.NODE:
1357         select_cells_with_data(proxy, on_points=[field_name])
1358     else:
1359         select_cells_with_data(proxy, on_cells=[field_name])
1360
1361     # Check vector mode
1362     nb_components = get_nb_components(proxy, entity, field_name)
1363     check_vector_mode(vector_mode, nb_components)
1364
1365     # Get time value
1366     time_value = get_time(proxy, timestamp_nb)
1367
1368     # Set timestamp
1369     pvs.GetRenderView().ViewTime = time_value
1370     pvs.UpdatePipeline(time_value, proxy)
1371
1372     # Create plot over line filter
1373     pol = pvs.PlotOverLine(proxy, Source="High Resolution Line Source")
1374     pol.Source.Point1 = point1
1375     pol.Source.Point2 = point2
1376     pol.UpdatePipeline()
1377
1378     # Get Cut Segment representation object
1379     cut_segment = pvs.GetRepresentation(pol)
1380
1381     # Get lookup table
1382     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1383
1384     # Set field range if necessary
1385     data_range = get_data_range(proxy, entity,
1386                                 field_name, vector_mode)
1387     lookup_table.LockScalarRange = 1
1388     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1389
1390     # Set properties
1391     cut_segment.ColorAttributeType = EntityType.get_pvtype(entity)
1392     cut_segment.ColorArrayName = field_name
1393     cut_segment.LookupTable = lookup_table
1394
1395     # Set wireframe represenatation mode
1396     cut_segment.Representation = 'Wireframe'
1397
1398     # Add scalar bar
1399     add_scalar_bar(field_name, nb_components,
1400                    vector_mode, lookup_table, time_value)
1401
1402     return cut_segment
1403
1404
1405 def VectorsOnField(proxy, entity, field_name, timestamp_nb,
1406                    scale_factor=None,
1407                    glyph_pos=GlyphPos.TAIL, glyph_type='2D Glyph',
1408                    is_colored=False, vector_mode='Magnitude'):
1409     """Creates Vectors presentation on the given field.
1410
1411     Arguments:
1412       proxy: the pipeline object, containig data
1413       entity: the entity type from PrsTypeEnum
1414       field_name: the field name
1415       timestamp_nb: the number of time step (1, 2, ...)
1416       scale_factor: scale factor
1417       glyph_pos: the position of glyphs
1418       glyph_type: the type of glyphs
1419       is_colored: this option allows to color the presentation according to
1420       the corresponding data array values
1421       vector_mode: the mode of transformation of vector values
1422       into scalar values, applicable only if the field contains vector values.
1423       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1424
1425     Returns:
1426       Vectors as representation object.
1427
1428     """
1429     if entity == EntityType.NODE:
1430         select_cells_with_data(proxy, on_points=[field_name])
1431     else:
1432         select_cells_with_data(proxy, on_cells=[field_name])
1433
1434     # Check vector mode
1435     nb_components = get_nb_components(proxy, entity, field_name)
1436     check_vector_mode(vector_mode, nb_components)
1437
1438     # Get time value
1439     time_value = get_time(proxy, timestamp_nb)
1440
1441     # Set timestamp
1442     pvs.GetRenderView().ViewTime = time_value
1443     pvs.UpdatePipeline(time_value, proxy)
1444
1445     # Extract only groups with data for the field
1446     source = proxy
1447
1448     # Cell centers
1449     if is_data_on_cells(proxy, field_name):
1450         cell_centers = pvs.CellCenters(source)
1451         cell_centers.VertexCells = 1
1452         source = cell_centers
1453
1454     vector_array = field_name
1455     # If the given vector array has only 2 components, add the third one
1456     if nb_components == 2:
1457         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1458         vector_array = calc.ResultArrayName
1459         source = calc
1460
1461     # Glyph
1462     glyph = pvs.Glyph(source)
1463     glyph.Vectors = vector_array
1464     glyph.ScaleMode = 'vector'
1465     glyph.MaskPoints = 0
1466
1467     # Set glyph type
1468     glyph.GlyphType = glyph_type
1469     if glyph_type == '2D Glyph':
1470         glyph.GlyphType.GlyphType = 'Arrow'
1471     elif glyph_type == 'Cone':
1472         glyph.GlyphType.Resolution = 7
1473         glyph.GlyphType.Height = 2
1474         glyph.GlyphType.Radius = 0.2
1475
1476     # Set glyph position if possible
1477     if glyph.GlyphType.GetProperty("Center"):
1478         if (glyph_pos == GlyphPos.TAIL):
1479             glyph.GlyphType.Center = [0.5, 0.0, 0.0]
1480         elif (glyph_pos == GlyphPos.HEAD):
1481             glyph.GlyphType.Center = [-0.5, 0.0, 0.0]
1482         elif (glyph_pos == GlyphPos.CENTER):
1483             glyph.GlyphType.Center = [0.0, 0.0, 0.0]
1484
1485     if scale_factor is not None:
1486         glyph.SetScaleFactor = scale_factor
1487     else:
1488         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1489                                       proxy, entity, field_name)
1490         glyph.SetScaleFactor = def_scale
1491
1492     glyph.UpdatePipeline()
1493
1494     # Get Vectors representation object
1495     vectors = pvs.GetRepresentation(glyph)
1496
1497     # Get lookup table
1498     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1499
1500     # Set field range if necessary
1501     data_range = get_data_range(proxy, entity,
1502                                 field_name, vector_mode)
1503     lookup_table.LockScalarRange = 1
1504     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1505
1506     # Set properties
1507     if (is_colored):
1508         vectors.ColorArrayName = 'GlyphVector'
1509     else:
1510         vectors.ColorArrayName = ''
1511     vectors.LookupTable = lookup_table
1512
1513     vectors.LineWidth = 1.0
1514
1515     # Set wireframe represenatation mode
1516     vectors.Representation = 'Wireframe'
1517
1518     # Add scalar bar
1519     add_scalar_bar(field_name, nb_components,
1520                    vector_mode, lookup_table, time_value)
1521
1522     return vectors
1523
1524
1525 def DeformedShapeOnField(proxy, entity, field_name,
1526                          timestamp_nb,
1527                          scale_factor=None, is_colored=False,
1528                          vector_mode='Magnitude'):
1529     """Creates Defromed Shape presentation on the given field.
1530
1531     Arguments:
1532       proxy: the pipeline object, containig data
1533       entity: the entity type from PrsTypeEnum
1534       field_name: the field name
1535       timestamp_nb: the number of time step (1, 2, ...)
1536       scale_factor: scale factor of the deformation
1537       is_colored: this option allows to color the presentation according to
1538       the corresponding data array values
1539       vector_mode: the mode of transformation of vector values
1540       into scalar values, applicable only if the field contains vector values.
1541       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1542
1543     Returns:
1544       Defromed Shape as representation object.
1545
1546     """
1547     # We don't need mesh parts with no data on them
1548     if entity == EntityType.NODE:
1549         select_cells_with_data(proxy, on_points=[field_name])
1550     else:
1551         select_cells_with_data(proxy, on_cells=[field_name])
1552
1553     # Check vector mode
1554     nb_components = get_nb_components(proxy, entity, field_name)
1555     check_vector_mode(vector_mode, nb_components)
1556
1557     # Get time value
1558     time_value = get_time(proxy, timestamp_nb)
1559
1560     # Set timestamp
1561     pvs.GetRenderView().ViewTime = time_value
1562     pvs.UpdatePipeline(time_value, proxy)
1563
1564     # Do merge
1565     source = pvs.MergeBlocks(proxy)
1566
1567     # Cell data to point data
1568     if is_data_on_cells(proxy, field_name):
1569         cell_to_point = pvs.CellDatatoPointData()
1570         cell_to_point.PassCellData = 1
1571         source = cell_to_point
1572
1573     vector_array = field_name
1574     # If the given vector array has only 2 components, add the third one
1575     if nb_components == 2:
1576         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1577         vector_array = calc.ResultArrayName
1578         source = calc
1579
1580     # Warp by vector
1581     warp_vector = pvs.WarpByVector(source)
1582     warp_vector.Vectors = [vector_array]
1583     if scale_factor is not None:
1584         warp_vector.ScaleFactor = scale_factor
1585     else:
1586         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1587                                       proxy, entity, field_name)
1588         warp_vector.ScaleFactor = def_scale
1589
1590     # Get Deformed Shape representation object
1591     defshape = pvs.GetRepresentation(warp_vector)
1592
1593     # Get lookup table
1594     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1595
1596     # Set field range if necessary
1597     data_range = get_data_range(proxy, entity,
1598                                 field_name, vector_mode)
1599     lookup_table.LockScalarRange = 1
1600     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1601
1602     # Set properties
1603     if is_colored:
1604         defshape.ColorAttributeType = EntityType.get_pvtype(entity)
1605         defshape.ColorArrayName = field_name
1606     else:
1607         defshape.ColorArrayName = ''
1608     defshape.LookupTable = lookup_table
1609
1610     # Set wireframe represenatation mode
1611     defshape.Representation = 'Wireframe'
1612
1613     # Add scalar bar
1614     add_scalar_bar(field_name, nb_components,
1615                    vector_mode, lookup_table, time_value)
1616
1617     return defshape
1618
1619
1620 def DeformedShapeAndScalarMapOnField(proxy, entity, field_name,
1621                                      timestamp_nb,
1622                                      scale_factor=None,
1623                                      scalar_entity=None,
1624                                      scalar_field_name=None,
1625                                      vector_mode='Magnitude'):
1626     """Creates Defromed Shape And Scalar Map presentation on the given field.
1627
1628     Arguments:
1629       proxy: the pipeline object, containig data
1630       entity: the entity type from PrsTypeEnum
1631       field_name: the field name
1632       timestamp_nb: the number of time step (1, 2, ...)
1633       scale_factor: scale factor of the deformation
1634       scalar_entity: scalar field entity
1635       scalar_field_name: scalar field, i.e. the field for coloring
1636       vector_mode: the mode of transformation of vector values
1637       into scalar values, applicable only if the field contains vector values.
1638       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1639
1640     Returns:
1641       Defromed Shape And Scalar Map as representation object.
1642
1643     """
1644     # We don't need mesh parts with no data on them
1645     on_points = []
1646     on_cells = []
1647
1648     if entity == EntityType.NODE:
1649         on_points.append(field_name)
1650     else:
1651         on_cells.append(field_name)
1652
1653     if scalar_entity and scalar_field_name:
1654         if scalar_entity == EntityType.NODE:
1655             on_points.append(scalar_field_name)
1656         else:
1657             on_cells.append(scalar_field_name)
1658
1659     select_cells_with_data(proxy, on_points, on_cells)
1660
1661     # Check vector mode
1662     nb_components = get_nb_components(proxy, entity, field_name)
1663     check_vector_mode(vector_mode, nb_components)
1664
1665     # Get time value
1666     time_value = get_time(proxy, timestamp_nb)
1667
1668     # Set timestamp
1669     pvs.GetRenderView().ViewTime = time_value
1670     pvs.UpdatePipeline(time_value, proxy)
1671
1672     # Set scalar field by default
1673     scalar_field_entity = scalar_entity
1674     scalar_field = scalar_field_name
1675     if (scalar_field_entity is None) or (scalar_field is None):
1676         scalar_field_entity = entity
1677         scalar_field = field_name
1678
1679     # Do merge
1680     source = pvs.MergeBlocks(proxy)
1681
1682     # Cell data to point data
1683     if is_data_on_cells(proxy, field_name):
1684         cell_to_point = pvs.CellDatatoPointData(source)
1685         cell_to_point.PassCellData = 1
1686         source = cell_to_point
1687
1688     vector_array = field_name
1689     # If the given vector array has only 2 components, add the third one
1690     if nb_components == 2:
1691         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1692         vector_array = calc.ResultArrayName
1693         source = calc
1694
1695     # Warp by vector
1696     warp_vector = pvs.WarpByVector(source)
1697     warp_vector.Vectors = [vector_array]
1698     if scale_factor is not None:
1699         warp_vector.ScaleFactor = scale_factor
1700     else:
1701         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1702                                       proxy, entity, field_name)
1703         warp_vector.ScaleFactor = def_scale
1704
1705     # Get Defromed Shape And Scalar Map representation object
1706     defshapemap = pvs.GetRepresentation(warp_vector)
1707
1708     # Get lookup table
1709     lookup_table = get_lookup_table(scalar_field, nb_components, vector_mode)
1710
1711     # Set field range if necessary
1712     data_range = get_data_range(proxy, scalar_field_entity,
1713                                 scalar_field, vector_mode)
1714     lookup_table.LockScalarRange = 1
1715     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1716
1717     # Set properties
1718     defshapemap.ColorArrayName = scalar_field
1719     defshapemap.LookupTable = lookup_table
1720     defshapemap.ColorAttributeType = EntityType.get_pvtype(scalar_field_entity)
1721
1722     # Add scalar bar
1723     add_scalar_bar(field_name, nb_components,
1724                    vector_mode, lookup_table, time_value)
1725
1726     return defshapemap
1727
1728
1729 def Plot3DOnField(proxy, entity, field_name, timestamp_nb,
1730                   orientation=Orientation.AUTO,
1731                   angle1=0, angle2=0,
1732                   position=0.5, is_relative=True,
1733                   scale_factor=None,
1734                   is_contour=False, nb_contours=32,
1735                   vector_mode='Magnitude'):
1736     """Creates Plot 3D presentation on the given field.
1737
1738     Arguments:
1739       proxy: the pipeline object, containig data
1740       entity: the entity type from PrsTypeEnum
1741       field_name: the field name
1742       timestamp_nb: the number of time step (1, 2, ...)
1743       orientation: the cut plane plane orientation in 3D space, if
1744       the input is planar - will not be taken into account
1745       angle1: rotation of the cut plane in 3d space around the first axis
1746       of the selected orientation (X axis for XY, Y axis for YZ,
1747       Z axis for ZX).
1748       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1749       angle2: rotation of the cut plane in 3d space around the second axis
1750       of the selected orientation. Acceptable range: [-45, 45].
1751       position: position of the cut plane in the object (ranging from 0 to 1).
1752       The value 0.5 corresponds to cutting by halves.
1753       is_relative: defines if the cut plane position is relative or absolute
1754       scale_factor: deformation scale factor
1755       is_contour: if True - Plot 3D will be represented with a set of contours,
1756       otherwise - Plot 3D will be represented with a smooth surface
1757       nb_contours: number of contours, applied if is_contour is True
1758       vector_mode: the mode of transformation of vector values
1759       into scalar values, applicable only if the field contains vector values.
1760       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1761
1762     Returns:
1763       Plot 3D as representation object.
1764
1765     """
1766     # We don't need mesh parts with no data on them
1767     if entity == EntityType.NODE:
1768         select_cells_with_data(proxy, on_points=[field_name])
1769     else:
1770         select_cells_with_data(proxy, on_cells=[field_name])
1771
1772     # Check vector mode
1773     nb_components = get_nb_components(proxy, entity, field_name)
1774     check_vector_mode(vector_mode, nb_components)
1775
1776     # Get time value
1777     time_value = get_time(proxy, timestamp_nb)
1778
1779     # Set timestamp
1780     pvs.GetRenderView().ViewTime = time_value
1781     pvs.UpdatePipeline(time_value, proxy)
1782
1783     # Do merge
1784     merge_blocks = pvs.MergeBlocks(proxy)
1785     merge_blocks.UpdatePipeline()
1786
1787     poly_data = None
1788
1789     # Cutting plane
1790
1791     # Define orientation if necessary (auto mode)
1792     plane_orientation = orientation
1793     if (orientation == Orientation.AUTO):
1794         plane_orientation = get_orientation(proxy)
1795
1796     # Get cutting plane normal
1797     normal = None
1798
1799     if (not is_planar_input(proxy)):
1800         normal = get_normal_by_orientation(plane_orientation,
1801                                            radians(angle1), radians(angle2))
1802
1803         # Create slice filter
1804         slice_filter = pvs.Slice(merge_blocks)
1805         slice_filter.SliceType = "Plane"
1806
1807         # Set cutting plane normal
1808         slice_filter.SliceType.Normal = normal
1809
1810         # Set cutting plane position
1811         if (is_relative):
1812             base_position = get_positions(1, normal,
1813                                           get_bounds(proxy), position)
1814             slice_filter.SliceOffsetValues = base_position
1815         else:
1816             slice_filter.SliceOffsetValues = position
1817
1818         slice_filter.UpdatePipeline()
1819         poly_data = slice_filter
1820     else:
1821         normal = get_normal_by_orientation(plane_orientation, 0, 0)
1822
1823     use_normal = 0
1824     # Geometry filter
1825     if not poly_data or poly_data.GetDataInformation().GetNumberOfCells() == 0:
1826         geometry_filter = pvs.GeometryFilter(merge_blocks)
1827         poly_data = geometry_filter
1828         use_normal = 1  # TODO(MZN): workaround
1829
1830     warp_scalar = None
1831     plot3d = None
1832     source = poly_data
1833
1834     if is_data_on_cells(poly_data, field_name):
1835         # Cell data to point data
1836         cell_to_point = pvs.CellDatatoPointData(poly_data)
1837         cell_to_point.PassCellData = 1
1838         source = cell_to_point
1839
1840     scalars = ['POINTS', field_name]
1841
1842     # Transform vector array to scalar array if necessary
1843     if (nb_components > 1):
1844         calc = get_calc_magnitude(source, EntityType.NODE, field_name)
1845         scalars = ['POINTS', calc.ResultArrayName]
1846         source = calc
1847
1848     # Warp by scalar
1849     warp_scalar = pvs.WarpByScalar(source)
1850     warp_scalar.Scalars = scalars
1851     warp_scalar.Normal = normal
1852     warp_scalar.UseNormal = use_normal
1853     if scale_factor is not None:
1854         warp_scalar.ScaleFactor = scale_factor
1855     else:
1856         def_scale = get_default_scale(PrsTypeEnum.PLOT3D,
1857                                       proxy, entity, field_name)
1858         warp_scalar.ScaleFactor = def_scale
1859
1860     warp_scalar.UpdatePipeline()
1861     source = warp_scalar
1862
1863     if (is_contour):
1864         # Contours
1865         contour = pvs.Contour(warp_scalar)
1866         contour.PointMergeMethod = "Uniform Binning"
1867         contour.ContourBy = ['POINTS', field_name]
1868         scalar_range = get_data_range(proxy, entity,
1869                                       field_name, vector_mode)
1870         contour.Isosurfaces = get_contours(scalar_range, nb_contours)
1871         contour.UpdatePipeline()
1872         source = contour
1873
1874     # Get Plot 3D representation object
1875     plot3d = pvs.GetRepresentation(source)
1876
1877     # Get lookup table
1878     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1879
1880     # Set field range if necessary
1881     data_range = get_data_range(proxy, entity,
1882                                 field_name, vector_mode)
1883     lookup_table.LockScalarRange = 1
1884     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1885
1886     # Set properties
1887     plot3d.ColorAttributeType = EntityType.get_pvtype(entity)
1888     plot3d.ColorArrayName = field_name
1889     plot3d.LookupTable = lookup_table
1890
1891     # Add scalar bar
1892     add_scalar_bar(field_name, nb_components,
1893                    vector_mode, lookup_table, time_value)
1894
1895     return plot3d
1896
1897
1898 def IsoSurfacesOnField(proxy, entity, field_name, timestamp_nb,
1899                        custom_range=None, nb_surfaces=10,
1900                        is_colored=True, color=None, vector_mode='Magnitude'):
1901     """Creates Iso Surfaces presentation on the given field.
1902
1903     Arguments:
1904       proxy: the pipeline object, containig data
1905       entity: the entity type from PrsTypeEnum
1906       field_name: the field name
1907       timestamp_nb: the number of time step (1, 2, ...)
1908       custom_range: scalar range, if undefined the source range will be applied
1909       nb_surfaces: number of surfaces, which will be generated
1910       is_colored: this option allows to color the presentation according to
1911       the corresponding data array values. If False - the presentation will
1912       be one-coloured.
1913       color: defines the presentation color as [R, G, B] triple. Taken into
1914       account only if is_colored is False.
1915       vector_mode: the mode of transformation of vector values
1916       into scalar values, applicable only if the field contains vector values.
1917       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1918
1919     Returns:
1920       Iso Surfaces as representation object.
1921
1922     """
1923     # We don't need mesh parts with no data on them
1924     if entity == EntityType.NODE:
1925         select_cells_with_data(proxy, on_points=[field_name])
1926     else:
1927         select_cells_with_data(proxy, on_cells=[field_name])
1928
1929     # Check vector mode
1930     nb_components = get_nb_components(proxy, entity, field_name)
1931     check_vector_mode(vector_mode, nb_components)
1932
1933     # Get time value
1934     time_value = get_time(proxy, timestamp_nb)
1935
1936     # Set timestamp
1937     pvs.GetRenderView().ViewTime = time_value
1938     pvs.UpdatePipeline(time_value, proxy)
1939
1940     # Do merge
1941     source = pvs.MergeBlocks(proxy)
1942
1943     # Transform cell data into point data if necessary
1944     if is_data_on_cells(proxy, field_name):
1945         cell_to_point = pvs.CellDatatoPointData(source)
1946         cell_to_point.PassCellData = 1
1947         source = cell_to_point
1948
1949     contour_by = ['POINTS', field_name]
1950
1951     # Transform vector array to scalar array if necessary
1952     if (nb_components > 1):
1953         calc = get_calc_magnitude(source, EntityType.NODE, field_name)
1954         contour_by = ['POINTS', calc.ResultArrayName]
1955         source = calc
1956
1957     # Contour filter settings
1958     contour = pvs.Contour(source)
1959     contour.ComputeScalars = 1
1960     contour.ContourBy = contour_by
1961
1962     # Specify the range
1963     scalar_range = custom_range
1964     if (scalar_range is None):
1965         scalar_range = get_data_range(proxy, entity,
1966                                       field_name, cut_off=True)
1967
1968     # Get contour values for the range
1969     surfaces = get_contours(scalar_range, nb_surfaces)
1970
1971     # Set contour values
1972     contour.Isosurfaces = surfaces
1973
1974     # Get Iso Surfaces representation object
1975     isosurfaces = pvs.GetRepresentation(contour)
1976
1977     # Get lookup table
1978     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1979
1980     # Set field range if necessary
1981     data_range = get_data_range(proxy, entity,
1982                                 field_name, vector_mode)
1983     lookup_table.LockScalarRange = 1
1984     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1985
1986     # Set display properties
1987     if (is_colored):
1988         isosurfaces.ColorAttributeType = EntityType.get_pvtype(entity)
1989         isosurfaces.ColorArrayName = field_name
1990     else:
1991         isosurfaces.ColorArrayName = ''
1992         if color:
1993             isosurfaces.DiffuseColor = color
1994     isosurfaces.LookupTable = lookup_table
1995
1996     # Add scalar bar
1997     add_scalar_bar(field_name, nb_components,
1998                    vector_mode, lookup_table, time_value)
1999
2000     return isosurfaces
2001
2002
2003 def GaussPointsOnField(proxy, entity, field_name,
2004                        timestamp_nb,
2005                        is_deformed=True, scale_factor=None,
2006                        is_colored=True, color=None,
2007                        primitive=GaussType.SPRITE,
2008                        is_proportional=True,
2009                        max_pixel_size=256,
2010                        multiplier=None, vector_mode='Magnitude'):
2011     """Creates Gauss Points on the given field.
2012
2013     Arguments:
2014
2015     proxy: the pipeline object, containig data
2016     entity: the field entity type from PrsTypeEnum
2017     field_name: the field name
2018     timestamp_nb: the number of time step (1, 2, ...)
2019     is_deformed: defines whether the Gauss Points will be deformed or not
2020     scale_factor -- the scale factor for deformation. Will be taken into
2021     account only if is_deformed is True.
2022     If not passed by user, default scale will be computed.
2023     is_colored -- defines whether the Gauss Points will be multicolored,
2024     using the corresponding data values
2025     color: defines the presentation color as [R, G, B] triple. Taken into
2026     account only if is_colored is False.
2027     primitive: primitive type from GaussType
2028     is_proportional: if True, the size of primitives will depends on
2029     the gauss point value
2030     max_pixel_size: the maximum sizr of the Gauss Points primitive in pixels
2031     multiplier: coefficient between data values and the size of primitives
2032     If not passed by user, default scale will be computed.
2033     vector_mode: the mode of transformation of vector values into
2034     scalar values, applicable only if the field contains vector values.
2035     Possible modes: 'Magnitude' - vector module;
2036     'X', 'Y', 'Z' - vector components.
2037
2038     Returns:
2039       Gauss Points as representation object.
2040
2041     """
2042     # We don't need mesh parts with no data on them
2043     on_gauss = select_cells_with_data(proxy, on_gauss=[field_name])
2044     if not on_gauss:
2045         if entity == EntityType.NODE:
2046             select_cells_with_data(proxy, on_points=[field_name])
2047         else:
2048             select_cells_with_data(proxy, on_cells=[field_name])
2049
2050     # Check vector mode
2051     nb_components = get_nb_components(proxy, entity, field_name)
2052     check_vector_mode(vector_mode, nb_components)
2053
2054     # Get time value
2055     time_value = get_time(proxy, timestamp_nb)
2056
2057     # Set timestamp
2058     pvs.GetRenderView().ViewTime = time_value
2059     proxy.UpdatePipeline(time=time_value)
2060
2061     source = proxy
2062
2063     # If no quadrature point array is passed, use cell centers
2064     if on_gauss:
2065         generate_qp = pvs.GenerateQuadraturePoints(source)
2066         generate_qp.QuadratureSchemeDef = ['CELLS', 'ELGA@0']
2067         source = generate_qp
2068     else:
2069         # Cell centers
2070         cell_centers = pvs.CellCenters(source)
2071         cell_centers.VertexCells = 1
2072         source = cell_centers
2073
2074     source.UpdatePipeline()
2075
2076     # Check if deformation enabled
2077     if is_deformed and nb_components > 1:
2078         vector_array = field_name
2079         # If the given vector array has only 2 components, add the third one
2080         if nb_components == 2:
2081             calc = get_add_component_calc(source,
2082                                           EntityType.NODE, field_name)
2083             vector_array = calc.ResultArrayName
2084             source = calc
2085
2086         # Warp by vector
2087         warp_vector = pvs.WarpByVector(source)
2088         warp_vector.Vectors = [vector_array]
2089         if scale_factor is not None:
2090             warp_vector.ScaleFactor = scale_factor
2091         else:
2092             def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE, proxy,
2093                                           entity, field_name)
2094             warp_vector.ScaleFactor = def_scale
2095         warp_vector.UpdatePipeline()
2096         source = warp_vector
2097
2098     # Get Gauss Points representation object
2099     gausspnt = pvs.GetRepresentation(source)
2100
2101     # Get lookup table
2102     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2103
2104     # Set field range if necessary
2105     data_range = get_data_range(proxy, entity,
2106                                 field_name, vector_mode)
2107     lookup_table.LockScalarRange = 1
2108     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2109
2110     # Set display properties
2111     if is_colored:
2112         gausspnt.ColorAttributeType = EntityType.get_pvtype(entity)
2113         gausspnt.ColorArrayName = field_name
2114     else:
2115         gausspnt.ColorArrayName = ''
2116         if color:
2117             gausspnt.DiffuseColor = color
2118
2119     gausspnt.LookupTable = lookup_table
2120
2121     # Add scalar bar
2122     add_scalar_bar(field_name, nb_components,
2123                    vector_mode, lookup_table, time_value)
2124
2125     # Set point sprite representation
2126     gausspnt.Representation = 'Point Sprite'
2127
2128     # Point sprite settings
2129     gausspnt.InterpolateScalarsBeforeMapping = 0
2130     gausspnt.MaxPixelSize = max_pixel_size
2131
2132     # Render mode
2133     gausspnt.RenderMode = GaussType.get_mode(primitive)
2134
2135     #if primitive == GaussType.SPRITE:
2136         # Set texture
2137         # TODO(MZN): replace with pvsimple high-level interface
2138     #    texture = sm.CreateProxy("textures", "SpriteTexture")
2139     #    alphamprop = texture.GetProperty("AlphaMethod")
2140     #    alphamprop.SetElement(0, 2)  # Clamp
2141     #    alphatprop = texture.GetProperty("AlphaThreshold")
2142     #    alphatprop.SetElement(0, 63)
2143     #    maxprop = texture.GetProperty("Maximum")
2144     #    maxprop.SetElement(0, 255)
2145     #    texture.UpdateVTKObjects()
2146
2147     #    gausspnt.Texture = texture
2148         #gausspnt.Texture.AlphaMethod = 'Clamp'
2149         #gausspnt.Texture.AlphaThreshold = 63
2150         #gausspnt.Texture.Maximum= 255
2151
2152     # Proportional radius
2153     gausspnt.RadiusUseScalarRange = 0
2154     gausspnt.RadiusIsProportional = 0
2155
2156     if is_proportional:
2157         mult = multiplier
2158         if mult is None:
2159             mult = abs(0.1 / data_range[1])
2160
2161         gausspnt.RadiusScalarRange = data_range
2162         gausspnt.RadiusTransferFunctionEnabled = 1
2163         gausspnt.RadiusMode = 'Scalar'
2164         gausspnt.RadiusArray = ['POINTS', field_name]
2165         if nb_components > 1:
2166             v_comp = get_vector_component(vector_mode)
2167             gausspnt.RadiusVectorComponent = v_comp
2168         gausspnt.RadiusTransferFunctionMode = 'Table'
2169         gausspnt.RadiusScalarRange = data_range
2170         gausspnt.RadiusUseScalarRange = 1
2171         gausspnt.RadiusIsProportional = 1
2172         gausspnt.RadiusProportionalFactor = mult
2173     else:
2174         gausspnt.RadiusTransferFunctionEnabled = 0
2175         gausspnt.RadiusMode = 'Constant'
2176         gausspnt.RadiusArray = ['POINTS', 'Constant Radius']
2177
2178     return gausspnt
2179
2180 def GaussPointsOnField1(proxy, entity, field_name,
2181                         timestamp_nb,
2182                         is_colored=True, color=None,
2183                         primitive=GaussType.SPHERE,
2184                         is_proportional=True,
2185                         max_pixel_size=256,
2186                         multiplier=None,
2187                         vector_mode='Magnitude'):
2188     """Creates Gauss Points on the given field. Use GaussPoints() Paraview interface.
2189
2190     Arguments:
2191     proxy: the pipeline object, containig data
2192     entity: the field entity type from PrsTypeEnum
2193     field_name: the field name
2194     timestamp_nb: the number of time step (1, 2, ...)
2195     is_colored -- defines whether the Gauss Points will be multicolored,
2196     using the corresponding data values
2197     color: defines the presentation color as [R, G, B] triple. Taken into
2198     account only if is_colored is False.
2199     primitive: primitive type from GaussType
2200     is_proportional: if True, the size of primitives will depends on
2201     the gauss point value
2202     max_pixel_size: the maximum sizr of the Gauss Points primitive in pixels
2203     multiplier: coefficient between data values and the size of primitives
2204     If not passed by user, default scale will be computed.
2205     vector_mode: the mode of transformation of vector values into
2206     scalar values, applicable only if the field contains vector values.
2207     Possible modes: 'Magnitude' - vector module;
2208     'X', 'Y', 'Z' - vector components.
2209
2210     Returns:
2211       Gauss Points as representation object.
2212
2213     """
2214     select_cells_with_data(proxy, on_gauss=[field_name])
2215
2216     nb_components = get_nb_components(proxy, entity, field_name)
2217
2218     # Get time value
2219     time_value = get_time(proxy, timestamp_nb)
2220
2221     # Set timestamp
2222     pvs.GetRenderView().ViewTime = time_value
2223     proxy.UpdatePipeline(time=time_value)
2224
2225     # Create Gauss Points object
2226     source = pvs.GaussPoints(proxy)
2227     source.UpdatePipeline()
2228   
2229     # Get Gauss Points representation object
2230     gausspnt = pvs.GetRepresentation(source)
2231
2232     # Get lookup table
2233     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2234
2235     # Set field range if necessary
2236     data_range = get_data_range(proxy, entity,
2237                                 field_name, vector_mode)
2238     lookup_table.LockScalarRange = 1
2239     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2240
2241     # Set display properties
2242     if is_colored:
2243         gausspnt.ColorAttributeType = EntityType.get_pvtype(entity)
2244         gausspnt.ColorArrayName = field_name
2245     else:
2246         gausspnt.ColorArrayName = ''
2247         if color:
2248             gausspnt.DiffuseColor = color
2249
2250     gausspnt.LookupTable = lookup_table
2251
2252     # Add scalar bar
2253     add_scalar_bar(field_name, nb_components,
2254                    vector_mode, lookup_table, time_value)
2255
2256     # Set point sprite representation
2257     gausspnt.Representation = 'Point Sprite'
2258
2259     # Point sprite settings
2260     gausspnt.InterpolateScalarsBeforeMapping = 0
2261     gausspnt.MaxPixelSize = max_pixel_size
2262
2263     # Render mode
2264     gausspnt.RenderMode = GaussType.get_mode(primitive)
2265
2266     #if primitive == GaussType.SPRITE:
2267         # Set texture
2268         # TODO(MZN): replace with pvsimple high-level interface
2269     #    texture = sm.CreateProxy("textures", "SpriteTexture")
2270     #    alphamprop = texture.GetProperty("AlphaMethod")
2271     #    alphamprop.SetElement(0, 2)  # Clamp
2272     #    alphatprop = texture.GetProperty("AlphaThreshold")
2273     #    alphatprop.SetElement(0, 63)
2274     #    maxprop = texture.GetProperty("Maximum")
2275     #    maxprop.SetElement(0, 255)
2276     #    texture.UpdateVTKObjects()
2277
2278     #    gausspnt.Texture = texture
2279         #gausspnt.Texture.AlphaMethod = 'Clamp'
2280         #gausspnt.Texture.AlphaThreshold = 63
2281         #gausspnt.Texture.Maximum= 255
2282
2283     # Proportional radius
2284     gausspnt.RadiusUseScalarRange = 0
2285     gausspnt.RadiusIsProportional = 0
2286
2287     if is_proportional:
2288         mult = multiplier
2289         if mult is None:
2290             mult = abs(0.1 / data_range[1])
2291
2292         gausspnt.RadiusScalarRange = data_range
2293         gausspnt.RadiusTransferFunctionEnabled = 1
2294         gausspnt.RadiusMode = 'Scalar'
2295         gausspnt.RadiusArray = ['POINTS', field_name]
2296         if nb_components > 1:
2297             v_comp = get_vector_component(vector_mode)
2298             gausspnt.RadiusVectorComponent = v_comp
2299         gausspnt.RadiusTransferFunctionMode = 'Table'
2300         gausspnt.RadiusScalarRange = data_range
2301         gausspnt.RadiusUseScalarRange = 1
2302         gausspnt.RadiusIsProportional = 1
2303         gausspnt.RadiusProportionalFactor = mult
2304     else:
2305         gausspnt.RadiusTransferFunctionEnabled = 0
2306         gausspnt.RadiusMode = 'Constant'
2307         gausspnt.RadiusArray = ['POINTS', 'Constant Radius']
2308
2309     return gausspnt
2310
2311 def StreamLinesOnField(proxy, entity, field_name, timestamp_nb,
2312                        direction='BOTH', is_colored=False, color=None,
2313                        vector_mode='Magnitude'):
2314     """Creates Stream Lines presentation on the given field.
2315
2316     Arguments:
2317       proxy: the pipeline object, containig data
2318       entity: the entity type from PrsTypeEnum
2319       field_name: the field name
2320       timestamp_nb: the number of time step (1, 2, ...)
2321       direction: the stream lines direction ('FORWARD', 'BACKWARD' or 'BOTH')
2322       is_colored: this option allows to color the presentation according to
2323       the corresponding data values. If False - the presentation will
2324       be one-coloured.
2325       color: defines the presentation color as [R, G, B] triple. Taken into
2326       account only if is_colored is False.
2327       vector_mode: the mode of transformation of vector values
2328       into scalar values, applicable only if the field contains vector values.
2329       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
2330
2331     Returns:
2332       Stream Lines as representation object.
2333
2334     """
2335     # We don't need mesh parts with no data on them
2336     if entity == EntityType.NODE:
2337         select_cells_with_data(proxy, on_points=[field_name])
2338     else:
2339         select_cells_with_data(proxy, on_cells=[field_name])
2340
2341     # Check vector mode
2342     nb_components = get_nb_components(proxy, entity, field_name)
2343     check_vector_mode(vector_mode, nb_components)
2344
2345     # Get time value
2346     time_value = get_time(proxy, timestamp_nb)
2347
2348     # Set timestamp
2349     pvs.GetRenderView().ViewTime = time_value
2350     pvs.UpdatePipeline(time_value, proxy)
2351
2352     # Do merge
2353     source = pvs.MergeBlocks(proxy)
2354
2355     # Cell data to point data
2356     if is_data_on_cells(proxy, field_name):
2357         cell_to_point = pvs.CellDatatoPointData(source)
2358         cell_to_point.PassCellData = 1
2359         cell_to_point.UpdatePipeline()
2360         source = cell_to_point
2361
2362     vector_array = field_name
2363     # If the given vector array has only 2 components, add the third one
2364     if nb_components == 2:
2365         calc = get_add_component_calc(source, EntityType.NODE, field_name)
2366         vector_array = calc.ResultArrayName
2367         calc.UpdatePipeline()
2368         source = calc
2369
2370     # Stream Tracer
2371     stream = pvs.StreamTracer(source)
2372     stream.SeedType = "Point Source"
2373     stream.Vectors = ['POINTS', vector_array]
2374     stream.SeedType = "Point Source"
2375     stream.IntegrationDirection = direction
2376     stream.IntegratorType = 'Runge-Kutta 2'
2377     stream.UpdatePipeline()
2378
2379     # Get Stream Lines representation object
2380     if is_empty(stream):
2381         return None
2382     streamlines = pvs.GetRepresentation(stream)
2383
2384     # Get lookup table
2385     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2386
2387     # Set field range if necessary
2388     data_range = get_data_range(proxy, entity,
2389                                 field_name, vector_mode)
2390     lookup_table.LockScalarRange = 1
2391     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2392
2393     # Set properties
2394     if is_colored:
2395         streamlines.ColorAttributeType = EntityType.get_pvtype(entity)
2396         streamlines.ColorArrayName = field_name
2397     else:
2398         streamlines.ColorArrayName = ''
2399         if color:
2400             streamlines.DiffuseColor = color
2401
2402     streamlines.LookupTable = lookup_table
2403
2404     # Add scalar bar
2405     add_scalar_bar(field_name, nb_components,
2406                    vector_mode, lookup_table, time_value)
2407
2408     return streamlines
2409
2410
2411 def MeshOnEntity(proxy, mesh_name, entity):
2412     """Creates submesh of the entity type for the mesh.
2413
2414     Arguments:
2415       proxy -- the pipeline object, containig data
2416       mesh_name -- the mesh name
2417       entity -- the entity type
2418
2419     Returns:
2420       Submesh as representation object of the given source.
2421
2422     """
2423     # Select all cell types
2424     select_all_cells(proxy)
2425
2426     # Get subset of groups on the given entity
2427     subset = get_group_names(proxy, mesh_name, entity)
2428
2429     # Select only groups of the given entity type
2430     proxy.Groups = subset
2431     proxy.UpdatePipeline()
2432
2433     # Get representation object if the submesh is not empty
2434     prs = None
2435     if (proxy.GetDataInformation().GetNumberOfPoints() or
2436         proxy.GetDataInformation().GetNumberOfCells()):
2437         prs = pvs.GetRepresentation(proxy)
2438         prs.ColorArrayName = ''
2439
2440     return prs
2441
2442
2443 def MeshOnGroup(proxy, group_name):
2444     """Creates submesh on the group.
2445
2446     Arguments:
2447       proxy -- the pipeline object, containig data
2448       group_name -- the full group name
2449
2450     Returns:
2451       Representation object of the given source with single group
2452       selected.
2453
2454     """
2455     # Select all cell types
2456     select_all_cells(proxy)
2457
2458     # Select only the group with the given name
2459     one_group = [group_name]
2460     proxy.Groups = one_group
2461     proxy.UpdatePipeline()
2462
2463     # Get representation object if the submesh is not empty
2464     prs = None
2465
2466     # Check if the group was set
2467     if proxy.Groups.GetData() == one_group:
2468         group_entity = get_group_entity(group_name)
2469         # Check if the submesh is not empty
2470         nb_items = 0
2471         if group_entity == EntityType.NODE:
2472             nb_items = proxy.GetDataInformation().GetNumberOfPoints()
2473         elif group_entity == EntityType.CELL:
2474             nb_items = proxy.GetDataInformation().GetNumberOfCells()
2475
2476         if nb_items:
2477             prs = pvs.GetRepresentation(proxy)
2478             prs.ColorArrayName = ''
2479
2480     return prs
2481
2482
2483 def CreatePrsForFile(paravis_instance, file_name, prs_types,
2484                      picture_dir, picture_ext):
2485     """Build presentations of the given types for the file.
2486
2487     Build presentations for all fields on all timestamps.
2488
2489     Arguments:
2490       paravis_instance: ParaVis module instance object
2491       file_name: full path to the MED file
2492       prs_types: the list of presentation types to build
2493       picture_dir: the directory path for saving snapshots
2494       picture_ext: graphics files extension (determines file type)
2495
2496     """
2497     # Import MED file
2498     print "Import " + file_name.split(os.sep)[-1] + "..."
2499
2500     try:
2501         proxy = pvs.MEDReader(FileName=file_name)
2502         if proxy is None:
2503             print "FAILED"
2504         else:
2505             proxy.UpdatePipeline()
2506             _med_field_sep = proxy.GetProperty("Separator")
2507             print "OK"
2508     except:
2509         print "FAILED"
2510     else:
2511         # Get view
2512         view = pvs.GetRenderView()
2513
2514         # Create required presentations for the proxy
2515         CreatePrsForProxy(proxy, view, prs_types,
2516                           picture_dir, picture_ext)
2517
2518 def CreatePrsForProxy(proxy, view, prs_types, picture_dir, picture_ext):
2519     """Build presentations of the given types for all fields of the proxy.
2520
2521     Save snapshots in graphics files (type depends on the given extension).
2522     Stores the files in the given directory.
2523
2524     Arguments:
2525       proxy: the pipeline object, containig data
2526       view: the render view
2527       prs_types: the list of presentation types to build
2528       picture_dir: the directory path for saving snapshots
2529       picture_ext: graphics files extension (determines file type)
2530
2531     """
2532     # List of the field names
2533     fields_info = proxy.GetProperty("FieldsTreeInfo")[::2]
2534     print fields_info
2535
2536     # Add path separator to the end of picture path if necessery
2537     if not picture_dir.endswith(os.sep):
2538         picture_dir += os.sep
2539
2540     # Mesh Presentation
2541     if PrsTypeEnum.MESH in prs_types:
2542         # Create Mesh presentation. Build all possible submeshes.
2543
2544         extGrp=pvs.ExtractGroup()
2545
2546         # Remember the current state
2547         groups = filter(lambda x:x[:4]=="GRP_",list(extGrp.GetProperty("GroupsFlagsInfo")[::2]))
2548
2549         # Iterate on meshes
2550         mesh_names = get_mesh_names(proxy)
2551         for mesh_name in mesh_names:
2552             # Build mesh on nodes and cells
2553             for entity in (EntityType.NODE, EntityType.CELL):
2554                 entity_name = EntityType.get_name(entity)
2555                 if if_possible(proxy, mesh_name, entity, PrsTypeEnum.MESH):
2556                     print "Creating submesh on " + entity_name + " for '" + mesh_name + "' mesh... "
2557                     prs = MeshOnEntity(proxy, mesh_name, entity)
2558                     if prs is None:
2559                         print "FAILED"
2560                         continue
2561                     else:
2562                         print "OK"
2563                     # Construct image file name
2564                     pic_name = picture_dir + mesh_name + "_" + entity_name + "." + picture_ext
2565
2566                     # Show and dump the presentation into a graphics file
2567                     process_prs_for_test(prs, view, pic_name, False)
2568
2569                 # Build submesh on all groups of the mesh
2570                 mesh_groups = get_group_names(proxy, mesh_name,
2571                                               entity, wo_nogroups=True)
2572                 for group in mesh_groups:
2573                     print "Creating submesh on group " + group + "... "
2574                     prs = MeshOnGroup(proxy, group)
2575                     if prs is None:
2576                         print "FAILED"
2577                         continue
2578                     else:
2579                         print "OK"
2580                     # Construct image file name
2581                     pic_name = picture_dir + group.replace('/', '_') + "." + picture_ext
2582
2583                     # Show and dump the presentation into a graphics file
2584                     process_prs_for_test(prs, view, pic_name, False)
2585
2586         # Restore the state
2587         extGrp.AllGroups = groups
2588         extGrp.UpdatePipelineInformation()
2589
2590     # Presentations on fields
2591     for field in fields_info:
2592         field_name = get_field_short_name(field)
2593         # Ignore mesh presentation
2594         if field_name == get_field_mesh_name(field):
2595             continue
2596         field_entity = get_field_entity(field)
2597         # Clear fields selection state
2598         proxy.AllArrays = []
2599         proxy.UpdatePipeline()
2600         # Select only the current field:
2601         # necessary for getting the right timestamps
2602         proxy.AllArrays = field
2603         proxy.UpdatePipeline()
2604
2605         # Get timestamps
2606         entity_data_info = proxy.GetCellDataInformation()
2607         timestamps = proxy.TimestepValues.GetData()
2608
2609         for prs_type in prs_types:
2610             # Ignore mesh presentation
2611             if prs_type == PrsTypeEnum.MESH:
2612                 continue
2613
2614             # Get name of presentation type
2615             prs_name = PrsTypeEnum.get_name(prs_type)
2616
2617             # Build the presentation if possible
2618             possible = if_possible(proxy, field_name,
2619                                    field_entity, prs_type)
2620             if possible:
2621                 # Presentation type for graphics file name
2622                 f_prs_type = prs_name.replace(' ', '').upper()
2623
2624                 for timestamp_nb in xrange(1, len(timestamps) + 1):
2625                     time = timestamps[timestamp_nb - 1]
2626                     print "Creating " + prs_name + " on " + field_name + ", time = " + str(time) + "... "
2627                     prs = create_prs(prs_type, proxy,
2628                                      field_entity, field_name, timestamp_nb)
2629                     if prs is None:
2630                         print "FAILED"
2631                         continue
2632                     else:
2633                         print "OK"
2634
2635                     # Construct image file name
2636                     pic_name = picture_dir + field_name + "_" + str(time) + "_" + f_prs_type + "." + picture_ext
2637
2638                     # Show and dump the presentation into a graphics file
2639                     process_prs_for_test(prs, view, pic_name)
2640     return