1 # -*- coding: iso-8859-1 -*-
2 # Copyright (C) 2016 CEA/DEN, EDF R&D
4 # This library is free software; you can redistribute it and/or
5 # modify it under the terms of the GNU Lesser General Public
6 # License as published by the Free Software Foundation; either
7 # version 2.1 of the License, or (at your option) any later version.
9 # This library is distributed in the hope that it will be useful,
10 # but WITHOUT ANY WARRANTY; without even the implied warranty of
11 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 # Lesser General Public License for more details.
14 # You should have received a copy of the GNU Lesser General Public
15 # License along with this library; if not, write to the Free Software
16 # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18 # See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 # Author : Anthony Geay (EDF R&D)
22 # non regression test that emulates https://ageay@git.salome-platform.org/gitpub/samples/datafiles.git Med/ResOK_0000.med
23 # This test point error during commit efd9331a9455785d0f04b75 in PARAVIS
25 from MEDLoader import *
26 fname="testMEDReader20.med"
27 png="testMEDReader20.png"
29 arrX=DataArrayDouble(nb+1) ; arrX.iota()
30 arrY=DataArrayDouble([0.,1.])
31 m=MEDCouplingCMesh() ; m.setCoords(arrX,arrY) ; m=m.buildUnstructured(); m.setName("mesh") ; m.simplexize(0)
32 mm=MEDFileUMesh() ; mm[0]=m
33 m1=m.computeSkin() ; mm[-1]=m1
35 f0=DataArrayInt(m1.getNumberOfCells()) ; f0.iota() ; mm.setFamilyFieldArr(-1,f0)
36 f1=DataArrayInt(m1.getNumberOfNodes()) ; f1.iota() ; mm.setFamilyFieldArr(1,f1) # <- very important the bug can be shown here
38 nbCells=m1.getNumberOfCells() ; nbNodes=m.getNumberOfNodes()
41 f=MEDCouplingFieldDouble(ON_CELLS) ; f.setMesh(m)
43 arr=DataArrayInt(2*nb) ; arr.iota(i) ; arr%=nb ; arr=arr.convertToDblArr()
44 f.setArray(arr) ; f.setTime(float(i),i,0)
45 WriteFieldUsingAlreadyWrittenMesh(fname,f)
47 f=MEDCouplingFieldDouble(ON_CELLS) ; f.setMesh(m1)
49 arr=DataArrayInt(nbCells) ; arr.iota(i) ; arr%=nbCells ; arr=arr.convertToDblArr()
50 f.setArray(arr) ; f.setTime(float(i),i,0)
51 WriteFieldUsingAlreadyWrittenMesh(fname,f)
53 f=MEDCouplingFieldDouble(ON_NODES) ; f.setMesh(m)
54 f.setName("FieldNode")
55 arr=DataArrayDouble(nbNodes) ; arr[:]=float(i)
56 f.setArray(arr) ; f.setTime(float(i),i,0)
57 WriteFieldUsingAlreadyWrittenMesh(fname,f)
60 from paraview.simple import *
61 #### disable automatic camera reset on 'Show'
62 paraview.simple._DisableFirstRenderCameraReset()
64 # create a new 'MED Reader'
65 testMEDReader20med = MEDReader(FileName=fname)
66 testMEDReader20med.AllArrays = ['TS0/mesh/ComSup0/Field@@][@@P0']
67 testMEDReader20med.AllTimeSteps = ['0000', '0001', '0002', '0003', '0004']
70 animationScene1 = GetAnimationScene()
72 # update animation scene based on data timesteps
73 animationScene1.UpdateAnimationUsingDataTimeSteps()
76 renderView1 = GetActiveViewOrCreate('RenderView')
77 # uncomment following to set a specific view size
78 # renderView1.ViewSize = [610, 477]
81 testMEDReader20medDisplay = Show(testMEDReader20med, renderView1)
82 # trace defaults for the display properties.
83 testMEDReader20medDisplay.ColorArrayName = [None, '']
84 testMEDReader20medDisplay.GlyphType = 'Arrow'
85 testMEDReader20medDisplay.ScalarOpacityUnitDistance = 4.664739046219201
86 testMEDReader20medDisplay.SelectUncertaintyArray = [None, '']
87 testMEDReader20medDisplay.UncertaintyTransferFunction = 'PiecewiseFunction'
88 testMEDReader20medDisplay.OpacityArray = [None, '']
89 testMEDReader20medDisplay.RadiusArray = [None, '']
90 testMEDReader20medDisplay.RadiusRange = [0.0, 10.0]
91 testMEDReader20medDisplay.ConstantRadius = 10.0
92 testMEDReader20medDisplay.PointSpriteDefaultsInitialized = 1
93 testMEDReader20medDisplay.SelectInputVectors = [None, '']
94 testMEDReader20medDisplay.WriteLog = ''
96 # reset view to fit data
97 renderView1.ResetCamera()
99 #changing interaction mode based on data extents
100 renderView1.InteractionMode = '2D'
101 renderView1.CameraPosition = [5.0, 0.5, 10000.0]
102 renderView1.CameraFocalPoint = [5.0, 0.5, 0.0]
104 # set scalar coloring
105 ColorBy(testMEDReader20medDisplay, ('FIELD', 'vtkBlockColors'))
107 # show color bar/color legend
108 testMEDReader20medDisplay.SetScalarBarVisibility(renderView1, True)
110 # get color transfer function/color map for 'vtkBlockColors'
111 vtkBlockColorsLUT = GetColorTransferFunction('vtkBlockColors')
113 # get opacity transfer function/opacity map for 'vtkBlockColors'
114 vtkBlockColorsPWF = GetOpacityTransferFunction('vtkBlockColors')
116 # set scalar coloring
117 ColorBy(testMEDReader20medDisplay, ('CELLS', 'Field'))
119 # rescale color and/or opacity maps used to include current data range
120 testMEDReader20medDisplay.RescaleTransferFunctionToDataRange(True)
122 # show color bar/color legend
123 testMEDReader20medDisplay.SetScalarBarVisibility(renderView1, True)
125 # get color transfer function/color map for 'Field'
126 fieldLUT = GetColorTransferFunction('Field')
128 # get opacity transfer function/opacity map for 'Field'
129 fieldPWF = GetOpacityTransferFunction('Field')
131 animationScene1.GoToNext() # <- very important to see the bug play with time steps...
132 animationScene1.GoToNext()
133 animationScene1.GoToNext()
134 animationScene1.GoToNext()
135 animationScene1.GoToPrevious()
136 animationScene1.GoToPrevious()
138 # current camera placement for renderView1
139 renderView1.InteractionMode = '2D'
140 renderView1.CameraPosition = [5.0, 0.5, 10000.0]
141 renderView1.CameraFocalPoint = [5.0, 0.5, 0.0]
142 renderView1.CameraParallelScale = 5.024937810560445
146 renderView1.ViewSize =[300,300]
150 #### saving camera placements for all active views
152 # current camera placement for renderView1
153 renderView1.InteractionMode = '2D'
154 renderView1.CameraPosition = [5.0, 0.5, 10000.0]
155 renderView1.CameraFocalPoint = [5.0, 0.5, 0.0]
156 renderView1.CameraParallelScale = 5.024937810560445
157 # compare with baseline image
161 baselineIndex = sys.argv.index('-B')+1
162 baselinePath = sys.argv[baselineIndex]
164 print "Could not get baseline directory. Test failed."
166 baseline_file = os.path.join(baselinePath,png)
167 import vtk.test.Testing
168 vtk.test.Testing.VTK_TEMP_DIR = vtk.util.misc.vtkGetTempDir()
169 vtk.test.Testing.compareImage(GetActiveView().GetRenderWindow(), baseline_file, threshold=25)
170 vtk.test.Testing.interact()