File: vtkImplicitArray.h

package info (click to toggle)
vtk9 9.5.2%2Bdfsg3-4
  • links: PTS, VCS
  • area: main
  • in suites: sid
  • size: 205,916 kB
  • sloc: cpp: 2,336,565; ansic: 327,116; python: 111,200; yacc: 4,104; java: 3,977; sh: 3,032; xml: 2,771; perl: 2,189; lex: 1,787; makefile: 178; javascript: 165; objc: 153; tcl: 59
file content (653 lines) | stat: -rw-r--r-- 26,771 bytes parent folder | download | duplicates (5)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
// Funded by CEA, DAM, DIF, F-91297 Arpajon, France
#ifndef vtkImplicitArray_h
#define vtkImplicitArray_h

#include "vtkCommonCoreModule.h" // for export macro
#include "vtkGenericDataArray.h"
#include "vtkImplicitArrayTraits.h" // for traits

#include <memory>
#include <type_traits>

/**
 * \class vtkImplicitArray
 * \brief A read only array class that wraps an implicit function from integers to any value type
 * supported by VTK
 *
 * This templated array class allows one to mimic the vtkDataArray interface using an implicit map
 * behind the scenes. The `BackendT` type can be a class or a struct that implements a const map
 * method that takes integers to any VTK value type. It can also be any type of Closure/Functor that
 * implements a const operator() method from integers to the value type of the array. If a void
 * mapTuple(vtkIdType, TupleType*) const method is also present, the array will use this method to
 * to populate the tuple instead of the map method. If a
 * ValueType mapComponent(vtkIdType, int) const method is also present, the array will use this
 * method to populate the GetTypedComponent function instead of the map method.
 *
 * The ordering of the array for tuples and components is implicitly AOS.
 *
 * The backend type can be trivially constructible, in which case the array gets initialized with a
 * default constructed instance of the BackendT, or not default constructible, in which case the
 * backend is initially a nullptr and must be set using the `SetBackend` method.
 *
 * Being a "read_only" array, any attempt to set a value in the array will result in a warning
 * message with no change to the backend itself. This may evolve in future versions of the class as
 * the needs of users become clearer.
 *
 * The `GetVoidPointer` method will create an internal vtkAOSDataArrayTemplate and populate it with
 * the values from the implicit array and can thus be very memory intensive. The `Squeeze` method
 * will destroy this internal memory array. Both deep and shallow copies to other types of arrays
 * will populate the other array with the implicit values contained in the implicit one. Deep and
 * shallow copies from other array into this one do not make sense and will result in undefined
 * behavior. Deep and shallow copies from implicit arrays of the same type will act the same way and
 * will transfer a shared backend object pointer. Deep and shallow copies from one type of implicit
 * array to a different type should result in a compile time error.
 *
 * Constraints on the backend type are enforced through `implicit_array_traits` found in the
 * vtkImplicitArrayTraits header file. These traits use metaprogramming to check the proposed
 * backend type at compile time. The decision to use this type of structure was taken for the
 * following reasons:
 * - static dispatch of the calls to the backend (when combined with the CRTP structure of
 * vtkGenericDataArray)
 * - flexibility regarding the complexity/simplicity/nature/type... of the backend one would like to
 * use
 *
 * Example for array that always returns 42:
 * @code
 * struct Const42
 * {
 *   int operator()(vtkIdType idx) const { return 42; }
 * };
 * vtkNew<vtkImplicitArray<Const42>> arr42;
 * @endcode
 *
 * Example for array that implements map and mapTuple
 * @code
 * struct ConstTupleStruct
 * {
 * int Tuple[3] = { 0, 0, 0 };
 * // constructor
 * ConstTupleStruct(int tuple[3])
 * {
 *  this->Tuple[0] = tuple[0];
 *  this->Tuple[1] = tuple[1];
 *  this->Tuple[2] = tuple[2];
 * }
 *
 * // used for GetValue
 * int map(int idx) const
 * {
 *   int tuple[3];
 *   this->mapTuple(idx / 3, tuple);
 *   return tuple[idx % 3];
 * }
 * // used for GetTypedTuple
 * void mapTuple(int vtkNotUsed(idx), int* tuple) const
 * {
 *   tuple[0] = this->tuple[0];
 *   tuple[1] = this->tuple[1];
 *   tuple[2] = this->tuple[2];
 * }
 * };
 * @endcode
 *
 * example for array that implements map and mapComponent
 * @code
 * struct ConstComponentStruct
 * {
 *   int Tuple[3] = { 0, 0, 0 };
 *
 * ConstComponentStruct(int tuple[3])
 * {
 *   this->Tuple[0] = tuple[0];
 *   this->Tuple[1] = tuple[1];
 *   this->Tuple[2] = tuple[2];
 * }
 *
 * // used for GetValue
 * int map(int idx) const { return this->mapComponent(idx / 3, idx % 3); }
 * // used for GetTypedComponent
 * int mapComponent(int vtkNotUsed(idx), int comp) const { return this->Tuple[comp]; }
 * };
 * @endcode
 *
 * A peculiarity of `vtkImplicitArray`s is that their `NewInstance` method no longer gives
 * an instance of the exact same array type. A `NewInstance` call on a `vtkImplicitArray`
 * will return a `vtkAOSDataArrayTemplate<ValueTypeT>` with the same value type as the
 * original implicit array. This is so that the following workflow (used extensively
 * throughout VTK) can work without issues:
 * @code
 * struct Const42
 * {
 *   int operator()(vtkIdType idx) const { return 42; }
 * };
 * vtkNew<vtkImplicitArray<Const42>> arr42;
 * arr42->SetNumberOfTuples(11);
 * vtkSmartPointer<vtkDataArray> arr43 = vtk::TakeSmartPointer(arr42->NewInstance());
 * arr43->SetNumberOfComponents(arr42->GetNumberOfComponents());
 * arr43->SetNumberOfTuples(arr42->GetNumberOfTuples());
 * arr43->Fill(43);
 * @endcode
 *
 * Optionally, `vtkImplicitArray`s backends can return their memory usage in KiB by defining
 * the function `getMemorySize` returning `unsigned long`. `vtkImplicitArray` then exposes this
 * function through the `GetActualMemorySize` function. If the backend does not define it,
 * `GetActualMemorySize` always returns 1.
 *
 * @sa
 * vtkGenericDataArray vtkImplicitArrayTraits vtkDataArray
 */

