TRACE User Guide  TRACE Version 9.6.1
Co-processing with ParaView Catalyst

How to build the Catalyst library for use with TRACE is described in Installing TRACE. Co-processing is activated by the control file command ActivateCatalyst.

To be able to use the filters distributed with TRACE, make sure to set the environment variable

PV_PLUGIN_PATH=/path/to/trace_suite/distrib/paraviewPlugins/
Note
All the examples below have been tested with Catalyst 5.11 and Python 3.9. In Catalyst versions before 5.10, the scripting language was different. We recommend using the latest version of Catalyst, but at least 5.8.

Writing VTK output

A generic script can be used to write the data passed to Catalyst in VTK format: allinputsgridwriter.py

from paraview.simple import *
from paraview import coprocessing
# a root directory under which all Catalyst output goes
rootDirectory = '../catalyst'
# the frequency to output everything
outputfrequency = 1
def getInputs(datadescription):
for i in range(datadescription.GetNumberOfInputDescriptions()):
yield datadescription.GetInputDescription(i), datadescription.GetInputDescriptionName(i)
# ----------------------- CoProcessor definition -----------------------
def CreateCoProcessor():
def _CreatePipeline(coprocessor, datadescription):
class Pipeline:
for inputdescription, name in getInputs(datadescription):
adaptorinput = coprocessor.CreateProducer( datadescription, name )
if adaptorinput is None:
continue
grid = adaptorinput.GetClientSideObject().GetOutputDataObject(0)
if grid.IsA('vtkImageData') or grid.IsA('vtkUniformGrid'):
writer = coprocessor.CreateWriter( XMLPImageDataWriter, f"{name}/{name}"+"_%t.pvti", outputfrequency )
elif grid.IsA('vtkRectilinearGrid'):
writer = coprocessor.CreateWriter( XMLPRectilinearGridWriter, f"{name}/{name}"+"_%t.pvtr", outputfrequency )
elif grid.IsA('vtkStructuredGrid'):
writer = coprocessor.CreateWriter( XMLPStructuredGridWriter, f"{name}/{name}"+"_%t.pvts", outputfrequency )
elif grid.IsA('vtkPolyData'):
writer = coprocessor.CreateWriter( XMLPPolyDataWriter, f"{name}/{name}"+"_%t.pvtp", outputfrequency )
elif grid.IsA('vtkUnstructuredGrid'):
writer = coprocessor.CreateWriter( XMLPUnstructuredGridWriter, f"{name}/{name}"+"_%t.pvtu", outputfrequency )
elif grid.IsA('vtkUniformGridAMR'):
writer = coprocessor.CreateWriter( XMLHierarchicalBoxDataWriter, f"{name}/{name}"+"_%t.vthb", outputfrequency )
elif grid.IsA('vtkMultiBlockDataSet'):
writer = coprocessor.CreateWriter( XMLMultiBlockDataWriter, f"{name}/{name}"+"_%t.vtm", outputfrequency )
else:
print("Don't know how to create a writer for a ", grid.GetClassName())
return Pipeline()
class CoProcessor(coprocessing.CoProcessor):
def CreatePipeline(self, datadescription):
self.Pipeline = _CreatePipeline(self, datadescription)
coprocessor = CoProcessor()
if rootDirectory:
coprocessor.SetRootDirectory(rootDirectory)
return coprocessor
#--------------------------------------------------------------
# Global variables that will hold the pipeline for each timestep
# Creating the CoProcessor object, doesn't actually create the ParaView pipeline.
# It will be automatically setup when coprocessor.UpdateProducers() is called the
# first time.
coprocessor = CreateCoProcessor()
#--------------------------------------------------------------
# Enable Live-Visualizaton with ParaView
coprocessor.EnableLiveVisualization(False)
# ---------------------- Data Selection method ----------------------
def RequestDataDescription(datadescription):
"Callback to populate the request for current timestep"
# We are just going to request all fields and meshes from the simulation
# code/adaptor.
for inputdescription, _ in getInputs(datadescription):
inputdescription.AllFieldsOn()
inputdescription.GenerateMeshOn()
# ------------------------ Processing method ------------------------
def DoCoProcessing(datadescription):
"Callback to do co-processing for current timestep"
global coprocessor
# Update the coprocessor by providing it the newly generated simulation data.
# If the pipeline hasn't been setup yet, this will setup the pipeline.
coprocessor.UpdateProducers(datadescription)
# Write output data, if appropriate.
coprocessor.WriteData(datadescription);
# Write image capture (Last arg: rescale lookup table), if appropriate.
coprocessor.WriteImages(datadescription, rescale_lookuptable=False)
# Live Visualization, if enabled.
coprocessor.DoLiveVisualization(datadescription, "localhost", 22222)

