PSCF v1.1
RFieldDft.h
1#ifndef PSPC_R_FIELD_DFT_H
2#define PSPC_R_FIELD_DFT_H
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 "Field.h"
12#include <pscf/math/IntVec.h>
13#include <util/global.h>
14
15#include <fftw3.h>
16
17namespace Pscf {
18namespace Pspc
19{
20
21 using namespace Util;
22 using namespace Pscf;
23
29 template <int D>
30 class RFieldDft : public Field<fftw_complex>
31 {
32
33 public:
34
38 RFieldDft();
39
47 RFieldDft(RFieldDft<D> const & other);
48
54 virtual ~RFieldDft();
55
67
68 using Field<fftw_complex>::allocate;
69
77 void allocate(IntVec<D> const & meshDimensions);
78
82 IntVec<D> const & meshDimensions() const;
83
91 IntVec<D> const & dftDimensions() const;
92
99 template <class Archive>
100 void serialize(Archive& ar, const unsigned int version);
101
102 private:
103
104 // Vector containing number of grid points in each direction.
105 IntVec<D> meshDimensions_;
106
107 // Vector containing dimensions of dft (Fourier) grid.
108 IntVec<D> dftDimensions_;
109
110 };
111
112 /*
113 * Allocate the underlying C array for an FFT grid.
114 */
115 template <int D>
116 void RFieldDft<D>::allocate(IntVec<D> const & meshDimensions)
117 {
118 int size = 1;
119 for (int i = 0; i < D; ++i) {
120 UTIL_CHECK(meshDimensions[i] > 0);
121 meshDimensions_[i] = meshDimensions[i];
122 if (i < D - 1) {
123 dftDimensions_[i] = meshDimensions[i];
124 size *= meshDimensions[i];
125 } else {
126 dftDimensions_[i] = (meshDimensions[i]/2 + 1);
127 size *= dftDimensions_[i];
128 }
129 }
131 }
132
133 /*
134 * Return mesh dimensions by constant reference.
135 */
136 template <int D>
138 { return meshDimensions_; }
139
140 /*
141 * Return dimensions of dft grid by constant reference.
142 */
143 template <int D>
145 { return dftDimensions_; }
146
147 /*
148 * Serialize a Field to/from an Archive.
149 */
150 template <int D>
151 template <class Archive>
152 void RFieldDft<D>::serialize(Archive& ar, const unsigned int version)
153 {
155 ar & meshDimensions_;
156 ar & dftDimensions_;
157 }
158
159 #ifndef PSPC_R_FIELD_DFT_TPP
160 extern template class RFieldDft<1>;
161 extern template class RFieldDft<2>;
162 extern template class RFieldDft<3>;
163 #endif
164
165}
166}
167#endif
Base class template for a field defined on a spatial grid.
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition: IntVec.h:27
void allocate(int capacity)
Allocate the underlying C array.
Definition: Field.tpp:51
void serialize(Archive &ar, const unsigned int version)
Serialize a Field to/from an Archive.
Fourier transform of a real field on an FFT mesh.
Definition: RFieldDft.h:31
IntVec< D > const & meshDimensions() const
Return vector of spatial mesh dimensions by constant reference.
Definition: RFieldDft.h:137
virtual ~RFieldDft()
Destructor.
Definition: RFieldDft.tpp:30
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array for an FFT grid.
Definition: RFieldDft.h:116
RFieldDft< D > & operator=(RFieldDft< D > const &other)
Assignment operator.
Definition: RFieldDft.tpp:68
void serialize(Archive &ar, const unsigned int version)
Serialize a Field to/from an Archive.
Definition: RFieldDft.h:152
RFieldDft()
Default constructor.
Definition: RFieldDft.tpp:22
IntVec< D > const & dftDimensions() const
Return vector of dft (Fourier) grid dimensions by constant reference.
Definition: RFieldDft.h:144
File containing preprocessor macros for error handling.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition: global.h:68
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1