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