ParaView
vtkSMTransferFunctionProxy.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: ParaView
4  Module: $RCSfile$
5 
6  Copyright (c) Kitware, Inc.
7  All rights reserved.
8  See Copyright.txt or http://www.paraview.org/HTML/Copyright.html for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
24 #ifndef vtkSMTransferFunctionProxy_h
25 #define vtkSMTransferFunctionProxy_h
26 
27 #include "vtkPVServerManagerRenderingModule.h" // needed for export macro
28 #include "vtkSMProxy.h"
29 
30 namespace Json
31 {
32 class Value;
33 }
34 
36 class VTKPVSERVERMANAGERRENDERING_EXPORT vtkSMTransferFunctionProxy : public vtkSMProxy
37 {
38 public:
39  static vtkSMTransferFunctionProxy* New();
41  void PrintSelf(ostream& os, vtkIndent indent);
42 
49  virtual bool RescaleTransferFunction(const double range[2], bool extend = false)
50  {
51  return this->RescaleTransferFunction(range[0], range[1], extend);
52  }
53  virtual bool RescaleTransferFunction(double rangeMin, double rangeMax, bool extend = false);
54 
56 
60  static bool RescaleTransferFunction(
61  vtkSMProxy* proxy, double rangeMin, double rangeMax, bool extend = false);
62  static bool RescaleTransferFunction(vtkSMProxy* proxy, const double range[2], bool extend = false)
63  {
64  return vtkSMTransferFunctionProxy::RescaleTransferFunction(proxy, range[0], range[1], extend);
65  }
67 
69 
74  virtual bool RescaleTransferFunctionToDataRange(bool extend = false);
75  static bool RescaleTransferFunctionToDataRange(vtkSMProxy* proxy, bool extend = false)
76  {
78  return self ? self->RescaleTransferFunctionToDataRange(extend) : false;
79  }
81 
85  virtual bool InvertTransferFunction();
86 
91  static bool InvertTransferFunction(vtkSMProxy*);
92 
100  virtual bool MapControlPointsToLogSpace(bool inverse = false);
101  virtual bool MapControlPointsToLinearSpace() { return this->MapControlPointsToLogSpace(true); }
102 
104 
108  static bool MapControlPointsToLogSpace(vtkSMProxy* proxy, bool inverse = false)
109  {
111  return self ? self->MapControlPointsToLogSpace(inverse) : false;
112  }
114 
120  {
122  }
123 
125 
133  virtual bool ApplyPreset(const Json::Value& value, bool rescale = true);
134  static bool ApplyPreset(vtkSMProxy* proxy, const Json::Value& value, bool rescale = true)
135  {
137  return self ? self->ApplyPreset(value, rescale) : false;
138  }
140 
141  virtual bool ApplyPreset(const char* presetname, bool rescale = true);
142  static bool ApplyPreset(vtkSMProxy* proxy, const char* presetname, bool rescale = true)
143  {
145  return self ? self->ApplyPreset(presetname, rescale) : false;
146  }
147 
149 
153  virtual Json::Value GetStateAsPreset();
154  static Json::Value GetStateAsPreset(vtkSMProxy* proxy);
156 
158 
163  virtual bool ApplyColorMap(const char* text);
164  virtual bool ApplyColorMap(vtkPVXMLElement* xml);
166 
168 
172  static bool ApplyColorMap(vtkSMProxy* proxy, const char* text)
173  {
175  return self ? self->ApplyColorMap(text) : false;
176  }
178 
180 
184  static bool ApplyColorMap(vtkSMProxy* proxy, vtkPVXMLElement* xml)
185  {
187  return self ? self->ApplyColorMap(xml) : false;
188  }
190 
195  virtual bool SaveColorMap(vtkPVXMLElement* xml);
196 
198 
202  static bool SaveColorMap(vtkSMProxy* proxy, vtkPVXMLElement* xml)
203  {
205  return self ? self->SaveColorMap(xml) : false;
206  }
208 
214  virtual bool IsScalarBarVisible(vtkSMProxy* view);
215 
217 
221  static bool IsScalarBarVisible(vtkSMProxy* proxy, vtkSMProxy* view)
222  {
224  return self ? self->IsScalarBarVisible(view) : false;
225  }
227 
233  virtual vtkSMProxy* FindScalarBarRepresentation(vtkSMProxy* view);
234 
236 
241  {
243  return self ? self->FindScalarBarRepresentation(view) : NULL;
244  }
246 
248 
253  virtual bool UpdateScalarBarsComponentTitle(vtkPVArrayInformation* arrayInfo);
255  {
257  return self ? self->UpdateScalarBarsComponentTitle(arrayInfo) : false;
258  }
260 
262 
267  virtual bool ComputeDataRange(double range[2]);
268  static bool ComputeDataRange(vtkSMProxy* proxy, double range[2])
269  {
271  return self ? self->ComputeDataRange(range) : false;
272  }
274 
275  // Helper method to compute the active annotated values in visible
276  // representations that use the transfer function.
277  virtual bool ComputeAvailableAnnotations(bool extend = false);
278  static bool ComputeAvailableAnnotations(vtkSMProxy* proxy, bool extend = false)
279  {
281  return self ? self->ComputeAvailableAnnotations(extend) : false;
282  }
283 
285 
290  virtual void ResetPropertiesToDefaults(const char* arrayName, bool preserve_range);
292  vtkSMProxy* proxy, const char* arrayName, bool preserve_range = false)
293  {
295  if (self)
296  {
297  self->ResetPropertiesToDefaults(arrayName, preserve_range);
298  }
299  }
300  using Superclass::ResetPropertiesToXMLDefaults;
302 
304 
309  static Json::Value ConvertLegacyColorMapXMLToJSON(vtkPVXMLElement* xml);
310  static Json::Value ConvertLegacyColorMapXMLToJSON(const char* xmlcontents);
312 
314 
318  static Json::Value ConvertMultipleLegacyColorMapXMLToJSON(vtkPVXMLElement* xml);
319  static Json::Value ConvertMultipleLegacyColorMapXMLToJSON(const char* xmlcontents);
321 
325  static bool ConvertLegacyColorMapsToJSON(const char* inxmlfile, const char* outjsonfile);
326 
328 
332  virtual bool GetRange(double range[2]);
333  static bool GetRange(vtkSMProxy* proxy, double range[2])
334  {
336  return self ? self->GetRange(range) : false;
337  }
339 
340 protected:
343 
344 private:
345  vtkSMTransferFunctionProxy(const vtkSMTransferFunctionProxy&) VTK_DELETE_FUNCTION;
346  void operator=(const vtkSMTransferFunctionProxy&) VTK_DELETE_FUNCTION;
347 };
348 
349 #endif
static vtkSMProxy * FindScalarBarRepresentation(vtkSMProxy *proxy, vtkSMProxy *view)
Safely call FindScalarBarRepresentation(..) after casting the proxy to the appropriate type...
virtual bool MapControlPointsToLogSpace(bool inverse=false)
Remaps control points by normalizing in linear-space and then interpolating in log-space.
virtual bool RescaleTransferFunction(const double range[2], bool extend=false)
Rescale the "RGBPoints" for the transfer function to match the new range.
static bool MapControlPointsToLogSpace(vtkSMProxy *proxy, bool inverse=false)
Safely call MapControlPointsToLogSpace() after casting the proxy to the appropriate type...
static bool MapControlPointsToLinearSpace(vtkSMProxy *proxy)
Safely call MapControlPointsToLinearSpace() after casting the proxy to the appropriate type...
static bool ApplyPreset(vtkSMProxy *proxy, const Json::Value &value, bool rescale=true)
Apply a preset.
vtkSMTransferFunctionProxy is the proxy used for "PVLookupTable", "ColorTransferFunction" and "Piecew...
static bool ApplyPreset(vtkSMProxy *proxy, const char *presetname, bool rescale=true)
static bool RescaleTransferFunctionToDataRange(vtkSMProxy *proxy, bool extend=false)
Locates all representations that are currently using this transfer function and then rescales the tra...
static void ResetPropertiesToDefaults(vtkSMProxy *proxy, const char *arrayName, bool preserve_range=false)
Helper method to reset a transfer function proxy to its defaults.
static bool SaveColorMap(vtkSMProxy *proxy, vtkPVXMLElement *xml)
Safely call ApplyColorMap(..) after casting the proxy to the appropriate type.
static bool ComputeDataRange(vtkSMProxy *proxy, double range[2])
Helper method used by RescaleTransferFunctionToDataRange() to compute range from all visible represen...
static bool ApplyColorMap(vtkSMProxy *proxy, vtkPVXMLElement *xml)
Safely call ApplyColorMap(..) after casting the proxy to the appropriate type.
static bool ComputeAvailableAnnotations(vtkSMProxy *proxy, bool extend=false)
proxy for a VTK object(s) on a server
Definition: vtkSMProxy.h:152
static bool IsScalarBarVisible(vtkSMProxy *proxy, vtkSMProxy *view)
Safely call IsScalarBarVisible(..) after casting the proxy to the appropriate type.
Data array information like type.
static bool GetRange(vtkSMProxy *proxy, double range[2])
Returns current transfer function data range.
static bool RescaleTransferFunction(vtkSMProxy *proxy, const double range[2], bool extend=false)
Safely call RescaleTransferFunction() after casting the proxy to appropriate type.
This is used by vtkPVXMLParser to represent an XML document starting at the root element.
static bool ApplyColorMap(vtkSMProxy *proxy, const char *text)
Safely call ApplyColorMap(..) after casting the proxy to the appropriate type.
static bool UpdateScalarBarsComponentTitle(vtkSMProxy *proxy, vtkPVArrayInformation *arrayInfo)
Update component titles for all scalar bars connected to this transfer function proxy.
static vtkSMTransferFunctionProxy * SafeDownCast(vtkObject *o)