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