The writer will generate one set of files per time step for each dataset passed to Catalyst. Each block group can be found as separate '<block group name>_*.vtm' and each panel group can be found as separate '<panel group name>_*.vtm'. You do not have to use all of the generated output for your Catalyst script.

Warning
When using allinputsgridwriter.py in a production run, be aware that this will trigger the pipeline in every time step irrespective of the specified output frequency. This overhead can be avoided with the sub-option --outputInterval for the command ActivateCatalyst.

Configuring the view with ParaView

  • Load the output generated in the step above ParaView to configure the view.
  • Apply the filters you need for your visualization.
  • Define exports using extractors for data or images.
    • Set up a writer to write the VTK data set if required.
    • Set up screen shots.
  • Export the state for the co-processor and set up rendering of images and/or live visualization using File->Save Catalyst state....
Note
In version 5.3.0 of ParaView the script generator results in Python files that contain commands not available in Catalyst. This seems to be a bug which will result in Python errors during a TRACE run. It aborts the co-processing but TRACE continues to run. Comment out the offending lines until you have a working script. This problem no longer occurs in version 5.8.0.
Example
Pipeline creating an image and a data extract using the DuplicatePeriodic filter. Each block group and panel group is a separate input channel.
# script-version: 2.0
# Catalyst state generated using paraview version 5.11.0
from paraview.simple import *
paraview.simple._DisableFirstRenderCameraReset()
# ----------------------------------------------------------------
# setup views used in the visualization
# ----------------------------------------------------------------
# get the material library
materialLibrary1 = GetMaterialLibrary()
# Create a new 'Render View'
renderView1 = CreateView('RenderView')
renderView1.ViewSize = [1424, 1184]
renderView1.AxesGrid = 'GridAxes3DActor'
renderView1.CenterOfRotation = [0.02664399892091751, 0.010776616632938385, 0.000697574985679239]
renderView1.StereoType = 'Crystal Eyes'
renderView1.CameraPosition = [0.02664399892091751, 0.010776616632938385, 0.5325977589873594]
renderView1.CameraFocalPoint = [0.02664399892091751, 0.010776616632938385, 0.000697574985679239]
renderView1.CameraFocalDisk = 1.0
renderView1.CameraParallelScale = 0.13766589771316995
renderView1.BackEnd = 'OSPRay raycaster'
renderView1.OSPRayMaterialLibrary = materialLibrary1
SetActiveView(None)
# ----------------------------------------------------------------
# setup view layouts
# ----------------------------------------------------------------
# create new layout object 'Layout #1'
layout1 = CreateLayout(name='Layout #1')
layout1.AssignView(0, renderView1)
layout1.SetSize(1424, 1184)
# ----------------------------------------------------------------
# restore active view
SetActiveView(renderView1)
# ----------------------------------------------------------------
# ----------------------------------------------------------------
# setup the data processing pipelines
# ----------------------------------------------------------------
# create a new 'XML MultiBlock Data Reader'
blade = XMLMultiBlockDataReader(registrationName='blade', FileName=['/tmp/translational_URANS/catalyst/blade/blade_10001.vtm', '/tmp/translational_URANS/catalyst/blade/blade_10002.vtm', '/tmp/translational_URANS/catalyst/blade/blade_10003.vtm', '/tmp/translational_URANS/catalyst/blade/blade_10004.vtm'])
blade.PointArrayStatus = ['Density', 'Velocity', 'Pressure']
# create a new 'XML MultiBlock Data Reader'
bladeWall = XMLMultiBlockDataReader(registrationName='bladeWall', FileName=['/tmp/translational_URANS/catalyst/bladeWall/bladeWall_10001.vtm', '/tmp/translational_URANS/catalyst/bladeWall/bladeWall_10002.vtm', '/tmp/translational_URANS/catalyst/bladeWall/bladeWall_10003.vtm', '/tmp/translational_URANS/catalyst/bladeWall/bladeWall_10004.vtm'])
bladeWall.CellArrayStatus = ['Pressure', 'Density', 'ShearStressWall']
# create a new 'XML MultiBlock Data Reader'
bar = XMLMultiBlockDataReader(registrationName='bar', FileName=['/tmp/translational_URANS/catalyst/bar/bar_10001.vtm', '/tmp/translational_URANS/catalyst/bar/bar_10002.vtm', '/tmp/translational_URANS/catalyst/bar/bar_10003.vtm', '/tmp/translational_URANS/catalyst/bar/bar_10004.vtm'])
bar.PointArrayStatus = ['Density', 'Velocity', 'Pressure']
# create a new 'Duplicate Periodic'
duplicatePeriodic1 = DuplicatePeriodic(registrationName='DuplicatePeriodic1', Input=bar)
duplicatePeriodic1.Pitchfactor = -1.0
# ----------------------------------------------------------------
# setup the visualization in view 'renderView1'
# ----------------------------------------------------------------
# show data from bar
barDisplay = Show(bar, renderView1, 'StructuredGridRepresentation')
# get 2D transfer function for 'Velocity'
velocityTF2D = GetTransferFunction2D('Velocity')
# get color transfer function/color map for 'Velocity'
velocityLUT = GetColorTransferFunction('Velocity')
velocityLUT.TransferFunction2D = velocityTF2D
velocityLUT.RGBPoints = [6.332992634021132e-09, 0.231373, 0.298039, 0.752941, 159.62285390953883, 0.865003, 0.865003, 0.865003, 319.2457078127447, 0.705882, 0.0156863, 0.14902]
velocityLUT.ScalarRangeInitialized = 1.0
# get opacity transfer function/opacity map for 'Velocity'
velocityPWF = GetOpacityTransferFunction('Velocity')
velocityPWF.Points = [6.332992634021132e-09, 0.0, 0.5, 0.0, 319.2457078127447, 1.0, 0.5, 0.0]
velocityPWF.ScalarRangeInitialized = 1
# trace defaults for the display properties.
barDisplay.Representation = 'Surface'
barDisplay.ColorArrayName = ['POINTS', 'Velocity']
barDisplay.LookupTable = velocityLUT
barDisplay.SelectTCoordArray = 'None'
barDisplay.SelectNormalArray = 'None'
barDisplay.SelectTangentArray = 'None'
barDisplay.OSPRayScaleArray = 'Density'
barDisplay.OSPRayScaleFunction = 'PiecewiseFunction'
barDisplay.SelectOrientationVectors = 'None'
barDisplay.ScaleFactor = 0.008835949446074664
barDisplay.SelectScaleArray = 'Density'
barDisplay.GlyphType = 'Arrow'
barDisplay.GlyphTableIndexArray = 'Density'
barDisplay.GaussianRadius = 0.00044179747230373324
barDisplay.SetScaleArray = ['POINTS', 'Density']
barDisplay.ScaleTransferFunction = 'PiecewiseFunction'
barDisplay.OpacityArray = ['POINTS', 'Density']
barDisplay.OpacityTransferFunction = 'PiecewiseFunction'
barDisplay.DataAxesGrid = 'GridAxesRepresentation'
barDisplay.PolarAxes = 'PolarAxesRepresentation'
barDisplay.ScalarOpacityFunction = velocityPWF
barDisplay.ScalarOpacityUnitDistance = 0.00490397951940153
barDisplay.SelectInputVectors = ['POINTS', 'Velocity']
barDisplay.WriteLog = ''
# init the 'PiecewiseFunction' selected for 'ScaleTransferFunction'
barDisplay.ScaleTransferFunction.Points = [0.06381958723068237, 0.0, 0.5, 0.0, 0.17131200432777405, 1.0, 0.5, 0.0]
# init the 'PiecewiseFunction' selected for 'OpacityTransferFunction'
barDisplay.OpacityTransferFunction.Points = [0.06381958723068237, 0.0, 0.5, 0.0, 0.17131200432777405, 1.0, 0.5, 0.0]
# show data from blade
bladeDisplay = Show(blade, renderView1, 'StructuredGridRepresentation')
# trace defaults for the display properties.
bladeDisplay.Representation = 'Surface'
bladeDisplay.ColorArrayName = ['POINTS', 'Velocity']
bladeDisplay.LookupTable = velocityLUT
bladeDisplay.SelectTCoordArray = 'None'
bladeDisplay.SelectNormalArray = 'None'
bladeDisplay.SelectTangentArray = 'None'
bladeDisplay.OSPRayScaleArray = 'Density'
bladeDisplay.OSPRayScaleFunction = 'PiecewiseFunction'
bladeDisplay.SelectOrientationVectors = 'None'
bladeDisplay.ScaleFactor = 0.020064663887023926
bladeDisplay.SelectScaleArray = 'Density'
bladeDisplay.GlyphType = 'Arrow'
bladeDisplay.GlyphTableIndexArray = 'Density'
bladeDisplay.GaussianRadius = 0.0010032331943511962
bladeDisplay.SetScaleArray = ['POINTS', 'Density']
bladeDisplay.ScaleTransferFunction = 'PiecewiseFunction'
bladeDisplay.OpacityArray = ['POINTS', 'Density']
bladeDisplay.OpacityTransferFunction = 'PiecewiseFunction'
bladeDisplay.DataAxesGrid = 'GridAxesRepresentation'
bladeDisplay.PolarAxes = 'PolarAxesRepresentation'
bladeDisplay.ScalarOpacityFunction = velocityPWF
bladeDisplay.ScalarOpacityUnitDistance = 0.00910452617438506
bladeDisplay.SelectInputVectors = ['POINTS', 'Velocity']
bladeDisplay.WriteLog = ''
# init the 'PiecewiseFunction' selected for 'ScaleTransferFunction'
bladeDisplay.ScaleTransferFunction.Points = [0.08122637122869492, 0.0, 0.5, 0.0, 0.14759118854999542, 1.0, 0.5, 0.0]
# init the 'PiecewiseFunction' selected for 'OpacityTransferFunction'
bladeDisplay.OpacityTransferFunction.Points = [0.08122637122869492, 0.0, 0.5, 0.0, 0.14759118854999542, 1.0, 0.5, 0.0]
# show data from duplicatePeriodic1
duplicatePeriodic1Display = Show(duplicatePeriodic1, renderView1, 'StructuredGridRepresentation')
# trace defaults for the display properties.
duplicatePeriodic1Display.Representation = 'Surface'
duplicatePeriodic1Display.ColorArrayName = ['POINTS', 'Velocity']
duplicatePeriodic1Display.LookupTable = velocityLUT
duplicatePeriodic1Display.SelectTCoordArray = 'None'
duplicatePeriodic1Display.SelectNormalArray = 'None'
duplicatePeriodic1Display.SelectTangentArray = 'None'
duplicatePeriodic1Display.OSPRayScaleArray = 'Density'
duplicatePeriodic1Display.OSPRayScaleFunction = 'PiecewiseFunction'
duplicatePeriodic1Display.SelectOrientationVectors = 'None'
duplicatePeriodic1Display.ScaleFactor = 0.008835949655622245
duplicatePeriodic1Display.SelectScaleArray = 'Density'
duplicatePeriodic1Display.GlyphType = 'Arrow'
duplicatePeriodic1Display.GlyphTableIndexArray = 'Density'
duplicatePeriodic1Display.GaussianRadius = 0.0004417974827811122
duplicatePeriodic1Display.SetScaleArray = ['POINTS', 'Density']
duplicatePeriodic1Display.ScaleTransferFunction = 'PiecewiseFunction'
duplicatePeriodic1Display.OpacityArray = ['POINTS', 'Density']
duplicatePeriodic1Display.OpacityTransferFunction = 'PiecewiseFunction'
duplicatePeriodic1Display.DataAxesGrid = 'GridAxesRepresentation'
duplicatePeriodic1Display.PolarAxes = 'PolarAxesRepresentation'
duplicatePeriodic1Display.ScalarOpacityFunction = velocityPWF
duplicatePeriodic1Display.ScalarOpacityUnitDistance = 0.0049039796043389445
duplicatePeriodic1Display.SelectInputVectors = ['POINTS', 'Velocity']
duplicatePeriodic1Display.WriteLog = ''
# init the 'PiecewiseFunction' selected for 'ScaleTransferFunction'
duplicatePeriodic1Display.ScaleTransferFunction.Points = [0.06381958723068237, 0.0, 0.5, 0.0, 0.17131200432777405, 1.0, 0.5, 0.0]
# init the 'PiecewiseFunction' selected for 'OpacityTransferFunction'
duplicatePeriodic1Display.OpacityTransferFunction.Points = [0.06381958723068237, 0.0, 0.5, 0.0, 0.17131200432777405, 1.0, 0.5, 0.0]
# setup the color legend parameters for each legend in this view
# get color legend/bar for velocityLUT in view renderView1
velocityLUTColorBar = GetScalarBar(velocityLUT, renderView1)
velocityLUTColorBar.WindowLocation = 'Upper Right Corner'
velocityLUTColorBar.Title = 'Velocity'
velocityLUTColorBar.ComponentTitle = 'Magnitude'
# set color bar visibility
velocityLUTColorBar.Visibility = 1
# show color legend
barDisplay.SetScalarBarVisibility(renderView1, True)
# show color legend
bladeDisplay.SetScalarBarVisibility(renderView1, True)
# show color legend
duplicatePeriodic1Display.SetScalarBarVisibility(renderView1, True)
# ----------------------------------------------------------------
# setup color maps and opacity mapes used in the visualization
# note: the Get..() functions create a new object, if needed
# ----------------------------------------------------------------
# ----------------------------------------------------------------
# setup extractors
# ----------------------------------------------------------------
# create extractor
pNG1 = CreateExtractor('PNG', renderView1, registrationName='PNG1')
# trace defaults for the extractor.
pNG1.Trigger = 'TimeStep'
# init the 'TimeStep' selected for 'Trigger'
pNG1.Trigger.Frequency = 20
# init the 'PNG' selected for 'Writer'
pNG1.Writer.FileName = 'RenderView1_{timestep:06d}{camera}.png'
pNG1.Writer.ImageResolution = [1424, 1184]
pNG1.Writer.Format = 'PNG'
# create extractor
vTM1 = CreateExtractor('VTM', bladeWall, registrationName='VTM1')
# trace defaults for the extractor.
vTM1.Trigger = 'TimeStep'
# init the 'TimeStep' selected for 'Trigger'
vTM1.Trigger.Frequency = 100
# init the 'VTM' selected for 'Writer'
vTM1.Writer.FileName = 'bladeWall_{timestep:06d}.vtm'
vTM1.Writer.CompressorType = 'ZLib'
# ----------------------------------------------------------------
# restore active source
SetActiveSource(bladeWall)
# ----------------------------------------------------------------
# ------------------------------------------------------------------------------
# Catalyst options
from paraview import catalyst
options = catalyst.Options()
options.ExtractsOutputDirectory = '../catalyst'
options.GlobalTrigger = 'TimeStep'
options.CatalystLiveTrigger = 'TimeStep'
# init the 'TimeStep' selected for 'GlobalTrigger'
options.GlobalTrigger.Frequency = 1
# ------------------------------------------------------------------------------
if __name__ == '__main__':
from paraview.simple import SaveExtractsUsingCatalystOptions
# Code for non in-situ environments; if executing in post-processing
# i.e. non-Catalyst mode, let's generate extracts using Catalyst options
SaveExtractsUsingCatalystOptions(options)

