PSCF v1.3.3
rpc/solvers/Block.tpp
1#ifndef RPC_BLOCK_TPP
2#define RPC_BLOCK_TPP
3
4/*
5* PSCF - Polymer Self-Consistent Field
6*
7* Copyright 2015 - 2025, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include "Block.h"
12#include "Propagator.h"
13
14#include <prdc/cpu/WaveList.h>
15#include <prdc/cpu/FFT.h>
16#include <prdc/crystal/UnitCell.h>
17#include <prdc/crystal/shiftToMinimum.h>
18
19#include <pscf/chem/Edge.h>
20#include <pscf/mesh/Mesh.h>
21#include <pscf/mesh/MeshIterator.h>
22#include <pscf/math/IntVec.h>
23#include <pscf/solvers/BlockTmpl.tpp>
24
25#include <util/containers/DMatrix.h>
26#include <util/containers/DArray.h>
27#include <util/containers/FSArray.h>
28
29namespace Pscf {
30namespace Rpc {
31
32 // Namespaces from which names may be used without qualification
33 using namespace Util;
34 using namespace Pscf::Prdc;
35 using namespace Pscf::Prdc::Cpu;
36
37 /*
38 * Constructor.
39 */
40 template <int D>
42 : meshPtr_(nullptr),
43 fftPtr_(nullptr),
44 unitCellPtr_(nullptr),
45 waveListPtr_(nullptr),
46 kMeshDimensions_(-1),
47 kSize_(-1),
48 ds_(-1.0),
49 dsTarget_(-1.0),
50 ns_(-1),
51 isAllocated_(false),
52 hasExpKsq_(false)
53 {
54 propagator(0).setBlock(*this);
55 propagator(1).setBlock(*this);
56 }
57
58 /*
59 * Destructor.
60 */
61 template <int D>
64
65 /*
66 * Store addresses of mesh, FFT and unit cell.
67 */
68 template <int D>
69 void Block<D>::associate(Mesh<D> const & mesh,
70 FFT<D> const& fft,
71 UnitCell<D> const& cell,
72 WaveList<D>& wavelist)
73 {
74 // Preconditions
75 UTIL_CHECK(!isAllocated_);
76
77 // Set pointers to mesh and fft
78 meshPtr_ = &mesh;
79 fftPtr_ = &fft;
80 unitCellPtr_ = &cell;
81 waveListPtr_ = &wavelist;
82
83 hasExpKsq_ = false;
84 }
85
86 /*
87 * Compute number of contour steps and allocate all memory.
88 */
89 template <int D>
91 {
92 UTIL_CHECK(ds > 0.0);
93 UTIL_CHECK(meshPtr_);
94 UTIL_CHECK(fftPtr_);
95 UTIL_CHECK(unitCellPtr_);
96 UTIL_CHECK(waveListPtr_);
97 UTIL_CHECK(mesh().size() > 1);
98 UTIL_CHECK(mesh().dimensions() == fft().meshDimensions());
99 UTIL_CHECK(!isAllocated_);
100
101 // Compute DFT grid dimensions kMeshDimensions_ and size kSize_
102 FFT<D>::computeKMesh(mesh().dimensions(), kMeshDimensions_, kSize_);
103
104 // Allocate work arrays for MDE solution
105 expKsq_.allocate(kMeshDimensions_);
106 expKsq2_.allocate(kMeshDimensions_);
107 expW_.allocate(mesh().dimensions());
108 qr_.allocate(mesh().dimensions());
109 qk_.allocate(mesh().dimensions());
110 qr2_.allocate(mesh().dimensions());
111 qk2_.allocate(mesh().dimensions());
113 expW2_.allocate(mesh().dimensions());
114 } else {
115 expWInv_.allocate(mesh().dimensions());
116 }
117
118 // Allocate block concentration field
119 cField().allocate(mesh().dimensions());
120
121 dsTarget_ = ds;
122
123 // Compute ns_
125 // Set contour length discretization for this block
126 UTIL_CHECK(length() > 0.0);
127 int tempNs;
128 tempNs = floor(length() / (2.0 *ds) + 0.5);
129 if (tempNs == 0) {
130 tempNs = 1; // ensure ns_ >= 3
131 }
132 ns_ = 2*tempNs + 1;
133 ds_ = length()/double(ns_-1);
134 } else
135 if (PolymerModel::isBead()) {
136 ds_ = 1.0;
137 ns_ = nBead() + 2;
138 }
139
140 // Allocate memory for solutions to MDE (requires ns_)
141 propagator(0).allocate(ns_, mesh());
142 propagator(1).allocate(ns_, mesh());
143
144 isAllocated_ = true;
145 hasExpKsq_ = false;
146 }
147
148 /*
149 * Set or reset the the block length.
150 */
151 template <int D>
152 void Block<D>::setLength(double newLength)
153 {
155 Edge::setLength(newLength);
156
157 if (isAllocated_) {
158
159 int oldNs = ns_;
160
161 // Reset contour length discretization
162 UTIL_CHECK(dsTarget_ > 0);
163 int tempNs;
164 tempNs = floor( length()/(2.0 *dsTarget_) + 0.5 );
165 if (tempNs == 0) {
166 tempNs = 1;
167 }
168 ns_ = 2*tempNs + 1;
169 ds_ = length()/double(ns_-1);
170
171 // Reallocate propagators if ns_ has changed
172 if (oldNs != ns_) {
173 propagator(0).reallocate(ns_);
174 propagator(1).reallocate(ns_);
175 }
176
177 }
178 hasExpKsq_ = false;
179 }
180
181 /*
182 * Set or reset the the block length.
183 */
184 template <int D>
186 {
188 hasExpKsq_ = false;
189 }
190
191 /*
192 * Clear all internal data that depends on the unit cell parameters.
193 */
194 template <int D>
196 {
197 hasExpKsq_ = false;
198 stress_.clear();
199 }
200
201 /*
202 * Compute all elements of expKsq_ and expKsq2_ arrays
203 */
204 template <int D>
205 void Block<D>::computeExpKsq()
206 {
207 UTIL_CHECK(isAllocated_);
208 UTIL_CHECK(waveListPtr_);
209
210 bool isThread = PolymerModel::isThread();
211 double bSqFactor = -1.0*kuhn()*kuhn() / 6.0;
212 if (isThread) {
213 bSqFactor *= ds_;
214 }
215
216 // Calculate KSq if necessary
217 if (!waveListPtr_->hasKSq()) {
218 waveListPtr_->computeKSq();
219 }
220 RField<D> const & kSq = waveListPtr_->kSq();
221
222 MeshIterator<D> iter;
223 iter.setDimensions(kMeshDimensions_);
224 double arg;
225 int i;
226 for (iter.begin(); !iter.atEnd(); ++iter) {
227 i = iter.rank();
228 arg = kSq[i]*bSqFactor;
229 expKsq_[i] = exp(arg);
230 expKsq2_[i] = exp(0.5*arg);
231 }
232
233 hasExpKsq_ = true;
234 }
235
236 /*
237 * Setup the the step algorithm for a specific field configuration.
238 */
239 template <int D>
240 void
242 {
243 // Preconditions
244 int nx = mesh().size();
245 UTIL_CHECK(nx > 0);
246 UTIL_CHECK(isAllocated_);
247
248 // Compute expW arrays
249 double arg;
251 double c = -0.5*ds_;
252 for (int i = 0; i < nx; ++i) {
253 arg = c*w[i];
254 expW_[i] = exp(arg);
255 expW2_[i] = exp(0.5*arg);
256 }
257 } else
258 if (PolymerModel::isBead()) {
259 for (int i = 0; i < nx; ++i) {
260 arg = -w[i];
261 expW_[i] = exp(arg);
262 expWInv_[i] = 1.0/expW_[i];
263 }
264 }
265
266 // Compute expKsq arrays if necessary
267 if (!hasExpKsq_) {
268 computeExpKsq();
269 }
270
271 }
272
273 /*
274 * Propagate solution by one step for the thread model.
275 */
276 template <int D>
277 void Block<D>::stepThread(RField<D> const & q, RField<D>& qout) const
278 {
280
281 // Preconditions on mesh and fft
282 int nx = mesh().size();
283 UTIL_CHECK(nx > 0);
284 UTIL_CHECK(fft().isSetup());
285 UTIL_CHECK(mesh().dimensions() == fft().meshDimensions());
286
287 // Internal preconditions
288 UTIL_CHECK(isAllocated_);
289 int nk = qk_.capacity();
290 UTIL_CHECK(nk > 0);
291 UTIL_CHECK(qr_.capacity() == nx);
292 UTIL_CHECK(expW_.capacity() == nx);
293 UTIL_CHECK(expW2_.capacity() == nx);
294 UTIL_CHECK(expKsq_.capacity() == nk);
295 UTIL_CHECK(expKsq2_.capacity() == nk);
296 UTIL_CHECK(hasExpKsq_);
297
298 // Preconditions on parameters
300 UTIL_CHECK(q.capacity() == nx);
301 UTIL_CHECK(qout.isAllocated());
302 UTIL_CHECK(qout.capacity() == nx);
303
304 // Apply pseudo-spectral algorithm
305
306 // Step by ds/2 for qr_, step by ds/4 for qr2_
307 int i;
308 for (i = 0; i < nx; ++i) {
309 qr_[i] = q[i]*expW_[i];
310 qr2_[i] = q[i]*expW2_[i];
311 }
312 fft().forwardTransform(qr_, qk_);
313 fft().forwardTransform(qr2_, qk2_);
314 for (i = 0; i < nk; ++i) {
315 qk_[i][0] *= expKsq_[i];
316 qk_[i][1] *= expKsq_[i];
317 qk2_[i][0] *= expKsq2_[i];
318 qk2_[i][1] *= expKsq2_[i];
319 }
320 fft().inverseTransformUnsafe(qk_, qr_); // overwrites qk_
321 fft().inverseTransformUnsafe(qk2_, qr2_); // overwrites qk2_
322 for (i = 0; i < nx; ++i) {
323 qr_[i] = qr_[i]*expW_[i];
324 qr2_[i] = qr2_[i]*expW_[i];
325 }
326
327 // Above, multiplying qr2_ by expW_, rather than by expW2_, combines
328 // required multiplications by expW2_ at the end of first half-step
329 // and at the beginning of the second.
330
331 // Finish second half-step of ds/2 for qr2_
332 fft().forwardTransform(qr2_, qk2_);
333 for (i = 0; i < nk; ++i) {
334 qk2_[i][0] *= expKsq2_[i];
335 qk2_[i][1] *= expKsq2_[i];
336 }
337 fft().inverseTransformUnsafe(qk2_, qr2_); // overwrites qk2_
338 for (i = 0; i < nx; ++i) {
339 qr2_[i] = qr2_[i]*expW2_[i];
340 }
341
342 // Richardson extrapolation
343 for (i = 0; i < nx; ++i) {
344 qout[i] = (4.0*qr2_[i] - qr_[i])/3.0;
345 }
346
347 }
348
349 /*
350 * Apply one step of MDE solution for the bead model.
351 */
352 template <int D>
353 void Block<D>::stepBead(RField<D> const & q, RField<D>& qout) const
354 {
356 stepBondBead(q, qout);
357 stepFieldBead(qout);
358 }
359
360 /*
361 * Apply the bond operator for the bead model.
362 */
363 template <int D>
364 void Block<D>::stepBondBead(RField<D> const & q, RField<D>& qout) const
365 {
366 // Prereconditions
367 UTIL_CHECK(isAllocated_);
368 UTIL_CHECK(hasExpKsq_);
369 int nx = mesh().size();
370 UTIL_CHECK(nx > 0);
371 UTIL_CHECK(q.capacity() == nx);
372 UTIL_CHECK(qout.capacity() == nx);
373 int nk = qk_.capacity();
374 UTIL_CHECK(nk > 0);
375 UTIL_CHECK(expKsq_.capacity() == nk);
376
377 // Apply bond operator
378 fft().forwardTransform(q, qk_);
379 for (int i = 0; i < nk; ++i) {
380 qk_[i][0] *= expKsq_[i];
381 qk_[i][1] *= expKsq_[i];
382 }
383 fft().inverseTransformUnsafe(qk_, qout); // destroys qk_
384 }
385
386 /*
387 * Apply the half-bond operator for the bead model.
388 */
389 template <int D>
390 void Block<D>::stepHalfBondBead(RField<D> const & q, RField<D>& qout) const
391 {
392 // Preconditions
393 UTIL_CHECK(isAllocated_);
394 UTIL_CHECK(hasExpKsq_);
395 int nx = mesh().size();
396 UTIL_CHECK(nx > 0);
397 UTIL_CHECK(q.capacity() == nx);
398 UTIL_CHECK(qout.capacity() == nx);
399 int nk = qk_.capacity();
400 UTIL_CHECK(nk > 0);
401 UTIL_CHECK(expKsq_.capacity() == nk);
402
403 // Apply bond operator
404 fft().forwardTransform(q, qk_);
405 for (int i = 0; i < nk; ++i) {
406 qk_[i][0] *= expKsq2_[i];
407 qk_[i][1] *= expKsq2_[i];
408 }
409 fft().inverseTransformUnsafe(qk_, qout); // destroys qk_
410 }
411
412 /*
413 * Apply the local field operator for the bead model.
414 */
415 template <int D>
417 {
418 // Preconditions
419 int nx = mesh().size();
420 UTIL_CHECK(nx > 0);
421 UTIL_CHECK(expW_.capacity() == nx);
422 UTIL_CHECK(q.capacity() == nx);
423
424 // Apply field operator
425 for (int i = 0; i < nx; ++i) {
426 q[i] *= expW_[i];
427 }
428 }
429
430 /*
431 * Integrate to calculate monomer concentration for this block
432 */
433 template <int D>
435 {
436 // Preconditions
437 UTIL_CHECK(isAllocated_);
438 int nx = mesh().size();
439 UTIL_CHECK(nx > 0);
440 UTIL_CHECK(ns_ > 0);
441 UTIL_CHECK(propagator(0).isSolved());
442 UTIL_CHECK(propagator(1).isSolved());
443 UTIL_CHECK(cField().capacity() == nx);
444
445 // Initialize cField to zero at all points
446 int i;
447 for (i = 0; i < nx; ++i) {
448 cField()[i] = 0.0;
449 }
450
451 // References to forward and reverse propagators
452 Propagator<D> const & p0 = propagator(0);
453 Propagator<D> const & p1 = propagator(1);
454
455 // Evaluate unnormalized integral
456
457 // Endpoint contributions
458 for (i = 0; i < nx; ++i) {
459 cField()[i] += p0.q(0)[i]*p1.q(ns_ - 1)[i];
460 cField()[i] += p0.q(ns_ -1)[i]*p1.q(0)[i];
461 }
462
463 // Odd indices
464 int j;
465 for (j = 1; j < (ns_ -1); j += 2) {
466 RField<D> const & qf = p0.q(j);
467 RField<D> const & qr = p1.q(ns_ - 1 - j);
468 for (i = 0; i < nx; ++i) {
469 //cField()[i] += p0.q(j)[i] * p1.q(ns_ - 1 - j)[i] * 4.0;
470 cField()[i] += qf[i] * qr[i] * 4.0;
471 }
472 }
473
474 // Even indices
475 for (j = 2; j < (ns_ -2); j += 2) {
476 RField<D> const & qf = p0.q(j);
477 RField<D> const & qr = p1.q(ns_ - 1 - j);
478 for (i = 0; i < nx; ++i) {
479 // cField()[i] += p0.q(j)[i] * p1.q(ns_ - 1 - j)[i] * 2.0;
480 cField()[i] += qf[i] * qr[i] * 2.0;
481 }
482 }
483
484 // Normalize the integral
485 prefactor *= ds_ / 3.0;
486 for (i = 0; i < nx; ++i) {
487 cField()[i] *= prefactor;
488 }
489
490 }
491
492 /*
493 * Calculate monomer concentration for this block, bead model.
494 */
495 template <int D>
497 {
498 // Preconditions
499 UTIL_CHECK(isAllocated_);
500 int nx = mesh().size();
501 UTIL_CHECK(nx > 0);
502 UTIL_CHECK(ns_ > 0);
503 UTIL_CHECK(propagator(0).isSolved());
504 UTIL_CHECK(propagator(1).isSolved());
505 UTIL_CHECK(cField().capacity() == nx);
506
507 // Initialize cField to zero at all points
508 int i;
509 for (i = 0; i < nx; ++i) {
510 cField()[i] = 0.0;
511 }
512
513 // References to forward and reverse propagators
514 Propagator<D> const & p0 = propagator(0);
515 Propagator<D> const & p1 = propagator(1);
516
517 // Sum over beads (j = 1, ... , ns_ -2)
518 int j;
519 for (j = 1; j < (ns_ -1); ++j) {
520 RField<D> const & qf = p0.q(j);
521 RField<D> const & qr = p1.q(ns_ - 1 - j);
522 for (i = 0; i < nx; ++i) {
523 cField()[i] += qf[i] * qr[i] * expWInv_[i];
524 }
525 }
526
527 // Normalize the integral
528 for (i = 0; i < nx; ++i) {
529 cField()[i] *= prefactor;
530 }
531 }
532
533 /*
534 * Integrate to compute stress exerted by this block (thread).
535 */
536 template <int D>
537 void Block<D>::computeStressThread(double prefactor)
538 {
539 // Preconditions
541 int nx = mesh().size();
542 UTIL_CHECK(nx > 0);
543 UTIL_CHECK(fft().isSetup());
544 UTIL_CHECK(mesh().dimensions() == fft().meshDimensions());
545 UTIL_CHECK(ns_ > 0);
546 UTIL_CHECK(ds_ > 0);
547
548 // References to forward and reverse propagators
549 Propagator<D> const & p0 = propagator(0);
550 Propagator<D> const & p1 = propagator(1);
553
554 // If necessary, update derivatives of |k|^2 in WaveList
555 UTIL_CHECK(waveListPtr_);
556 if (!waveListPtr_->hasdKSq()) {
557 waveListPtr_->computedKSq();
558 }
559
560 // Initialize dQ work array to zero
562 const int nParam = unitCell().nParameter();
563 for (int i = 0; i < nParam; ++i) {
564 dQ.append(0.0);
565 }
566
567 // Work variables
568 const double bSq = kuhn()*kuhn()/6.0;
569 double dels, prod, increment;
570 int m, n;
571
572 // Evaluate unnormalized integral over contour
573 for (int j = 0; j < ns_ ; ++j) {
574
575 qr_ = p0.q(j);
576 fft().forwardTransform(qr_, qk_);
577
578 qr2_ = p1.q(ns_ - 1 - j);
579 fft().forwardTransform(qr2_, qk2_);
580
581 // Compute prefactor dels for Simpson's rule
582 dels = ds_ / 3.0;
583 if (j != 0 && j != ns_ - 1) {
584 if (j % 2 == 0) {
585 dels *= 2.0;
586 } else {
587 dels *= 4.0;
588 }
589 }
590
591 // Loop over unit cell parameters
592 for (n = 0; n < nParam ; ++n) {
593 RField<D> dKSq = waveListPtr_->dKSq(n);
594 increment = 0.0;
595 // Loop over wavevectors
596 for (m = 0; m < kSize_ ; ++m) {
597 prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
598 prod *= dKSq[m];
599 increment += prod;
600 }
601 increment *= bSq * dels;
602 dQ[n] = dQ[n] - increment;
603 }
604
605 }
606
607 // Compute stress_ from dQ
608 stress_.clear();
609 for (int n = 0; n < nParam; ++n) {
610 stress_.append(-1.0*prefactor*dQ[n]);
611 }
612
613 }
614
615 /*
616 * Compute contribution of this block to stress for bead model.
617 */
618 template <int D>
619 void Block<D>::computeStressBead(double prefactor)
620 {
621 // Preconditions
623 int nx = mesh().size();
624 UTIL_CHECK(nx > 0);
625 UTIL_CHECK(fft().isSetup());
626 UTIL_CHECK(mesh().dimensions() == fft().meshDimensions());
627 UTIL_CHECK(ns_ > 0);
628
629 // References to forward to reverse propagators
630 Propagator<D> const & p0 = propagator(0);
631 Propagator<D> const & p1 = propagator(1);
632 UTIL_CHECK(p0.isSolved());
633 UTIL_CHECK(p1.isSolved());
634
635
636 // If necessary, update derivatives of |k|^2 in WaveList
637 UTIL_CHECK(waveListPtr_);
638 if (!waveListPtr_->hasdKSq()) {
639 waveListPtr_->computedKSq();
640 }
641
642 // Initialize dQ work array to zero
644 int nParam = unitCell().nParameter();
645 for (int i = 0; i < nParam; ++i) {
646 dQ.append(0.0);
647 }
648
649 const double bSq = kuhn()*kuhn()/6.0;
650 double increment, prod;
651
652 // Half-bond from head junction to first bead, if not a chain end
653 if (!p0.isHeadEnd()) {
654
655 // Bead j = 0, forward propagator
656 qr_ = p0.q(0);
657 fft().forwardTransform(qr_, qk_);
658
659 // Next bead (ns_ - 2), reverse propagator
660 qr2_ = p1.q(ns_ - 2);
661 fft().forwardTransform(qr2_, qk2_);
662
663 // Loop over unit cell parameters
664 for (int n = 0; n < nParam ; ++n) {
665 RField<D> dKSq = waveListPtr_->dKSq(n);
666 increment = 0.0;
667 // Loop over wavevectors
668 for (int m = 0; m < kSize_ ; ++m) {
669 prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
670 prod *= dKSq[m]*expKsq2_[m];
671 increment += prod;
672 }
673 increment *= 0.5*bSq;
674 dQ[n] = dQ[n] - increment;
675 }
676 }
677
678 // Loop over nBead - 1 full bonds within this block
679 for (int j = 1; j < ns_ - 2 ; ++j) {
680
681 // Bead j, forward propagator
682 qr_ = p0.q(j);
683 fft().forwardTransform(qr_, qk_);
684
685 // Next bead (ns_- 2 - j), reverse propagator
686 qr2_ = p1.q(ns_ - 2 - j);
687 fft().forwardTransform(qr2_, qk2_);
688
689 // Loop over unit cell parameters
690 for (int n = 0; n < nParam ; ++n) {
691 RField<D> dKSq = waveListPtr_->dKSq(n);
692 increment = 0.0;
693 // Loop over wavevectors
694 for (int m = 0; m < kSize_ ; ++m) {
695 prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
696 prod *= dKSq[m]*expKsq_[m];
697 increment += prod;
698 }
699 increment *= bSq;
700 dQ[n] = dQ[n] - increment;
701 }
702
703 }
704
705 // Half-bond from last bead to tail junction, if not a chain end
706 if (!p0.isTailEnd()) {
707
708 // Second-to-last bead, j = ns_ - 2, forward propagator
709 qr_ = p0.q(ns_-2);
710 fft().forwardTransform(qr_, qk_);
711
712 // Last bead, j=0, reverse propagator
713 qr2_ = p1.q(0);
714 fft().forwardTransform(qr2_, qk2_);
715
716 // Loop over unit cell parameters
717 for (int n = 0; n < nParam ; ++n) {
718 RField<D> dKSq = waveListPtr_->dKSq(n);
719 increment = 0.0;
720 // Loop over wavevectors
721 for (int m = 0; m < kSize_ ; ++m) {
722 prod = (qk2_[m][0] * qk_[m][0]) + (qk2_[m][1] * qk_[m][1]);
723 prod *= dKSq[m]*expKsq2_[m];
724 increment += prod;
725 }
726 increment *= 0.5*bSq;
727 dQ[n] = dQ[n] - increment;
728 }
729
730 }
731
732 // Compute stress_ from dQ
733 stress_.clear();
734 for (int i = 0; i < nParam; ++i) {
735 stress_.append(-1.0*prefactor*dQ[i]);
736 }
737
738 }
739
740}
741}
742#endif
virtual void setLength(double length)
Set the length of this block (only valid for thread model).
Definition Edge.cpp:68
Iterator over points in a Mesh<D>.
int rank() const
Get the rank of current element.
void begin()
Set iterator to the first point in the mesh.
bool atEnd() const
Is this the end (i.e., one past the last point)?
void setDimensions(const IntVec< D > &dimensions)
Set the grid dimensions in all directions.
Description of a regular grid of points in a periodic domain.
Definition Mesh.h:61
Fourier transform wrapper.
Definition cpu/FFT.h:38
static void computeKMesh(IntVec< D > const &rMeshDimensions, IntVec< D > &kMeshDimensions, int &kSize)
Compute dimensions and size of k-space mesh for DFT of real data.
Definition cpu/FFT.tpp:262
bool isAllocated() const
Return true if the FftwDArray has been allocated, false otherwise.
Definition FftwDArray.h:102
Field of real double precision values on an FFT mesh.
Definition cpu/RField.h:29
Class to calculate and store properties of wavevectors.
Base template for UnitCell<D> classes, D=1, 2 or 3.
Definition UnitCell.h:56
void stepBondBead(RField< D > const &qin, RField< D > &qout) const
Apply a bond operator for the bead model.
void setLength(double newLength)
Set or reset block length (only used in thread model).
void computeConcentrationBead(double prefactor)
Compute the concentration for this block, using the bead model.
void setKuhn(double kuhn)
Set or reset monomer statistical segment length.
void clearUnitCellData()
Clear all internal data that depends on the unit cell parameters.
void setupSolver(RField< D > const &w)
Set up the MDE solver for this block.
double kuhn() const
Get monomer statistical segment length.
void allocate(double ds)
Allocate memory and set contour step size.
void computeStressBead(double prefactor)
Compute stress contribution for this block, using bead model.
int nBead() const
Get the number of beads in this block, in the bead model.
Definition Edge.h:302
void associate(Mesh< D > const &mesh, FFT< D > const &fft, UnitCell< D > const &cell, WaveList< D > &wavelist)
Create permanent associations with related objects.
Propagator< D > & propagator(int directionId)
Get a Propagator for a specified direction.
void computeConcentrationThread(double prefactor)
Compute the concentration for this block, for the thread model.
void stepThread(RField< D > const &qin, RField< D > &qout) const
Compute one step of solution of MDE for the thread model.
void computeStressThread(double prefactor)
Compute stress contribution for this block, using thread model.
double ds() const
Get contour length step size.
void stepFieldBead(RField< D > &q) const
Apply the exponential field operator for the bead model.
void stepHalfBondBead(RField< D > const &qin, RField< D > &qout) const
Apply a half-bond operator for the bead model.
RField< D > & cField()
Get the associated monomer concentration field.
void stepBead(RField< D > const &qin, RField< D > &qout) const
Compute one step of solution of MDE for the bead model.
double length() const
Get the length of this block, in the thread model.
Definition Edge.h:311
MDE solver for one direction of one block.
bool isHeadEnd() const
Is the head vertex a chain end?
bool isSolved() const
Has the modified diffusion equation been solved?
bool isTailEnd() const
Is the tail vertex a chain end?
bool isAllocated() const
Has memory been allocated for this propagator?
const FieldT & q(int i) const
Return q-field at specified step.
int capacity() const
Return allocated size.
Definition Array.h:159
A fixed capacity (static) contiguous array with a variable logical size.
Definition FSArray.h:38
void append(Data const &data)
Append data to the end of the array.
Definition FSArray.h:258
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
bool isThread()
Is the thread model in use ?
bool isBead()
Is the bead model in use ?
Fields and FFTs for periodic boundary conditions (CPU)
Definition CField.cpp:12
Periodic fields and crystallography.
Definition CField.cpp:11
Real periodic fields, SCFT and PS-FTS (CPU).
Definition param_pc.dox:2
PSCF package top-level namespace.