PSCF v1.1
RField.tpp
1#ifndef PSPC_R_FIELD_TPP
2#define PSPC_R_FIELD_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 "RField.h"
12
13namespace Pscf {
14namespace Pspc
15{
16
17 using namespace Util;
18
22 template <int D>
24 : Field<double>()
25 {}
26
27 /*
28 * Destructor.
29 */
30 template <int D>
32 {}
33
34 /*
35 * Copy constructor.
36 *
37 * Allocates new memory and copies all elements by value.
38 *
39 *\param other the Field to be copied.
40 */
41 template <int D>
43 : Field<double>(),
44 meshDimensions_(0)
45 {
46 if (!other.isAllocated()) {
47 UTIL_THROW("Other Field must be allocated.");
48 }
49 data_ = (double*) fftw_malloc(sizeof(double)*other.capacity_);
50 capacity_ = other.capacity_;
51 for (int i = 0; i < capacity_; ++i) {
52 data_[i] = other.data_[i];
53 }
54 meshDimensions_ = other.meshDimensions_;
55 }
56
57 /*
58 * Assignment, element-by-element.
59 *
60 * This operator will allocate memory if not allocated previously.
61 *
62 * \throw Exception if other Field is not allocated.
63 * \throw Exception if both Fields are allocated with unequal capacities.
64 *
65 * \param other the rhs Field
66 */
67 template <int D>
69 {
70 // Check for self assignment
71 if (this == &other) return *this;
72
73 // Precondition
74 if (!other.isAllocated()) {
75 UTIL_THROW("Other Field must be allocated.");
76 }
77
78 if (!isAllocated()) {
79 allocate(other.capacity());
80 } else if (capacity_ != other.capacity_) {
81 UTIL_THROW("Cannot assign Fields of unequal capacity");
82 }
83
84 // Copy elements
85 for (int i = 0; i < capacity_; ++i) {
86 data_[i] = other[i];
87 }
88 meshDimensions_ = other.meshDimensions_;
89
90 return *this;
91 }
92
93 /*
94 * Allocate the underlying C array for an FFT grid.
95 */
96 template <int D>
97 void RField<D>::allocate(const IntVec<D>& meshDimensions)
98 {
99 int size = 1;
100 for (int i = 0; i < D; ++i) {
101 UTIL_CHECK(meshDimensions[i] > 0);
102 meshDimensions_[i] = meshDimensions[i];
103 size *= meshDimensions[i];
104 }
106 }
107
108}
109}
110#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
Field of real double precision values on an FFT mesh.
Definition: RField.h:29
RField & operator=(const RField &other)
Assignment operator.
Definition: RField.tpp:68
virtual ~RField()
Destructor.
Definition: RField.tpp:31
RField()
Default constructor.
Definition: RField.tpp:23
void allocate(const IntVec< D > &meshDimensions)
Allocate the underlying C array for an FFT grid.
Definition: RField.tpp:97
Data * data_
Pointer to an array of Data elements.
Definition: Array.h:109
int capacity() const
Return allocated size.
Definition: Array.h:159
int capacity_
Allocated size of the data_ array.
Definition: Array.h:112
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
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1