PSCF v1.3
cpu/RFieldDft.tpp
1#ifndef PRDC_CPU_R_FIELD_DFT_TPP
2#define PRDC_CPU_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 <prdc/cpu/FFT.h>
13
14namespace Pscf {
15namespace Prdc {
16namespace Cpu {
17
18 using namespace Util;
19
20 /*
21 * Default constructor.
22 */
23 template <int D>
25 : FftwDArray<fftw_complex>(),
26 meshDimensions_(0),
27 dftDimensions_(0)
28 {}
29
30 /*
31 * Destructor.
32 */
33 template <int D>
36
37 /*
38 * Copy constructor.
39 *
40 * Allocates new memory and copies all elements by value.
41 */
42 template <int D>
44 : FftwDArray<fftw_complex>(),
45 meshDimensions_(0),
46 dftDimensions_(0)
47 {
48 if (!other.isAllocated()) {
49 UTIL_THROW("Other Field must be allocated.");
50 }
51 allocate(other.meshDimensions_);
52 for (int i = 0; i < capacity_; ++i) {
53 data_[i][0] = other.data_[i][0];
54 data_[i][1] = other.data_[i][1];
55 }
56 }
57
58 /*
59 * Assignment, element-by-element.
60 *
61 * This operator will allocate memory if not allocated previously.
62 *
63 * \throw Exception if other Field is not allocated.
64 * \throw Exception if both Fields are allocated with unequal capacities.
65 *
66 * \param other the rhs Field
67 */
68 template <int D>
70 {
71 // Check for self assignment
72 if (this == &other) return *this;
73
74 // Precondition
75 if (!other.isAllocated()) {
76 UTIL_THROW("Other RFieldDft must be allocated.");
77 }
78
79 // Allocate if necessary, check dimensions
80 if (!isAllocated()) {
81 allocate(other.meshDimensions_);
82 }
84 UTIL_CHECK(meshDimensions_ == other.meshDimensions_);
85 UTIL_CHECK(dftDimensions_ == other.dftDimensions_);
86
87 // Copy elements
88 for (int i = 0; i < capacity_; ++i) {
89 data_[i][0] = other.data_[i][0];
90 data_[i][1] = other.data_[i][1];
91 }
92
93 return *this;
94 }
95
96 /*
97 * Allocate the underlying array for an FFT grid.
98 */
99 template <int D>
101 {
102 // Copy real space grid dimensions
103 for (int i = 0; i < D; ++i) {
105 meshDimensions_[i] = meshDimensions[i];
106 }
107
108 // Compute dimensions and size for k-space grid
109 int size;
110 FFT<D>::computeKMesh(meshDimensions, dftDimensions_, size);
111
112 // Allocate memory
114 }
115
116 /*
117 * Dellocate the underlying C array and clear dimensions.
118 */
119 template <int D>
121 {
123 for (int i = 0; i < D; ++i) {
124 meshDimensions_[i] = 0;
125 dftDimensions_[i] = 0;
126 }
127 }
128
129}
130}
131}
132#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
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.
Definition cpu/FFT.tpp:262
void allocate(int capacity)
Allocate the underlying C array.
virtual void deallocate()
Dellocate the underlying C array.
bool isAllocated() const
Return true if the FftwDArray has been allocated, false otherwise.
Definition FftwDArray.h:102
virtual void deallocate()
Deallocate underlying C array and clear mesh dimensions.
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array and set mesh dimensions.
virtual ~RFieldDft()
Destructor.
RFieldDft< D > & operator=(RFieldDft< D > const &other)
Assignment operator.
RFieldDft()
Default constructor.
IntVec< D > const & meshDimensions() const
Return vector of spatial mesh dimensions by constant reference.
fftw_complex * data_
Definition Array.h:109
#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 and FFTs for periodic boundary conditions (CPU)
Definition CField.cpp:12
Periodic fields and crystallography.
Definition CField.cpp:11
PSCF package top-level namespace.
Definition param_pc.dox:1