Salome HOME
Update dump trace/study tests.
[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     if hasattr(lookup_table,"LockDataRange"):
959         lookup_table.LockDataRange = 0
960     elif hasattr(lookup_table,"LockScalarRange"):
961         lookup_table.LockScalarRange = 0
962     else:
963         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
964
965     return lookup_table
966
967
968 def get_group_mesh_name(full_group_name):
969     """Return mesh name of the group by its full name."""
970     aList = full_group_name.split('/')
971     if len(aList) >= 2 :
972         group_name = full_group_name.split('/')[1]
973         return group_name
974
975 def get_group_entity(full_group_name):
976     """Return entity type of the group by its full name."""
977     aList = full_group_name.split('/')
978     if len(aList) >= 3 :
979         entity_name = full_group_name.split('/')[2]
980         entity = EntityType.get_type(entity_name)
981         return entity
982
983
984 def get_group_short_name(full_group_name):
985     """Return short name of the group by its full name."""
986     short_name = re.sub('^GRP_', '', full_group_name)
987     return short_name
988
989
990 def get_mesh_full_names(proxy):
991     """Return all mesh names in the given proxy as a set."""
992     proxy.UpdatePipeline()
993     fields = proxy.GetProperty("FieldsTreeInfo")[::2]
994     mesh_full_names = set([item for item in fields if get_field_mesh_name(item) == get_field_short_name(item)])
995     return mesh_full_names
996
997
998 def get_group_names(extrGrps):
999     """Return full names of all groups of the given 'ExtractGroup' filter object.
1000     """
1001     group_names = filter(lambda x:x[:4]=="GRP_",list(extrGrps.GetProperty("GroupsFlagsInfo")[::2]))
1002     return group_names
1003
1004
1005 def get_time(proxy, timestamp_nb):
1006     """Get time value by timestamp number."""
1007     #proxy.UpdatePipeline()
1008     # Check timestamp number
1009     timestamps = []
1010
1011     if (hasattr(proxy, 'TimestepValues')):
1012         timestamps = proxy.TimestepValues.GetData()
1013     elif (hasattr(proxy.Input, 'TimestepValues')):
1014         timestamps = proxy.Input.TimestepValues.GetData()
1015
1016     length = len(timestamps)
1017     if (timestamp_nb > 0 and (timestamp_nb - 1) not in xrange(length) ) or (timestamp_nb < 0 and -timestamp_nb > length):
1018         raise ValueError("Timestamp number is out of range: " + str(timestamp_nb))
1019
1020     if not timestamps:
1021         return 0.0
1022
1023     # Return time value
1024     if timestamp_nb > 0:
1025         return timestamps[timestamp_nb - 1]
1026     else:
1027         return timestamps[timestamp_nb]
1028
1029 def create_prs(prs_type, proxy, field_entity, field_name, timestamp_nb):
1030     """Auxiliary function.
1031
1032     Build presentation of the given type on the given field and
1033     timestamp number.
1034     Set the presentation properties like visu.CreatePrsForResult() do.
1035
1036     """
1037     proxy.UpdatePipeline()
1038     prs = None
1039
1040     if prs_type == PrsTypeEnum.SCALARMAP:
1041         prs = ScalarMapOnField(proxy, field_entity, field_name, timestamp_nb)
1042     elif prs_type == PrsTypeEnum.CUTPLANES:
1043         prs = CutPlanesOnField(proxy, field_entity, field_name, timestamp_nb,
1044                                orientation=Orientation.ZX)
1045     elif prs_type == PrsTypeEnum.CUTLINES:
1046         prs = CutLinesOnField(proxy, field_entity, field_name, timestamp_nb,
1047                               orientation1=Orientation.XY,
1048                               orientation2=Orientation.ZX)
1049     elif prs_type == PrsTypeEnum.DEFORMEDSHAPE:
1050         prs = DeformedShapeOnField(proxy, field_entity,
1051                                    field_name, timestamp_nb)
1052     elif prs_type == PrsTypeEnum.DEFORMEDSHAPESCALARMAP:
1053         prs = DeformedShapeAndScalarMapOnField(proxy, field_entity,
1054                                                field_name, timestamp_nb)
1055     elif prs_type == PrsTypeEnum.VECTORS:
1056         prs = VectorsOnField(proxy, field_entity, field_name, timestamp_nb)
1057     elif prs_type == PrsTypeEnum.PLOT3D:
1058         prs = Plot3DOnField(proxy, field_entity, field_name, timestamp_nb)
1059     elif prs_type == PrsTypeEnum.ISOSURFACES:
1060         prs = IsoSurfacesOnField(proxy, field_entity, field_name, timestamp_nb)
1061     elif prs_type == PrsTypeEnum.GAUSSPOINTS:
1062         prs = GaussPointsOnField(proxy, field_entity, field_name, timestamp_nb)
1063     elif prs_type == PrsTypeEnum.STREAMLINES:
1064         prs = StreamLinesOnField(proxy, field_entity, field_name, timestamp_nb)
1065     else:
1066         raise ValueError("Unexistent presentation type.")
1067
1068     return prs
1069
1070
1071 # Functions for building Post-Pro presentations
1072 def ScalarMapOnField(proxy, entity, field_name, timestamp_nb,
1073                      vector_mode='Magnitude'):
1074     """Creates Scalar Map presentation on the given field.
1075
1076     Arguments:
1077       proxy: the pipeline object, containig data
1078       entity: the entity type from PrsTypeEnum
1079       field_name: the field name
1080       timestamp_nb: the number of time step (1, 2, ...)
1081       vector_mode: the mode of transformation of vector values
1082       into scalar values, applicable only if the field contains vector values.
1083       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1084
1085     Returns:
1086       Scalar Map as representation object.
1087
1088     """
1089     proxy.UpdatePipeline()
1090     # We don't need mesh parts with no data on them
1091     if entity == EntityType.NODE:
1092         select_cells_with_data(proxy, on_points=[field_name])
1093     else:
1094         select_cells_with_data(proxy, on_cells=[field_name])
1095
1096     # Check vector mode
1097     nb_components = get_nb_components(proxy, entity, field_name)
1098     check_vector_mode(vector_mode, nb_components)
1099
1100     # Get time value
1101     time_value = get_time(proxy, timestamp_nb)
1102
1103     # Set timestamp
1104     pvs.GetRenderView().ViewTime = time_value
1105     pvs.UpdatePipeline(time_value, proxy)
1106
1107     # Get Scalar Map representation object
1108     scalarmap = pvs.GetRepresentation(proxy)
1109
1110     # Get lookup table
1111     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1112
1113     # Set field range if necessary
1114     data_range = get_data_range(proxy, entity,
1115                                 field_name, vector_mode)
1116     if hasattr(lookup_table,"LockDataRange"):
1117         lookup_table.LockDataRange = 1
1118     elif hasattr(lookup_table,"LockScalarRange"):
1119         lookup_table.LockScalarRange = 1
1120     else:
1121         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1122
1123     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1124     # Set properties
1125     pvs.ColorBy(scalarmap, (EntityType.get_pvtype(entity), field_name))
1126     scalarmap.LookupTable = lookup_table
1127
1128     # Add scalar bar
1129     bar_title = field_name + ", " + str(time_value)
1130     if (nb_components > 1):
1131         bar_title += "\n" + vector_mode
1132     add_scalar_bar(field_name, nb_components, vector_mode,
1133                    lookup_table, time_value)
1134
1135     return scalarmap
1136
1137
1138 def CutPlanesOnField(proxy, entity, field_name, timestamp_nb,
1139                      nb_planes=10, orientation=Orientation.YZ,
1140                      angle1=0, angle2=0,
1141                      displacement=0.5, vector_mode='Magnitude'):
1142     """Creates Cut Planes presentation on the given field.
1143
1144     Arguments:
1145       proxy: the pipeline object, containig data
1146       entity: the entity type from PrsTypeEnum
1147       field_name: the field name
1148       timestamp_nb: the number of time step (1, 2, ...)
1149       nb_planes: number of cutting planes
1150       orientation: cutting planes orientation in 3D space
1151       angle1: rotation of the planes in 3d space around the first axis of the
1152       selected orientation (X axis for XY, Y axis for YZ, Z axis for ZX).
1153       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1154       angle2: rotation of the planes in 3d space around the second axis of the
1155       selected orientation. Acceptable range: [-45, 45].
1156       displacement: the displacement of the planes into one or another side
1157       vector_mode: the mode of transformation of vector values
1158       into scalar values, applicable only if the field contains vector values.
1159       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1160
1161     Returns:
1162       Cut Planes as representation object.
1163
1164     """
1165     proxy.UpdatePipeline()
1166     if entity == EntityType.NODE:
1167         select_cells_with_data(proxy, on_points=[field_name])
1168     else:
1169         select_cells_with_data(proxy, on_cells=[field_name])
1170
1171     # Check vector mode
1172     nb_components = get_nb_components(proxy, entity, field_name)
1173     check_vector_mode(vector_mode, nb_components)
1174
1175     # Get time value
1176     time_value = get_time(proxy, timestamp_nb)
1177
1178     # Set timestamp
1179     pvs.GetRenderView().ViewTime = time_value
1180     pvs.UpdatePipeline(time_value, proxy)
1181
1182     # Create slice filter
1183     slice_filter = pvs.Slice(proxy)
1184     slice_filter.SliceType = "Plane"
1185
1186     # Set cut planes normal
1187     normal = get_normal_by_orientation(orientation,
1188                                        radians(angle1), radians(angle2))
1189     slice_filter.SliceType.Normal = normal
1190
1191     # Set cut planes positions
1192     positions = get_positions(nb_planes, normal,
1193                               get_bounds(proxy), displacement)
1194     slice_filter.SliceOffsetValues = positions
1195
1196     # Get Cut Planes representation object
1197     cut_planes = pvs.GetRepresentation(slice_filter)
1198
1199     # Get lookup table
1200     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1201
1202     # Set field range if necessary
1203     data_range = get_data_range(proxy, entity,
1204                                 field_name, vector_mode)
1205
1206     if hasattr(lookup_table,"LockDataRange"):
1207         lookup_table.LockDataRange = 1
1208     elif hasattr(lookup_table,"LockScalarRange"):
1209         lookup_table.LockScalarRange = 1
1210     else:
1211         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1212
1213     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1214
1215     # Set properties
1216     pvs.ColorBy(cut_planes, (EntityType.get_pvtype(entity), field_name))
1217     cut_planes.LookupTable = lookup_table
1218
1219     # Add scalar bar
1220     add_scalar_bar(field_name, nb_components,
1221                    vector_mode, lookup_table, time_value)
1222
1223     return cut_planes
1224
1225
1226 def CutLinesOnField(proxy, entity, field_name, timestamp_nb,
1227                     nb_lines=10,
1228                     orientation1=Orientation.XY,
1229                     base_angle1=0, base_angle2=0,
1230                     orientation2=Orientation.YZ,
1231                     cut_angle1=0, cut_angle2=0,
1232                     displacement1=0.5, displacement2=0.5,
1233                     generate_curves=False,
1234                     vector_mode='Magnitude'):
1235     """Creates Cut Lines presentation on the given field.
1236
1237     Arguments:
1238       proxy: the pipeline object, containig data
1239       entity: the entity type from PrsTypeEnum
1240       field_name: the field name
1241       timestamp_nb: the number of time step (1, 2, ...)
1242       nb_lines: number of lines
1243       orientation1: base plane orientation in 3D space
1244       base_angle1: rotation of the base plane in 3d space around the first
1245       axis of the orientation1 (X axis for XY, Y axis for YZ, Z axis for ZX).
1246       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1247       base_angle2: rotation of the base plane in 3d space around the second
1248       axis of the orientation1. Acceptable range: [-45, 45].
1249       orientation2: cutting planes orientation in 3D space
1250       cut_angle1: rotation of the cut planes in 3d space around the first
1251       axis of the orientation2. Acceptable range: [-45, 45].
1252       cut_angle2: rotation of the cuting planes in 3d space around the second
1253       axis of the orientation2. Acceptable range: [-45, 45].
1254       displacement1: base plane displacement
1255       displacement2: cutting planes displacement
1256       generate_curves: if true, 'PlotOverLine' filter will be created
1257       for each cut line
1258       vector_mode: the mode of transformation of vector values
1259       into scalar values, applicable only if the field contains vector values.
1260       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1261
1262     Returns:
1263       Cut Lines as representation object if generate_curves == False,
1264       (Cut Lines as representation object, list of 'PlotOverLine') otherwise
1265
1266     """
1267     proxy.UpdatePipeline()
1268     if entity == EntityType.NODE:
1269         select_cells_with_data(proxy, on_points=[field_name])
1270     else:
1271         select_cells_with_data(proxy, on_cells=[field_name])
1272
1273     # Check vector mode
1274     nb_components = get_nb_components(proxy, entity, field_name)
1275     check_vector_mode(vector_mode, nb_components)
1276
1277     # Get time value
1278     time_value = get_time(proxy, timestamp_nb)
1279
1280     # Set timestamp
1281     pvs.GetRenderView().ViewTime = time_value
1282     pvs.UpdatePipeline(time_value, proxy)
1283
1284     # Create base plane
1285     base_plane = pvs.Slice(proxy)
1286     base_plane.SliceType = "Plane"
1287
1288     # Set base plane normal
1289     base_normal = get_normal_by_orientation(orientation1,
1290                                             radians(base_angle1),
1291                                             radians(base_angle2))
1292     base_plane.SliceType.Normal = base_normal
1293
1294     # Set base plane position
1295     base_position = get_positions(1, base_normal,
1296                                   get_bounds(proxy), displacement1)
1297     base_plane.SliceOffsetValues = base_position
1298
1299     # Check base plane
1300     base_plane.UpdatePipeline()
1301     if (base_plane.GetDataInformation().GetNumberOfCells() == 0):
1302         base_plane = proxy
1303
1304     # Create cutting planes
1305     cut_planes = pvs.Slice(base_plane)
1306     cut_planes.SliceType = "Plane"
1307
1308     # Set cutting planes normal and get positions
1309     cut_normal = get_normal_by_orientation(orientation2,
1310                                            radians(cut_angle1),
1311                                            radians(cut_angle2))
1312     cut_planes.SliceType.Normal = cut_normal
1313
1314     # Set cutting planes position
1315     cut_positions = get_positions(nb_lines, cut_normal,
1316                                   get_bounds(base_plane), displacement2)
1317
1318     # Generate curves
1319     curves = []
1320     if generate_curves:
1321         index = 0
1322         for pos in cut_positions:
1323             # Get points for plot over line objects
1324             cut_planes.SliceOffsetValues = pos
1325             cut_planes.UpdatePipeline()
1326             bounds = get_bounds(cut_planes)
1327             point1 = [bounds[0], bounds[2], bounds[4]]
1328             point2 = [bounds[1], bounds[3], bounds[5]]
1329
1330             # Create plot over line filter
1331             pol = pvs.PlotOverLine(cut_planes,
1332                                   Source="High Resolution Line Source")
1333             pvs.RenameSource('Y' + str(index), pol)
1334             pol.Source.Point1 = point1
1335             pol.Source.Point2 = point2
1336             pol.UpdatePipeline()
1337             curves.append(pol)
1338
1339             index += 1
1340
1341     cut_planes.SliceOffsetValues = cut_positions
1342     cut_planes.UpdatePipeline()
1343
1344     # Get Cut Lines representation object
1345     cut_lines = pvs.GetRepresentation(cut_planes)
1346
1347     # Get lookup table
1348     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1349
1350     # Set field range if necessary
1351     data_range = get_data_range(proxy, entity,
1352                                 field_name, vector_mode)
1353     if hasattr(lookup_table,"LockDataRange"):
1354         lookup_table.LockDataRange = 1
1355     elif hasattr(lookup_table,"LockScalarRange"):
1356         lookup_table.LockScalarRange = 1
1357     else:
1358         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1359
1360     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1361
1362     # Set properties
1363     pvs.ColorBy(cut_lines, (EntityType.get_pvtype(entity), field_name))
1364     cut_lines.LookupTable = lookup_table
1365
1366     # Set wireframe represenatation mode
1367     cut_lines.Representation = 'Wireframe'
1368
1369     # Add scalar bar
1370     add_scalar_bar(field_name, nb_components,
1371                    vector_mode, lookup_table, time_value)
1372
1373     result = cut_lines
1374     # If curves were generated return tuple (cut lines, list of curves)
1375     if curves:
1376         result = cut_lines, curves
1377
1378     return result
1379
1380
1381 def CutSegmentOnField(proxy, entity, field_name, timestamp_nb,
1382                       point1, point2, vector_mode='Magnitude'):
1383     """Creates Cut Segment presentation on the given field.
1384
1385     Arguments:
1386       proxy: the pipeline object, containig data
1387       entity: the entity type from PrsTypeEnum
1388       field_name: the field name
1389       timestamp_nb: the number of time step (1, 2, ...)
1390       point1: set the first point of the segment (as [x, y, z])
1391       point1: set the second point of the segment (as [x, y, z])
1392       vector_mode: the mode of transformation of vector values
1393       into scalar values, applicable only if the field contains vector values.
1394       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1395
1396     Returns:
1397       Cut Segment as 3D representation object.
1398
1399     """
1400     proxy.UpdatePipeline()
1401     if entity == EntityType.NODE:
1402         select_cells_with_data(proxy, on_points=[field_name])
1403     else:
1404         select_cells_with_data(proxy, on_cells=[field_name])
1405
1406     # Check vector mode
1407     nb_components = get_nb_components(proxy, entity, field_name)
1408     check_vector_mode(vector_mode, nb_components)
1409
1410     # Get time value
1411     time_value = get_time(proxy, timestamp_nb)
1412
1413     # Set timestamp
1414     pvs.GetRenderView().ViewTime = time_value
1415     pvs.UpdatePipeline(time_value, proxy)
1416
1417     # Create plot over line filter
1418     pol = pvs.PlotOverLine(proxy, Source="High Resolution Line Source")
1419     pol.Source.Point1 = point1
1420     pol.Source.Point2 = point2
1421     pol.UpdatePipeline()
1422
1423     # Get Cut Segment representation object
1424     cut_segment = pvs.GetRepresentation(pol)
1425
1426     # Get lookup table
1427     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1428
1429     # Set field range if necessary
1430     data_range = get_data_range(proxy, entity,
1431                                 field_name, vector_mode)
1432     if hasattr(lookup_table,"LockDataRange"):
1433         lookup_table.LockDataRange = 1
1434     elif hasattr(lookup_table,"LockScalarRange"):
1435         lookup_table.LockScalarRange = 1
1436     else:
1437         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1438
1439     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1440
1441     # Set properties
1442     pvs.ColorBy(cut_segment, (EntityType.get_pvtype(entity), field_name))
1443     cut_segment.LookupTable = lookup_table
1444
1445     # Set wireframe represenatation mode
1446     cut_segment.Representation = 'Wireframe'
1447
1448     # Add scalar bar
1449     add_scalar_bar(field_name, nb_components,
1450                    vector_mode, lookup_table, time_value)
1451
1452     return cut_segment
1453
1454
1455 def VectorsOnField(proxy, entity, field_name, timestamp_nb,
1456                    scale_factor=None,
1457                    glyph_pos=GlyphPos.TAIL, glyph_type='2D Glyph',
1458                    is_colored=False, vector_mode='Magnitude'):
1459     """Creates Vectors presentation on the given field.
1460
1461     Arguments:
1462       proxy: the pipeline object, containig data
1463       entity: the entity type from PrsTypeEnum
1464       field_name: the field name
1465       timestamp_nb: the number of time step (1, 2, ...)
1466       scale_factor: scale factor
1467       glyph_pos: the position of glyphs
1468       glyph_type: the type of glyphs
1469       is_colored: this option allows to color the presentation according to
1470       the corresponding data array values
1471       vector_mode: the mode of transformation of vector values
1472       into scalar values, applicable only if the field contains vector values.
1473       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1474
1475     Returns:
1476       Vectors as representation object.
1477
1478     """
1479     proxy.UpdatePipeline()
1480     if entity == EntityType.NODE:
1481         select_cells_with_data(proxy, on_points=[field_name])
1482     else:
1483         select_cells_with_data(proxy, on_cells=[field_name])
1484
1485     # Check vector mode
1486     nb_components = get_nb_components(proxy, entity, field_name)
1487     check_vector_mode(vector_mode, nb_components)
1488
1489     # Get time value
1490     time_value = get_time(proxy, timestamp_nb)
1491
1492     # Set timestamp
1493     pvs.GetRenderView().ViewTime = time_value
1494     pvs.UpdatePipeline(time_value, proxy)
1495
1496     # Extract only groups with data for the field
1497     source = proxy
1498
1499     # Cell centers
1500     if is_data_on_cells(proxy, field_name):
1501         cell_centers = pvs.CellCenters(source)
1502         cell_centers.VertexCells = 1
1503         source = cell_centers
1504
1505     vector_array = field_name
1506     # If the given vector array has only 2 components, add the third one
1507     if nb_components == 2:
1508         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1509         vector_array = calc.ResultArrayName
1510         source = calc
1511
1512     # Glyph
1513     glyph = pvs.Glyph(source)
1514     glyph.Vectors = vector_array
1515     glyph.ScaleMode = 'vector'
1516     #glyph.MaskPoints = 0
1517
1518     # Set glyph type
1519     glyph.GlyphType = glyph_type
1520     if glyph_type == '2D Glyph':
1521         glyph.GlyphType.GlyphType = 'Arrow'
1522     elif glyph_type == 'Cone':
1523         glyph.GlyphType.Resolution = 7
1524         glyph.GlyphType.Height = 2
1525         glyph.GlyphType.Radius = 0.2
1526
1527     # Set glyph position if possible
1528     if glyph.GlyphType.GetProperty("Center"):
1529         if (glyph_pos == GlyphPos.TAIL):
1530             glyph.GlyphType.Center = [0.5, 0.0, 0.0]
1531         elif (glyph_pos == GlyphPos.HEAD):
1532             glyph.GlyphType.Center = [-0.5, 0.0, 0.0]
1533         elif (glyph_pos == GlyphPos.CENTER):
1534             glyph.GlyphType.Center = [0.0, 0.0, 0.0]
1535
1536     if scale_factor is not None:
1537         glyph.ScaleFactor = scale_factor
1538     else:
1539         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1540                                       proxy, entity, field_name)
1541         glyph.ScaleFactor = def_scale
1542
1543     glyph.UpdatePipeline()
1544
1545     # Get Vectors representation object
1546     vectors = pvs.GetRepresentation(glyph)
1547
1548     # Get lookup table
1549     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1550
1551     # Set field range if necessary
1552     data_range = get_data_range(proxy, entity,
1553                                 field_name, vector_mode)
1554     if hasattr(lookup_table,"LockDataRange"):
1555         lookup_table.LockDataRange = 1
1556     elif hasattr(lookup_table,"LockScalarRange"):
1557         lookup_table.LockScalarRange = 1
1558     else:
1559         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1560     
1561     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1562
1563     # Set properties
1564     if (is_colored):
1565         pvs.ColorBy(vectors, (EntityType.get_pvtype(entity), 'GlyphVector'))
1566     else:
1567         pvs.ColorBy(vectors, (EntityType.get_pvtype(entity), None))
1568     vectors.LookupTable = lookup_table
1569
1570     vectors.LineWidth = 1.0
1571
1572     # Set wireframe represenatation mode
1573     vectors.Representation = 'Wireframe'
1574
1575     # Add scalar bar
1576     add_scalar_bar(field_name, nb_components,
1577                    vector_mode, lookup_table, time_value)
1578
1579     return vectors
1580
1581
1582 def DeformedShapeOnField(proxy, entity, field_name,
1583                          timestamp_nb,
1584                          scale_factor=None, is_colored=False,
1585                          vector_mode='Magnitude'):
1586     """Creates Defromed Shape presentation on the given field.
1587
1588     Arguments:
1589       proxy: the pipeline object, containig data
1590       entity: the entity type from PrsTypeEnum
1591       field_name: the field name
1592       timestamp_nb: the number of time step (1, 2, ...)
1593       scale_factor: scale factor of the deformation
1594       is_colored: this option allows to color the presentation according to
1595       the corresponding data array values
1596       vector_mode: the mode of transformation of vector values
1597       into scalar values, applicable only if the field contains vector values.
1598       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1599
1600     Returns:
1601       Defromed Shape as representation object.
1602
1603     """
1604     proxy.UpdatePipeline()
1605     # We don't need mesh parts with no data on them
1606     if entity == EntityType.NODE:
1607         select_cells_with_data(proxy, on_points=[field_name])
1608     else:
1609         select_cells_with_data(proxy, on_cells=[field_name])
1610
1611     # Check vector mode
1612     nb_components = get_nb_components(proxy, entity, field_name)
1613     check_vector_mode(vector_mode, nb_components)
1614
1615     # Get time value
1616     time_value = get_time(proxy, timestamp_nb)
1617
1618     # Set timestamp
1619     pvs.GetRenderView().ViewTime = time_value
1620     pvs.UpdatePipeline(time_value, proxy)
1621
1622     # Do merge
1623     source = pvs.MergeBlocks(proxy)
1624     pvs.UpdatePipeline()
1625
1626     # Cell data to point data
1627     if is_data_on_cells(proxy, field_name):
1628         cell_to_point = pvs.CellDatatoPointData()
1629         cell_to_point.PassCellData = 1
1630         source = cell_to_point
1631
1632     vector_array = field_name
1633     # If the given vector array has only 2 components, add the third one
1634     if nb_components == 2:
1635         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1636         vector_array = calc.ResultArrayName
1637         source = calc
1638
1639     # Warp by vector
1640     warp_vector = pvs.WarpByVector(source)
1641     warp_vector.Vectors = [vector_array]
1642     if scale_factor is not None:
1643         warp_vector.ScaleFactor = scale_factor
1644     else:
1645         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1646                                       proxy, entity, field_name)
1647         warp_vector.ScaleFactor = def_scale
1648
1649     # Get Deformed Shape representation object
1650     defshape = pvs.GetRepresentation(warp_vector)
1651
1652     # Get lookup table
1653     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1654
1655     # Set field range if necessary
1656     data_range = get_data_range(proxy, entity,
1657                                 field_name, vector_mode)
1658     if hasattr(lookup_table,"LockDataRange"):
1659         lookup_table.LockDataRange = 1
1660     elif hasattr(lookup_table,"LockScalarRange"):
1661         lookup_table.LockScalarRange = 1
1662     else:
1663         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1664
1665     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1666
1667     # Set properties
1668     if is_colored:
1669         pvs.ColorBy(defshape, (EntityType.get_pvtype(entity), field_name))
1670     else:
1671         pvs.ColorBy(defshape, (EntityType.get_pvtype(entity), None))
1672     defshape.LookupTable = lookup_table
1673
1674     # Set wireframe represenatation mode
1675     defshape.Representation = 'Wireframe'
1676
1677     # Add scalar bar
1678     add_scalar_bar(field_name, nb_components,
1679                    vector_mode, lookup_table, time_value)
1680
1681     return defshape
1682
1683
1684 def DeformedShapeAndScalarMapOnField(proxy, entity, field_name,
1685                                      timestamp_nb,
1686                                      scale_factor=None,
1687                                      scalar_entity=None,
1688                                      scalar_field_name=None,
1689                                      vector_mode='Magnitude'):
1690     """Creates Defromed Shape And Scalar Map presentation on the given field.
1691
1692     Arguments:
1693       proxy: the pipeline object, containig data
1694       entity: the entity type from PrsTypeEnum
1695       field_name: the field name
1696       timestamp_nb: the number of time step (1, 2, ...)
1697       scale_factor: scale factor of the deformation
1698       scalar_entity: scalar field entity
1699       scalar_field_name: scalar field, i.e. the field for coloring
1700       vector_mode: the mode of transformation of vector values
1701       into scalar values, applicable only if the field contains vector values.
1702       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1703
1704     Returns:
1705       Defromed Shape And Scalar Map as representation object.
1706
1707     """
1708     proxy.UpdatePipeline()
1709     # We don't need mesh parts with no data on them
1710     on_points = []
1711     on_cells = []
1712
1713     if entity == EntityType.NODE:
1714         on_points.append(field_name)
1715     else:
1716         on_cells.append(field_name)
1717
1718     if scalar_entity and scalar_field_name:
1719         if scalar_entity == EntityType.NODE:
1720             on_points.append(scalar_field_name)
1721         else:
1722             on_cells.append(scalar_field_name)
1723
1724     nb_components = get_nb_components(proxy, entity, field_name)
1725
1726     # Select fields
1727     select_cells_with_data(proxy, on_points, on_cells)
1728
1729     # Check vector mode
1730     check_vector_mode(vector_mode, nb_components)
1731
1732     # Get time value
1733     time_value = get_time(proxy, timestamp_nb)
1734
1735     # Set timestamp
1736     pvs.GetRenderView().ViewTime = time_value
1737     pvs.UpdatePipeline(time_value, proxy)
1738
1739     # Set scalar field by default
1740     scalar_field_entity = scalar_entity
1741     scalar_field = scalar_field_name
1742     if (scalar_field_entity is None) or (scalar_field is None):
1743         scalar_field_entity = entity
1744         scalar_field = field_name
1745
1746     # Do merge
1747     source = pvs.MergeBlocks(proxy)
1748     pvs.UpdatePipeline()
1749
1750     # Cell data to point data
1751     if is_data_on_cells(proxy, field_name):
1752         cell_to_point = pvs.CellDatatoPointData(source)
1753         cell_to_point.PassCellData = 1
1754         source = cell_to_point
1755
1756     vector_array = field_name
1757     # If the given vector array has only 2 components, add the third one
1758     if nb_components == 2:
1759         calc = get_add_component_calc(source, EntityType.NODE, field_name)
1760         vector_array = calc.ResultArrayName
1761         source = calc
1762
1763     # Warp by vector
1764     warp_vector = pvs.WarpByVector(source)
1765     warp_vector.Vectors = [vector_array]
1766     if scale_factor is not None:
1767         warp_vector.ScaleFactor = scale_factor
1768     else:
1769         def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE,
1770                                       proxy, entity, field_name)
1771         warp_vector.ScaleFactor = def_scale
1772
1773     # Get Defromed Shape And Scalar Map representation object
1774     defshapemap = pvs.GetRepresentation(warp_vector)
1775
1776     # Get lookup table
1777     lookup_table = get_lookup_table(scalar_field, nb_components, vector_mode)
1778
1779     # Set field range if necessary
1780     data_range = get_data_range(proxy, scalar_field_entity,
1781                                 scalar_field, vector_mode)
1782     if hasattr(lookup_table,"LockDataRange"):
1783         lookup_table.LockDataRange = 1
1784     elif hasattr(lookup_table,"LockScalarRange"):
1785         lookup_table.LockScalarRange = 1
1786     else:
1787         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1788
1789     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1790
1791     # Set properties
1792     pvs.ColorBy(defshapemap, (EntityType.get_pvtype(scalar_field_entity), scalar_field))
1793     defshapemap.LookupTable = lookup_table
1794
1795     # Add scalar bar
1796     add_scalar_bar(field_name, nb_components,
1797                    vector_mode, lookup_table, time_value)
1798
1799     return defshapemap
1800
1801
1802 def Plot3DOnField(proxy, entity, field_name, timestamp_nb,
1803                   orientation=Orientation.AUTO,
1804                   angle1=0, angle2=0,
1805                   position=0.5, is_relative=True,
1806                   scale_factor=None,
1807                   is_contour=False, nb_contours=32,
1808                   vector_mode='Magnitude'):
1809     """Creates Plot 3D presentation on the given field.
1810
1811     Arguments:
1812       proxy: the pipeline object, containig data
1813       entity: the entity type from PrsTypeEnum
1814       field_name: the field name
1815       timestamp_nb: the number of time step (1, 2, ...)
1816       orientation: the cut plane plane orientation in 3D space, if
1817       the input is planar - will not be taken into account
1818       angle1: rotation of the cut plane in 3d space around the first axis
1819       of the selected orientation (X axis for XY, Y axis for YZ,
1820       Z axis for ZX).
1821       The angle of rotation is set in degrees. Acceptable range: [-45, 45].
1822       angle2: rotation of the cut plane in 3d space around the second axis
1823       of the selected orientation. Acceptable range: [-45, 45].
1824       position: position of the cut plane in the object (ranging from 0 to 1).
1825       The value 0.5 corresponds to cutting by halves.
1826       is_relative: defines if the cut plane position is relative or absolute
1827       scale_factor: deformation scale factor
1828       is_contour: if True - Plot 3D will be represented with a set of contours,
1829       otherwise - Plot 3D will be represented with a smooth surface
1830       nb_contours: number of contours, applied if is_contour is True
1831       vector_mode: the mode of transformation of vector values
1832       into scalar values, applicable only if the field contains vector values.
1833       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1834
1835     Returns:
1836       Plot 3D as representation object.
1837
1838     """
1839     proxy.UpdatePipeline()
1840     # We don't need mesh parts with no data on them
1841     if entity == EntityType.NODE:
1842         select_cells_with_data(proxy, on_points=[field_name])
1843     else:
1844         select_cells_with_data(proxy, on_cells=[field_name])
1845
1846     # Check vector mode
1847     nb_components = get_nb_components(proxy, entity, field_name)
1848     check_vector_mode(vector_mode, nb_components)
1849
1850     # Get time value
1851     time_value = get_time(proxy, timestamp_nb)
1852
1853     # Set timestamp
1854     pvs.GetRenderView().ViewTime = time_value
1855     pvs.UpdatePipeline(time_value, proxy)
1856
1857     # Do merge
1858     merge_blocks = pvs.MergeBlocks(proxy)
1859     merge_blocks.UpdatePipeline()
1860
1861     poly_data = None
1862
1863     # Cutting plane
1864
1865     # Define orientation if necessary (auto mode)
1866     plane_orientation = orientation
1867     if (orientation == Orientation.AUTO):
1868         plane_orientation = get_orientation(proxy)
1869
1870     # Get cutting plane normal
1871     normal = None
1872
1873     if (not is_planar_input(proxy)):
1874         normal = get_normal_by_orientation(plane_orientation,
1875                                            radians(angle1), radians(angle2))
1876
1877         # Create slice filter
1878         slice_filter = pvs.Slice(merge_blocks)
1879         slice_filter.SliceType = "Plane"
1880
1881         # Set cutting plane normal
1882         slice_filter.SliceType.Normal = normal
1883
1884         # Set cutting plane position
1885         if (is_relative):
1886             base_position = get_positions(1, normal,
1887                                           get_bounds(proxy), position)
1888             slice_filter.SliceOffsetValues = base_position
1889         else:
1890             slice_filter.SliceOffsetValues = position
1891
1892         slice_filter.UpdatePipeline()
1893         poly_data = slice_filter
1894     else:
1895         normal = get_normal_by_orientation(plane_orientation, 0, 0)
1896
1897     use_normal = 0
1898     # Geometry filter
1899     if not poly_data or poly_data.GetDataInformation().GetNumberOfCells() == 0:
1900         geometry_filter = pvs.GeometryFilter(merge_blocks)
1901         poly_data = geometry_filter
1902         use_normal = 1  # TODO(MZN): workaround
1903
1904     warp_scalar = None
1905     plot3d = None
1906     source = poly_data
1907
1908     if is_data_on_cells(poly_data, field_name):
1909         # Cell data to point data
1910         cell_to_point = pvs.CellDatatoPointData(poly_data)
1911         cell_to_point.PassCellData = 1
1912         source = cell_to_point
1913
1914     scalars = ['POINTS', field_name]
1915
1916     # Transform vector array to scalar array if necessary
1917     if (nb_components > 1):
1918         calc = get_calc_magnitude(source, EntityType.NODE, field_name)
1919         scalars = ['POINTS', calc.ResultArrayName]
1920         source = calc
1921
1922     # Warp by scalar
1923     warp_scalar = pvs.WarpByScalar(source)
1924     warp_scalar.Scalars = scalars
1925     warp_scalar.Normal = normal
1926     warp_scalar.UseNormal = use_normal
1927     if scale_factor is not None:
1928         warp_scalar.ScaleFactor = scale_factor
1929     else:
1930         def_scale = get_default_scale(PrsTypeEnum.PLOT3D,
1931                                       proxy, entity, field_name)
1932         warp_scalar.ScaleFactor = def_scale
1933
1934     warp_scalar.UpdatePipeline()
1935     source = warp_scalar
1936
1937     if (is_contour):
1938         # Contours
1939         contour = pvs.Contour(warp_scalar)
1940         contour.PointMergeMethod = "Uniform Binning"
1941         contour.ContourBy = ['POINTS', field_name]
1942         scalar_range = get_data_range(proxy, entity,
1943                                       field_name, vector_mode)
1944         contour.Isosurfaces = get_contours(scalar_range, nb_contours)
1945         contour.UpdatePipeline()
1946         source = contour
1947
1948     # Get Plot 3D representation object
1949     plot3d = pvs.GetRepresentation(source)
1950
1951     # Get lookup table
1952     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
1953
1954     # Set field range if necessary
1955     data_range = get_data_range(proxy, entity,
1956                                 field_name, vector_mode)
1957     if hasattr(lookup_table,"LockDataRange"):
1958         lookup_table.LockDataRange = 1
1959     elif hasattr(lookup_table,"LockScalarRange"):
1960         lookup_table.LockScalarRange = 1
1961     else:
1962         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
1963
1964     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
1965
1966     # Set properties
1967     pvs.ColorBy(plot3d, (EntityType.get_pvtype(entity), field_name))
1968     plot3d.LookupTable = lookup_table
1969
1970     # Add scalar bar
1971     add_scalar_bar(field_name, nb_components,
1972                    vector_mode, lookup_table, time_value)
1973
1974     return plot3d
1975
1976
1977 def IsoSurfacesOnField(proxy, entity, field_name, timestamp_nb,
1978                        custom_range=None, nb_surfaces=10,
1979                        is_colored=True, color=None, vector_mode='Magnitude'):
1980     """Creates Iso Surfaces presentation on the given field.
1981
1982     Arguments:
1983       proxy: the pipeline object, containig data
1984       entity: the entity type from PrsTypeEnum
1985       field_name: the field name
1986       timestamp_nb: the number of time step (1, 2, ...)
1987       custom_range: scalar range, if undefined the source range will be applied
1988       nb_surfaces: number of surfaces, which will be generated
1989       is_colored: this option allows to color the presentation according to
1990       the corresponding data array values. If False - the presentation will
1991       be one-coloured.
1992       color: defines the presentation color as [R, G, B] triple. Taken into
1993       account only if is_colored is False.
1994       vector_mode: the mode of transformation of vector values
1995       into scalar values, applicable only if the field contains vector values.
1996       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
1997
1998     Returns:
1999       Iso Surfaces as representation object.
2000
2001     """
2002     proxy.UpdatePipeline()
2003     # We don't need mesh parts with no data on them
2004     if entity == EntityType.NODE:
2005         select_cells_with_data(proxy, on_points=[field_name])
2006     else:
2007         select_cells_with_data(proxy, on_cells=[field_name])
2008
2009     # Check vector mode
2010     nb_components = get_nb_components(proxy, entity, field_name)
2011     check_vector_mode(vector_mode, nb_components)
2012
2013     # Get time value
2014     time_value = get_time(proxy, timestamp_nb)
2015
2016     # Set timestamp
2017     pvs.GetRenderView().ViewTime = time_value
2018     pvs.UpdatePipeline(time_value, proxy)
2019
2020     # Do merge
2021     source = pvs.MergeBlocks(proxy)
2022     pvs.UpdatePipeline()
2023
2024     # Transform cell data into point data if necessary
2025     if is_data_on_cells(proxy, field_name):
2026         cell_to_point = pvs.CellDatatoPointData(source)
2027         cell_to_point.PassCellData = 1
2028         source = cell_to_point
2029
2030     contour_by = ['POINTS', field_name]
2031
2032     # Transform vector array to scalar array if necessary
2033     if (nb_components > 1):
2034         calc = get_calc_magnitude(source, EntityType.NODE, field_name)
2035         contour_by = ['POINTS', calc.ResultArrayName]
2036         source = calc
2037
2038     # Contour filter settings
2039     contour = pvs.Contour(source)
2040     contour.ComputeScalars = 1
2041     contour.ContourBy = contour_by
2042
2043     # Specify the range
2044     scalar_range = custom_range
2045     if (scalar_range is None):
2046         scalar_range = get_data_range(proxy, entity,
2047                                       field_name, cut_off=True)
2048         
2049     # Get contour values for the range
2050     surfaces = get_contours(scalar_range, nb_surfaces)
2051
2052     # Set contour values
2053     contour.Isosurfaces = surfaces
2054
2055     # Get Iso Surfaces representation object
2056     isosurfaces = pvs.GetRepresentation(contour)
2057
2058     # Get lookup table
2059     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2060
2061     # Set field range if necessary
2062     data_range = get_data_range(proxy, entity,
2063                                 field_name, vector_mode)
2064     if hasattr(lookup_table,"LockDataRange"):
2065         lookup_table.LockDataRange = 1
2066     elif hasattr(lookup_table,"LockScalarRange"):
2067         lookup_table.LockScalarRange = 1
2068     else:
2069         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
2070
2071     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2072
2073     # Set display properties
2074     if (is_colored):
2075         pvs.ColorBy(isosurfaces, (EntityType.get_pvtype(entity), field_name))
2076     else:
2077         pvs.ColorBy(isosurfaces, (EntityType.get_pvtype(entity), None))
2078         if color:
2079             isosurfaces.DiffuseColor = color
2080     isosurfaces.LookupTable = lookup_table
2081
2082     # Add scalar bar
2083     add_scalar_bar(field_name, nb_components,
2084                    vector_mode, lookup_table, time_value)
2085
2086     return isosurfaces
2087
2088
2089 def GaussPointsOnField(proxy, entity, field_name,
2090                        timestamp_nb,
2091                        is_deformed=True, scale_factor=None,
2092                        is_colored=True, color=None,
2093                        primitive=GaussType.SPRITE,
2094                        is_proportional=True,
2095                        max_pixel_size=256,
2096                        multiplier=None, vector_mode='Magnitude'):
2097     """Creates Gauss Points on the given field.
2098
2099     Arguments:
2100
2101     proxy: the pipeline object, containig data
2102     entity: the field entity type from PrsTypeEnum
2103     field_name: the field name
2104     timestamp_nb: the number of time step (1, 2, ...)
2105     is_deformed: defines whether the Gauss Points will be deformed or not
2106     scale_factor -- the scale factor for deformation. Will be taken into
2107     account only if is_deformed is True.
2108     If not passed by user, default scale will be computed.
2109     is_colored -- defines whether the Gauss Points will be multicolored,
2110     using the corresponding data values
2111     color: defines the presentation color as [R, G, B] triple. Taken into
2112     account only if is_colored is False.
2113     primitive: primitive type from GaussType
2114     is_proportional: if True, the size of primitives will depends on
2115     the gauss point value
2116     max_pixel_size: the maximum sizr of the Gauss Points primitive in pixels
2117     multiplier: coefficient between data values and the size of primitives
2118     If not passed by user, default scale will be computed.
2119     vector_mode: the mode of transformation of vector values into
2120     scalar values, applicable only if the field contains vector values.
2121     Possible modes: 'Magnitude' - vector module;
2122     'X', 'Y', 'Z' - vector components.
2123
2124     Returns:
2125       Gauss Points as representation object.
2126
2127     """
2128     proxy.UpdatePipeline()
2129     # We don't need mesh parts with no data on them
2130     on_gauss = select_cells_with_data(proxy, on_gauss=[field_name])
2131     if not on_gauss:
2132         if entity == EntityType.NODE:
2133             select_cells_with_data(proxy, on_points=[field_name])
2134         else:
2135             select_cells_with_data(proxy, on_cells=[field_name])
2136
2137     # Check vector mode
2138     nb_components = get_nb_components(proxy, entity, field_name)
2139     check_vector_mode(vector_mode, nb_components)
2140
2141     # Get time value
2142     time_value = get_time(proxy, timestamp_nb)
2143
2144     # Set timestamp
2145     pvs.GetRenderView().ViewTime = time_value
2146     pvs.UpdatePipeline(time_value, proxy)
2147
2148     source = proxy
2149
2150     # If no quadrature point array is passed, use cell centers
2151     if on_gauss:
2152         generate_qp = pvs.GenerateQuadraturePoints(source)
2153         generate_qp.QuadratureSchemeDef = ['CELLS', 'ELGA@0']
2154         source = generate_qp
2155     else:
2156         # Cell centers
2157         cell_centers = pvs.CellCenters(source)
2158         cell_centers.VertexCells = 1
2159         source = cell_centers
2160
2161     source.UpdatePipeline()
2162
2163     # Check if deformation enabled
2164     if is_deformed and nb_components > 1:
2165         vector_array = field_name
2166         # If the given vector array has only 2 components, add the third one
2167         if nb_components == 2:
2168             calc = get_add_component_calc(source, EntityType.NODE, field_name)
2169             vector_array = calc.ResultArrayName
2170             source = calc
2171
2172         # Warp by vector
2173         warp_vector = pvs.WarpByVector(source)
2174         warp_vector.Vectors = [vector_array]
2175         if scale_factor is not None:
2176             warp_vector.ScaleFactor = scale_factor
2177         else:
2178             def_scale = get_default_scale(PrsTypeEnum.DEFORMEDSHAPE, proxy,
2179                                           entity, field_name)
2180             warp_vector.ScaleFactor = def_scale
2181         warp_vector.UpdatePipeline()
2182         source = warp_vector
2183
2184     # Get Gauss Points representation object
2185     gausspnt = pvs.GetRepresentation(source)
2186
2187     # Get lookup table
2188     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2189
2190     # Set field range if necessary
2191     data_range = get_data_range(proxy, entity,
2192                                 field_name, vector_mode)
2193     if hasattr(lookup_table,"LockDataRange"):
2194         lookup_table.LockDataRange = 1
2195     elif hasattr(lookup_table,"LockScalarRange"):
2196         lookup_table.LockScalarRange = 1
2197     else:
2198         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
2199
2200     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2201
2202     # Set display properties
2203     if is_colored:
2204         pvs.ColorBy(gausspnt, (EntityType.get_pvtype(entity), field_name))
2205     else:
2206         pvs.ColorBy(gausspnt, (EntityType.get_pvtype(entity), None))
2207         if color:
2208             gausspnt.DiffuseColor = color
2209
2210     gausspnt.LookupTable = lookup_table
2211
2212     # Add scalar bar
2213     add_scalar_bar(field_name, nb_components,
2214                    vector_mode, lookup_table, time_value)
2215
2216     # Set point sprite representation
2217     gausspnt.Representation = 'Point Sprite'
2218
2219     # Point sprite settings
2220     gausspnt.InterpolateScalarsBeforeMapping = 0
2221     gausspnt.MaxPixelSize = max_pixel_size
2222
2223     # Render mode
2224     gausspnt.RenderMode = GaussType.get_mode(primitive)
2225
2226     #if primitive == GaussType.SPRITE:
2227         # Set texture
2228         # TODO(MZN): replace with pvsimple high-level interface
2229     #    texture = sm.CreateProxy("textures", "SpriteTexture")
2230     #    alphamprop = texture.GetProperty("AlphaMethod")
2231     #    alphamprop.SetElement(0, 2)  # Clamp
2232     #    alphatprop = texture.GetProperty("AlphaThreshold")
2233     #    alphatprop.SetElement(0, 63)
2234     #    maxprop = texture.GetProperty("Maximum")
2235     #    maxprop.SetElement(0, 255)
2236     #    texture.UpdateVTKObjects()
2237
2238     #    gausspnt.Texture = texture
2239         #gausspnt.Texture.AlphaMethod = 'Clamp'
2240         #gausspnt.Texture.AlphaThreshold = 63
2241         #gausspnt.Texture.Maximum= 255
2242
2243     # Proportional radius
2244     gausspnt.RadiusUseScalarRange = 0
2245     gausspnt.RadiusIsProportional = 0
2246
2247     if is_proportional:
2248         mult = multiplier
2249         if mult is None and data_range[1] != 0:
2250             mult = abs(0.1 / data_range[1])
2251
2252         gausspnt.RadiusScalarRange = data_range
2253         gausspnt.RadiusTransferFunctionEnabled = 1
2254         gausspnt.RadiusMode = 'Scalar'
2255         gausspnt.RadiusArray = ['POINTS', field_name]
2256         if nb_components > 1:
2257             v_comp = get_vector_component(vector_mode)
2258             gausspnt.RadiusVectorComponent = v_comp
2259         gausspnt.RadiusTransferFunctionMode = 'Table'
2260         gausspnt.RadiusScalarRange = data_range
2261         gausspnt.RadiusUseScalarRange = 1
2262         if mult is not None:
2263             gausspnt.RadiusIsProportional = 1
2264             gausspnt.RadiusProportionalFactor = mult
2265     else:
2266         gausspnt.RadiusTransferFunctionEnabled = 0
2267         gausspnt.RadiusMode = 'Constant'
2268         gausspnt.RadiusArray = ['POINTS', 'Constant Radius']
2269
2270     return gausspnt
2271
2272 def GaussPointsOnField1(proxy, entity, field_name,
2273                         timestamp_nb,
2274                         is_colored=True, color=None,
2275                         primitive=GaussType.SPHERE,
2276                         is_proportional=True,
2277                         max_pixel_size=256,
2278                         multiplier=None,
2279                         vector_mode='Magnitude'):
2280     """Creates Gauss Points on the given field. Use GaussPoints() Paraview interface.
2281
2282     Arguments:
2283     proxy: the pipeline object, containig data
2284     entity: the field entity type from PrsTypeEnum
2285     field_name: the field name
2286     timestamp_nb: the number of time step (1, 2, ...)
2287     is_colored -- defines whether the Gauss Points will be multicolored,
2288     using the corresponding data values
2289     color: defines the presentation color as [R, G, B] triple. Taken into
2290     account only if is_colored is False.
2291     primitive: primitive type from GaussType
2292     is_proportional: if True, the size of primitives will depends on
2293     the gauss point value
2294     max_pixel_size: the maximum sizr of the Gauss Points primitive in pixels
2295     multiplier: coefficient between data values and the size of primitives
2296     If not passed by user, default scale will be computed.
2297     vector_mode: the mode of transformation of vector values into
2298     scalar values, applicable only if the field contains vector values.
2299     Possible modes: 'Magnitude' - vector module;
2300     'X', 'Y', 'Z' - vector components.
2301
2302     Returns:
2303       Gauss Points as representation object.
2304
2305     """
2306     proxy.UpdatePipeline()
2307     select_cells_with_data(proxy, on_gauss=[field_name])
2308
2309     nb_components = get_nb_components(proxy, entity, field_name)
2310
2311     # Get time value
2312     time_value = get_time(proxy, timestamp_nb)
2313
2314     # Set timestamp
2315     pvs.GetRenderView().ViewTime = time_value
2316     proxy.UpdatePipeline(time=time_value)
2317
2318     # Create Gauss Points object
2319     source = pvs.GaussPoints(proxy)
2320     source.UpdatePipeline()
2321
2322     # Get Gauss Points representation object
2323     gausspnt = pvs.GetRepresentation(source)
2324
2325     # Get lookup table
2326     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2327
2328     # Set field range if necessary
2329     data_range = get_data_range(proxy, entity,
2330                                 field_name, vector_mode)
2331     if hasattr(lookup_table,"LockDataRange"):
2332         lookup_table.LockDataRange = 1
2333     elif hasattr(lookup_table,"LockScalarRange"):
2334         lookup_table.LockScalarRange = 1
2335     else:
2336         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
2337
2338     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2339
2340     # Set display properties
2341     if is_colored:
2342         pvs.ColorBy(gausspnt, (EntityType.get_pvtype(entity), field_name))
2343     else:
2344         pvs.ColorBy(gausspnt, (EntityType.get_pvtype(entity), None))
2345         if color:
2346             gausspnt.DiffuseColor = color
2347
2348     gausspnt.LookupTable = lookup_table
2349
2350     # Add scalar bar
2351     add_scalar_bar(field_name, nb_components,
2352                    vector_mode, lookup_table, time_value)
2353
2354     # Set point sprite representation
2355     gausspnt.Representation = 'Point Sprite'
2356
2357     # Point sprite settings
2358     gausspnt.InterpolateScalarsBeforeMapping = 0
2359     gausspnt.MaxPixelSize = max_pixel_size
2360
2361     # Render mode
2362     gausspnt.RenderMode = GaussType.get_mode(primitive)
2363
2364     #if primitive == GaussType.SPRITE:
2365         # Set texture
2366         # TODO(MZN): replace with pvsimple high-level interface
2367     #    texture = sm.CreateProxy("textures", "SpriteTexture")
2368     #    alphamprop = texture.GetProperty("AlphaMethod")
2369     #    alphamprop.SetElement(0, 2)  # Clamp
2370     #    alphatprop = texture.GetProperty("AlphaThreshold")
2371     #    alphatprop.SetElement(0, 63)
2372     #    maxprop = texture.GetProperty("Maximum")
2373     #    maxprop.SetElement(0, 255)
2374     #    texture.UpdateVTKObjects()
2375
2376     #    gausspnt.Texture = texture
2377         #gausspnt.Texture.AlphaMethod = 'Clamp'
2378         #gausspnt.Texture.AlphaThreshold = 63
2379         #gausspnt.Texture.Maximum= 255
2380
2381     # Proportional radius
2382     gausspnt.RadiusUseScalarRange = 0
2383     gausspnt.RadiusIsProportional = 0
2384
2385     if is_proportional:
2386         mult = multiplier
2387         if mult is None and data_range[1] != 0:
2388             mult = abs(0.1 / data_range[1])
2389
2390         gausspnt.RadiusScalarRange = data_range
2391         gausspnt.RadiusTransferFunctionEnabled = 1
2392         gausspnt.RadiusMode = 'Scalar'
2393         gausspnt.RadiusArray = ['POINTS', field_name]
2394         if nb_components > 1:
2395             v_comp = get_vector_component(vector_mode)
2396             gausspnt.RadiusVectorComponent = v_comp
2397         gausspnt.RadiusTransferFunctionMode = 'Table'
2398         gausspnt.RadiusScalarRange = data_range
2399         gausspnt.RadiusUseScalarRange = 1
2400         if mult is not None:
2401             gausspnt.RadiusIsProportional = 1
2402             gausspnt.RadiusProportionalFactor = mult
2403     else:
2404         gausspnt.RadiusTransferFunctionEnabled = 0
2405         gausspnt.RadiusMode = 'Constant'
2406         gausspnt.RadiusArray = ['POINTS', 'Constant Radius']
2407
2408     return gausspnt
2409
2410 def StreamLinesOnField(proxy, entity, field_name, timestamp_nb,
2411                        direction='BOTH', is_colored=False, color=None,
2412                        vector_mode='Magnitude'):
2413     """Creates Stream Lines presentation on the given field.
2414
2415     Arguments:
2416       proxy: the pipeline object, containig data
2417       entity: the entity type from PrsTypeEnum
2418       field_name: the field name
2419       timestamp_nb: the number of time step (1, 2, ...)
2420       direction: the stream lines direction ('FORWARD', 'BACKWARD' or 'BOTH')
2421       is_colored: this option allows to color the presentation according to
2422       the corresponding data values. If False - the presentation will
2423       be one-coloured.
2424       color: defines the presentation color as [R, G, B] triple. Taken into
2425       account only if is_colored is False.
2426       vector_mode: the mode of transformation of vector values
2427       into scalar values, applicable only if the field contains vector values.
2428       Possible modes: 'Magnitude', 'X', 'Y' or 'Z'.
2429
2430     Returns:
2431       Stream Lines as representation object.
2432
2433     """
2434     proxy.UpdatePipeline()
2435     # We don't need mesh parts with no data on them
2436     if entity == EntityType.NODE:
2437         select_cells_with_data(proxy, on_points=[field_name])
2438     else:
2439         select_cells_with_data(proxy, on_cells=[field_name])
2440
2441     # Check vector mode
2442     nb_components = get_nb_components(proxy, entity, field_name)
2443     check_vector_mode(vector_mode, nb_components)
2444
2445     # Get time value
2446     time_value = get_time(proxy, timestamp_nb)
2447
2448     # Set timestamp
2449     pvs.GetRenderView().ViewTime = time_value
2450     pvs.UpdatePipeline(time_value, proxy)
2451
2452     # Do merge
2453     source = pvs.MergeBlocks(proxy)
2454     pvs.UpdatePipeline()
2455
2456     # Cell data to point data
2457     if is_data_on_cells(proxy, field_name):
2458         cell_to_point = pvs.CellDatatoPointData(source)
2459         cell_to_point.PassCellData = 1
2460         pvs.UpdatePipeline()
2461         source = cell_to_point
2462
2463     vector_array = field_name
2464     # If the given vector array has only 2 components, add the third one
2465     if nb_components == 2:
2466         calc = get_add_component_calc(source, EntityType.NODE, field_name)
2467         vector_array = calc.ResultArrayName
2468         pvs.UpdatePipeline()
2469         source = calc
2470
2471     # Stream Tracer
2472     stream = pvs.StreamTracer(source)
2473     stream.SeedType = "Point Source"
2474     stream.Vectors = ['POINTS', vector_array]
2475     stream.IntegrationDirection = direction
2476     stream.IntegratorType = 'Runge-Kutta 2'
2477     stream.SeedType = 'High Resolution Line Source'
2478     stream.UpdatePipeline()
2479
2480     # Get Stream Lines representation object
2481     if is_empty(stream):
2482         return None
2483     streamlines = pvs.GetRepresentation(stream)
2484
2485     # Get lookup table
2486     lookup_table = get_lookup_table(field_name, nb_components, vector_mode)
2487
2488     # Set field range if necessary
2489     data_range = get_data_range(proxy, entity,
2490                                 field_name, vector_mode)
2491     if hasattr(lookup_table,"LockDataRange"):
2492         lookup_table.LockDataRange = 1
2493     elif hasattr(lookup_table,"LockScalarRange"):
2494         lookup_table.LockScalarRange = 1
2495     else:
2496         raise RuntimeError("Object %s has no 'LockDataRange' or 'LockScalarRange' attribute!"%(lookup_table))
2497
2498     lookup_table.RGBPoints = [data_range[0], 0, 0, 1, data_range[1], 1, 0, 0]
2499
2500     # Set properties
2501     if is_colored:
2502         pvs.ColorBy(streamlines, (EntityType.get_pvtype(entity), field_name))
2503     else:
2504         pvs.ColorBy(streamlines, (EntityType.get_pvtype(entity), None))
2505         if color:
2506             streamlines.DiffuseColor = color
2507
2508     streamlines.LookupTable = lookup_table
2509
2510     # Add scalar bar
2511     add_scalar_bar(field_name, nb_components,
2512                    vector_mode, lookup_table, time_value)
2513
2514     return streamlines
2515
2516
2517 def MeshOnEntity(proxy, mesh_name, entity):
2518     """Creates submesh of the entity type for the mesh.
2519
2520     Arguments:
2521       proxy -- the pipeline object, containig data
2522       mesh_name -- the full or short name of mesh field
2523
2524     Returns:
2525       Submesh as representation object of the given source.
2526
2527     """
2528     proxy.UpdatePipeline()
2529     mesh_full_name = None
2530     aList = mesh_name.split('/')
2531     if len(aList) >= 2:
2532         mesh_full_name = mesh_name
2533     else:
2534         mesh_full_name = find_mesh_full_name(proxy, mesh_name)
2535     if not mesh_full_name:
2536         raise RuntimeError, "The given mesh name was not found"
2537     # Select only the given mesh
2538     proxy.AllArrays = [mesh_full_name]
2539     proxy.UpdatePipeline()
2540
2541     # Get representation object if the submesh is not empty
2542     prs = None
2543     if (proxy.GetDataInformation().GetNumberOfPoints() or
2544         proxy.GetDataInformation().GetNumberOfCells()):
2545         my_view = pvs.GetRenderView()
2546         prs = pvs.GetRepresentation(proxy, view=my_view)
2547         pvs.ColorBy(prs, (EntityType.get_pvtype(entity), None))
2548
2549     return prs
2550
2551
2552 def MeshOnGroup(proxy, extrGroups, group_name):
2553     """Creates submesh on the group.
2554
2555     Arguments:
2556       proxy -- the pipeline object, containig data
2557       group_name -- the full group name
2558       extrGroups -- all extracted groups object
2559
2560     Returns:
2561       Representation object of the given source with single group
2562       selected.
2563
2564     """
2565     proxy.UpdatePipeline()
2566     # Deselect all groups
2567     extrGroups.AllGroups = []
2568     extrGroups.UpdatePipelineInformation()
2569     # Select only the group with the given name
2570     extrGroups.AllGroups = [group_name]
2571     extrGroups.UpdatePipelineInformation()
2572
2573     # Get representation object if the submesh is not empty
2574     prs = None
2575
2576     # Check if the group was set
2577     if len(extrGroups.AllGroups) == 1 and \
2578        extrGroups.AllGroups[0] == group_name:
2579         # Check if the submesh is not empty
2580         nb_points = proxy.GetDataInformation().GetNumberOfPoints()
2581         nb_cells = proxy.GetDataInformation().GetNumberOfCells()
2582
2583         if nb_points or nb_cells:
2584 #            prs = pvs.GetRepresentation(proxy)
2585             prs = pvs.Show()
2586             prs.ColorArrayName = (None, '')
2587             display_only(prs)
2588
2589     return prs
2590
2591
2592 def CreatePrsForFile(file_name, prs_types,
2593                      picture_dir, picture_ext):
2594     """Build presentations of the given types for the file.
2595
2596     Build presentations for all fields on all timestamps.
2597
2598     Arguments:
2599       file_name: full path to the MED file
2600       prs_types: the list of presentation types to build
2601       picture_dir: the directory path for saving snapshots
2602       picture_ext: graphics files extension (determines file type)
2603
2604     """
2605     # Import MED file
2606     print "Import " + file_name.split(os.sep)[-1] + "..."
2607
2608     try:
2609         proxy = pvs.MEDReader(FileName=file_name)
2610         if proxy is None:
2611             print "FAILED"
2612         else:
2613             #proxy.UpdatePipeline()
2614             print "OK"
2615     except:
2616         print "FAILED"
2617     else:
2618         # Get view
2619         view = pvs.GetRenderView()
2620         time_value = get_time(proxy, 0)
2621         view.ViewTime = time_value
2622         pvs.UpdatePipeline(time=time_value, proxy=proxy)
2623
2624         # Create required presentations for the proxy
2625         CreatePrsForProxy(proxy, view, prs_types,
2626                           picture_dir, picture_ext)
2627
2628 def CreatePrsForProxy(proxy, view, prs_types, picture_dir, picture_ext):
2629     """Build presentations of the given types for all fields of the proxy.
2630
2631     Save snapshots in graphics files (type depends on the given extension).
2632     Stores the files in the given directory.
2633
2634     Arguments:
2635       proxy: the pipeline object, containig data
2636       view: the render view
2637       prs_types: the list of presentation types to build
2638       picture_dir: the directory path for saving snapshots
2639       picture_ext: graphics files extension (determines file type)
2640
2641     """
2642     proxy.UpdatePipeline()
2643     # List of the field names
2644     fields_info = proxy.GetProperty("FieldsTreeInfo")[::2]
2645
2646     # Add path separator to the end of picture path if necessery
2647     if not picture_dir.endswith(os.sep):
2648         picture_dir += os.sep
2649
2650     # Mesh Presentation
2651     if PrsTypeEnum.MESH in prs_types:
2652         # Iterate on meshes
2653         mesh_names = get_mesh_full_names(proxy)
2654         for mesh_name in mesh_names:
2655             # Build mesh field presentation
2656             print "Creating submesh for '" + get_field_short_name(mesh_name) + "' mesh... "
2657             prs = MeshOnEntity(proxy, mesh_name, None)
2658             if prs is None:
2659                 print "FAILED"
2660                 continue
2661             else:
2662                 print "OK"
2663             # Construct image file name
2664             pic_name = picture_dir + get_field_short_name(mesh_name) + "." + picture_ext
2665
2666             # Show and dump the presentation into a graphics file
2667             process_prs_for_test(prs, view, pic_name, False)
2668
2669             # Create Mesh presentation. Build all groups.
2670             extGrp = pvs.ExtractGroup()
2671             extGrp.UpdatePipelineInformation()
2672             if if_possible(proxy, None, None, PrsTypeEnum.MESH, extGrp):
2673                 for group in get_group_names(extGrp):
2674                     print "Creating submesh on group " + get_group_short_name(group) + "... "
2675                     prs = MeshOnGroup(proxy, extGrp, group)
2676                     if prs is None:
2677                         print "FAILED"
2678                         continue
2679                     else:
2680                         print "OK"
2681                     # Construct image file name
2682                     pic_name = picture_dir + get_group_short_name(group) + "." + picture_ext
2683
2684                     # Show and dump the presentation into a graphics file
2685                     process_prs_for_test(prs, view, pic_name, False)
2686
2687     # Presentations on fields
2688     for field in fields_info:
2689         field_name = get_field_short_name(field)
2690         # Ignore mesh presentation
2691         if field_name == get_field_mesh_name(field):
2692             continue
2693         field_entity = get_field_entity(field)
2694         # Select only the current field:
2695         # necessary for getting the right timestamps
2696         proxy.AllArrays = [field]
2697         proxy.UpdatePipeline()
2698
2699         # Get timestamps
2700         timestamps = proxy.TimestepValues.GetData()
2701
2702         for prs_type in prs_types:
2703             # Ignore mesh presentation
2704             if prs_type == PrsTypeEnum.MESH:
2705                 continue
2706
2707             # Get name of presentation type
2708             prs_name = PrsTypeEnum.get_name(prs_type)
2709
2710             # Build the presentation if possible
2711             possible = if_possible(proxy, field_name,
2712                                    field_entity, prs_type)
2713             if possible:
2714                 # Presentation type for graphics file name
2715                 f_prs_type = prs_name.replace(' ', '').upper()
2716
2717                 for timestamp_nb in xrange(1, len(timestamps) + 1):
2718                     time = timestamps[timestamp_nb - 1]
2719                     if (time == 0.0):
2720                         scalar_range = get_data_range(proxy, field_entity,
2721                                                       field_name, cut_off=True)
2722                         # exclude time stamps with null lenght of scalar range
2723                         if (scalar_range[0] == scalar_range[1]):
2724                             continue
2725                     print "Creating " + prs_name + " on " + field_name + ", time = " + str(time) + "... "
2726                     try:
2727                         prs = create_prs(prs_type, proxy,
2728                                          field_entity, field_name, timestamp_nb)
2729                     except ValueError:
2730                         """ This exception comes from get_nb_components(...) function.
2731                             The reason of exception is an implementation of MEDReader
2732                             activating the first leaf when reading MED file (refer to
2733                             MEDFileFieldRepresentationTree::activateTheFirst() and
2734                             MEDFileFieldRepresentationTree::getTheSingleActivated(...) methods).
2735                         """
2736                         print "ValueError exception is catched"
2737                         continue
2738                     if prs is None:
2739                         print "FAILED"
2740                         continue
2741                     else:
2742                         print "OK"
2743
2744                     # Construct image file name
2745                     pic_name = picture_dir + field_name + "_" + str(time) + "_" + f_prs_type + "." + picture_ext
2746
2747                     # Show and dump the presentation into a graphics file
2748                     process_prs_for_test(prs, view, pic_name)
2749     return