ABLATE Source Documentation  0.12.33
pointData.hpp
1 #ifndef ABLATELIBRARY_POINTDATA_HPP
2 #define ABLATELIBRARY_POINTDATA_HPP
3 
4 #include <petsc.h>
5 #include <map>
6 #include "particles/field.hpp"
7 #include "utilities/petscUtilities.hpp"
8 
9 namespace ablate::particles::accessors {
13 template <class DataType>
14 struct Data {
16  DataType* values = nullptr;
17 
19  PetscInt numberComponents = 0;
20 
22  PetscInt dataSize = 0;
23 
25  PetscInt offset = 0;
26 
30  Data() = default;
31 
39  Data(DataType* values, PetscInt numberComponents, PetscInt dataSizeIn = 0, PetscInt offset = 0)
41 
50 
52  template <class IndexType>
53  inline DataType* operator[](IndexType particle) const {
54  return values + (particle * dataSize + offset);
55  }
56 
63  template <class IndexType>
64  inline DataType& operator()(IndexType particle) {
65  return values[particle * dataSize + offset];
66  }
67 
74  template <class IndexType, class DimType>
75  inline DataType& operator()(IndexType particle, DimType dim) {
76  return values[particle * dataSize + offset + dim];
77  }
78 
85  template <class DestinationDataType, class IndexType>
86  inline void CopyFrom(DestinationDataType* source, IndexType p) const {
87  for (PetscInt d = 0; d < numberComponents; d++) {
88  values[p * dataSize + offset + d] = source[d];
89  }
90  }
91 
98  template <class DestinationDataType, class IndexType>
99  inline void AddFrom(DestinationDataType* source, IndexType p) const {
100  for (PetscInt d = 0; d < numberComponents; d++) {
101  values[p * dataSize + offset + d] += source[d];
102  }
103  }
104 
111  template <class DestinationDataType, class IndexType>
112  inline void CopyAll(DestinationDataType* destination, IndexType np) const {
113  for (IndexType p = 0; p < np; p++) {
114  for (PetscInt d = 0; d < numberComponents; d++) {
115  destination[p * numberComponents + d] = values[p * dataSize + offset + d];
116  }
117  }
118  }
119 };
120 
124 using ConstPointData = Data<const PetscReal>;
125 using PointData = Data<PetscReal>;
126 
127 } // namespace ablate::particles::accessors
128 #endif // ABLATELIBRARY_POINTDATA_HPP
Definition: field.hpp:11
Definition: pointData.hpp:14
DataType & operator()(IndexType particle)
Definition: pointData.hpp:64
PetscInt numberComponents
The number of the components.
Definition: pointData.hpp:19
void CopyFrom(DestinationDataType *source, IndexType p) const
Definition: pointData.hpp:86
void CopyAll(DestinationDataType *destination, IndexType np) const
Definition: pointData.hpp:112
Data(DataType *values, PetscInt numberComponents, PetscInt dataSizeIn=0, PetscInt offset=0)
Definition: pointData.hpp:39
PetscInt dataSize
The size of the component for this data.
Definition: pointData.hpp:22
DataType * values
the array for the solution values
Definition: pointData.hpp:16
DataType & operator()(IndexType particle, DimType dim)
Definition: pointData.hpp:75
Data(DataType *values, const ablate::particles::Field &field)
Definition: pointData.hpp:49
PetscInt offset
The offset in the local array, 0 for aux, computed for sol.
Definition: pointData.hpp:25
void AddFrom(DestinationDataType *source, IndexType p) const
Definition: pointData.hpp:99
DataType * operator[](IndexType particle) const
Inline function to compute the memory address at this particle.
Definition: pointData.hpp:53