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