PSCF v1.2
rpg/solvers/Block.h
1#ifndef RPG_BLOCK_H
2#define RPG_BLOCK_H
3
4/*
5* PSCF - Polymer Self-Consistent Field Theory
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 "Propagator.h" // base class argument
12#include <prdc/cuda/RField.h> // member
13#include <prdc/cuda/RFieldDft.h> // member
14#include <prdc/cuda/FFT.h> // member
15#include <prdc/cuda/FFTBatched.h> // member
16#include <rpg/solvers/WaveList.h>
17#include <prdc/crystal/UnitCell.h>
18#include <pscf/solvers/BlockTmpl.h> // base class template
19#include <util/containers/FSArray.h>
20
21namespace Pscf {
22
23 template <int D> class Mesh;
24
25namespace Rpg {
26
27 using namespace Util;
28 using namespace Pscf::Prdc;
29 using namespace Pscf::Prdc::Cuda;
30
31 // Wrappers for the CUDA kernels used internally
32
44 DeviceArray<cudaReal> const & d);
45
55 DeviceArray<cudaReal> const & qr,
56 DeviceArray<cudaReal> const & qr2,
57 DeviceArray<cudaReal> const & expW2);
58
68 DeviceArray<cudaReal> const & c, cudaReal const d);
69
78 template <int D>
79 class Block : public BlockTmpl< Propagator<D> >
80 {
81
82 public:
83
87 Block();
88
92 ~Block();
93
104 void associate(Mesh<D> const & mesh, FFT<D> const & fft,
105 UnitCell<D> const & cell, WaveList<D>& wavelist);
106
115 void allocate(double ds, bool useBatchedFFT = true);
116
124 void clearUnitCellData();
125
131 void setLength(double newLength);
132
138 void setKuhn(double kuhn);
139
145 void setupSolver(RField<D> const & w);
146
153 void step(RField<D> const & q, RField<D>& qNew);
154
167 void computeConcentration(double prefactor);
168
176 void computeStress(double prefactor);
177
183 double stress(int n);
184
188 Mesh<D> const & mesh() const;
189
193 FFT<D> const & fft() const;
194
198 double ds() const;
199
203 int ns() const;
204
205 // Functions with non-dependent names from BlockTmpl<Propagator<D>>
211
212 // Functions with non-dependent names from BlockDescriptor
222
223 private:
224
226 FFTBatched<D> fftBatchedAll_;
227
229 FFTBatched<D> fftBatchedPair_;
230
232 FSArray<double, 6> stress_;
233
235 RField<D> expKsq_;
236
238 RField<D> expKsq2_;
239
241 RField<D> expW_;
242
244 RField<D> expW2_;
245
253 DeviceArray<cudaReal> qrPair_;
254
263
265 DeviceArray<cudaComplex> qkBatched_;
266
268 DeviceArray<cudaComplex> qk2Batched_;
269
270 // Propagators on r-grid
271 RField<D> q1_;
272 RField<D> q2_;
273
275 HostDArray<cudaReal> expKsq_h_;
276
278 HostDArray<cudaReal> expKsq2_h_;
279
281 Mesh<D> const * meshPtr_;
282
284 FFT<D> const * fftPtr_;
285
287 UnitCell<D> const * unitCellPtr_;
288
290 WaveList<D> * waveListPtr_;
291
293 IntVec<D> kMeshDimensions_;
294
296 int kSize_;
297
299 double ds_;
300
302 double dsTarget_;
303
305 int ns_;
306
308 bool isAllocated_;
309
311 bool hasExpKsq_;
312
314 bool useBatchedFFT_;
315
317 UnitCell<D> const & unitCell() const
318 { return *unitCellPtr_; }
319
321 WaveList<D> const & wavelist() const
322 { return *waveListPtr_; }
323
325 int nParams_;
326
328 void computeExpKsq();
329
330 };
331
332 // Inline member functions
333
334 // Get number of contour steps.
335 template <int D>
336 inline int Block<D>::ns() const
337 { return ns_; }
338
339 // Get contour length step size.
340 template <int D>
341 inline double Block<D>::ds() const
342 { return ds_; }
343
344 // Get derivative of free energy w/ respect to a unit cell parameter.
345 template <int D>
346 inline double Block<D>::stress(int n)
347 { return stress_[n]; }
348
349 // Get Mesh by reference.
350 template <int D>
351 inline Mesh<D> const & Block<D>::mesh() const
352 {
353 UTIL_ASSERT(meshPtr_);
354 return *meshPtr_;
355 }
356
357 // Get FFT by reference.
358 template <int D>
359 inline FFT<D> const & Block<D>::fft() const
360 {
361 UTIL_ASSERT(fftPtr_);
362 return *fftPtr_;
363 }
364
365 #ifndef RPG_BLOCK_TPP
366 // Suppress implicit instantiation
367 extern template class Block<1>;
368 extern template class Block<2>;
369 extern template class Block<3>;
370 #endif
371
372}
373}
374#endif
double length() const
Get the length (number of monomers) in this block.
int vertexId(int i) const
Get id of an associated vertex.
void setId(int id)
Set the id for this block.
int monomerId() const
Get the monomer type id.
const Pair< int > & vertexIds() const
Get the pair of associated vertex ids.
virtual void setLength(double length)
Set the length of this block.
void setVertexIds(int vertexAId, int vertexBId)
Set indices of associated vertices.
int id() const
Get the id of this block.
void setMonomerId(int monomerId)
Set the monomer type id.
Class template for a block in a block copolymer.
Definition BlockTmpl.h:106
Propagator< D > & propagator(int directionId)
Definition BlockTmpl.h:191
Propagator< D >::CField & cField()
Definition BlockTmpl.h:207
Dynamic array on the GPU device with aligned data.
Definition rpg/System.h:32
Template for dynamic array stored in host CPU memory.
Definition HostDArray.h:43
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Description of a regular grid of points in a periodic domain.
Fourier transform wrapper.
Field of real double precision values on an FFT mesh.
Batched Fourier transform wrapper for real data.
Definition FFTBatched.h:32
Base template for UnitCell<D> classes, D=1, 2 or 3.
Definition rpg/System.h:34
Block within a branched polymer.
void setLength(double newLength)
Set or reset block length.
double ds() const
Contour length step size.
void computeConcentration(double prefactor)
Compute unnormalized concentration for block by integration.
double length() const
Get the length (number of monomers) in this block.
void allocate(double ds, bool useBatchedFFT=true)
Allocate internal data containers.
Mesh< D > const & mesh() const
Return associated spatial Mesh by const reference.
void setKuhn(double kuhn)
Set or reset monomer statistical segment length.
int ns() const
Number of contour length steps.
void associate(Mesh< D > const &mesh, FFT< D > const &fft, UnitCell< D > const &cell, WaveList< D > &wavelist)
Associate this object with a mesh, FFT, UnitCell, and WaveList object.
void computeStress(double prefactor)
Compute derivatives of free energy w/ respect to cell parameters.
void setupSolver(RField< D > const &w)
Set solver for this block.
double stress(int n)
Get derivative of free energy w/ respect to a unit cell parameter.
void step(RField< D > const &q, RField< D > &qNew)
Compute step of integration loop, from i to i+1.
void clearUnitCellData()
Clear all internal data that depends on lattice parameters.
FFT< D > const & fft() const
Return associated FFT<D> object by const reference.
MDE solver for one-direction of one block.
A fixed capacity (static) contiguous array with a variable logical size.
Definition rpg/System.h:28
#define UTIL_ASSERT(condition)
Assertion macro suitable for debugging serial or parallel code.
Definition global.h:75
Fields, FFTs, and utilities for periodic boundary conditions (CUDA)
Definition CField.cu:12
Periodic fields and crystallography.
Definition CField.cpp:11
void richardsonEx(DeviceArray< cudaReal > &qNew, DeviceArray< cudaReal > const &qr, DeviceArray< cudaReal > const &qr2, DeviceArray< cudaReal > const &expW2)
Performs qNew = (4 * (qr2 * expW2) - qr) / 3 elementwise, kernel wrapper.
void realMulVConjVV(DeviceArray< cudaReal > &a, DeviceArray< cudaComplex > const &b, DeviceArray< cudaComplex > const &c, DeviceArray< cudaReal > const &d)
Element-wise calculation of a = real(b * conj(c) * d), kernel wrapper.
void addEqMulVVc(DeviceArray< cudaReal > &a, DeviceArray< cudaReal > const &b, DeviceArray< cudaReal > const &c, cudaReal const d)
Performs a[i] += b[i] * c[i] * d, kernel wrapper.
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.