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