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