PSCF v1.2
cuda/RFieldDft.tpp
1#ifndef PRDC_CUDA_R_FIELD_DFT_TPP
2#define PRDC_CUDA_R_FIELD_DFT_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 "RFieldDft.h"
12
13namespace Pscf {
14namespace Prdc {
15namespace Cuda {
16
17 using namespace Util;
18
19 /*
20 * Default constructor.
21 */
22 template <int D>
24 : DeviceArray<cudaComplex>()
25 {}
26
27 /*
28 * Allocating constructor (calls allocate).
29 */
30 template <int D>
31 RFieldDft<D>::RFieldDft(IntVec<D> const & meshDimensions)
32 : DeviceArray<cudaComplex>()
34
35 /*
36 * Destructor.
37 */
38 template <int D>
40 {}
41
42 /*
43 * Copy constructor.
44 */
45 template <int D>
47 : DeviceArray<cudaComplex>(other)
48 {
49 meshDimensions_ = other.meshDimensions_;
50 dftDimensions_ = other.dftDimensions_;
51 }
52
53 /*
54 * Assignment, element-by-element from another RFieldDft<D>.
55 */
56 template <int D>
57 RFieldDft<D>& RFieldDft<D>::operator = (const RFieldDft<D>& other)
58 {
59
61 meshDimensions_ = other.meshDimensions_;
62 dftDimensions_ = other.dftDimensions_;
63
64 return *this;
65 }
66
67 /*
68 * Assignment of RFieldDft<D> from RHS HostDArray<Data> host array.
69 */
70 template <int D>
71 RFieldDft<D>&
73 {
74 // Preconditions: both arrays must be allocated with equal capacities
75 if (!other.isAllocated()) {
76 UTIL_THROW("Error: RHS HostDArray<cudaComplex> is not allocated.");
77 }
78 if (!isAllocated()) {
79 UTIL_THROW("Error: LHS RFieldDft<D> is not allocated.");
80 }
81 if (capacity_ != other.capacity()) {
82 UTIL_THROW("Cannot assign Fields of unequal capacity");
83 }
84
85 // Use base class assignment operator to copy elements
87
88 return *this;
89 }
90
91 /*
92 * Allocate the underlying C array for the dftDimensions mesh.
93 */
94 template <int D>
95 void RFieldDft<D>::allocate(const IntVec<D>& meshDimensions)
96 {
97 int size = 1;
98 for (int i = 0; i < D; ++i) {
99 UTIL_CHECK(meshDimensions[i] > 0);
100 meshDimensions_[i] = meshDimensions[i];
101 if (i < D - 1) {
102 dftDimensions_[i] = meshDimensions[i];
103 size *= meshDimensions[i];
104 } else {
105 dftDimensions_[i] = (meshDimensions[i]/2 + 1);
106 size *= (meshDimensions[i]/2 + 1);
107 }
108 }
109 // Note: size denotes size of mesh with dftDimensions
111 }
112
113 /*
114 * Associate this object with a slice of another DeviceArray.
115 */
116 template <int D>
118 IntVec<D> const & meshDimensions)
119 {
120 int size = 1;
121 for (int i = 0; i < D; ++i) {
122 UTIL_CHECK(meshDimensions[i] > 0);
123 meshDimensions_[i] = meshDimensions[i];
124 if (i < D - 1) {
125 dftDimensions_[i] = meshDimensions[i];
126 size *= meshDimensions[i];
127 } else {
128 dftDimensions_[i] = (meshDimensions[i]/2 + 1);
129 size *= (meshDimensions[i]/2 + 1);
130 }
131 }
132 // Note: size denotes size of mesh with dftDimensions
133 DeviceArray<cudaComplex>::associate(arr, beginId, size);
134 }
135
136}
137}
138}
139#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
Fourier transform of a real field on an FFT mesh.
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array and set mesh dimensions.
IntVec< D > const & meshDimensions() const
Return vector of spatial mesh dimensions by constant reference.
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array for an FFT grid.
virtual ~RFieldDft()
Destructor.
void associate(DeviceArray< cudaComplex > &arr, int beginId, IntVec< D > const &meshDimensions)
Associate this object with a slice of another DeviceArray.
RFieldDft< D > & operator=(RFieldDft< D > const &other)
Assignment operator, assignment from another RFieldDft<D>.
RFieldDft()
Default constructor.
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.