//-------------------------------------------------------------------------------------------------
// Special macro for implicit array types modifying the behavior of NewInstance to provide writable
// AOS arrays instead of empty implicit arrays
#define vtkImplicitArrayTypeMacro(thisClass, superclass)                                           \
  vtkAbstractTypeMacroWithNewInstanceType(thisClass, superclass,                                   \
    vtkAOSDataArrayTemplate<typename thisClass::ValueType>, typeid(thisClass).name());             \
                                                                                                   \
protected:                                                                                         \
  vtkObjectBase* NewInstanceInternal() const override                                              \
  {                                                                                                \
    return vtkAOSDataArrayTemplate<typename thisClass::ValueType>::New();                          \
  }                                                                                                \
                                                                                                   \
public:
//-------------------------------------------------------------------------------------------------

VTK_ABI_NAMESPACE_BEGIN
template <class BackendT>
class vtkImplicitArray
  : public vtkGenericDataArray<vtkImplicitArray<BackendT>,
      typename vtk::detail::implicit_array_traits<BackendT>::rtype>
{
  using trait = vtk::detail::implicit_array_traits<BackendT>;
  static_assert(trait::can_read,
    "Supplied backend type does not have mandatory read trait. Must implement either map() const "
    "or operator() const.");
  using ValueTypeT = typename trait::rtype;
  using GenericDataArrayType = vtkGenericDataArray<vtkImplicitArray<BackendT>, ValueTypeT>;

public:
  using SelfType = vtkImplicitArray<BackendT>;
  vtkImplicitArrayTypeMacro(SelfType, GenericDataArrayType);
  using ValueType = typename GenericDataArrayType::ValueType;
  using BackendType = BackendT;

  static vtkImplicitArray* New();

  ///@{
  /**
   * Implementation of vtkGDAConceptMethods
   */
  /**
   * Get the value at @a idx. @a idx assumes AOS ordering.
   */
  ValueType GetValue(vtkIdType idx) const { return this->GetValueImpl<BackendT>(idx); }

  /**
   * Will not do anything for these read only arrays!
   */
  void SetValue(vtkIdType idx, ValueType value);

  /**
   * Copy the tuple at @a idx into @a tuple.
   */
  void GetTypedTuple(vtkIdType idx, ValueType* tuple) const
  {
    this->GetTypedTupleImpl<BackendT>(idx, tuple);
  }

  /**
   * Will not do anything for these read only arrays!
   */
  void SetTypedTuple(vtkIdType tupleIdx, const ValueType* tuple);

  /**
   * Get component @a comp of the tuple at @a idx.
   */
  ValueType GetTypedComponent(vtkIdType idx, int comp) const
  {
    return this->GetTypedComponentImpl<BackendT>(idx, comp);
  }

  /**
   * Will not do anything for these read only arrays!
   */
  void SetTypedComponent(vtkIdType tupleIdx, int comp, ValueType value);
  ///@}

  ///@{
  /**
   * Setter/Getter for Backend
   */
  void SetBackend(std::shared_ptr<BackendT> newBackend)
  {
    this->Backend = newBackend;
    this->Modified();
  }
  std::shared_ptr<BackendT> GetBackend() { return this->Backend; }
  ///@}

  /**
   * Utility method for setting backend parameterization directly
   */
  template <typename... Params>
  void ConstructBackend(Params&&... params)
  {
    this->SetBackend(std::make_shared<BackendT>(std::forward<Params>(params)...));
  }

  /**
   * Use of this method is discouraged, it creates a memory copy of the data into
   * a contiguous AoS-ordered buffer internally.
   *
   * Implicit array aims to limit memory consumption. Calling this method breaks
   * this paradigm and can cause unexpected memory consumption,
   * specially when called indirectly by some implementation details.
   * E.g. when using the numpy wrapping, see #19304.
   */
  void* GetVoidPointer(vtkIdType valueIdx) override;

  /**
   * Release all extraneous internal memory including the void pointer used by `GetVoidPointer`
   */
  void Squeeze() override;

  /**
   * Get the type of array this is when down casting
   */
  int GetArrayType() const override { return vtkAbstractArray::ImplicitArray; }

  /**
   * Reset the array to default construction
   */
  void Initialize() override
  {
    this->Initialize<BackendT>();
    this->Squeeze();
  }

  /**
   * Return the memory in kibibytes (1024 bytes) consumed by this implicit data array.
   *
   * The value returned is guaranteed to be greater than or equal to the memory required to
   * actually represent the data represented by this object.
   *
   * Implicit array backends can implement the `getMemorySize` function to override the default
   * implementation, which always returns 1.
   */
  unsigned long GetActualMemorySize() const override
  {
    return this->GetActualMemorySizeImpl<BackendT>();
  }

  /**
   * Specific DeepCopy for implicit arrays
   *
   * This method should be preferred for two implicit arrays having the same backend. We cannot call
   * the method `DeepCopy` since that conflicts with the virtual function of the same name that
   * cannot be templated. The non-interopability of templates and virtual functions is a language
   * limitation at the time of writing this documentation.  We can call this from the dispatched
   * version of the `DeepCopy` in `vtkDataArray`. However, the implicit array needs to be
   * dispatchable in order to to not enter into the Generic implementation of the deep copy. This
   * dispatch is not always the case for all implicit arrays.
   */
  template <typename OtherBackend>
  void ImplicitDeepCopy(vtkImplicitArray<OtherBackend>* other)
  {
    static_assert(std::is_same<BackendT, OtherBackend>::value,
      "Cannot copy implicit array with one type of backend to an implicit array with a different "
      "type of backend");
    this->SetNumberOfComponents(other->GetNumberOfComponents());
    this->SetNumberOfTuples(other->GetNumberOfTuples());
    this->SetBackend(other->GetBackend());
  }

  ///@{
  /**
   * Perform a fast, safe cast from a vtkAbstractArray to a vtkDataArray.
   * This method checks if source->GetArrayType() returns DataArray
   * or a more derived type, and performs a static_cast to return
   * source as a vtkDataArray pointer. Otherwise, nullptr is returned.
   */
  static vtkImplicitArray<BackendT>* FastDownCast(vtkAbstractArray* source);
  ///@}

protected:
  vtkImplicitArray();
  ~vtkImplicitArray() override;

  ///@{
  /**
   * No allocation necessary
   */
  bool AllocateTuples(vtkIdType vtkNotUsed(numTuples)) { return true; }
  bool ReallocateTuples(vtkIdType vtkNotUsed(numTuples)) { return true; }
  ///@}

  struct vtkInternals;
  std::unique_ptr<vtkInternals> Internals;

  /**
   * The backend object actually mapping the indexes
   */
  std::shared_ptr<BackendT> Backend;

private:
  vtkImplicitArray(const vtkImplicitArray&) = delete;
  void operator=(const vtkImplicitArray&) = delete;

  ///@{
  /**
   * Methods for static dispatch towards map trait
   */
  template <typename U>
  typename std::enable_if<vtk::detail::has_map_trait<U>::value, ValueType>::type GetValueImpl(
    vtkIdType idx) const
  {
    return this->Backend->map(idx);
  }
  ///@}

  ///@{
  /**
   * Methods for static dispatch towards closure trait
   */
  template <typename U>
  typename std::enable_if<vtk::detail::is_closure_trait<U>::value, ValueType>::type GetValueImpl(
    vtkIdType idx) const
  {
    return (*this->Backend)(idx);
  }
  ///@}

  ///@{
  /**
   * Static dispatch Initialize for default constructible things
   */
  template <typename U>
  typename std::enable_if<vtk::detail::implicit_array_traits<U>::default_constructible, void>::type
  Initialize()
  {
    this->Backend = std::make_shared<BackendT>();
  }
  ///@}

  ///@{
  /**
   * Static dispatch Initialize for non-default constructible things
   */
  template <typename U>
  typename std::enable_if<!vtk::detail::implicit_array_traits<U>::default_constructible, void>::type
  Initialize()
  {
    this->Backend = nullptr;
  }
  ///@}

  ///@{
  /**
   * Static dispatch tuple mapping for compatible backends
   */
  template <typename U>
  typename std::enable_if<vtk::detail::implicit_array_traits<U>::can_direct_read_tuple, void>::type
  GetTypedTupleImpl(vtkIdType idx, ValueType* tuple) const
  {
    static_assert(
      std::is_same<typename vtk::detail::can_map_tuple_trait<U>::rtype, ValueType>::value,
      "Tuple type should be the same as the return type of the mapTuple");
    this->Backend->mapTuple(idx, tuple);
  }
  ///@}

  ///@{
  /**
   * Static dispatch tuple mapping using component mapping for compatible backends
   */
  template <typename U>
  typename std::enable_if<!vtk::detail::implicit_array_traits<U>::can_direct_read_tuple &&
      vtk::detail::implicit_array_traits<U>::can_direct_read_component,
    void>::type
  GetTypedTupleImpl(vtkIdType idx, ValueType* tuple) const
  {
    for (vtkIdType comp = 0; comp < this->NumberOfComponents; comp++)
    {
      tuple[comp] = this->GetTypedComponent(idx, comp);
    }
  }

  /**
   * Static dispatch tuple mapping for incompatible backends
   */
  template <typename U>
  typename std::enable_if<!vtk::detail::implicit_array_traits<U>::can_direct_read_tuple &&
      !vtk::detail::implicit_array_traits<U>::can_direct_read_component,
    void>::type
  GetTypedTupleImpl(vtkIdType idx, ValueType* tuple) const
  {
    const vtkIdType tupIdx = idx * this->NumberOfComponents;
    for (vtkIdType comp = 0; comp < this->NumberOfComponents; comp++)
    {
      tuple[comp] = this->GetValue(tupIdx + comp);
    }
  }
  ///@}

  ///@{
  /**
   * Static dispatch component mapping for compatible backends
   */
  template <typename U>
  typename std::enable_if<vtk::detail::implicit_array_traits<U>::can_direct_read_component,
    ValueType>::type
  GetTypedComponentImpl(vtkIdType idx, int comp) const
  {
    static_assert(
      std::is_same<typename vtk::detail::can_map_component_trait<U>::rtype, ValueType>::value,
      "Component return type should be the same as the return type of the mapComponent");
    return this->Backend->mapComponent(idx, comp);
  }
  ///@}

  ///@{
  /**
   * Static dispatch component mapping for incompatible backends
   */
  template <typename U>
  typename std::enable_if<!vtk::detail::implicit_array_traits<U>::can_direct_read_component,
    ValueType>::type
  GetTypedComponentImpl(vtkIdType idx, int comp) const
  {
    return this->GetValue(idx * this->NumberOfComponents + comp);
  }
  ///@}

  ///@{
  /**
   * Static call to get memory size for compatible backends
   */
  template <typename U>
  typename std::enable_if<vtk::detail::implicit_array_traits<U>::can_get_memory_size,
    unsigned long>::type
  GetActualMemorySizeImpl() const
  {
    return this->Backend->getMemorySize();
  }

  /**
   * Static call to get memory size for incompatible backends
   * For those backends, dhe default memory size is 1KiB.
   */
  template <typename U>
  typename std::enable_if<!vtk::detail::implicit_array_traits<U>::can_get_memory_size,
    unsigned long>::type
  GetActualMemorySizeImpl() const
  {
    return 1;
  }
  ///@}

  friend class vtkGenericDataArray<vtkImplicitArray<BackendT>, ValueTypeT>;
};

