PSCF v1.1
RDField.h
1#ifndef PSPG_R_DFIELD_H
2#define PSPG_R_DFIELD_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 "DField.h"
12#include <pscf/math/IntVec.h>
13#include <util/global.h>
14#include <pspg/math/GpuResources.h>
15
16#include <cufft.h>
17
18namespace Pscf {
19namespace Pspg
20{
21
22 using namespace Util;
23 using namespace Pscf;
24
32 template <int D>
33 class RDField : public DField<cudaReal>
34 {
35
36 public:
37
41 RDField();
42
50 RDField(const RDField& other);
51
57 virtual ~RDField();
58
69 RDField& operator = (const RDField& other);
70
79
83 const IntVec<D>& meshDimensions() const;
84
91 template <class Archive>
92 void serialize(Archive& ar, const unsigned int version);
93
94 using DField<cudaReal>::allocate;
95 using DField<cudaReal>::operator=;
96
97 private:
98
99 // Vector containing number of grid points in each direction.
100 IntVec<D> meshDimensions_;
101
102 };
103
104 /*
105 * Allocate the underlying C array for an FFT grid.
106 */
107 template <int D>
108 void RDField<D>::allocate(const IntVec<D>& meshDimensions)
109 {
110 int size = 1;
111 for (int i = 0; i < D; ++i) {
112 UTIL_CHECK(meshDimensions[i] > 0);
113 meshDimensions_[i] = meshDimensions[i];
114 size *= meshDimensions[i];
115 }
117 }
118
119 /*
120 * Return mesh dimensions by constant reference.
121 */
122 template <int D>
124 { return meshDimensions_; }
125
126 /*
127 * Serialize a Field to/from an Archive.
128 */
129 template <int D>
130 template <class Archive>
131 void RDField<D>::serialize(Archive& ar, const unsigned int version)
132 {
133 int capacity;
134 if (Archive::is_saving()) {
135 capacity = capacity_;
136 }
137 ar & capacity;
138 if (Archive::is_loading()) {
139 if (!isAllocated()) {
140 if (capacity > 0) {
141 allocate(capacity);
142 }
143 } else {
144 if (capacity != capacity_) {
145 UTIL_THROW("Inconsistent Field capacities");
146 }
147 }
148 }
149
150 if (isAllocated()) {
151 double* tempData = new double[capacity];
152 cudaMemcpy(tempData, data_, capacity * sizeof(cudaReal),
153 cudaMemcpyDeviceToHost);
154 for (int i = 0; i < capacity_; ++i) {
155 ar & tempData[i];
156 }
157 delete[] tempData;
158 }
159 ar & meshDimensions_;
160 }
161
162
163
164
165
166}
167}
168#include "RDField.tpp"
169#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition: IntVec.h:27
Dynamic array on the GPU with alligned data.
Definition: DField.h:30
void allocate(int capacity)
Allocate the underlying C array on the device.
Definition: DField.tpp:52
Field of real single precision values on an FFT mesh on a device.
Definition: RDField.h:34
void allocate(const IntVec< D > &meshDimensions)
Allocate the underlying C array for an FFT grid.
Definition: RDField.h:108
RDField()
Default constructor.
Definition: RDField.tpp:23
void serialize(Archive &ar, const unsigned int version)
Serialize a Field to/from an Archive.
Definition: RDField.h:131
RDField & operator=(const RDField &other)
Assignment operator.
Definition: RDField.tpp:60
virtual ~RDField()
Destructor.
Definition: RDField.tpp:31
const IntVec< D > & meshDimensions() const
Return mesh dimensions by constant reference.
Definition: RDField.h:123
File containing preprocessor macros for error handling.
#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
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1