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