// Declare vtkArrayDownCast implementations for implicit containers:
vtkArrayDownCast_TemplateFastCastMacro(vtkImplicitArray);
VTK_ABI_NAMESPACE_END

#include "vtkImplicitArray.txx"

#endif // vtkImplicitArray_h

// See vtkGenericDataArray for similar section
#ifdef VTK_IMPLICIT_VALUERANGE_INSTANTIATING
VTK_ABI_NAMESPACE_BEGIN
template <typename ValueType>
struct vtkAffineImplicitBackend;
template <typename ValueType>
class vtkCompositeImplicitBackend;
template <typename ValueType>
struct vtkConstantImplicitBackend;
template <typename ValueType>
class vtkStructuredPointBackend;
template <typename ValueType>
class vtkIndexedImplicitBackend;
VTK_ABI_NAMESPACE_END
#include <functional>

// Needed to export for this module and not CommonCore
#define VTK_INSTANTIATE_VALUERANGE_ARRAYTYPE(ArrayType, ValueType)                                 \
  template VTKCOMMONCORE_EXPORT bool DoComputeScalarRange(                                         \
    ArrayType*, ValueType*, vtkDataArrayPrivate::AllValues, const unsigned char*, unsigned char);  \
  template VTKCOMMONCORE_EXPORT bool DoComputeScalarRange(ArrayType*, ValueType*,                  \
    vtkDataArrayPrivate::FiniteValues, const unsigned char*, unsigned char);                       \
  template VTKCOMMONCORE_EXPORT bool DoComputeVectorRange(ArrayType*, ValueType[2],                \
    vtkDataArrayPrivate::AllValues, const unsigned char*, unsigned char);                          \
  template VTKCOMMONCORE_EXPORT bool DoComputeVectorRange(ArrayType*, ValueType[2],                \
    vtkDataArrayPrivate::FiniteValues, const unsigned char*, unsigned char);

