BioDynaMo  v1.05.120-25dc9790
mapped_data_array.h
Go to the documentation of this file.
1 // -----------------------------------------------------------------------------
2 //
3 // Copyright (C) 2021 CERN & University of Surrey for the benefit of the
4 // BioDynaMo collaboration. All Rights Reserved.
5 //
6 // Licensed under the Apache License, Version 2.0 (the "License");
7 // you may not use this file except in compliance with the License.
8 //
9 // See the LICENSE file distributed with this work for details.
10 // See the NOTICE file distributed with this work for additional information
11 // regarding copyright ownership.
12 //
13 // -----------------------------------------------------------------------------
14 
15 #ifndef CORE_VISUALIZATION_PARAVIEW_MAPPED_DATA_ARRAY_H_
16 #define CORE_VISUALIZATION_PARAVIEW_MAPPED_DATA_ARRAY_H_
17 
18 #include <algorithm>
19 #include <string>
20 
21 #include <vtkMappedDataArray.h>
22 #include <vtkObjectFactory.h>
23 #include <thread>
24 #include <vector>
25 
26 #include "core/agent/agent.h"
28 #include "core/functor.h"
29 #include "core/param/param.h"
30 #include "core/util/thread_info.h"
31 #include "core/util/type.h"
32 
33 namespace bdm {
34 
35 // -----------------------------------------------------------------------------
38 template <typename TReturn, typename TClass, typename TDataMember>
40  uint64_t dm_offset_;
41  using TempValueType = typename std::remove_pointer<TReturn>::type;
43  mutable std::vector<std::array<TempValueType, 64 / sizeof(TempValueType)>>
45 
46  GetDataMemberForVis() { temp_values_.resize(256); }
47 
48  void Update() {
49  // We assume that number of threads won't real_t within one iteration
50  temp_values_.resize(
51  2 * std::max<unsigned long long>(
52  ThreadInfo::GetInstance()->GetMaxUniversalThreadId(), 256ULL));
53  }
54 
56 
57  template <typename T>
58  static constexpr DataType GetDataType() {
59  if (IsArray<T>::value) {
60  return DataType::kArray;
61  } else if (std::is_same<T, AgentUid>::value) {
62  return DataType::kAgentUid;
63  } else if (is_agent_ptr<T>::value) {
64  return DataType::kSoPointer;
65  }
66  return DataType::kDefault;
67  }
68 
69  template <typename TTDataMember = TDataMember>
70  typename std::enable_if<GetDataType<TTDataMember>() == DataType::kDefault,
71  TReturn>::type
72  operator()(Agent* agent) const {
73  auto* casted_agent = static_cast<TClass*>(agent);
74  return reinterpret_cast<TDataMember*>(
75  reinterpret_cast<char*>(casted_agent) + dm_offset_);
76  }
77 
78  template <typename TTDataMember = TDataMember>
79  typename std::enable_if<GetDataType<TTDataMember>() == DataType::kArray,
80  TReturn>::type
81  operator()(Agent* agent) const {
82  auto* casted_agent = static_cast<TClass*>(agent);
83  auto* data = reinterpret_cast<TDataMember*>(
84  reinterpret_cast<char*>(casted_agent) + dm_offset_)
85  ->data();
86  return const_cast<TReturn>(data);
87  }
88 
89  template <typename TTDataMember = TDataMember>
90  typename std::enable_if<GetDataType<TTDataMember>() == DataType::kAgentUid,
91  TReturn>::type
92  operator()(Agent* agent) const {
93  auto* casted_agent = static_cast<TClass*>(agent);
94  auto* data = reinterpret_cast<TDataMember*>(
95  reinterpret_cast<char*>(casted_agent) + dm_offset_);
96  uint64_t uid = *data;
98  assert(temp_values_.size() > tid);
99  temp_values_[tid][0] = uid;
100  return &temp_values_[tid][0];
101  }
102 
103  template <typename TTDataMember = TDataMember>
104  typename std::enable_if<GetDataType<TTDataMember>() == DataType::kSoPointer,
105  TReturn>::type
106  operator()(Agent* agent) const {
107  auto* casted_agent = static_cast<TClass*>(agent);
108  auto* data = reinterpret_cast<TDataMember*>(
109  reinterpret_cast<char*>(casted_agent) + dm_offset_);
110  uint64_t uid = data->GetUidAsUint64();
112  assert(temp_values_.size() > tid);
113  temp_values_[tid][0] = uid;
114  return &temp_values_[tid][0];
115  }
116 };
117 
118 // -----------------------------------------------------------------------------
120  MappedDataArrayInterface() = default;
121  virtual ~MappedDataArrayInterface() = default;
122  virtual void Update(const std::vector<Agent*>* agents, uint64_t start,
123  uint64_t end) = 0;
124 };
125 
126 // -----------------------------------------------------------------------------
127 template <typename TScalar, typename TClass, typename TDataMember>
128 class MappedDataArray : public vtkMappedDataArray<TScalar>,
129  public MappedDataArrayInterface {
130  public:
131  vtkAbstractTemplateTypeMacro(MappedDataArray, vtkMappedDataArray<TScalar>)
132  vtkMappedDataArrayNewInstanceMacro(
133  MappedDataArray) static MappedDataArray* New();
134  typedef typename Superclass::ValueType ValueType;
135 
136  void Initialize(Param::MappedDataArrayMode mode, const std::string& name,
137  uint64_t num_components, uint64_t dm_offset);
138  void Update(const std::vector<Agent*>* agents, uint64_t start,
139  uint64_t end) final;
140 
141  // Reimplemented virtuals -- see superclasses for descriptions:
142  void PrintSelf(ostream& os, vtkIndent indent) final;
143  void Initialize() final;
144  void GetTuples(vtkIdList* pt_ids, vtkAbstractArray* output) final;
145  void GetTuples(vtkIdType p1, vtkIdType p2, vtkAbstractArray* output) final;
146  void Squeeze() final;
147  vtkArrayIterator* NewIterator() final;
148  vtkIdType LookupValue(vtkVariant value) final;
149  void LookupValue(vtkVariant value, vtkIdList* ids) final;
150  vtkVariant GetVariantValue(vtkIdType idx) final;
151  void ClearLookup() final;
152  double* GetTuple(vtkIdType i) final;
153  void GetTuple(vtkIdType i, double* tuple) final;
154  vtkIdType LookupTypedValue(TScalar value) final;
155  void LookupTypedValue(TScalar value, vtkIdList* ids) final;
156  ValueType GetValue(vtkIdType idx) const final;
157  TScalar& GetValueReference(vtkIdType idx) final;
158  void GetTypedTuple(vtkIdType idx, TScalar* t) const final;
159 
160  // Description:
161  // This container is read only -- this method does nothing but print a
162  // warning.
163  int Allocate(vtkIdType sz, vtkIdType ext) final;
164  int Resize(vtkIdType num_tuples) final;
165  void SetNumberOfTuples(vtkIdType number) final;
166  void SetTuple(vtkIdType i, vtkIdType j, vtkAbstractArray* source) final;
167  void SetTuple(vtkIdType i, const float* source) final;
168  void SetTuple(vtkIdType i, const double* source) final;
169  void InsertTuple(vtkIdType i, vtkIdType j, vtkAbstractArray* source) final;
170  void InsertTuple(vtkIdType i, const float* source) final;
171  void InsertTuple(vtkIdType i, const double* source) final;
172  void InsertTuples(vtkIdList* dst_ids, vtkIdList* src_ids,
173  vtkAbstractArray* source) final;
174  void InsertTuples(vtkIdType dst_start, vtkIdType n, vtkIdType src_start,
175  vtkAbstractArray* source) final;
176  vtkIdType InsertNextTuple(vtkIdType j, vtkAbstractArray* source) final;
177  vtkIdType InsertNextTuple(const float* source) final;
178  vtkIdType InsertNextTuple(const double* source) final;
179  void DeepCopy(vtkAbstractArray* aa) final;
180  void DeepCopy(vtkDataArray* da) final;
181  void InterpolateTuple(vtkIdType i, vtkIdList* pt_indices,
182  vtkAbstractArray* source, double* weights) final;
183  void InterpolateTuple(vtkIdType i, vtkIdType id1, vtkAbstractArray* source1,
184  vtkIdType id2, vtkAbstractArray* source2,
185  double t) final;
186  void SetVariantValue(vtkIdType idx, vtkVariant value) final;
187  void InsertVariantValue(vtkIdType idx, vtkVariant value) final;
188  void RemoveTuple(vtkIdType id) final;
189  void RemoveFirstTuple() final;
190  void RemoveLastTuple() final;
191  void SetTypedTuple(vtkIdType i, const TScalar* t) final;
192  void InsertTypedTuple(vtkIdType i, const TScalar* t) final;
193  vtkIdType InsertNextTypedTuple(const TScalar* t) final;
194  void SetValue(vtkIdType idx, TScalar value) final;
195  vtkIdType InsertNextValue(TScalar v) final;
196  void InsertValue(vtkIdType idx, TScalar v) final;
197 
198  protected:
199  MappedDataArray();
200  ~MappedDataArray() override;
201 
203  GetDataMemberForVis<TScalar*, TClass, TDataMember> get_dm_;
204  const std::vector<Agent*>* agents_ = nullptr;
205  uint64_t start_ = 0;
206  uint64_t end_ = 0;
207  double* temp_array_ = nullptr;
208 
209  Param::MappedDataArrayMode mode_;
212  mutable uint64_t match_value_ = 0;
216  mutable std::vector<uint64_t> is_matching_;
218  mutable std::vector<TScalar> data_;
219 
220  private:
221  MappedDataArray(const MappedDataArray&) = delete;
222  void operator=(const MappedDataArray&) = delete;
223 
224  vtkIdType Lookup(const TScalar& val, vtkIdType start_index);
225 };
226 
227 // ----------------------------------------------------------------------------
228 // Implementation
229 template <typename TScalar, typename TClass, typename TDataMember>
230 void MappedDataArray<TScalar, TClass, TDataMember>::Initialize(
231  Param::MappedDataArrayMode mode, const std::string& name,
232  uint64_t num_components, uint64_t dm_offset) {
233  get_dm_.dm_offset_ = dm_offset;
234  mode_ = mode;
235  this->NumberOfComponents = num_components;
236  this->SetName(name.c_str());
237 }
238 
239 template <typename TScalar, typename TClass, typename TDataMember>
241  const std::vector<Agent*>* agents, uint64_t start, uint64_t end) {
242  agents_ = agents;
243  start_ = start;
244  end_ = end;
245  get_dm_.Update();
246  this->Size = this->NumberOfComponents * (end - start);
247  this->MaxId = this->Size - 1;
248 
249  if (this->Size <= 0) {
250  this->Size = 0;
251  this->MaxId = -1;
252  this->Modified();
253  return;
254  }
255 
256  this->Modified();
257 
258  if (mode_ != Param::MappedDataArrayMode::kZeroCopy) {
259  if (data_.capacity() < this->Size) {
260  data_.reserve(this->Size * 1.5);
261  }
262  if (mode_ == Param::MappedDataArrayMode::kCopy) {
263  uint64_t counter = 0;
264  for (uint64_t i = start; i < end; ++i) {
265  auto* data = get_dm_((*agents_)[i]);
266  for (uint64_t c = 0; c < this->NumberOfComponents; ++c) {
267  data_[counter++] = data[c];
268  }
269  }
270  } else {
271  // mode_ == Param::MappedDataArrayMode::kCache
272  auto cap = is_matching_.capacity();
273  if (cap < this->Size) {
274  is_matching_.reserve(this->Size * 1.5);
275  for (uint64_t i = cap; i < is_matching_.capacity(); ++i) {
277  }
278  }
279  match_value_++;
280  }
281  }
282 }
283 
284 //------------------------------------------------------------------------------
285 // Can't use vtkStandardNewMacro with a template.
286 template <typename TScalar, typename TClass, typename TDataMember>
289  VTK_STANDARD_NEW_BODY(MappedDataArray);
290 }
291 
292 //------------------------------------------------------------------------------
293 template <typename TScalar, typename TClass, typename TDataMember>
295  ostream& os, vtkIndent indent) {
297  os, indent);
298  os << indent << "agents_: " << this->agents_ << std::endl;
299  os << indent << "start_: " << this->start_ << std::endl;
300  os << indent << "end_: " << this->end_ << std::endl;
301  os << indent << "temp_array_: " << this->temp_array_ << std::endl;
302  os << indent << "mode_: " << this->mode_ << std::endl;
303  os << indent << "match_value_: " << this->match_value_ << std::endl;
304  os << indent << "is_matching_.size(): " << this->is_matching_.size()
305  << std::endl;
306  os << indent << "data_.size(): " << this->data_.size() << std::endl;
307 }
308 
309 //------------------------------------------------------------------------------
310 template <typename TScalar, typename TClass, typename TDataMember>
312  this->MaxId = -1;
313  this->Size = 0;
314  this->NumberOfComponents = 1;
315 }
316 
317 //------------------------------------------------------------------------------
318 template <typename TScalar, typename TClass, typename TDataMember>
320  vtkIdList* pt_ids, vtkAbstractArray* output) {
321  vtkDataArray* out_array = vtkDataArray::FastDownCast(output);
322  if (!out_array) {
323  vtkWarningMacro(<< "Input is not a vtkDataArray");
324  return;
325  }
326 
327  vtkIdType num_tuples = pt_ids->GetNumberOfIds();
328 
329  out_array->SetNumberOfComponents(this->NumberOfComponents);
330  out_array->SetNumberOfTuples(num_tuples);
331 
332  for (vtkIdType i = 0; i < pt_ids->GetNumberOfIds(); ++i) {
333  out_array->SetTuple(i, this->GetTuple(pt_ids->GetId(i)));
334  }
335 }
336 
337 //------------------------------------------------------------------------------
338 template <typename TScalar, typename TClass, typename TDataMember>
340  vtkIdType p1, vtkIdType p2, vtkAbstractArray* output) {
341  vtkDataArray* da = vtkDataArray::FastDownCast(output);
342  if (!da) {
343  vtkErrorMacro(<< "Input is not a vtkDataArray");
344  return;
345  }
346 
347  if (da->GetNumberOfComponents() != this->GetNumberOfComponents()) {
348  vtkErrorMacro(<< "Incorrect number of components in input array.");
349  return;
350  }
351 
352  for (vtkIdType id = 0; p1 <= p2; ++p1) {
353  da->SetTuple(id++, this->GetTuple(p1));
354  }
355 }
356 
357 //------------------------------------------------------------------------------
358 template <typename TScalar, typename TClass, typename TDataMember>
360  // noop
361 }
362 
363 //------------------------------------------------------------------------------
364 template <typename TScalar, typename TClass, typename TDataMember>
366  vtkErrorMacro(<< "Not implemented.");
367  return NULL;
368 }
369 
370 //------------------------------------------------------------------------------
371 template <typename TScalar, typename TClass, typename TDataMember>
373  vtkVariant value) {
374  bool valid = true;
375  TScalar val = vtkVariantCast<TScalar>(value, &valid);
376  if (valid) {
377  return this->Lookup(val, 0);
378  }
379  return -1;
380 }
381 
382 //------------------------------------------------------------------------------
383 template <typename TScalar, typename TClass, typename TDataMember>
385  vtkVariant value, vtkIdList* ids) {
386  bool valid = true;
387  TScalar val = vtkVariantCast<TScalar>(value, &valid);
388  ids->Reset();
389  if (valid) {
390  vtkIdType index = 0;
391  while ((index = this->Lookup(val, index)) >= 0) {
392  ids->InsertNextId(index++);
393  }
394  }
395 }
396 
397 //------------------------------------------------------------------------------
398 template <typename TScalar, typename TClass, typename TDataMember>
400  vtkIdType idx) {
401  return vtkVariant(this->GetValueReference(idx));
402 }
403 
404 //------------------------------------------------------------------------------
405 template <typename TScalar, typename TClass, typename TDataMember>
407  // no-op, no fast lookup implemented.
408 }
409 
410 //------------------------------------------------------------------------------
411 template <typename TScalar, typename TClass, typename TDataMember>
413  this->GetTuple(i, this->temp_array_);
414  return this->temp_array_;
415 }
416 
417 //------------------------------------------------------------------------------
418 template <typename TScalar, typename TClass, typename TDataMember>
420  double* tuple) {
421  uint64_t idx = tuple_id * this->NumberOfComponents;
422  TScalar* data;
423  switch (mode_) {
424  case Param::MappedDataArrayMode::kCache: {
425  bool all_matches = true;
426  for (uint64_t i = idx;
427  i < idx + static_cast<uint64_t>(this->NumberOfComponents); ++i) {
428  all_matches |= (is_matching_[i] == match_value_);
429  }
430  if (all_matches) {
431  uint64_t cnt = 0;
432  for (uint64_t i = idx;
433  i < idx + static_cast<uint64_t>(this->NumberOfComponents); ++i) {
434  tuple[cnt++] = static_cast<real_t>(data_[i]);
435  }
436  return;
437  }
438  }
439  case Param::MappedDataArrayMode::kZeroCopy: {
440  auto* temporary_data = get_dm_((*agents_)[start_ + tuple_id]);
441  for (uint64_t i = 0; i < static_cast<uint64_t>(this->NumberOfComponents);
442  ++i) {
443  tuple[i] = static_cast<real_t>(temporary_data[i]);
444  if (mode_ == Param::MappedDataArrayMode::kCache) {
445  data_[idx + i] = temporary_data[i];
446  is_matching_[idx + i] = match_value_;
447  }
448  }
449  return;
450  } break;
451  case Param::MappedDataArrayMode::kCopy:
452  uint64_t cnt = 0;
453  for (uint64_t i = idx;
454  i < idx + static_cast<uint64_t>(this->NumberOfComponents); ++i) {
455  tuple[cnt++] = static_cast<real_t>(data_[i]);
456  }
457  }
458 }
459 
460 //------------------------------------------------------------------------------
461 template <typename TScalar, typename TClass, typename TDataMember>
463  TScalar value) {
464  return this->Lookup(value, 0);
465 }
466 
467 //------------------------------------------------------------------------------
468 template <typename TScalar, typename TClass, typename TDataMember>
470  TScalar value, vtkIdList* ids) {
471  ids->Reset();
472  vtkIdType index = 0;
473  while ((index = this->Lookup(value, index)) >= 0) {
474  ids->InsertNextId(index++);
475  }
476 }
477 
478 //------------------------------------------------------------------------------
479 template <typename TScalar, typename TClass, typename TDataMember>
482  // TODO(lukas) resolve code duplication with GetValueReference
483  // make non virtual private helper function and call it from here and
484  // GetValueReference
485  switch (mode_) {
486  case Param::MappedDataArrayMode::kCache:
487  if (is_matching_[idx] == match_value_) {
488  return data_[idx];
489  }
490  case Param::MappedDataArrayMode::kZeroCopy: {
491  if (this->NumberOfComponents == 1) {
492  auto* data = get_dm_((*agents_)[start_ + idx]);
493  if (mode_ == Param::MappedDataArrayMode::kCache) {
494  data_[idx] = *data;
495  is_matching_[idx] = match_value_;
496  }
497  return *data;
498  }
499  const vtkIdType tuple = idx / this->NumberOfComponents;
500  const vtkIdType comp = idx % this->NumberOfComponents;
501  auto* data = get_dm_((*agents_)[start_ + tuple]);
502  if (mode_ == Param::MappedDataArrayMode::kCache) {
503  data_[idx] = data[comp];
504  is_matching_[idx] = match_value_;
505  }
506  return data[comp];
507  } break;
508  case Param::MappedDataArrayMode::kCopy:
509  return data_[idx];
510  }
511 }
512 
513 //------------------------------------------------------------------------------
514 template <typename TScalar, typename TClass, typename TDataMember>
516  vtkIdType idx) {
517  switch (mode_) {
518  case Param::MappedDataArrayMode::kCache:
519  if (is_matching_[idx] == match_value_) {
520  return data_[idx];
521  }
522  case Param::MappedDataArrayMode::kZeroCopy: {
523  if (this->NumberOfComponents == 1) {
524  auto* data = get_dm_((*agents_)[start_ + idx]);
525  if (mode_ == Param::MappedDataArrayMode::kCache) {
526  data_[idx] = *data;
527  is_matching_[idx] = match_value_;
528  }
529  return *data;
530  }
531  const vtkIdType tuple = idx / this->NumberOfComponents;
532  const vtkIdType comp = idx % this->NumberOfComponents;
533  auto* data = get_dm_((*agents_)[start_ + tuple]);
534  if (mode_ == Param::MappedDataArrayMode::kCache) {
535  data_[idx] = data[comp];
536  is_matching_[idx] = match_value_;
537  }
538  return data[comp];
539  } break;
540  case Param::MappedDataArrayMode::kCopy:
541  return data_[idx];
542  }
543 }
544 
545 //------------------------------------------------------------------------------
546 template <typename TScalar, typename TClass, typename TDataMember>
548  vtkIdType tuple_id, TScalar* tuple) const {
549  auto* data = get_dm_((*agents_)[start_ + tuple_id]);
550  for (uint64_t i = 0; i < static_cast<uint64_t>(this->NumberOfComponents);
551  ++i) {
552  tuple[i] = data[i];
553  }
554 }
555 
556 //------------------------------------------------------------------------------
557 template <typename TScalar, typename TClass, typename TDataMember>
559  vtkIdType) {
560  vtkErrorMacro("Read only container.");
561  return 0;
562 }
563 
564 //------------------------------------------------------------------------------
565 template <typename TScalar, typename TClass, typename TDataMember>
567  vtkErrorMacro("Read only container.");
568  return 0;
569 }
570 
571 //------------------------------------------------------------------------------
572 template <typename TScalar, typename TClass, typename TDataMember>
574  vtkIdType) {
575  vtkErrorMacro("Read only container.");
576  return;
577 }
578 
579 //------------------------------------------------------------------------------
580 template <typename TScalar, typename TClass, typename TDataMember>
582  vtkIdType, vtkIdType, vtkAbstractArray*) {
583  vtkErrorMacro("Read only container.");
584  return;
585 }
586 
587 //------------------------------------------------------------------------------
588 template <typename TScalar, typename TClass, typename TDataMember>
590  const float*) {
591  vtkErrorMacro("Read only container.");
592  return;
593 }
594 
595 //------------------------------------------------------------------------------
596 template <typename TScalar, typename TClass, typename TDataMember>
598  const double*) {
599  vtkErrorMacro("Read only container.");
600  return;
601 }
602 
603 //------------------------------------------------------------------------------
604 template <typename TScalar, typename TClass, typename TDataMember>
606  vtkIdType, vtkIdType, vtkAbstractArray*) {
607  vtkErrorMacro("Read only container.");
608  return;
609 }
610 
611 //------------------------------------------------------------------------------
612 template <typename TScalar, typename TClass, typename TDataMember>
614  const float*) {
615  vtkErrorMacro("Read only container.");
616  return;
617 }
618 
619 //------------------------------------------------------------------------------
620 template <typename TScalar, typename TClass, typename TDataMember>
622  const double*) {
623  vtkErrorMacro("Read only container.");
624  return;
625 }
626 
627 //------------------------------------------------------------------------------
628 template <typename TScalar, typename TClass, typename TDataMember>
630  vtkIdList*, vtkIdList*, vtkAbstractArray*) {
631  vtkErrorMacro("Read only container.");
632  return;
633 }
634 
635 //------------------------------------------------------------------------------
636 template <typename TScalar, typename TClass, typename TDataMember>
638  vtkIdType, vtkIdType, vtkIdType, vtkAbstractArray*) {
639  vtkErrorMacro("Read only container.");
640  return;
641 }
642 
643 //------------------------------------------------------------------------------
644 template <typename TScalar, typename TClass, typename TDataMember>
646  vtkIdType, vtkAbstractArray*) {
647  vtkErrorMacro("Read only container.");
648  return -1;
649 }
650 
651 //------------------------------------------------------------------------------
652 template <typename TScalar, typename TClass, typename TDataMember>
654  const float*) {
655  vtkErrorMacro("Read only container.");
656  return -1;
657 }
658 
659 //------------------------------------------------------------------------------
660 template <typename TScalar, typename TClass, typename TDataMember>
662  const double*) {
663  vtkErrorMacro("Read only container.");
664  return -1;
665 }
666 
667 //------------------------------------------------------------------------------
668 template <typename TScalar, typename TClass, typename TDataMember>
670  vtkAbstractArray*) {
671  vtkErrorMacro("Read only container.");
672  return;
673 }
674 
675 //------------------------------------------------------------------------------
676 template <typename TScalar, typename TClass, typename TDataMember>
678  vtkErrorMacro("Read only container.");
679  return;
680 }
681 
682 //------------------------------------------------------------------------------
683 template <typename TScalar, typename TClass, typename TDataMember>
685  vtkIdType, vtkIdList*, vtkAbstractArray*, double*) {
686  vtkErrorMacro("Read only container.");
687  return;
688 }
689 
690 //------------------------------------------------------------------------------
691 template <typename TScalar, typename TClass, typename TDataMember>
693  vtkIdType, vtkIdType, vtkAbstractArray*, vtkIdType, vtkAbstractArray*,
694  double) {
695  vtkErrorMacro("Read only container.");
696  return;
697 }
698 
699 //------------------------------------------------------------------------------
700 template <typename TScalar, typename TClass, typename TDataMember>
702  vtkIdType, vtkVariant) {
703  vtkErrorMacro("Read only container.");
704  return;
705 }
706 
707 //------------------------------------------------------------------------------
708 template <typename TScalar, typename TClass, typename TDataMember>
710  vtkIdType, vtkVariant) {
711  vtkErrorMacro("Read only container.");
712  return;
713 }
714 
715 //------------------------------------------------------------------------------
716 template <typename TScalar, typename TClass, typename TDataMember>
718  vtkErrorMacro("Read only container.");
719  return;
720 }
721 
722 //------------------------------------------------------------------------------
723 template <typename TScalar, typename TClass, typename TDataMember>
725  vtkErrorMacro("Read only container.");
726  return;
727 }
728 
729 //------------------------------------------------------------------------------
730 template <typename TScalar, typename TClass, typename TDataMember>
732  vtkErrorMacro("Read only container.");
733  return;
734 }
735 
736 //------------------------------------------------------------------------------
737 template <typename TScalar, typename TClass, typename TDataMember>
739  vtkIdType, const TScalar*) {
740  vtkErrorMacro("Read only container.");
741  return;
742 }
743 
744 //------------------------------------------------------------------------------
745 template <typename TScalar, typename TClass, typename TDataMember>
747  vtkIdType, const TScalar*) {
748  vtkErrorMacro("Read only container.");
749  return;
750 }
751 
752 //------------------------------------------------------------------------------
753 template <typename TScalar, typename TClass, typename TDataMember>
755  const TScalar*) {
756  vtkErrorMacro("Read only container.");
757  return -1;
758 }
759 
760 //------------------------------------------------------------------------------
761 template <typename TScalar, typename TClass, typename TDataMember>
763  TScalar) {
764  vtkErrorMacro("Read only container.");
765  return;
766 }
767 
768 //------------------------------------------------------------------------------
769 template <typename TScalar, typename TClass, typename TDataMember>
771  TScalar) {
772  vtkErrorMacro("Read only container.");
773  return -1;
774 }
775 
776 //------------------------------------------------------------------------------
777 template <typename TScalar, typename TClass, typename TDataMember>
779  TScalar) {
780  vtkErrorMacro("Read only container.");
781  return;
782 }
783 
784 //------------------------------------------------------------------------------
785 template <typename TScalar, typename TClass, typename TDataMember>
787  this->temp_array_ = new double[this->NumberOfComponents];
788 }
789 
790 //------------------------------------------------------------------------------
791 template <typename TScalar, typename TClass, typename TDataMember>
793  delete[] this->temp_array_;
794 }
795 
796 //------------------------------------------------------------------------------
797 template <typename TScalar, typename TClass, typename TDataMember>
799  const TScalar& val, vtkIdType index) {
800  while (index <= this->MaxId) {
801  if (this->GetValueReference(index++) == val) {
802  return index;
803  }
804  }
805  return -1;
806 }
807 
808 } // namespace bdm
809 
810 #endif // CORE_VISUALIZATION_PARAVIEW_MAPPED_DATA_ARRAY_H_
bdm::MappedDataArrayInterface::Update
virtual void Update(const std::vector< Agent * > *agents, uint64_t start, uint64_t end)=0
bdm::MappedDataArray::RemoveFirstTuple
void RemoveFirstTuple() final
Definition: mapped_data_array.h:724
bdm::MappedDataArray::Update
void Update(const std::vector< Agent * > *agents, uint64_t start, uint64_t end) final
Definition: mapped_data_array.h:240
bdm::GetDataMemberForVis< TScalar *, TClass, TDataMember >::TempValueType
typename std::remove_pointer< TScalar * >::type TempValueType
Definition: mapped_data_array.h:41
bdm::GetDataMemberForVis
Definition: mapped_data_array.h:39
bdm::MappedDataArray::get_dm_
GetDataMemberForVis< TScalar *, TClass, TDataMember > get_dm_
Access agent data member functor.
Definition: mapped_data_array.h:203
bdm::GetDataMemberForVis::kArray
@ kArray
Definition: mapped_data_array.h:55
bdm::MappedDataArray::agents_
const std::vector< Agent * > * agents_
Definition: mapped_data_array.h:204
bdm::MappedDataArray::SetVariantValue
void SetVariantValue(vtkIdType idx, vtkVariant value) final
Definition: mapped_data_array.h:701
bdm::ThreadInfo::GetInstance
static ThreadInfo * GetInstance()
Definition: thread_info.cc:21
bdm::GetDataMemberForVis::kSoPointer
@ kSoPointer
Definition: mapped_data_array.h:55
bdm::MappedDataArray::vtkAbstractTemplateTypeMacro
vtkAbstractTemplateTypeMacro(MappedDataArray, vtkMappedDataArray< TScalar >) vtkMappedDataArrayNewInstanceMacro(MappedDataArray) static MappedDataArray *New()
bdm::MappedDataArray::InsertVariantValue
void InsertVariantValue(vtkIdType idx, vtkVariant value) final
Definition: mapped_data_array.h:709
bdm
Definition: agent.cc:39
bdm::MappedDataArray::mode_
Param::MappedDataArrayMode mode_
Definition: mapped_data_array.h:209
thread_info.h
bdm::GetDataMemberForVis::operator()
std::enable_if< GetDataType< TTDataMember >)==DataType::kAgentUid, TReturn >::type operator()(Agent *agent) const
Definition: mapped_data_array.h:92
bdm::MappedDataArray::Allocate
int Allocate(vtkIdType sz, vtkIdType ext) final
Definition: mapped_data_array.h:558
bdm::GetDataMemberForVis< TScalar *, TClass, TDataMember >::DataType
DataType
Definition: mapped_data_array.h:55
bdm::MappedDataArrayInterface
Definition: mapped_data_array.h:119
bdm::MappedDataArray::temp_array_
double * temp_array_
Definition: mapped_data_array.h:207
bdm::GetDataMemberForVis::operator()
std::enable_if< GetDataType< TTDataMember >)==DataType::kDefault, TReturn >::type operator()(Agent *agent) const
Definition: mapped_data_array.h:72
bdm::GetDataMemberForVis::Update
void Update()
Definition: mapped_data_array.h:48
bdm::MappedDataArray::InsertNextTypedTuple
vtkIdType InsertNextTypedTuple(const TScalar *t) final
Definition: mapped_data_array.h:754
bdm::MappedDataArray
Definition: mapped_data_array.h:128
bdm::MappedDataArray::PrintSelf
void PrintSelf(ostream &os, vtkIndent indent) final
Definition: mapped_data_array.h:294
bdm::real_t
double real_t
Definition: real_t.h:21
bdm::MappedDataArray::Lookup
vtkIdType Lookup(const TScalar &val, vtkIdType start_index)
Definition: mapped_data_array.h:798
agent_pointer.h
bdm::GetDataMemberForVis::GetDataType
static constexpr DataType GetDataType()
Definition: mapped_data_array.h:58
type.h
bdm::MappedDataArray::InterpolateTuple
void InterpolateTuple(vtkIdType i, vtkIdList *pt_indices, vtkAbstractArray *source, double *weights) final
Definition: mapped_data_array.h:684
bdm::GetDataMemberForVis::dm_offset_
uint64_t dm_offset_
Definition: mapped_data_array.h:40
bdm::MappedDataArray::ClearLookup
void ClearLookup() final
Definition: mapped_data_array.h:406
bdm::Agent
Contains code required by all agents.
Definition: agent.h:79
bdm::MappedDataArrayInterface::MappedDataArrayInterface
MappedDataArrayInterface()=default
bdm::MappedDataArray::GetValueReference
TScalar & GetValueReference(vtkIdType idx) final
Definition: mapped_data_array.h:515
bdm::Param::MappedDataArrayMode
MappedDataArrayMode
Definition: param.h:545
bdm::MappedDataArray::LookupValue
vtkIdType LookupValue(vtkVariant value) final
Definition: mapped_data_array.h:372
bdm::MappedDataArray::start_
uint64_t start_
Definition: mapped_data_array.h:205
bdm::MappedDataArray::NewIterator
vtkArrayIterator * NewIterator() final
Definition: mapped_data_array.h:365
bdm::MappedDataArray::match_value_
uint64_t match_value_
Definition: mapped_data_array.h:212
bdm::MappedDataArray::RemoveTuple
void RemoveTuple(vtkIdType id) final
Definition: mapped_data_array.h:717
bdm::MappedDataArray::InsertNextTuple
vtkIdType InsertNextTuple(vtkIdType j, vtkAbstractArray *source) final
Definition: mapped_data_array.h:645
bdm::MappedDataArray::LookupTypedValue
vtkIdType LookupTypedValue(TScalar value) final
Definition: mapped_data_array.h:462
bdm::MappedDataArray::DeepCopy
void DeepCopy(vtkAbstractArray *aa) final
Definition: mapped_data_array.h:669
bdm::MappedDataArray::MappedDataArray
MappedDataArray()
Definition: mapped_data_array.h:786
bdm::GetDataMemberForVis::operator()
std::enable_if< GetDataType< TTDataMember >)==DataType::kSoPointer, TReturn >::type operator()(Agent *agent) const
Definition: mapped_data_array.h:106
bdm::MappedDataArray::InsertTuple
void InsertTuple(vtkIdType i, vtkIdType j, vtkAbstractArray *source) final
Definition: mapped_data_array.h:605
bdm::MappedDataArray::~MappedDataArray
~MappedDataArray() override
Definition: mapped_data_array.h:792
bdm::GetDataMemberForVis::operator()
std::enable_if< GetDataType< TTDataMember >)==DataType::kArray, TReturn >::type operator()(Agent *agent) const
Definition: mapped_data_array.h:81
bdm::MappedDataArray::InsertNextValue
vtkIdType InsertNextValue(TScalar v) final
Definition: mapped_data_array.h:770
bdm::GetDataMemberForVis::kDefault
@ kDefault
Definition: mapped_data_array.h:55
bdm::MappedDataArray::data_
std::vector< TScalar > data_
If mode_ is kCopy or kCache, data from Agents is written to data_.
Definition: mapped_data_array.h:218
agent.h
bdm::MappedDataArray::Initialize
void Initialize() final
Definition: mapped_data_array.h:311
bdm::MappedDataArray::SetNumberOfTuples
void SetNumberOfTuples(vtkIdType number) final
Definition: mapped_data_array.h:573
bdm::MappedDataArray::GetTypedTuple
void GetTypedTuple(vtkIdType idx, TScalar *t) const final
Definition: mapped_data_array.h:547
bdm::MappedDataArrayInterface::~MappedDataArrayInterface
virtual ~MappedDataArrayInterface()=default
bdm::MappedDataArray::GetVariantValue
vtkVariant GetVariantValue(vtkIdType idx) final
Definition: mapped_data_array.h:399
bdm::MappedDataArray::is_matching_
std::vector< uint64_t > is_matching_
Definition: mapped_data_array.h:216
bdm::MappedDataArray::Squeeze
void Squeeze() final
Definition: mapped_data_array.h:359
std
Definition: agent_uid.h:117
bdm::GetDataMemberForVis::kAgentUid
@ kAgentUid
Definition: mapped_data_array.h:55
bdm::MappedDataArray::ValueType
Superclass::ValueType ValueType
Definition: mapped_data_array.h:134
bdm::MappedDataArray::GetTuples
void GetTuples(vtkIdList *pt_ids, vtkAbstractArray *output) final
Definition: mapped_data_array.h:319
bdm::GetDataMemberForVis::temp_values_
std::vector< std::array< TempValueType, 64/sizeof(TempValueType)> > temp_values_
Thread local storage for temporary values.
Definition: mapped_data_array.h:44
bdm::MappedDataArray::InsertTypedTuple
void InsertTypedTuple(vtkIdType i, const TScalar *t) final
Definition: mapped_data_array.h:746
bdm::IsArray
Checks whether T is std::array or bdm::MathArray.
Definition: type.h:49
bdm::MappedDataArray::RemoveLastTuple
void RemoveLastTuple() final
Definition: mapped_data_array.h:731
param.h
bdm::MappedDataArray::GetTuple
double * GetTuple(vtkIdType i) final
Definition: mapped_data_array.h:412
bdm::MappedDataArray::GetValue
ValueType GetValue(vtkIdType idx) const final
Definition: mapped_data_array.h:481
bdm::MappedDataArray::SetValue
void SetValue(vtkIdType idx, TScalar value) final
Definition: mapped_data_array.h:762
bdm::MappedDataArray::SetTypedTuple
void SetTypedTuple(vtkIdType i, const TScalar *t) final
Definition: mapped_data_array.h:738
bdm::Param
Definition: param.h:35
bdm::is_agent_ptr
Definition: agent_pointer.h:294
bdm::MappedDataArray::Resize
int Resize(vtkIdType num_tuples) final
Definition: mapped_data_array.h:566
bdm::MappedDataArray::SetTuple
void SetTuple(vtkIdType i, vtkIdType j, vtkAbstractArray *source) final
Definition: mapped_data_array.h:581
bdm::MappedDataArray::InsertTuples
void InsertTuples(vtkIdList *dst_ids, vtkIdList *src_ids, vtkAbstractArray *source) final
Definition: mapped_data_array.h:629
functor.h
bdm::ThreadInfo::GetUniversalThreadId
uint64_t GetUniversalThreadId() const
Definition: thread_info.cc:26
bdm::GetDataMemberForVis::GetDataMemberForVis
GetDataMemberForVis()
Definition: mapped_data_array.h:46
bdm::MappedDataArray::end_
uint64_t end_
Definition: mapped_data_array.h:206
bdm::MappedDataArray::InsertValue
void InsertValue(vtkIdType idx, TScalar v) final
Definition: mapped_data_array.h:778