PSCF v1.1
DField.tpp
1#ifndef PSPG_DFIELD_TPP
2#define PSPG_DFIELD_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 "DField.h"
12#include <pspg/math/GpuResources.h>
13#include <cuda_runtime.h>
14
15
16
17namespace Pscf {
18namespace Pspg
19{
20
21 using namespace Util;
22
23 /*
24 * Default constructor.
25 */
26 template <typename Data>
28 : data_(0),
29 capacity_(0)
30 {}
31
32 /*
33 * Destructor.
34 */
35 template <typename Data>
37 {
38 if (isAllocated()) {
39 cudaFree(data_);
40 capacity_ = 0;
41 }
42 }
43
44 /*
45 * Allocate the underlying C array.
46 *
47 * Throw an Exception if the Field has already allocated.
48 *
49 * \param capacity number of elements to allocate.
50 */
51 template <typename Data>
52 void DField<Data>::allocate(int capacity)
53 {
54 if (isAllocated()) {
55 UTIL_THROW("Attempt to re-allocate a DField");
56 }
57 if (capacity <= 0) {
58 UTIL_THROW("Attempt to allocate with capacity <= 0");
59 }
60 gpuErrchk(cudaMalloc((void**) &data_, capacity * sizeof(Data)));
61 capacity_ = capacity;
62 }
63
64 /*
65 * Deallocate the underlying C array.
66 *
67 * Throw an Exception if this Field is not allocated.
68 */
69 template <typename Data>
71 {
72 if (!isAllocated()) {
73 UTIL_THROW("Array is not allocated");
74 }
75 cudaFree(data_);
76 capacity_ = 0;
77 }
78
79 /*
80 * Copy constructor.
81 *
82 * Allocates new memory and copies all elements by value.
83 *
84 *\param other the Field to be copied.
85 */
86 template <typename Data>
88 {
89 if (!other.isAllocated()) {
90 UTIL_THROW("Other Field must be allocated.");
91 }
92
93 allocate(other.capacity_);
94 cudaMemcpy(data_, other.cDField(),
95 capacity_ * sizeof(Data), cudaMemcpyDeviceToDevice);
97 }
98
99 /*
100 * Assignment.
101 *
102 * This operator will allocate memory if not allocated previously.
103 *
104 * \throw Exception if other Field is not allocated.
105 * \throw Exception if both Fields are allocated with unequal capacities.
106 *
107 * \param other the rhs Field
108 */
109 template <typename Data>
111 {
112 // Check for self assignment
113 if (this == &other) return *this;
114
115 // Precondition
116 if (!other.isAllocated()) {
117 UTIL_THROW("Other Field must be allocated.");
118 }
119
120 if (!isAllocated()) {
121 allocate(other.capacity());
122 } else if (capacity_ != other.capacity_) {
123 UTIL_THROW("Cannot assign Fields of unequal capacity");
124 }
125
126 // Copy elements
127 cudaMemcpy(data_, other.cDField(),
128 capacity_ * sizeof(Data), cudaMemcpyDeviceToDevice);
129
130 return *this;
131 }
132
133}
134}
135#endif
Dynamic array on the GPU with alligned data.
Definition: DField.h:30
virtual DField< Data > & operator=(const DField< Data > &other)
Assignment operator.
Definition: DField.tpp:110
Data * cDField()
Return pointer to underlying C array.
Definition: DField.h:119
void allocate(int capacity)
Allocate the underlying C array on the device.
Definition: DField.tpp:52
int capacity_
Allocated size of the data_ array.
Definition: DField.h:104
void deallocate()
Dellocate the underlying C array.
Definition: DField.tpp:70
virtual ~DField()
Destructor.
Definition: DField.tpp:36
DField()
Default constructor.
Definition: DField.tpp:27
int capacity() const
Return allocated size.
Definition: DField.h:112
bool isAllocated() const
Return true if the Field has been allocated, false otherwise.
Definition: DField.h:133
#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