#define VTK_INSTANTIATE_VALUERANGE_VALUETYPE(ValueType)                                            \
  VTK_INSTANTIATE_VALUERANGE_ARRAYTYPE(                                                            \
    vtkImplicitArray<vtkAffineImplicitBackend<ValueType>>, ValueType)                              \
  VTK_INSTANTIATE_VALUERANGE_ARRAYTYPE(                                                            \
    vtkImplicitArray<vtkCompositeImplicitBackend<ValueType>>, ValueType)                           \
  VTK_INSTANTIATE_VALUERANGE_ARRAYTYPE(                                                            \
    vtkImplicitArray<vtkConstantImplicitBackend<ValueType>>, ValueType)                            \
  VTK_INSTANTIATE_VALUERANGE_ARRAYTYPE(                                                            \
    vtkImplicitArray<vtkStructuredPointBackend<ValueType>>, ValueType)                             \
  VTK_INSTANTIATE_VALUERANGE_ARRAYTYPE(                                                            \
    vtkImplicitArray<vtkIndexedImplicitBackend<ValueType>>, ValueType)                             \
  VTK_INSTANTIATE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<ValueType(int)>>, ValueType)

#elif defined(VTK_USE_EXTERN_TEMPLATE) // VTK_IMPLICIT_VALUERANGE_INSTANTIATING

