PSCF v1.2
rpg/System.tpp
1#ifndef RPG_SYSTEM_TPP
2#define RPG_SYSTEM_TPP
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 "System.h"
12
13#include <rpg/fts/simulator/Simulator.h>
14#include <rpg/fts/simulator/SimulatorFactory.h>
15#include <rpg/fts/compressor/Compressor.h>
16#include <rpg/scft/sweep/Sweep.h>
17#include <rpg/scft/sweep/SweepFactory.h>
18#include <rpg/scft/iterator/Iterator.h>
19#include <rpg/scft/iterator/IteratorFactory.h>
20#include <rpg/solvers/Polymer.h>
21#include <rpg/solvers/Solvent.h>
22
23#include <prdc/cuda/resources.h>
24#include <prdc/cuda/RField.h>
25#include <prdc/cuda/RFieldComparison.h>
26#include <prdc/crystal/BFieldComparison.h>
27
28#include <pscf/inter/Interaction.h>
29#include <pscf/math/IntVec.h>
30#include <pscf/homogeneous/Clump.h>
31
32#include <util/containers/FSArray.h>
33#include <util/param/BracketPolicy.h>
34#include <util/param/ParamComponent.h>
35#include <util/format/Str.h>
36#include <util/format/Int.h>
37#include <util/format/Dbl.h>
38#include <util/misc/ioUtil.h>
39
40#include <string>
41#include <unistd.h>
42
43namespace Pscf {
44namespace Rpg {
45
46 using namespace Util;
47 using namespace Pscf::Prdc;
48 using namespace Pscf::Prdc::Cuda;
49
50 /*
51 * Constructor.
52 */
53 template <int D>
55 : mixture_(),
56 domain_(),
57 fileMaster_(),
58 homogeneous_(),
59 interactionPtr_(nullptr),
60 iteratorPtr_(nullptr),
61 iteratorFactoryPtr_(nullptr),
62 sweepPtr_(nullptr),
63 sweepFactoryPtr_(nullptr),
64 simulatorPtr_(nullptr),
65 simulatorFactoryPtr_(nullptr),
66 w_(),
67 c_(),
68 h_(),
69 mask_(),
70 fHelmholtz_(0.0),
71 fIdeal_(0.0),
72 fInter_(0.0),
73 fExt_(0.0),
74 pressure_(0.0),
75 isAllocatedGrid_(false),
76 isAllocatedBasis_(false),
77 hasMixture_(false),
78 hasCFields_(false),
79 hasFreeEnergy_(false)
80 {
81 setClassName("System");
82 domain_.setFileMaster(fileMaster_);
83 w_.setFieldIo(domain_.fieldIo());
84 h_.setFieldIo(domain_.fieldIo());
85 mask_.setFieldIo(domain_.fieldIo());
86
87 interactionPtr_ = new Interaction();
88 iteratorFactoryPtr_ = new IteratorFactory<D>(*this);
89 sweepFactoryPtr_ = new SweepFactory<D>(*this);
90 simulatorFactoryPtr_ = new SimulatorFactory<D>(*this);
91 BracketPolicy::set(BracketPolicy::Optional);
93 }
94
95 /*
96 * Destructor.
97 */
98 template <int D>
100 {
101 if (interactionPtr_) {
102 delete interactionPtr_;
103 }
104 if (iteratorPtr_) {
105 delete iteratorPtr_;
106 }
107 if (iteratorFactoryPtr_) {
108 delete iteratorFactoryPtr_;
109 }
110 if (sweepPtr_) {
111 delete sweepPtr_;
112 }
113 if (sweepFactoryPtr_) {
114 delete sweepFactoryPtr_;
115 }
116 if (simulatorPtr_) {
117 delete simulatorPtr_;
118 }
119 if (simulatorFactoryPtr_) {
120 delete simulatorFactoryPtr_;
121 }
122 }
123
124 /*
125 * Process command line options.
126 */
127 template <int D>
128 void System<D>::setOptions(int argc, char **argv)
129 {
130 bool eFlag = false; // echo
131 bool dFlag = false; // spatial dimension (1, 2, or 3)
132 bool pFlag = false; // param file
133 bool cFlag = false; // command file
134 bool iFlag = false; // input prefix
135 bool oFlag = false; // output prefix
136 bool tFlag = false; // GPU threads (max # threads per block)
137 char* pArg = 0;
138 char* cArg = 0;
139 char* iArg = 0;
140 char* oArg = 0;
141 int dArg = 0;
142 int tArg = 0;
143
144 // Read program arguments
145 int c;
146 opterr = 0;
147 while ((c = getopt(argc, argv, "ed:p:c:i:o:t:")) != -1) {
148 switch (c) {
149 case 'e':
150 eFlag = true;
151 break;
152 case 'd':
153 dFlag = true;
154 dArg = atoi(optarg);
155 break;
156 case 'p': // parameter file
157 pFlag = true;
158 pArg = optarg;
159 break;
160 case 'c': // command file
161 cFlag = true;
162 cArg = optarg;
163 break;
164 case 'i': // input prefix
165 iFlag = true;
166 iArg = optarg;
167 break;
168 case 'o': // output prefix
169 oFlag = true;
170 oArg = optarg;
171 break;
172 case 't': //number of threads per block, user set
173 tFlag = true;
174 tArg = atoi(optarg);
175 break;
176 case '?':
177 Log::file() << "Unknown option -" << optopt << std::endl;
178 UTIL_THROW("Invalid command line option");
179 }
180 }
181
182 // Set flag to echo parameters as they are read.
183 if (eFlag) {
185 }
186
187 // Check -d flag
188 if (dFlag) {
189 UTIL_CHECK(D == dArg);
190 } else {
191 UTIL_THROW("Missing required -d option");
192 }
193
194 // Check option -p, set parameter file name
195 if (pFlag) {
196 fileMaster_.setParamFileName(std::string(pArg));
197 } else {
198 UTIL_THROW("Missing required -p option - no parameter file");
199 }
200
201 // Check option -c, set command file name
202 if (cFlag) {
203 fileMaster_.setCommandFileName(std::string(cArg));
204 } else {
205 UTIL_THROW("Missing required -c option - no command file");
206 }
207
208 // If option -i, set path prefix for input files
209 if (iFlag) {
210 fileMaster_.setInputPrefix(std::string(iArg));
211 }
212
213 // If option -o, set path prefix for output files
214 if (oFlag) {
215 fileMaster_.setOutputPrefix(std::string(oArg));
216 }
217
218 // If option -t, set the threads per block.
219 if (tFlag) {
222 }
223
224 }
225
226 /*
227 * Read parameters and initialize.
228 */
229 template <int D>
230 void System<D>::readParameters(std::istream& in)
231 {
232 // Read the Mixture{ ... } block
233 readParamComposite(in, mixture_);
234 hasMixture_ = true;
235
236 int nm = mixture_.nMonomer();
237 int np = mixture_.nPolymer();
238 int ns = mixture_.nSolvent();
239 UTIL_CHECK(nm > 0);
240 UTIL_CHECK(np >= 0);
241 UTIL_CHECK(ns >= 0);
242 UTIL_CHECK(np + ns > 0);
243
244 // Read the Interaction{ ... } block
245 interaction().setNMonomer(nm);
246 readParamComposite(in, interaction());
247
248 // Read the Domain{ ... } block
249 readParamComposite(in, domain_);
250 UTIL_CHECK(domain_.mesh().size() > 0);
251 UTIL_CHECK(domain_.unitCell().nParameter() > 0);
252 UTIL_CHECK(domain_.unitCell().lattice() != UnitCell<D>::Null);
253
254 // Setup the mixture
255 mixture_.associate(domain_.mesh(), fft(),
256 domain_.unitCell(), domain_.waveList());
257 mixture_.allocate();
258 // mixture_.clearUnitCellData();
259
260 // Allocate memory for r-grid w and c fields
261 allocateFieldsGrid();
262
263 // Optionally instantiate an Iterator object
264 std::string className;
265 bool isEnd;
266 iteratorPtr_ =
267 iteratorFactoryPtr_->readObjectOptional(in, *this, className,
268 isEnd);
269 if (!iteratorPtr_ && ParamComponent::echo()) {
270 Log::file() << indent() << " Iterator{ [absent] }\n";
271 }
272
273 // Optionally instantiate a Sweep object (if an iterator exists)
274 if (iteratorPtr_ && !isEnd) {
275 sweepPtr_ =
276 sweepFactoryPtr_->readObjectOptional(in, *this, className,
277 isEnd);
278 if (!sweepPtr_ && ParamComponent::echo()) {
279 Log::file() << indent() << " Sweep{ [absent] }\n";
280 }
281 }
282
283 // Optionally instantiate a Simulator object
284 if (!isEnd) {
285 simulatorPtr_ =
286 simulatorFactoryPtr_->readObjectOptional(in, *this,
287 className, isEnd);
288 if (!simulatorPtr_ && ParamComponent::echo()) {
289 Log::file() << indent() << " Simulator{ [absent] }\n";
290 }
291 }
292
293 // Initialize homogeneous object
294 // NOTE: THIS OBJECT IS NOT USED AT ALL.
295 homogeneous_.setNMolecule(np+ns);
296 homogeneous_.setNMonomer(nm);
297 initHomogeneous();
298
299 }
300
301 /*
302 * Read parameter file (including open and closing brackets).
303 */
304 template <int D>
305 void System<D>::readParam(std::istream& in)
306 {
307 readBegin(in, className().c_str());
308 readParameters(in);
309 readEnd(in);
310 }
311
312 /*
313 * Read default parameter file.
314 */
315 template <int D>
317 { readParam(fileMaster_.paramFile()); }
318
319 /*
320 * Read and execute commands from a specified command file.
321 */
322 template <int D>
323 void System<D>::readCommands(std::istream &in)
324 {
325 UTIL_CHECK(isAllocatedGrid_);
326 std::string command, filename, inFileName, outFileName;
327
328 bool readNext = true;
329 while (readNext) {
330
331 in >> command;
332
333 if (in.eof()) {
334 break;
335 } else {
336 Log::file() << command << std::endl;
337 }
338
339 if (command == "FINISH") {
340 Log::file() << std::endl;
341 readNext = false;
342 } else
343 if (command == "READ_W_BASIS") {
344 readEcho(in, filename);
345 readWBasis(filename);
346 } else
347 if (command == "READ_W_RGRID") {
348 readEcho(in, filename);
349 readWRGrid(filename);
350 } else
351 if (command == "ESTIMATE_W_FROM_C") {
352 readEcho(in, inFileName);
353 estimateWfromC(inFileName);
354 } else
355 if (command == "SET_UNIT_CELL") {
356 UnitCell<D> unitCell;
357 in >> unitCell;
358 Log::file() << " " << unitCell << std::endl;
359 setUnitCell(unitCell);
360 } else
361 if (command == "COMPUTE") {
362 // Solve the modified diffusion equation, without iteration
363 compute();
364 } else
365 if (command == "ITERATE") {
366 // Attempt to iteratively solve a single SCFT problem
367 bool isContinuation = false;
368 int fail = iterate(isContinuation);
369 if (fail) {
370 readNext = false;
371 }
372 } else
373 if (command == "SWEEP") {
374 // Attempt to solve a sequence of SCFT problems along a path
375 // through parameter space
376 sweep();
377 } else
378 if (command == "COMPRESS") {
379 // Impose incompressibility
380 UTIL_CHECK(hasSimulator());
381 simulator().compressor().compress();
382 } else
383 if (command == "SIMULATE") {
384 // Perform a field theoretic simulation
385 int nStep;
386 in >> nStep;
387 Log::file() << " " << nStep << "\n";
388 simulate(nStep);
389 } else
390 if (command == "ANALYZE" || command == "ANALYZE_TRAJECTORY") {
391 // Read and analyze a field trajectory file
392 int min, max;
393 in >> min;
394 in >> max;
395 Log::file() << " " << min ;
396 Log::file() << " " << max << "\n";
397 std::string classname;
398 readEcho(in, classname);
399 readEcho(in, filename);
400 simulator().analyze(min, max, classname, filename);
401 } else
402 if (command == "WRITE_TIMERS") {
403 readEcho(in, filename);
404 std::ofstream file;
405 fileMaster_.openOutputFile(filename, file);
406 writeTimers(file);
407 file.close();
408 } else
409 if (command == "CLEAR_TIMERS") {
410 clearTimers();
411 } else
412 if (command == "WRITE_PARAM") {
413 readEcho(in, filename);
414 std::ofstream file;
415 fileMaster_.openOutputFile(filename, file);
416 writeParamNoSweep(file);
417 file.close();
418 } else
419 if (command == "WRITE_THERMO") {
420 readEcho(in, filename);
421 std::ofstream file;
422 fileMaster_.openOutputFile(filename, file,
423 std::ios_base::app);
424 writeThermo(file);
425 file.close();
426 } else
427 if (command == "WRITE_STRESS") {
428 readEcho(in, filename);
429 std::ofstream file;
430 fileMaster_.openOutputFile(filename, file,
431 std::ios_base::app);
432 writeStress(file);
433 file.close();
434 } else
435 if (command == "WRITE_W_BASIS") {
436 readEcho(in, filename);
437 writeWBasis(filename);
438 } else
439 if (command == "WRITE_W_RGRID") {
440 readEcho(in, filename);
441 writeWRGrid(filename);
442 } else
443 if (command == "WRITE_C_BASIS") {
444 readEcho(in, filename);
445 writeCBasis(filename);
446 } else
447 if (command == "WRITE_C_RGRID") {
448 readEcho(in, filename);
449 writeCRGrid(filename);
450 } else
451 if (command == "WRITE_BLOCK_C_RGRID") {
452 readEcho(in, filename);
453 writeBlockCRGrid(filename);
454 } else
455 if (command == "WRITE_Q_SLICE") {
456 int polymerId, blockId, directionId, segmentId;
457 readEcho(in, filename);
458 in >> polymerId;
459 in >> blockId;
460 in >> directionId;
461 in >> segmentId;
462 Log::file() << Str("polymer ID ", 21) << polymerId << "\n"
463 << Str("block ID ", 21) << blockId << "\n"
464 << Str("direction ID ", 21) << directionId
465 << "\n"
466 << Str("segment ID ", 21) << segmentId
467 << std::endl;
468 writeQSlice(filename, polymerId, blockId, directionId,
469 segmentId);
470 } else
471 if (command == "WRITE_Q_TAIL") {
472 int polymerId, blockId, directionId;
473 readEcho(in, filename);
474 in >> polymerId;
475 in >> blockId;
476 in >> directionId;
477 Log::file() << Str("polymer ID ", 21) << polymerId << "\n"
478 << Str("block ID ", 21) << blockId << "\n"
479 << Str("direction ID ", 21) << directionId
480 << "\n";
481 writeQTail(filename, polymerId, blockId, directionId);
482 } else
483 if (command == "WRITE_Q") {
484 int polymerId, blockId, directionId;
485 readEcho(in, filename);
486 in >> polymerId;
487 in >> blockId;
488 in >> directionId;
489 Log::file() << Str("polymer ID ", 21) << polymerId << "\n"
490 << Str("block ID ", 21) << blockId << "\n"
491 << Str("direction ID ", 21) << directionId << "\n";
492 writeQ(filename, polymerId, blockId, directionId);
493 } else
494 if (command == "WRITE_Q_ALL") {
495 readEcho(in, filename);
496 writeQAll(filename);
497 } else
498 if (command == "WRITE_STARS") {
499 readEcho(in, filename);
500 writeStars(filename);
501 } else
502 if (command == "WRITE_WAVES") {
503 readEcho(in, filename);
504 writeWaves(filename);
505 } else
506 if (command == "WRITE_GROUP") {
507 readEcho(in, filename);
508 writeGroup(filename);
509 //writeGroup(filename, domain_.group());
510 } else
511 if (command == "BASIS_TO_RGRID") {
512 readEcho(in, inFileName);
513 readEcho(in, outFileName);
514 basisToRGrid(inFileName, outFileName);
515 } else
516 if (command == "RGRID_TO_BASIS") {
517 readEcho(in, inFileName);
518 readEcho(in, outFileName);
519 rGridToBasis(inFileName, outFileName);
520 } else
521 if (command == "KGRID_TO_RGRID") {
522 readEcho(in, inFileName);
523 readEcho(in, outFileName);
524 kGridToRGrid(inFileName, outFileName);
525 } else
526 if (command == "RGRID_TO_KGRID") {
527 readEcho(in, inFileName);
528 readEcho(in, outFileName);
529 rGridToKGrid(inFileName, outFileName);
530 } else
531 if (command == "KGRID_TO_BASIS") {
532 readEcho(in, inFileName);
533 readEcho(in, outFileName);
534 kGridToBasis(inFileName, outFileName);
535 } else
536 if (command == "CHECK_RGRID_SYMMETRY") {
537 double epsilon;
538 readEcho(in, inFileName);
539 readEcho(in, epsilon);
540 bool hasSymmetry;
541 hasSymmetry = checkRGridFieldSymmetry(inFileName, epsilon);
542 if (hasSymmetry) {
543 Log::file() << std::endl
544 << "Symmetry of r-grid file matches this space group."
545 << std::endl << std::endl;
546 } else {
547 Log::file() << std::endl
548 << "Symmetry of r-grid file does not match this\n"
549 << "space group to within error threshold of "
550 << Dbl(epsilon) << " ." << std::endl << std::endl;
551 }
552 } else
553 if (command == "COMPARE_BASIS") {
554
555 // Get two filenames for comparison
556 std::string filecompare1, filecompare2;
557 readEcho(in, filecompare1);
558 readEcho(in, filecompare2);
559
560 DArray< DArray<double> > Bfield1, Bfield2;
561 UnitCell<D> tmpUnitCell;
562 domain().fieldIo().readFieldsBasis(filecompare1, Bfield1,
563 tmpUnitCell);
564 domain().fieldIo().readFieldsBasis(filecompare2, Bfield2,
565 tmpUnitCell);
566 // Note: Bfield1 & Bfield2 are allocated by readFieldsBasis
567
568 // Compare and output report
569 compare(Bfield1, Bfield2);
570
571 } else
572 if (command == "COMPARE_RGRID") {
573 // Get two filenames for comparison
574 std::string filecompare1, filecompare2;
575 readEcho(in, filecompare1);
576 readEcho(in, filecompare2);
577
578 DArray< RField<D> > Rfield1, Rfield2;
579 UnitCell<D> tmpUnitCell;
580 domain().fieldIo().readFieldsRGrid(filecompare1, Rfield1,
581 tmpUnitCell);
582 domain().fieldIo().readFieldsRGrid(filecompare2, Rfield2,
583 tmpUnitCell);
584 // Note: Rfield1, Rfield2 will be allocated by readFieldsRGrid
585
586 // Compare and output report
587 compare(Rfield1, Rfield2);
588
589 } else
590 if (command == "SCALE_BASIS") {
591 double factor;
592 readEcho(in, inFileName);
593 readEcho(in, outFileName);
594 readEcho(in, factor);
595 scaleFieldsBasis(inFileName, outFileName, factor);
596 } else
597 if (command == "SCALE_RGRID") {
598 double factor;
599 readEcho(in, inFileName);
600 readEcho(in, outFileName);
601 readEcho(in, factor);
602 scaleFieldsRGrid(inFileName, outFileName, factor);
603 } else
604 if (command == "EXPAND_RGRID_DIMENSION") {
605 readEcho(in, inFileName);
606 readEcho(in, outFileName);
607
608 // Read expanded spatial dimension d
609 int d;
610 in >> d;
611 UTIL_CHECK(d > D);
612 Log::file() << Str("Expand dimension to: ") << d << "\n";
613
614 // Read numbers of grid points along added dimensions
615 DArray<int> newGridDimensions;
616 newGridDimensions.allocate(d-D);
617 for (int i = 0; i < d-D; i++) {
618 in >> newGridDimensions[i];
619 }
620 Log::file()
621 << Str("Numbers of grid points in added dimensions : ");
622 for (int i = 0; i < d-D; i++) {
623 Log::file() << newGridDimensions[i] << " ";
624 }
625 Log::file() << "\n";
626
627 expandRGridDimension(inFileName, outFileName,
628 d, newGridDimensions);
629
630 } else
631 if (command == "REPLICATE_UNIT_CELL") {
632 readEcho(in, inFileName);
633 readEcho(in, outFileName);
634
635 // Read numbers of replicas along each direction
636 IntVec<D> replicas;
637 for (int i = 0; i < D; i++){
638 in >> replicas[i];
639 }
640 for (int i = 0; i < D; i++){
641 Log::file() << "Replicate unit cell in direction "
642 << i << " : ";
643 Log::file() << replicas[i] << " times ";
644 Log::file() << "\n";
645 }
646
647 replicateUnitCell(inFileName, outFileName, replicas);
648 } else
649 if (command == "BASIS_TO_KGRID") {
650 readEcho(in, inFileName);
651 readEcho(in, outFileName);
652 basisToKGrid(inFileName, outFileName);
653 } else
654 if (command == "READ_H_BASIS") {
655 readEcho(in, filename);
656 if (!h_.isAllocatedBasis()) {
657 h_.allocateBasis(domain().basis().nBasis());
658 }
659 if (!h_.isAllocatedRGrid()) {
660 h_.allocateRGrid(domain().mesh().dimensions());
661 }
662 UnitCell<D> tmpUnitCell;
663 h_.readBasis(filename, tmpUnitCell);
664 } else
665 if (command == "READ_H_RGRID") {
666 readEcho(in, filename);
667 if (!h_.isAllocatedRGrid()) {
668 h_.allocateRGrid(domain().mesh().dimensions());
669 }
670 if (iterator().isSymmetric() && !h_.isAllocatedBasis()) {
671 mask_.allocateBasis(domain().basis().nBasis());
672 }
673 UnitCell<D> tmpUnitCell;
674 h_.readRGrid(filename, tmpUnitCell);
675 } else
676 if (command == "WRITE_H_BASIS") {
677 readEcho(in, filename);
678 UTIL_CHECK(h_.hasData());
679 UTIL_CHECK(h_.isSymmetric());
680 domain().fieldIo().writeFieldsBasis(filename, h_.basis(),
681 domain().unitCell());
682 } else
683 if (command == "WRITE_H_RGRID") {
684 readEcho(in, filename);
685 UTIL_CHECK(h_.hasData());
686 domain().fieldIo().writeFieldsRGrid(filename, h_.rgrid(),
687 domain().unitCell());
688 } else
689 if (command == "READ_MASK_BASIS") {
690 readEcho(in, filename);
691 UTIL_CHECK(domain().basis().isInitialized());
692 if (!mask_.isAllocatedBasis()) {
693 mask_.allocateBasis(domain().basis().nBasis());
694 }
695 if (!mask_.isAllocatedRGrid()) {
696 mask_.allocateRGrid(domain().mesh().dimensions());
697 }
698 UnitCell<D> tmpUnitCell;
699 mask_.readBasis(filename, tmpUnitCell);
700 } else
701 if (command == "READ_MASK_RGRID") {
702 readEcho(in, filename);
703 if (!mask_.isAllocatedRGrid()) {
704 mask_.allocateRGrid(domain().mesh().dimensions());
705 }
706 if (iterator().isSymmetric() && !mask_.isAllocatedBasis()) {
707 UTIL_CHECK(domain().basis().isInitialized());
708 mask_.allocateBasis(domain().basis().nBasis());
709 }
710 UnitCell<D> tmpUnitCell;
711 mask_.readRGrid(filename, tmpUnitCell);
712 } else
713 if (command == "WRITE_MASK_BASIS") {
714 readEcho(in, filename);
715 UTIL_CHECK(mask_.hasData());
716 UTIL_CHECK(mask_.isSymmetric());
717 domain().fieldIo().writeFieldBasis(filename, mask_.basis(),
718 domain().unitCell());
719 } else
720 if (command == "WRITE_MASK_RGRID") {
721 readEcho(in, filename);
722 UTIL_CHECK(mask_.hasData());
723 domain().fieldIo().writeFieldRGrid(filename, mask_.rgrid(),
724 domain().unitCell(),
725 mask_.isSymmetric());
726 } else {
727 Log::file() << "Error: Unknown command "
728 << command << std::endl;
729 readNext = false;
730 }
731 }
732 }
733
734 /*
735 * Read and execute commands from the default command file.
736 */
737 template <int D>
739 {
740 if (fileMaster_.commandFileName().empty()) {
741 UTIL_THROW("Empty command file name");
742 }
743 readCommands(fileMaster_.commandFile());
744 }
745
746 // W-Field Modifier Functions
747
748 /*
749 * Read w-field in symmetry adapted basis format.
750 */
751 template <int D>
752 void System<D>::readWBasis(const std::string & filename)
753 {
754 // Precondition
755 UTIL_CHECK(domain_.hasGroup());
756
757 if (!domain_.unitCell().isInitialized()) {
758 readFieldHeader(filename);
759 }
760 UTIL_CHECK(domain_.unitCell().isInitialized());
761 UTIL_CHECK(domain_.basis().isInitialized());
762 UTIL_CHECK(domain_.basis().nBasis() > 0);
763 if (!isAllocatedBasis_) {
764 allocateFieldsBasis();
765 }
766
767 // Read w fields
768 w_.readBasis(filename, domain_.unitCell());
769 hasCFields_ = false;
770 hasFreeEnergy_ = false;
771
772 // Clear unit cell data from waveList and mixture
773 domain_.waveList().clearUnitCellData();
774 mixture_.clearUnitCellData();
775 }
776
777 template <int D>
778 void System<D>::readWRGrid(const std::string & filename)
779 {
780 UTIL_CHECK(isAllocatedGrid_);
781
782 // If necessary, peek at header to initialize unit cell
783 if (!domain_.unitCell().isInitialized()) {
784 readFieldHeader(filename);
785 }
786 UTIL_CHECK(domain_.unitCell().isInitialized());
787 if (domain_.hasGroup() && !isAllocatedBasis_) {
788 UTIL_CHECK(domain_.basis().isInitialized());
789 UTIL_CHECK(domain_.basis().nBasis() > 0);
790 allocateFieldsBasis();
791 }
792
793 // Read w fields
794 w_.readRGrid(filename, domain_.unitCell());
795 hasCFields_ = false;
796 hasFreeEnergy_ = false;
797
798 // Clear unit cell data from waveList and mixture
799 domain_.waveList().clearUnitCellData();
800 mixture_.clearUnitCellData();
801 }
802
803 /*
804 * Construct estimate for w-fields from c-fields.
805 *
806 * Modifies wFields and wFieldsRGrid.
807 */
808 template <int D>
809 void System<D>::estimateWfromC(std::string const & filename)
810 {
811 // Preconditions
812 UTIL_CHECK(hasMixture_);
813 const int nm = mixture_.nMonomer();
814 UTIL_CHECK(nm > 0);
815 UTIL_CHECK(isAllocatedGrid_);
816 UTIL_CHECK(domain_.hasGroup());
817
818 // Check required initializations
819 if (!domain_.unitCell().isInitialized()) {
820 readFieldHeader(filename);
821 }
822 UTIL_CHECK(domain_.unitCell().isInitialized());
823 UTIL_CHECK(domain_.basis().isInitialized());
824 if (!isAllocatedBasis_) {
825 allocateFieldsBasis();
826 }
827
828 // Allocate temporary storage
829 const int nb = domain_.basis().nBasis();
830 UTIL_CHECK(nb > 0);
831 DArray< DArray<double> > tmpCFieldsBasis;
832 tmpCFieldsBasis.allocate(nm);
833 for (int i = 0; i < nm; ++i) {
834 tmpCFieldsBasis[i].allocate(nb);
835 }
836
837 // Read c fields from input file
838 fieldIo().readFieldsBasis(filename, tmpCFieldsBasis,
839 domain_.unitCell());
840
841 // Compute w fields from c fields
842 for (int i = 0; i < nb; ++i) {
843 for (int j = 0; j < nm; ++j) {
844 tmpFieldsBasis_[j][i] = 0.0;
845 for (int k = 0; k < nm; ++k) {
846 tmpFieldsBasis_[j][i] += interaction().chi(j,k)
847 * tmpCFieldsBasis[k][i];
848 }
849 }
850 }
851
852 // Store estimated w fields in System w container
853 w_.setBasis(tmpFieldsBasis_);
854 hasCFields_ = false;
855 hasFreeEnergy_ = false;
856
857 // Clear unit cell data from waveList and mixture
858 domain_.waveList().clearUnitCellData();
859 mixture_.clearUnitCellData();
860 }
861
862 /*
863 * Set new w-field values.
864 */
865 template <int D>
867 {
868 UTIL_CHECK(domain_.hasGroup());
869 UTIL_CHECK(domain_.basis().isInitialized());
870 UTIL_CHECK(isAllocatedGrid_);
871 UTIL_CHECK(isAllocatedBasis_);
872 w_.setBasis(fields);
873 hasCFields_ = false;
874 hasFreeEnergy_ = false;
875 }
876
877 /*
878 * Set new w-field values, using r-grid fields as inputs.
879 */
880 template <int D>
881 void System<D>::setWRGrid(DArray< RField<D> > const & fields)
882 {
883 UTIL_CHECK(isAllocatedGrid_);
884 w_.setRGrid(fields);
885 hasCFields_ = false;
886 hasFreeEnergy_ = false;
887 }
888
889 /*
890 * Set new w-field values, using unfoldeded array of r-grid fields.
891 */
892 template <int D>
894 {
895 UTIL_CHECK(isAllocatedGrid_);
896 w_.setRGrid(fields);
897 hasCFields_ = false;
898 hasFreeEnergy_ = false;
899 }
900
901 /*
902 * Set new w-field values, using array of r-grid fields as input.
903 */
904 template <int D>
906 { w_.symmetrize(); }
907
908 // Unit Cell Modifiers
909
910 /*
911 * Set the associated unit cell.
912 */
913 template <int D>
914 void System<D>::setUnitCell(UnitCell<D> const & unitCell)
915 {
916 domain_.setUnitCell(unitCell);
917 // Note - Domain::setUnitCell updates the WaveList
918 mixture_.clearUnitCellData();
919 if (domain_.hasGroup() && !isAllocatedBasis_) {
920 UTIL_CHECK(domain_.basis().isInitialized());
921 allocateFieldsBasis();
922 }
923 }
924
925 /*
926 * Set the associated unit cell.
927 */
928 template <int D>
929 void
931 FSArray<double, 6> const & parameters)
932 {
933 domain_.setUnitCell(lattice, parameters);
934 // Note - Domain::setUnitCell updates the WaveList
935 mixture_.clearUnitCellData();
936 if (domain_.hasGroup() && !isAllocatedBasis_) {
937 UTIL_CHECK(domain_.basis().isInitialized());
938 allocateFieldsBasis();
939 }
940 }
941
942 /*
943 * Set parameters of the associated unit cell.
944 */
945 template <int D>
947 {
948 domain_.setUnitCell(parameters);
949 // Note - Domain::setUnitCell updates the WaveList
950 mixture_.clearUnitCellData();
951 if (domain_.hasGroup() && !isAllocatedBasis_) {
952 UTIL_CHECK(domain_.basis().isInitialized());
953 allocateFieldsBasis();
954 }
955 }
956
957 // Primary Field Theory Computations
958
959 /*
960 * Solve MDE for current w-fields, without iteration.
961 */
962 template <int D>
963 void System<D>::compute(bool needStress)
964 {
965 UTIL_CHECK(w_.isAllocatedRGrid());
966 UTIL_CHECK(c_.isAllocatedRGrid());
967 UTIL_CHECK(w_.hasData());
968
969 // Solve the modified diffusion equation (without iteration)
970 mixture_.compute(w_.rgrid(), c_.rgrid(), mask_.phiTot());
971 hasCFields_ = true;
972 hasFreeEnergy_ = false;
973
974 // Convert c fields from r-grid to basis format
975 if (w_.isSymmetric()) {
976 UTIL_CHECK(c_.isAllocatedBasis());
977 fieldIo().convertRGridToBasis(c_.rgrid(), c_.basis());
978 }
979
980 // Compute stress if needed
981 if (needStress) {
982 mixture_.computeStress(mask().phiTot());
983 }
984 }
985
986 /*
987 * Iteratively solve a SCFT problem for specified parameters.
988 */
989 template <int D>
990 int System<D>::iterate(bool isContinuation)
991 {
992 UTIL_CHECK(hasIterator());
993 UTIL_CHECK(w_.hasData());
994 if (iterator().isSymmetric()) {
995 UTIL_CHECK(w_.isSymmetric());
996 UTIL_CHECK(isAllocatedBasis_);
997 }
998 hasCFields_ = false;
999 hasFreeEnergy_ = false;
1000
1001 Log::file() << std::endl;
1002 Log::file() << std::endl;
1003
1004 // Call iterator (return 0 for convergence, 1 for failure)
1005 int error = iterator().solve(isContinuation);
1006
1007 hasCFields_ = true;
1008
1009 if (w_.isSymmetric()) {
1010 fieldIo().convertRGridToBasis(c_.rgrid(), c_.basis());
1011 }
1012
1013 // If converged, compute related thermodynamic properties
1014 if (!error) {
1015 if (!iterator().isFlexible()) {
1016 mixture_.computeStress(mask().phiTot());
1017 }
1018 computeFreeEnergy(); // Sets hasFreeEnergy_ = true
1019 writeThermo(Log::file());
1020 if (!iterator().isFlexible()) {
1021 writeStress(Log::file());
1022 }
1023 }
1024 return error;
1025 }
1026
1027 /*
1028 * Perform sweep along a line in parameter space.
1029 */
1030 template <int D>
1032 {
1033 UTIL_CHECK(hasIterator());
1034 UTIL_CHECK(w_.hasData());
1035 if (iterator().isSymmetric()) {
1036 UTIL_CHECK(w_.isSymmetric());
1037 UTIL_CHECK(isAllocatedBasis_);
1038 }
1039 UTIL_CHECK(hasSweep());
1040 hasCFields_ = false;
1041 hasFreeEnergy_ = false;
1042
1043 Log::file() << std::endl;
1044 Log::file() << std::endl;
1045
1046 // Perform sweep
1047 sweepPtr_->sweep();
1048 }
1049
1050 /*
1051 * Perform a field theoretic simulation of nStep steps.
1052 */
1053 template <int D>
1054 void System<D>::simulate(int nStep)
1055 {
1056 UTIL_CHECK(nStep > 0);
1057 UTIL_CHECK(hasSimulator());
1058 hasCFields_ = false;
1059 hasFreeEnergy_ = false;
1060
1061 simulator().simulate(nStep);
1062 hasCFields_ = true;
1063 }
1064
1065 /*
1066 * Expand the number of spatial dimensions of an RField.
1067 */
1068 template <int D>
1069 void System<D>::expandRGridDimension(std::string const & inFileName,
1070 std::string const & outFileName,
1071 int d,
1072 DArray<int> newGridDimensions)
1073 {
1074 UTIL_CHECK(d > D);
1075
1076 // Read fields
1077 UnitCell<D> tmpUnitCell;
1078 domain().fieldIo().readFieldsRGrid(inFileName,
1079 tmpFieldsRGrid_,
1080 tmpUnitCell);
1081
1082 // Expand Fields
1083 domain().fieldIo().expandRGridDimension(outFileName,
1084 tmpFieldsRGrid_,
1085 tmpUnitCell,
1086 d, newGridDimensions);
1087 }
1088
1089 /*
1090 * Replicate unit cell a specified number of times in each direction.
1091 */
1092 template <int D>
1093 void System<D>::replicateUnitCell(std::string const & inFileName,
1094 std::string const & outFileName,
1095 IntVec<D> const & replicas)
1096 {
1097 // Read fields
1098 UnitCell<D> tmpUnitCell;
1099 domain().fieldIo().readFieldsRGrid(inFileName, tmpFieldsRGrid_,
1100 tmpUnitCell);
1101
1102 // Replicate fields
1103 domain().fieldIo().replicateUnitCell(outFileName,
1104 tmpFieldsRGrid_,
1105 tmpUnitCell,
1106 replicas);
1107 }
1108
1109 // Thermodynamic Properties
1110
1111 /*
1112 * Compute Helmholtz free energy and pressure
1113 */
1114 template <int D>
1116 {
1117 if (hasFreeEnergy_) return;
1118
1119 UTIL_CHECK(w_.hasData());
1120 UTIL_CHECK(hasCFields_);
1121
1122 int nm = mixture().nMonomer(); // number of monomer types
1123 int np = mixture().nPolymer(); // number of polymer species
1124 int ns = mixture().nSolvent(); // number of solvent species
1125
1126 // Initialize all free energy contributions to zero
1127 fHelmholtz_ = 0.0;
1128 fIdeal_ = 0.0;
1129 fInter_ = 0.0;
1130 fExt_ = 0.0;
1131
1132 double phi, mu;
1133
1134 // Compute polymer ideal gas contributions to fIdeal_
1135 if (np > 0) {
1136 Polymer<D>* polymerPtr;
1137 double length;
1138 for (int i = 0; i < np; ++i) {
1139 polymerPtr = &mixture_.polymer(i);
1140 phi = polymerPtr->phi();
1141 mu = polymerPtr->mu();
1142 // Recall: mu = ln(phi/q)
1143 length = polymerPtr->length();
1144 if (phi > 1.0E-08) {
1145 fIdeal_ += phi*( mu - 1.0 )/length;
1146 }
1147 }
1148 }
1149
1150 // Compute solvent ideal gas contributions to fHelhmoltz_
1151 if (ns > 0) {
1152 Solvent<D>* solventPtr;
1153 double size;
1154 for (int i = 0; i < ns; ++i) {
1155 solventPtr = &mixture_.solvent(i);
1156 phi = solventPtr->phi();
1157 mu = solventPtr->mu();
1158 size = solventPtr->size();
1159 if (phi > 1.0E-08) {
1160 fIdeal_ += phi*( mu - 1.0 )/size;
1161 }
1162 }
1163 }
1164
1165 // Volume integrals with a mask: If the system has a mask, then the
1166 // volume that should be used in calculating free energy/pressure
1167 // is the volume available to the material, not the total unit cell
1168 // volume. We thus divide all terms that involve integrating over
1169 // the unit cell volume by quantity mask().phiTot(), which is the
1170 // volume fraction of the unit cell that is occupied by material.
1171 // This properly scales them to the correct value. fExt_, fInter_,
1172 // and the Legendre transform component of fIdeal_ all require
1173 // this scaling. If no mask is present, mask.phiTot() = 1 and no
1174 // scaling occurs.
1175 const double phiTot = mask().phiTot();
1176
1177 // Compute Legendre transform subtraction from fIdeal_
1178 double temp = 0.0;
1179 if (w_.isSymmetric()) {
1180 // Use expansion in symmetry-adapted orthonormal basis
1181 const int nBasis = domain_.basis().nBasis();
1182 for (int i = 0; i < nm; ++i) {
1183 for (int k = 0; k < nBasis; ++k) {
1184 temp -= w_.basis(i)[k] * c_.basis(i)[k];
1185 }
1186 }
1187 } else {
1188 // Use summation over grid points
1189 const int meshSize = domain().mesh().size();
1190 for (int i = 0; i < nm; ++i) {
1191 temp -= Reduce::innerProduct(w_.rgrid(i), c_.rgrid(i));
1192 }
1193 temp /= double(meshSize);
1194 }
1195 temp /= phiTot;
1196 fIdeal_ += temp;
1197 fHelmholtz_ += fIdeal_;
1198
1199 // Compute contribution from external fields, if they exist
1200 if (hasExternalFields()) {
1201 if (iterator().isSymmetric()) {
1202 // Use expansion in symmetry-adapted orthonormal basis
1203 const int nBasis = domain_.basis().nBasis();
1204 for (int i = 0; i < nm; ++i) {
1205 for (int k = 0; k < nBasis; ++k) {
1206 fExt_ += h_.basis(i)[k] * c_.basis(i)[k];
1207 }
1208 }
1209 } else {
1210 // Use summation over grid points
1211 const int meshSize = domain().mesh().size();
1212 for (int i = 0; i < nm; ++i) {
1213 fExt_ += Reduce::innerProduct(h_.rgrid(i), c_.rgrid(i));
1214 }
1215 fExt_ /= double(meshSize);
1216 }
1217 fExt_ /= phiTot;
1218 fHelmholtz_ += fExt_;
1219 }
1220
1221 // Compute excess interaction free energy [ phi^{T}*chi*phi/2 ]
1222 if (w_.isSymmetric()) {
1223 const int nBasis = domain_.basis().nBasis();
1224 for (int i = 0; i < nm; ++i) {
1225 for (int j = i; j < nm; ++j) {
1226 const double chi = interaction().chi(i,j);
1227 if (std::abs(chi) > 1.0E-9) {
1228 double temp = 0.0;
1229 for (int k = 0; k < nBasis; ++k) {
1230 temp += c_.basis(i)[k] * c_.basis(j)[k];
1231 }
1232 if (i == j) {
1233 fInter_ += 0.5*chi*temp;
1234 } else {
1235 fInter_ += chi*temp;
1236 }
1237 }
1238 }
1239 }
1240 } else {
1241 const int meshSize = domain().mesh().size();
1242 for (int i = 0; i < nm; ++i) {
1243 for (int j = i; j < nm; ++j) {
1244 const double chi = interaction().chi(i,j);
1245 if (std::abs(chi) > 1.0E-9) {
1246 double temp = Reduce::innerProduct(c_.rgrid(i), c_.rgrid(j));
1247 if (i == j) {
1248 fInter_ += 0.5*chi*temp;
1249 } else {
1250 fInter_ += chi*temp;
1251 }
1252 }
1253 }
1254 }
1255 fInter_ /= double(meshSize);
1256 }
1257 fInter_ /= phiTot;
1258 fHelmholtz_ += fInter_;
1259
1260 // Initialize pressure
1261 pressure_ = -fHelmholtz_;
1262
1263 // Polymer corrections to pressure
1264 if (np > 0) {
1265 Polymer<D>* polymerPtr;
1266 double length;
1267 for (int i = 0; i < np; ++i) {
1268 polymerPtr = &mixture_.polymer(i);
1269 phi = polymerPtr->phi();
1270 mu = polymerPtr->mu();
1271 length = polymerPtr->length();
1272 if (phi > 1.0E-08) {
1273 pressure_ += mu * phi / length;
1274 }
1275 }
1276 }
1277
1278 // Solvent corrections to pressure
1279 if (ns > 0) {
1280 Solvent<D>* solventPtr;
1281 double size;
1282 for (int i = 0; i < ns; ++i) {
1283 solventPtr = &mixture_.solvent(i);
1284 phi = solventPtr->phi();
1285 mu = solventPtr->mu();
1286 size = solventPtr->size();
1287 if (phi > 1.0E-08) {
1288 pressure_ += mu * phi / size;
1289 }
1290 }
1291 }
1292
1293 hasFreeEnergy_ = true;
1294 }
1295
1296 // Output Operations
1297 /*
1298 * Write time cost to file.
1299 */
1300 template <int D>
1301 void System<D>::writeTimers(std::ostream& out)
1302 {
1303 if (hasIterator()) {
1304 iterator().outputTimers(Log::file());
1305 iterator().outputTimers(out);
1306 }
1307 if (hasSimulator()){
1308 simulator().outputTimers(Log::file());
1309 simulator().outputTimers(out);
1310 }
1311 }
1312
1313 /*
1314 * Clear timers.
1315 */
1316 template <int D>
1318 {
1319 if (hasIterator()) {
1320 iterator().clearTimers();
1321 }
1322 if (hasSimulator()){
1323 simulator().clearTimers();
1324 }
1325 }
1326 /*
1327 * Write parameter file, omitting the sweep block.
1328 */
1329 template <int D>
1330 void System<D>::writeParamNoSweep(std::ostream& out) const
1331 {
1332 out << "System{" << std::endl;
1333 mixture_.writeParam(out);
1334 interaction().writeParam(out);
1335 domain_.writeParam(out);
1336 if (hasIterator()) {
1337 iteratorPtr_->writeParam(out);
1338 }
1339 out << "}" << std::endl;
1340 }
1341
1342 /*
1343 * Write thermodynamic properties to file.
1344 */
1345 template <int D>
1346 void System<D>::writeThermo(std::ostream& out)
1347 {
1348 if (!hasFreeEnergy_) {
1349 computeFreeEnergy();
1350 }
1351
1352 out << std::endl;
1353 out << "fHelmholtz " << Dbl(fHelmholtz(), 18, 11) << std::endl;
1354 out << "pressure " << Dbl(pressure(), 18, 11) << std::endl;
1355 out << std::endl;
1356 out << "fIdeal " << Dbl(fIdeal_, 18, 11) << std::endl;
1357 out << "fInter " << Dbl(fInter_, 18, 11) << std::endl;
1358 if (hasExternalFields()) {
1359 out << "fExt " << Dbl(fExt_, 18, 11) << std::endl;
1360 }
1361 out << std::endl;
1362
1363 int np = mixture_.nPolymer();
1364 int ns = mixture_.nSolvent();
1365
1366 if (np > 0) {
1367 out << "polymers:" << std::endl;
1368 out << " "
1369 << " phi "
1370 << " mu "
1371 << std::endl;
1372 for (int i = 0; i < np; ++i) {
1373 out << Int(i, 5)
1374 << " " << Dbl(mixture_.polymer(i).phi(),18, 11)
1375 << " " << Dbl(mixture_.polymer(i).mu(), 18, 11)
1376 << std::endl;
1377 }
1378 out << std::endl;
1379 }
1380
1381 if (ns > 0) {
1382 out << "solvents:" << std::endl;
1383 out << " "
1384 << " phi "
1385 << " mu "
1386 << std::endl;
1387 for (int i = 0; i < ns; ++i) {
1388 out << Int(i, 5)
1389 << " " << Dbl(mixture_.solvent(i).phi(),18, 11)
1390 << " " << Dbl(mixture_.solvent(i).mu(), 18, 11)
1391 << std::endl;
1392 }
1393 out << std::endl;
1394 }
1395
1396 out << "cellParams:" << std::endl;
1397 for (int i = 0; i < domain_.unitCell().nParameter(); ++i) {
1398 out << Int(i, 5)
1399 << " "
1400 << Dbl(domain_.unitCell().parameter(i), 18, 11)
1401 << std::endl;
1402 }
1403 out << std::endl;
1404 }
1405
1406 /*
1407 * Write stress properties to file.
1408 */
1409 template <int D>
1410 void System<D>::writeStress(std::ostream& out)
1411 {
1412 out << "stress:" << std::endl;
1413 for (int i = 0; i < domain().unitCell().nParameter(); ++i) {
1414 out << Int(i, 5)
1415 << " "
1416 << Dbl(mixture_.stress(i), 18, 11)
1417 << std::endl;
1418 }
1419 out << std::endl;
1420 }
1421
1422 /*
1423 * Write w-fields in symmetry-adapted basis format.
1424 */
1425 template <int D>
1426 void System<D>::writeWBasis(const std::string & filename)
1427 {
1428 UTIL_CHECK(w_.hasData());
1429 UTIL_CHECK(w_.isSymmetric());
1430 fieldIo().writeFieldsBasis(filename, w_.basis(), domain_.unitCell());
1431 }
1432
1433 /*
1434 * Write w-fields in real space grid file format.
1435 */
1436 template <int D>
1437 void System<D>::writeWRGrid(const std::string & filename) const
1438 {
1439 UTIL_CHECK(w_.hasData());
1440 fieldIo().writeFieldsRGrid(filename, w_.rgrid(),
1441 domain_.unitCell(),
1442 w_.isSymmetric());
1443 }
1444
1445 /*
1446 * Write all concentration fields in symmetry-adapted basis format.
1447 */
1448 template <int D>
1449 void System<D>::writeCBasis(const std::string & filename)
1450 {
1451 UTIL_CHECK(hasCFields_);
1452 UTIL_CHECK(w_.isSymmetric());
1453 fieldIo().writeFieldsBasis(filename, c_.basis(), domain_.unitCell());
1454 }
1455
1456 /*
1457 * Write all concentration fields in real space (r-grid) format.
1458 */
1459 template <int D>
1460 void System<D>::writeCRGrid(const std::string & filename) const
1461 {
1462 UTIL_CHECK(hasCFields_);
1463 fieldIo().writeFieldsRGrid(filename, c_.rgrid(),
1464 domain_.unitCell(),
1465 w_.isSymmetric());
1466 }
1467
1468 /*
1469 * Write all concentration fields in real space (r-grid) format, for
1470 * each block (or solvent) individually rather than for each species.
1471 */
1472 template <int D>
1473 void System<D>::writeBlockCRGrid(const std::string & filename) const
1474 {
1475 UTIL_CHECK(hasCFields_);
1476
1477 // Create and allocate the DArray of fields to be written
1478 DArray< RField<D> > blockCFields;
1479 blockCFields.allocate(mixture_.nSolvent() + mixture_.nBlock());
1480 int n = blockCFields.capacity();
1481 for (int i = 0; i < n; i++) {
1482 blockCFields[i].allocate(domain_.mesh().dimensions());
1483 }
1484
1485 // Get data from Mixture and write to file
1486 mixture_.createBlockCRGrid(blockCFields);
1487 fieldIo().writeFieldsRGrid(filename, blockCFields,
1488 domain_.unitCell(),
1489 w_.isSymmetric());
1490 }
1491
1492 /*
1493 * Write the last time slice of the propagator in r-grid format.
1494 */
1495 template <int D>
1496 void System<D>::writeQSlice(const std::string & filename,
1497 int polymerId, int blockId,
1498 int directionId, int segmentId)
1499 const
1500 {
1501 UTIL_CHECK(polymerId >= 0);
1502 UTIL_CHECK(polymerId < mixture_.nPolymer());
1503 Polymer<D> const& polymer = mixture_.polymer(polymerId);
1504 UTIL_CHECK(blockId >= 0);
1505 UTIL_CHECK(blockId < polymer.nBlock());
1506 UTIL_CHECK(directionId >= 0);
1507 UTIL_CHECK(directionId <= 1);
1508 Propagator<D> const&
1509 propagator = polymer.propagator(blockId, directionId);
1510 fieldIo().writeFieldRGrid(filename, propagator.q(segmentId),
1511 domain_.unitCell(), w_.isSymmetric());
1512 }
1513
1514 /*
1515 * Write the last time slice of the propagator in r-grid format.
1516 */
1517 template <int D>
1518 void System<D>::writeQTail(const std::string & filename,
1519 int polymerId, int blockId, int directionId)
1520 const
1521 {
1522 UTIL_CHECK(polymerId >= 0);
1523 UTIL_CHECK(polymerId < mixture_.nPolymer());
1524 Polymer<D> const& polymer = mixture_.polymer(polymerId);
1525 UTIL_CHECK(blockId >= 0);
1526 UTIL_CHECK(blockId < polymer.nBlock());
1527 UTIL_CHECK(directionId >= 0);
1528 UTIL_CHECK(directionId <= 1);
1529 Propagator<D> const&
1530 propagator = polymer.propagator(blockId, directionId);
1531 fieldIo().writeFieldRGrid(filename, propagator.tail(),
1532 domain_.unitCell(), w_.isSymmetric());
1533 }
1534
1535 /*
1536 * Write the propagator for a block and direction.
1537 */
1538 template <int D>
1539 void System<D>::writeQ(const std::string & filename,
1540 int polymerId, int blockId, int directionId)
1541 const
1542 {
1543 UTIL_CHECK(polymerId >= 0);
1544 UTIL_CHECK(polymerId < mixture_.nPolymer());
1545 Polymer<D> const& polymer = mixture_.polymer(polymerId);
1546 UTIL_CHECK(blockId >= 0);
1547 UTIL_CHECK(blockId < polymer.nBlock());
1548 UTIL_CHECK(directionId >= 0);
1549 UTIL_CHECK(directionId <= 1);
1550 Propagator<D> const&
1551 propagator = polymer.propagator(blockId, directionId);
1552 int ns = propagator.ns();
1553
1554 // Open file
1555 std::ofstream file;
1556 fileMaster_.openOutputFile(filename, file);
1557
1558 // Write header
1559 fieldIo().writeFieldHeader(file, 1, domain_.unitCell(),
1560 w_.isSymmetric());
1561 file << "ngrid" << std::endl
1562 << " " << domain_.mesh().dimensions() << std::endl
1563 << "nslice" << std::endl
1564 << " " << ns << std::endl;
1565
1566 // Write data
1567 bool hasHeader = false;
1568 for (int i = 0; i < ns; ++i) {
1569 file << "slice " << i << std::endl;
1570 fieldIo().writeFieldRGrid(file, propagator.q(i),
1571 domain_.unitCell(), hasHeader);
1572 }
1573 }
1574
1575 /*
1576 * Write propagators for all blocks of all polymers to files.
1577 */
1578 template <int D>
1579 void System<D>::writeQAll(std::string const & basename)
1580 {
1581 std::string filename;
1582 int np, nb, ip, ib, id;
1583 np = mixture_.nPolymer();
1584 for (ip = 0; ip < np; ++ip) {
1585 //Polymer<D> const * polymerPtr = &mixture_.polymer(ip);
1586 //nb = polymerPtr->nBlock();
1587 nb = mixture_.polymer(ip).nBlock();
1588 for (ib = 0; ib < nb; ++ib) {
1589 for (id = 0; id < 2; ++id) {
1590 filename = basename;
1591 filename += "_";
1592 filename += toString(ip);
1593 filename += "_";
1594 filename += toString(ib);
1595 filename += "_";
1596 filename += toString(id);
1597 filename += ".rq";
1598 writeQ(filename, ip, ib, id);
1599 }
1600 }
1601 }
1602 }
1603
1604 /*
1605 * Write description of symmetry-adapted stars and basis to file.
1606 */
1607 template <int D>
1608 void System<D>::writeStars(const std::string & filename) const
1609 {
1610 UTIL_CHECK(domain_.hasGroup());
1611 UTIL_CHECK(domain_.basis().isInitialized());
1612 std::ofstream outFile;
1613 fileMaster_.openOutputFile(filename, outFile);
1614 bool isSymmetric = true;
1615 fieldIo().writeFieldHeader(outFile, mixture_.nMonomer(),
1616 domain_.unitCell(), isSymmetric);
1617 basis().outputStars(outFile);
1618 }
1619
1620 /*
1621 * Write a list of waves and associated stars to file.
1622 */
1623 template <int D>
1624 void System<D>::writeWaves(const std::string & filename) const
1625 {
1626 UTIL_CHECK(domain_.hasGroup());
1627 UTIL_CHECK(domain_.basis().isInitialized());
1628 std::ofstream outFile;
1629 fileMaster_.openOutputFile(filename, outFile);
1630 bool isSymmetric = true;
1631 fieldIo().writeFieldHeader(outFile, mixture_.nMonomer(),
1632 domain_.unitCell(), isSymmetric);
1633 basis().outputWaves(outFile);
1634 }
1635
1636 /*
1637 * Write all elements of the space group to a file.
1638 */
1639 template <int D>
1640 void System<D>::writeGroup(const std::string & filename) const
1641 {
1642 UTIL_CHECK(domain_.hasGroup());
1643 Pscf::Prdc::writeGroup(filename, domain_.group());
1644 }
1645
1646 // Field File Operations
1647
1648 /*
1649 * Convert fields from symmetry-adpated basis to real-space grid format.
1650 */
1651 template <int D>
1652 void System<D>::basisToRGrid(const std::string & inFileName,
1653 const std::string & outFileName)
1654 {
1655 // Precondition
1656 UTIL_CHECK(domain_.hasGroup());
1657
1658 // Check required initializations
1659 if (!domain_.unitCell().isInitialized()) {
1660 readFieldHeader(inFileName);
1661 }
1662 UTIL_CHECK(domain_.unitCell().isInitialized());
1663 UTIL_CHECK(domain_.basis().isInitialized());
1664 if (!isAllocatedBasis_) {
1665 allocateFieldsBasis();
1666 }
1667
1668 // Read, convert and write fields
1669 UnitCell<D> tmpUnitCell;
1670 fieldIo().readFieldsBasis(inFileName, tmpFieldsBasis_, tmpUnitCell);
1671 fieldIo().convertBasisToRGrid(tmpFieldsBasis_, tmpFieldsRGrid_);
1672 fieldIo().writeFieldsRGrid(outFileName, tmpFieldsRGrid_,
1673 tmpUnitCell);
1674 }
1675
1676 /*
1677 * Convert fields from real-space grid to symmetry-adapted basis format.
1678 */
1679 template <int D>
1680 void System<D>::rGridToBasis(const std::string & inFileName,
1681 const std::string & outFileName)
1682 {
1683 // Precondition
1684 UTIL_CHECK(domain_.hasGroup());
1685
1686 // Check required initializations
1687 if (!domain_.unitCell().isInitialized()) {
1688 readFieldHeader(inFileName);
1689 }
1690 UTIL_CHECK(domain_.unitCell().isInitialized());
1691 UTIL_CHECK(domain_.basis().isInitialized());
1692 if (!isAllocatedBasis_) {
1693 allocateFieldsBasis();
1694 }
1695
1696 // Read, convert and write fields
1697 UnitCell<D> tmpUnitCell;
1698 fieldIo().readFieldsRGrid(inFileName, tmpFieldsRGrid_, tmpUnitCell);
1699 fieldIo().convertRGridToBasis(tmpFieldsRGrid_, tmpFieldsBasis_);
1700 fieldIo().writeFieldsBasis(outFileName, tmpFieldsBasis_,
1701 tmpUnitCell);
1702 }
1703
1704 /*
1705 * Convert fields from Fourier (k-grid) to real-space (r-grid) format.
1706 */
1707 template <int D>
1708 void System<D>::kGridToRGrid(const std::string & inFileName,
1709 const std::string& outFileName)
1710 {
1711 // Check required initializations
1712 if (!domain_.unitCell().isInitialized()) {
1713 readFieldHeader(inFileName);
1714 }
1715 UTIL_CHECK(domain_.unitCell().isInitialized());
1716
1717 // Read, convert and write fields
1718 UnitCell<D> tmpUnitCell;
1719 fieldIo().readFieldsKGrid(inFileName, tmpFieldsKGrid_, tmpUnitCell);
1720 for (int i = 0; i < mixture_.nMonomer(); ++i) {
1721 fft().inverseTransformUnsafe(tmpFieldsKGrid_[i], tmpFieldsRGrid_[i]);
1722 }
1723 fieldIo().writeFieldsRGrid(outFileName, tmpFieldsRGrid_,
1724 tmpUnitCell);
1725 }
1726
1727 /*
1728 * Convert fields from real-space (r-grid) to Fourier (k-grid) format.
1729 */
1730 template <int D>
1731 void System<D>::rGridToKGrid(const std::string & inFileName,
1732 const std::string & outFileName)
1733 {
1734 // Check required initializations
1735 if (!domain_.unitCell().isInitialized()) {
1736 readFieldHeader(inFileName);
1737 }
1738 UTIL_CHECK(domain_.unitCell().isInitialized());
1739
1740 // Read, convert and write fields
1741 UnitCell<D> tmpUnitCell;
1742 fieldIo().readFieldsRGrid(inFileName, tmpFieldsRGrid_,
1743 tmpUnitCell);
1744 for (int i = 0; i < mixture_.nMonomer(); ++i) {
1745 fft().forwardTransform(tmpFieldsRGrid_[i], tmpFieldsKGrid_[i]);
1746 }
1747 fieldIo().writeFieldsKGrid(outFileName, tmpFieldsKGrid_,
1748 tmpUnitCell);
1749 }
1750
1751 /*
1752 * Convert fields from Fourier (k-grid) to symmetry-adapted basis format.
1753 */
1754 template <int D>
1755 void System<D>::kGridToBasis(const std::string & inFileName,
1756 const std::string& outFileName)
1757 {
1758 // Precondition
1759 UTIL_CHECK(domain_.hasGroup());
1760
1761 // Check for required initializations
1762 if (!domain_.unitCell().isInitialized()) {
1763 readFieldHeader(inFileName);
1764 }
1765 UTIL_CHECK(domain_.unitCell().isInitialized());
1766 UTIL_CHECK(domain_.basis().isInitialized());
1767 if (!isAllocatedBasis_) {
1768 allocateFieldsBasis();
1769 }
1770
1771 // Read, convert and write fields
1772 UnitCell<D> tmpUnitCell;
1773 fieldIo().readFieldsKGrid(inFileName, tmpFieldsKGrid_, tmpUnitCell);
1774 fieldIo().convertKGridToBasis(tmpFieldsKGrid_, tmpFieldsBasis_);
1775 fieldIo().writeFieldsBasis(outFileName, tmpFieldsBasis_,
1776 tmpUnitCell);
1777 }
1778
1779 /*
1780 * Convert fields from symmetry-adapted basis to Fourier (k-grid) format.
1781 */
1782 template <int D>
1783 void System<D>::basisToKGrid(const std::string & inFileName,
1784 const std::string & outFileName)
1785 {
1786 // Precondition
1787 UTIL_CHECK(domain_.hasGroup());
1788
1789 // Check for required initializations
1790 if (!domain_.unitCell().isInitialized()) {
1791 readFieldHeader(inFileName);
1792 }
1793 UTIL_CHECK(domain_.unitCell().isInitialized());
1794 UTIL_CHECK(domain_.basis().isInitialized());
1795 if (!isAllocatedBasis_) {
1796 allocateFieldsBasis();
1797 }
1798
1799 // Read, convert and write fields
1800 UnitCell<D> tmpUnitCell;
1801 fieldIo().readFieldsBasis(inFileName, tmpFieldsBasis_, tmpUnitCell);
1802 fieldIo().convertBasisToKGrid(tmpFieldsBasis_, tmpFieldsKGrid_);
1803 fieldIo().writeFieldsKGrid(outFileName, tmpFieldsKGrid_,
1804 tmpUnitCell);
1805 }
1806
1807 /*
1808 * Check if r-grid fields have declared space group symmetry.
1809 */
1810 template <int D>
1811 bool
1812 System<D>::checkRGridFieldSymmetry(std::string const & inFileName,
1813 double epsilon)
1814 {
1815 UTIL_CHECK(domain_.hasGroup());
1816
1817 // If basis fields are not allocated, peek at field file header to
1818 // get unit cell parameters, initialize basis and allocate fields.
1819 if (!isAllocatedBasis_) {
1820 readFieldHeader(inFileName);
1821 allocateFieldsBasis();
1822 }
1823
1824 // Read fields
1825 UnitCell<D> tmpUnitCell;
1826 domain_.fieldIo().readFieldsRGrid(inFileName,
1827 tmpFieldsRGrid_, tmpUnitCell);
1828
1829 // Check symmetry for all fields
1830 for (int i = 0; i < mixture().nMonomer(); ++i) {
1831 bool symmetric;
1832 symmetric = domain_.fieldIo().hasSymmetry(tmpFieldsRGrid_[i],
1833 epsilon);
1834 if (!symmetric) {
1835 return false;
1836 }
1837 }
1838 return true;
1839
1840 }
1841
1842 /*
1843 * Rescale a field by a constant, write rescaled field to a file.
1844 */
1845 template <int D>
1846 void System<D>::scaleFieldsBasis(std::string const & inFileName,
1847 std::string const & outFileName,
1848 double factor)
1849 {
1850 // If basis fields are not allocated, peek at field file header to
1851 // get unit cell parameters, initialize basis and allocate fields.
1852 if (!isAllocatedBasis_) {
1853 readFieldHeader(inFileName);
1854 allocateFieldsBasis();
1855 }
1856
1857 UnitCell<D> tmpUnitCell;
1858 FieldIo<D> const & fieldIo = domain().fieldIo();
1859 fieldIo.readFieldsBasis(inFileName, tmpFieldsBasis_, tmpUnitCell);
1860 fieldIo.scaleFieldsBasis(tmpFieldsBasis_, factor);
1861 fieldIo.writeFieldsBasis(outFileName, tmpFieldsBasis_,
1862 tmpUnitCell);
1863 }
1864
1865 /*
1866 * Rescale a field by a constant, write rescaled field to a file.
1867 */
1868 template <int D>
1869 void System<D>::scaleFieldsRGrid(std::string const & inFileName,
1870 std::string const & outFileName,
1871 double factor) const
1872 {
1873 UnitCell<D> tmpUnitCell;
1874 FieldIo<D> const & fieldIo = domain().fieldIo();
1875 fieldIo.readFieldsRGrid(inFileName, tmpFieldsRGrid_,
1876 tmpUnitCell);
1877 fieldIo.scaleFieldsRGrid(tmpFieldsRGrid_, factor);
1878 fieldIo.writeFieldsRGrid(outFileName, tmpFieldsRGrid_,
1879 tmpUnitCell);
1880 }
1881
1882 /*
1883 * Compare two fields in basis format, write report to Log file.
1884 */
1885 template <int D>
1887 const DArray< DArray<double> > field2)
1888 {
1889 UTIL_CHECK(domain_.hasGroup());
1890
1891 BFieldComparison comparison(1);
1892 comparison.compare(field1,field2);
1893
1894 Log::file() << "\n Basis expansion field comparison results"
1895 << std::endl;
1896 Log::file() << " Maximum Absolute Difference: "
1897 << comparison.maxDiff() << std::endl;
1898 Log::file() << " Root-Mean-Square Difference: "
1899 << comparison.rmsDiff() << "\n" << std::endl;
1900 }
1901
1902 /*
1903 * Compare two fields in r-grid format, output report to Log file.
1904 */
1905 template <int D>
1907 const DArray< RField<D> > field2)
1908 {
1909 RFieldComparison<D> comparison;
1910 comparison.compare(field1, field2);
1911
1912 Log::file() << "\n Real-space field comparison results"
1913 << std::endl;
1914 Log::file() << " Maximum Absolute Difference: "
1915 << comparison.maxDiff() << std::endl;
1916 Log::file() << " Root-Mean-Square Difference: "
1917 << comparison.rmsDiff() << "\n" << std::endl;
1918 }
1919
1920 // Private member functions
1921
1922 /*
1923 * Allocate memory for fields in grid format.
1924 */
1925 template <int D>
1927 {
1928 // Preconditions
1929 UTIL_CHECK(hasMixture_);
1930 int nMonomer = mixture_.nMonomer();
1931 UTIL_CHECK(nMonomer > 0);
1932 UTIL_CHECK(domain_.mesh().size() > 0);
1933 UTIL_CHECK(!isAllocatedGrid_);
1934
1935 // Alias for mesh dimensions
1936 IntVec<D> const & dimensions = domain_.mesh().dimensions();
1937
1938 // Allocate c (chemical potential) fields
1939 w_.setNMonomer(nMonomer);
1940 w_.allocateRGrid(dimensions);
1941
1942 // Allocate c (monomer concentration) fields
1943 c_.setNMonomer(nMonomer);
1944 c_.allocateRGrid(dimensions);
1945
1946 h_.setNMonomer(nMonomer);
1947
1948 // Allocate work space field arrays
1949 tmpFieldsRGrid_.allocate(nMonomer);
1950 tmpFieldsKGrid_.allocate(nMonomer);
1951 for (int i = 0; i < nMonomer; ++i) {
1952 tmpFieldsRGrid_[i].allocate(dimensions);
1953 tmpFieldsKGrid_[i].allocate(dimensions);
1954 }
1955 workArray_.allocate(dimensions);
1956
1957 isAllocatedGrid_ = true;
1958 }
1959
1960 /*
1961 * Allocate memory for fields in basis format.
1962 */
1963 template <int D>
1964 void System<D>::allocateFieldsBasis()
1965 {
1966 // Preconditions and constants
1967 UTIL_CHECK(hasMixture_);
1968 const int nMonomer = mixture_.nMonomer();
1969 UTIL_CHECK(nMonomer > 0);
1970 UTIL_CHECK(domain_.unitCell().isInitialized());
1971 UTIL_CHECK(domain_.unitCell().nParameter() > 0);
1972 UTIL_CHECK(domain_.mesh().size() > 0);
1973 UTIL_CHECK(domain_.basis().isInitialized());
1974 UTIL_CHECK(isAllocatedGrid_);
1975 UTIL_CHECK(!isAllocatedBasis_);
1976 const int nBasis = basis().nBasis();
1977 UTIL_CHECK(nBasis > 0);
1978
1979 // Allocate basis field containers
1980 w_.allocateBasis(nBasis);
1981 c_.allocateBasis(nBasis);
1982
1983 // Temporary work space
1984 tmpFieldsBasis_.allocate(nMonomer);
1985 for (int i = 0; i < nMonomer; ++i) {
1986 tmpFieldsBasis_[i].allocate(nBasis);
1987 }
1988
1989 isAllocatedBasis_ = true;
1990 }
1991
1992 /*
1993 * Peek at field file header, initialize unit cell parameters and basis.
1994 */
1995 template <int D>
1996 void System<D>::readFieldHeader(std::string filename)
1997 {
1998 UTIL_CHECK(hasMixture_);
1999 UTIL_CHECK(mixture_.nMonomer() > 0);
2000 UTIL_CHECK(isAllocatedGrid_);
2001
2002 // Open field file
2003 std::ifstream file;
2004 fileMaster_.openInputFile(filename, file);
2005
2006 // Read field file header.
2007 int nMonomer;
2008 bool isSymmetric;
2009 domain_.fieldIo().readFieldHeader(file, nMonomer,
2010 domain_.unitCell(), isSymmetric);
2011 // Note: FieldIo<D>::readFieldHeader computes WaveList minimum
2012 // images and initializes the Basis if needed.
2013 file.close();
2014
2015 // Postconditions
2016 UTIL_CHECK(mixture_.nMonomer() == nMonomer);
2017 UTIL_CHECK(domain_.unitCell().nParameter() > 0);
2018 UTIL_CHECK(domain_.unitCell().lattice() != UnitCell<D>::Null);
2019 UTIL_CHECK(domain_.unitCell().isInitialized());
2020 if (domain_.hasGroup()) {
2021 UTIL_CHECK(domain_.basis().isInitialized());
2022 UTIL_CHECK(domain_.basis().nBasis() > 0);
2023 }
2024 }
2025
2026 /*
2027 * Read a filename string and echo to log file (used in readCommands).
2028 */
2029 template <int D>
2030 void System<D>::readEcho(std::istream& in, std::string& string) const
2031 {
2032 in >> string;
2033 Log::file() << " " << Str(string, 20) << std::endl;
2034 }
2035
2036 /*
2037 * Read floating point number, echo to log file (used in readCommands).
2038 */
2039 template <int D>
2040 void System<D>::readEcho(std::istream& in, double& value) const
2041 {
2042 in >> value;
2043 if (in.fail()) {
2044 UTIL_THROW("Unable to read floating point parameter.");
2045 }
2046 Log::file() << " " << Dbl(value, 20) << std::endl;
2047 }
2048
2049 /*
2050 * Initialize Pscf::Homogeneous::Mixture homogeneous_ member.
2051 */
2052 template <int D>
2053 void System<D>::initHomogeneous()
2054 {
2055
2056 // Set number of molecular species and monomers
2057 int nm = mixture_.nMonomer();
2058 int np = mixture_.nPolymer();
2059 int ns = mixture_.nSolvent();
2060
2061 UTIL_CHECK(homogeneous_.nMolecule() == np + ns);
2062 UTIL_CHECK(homogeneous_.nMonomer() == nm);
2063
2064 int i; // molecule index
2065 int j; // monomer index
2066
2067 // Loop over polymer molecule species
2068 if (np > 0) {
2069
2070 DArray<double> cTmp;
2071 cTmp.allocate(nm);
2072
2073 int k; // block or clump index
2074 int nb; // number of blocks
2075 int nc; // number of clumps
2076 for (i = 0; i < np; ++i) {
2077
2078 // Initial array of clump sizes
2079 for (j = 0; j < nm; ++j) {
2080 cTmp[j] = 0.0;
2081 }
2082
2083 // Compute clump sizes for all monomer types.
2084 nb = mixture_.polymer(i).nBlock();
2085 for (k = 0; k < nb; ++k) {
2086 Block<D>& block = mixture_.polymer(i).block(k);
2087 j = block.monomerId();
2088 cTmp[j] += block.length();
2089 }
2090
2091 // Count the number of clumps of nonzero size
2092 nc = 0;
2093 for (j = 0; j < nm; ++j) {
2094 if (cTmp[j] > 1.0E-8) {
2095 ++nc;
2096 }
2097 }
2098 homogeneous_.molecule(i).setNClump(nc);
2099
2100 // Set clump properties for this Homogeneous::Molecule
2101 k = 0; // Clump index
2102 for (j = 0; j < nm; ++j) {
2103 if (cTmp[j] > 1.0E-8) {
2104 homogeneous_.molecule(i).clump(k).setMonomerId(j);
2105 homogeneous_.molecule(i).clump(k).setSize(cTmp[j]);
2106 ++k;
2107 }
2108 }
2109 homogeneous_.molecule(i).computeSize();
2110
2111 }
2112 }
2113
2114 // Add solvent contributions
2115 if (ns > 0) {
2116 double size;
2117 int monomerId;
2118 for (int is = 0; is < ns; ++is) {
2119 i = is + np;
2120 monomerId = mixture_.solvent(is).monomerId();
2121 size = mixture_.solvent(is).size();
2122 homogeneous_.molecule(i).setNClump(1);
2123 homogeneous_.molecule(i).clump(0).setMonomerId(monomerId);
2124 homogeneous_.molecule(i).clump(0).setSize(size);
2125 homogeneous_.molecule(i).computeSize();
2126 }
2127 }
2128 }
2129
2130} // namespace Rpg
2131} // namespace Pscf
2132#endif
Dynamic array on the GPU device with aligned data.
Definition rpg/System.h:32
double rmsDiff() const
Return the precomputed root-mean-squared difference.
double compare(FT const &a, FT const &b)
Compare individual fields.
double maxDiff() const
Return the precomputed maximum element-by-element difference.
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Flory-Huggins excess free energy model.
Definition Interaction.h:26
int nBlock() const
Number of blocks.
Propagator & propagator(int blockId, int directionId)
Get propagator for a specific block and direction.
double length() const
Sum of the lengths of all blocks in the polymer.
Comparator for fields in symmetry-adapted basis format.
Comparator for fields in real-space (r-grid) format.
Field of real double precision values on an FFT mesh.
void scaleFieldsRGrid(DArray< RFRT > &fields, double factor) const
Scale an array of r-grid fields by a real scalar.
void readFieldsBasis(std::istream &in, DArray< DArray< double > > &fields, UnitCell< D > &unitCell) const
Read concentration or chemical potential fields from file.
void scaleFieldsBasis(DArray< DArray< double > > &fields, double factor) const
Scale an array of fields in basis format by a real scalar.
void writeFieldsBasis(std::ostream &out, DArray< DArray< double > > const &fields, UnitCell< D > const &unitCell) const
Write concentration or chemical potential field components to file.
Base template for UnitCell<D> classes, D=1, 2 or 3.
Definition rpg/System.h:34
File input/output operations and format conversions for fields.
void writeFieldsRGrid(std::ostream &out, DArray< RField< D > > const &fields, UnitCell< D > const &unitCell, bool writeHeader=true, bool isSymmetric=true, bool writeMeshSize=true) const
Write array of RField objects (fields on r-space grid) to a stream.
void readFieldsRGrid(std::istream &in, DArray< RField< D > > &fields, UnitCell< D > &unitCell) const
Read array of RField objects (r-grid fields) from a stream.
Factory for subclasses of Iterator.
Definition rpg/System.h:38
Descriptor and solver for a branched polymer species.
MDE solver for one-direction of one block.
int ns() const
Get the number of chain contour points.
RField< D > const & tail() const
Return q-field at end of block (after propagator is solved).
RField< D > const & q(int i) const
Return const q-field at specified step by reference (after solving).
Factory for subclasses of Simulator.
Definition rpg/System.h:42
Solver and descriptor for a solvent species.
Default Factory for subclasses of Sweep.
Definition rpg/System.h:40
Main class for calculations that represent one system.
Definition rpg/System.h:107
void expandRGridDimension(const std::string &inFileName, const std::string &outFileName, int d, DArray< int > newGridDimensions)
Expand the number of spatial dimensions of an r-grid field.
void kGridToBasis(const std::string &inFileName, const std::string &outFileName)
Convert fields from Fourier (k-grid) to symmetrized basis format.
void writeCBasis(const std::string &filename)
Write concentrations in symmetry-adapted basis format.
void estimateWfromC(const std::string &filename)
Construct trial w-fields from c-fields.
void writeQTail(std::string const &filename, int polymerId, int blockId, int directionId) const
Write the final slice of a propagator in r-grid format.
void readParam()
Read input parameters from default param file.
void computeFreeEnergy()
Compute free energy density and pressure for current fields.
void writeBlockCRGrid(const std::string &filename) const
Write c fields for all blocks and solvents in r-grid format.
void compute(bool needStress=false)
Solve the modified diffusion equation once, without iteration.
void setWRGrid(DArray< RField< D > > const &fields)
Set new w fields, in real-space (r-grid) format.
void writeQ(std::string const &filename, int polymerId, int blockId, int directionId) const
Write one propagator for one block, in r-grid format.
void writeTimers(std::ostream &out)
Write timer file to an ostream.
void rGridToKGrid(const std::string &inFileName, const std::string &outFileName)
Convert fields from real-space (r-grid) to Fourier (k-grid) format.
void rGridToBasis(const std::string &inFileName, const std::string &outFileName)
Convert a field from real-space grid to symmetrized basis format.
void setWBasis(DArray< DArray< double > > const &fields)
Set chemical potential fields, in symmetry-adapted basis format.
void kGridToRGrid(const std::string &inFileName, const std::string &outFileName)
Convert fields from Fourier (k-grid) to real-space (r-grid) format.
void writeParamNoSweep(std::ostream &out) const
Write parameter file to an ostream, omitting the sweep block.
void readCommands()
Read and process commands from the default command file.
void scaleFieldsRGrid(const std::string &inFileName, const std::string &outFileName, double factor) const
Multiply all elements of an array of r-grid fields by a scalar.
void writeWBasis(const std::string &filename)
Write chemical potential fields in symmetry adapted basis format.
void writeGroup(std::string const &filename) const
Output all elements of the space group.
void writeWRGrid(const std::string &filename) const
Write chemical potential fields in real space grid (r-grid) format.
void symmetrizeWFields()
Symmetrize r-grid w-fields, compute basis components.
void readWBasis(const std::string &filename)
Read chemical potential fields in symmetry adapted basis format.
void writeStars(const std::string &filename) const
Output information about stars and symmetrized basis functions.
virtual void readParameters(std::istream &in)
Read body of parameter file (without opening, closing lines).
void sweep()
Sweep in parameter space, solving an SCF problem at each point.
void basisToKGrid(const std::string &inFileName, const std::string &outFileName)
Convert fields from symmetrized basis to Fourier (k-grid) format.
~System()
Destructor.
int iterate(bool isContinuation=false)
Iteratively solve a SCFT problem.
void basisToRGrid(const std::string &inFileName, const std::string &outFileName)
Convert a field from symmetry-adapted basis to r-grid format.
System()
Constructor.
void writeStress(std::ostream &out)
Write stress properties to a file.
void readWRGrid(const std::string &filename)
Read chemical potential fields in real space grid (r-grid) format.
void setOptions(int argc, char **argv)
Process command line options.
void clearTimers()
Clear timers.
void writeWaves(const std::string &filename) const
Output information about waves.
void replicateUnitCell(const std::string &inFileName, const std::string &outFileName, IntVec< D > const &replicas)
Replicate the crystal unit cell to create a larger cell.
void writeQSlice(std::string const &filename, int polymerId, int blockId, int directionId, int segmentId) const
Write specified slice of a propagator at fixed s in r-grid format.
void writeCRGrid(const std::string &filename) const
Write concentration fields in real space grid (r-grid) format.
void writeQAll(std::string const &basename)
Write all propagators of all blocks, each to a separate file.
void compare(const DArray< DArray< double > > field1, const DArray< DArray< double > > field2)
Compare arrays of fields in basis format, output a report.
bool checkRGridFieldSymmetry(const std::string &inFileName, double epsilon=1.0E-8)
Check if r-grid fields have the declared space group symmetry.
void setUnitCell(UnitCell< D > const &unitCell)
Set parameters of the associated unit cell.
void writeThermo(std::ostream &out)
Write thermodynamic properties to a file.
void simulate(int nStep)
Perform a field theoretic simulation.
void scaleFieldsBasis(const std::string &inFileName, const std::string &outFileName, double factor)
Multiply all components of an array of basis fields by a scalar.
double size() const
Get the size (number of monomers) in this solvent.
double phi() const
Get the overall volume fraction for this species.
Definition Species.h:90
double mu() const
Get the chemical potential for this species (units kT=1).
Definition Species.h:96
int capacity() const
Return allocated size.
Definition Array.h:159
Dynamically allocatable contiguous array template.
void allocate(int capacity)
Allocate the underlying C array.
Definition DArray.h:199
Wrapper for a double precision number, for formatted ostream output.
Definition Dbl.h:40
A fixed capacity (static) contiguous array with a variable logical size.
Definition rpg/System.h:28
Wrapper for an int, for formatted ostream output.
Definition Int.h:37
static std::ostream & file()
Get log ostream by reference.
Definition Log.cpp:57
static bool echo()
Get echo parameter.
static void setEcho(bool echo=true)
Enable or disable echoing for all subclasses of ParamComponent.
void setClassName(const char *className)
Set class name string.
Wrapper for a std::string, for formatted ostream output.
Definition Str.h:37
#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
std::string toString(int n)
Return string representation of an integer.
Definition ioUtil.cpp:52
void readFieldHeader(std::istream &in, int &ver1, int &ver2, UnitCell< D > &cell, std::string &groupName, int &nMonomer)
Read common part of field header (fortran PSCF format).
void replicateUnitCell(std::ostream &out, DArray< AT > const &fields, IntVec< D > const &meshDimensions, UnitCell< D > const &unitCell, IntVec< D > const &replicas)
Write r-grid fields in a replicated unit cell to std::ostream.
bool hasSymmetry(AT const &in, Basis< D > const &basis, IntVec< D > const &dftDimensions, double epsilon=1.0e-8, bool verbose=true)
Check if a k-grid field has the declared space group symmetry.
void expandRGridDimension(std::ostream &out, DArray< AT > const &fields, IntVec< D > const &meshDimensions, UnitCell< D > const &unitCell, int d, DArray< int > newGridDimensions)
Expand the dimensionality of space from D to d.
void init()
Initialize static variables in Pscf::ThreadArray namespace.
void setThreadsPerBlock()
Set the number of threads per block to a default value.
void setThreadsPerBlock(int blockSize)
Manually set the block size that should be used by default.
cudaReal innerProduct(DeviceArray< cudaReal > const &a, DeviceArray< cudaReal > const &b)
Compute inner product of two real arrays (GPU kernel wrapper).
Definition Reduce.cu:891
Fields, FFTs, and utilities for periodic boundary conditions (CUDA)
Definition CField.cu:12
Periodic fields and crystallography.
Definition CField.cpp:11
PSCF package top-level namespace.
Definition param_pc.dox:1
void set(BracketPolicy::Type policy)
Set policy regarding use of bracket delimiters on arrays.
Utility classes for scientific computation.