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