#ifndef VTK_IMPLICIT_TEMPLATE_EXTERN
#define VTK_IMPLICIT_TEMPLATE_EXTERN
#ifdef _MSC_VER
#pragma warning(push)
// The following is needed when the following is declared
// dllexport and is used from another class in vtkCommonCore
#pragma warning(disable : 4910) // extern and dllexport incompatible
#endif

VTK_ABI_NAMESPACE_BEGIN
template <typename ValueType>
struct vtkAffineImplicitBackend;
template <typename ValueType>
class vtkCompositeImplicitBackend;
template <typename ValueType>
struct vtkConstantImplicitBackend;
template <typename ValueType>
class vtkStructuredPointBackend;
template <typename ValueType>
class vtkIndexedImplicitBackend;
VTK_ABI_NAMESPACE_END
#include <functional>

namespace vtkDataArrayPrivate
{
VTK_ABI_NAMESPACE_BEGIN
template <typename A, typename R, typename T>
VTKCOMMONCORE_EXPORT bool DoComputeScalarRange(
  A*, R*, T, const unsigned char* ghosts, unsigned char ghostsToSkip);
template <typename A, typename R>
VTKCOMMONCORE_EXPORT bool DoComputeVectorRange(
  A*, R[2], AllValues, const unsigned char* ghosts, unsigned char ghostsToSkip);
template <typename A, typename R>
VTKCOMMONCORE_EXPORT bool DoComputeVectorRange(
  A*, R[2], FiniteValues, const unsigned char* ghosts, unsigned char ghostsToSkip);
VTK_ABI_NAMESPACE_END
} // namespace vtkDataArrayPrivate

