PSCF v1.3
cuda/RFieldDft.tpp
1#ifndef PRDC_CUDA_R_FIELD_DFT_TPP
2#define PRDC_CUDA_R_FIELD_DFT_TPP
3
4/*
5* PSCF - Polymer Self-Consistent Field
6*
7* Copyright 2015 - 2025, 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#include "FFT.h"
13
14namespace Pscf {
15namespace Prdc {
16namespace Cuda {
17
18 using namespace Util;
19
20 /*
21 * Default constructor.
22 */
23 template <int D>
25 : DeviceArray<cudaComplex>()
26 {}
27
28 /*
29 * Allocating constructor (calls allocate).
30 */
31 template <int D>
35
36 /*
37 * Destructor.
38 */
39 template <int D>
41 {}
42
43 /*
44 * Copy constructor.
45 */
46 template <int D>
48 : DeviceArray<cudaComplex>(other)
49 {
50 meshDimensions_ = other.meshDimensions_;
51 dftDimensions_ = other.dftDimensions_;
52 }
53
54 /*
55 * Assignment, element-by-element from another RFieldDft<D>.
56 */
57 template <int D>
58 RFieldDft<D>& RFieldDft<D>::operator = (const RFieldDft<D>& other)
59 {
60 // Assign data and size of underlying array
62
63 // Assign more specialized data members
64 meshDimensions_ = other.meshDimensions_;
65 dftDimensions_ = other.dftDimensions_;
66
67 return *this;
68 }
69
70 /*
71 * Assignment of RFieldDft<D> from RHS HostDArray<Data> host array.
72 */
73 template <int D>
74 RFieldDft<D>&
76 {
77 // Preconditions: both arrays must be allocated with equal capacities
78 if (!other.isAllocated()) {
79 UTIL_THROW("Error: RHS HostDArray<cudaComplex> is not allocated.");
80 }
81 if (!isAllocated()) {
82 UTIL_THROW("Error: LHS RFieldDft<D> is not allocated.");
83 }
84 if (capacity_ != other.capacity()) {
85 UTIL_THROW("Cannot assign Fields of unequal capacity");
86 }
87
88 // Use base class assignment operator to copy elements
90
91 return *this;
92 }
93
94 /*
95 * Allocate underlying DeviceArray<cudaComplex> for the DFT mesh.
96 */
97 template <int D>
98 void RFieldDft<D>::allocate(const IntVec<D>& meshDimensions)
99 {
100 // Copy and validate dimensions of real space grid
101 for (int i = 0; i < D; ++i) {
102 UTIL_CHECK(meshDimensions[i] > 0);
103 meshDimensions_[i] = meshDimensions[i];
104 }
105
106 // Compute dimensions and size of Fourier space mesh
107 int size;
108 FFT<D>::computeKMesh(meshDimensions, dftDimensions_, size);
109
110 // Allocate complex array on the GPU with size of DFT mesh
112 }
113
114 /*
115 * Associate this object with a slice of another DeviceArray<cudaComplex>.
116 */
117 template <int D>
119 int beginId,
121 {
122 // Copy and validate dimensions of real space grid
123 for (int i = 0; i < D; ++i) {
125 meshDimensions_[i] = meshDimensions[i];
126 }
127
128 // Compute dimensions and size of Fourier space mesh
129 int size;
130 FFT<D>::computeKMesh(meshDimensions, dftDimensions_, size);
131
132 // Associate data with a slice of input array arr
133 DeviceArray<cudaComplex>::associate(arr, beginId, size);
134 }
135
136} // namespace Cuda
137} // namespace Prdc
138} // namespace Pscf
139#endif
Dynamic array on the GPU device with aligned data.
Definition DeviceArray.h:43
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.
static void computeKMesh(IntVec< D > const &rMeshDimensions, IntVec< D > &kMeshDimensions, int &kSize)
Compute dimensions and size of k-space mesh for DFT of real data.
IntVec< D > const & meshDimensions() const
Return vector of real-space 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:49
Fields, FFTs, and utilities for periodic boundary conditions (CUDA)
Definition Reduce.cpp:14
Periodic fields and crystallography.
Definition CField.cpp:11
PSCF package top-level namespace.
Definition param_pc.dox:1