PSCF v1.2
cuda/CField.tpp
1#ifndef PRDC_CUDA_C_FIELD_TPP
2#define PRDC_CUDA_C_FIELD_TPP
3
4/*
5* PSCF Package
6*
7* Copyright 2016 - 2022, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include "CField.h"
12
13namespace Pscf {
14namespace Prdc {
15namespace Cuda {
16
17 using namespace Util;
18 using namespace Pscf;
19
23 template <int D>
25 : DeviceArray<cudaComplex>()
26 {}
27
31 template <int D>
32 CField<D>::CField(IntVec<D> const & meshDimensions)
33 : DeviceArray<cudaComplex>()
35
36 /*
37 * Destructor.
38 */
39 template <int D>
41 {}
42
43 /*
44 * Copy constructor.
45 */
46 template <int D>
47 CField<D>::CField(const CField<D>& other)
48 : DeviceArray<cudaComplex>(other),
49 meshDimensions_(0)
50 {
51 meshDimensions_ = other.meshDimensions_;
52 }
53
54 /*
55 * Assignment from another RField<D>.
56 */
57 template <int D>
58 CField<D>& CField<D>::operator = (const CField<D>& other)
59 {
61 meshDimensions_ = other.meshDimensions_;
62
63 return *this;
64 }
65
66 /*
67 * Assignment from RHS HostDArray<Data> host array.
68 */
69 template <int D>
71 {
72 // Preconditions: both arrays must be allocated with equal capacities
73 if (!other.isAllocated()) {
74 UTIL_THROW("Error: RHS HostDArray<cudaComplex> is not allocated.");
75 }
76 if (!isAllocated()) {
77 UTIL_THROW("Error: LHS CField<D> is not allocated.");
78 }
79 if (capacity_ != other.capacity()) {
80 UTIL_THROW("Cannot assign Fields of unequal capacity");
81 }
82
83 // Use base class assignment operator to copy elements
85
86 return *this;
87 }
88
89 /*
90 * Allocate the underlying C array sized for an associated mesh.
91 */
92 template <int D>
93 void CField<D>::allocate(IntVec<D> const & meshDimensions)
94 {
95 int size = 1;
96 for (int i = 0; i < D; ++i) {
97 UTIL_CHECK(meshDimensions[i] > 0);
98 meshDimensions_[i] = meshDimensions[i];
99 size *= meshDimensions[i];
100 }
102 }
103
104 /*
105 * Associate this object with a slice of another DeviceArray.
106 */
107 template <int D>
109 IntVec<D> const & meshDimensions)
110 {
111 int size = 1;
112 for (int i = 0; i < D; ++i) {
113 UTIL_CHECK(meshDimensions[i] > 0);
114 meshDimensions_[i] = meshDimensions[i];
115 size *= meshDimensions[i];
116 }
117 DeviceArray<cudaComplex>::associate(arr, beginId, size);
118 }
119
120}
121}
122}
123#endif
Dynamic array on the GPU device with aligned data.
Definition rpg/System.h:32
virtual DeviceArray< Data > & operator=(const DeviceArray< Data > &other)
Assignment operator, assign from another DeviceArray<Data> array.
void associate(DeviceArray< Data > &arr, int beginId, int capacity)
Associate this object with a slice of a different DeviceArray.
void allocate(int capacity)
Allocate the underlying C array on the device.
Template for dynamic array stored in host CPU memory.
Definition HostDArray.h:43
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Field of complex double precision values on an FFT mesh.
Definition cpu/CField.h:30
const IntVec< D > & meshDimensions() const
Return mesh dimensions by constant reference.
Definition cpu/CField.h:127
void allocate(const IntVec< D > &meshDimensions)
Allocate the underlying C array for an FFT grid.
CField< D > & operator=(CField< D > const &other)
Assignment operator, assignment from another CField<D>.
void associate(DeviceArray< cudaComplex > &arr, int beginId, IntVec< D > const &meshDimensions)
Associate this object with a slice of another DeviceArray.
virtual ~CField()
Destructor.
CField()
Default constructor.
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array for data on a regular mesh.
int capacity() const
Return allocated size.
Definition Array.h:159
bool isAllocated() const
Return true if this DArray has been allocated, false otherwise.
Definition DArray.h:247
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
#define UTIL_THROW(msg)
Macro for throwing an Exception, reporting function, file and line number.
Definition global.h:51
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.