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