#define VTK_DECLARE_VALUERANGE_ARRAYTYPE(ArrayType, ValueType)                                     \
  extern template VTKCOMMONCORE_EXPORT bool DoComputeScalarRange(                                  \
    ArrayType*, ValueType*, vtkDataArrayPrivate::AllValues, const unsigned char*, unsigned char);  \
  extern template VTKCOMMONCORE_EXPORT bool DoComputeScalarRange(ArrayType*, ValueType*,           \
    vtkDataArrayPrivate::FiniteValues, const unsigned char*, unsigned char);                       \
  extern template VTKCOMMONCORE_EXPORT bool DoComputeVectorRange(ArrayType*, ValueType[2],         \
    vtkDataArrayPrivate::AllValues, const unsigned char*, unsigned char);                          \
  extern template VTKCOMMONCORE_EXPORT bool DoComputeVectorRange(ArrayType*, ValueType[2],         \
    vtkDataArrayPrivate::FiniteValues, const unsigned char*, unsigned char);

#define VTK_DECLARE_VALUERANGE_VALUETYPE(ValueType)                                                \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(                                                                \
    vtkImplicitArray<vtkAffineImplicitBackend<ValueType>>, ValueType)                              \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(                                                                \
    vtkImplicitArray<vtkCompositeImplicitBackend<ValueType>>, ValueType)                           \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(                                                                \
    vtkImplicitArray<vtkConstantImplicitBackend<ValueType>>, ValueType)                            \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(                                                                \
    vtkImplicitArray<vtkStructuredPointBackend<ValueType>>, ValueType)                             \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(                                                                \
    vtkImplicitArray<vtkIndexedImplicitBackend<ValueType>>, ValueType)                             \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<ValueType(int)>>, ValueType)

#define VTK_DECLARE_VALUERANGE_IMPLICIT_BACKENDTYPE(BackendT)                                      \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<float>>, double)                      \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<double>>, double)                     \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<char>>, double)                       \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<signed char>>, double)                \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<unsigned char>>, double)              \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<short>>, double)                      \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<unsigned short>>, double)             \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<int>>, double)                        \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<unsigned int>>, double)               \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<long>>, double)                       \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<unsigned long>>, double)              \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<long long>>, double)                  \
  VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<BackendT<unsigned long long>>, double)

namespace vtkDataArrayPrivate
{
VTK_ABI_NAMESPACE_BEGIN
VTK_DECLARE_VALUERANGE_IMPLICIT_BACKENDTYPE(vtkAffineImplicitBackend)
VTK_DECLARE_VALUERANGE_IMPLICIT_BACKENDTYPE(vtkConstantImplicitBackend)
VTK_DECLARE_VALUERANGE_IMPLICIT_BACKENDTYPE(vtkCompositeImplicitBackend)
VTK_DECLARE_VALUERANGE_IMPLICIT_BACKENDTYPE(vtkStructuredPointBackend)
VTK_DECLARE_VALUERANGE_IMPLICIT_BACKENDTYPE(vtkIndexedImplicitBackend)

VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<float(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<double(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<char(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<signed char>(int)>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<unsigned char(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<short(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<unsigned short(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<int(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<unsigned int(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<long(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<unsigned long(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<long long(int)>>, double)
VTK_DECLARE_VALUERANGE_ARRAYTYPE(vtkImplicitArray<std::function<unsigned long long(int)>>, double)
VTK_ABI_NAMESPACE_END
}

#undef VTK_DECLARE_VALUERANGE_ARRAYTYPE
#undef VTK_DECLARE_VALUERANGE_VALUETYPE

#ifdef _MSC_VER
#pragma warning(pop)
#endif
#endif // VTK_IMPLICIT_TEMPLATE_EXTERN

#endif // VTK_IMPLICIT_VALUERANGE_INSTANTIATING