Live visualization

It is possible to connect ParaView directly to a running TRACE computation. For this to work, the live visualization must be enabled in the co-processing script and a host and port to which data should be sent must be specified. A minimal example which does only live visualization and sends data to the localhost on port 22222 every 10 time steps is given in live.py:

#--------------------------------------------------------------
# Global timestep output options
timeStepToStartOutputAt=0
forceOutputAtFirstCall=False
# Global screenshot output options
imageFileNamePadding=0
rescale_lookuptable=False
# Whether or not to request specific arrays from the adaptor.
requestSpecificArrays=False
# a root directory under which all Catalyst output goes
rootDirectory=''
# makes a cinema D index table
make_cinema_table=False
#--------------------------------------------------------------
# Code generated from cpstate.py to create the CoProcessor.
# paraview version 5.8.0
#--------------------------------------------------------------
from paraview.simple import *
from paraview import coprocessing
# ----------------------- CoProcessor definition -----------------------
def CreateCoProcessor():
def _CreatePipeline(coprocessor, datadescription):
class Pipeline:
input = coprocessor.CreateProducer(datadescription, 'input')
return Pipeline()
class CoProcessor(coprocessing.CoProcessor):
def CreatePipeline(self, datadescription):
self.Pipeline = _CreatePipeline(self, datadescription)
coprocessor = CoProcessor()
# these are the frequencies at which the coprocessor updates.
freqs = {'input': [10]}
coprocessor.SetUpdateFrequencies(freqs)
if requestSpecificArrays:
arrays = [['Density', 0], ['DensityGradientMagnitude', 0], ['Lambda2', 0], ['Pressure', 0], ['VelocityX', 0], ['VelocityY', 0], ['VelocityZ', 0], ['VorticityStreamwise', 0]]
coprocessor.SetRequestedArrays('input', arrays)
coprocessor.SetInitialOutputOptions(timeStepToStartOutputAt,forceOutputAtFirstCall)
if rootDirectory:
coprocessor.SetRootDirectory(rootDirectory)
if make_cinema_table:
coprocessor.EnableCinemaDTable()
return coprocessor
#--------------------------------------------------------------
# Global variable that will hold the pipeline for each timestep
# Creating the CoProcessor object, doesn't actually create the ParaView pipeline.
# It will be automatically setup when coprocessor.UpdateProducers() is called the
# first time.
coprocessor = CreateCoProcessor()
#--------------------------------------------------------------
# Enable Live-Visualizaton with ParaView and the update frequency
coprocessor.EnableLiveVisualization(True, 10)
# ---------------------- Data Selection method ----------------------
def RequestDataDescription(datadescription):
"Callback to populate the request for current timestep"
global coprocessor
# setup requests for all inputs based on the requirements of the
# pipeline.
coprocessor.LoadRequestedData(datadescription)
# ------------------------ Processing method ------------------------
def DoCoProcessing(datadescription):
"Callback to do co-processing for current timestep"
global coprocessor
# Update the coprocessor by providing it the newly generated simulation data.
# If the pipeline hasn't been setup yet, this will setup the pipeline.
coprocessor.UpdateProducers(datadescription)
# Write output data, if appropriate.
coprocessor.WriteData(datadescription);
# Write image capture (Last arg: rescale lookup table), if appropriate.
coprocessor.WriteImages(datadescription, rescale_lookuptable=rescale_lookuptable,
image_quality=0, padding_amount=imageFileNamePadding)
# Live Visualization, if enabled.
coprocessor.DoLiveVisualization(datadescription, "localhost", 22222)

