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