Adding static data

If you want to have static data in your view, you can manually replace the source in the catalyst script by an appropriate reader. Make sure to remove the source from the coprocessor settings below. This can be used to display solid walls.

@@ -107,7 +107,8 @@ def CreateCoProcessor():
# create a new 'XML Partitioned Unstructured Grid Reader'
# create a producer from a simulation input
- geometrypvtu = coprocessor.CreateProducer(datadescription, 'geometry.pvtu')
+ geometrypvtu = XMLPartitionedUnstructuredGridReader(FileName=['geometry.pvtu'])
+ geometrypvtu.CellArrayStatus = []
# ----------------------------------------------------------------
# setup the visualization in view 'renderView1'
@@ -196,13 +197,11 @@ def CreateCoProcessor():
coprocessor = CoProcessor()
# these are the frequencies at which the coprocessor updates.
- freqs = {'input': [1], 'geometry.pvtu': [1]}
+ freqs = {'input': [1]}
coprocessor.SetUpdateFrequencies(freqs)
if requestSpecificArrays:
arrays = [['Density', 0], ['DensityGradientMagnitude', 0], ['Lambda2', 0], ['Pressure', 0], ['Velocity', 0], ['VorticityStreamwise', 0]]
coprocessor.SetRequestedArrays('input', arrays)
- arrays = []
- coprocessor.SetRequestedArrays('geometry.pvtu', arrays)
coprocessor.SetInitialOutputOptions(timeStepToStartOutputAt,forceOutputAtFirstCall)
if rootDirectory: