PSCF v1.2
rpc/System.tpp
1#ifndef RPC_SYSTEM_TPP
2#define RPC_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 <rpc/fts/simulator/Simulator.h>
14#include <rpc/fts/simulator/SimulatorFactory.h>
15#include <rpc/fts/compressor/Compressor.h>
16#include <rpc/scft/sweep/Sweep.h>
17#include <rpc/scft/sweep/SweepFactory.h>
18#include <rpc/scft/iterator/Iterator.h>
19#include <rpc/scft/iterator/IteratorFactory.h>
20#include <rpc/solvers/Polymer.h>
21#include <rpc/solvers/Solvent.h>
22
23
24#include <prdc/cpu/RField.h>
25#include <prdc/cpu/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 Rpc {
45
46 using namespace Util;
47 using namespace Pscf::Prdc;
48 using namespace Pscf::Prdc::Cpu;
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);
92
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; // nThread
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':
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 existence and consistency of -d option
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 (tFlag) {
219 if (tArg <= 0) {
220 UTIL_THROW("Error: Non-positive thread count -t option");
221 }
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(), domain().fft(),
256 domain().unitCell());
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 << "\n"
465 << Str("segment ID ", 21) << segmentId
466 << std::endl;
467 writeQSlice(filename, polymerId, blockId, directionId,
468 segmentId);
469 } else
470 if (command == "WRITE_Q_TAIL") {
471 readEcho(in, filename);
472 int polymerId, blockId, directionId;
473 in >> polymerId;
474 in >> blockId;
475 in >> directionId;
476 Log::file() << Str("polymer ID ", 21) << polymerId << "\n"
477 << Str("block ID ", 21) << blockId << "\n"
478 << Str("direction ID ", 21) << directionId
479 << "\n";
480 writeQTail(filename, polymerId, blockId, directionId);
481 } else
482 if (command == "WRITE_Q") {
483 readEcho(in, filename);
484 int polymerId, blockId, directionId;
485 in >> polymerId;
486 in >> blockId;
487 in >> directionId;
488 Log::file() << Str("polymer ID ", 21) << polymerId << "\n"
489 << Str("block ID ", 21) << blockId << "\n"
490 << Str("direction ID ", 21) << directionId
491 << "\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 } else
510 if (command == "BASIS_TO_RGRID") {
511 readEcho(in, inFileName);
512 readEcho(in, outFileName);
513 basisToRGrid(inFileName, outFileName);
514 } else
515 if (command == "RGRID_TO_BASIS") {
516 readEcho(in, inFileName);
517 readEcho(in, outFileName);
518 rGridToBasis(inFileName, outFileName);
519 } else
520 if (command == "KGRID_TO_RGRID") {
521 readEcho(in, inFileName);
522 readEcho(in, outFileName);
523 kGridToRGrid(inFileName, outFileName);
524 } else
525 if (command == "RGRID_TO_KGRID") {
526 readEcho(in, inFileName);
527 readEcho(in, outFileName);
528 rGridToKGrid(inFileName, outFileName);
529 } else
530 if (command == "BASIS_TO_KGRID") {
531 readEcho(in, inFileName);
532 readEcho(in, outFileName);
533 basisToKGrid(inFileName, outFileName);
534 } else
535 if (command == "KGRID_TO_BASIS") {
536 readEcho(in, inFileName);
537 readEcho(in, outFileName);
538 kGridToBasis(inFileName, outFileName);
539 } else
540 if (command == "CHECK_RGRID_SYMMETRY") {
541 double epsilon;
542 readEcho(in, inFileName);
543 readEcho(in, epsilon);
544 bool hasSymmetry;
545 hasSymmetry = checkRGridFieldSymmetry(inFileName, epsilon);
546 if (hasSymmetry) {
547 Log::file() << std::endl
548 << "Symmetry of r-grid file matches this space group."
549 << std::endl << std::endl;
550 } else {
551 Log::file() << std::endl
552 << "Symmetry of r-grid file does not match this\n"
553 << "space group to within error threshold of "
554 << Dbl(epsilon) << " ." << std::endl << std::endl;
555 }
556 } else
557 if (command == "COMPARE_BASIS") {
558
559 // Get two filenames for comparison
560 std::string filecompare1, filecompare2;
561 readEcho(in, filecompare1);
562 readEcho(in, filecompare2);
563
564 DArray< DArray<double> > Bfield1, Bfield2;
565 UnitCell<D> tmpUnitCell;
566 domain().fieldIo().readFieldsBasis(filecompare1, Bfield1,
567 tmpUnitCell);
568 domain().fieldIo().readFieldsBasis(filecompare2, Bfield2,
569 tmpUnitCell);
570 // Note: Bfield1 & Bfield2 are allocated by readFieldsBasis
571
572 // Compare and output report
573 compare(Bfield1, Bfield2);
574
575 } else
576 if (command == "COMPARE_RGRID") {
577 // Get two filenames for comparison
578 std::string filecompare1, filecompare2;
579 readEcho(in, filecompare1);
580 readEcho(in, filecompare2);
581
582 DArray< RField<D> > Rfield1, Rfield2;
583 UnitCell<D> tmpUnitCell;
584 domain().fieldIo().readFieldsRGrid(filecompare1, Rfield1,
585 tmpUnitCell);
586 domain().fieldIo().readFieldsRGrid(filecompare2, Rfield2,
587 tmpUnitCell);
588 // Note: Rfield1, Rfield2 will be allocated by readFieldsRGrid
589
590 // Compare and output report
591 compare(Rfield1, Rfield2);
592
593 } else
594 if (command == "SCALE_BASIS") {
595 double factor;
596 readEcho(in, inFileName);
597 readEcho(in, outFileName);
598 readEcho(in, factor);
599 scaleFieldsBasis(inFileName, outFileName, factor);
600 } else
601 if (command == "SCALE_RGRID") {
602 double factor;
603 readEcho(in, inFileName);
604 readEcho(in, outFileName);
605 readEcho(in, factor);
606 scaleFieldsRGrid(inFileName, outFileName, factor);
607 } else
608 if (command == "EXPAND_RGRID_DIMENSION") {
609 readEcho(in, inFileName);
610 readEcho(in, outFileName);
611
612 // Read expanded spatial dimension d
613 int d;
614 in >> d;
615 UTIL_CHECK(d > D);
616 Log::file() << Str("Expand dimension to: ") << d << "\n";
617
618 // Read numbers of grid points along added dimensions
619 DArray<int> newGridDimensions;
620 newGridDimensions.allocate(d-D);
621 for (int i = 0; i < d-D; i++) {
622 in >> newGridDimensions[i];
623 }
624 Log::file()
625 << Str("Numbers of grid points in added dimensions : ");
626 for (int i = 0; i < d-D; i++) {
627 Log::file() << newGridDimensions[i] << " ";
628 }
629 Log::file() << "\n";
630
631 expandRGridDimension(inFileName, outFileName,
632 d, newGridDimensions);
633
634 } else
635 if (command == "REPLICATE_UNIT_CELL") {
636 readEcho(in, inFileName);
637 readEcho(in, outFileName);
638
639 // Read numbers of replicas along each direction
640 IntVec<D> replicas;
641 for (int i = 0; i < D; i++){
642 in >> replicas[i];
643 }
644 for (int i = 0; i < D; i++){
645 Log::file() << "Replicate unit cell in direction "
646 << i << " : ";
647 Log::file() << replicas[i] << " times ";
648 Log::file() << "\n";
649 }
650
651 replicateUnitCell(inFileName, outFileName, replicas);
652 } else
653 if (command == "READ_H_BASIS") {
654 readEcho(in, filename);
655 if (!h_.isAllocatedBasis()) {
656 h_.allocateBasis(domain().basis().nBasis());
657 }
658 if (!h_.isAllocatedRGrid()) {
659 h_.allocateRGrid(domain().mesh().dimensions());
660 }
661 UnitCell<D> tmpUnitCell;
662 h_.readBasis(filename, tmpUnitCell);
663 } else
664 if (command == "READ_H_RGRID") {
665 readEcho(in, filename);
666 if (!h_.isAllocatedRGrid()) {
667 h_.allocateRGrid(domain().mesh().dimensions());
668 }
669 if (iterator().isSymmetric() && !h_.isAllocatedBasis()) {
670 mask_.allocateBasis(domain().basis().nBasis());
671 }
672 UnitCell<D> tmpUnitCell;
673 h_.readRGrid(filename, tmpUnitCell);
674 } else
675 if (command == "WRITE_H_BASIS") {
676 readEcho(in, filename);
677 UTIL_CHECK(h_.hasData());
678 UTIL_CHECK(h_.isSymmetric());
679 domain().fieldIo().writeFieldsBasis(filename, h_.basis(),
680 domain().unitCell());
681 } else
682 if (command == "WRITE_H_RGRID") {
683 readEcho(in, filename);
684 UTIL_CHECK(h_.hasData());
685 domain().fieldIo().writeFieldsRGrid(filename, h_.rgrid(),
686 domain().unitCell());
687 } else
688 if (command == "READ_MASK_BASIS") {
689 readEcho(in, filename);
690 UTIL_CHECK(domain().basis().isInitialized());
691 if (!mask_.isAllocatedBasis()) {
692 mask_.allocateBasis(domain().basis().nBasis());
693 }
694 if (!mask_.isAllocatedRGrid()) {
695 mask_.allocateRGrid(domain().mesh().dimensions());
696 }
697 UnitCell<D> tmpUnitCell;
698 mask_.readBasis(filename, tmpUnitCell);
699 } else
700 if (command == "READ_MASK_RGRID") {
701 readEcho(in, filename);
702 if (!mask_.isAllocatedRGrid()) {
703 mask_.allocateRGrid(domain().mesh().dimensions());
704 }
705 if (iterator().isSymmetric() && !mask_.isAllocatedBasis()) {
706 UTIL_CHECK(domain().basis().isInitialized());
707 mask_.allocateBasis(domain().basis().nBasis());
708 }
709 UnitCell<D> tmpUnitCell;
710 mask_.readRGrid(filename, tmpUnitCell);
711 } else
712 if (command == "WRITE_MASK_BASIS") {
713 readEcho(in, filename);
714 UTIL_CHECK(mask_.hasData());
715 UTIL_CHECK(mask_.isSymmetric());
716 domain().fieldIo().writeFieldBasis(filename, mask_.basis(),
717 domain().unitCell());
718 } else
719 if (command == "WRITE_MASK_RGRID") {
720 readEcho(in, filename);
721 UTIL_CHECK(mask_.hasData());
722 domain().fieldIo().writeFieldRGrid(filename, mask_.rgrid(),
723 domain().unitCell(),
724 mask_.isSymmetric());
725 } else {
726 Log::file() << "Error: Unknown command "
727 << command << std::endl;
728 readNext = false;
729 }
730 }
731 }
732
733 /*
734 * Read and execute commands from the default command file.
735 */
736 template <int D>
738 {
739 if (fileMaster_.commandFileName().empty()) {
740 UTIL_THROW("Empty command file name");
741 }
742 readCommands(fileMaster_.commandFile());
743 }
744
745 // W Field Modifier Functions
746
747 /*
748 * Read w-field in symmetry adapted basis format.
749 */
750 template <int D>
751 void System<D>::readWBasis(std::string const & filename)
752 {
753 // Precondition
754 UTIL_CHECK(domain().hasGroup());
755
756 if (!domain().unitCell().isInitialized()) {
757 readFieldHeader(filename);
758 }
759 UTIL_CHECK(domain().unitCell().isInitialized());
760 UTIL_CHECK(domain().basis().isInitialized());
761 UTIL_CHECK(domain().basis().nBasis() > 0);
762 if (!isAllocatedBasis_) {
763 allocateFieldsBasis();
764 }
765
766 // Read w fields
767 w_.readBasis(filename, domain_.unitCell());
768
769 // Update UnitCell in Mixture
770 mixture_.clearUnitCellData();
771
772 hasCFields_ = false;
773 hasFreeEnergy_ = false;
774
775 }
776
777 /*
778 * Read w-fields in real-space grid (r-grid) format.
779 */
780 template <int D>
781 void System<D>::readWRGrid(std::string const & filename)
782 {
783 UTIL_CHECK(isAllocatedGrid_);
784
785 // If necessary, peek at header to initialize unit cell
786 if (!domain().unitCell().isInitialized()) {
787 readFieldHeader(filename);
788 }
789 UTIL_CHECK(domain().unitCell().isInitialized());
790 if (domain().hasGroup() && !isAllocatedBasis_) {
791 UTIL_CHECK(domain().basis().isInitialized());
792 UTIL_CHECK(domain().basis().nBasis() > 0);
793 allocateFieldsBasis();
794 }
795
796 // Read w fields
797 w_.readRGrid(filename, domain_.unitCell());
798
799 // Update UnitCell in Mixture
800 mixture_.clearUnitCellData();
801
802 hasCFields_ = false;
803 hasFreeEnergy_ = false;
804 }
805
806 /*
807 * Construct estimate for w fields from c fields, by setting xi=0.
808 *
809 * Modifies wFields and wFieldsRGrid.
810 */
811 template <int D>
812 void System<D>::estimateWfromC(std::string const & filename)
813 {
814 UTIL_CHECK(hasMixture_);
815 const int nm = mixture().nMonomer();
816 UTIL_CHECK(nm > 0);
817 UTIL_CHECK(isAllocatedGrid_);
818 UTIL_CHECK(domain().hasGroup());
819
820 if (!domain().unitCell().isInitialized()) {
821 readFieldHeader(filename);
822 }
823 UTIL_CHECK(domain().unitCell().isInitialized());
824 UTIL_CHECK(domain().basis().isInitialized());
825 const int nb = domain().basis().nBasis();
826 UTIL_CHECK(nb > 0);
827 if (!isAllocatedBasis_) {
828 allocateFieldsBasis();
829 }
830
831 // Read c fields into temporary array and set unit cell
832 domain().fieldIo().readFieldsBasis(filename, tmpFieldsBasis_,
833 domain_.unitCell());
834
835 // Update UnitCell in Mixture
836 mixture_.clearUnitCellData();
837
838 // Allocate work space array
839 DArray<double> wtmp;
840 wtmp.allocate(nm);
841
842 // Compute estimated w fields from c fields
843 int i, j, k;
844 for (i = 0; i < nb; ++i) {
845 for (j = 0; j < nm; ++j) {
846 wtmp[j] = 0.0;
847 for (k = 0; k < nm; ++k) {
848 wtmp[j] += interaction().chi(j,k)*tmpFieldsBasis_[k][i];
849 }
850 }
851 for (j = 0; j < nm; ++j) {
852 tmpFieldsBasis_[j][i] = wtmp[j];
853 }
854 }
855
856 // Set estimated w fields in system w-field container
857 w_.setBasis(tmpFieldsBasis_);
858
859 hasCFields_ = false;
860 hasFreeEnergy_ = false;
861 }
862
863 /*
864 * Set new w-field values.
865 */
866 template <int D>
868 {
869 UTIL_CHECK(domain().hasGroup());
870 UTIL_CHECK(domain().basis().isInitialized());
871 UTIL_CHECK(isAllocatedGrid_);
872 UTIL_CHECK(isAllocatedBasis_);
873 w_.setBasis(fields);
874 hasCFields_ = false;
875 hasFreeEnergy_ = false;
876 }
877
878 /*
879 * Set new w-field values, using r-grid fields as inputs.
880 */
881 template <int D>
882 void System<D>::setWRGrid(DArray< RField<D> > const & fields)
883 {
884 UTIL_CHECK(isAllocatedGrid_);
885 w_.setRGrid(fields);
886 hasCFields_ = false;
887 hasFreeEnergy_ = false;
888 }
889
890 // Unit Cell Modifiers
891
892 /*
893 * Set the system unit cell.
894 */
895 template <int D>
896 void System<D>::setUnitCell(UnitCell<D> const & unitCell)
897 {
898 domain_.setUnitCell(unitCell);
899 // Note: Domain::setUnitCell checks agreement of lattice system
900 mixture_.clearUnitCellData();
901 if (domain().hasGroup() && !isAllocatedBasis_) {
902 allocateFieldsBasis();
903 }
904 }
905
906 /*
907 * Set state of the system unit cell.
908 */
909 template <int D>
910 void
912 FSArray<double, 6> const & parameters)
913 {
914 domain_.setUnitCell(lattice, parameters);
915 // Note: Domain::setUnitCell checks agreement of lattice system
916 mixture_.clearUnitCellData();
917 if (domain().hasGroup() && !isAllocatedBasis_) {
918 allocateFieldsBasis();
919 }
920 }
921
922 /*
923 * Set parameters of the system unit cell.
924 */
925 template <int D>
927 {
928 domain_.setUnitCell(parameters);
929 // Note: Domain::setUnitCell requires lattice system is set on entry
930 mixture_.clearUnitCellData();
931 if (domain().hasGroup() && !isAllocatedBasis_) {
932 allocateFieldsBasis();
933 }
934 }
935
936 // Primary Field Theory Computations
937
938 /*
939 * Solve MDE for current w-fields, without iteration.
940 */
941 template <int D>
942 void System<D>::compute(bool needStress)
943 {
944 UTIL_CHECK(w_.isAllocatedRGrid());
945 UTIL_CHECK(c_.isAllocatedRGrid());
946 UTIL_CHECK(w_.hasData());
947
948 // Solve the modified diffusion equation (without iteration)
949 mixture_.compute(w_.rgrid(), c_.rgrid(), mask_.phiTot());
950 hasCFields_ = true;
951 hasFreeEnergy_ = false;
952
953 // If w fields are symmetric, compute basis components for c-fields
954 if (w_.isSymmetric()) {
955 UTIL_CHECK(c_.isAllocatedBasis());
956 domain().fieldIo().convertRGridToBasis(c_.rgrid(), c_.basis(),
957 false);
958 }
959
960 // Compute stress if it is needed
961 if (needStress) {
962 mixture_.computeStress(mask().phiTot());
963 }
964
965 }
966
967 /*
968 * Iteratively solve a SCFT problem for specified parameters.
969 */
970 template <int D>
971 int System<D>::iterate(bool isContinuation)
972 {
973 UTIL_CHECK(iteratorPtr_);
974 UTIL_CHECK(w_.hasData());
975 if (iterator().isSymmetric()) {
976 UTIL_CHECK(w_.isSymmetric());
977 }
978 hasCFields_ = false;
979 hasFreeEnergy_ = false;
980
981 Log::file() << std::endl;
982 Log::file() << std::endl;
983
984 // Call iterator (return 0 for convergence, 1 for failure)
985 int error = iterator().solve(isContinuation);
986 hasCFields_ = true;
987
988 // If converged, compute related thermodynamic properties
989 if (!error) {
990 computeFreeEnergy(); // Sets hasFreeEnergy_ = true
991 if (!iterator().isFlexible()) {
992 mixture_.computeStress(mask().phiTot());
993 }
994 writeThermo(Log::file());
995 if (!iterator().isFlexible()) {
996 writeStress(Log::file());
997 }
998 }
999
1000 return error;
1001 }
1002
1003 /*
1004 * Perform sweep of SCFT calculations along a path in parameter space.
1005 */
1006 template <int D>
1008 {
1009 UTIL_CHECK(w_.hasData());
1010 UTIL_CHECK(hasSweep());
1011 Log::file() << std::endl;
1012 Log::file() << std::endl;
1013
1014 // Perform SCFT sweep
1015 sweepPtr_->sweep();
1016 }
1017
1018 /*
1019 * Perform a stochast field theoretic simulation of nStep steps.
1020 */
1021 template <int D>
1022 void System<D>::simulate(int nStep)
1023 {
1024 UTIL_CHECK(nStep > 0);
1025 UTIL_CHECK(hasSimulator());
1026 hasCFields_ = false;
1027 hasFreeEnergy_ = false;
1028
1029 simulator().simulate(nStep);
1030 hasCFields_ = true;
1031 }
1032
1033 /*
1034 * Expand the number of spatial dimensions of an RField.
1035 */
1036 template <int D>
1037 void System<D>::expandRGridDimension(std::string const & inFileName,
1038 std::string const & outFileName,
1039 int d,
1040 DArray<int> newGridDimensions)
1041 {
1042 UTIL_CHECK(d > D);
1043
1044 // Read fields
1045 UnitCell<D> tmpUnitCell;
1046 domain().fieldIo().readFieldsRGrid(inFileName,
1047 tmpFieldsRGrid_,
1048 tmpUnitCell);
1049
1050 // Expand Fields
1051 domain().fieldIo().expandRGridDimension(outFileName,
1052 tmpFieldsRGrid_,
1053 tmpUnitCell,
1054 d, newGridDimensions);
1055 }
1056
1057 /*
1058 * Replicate unit cell a specified number of times in each direction.
1059 */
1060 template <int D>
1061 void System<D>::replicateUnitCell(std::string const & inFileName,
1062 std::string const & outFileName,
1063 IntVec<D> const & replicas)
1064 {
1065 // Read fields
1066 UnitCell<D> tmpUnitCell;
1067 domain().fieldIo().readFieldsRGrid(inFileName, tmpFieldsRGrid_,
1068 tmpUnitCell);
1069
1070 // Replicate fields
1071 domain().fieldIo().replicateUnitCell(outFileName,
1072 tmpFieldsRGrid_,
1073 tmpUnitCell,
1074 replicas);
1075 }
1076
1077 // Thermodynamic Properties
1078
1079 /*
1080 * Compute Helmholtz free energy and pressure.
1081 */
1082 template <int D>
1084 {
1085 if (hasFreeEnergy_) return;
1086
1087 UTIL_CHECK(w_.hasData());
1088 UTIL_CHECK(hasCFields_);
1089
1090 int nm = mixture().nMonomer(); // number of monomer types
1091 int np = mixture().nPolymer(); // number of polymer species
1092 int ns = mixture().nSolvent(); // number of solvent species
1093
1094 // Initialize all free energy contributions to zero
1095 fHelmholtz_ = 0.0;
1096 fIdeal_ = 0.0;
1097 fInter_ = 0.0;
1098 fExt_ = 0.0;
1099
1100 double phi, mu;
1101
1102 // Compute polymer ideal gas contribution to fIdeal_
1103 if (np > 0) {
1104 Polymer<D>* polymerPtr;
1105 double length;
1106 for (int i = 0; i < np; ++i) {
1107 polymerPtr = &mixture_.polymer(i);
1108 phi = polymerPtr->phi();
1109 mu = polymerPtr->mu();
1110 length = polymerPtr->length();
1111 // Recall: mu = ln(phi/q)
1112 if (phi > 1.0E-08) {
1113 fIdeal_ += phi*( mu - 1.0 )/length;
1114 }
1115 }
1116 }
1117
1118 // Compute solvent ideal gas contributions to fIdeal_
1119 if (ns > 0) {
1120 Solvent<D>* solventPtr;
1121 double size;
1122 for (int i = 0; i < ns; ++i) {
1123 solventPtr = &mixture_.solvent(i);
1124 phi = solventPtr->phi();
1125 mu = solventPtr->mu();
1126 size = solventPtr->size();
1127 if (phi > 1.0E-08) {
1128 fIdeal_ += phi*( mu - 1.0 )/size;
1129 }
1130 }
1131 }
1132
1133 // Volume integrals with a mask: If the system has a mask, then the
1134 // volume that should be used in calculating free energy/pressure
1135 // is the volume available to the material, not the total unit cell
1136 // volume. We thus divide all terms that involve integrating over
1137 // the unit cell volume by quantity mask().phiTot(), which is the
1138 // volume fraction of the unit cell that is occupied by material.
1139 // This properly scales them to the correct value. fExt_, fInter_,
1140 // and the Legendre transform component of fIdeal_ all require
1141 // this scaling. If no mask is present, mask.phiTot() = 1 and no
1142 // scaling occurs.
1143 const double phiTot = mask().phiTot();
1144
1145 // Compute Legendre transform subtraction from fIdeal_
1146 double temp = 0.0;
1147 if (w_.isSymmetric()) {
1148 // Use expansion in symmetry-adapted orthonormal basis
1149 const int nBasis = domain_.basis().nBasis();
1150 for (int i = 0; i < nm; ++i) {
1151 for (int k = 0; k < nBasis; ++k) {
1152 temp -= w_.basis(i)[k] * c_.basis(i)[k];
1153 }
1154 }
1155 } else {
1156 // Use summation over grid points
1157 const int meshSize = domain().mesh().size();
1158 for (int i = 0; i < nm; ++i) {
1159 for (int k = 0; k < meshSize; ++k) {
1160 temp -= w_.rgrid(i)[k] * c_.rgrid(i)[k];
1161 }
1162 }
1163 temp /= double(meshSize);
1164 }
1165 temp /= phiTot;
1166 fIdeal_ += temp;
1167 fHelmholtz_ += fIdeal_;
1168
1169 // Compute contribution from external fields, if they exist
1170 if (hasExternalFields()) {
1171 if (w_.isSymmetric()) {
1172 // Use expansion in symmetry-adapted orthonormal basis
1173 const int nBasis = domain_.basis().nBasis();
1174 for (int i = 0; i < nm; ++i) {
1175 for (int k = 0; k < nBasis; ++k) {
1176 fExt_ += h_.basis(i)[k] * c_.basis(i)[k];
1177 }
1178 }
1179 } else {
1180 // Use summation over grid points
1181 const int meshSize = domain().mesh().size();
1182 for (int i = 0; i < nm; ++i) {
1183 for (int k = 0; k < meshSize; ++k) {
1184 fExt_ += h_.rgrid(i)[k] * c_.rgrid(i)[k];
1185 }
1186 }
1187 fExt_ /= double(meshSize);
1188 }
1189 fExt_ /= phiTot;
1190 fHelmholtz_ += fExt_;
1191 }
1192
1193 // Compute excess interaction free energy [ phi^{T}*chi*phi/2 ]
1194 if (w_.isSymmetric()) {
1195 const int nBasis = domain_.basis().nBasis();
1196 for (int i = 0; i < nm; ++i) {
1197 for (int j = i; j < nm; ++j) {
1198 const double chi = interaction().chi(i,j);
1199 if (std::abs(chi) > 1.0E-9) {
1200 double temp = 0.0;
1201 for (int k = 0; k < nBasis; ++k) {
1202 temp += c_.basis(i)[k] * c_.basis(j)[k];
1203 }
1204 if (i == j) {
1205 fInter_ += 0.5*chi*temp;
1206 } else {
1207 fInter_ += chi*temp;
1208 }
1209 }
1210 }
1211 }
1212 } else {
1213 const int meshSize = domain().mesh().size();
1214 for (int i = 0; i < nm; ++i) {
1215 for (int j = i; j < nm; ++j) {
1216 const double chi = interaction().chi(i,j);
1217 if (std::abs(chi) > 1.0E-9) {
1218 double temp = 0.0;
1219 for (int k = 0; k < meshSize; ++k) {
1220 temp += c_.rgrid(i)[k] * c_.rgrid(j)[k];
1221 }
1222 if (i == j) {
1223 fInter_ += 0.5*chi*temp;
1224 } else {
1225 fInter_ += chi*temp;
1226 }
1227 }
1228 }
1229 }
1230 fInter_ /= double(meshSize);
1231 }
1232 fInter_ /= phiTot;
1233 fHelmholtz_ += fInter_;
1234
1235 // Initialize pressure (-1 x grand-canonical free energy / monomer)
1236 pressure_ = -fHelmholtz_;
1237
1238 // Polymer chemical potential corrections to pressure
1239 if (np > 0) {
1240 Polymer<D>* polymerPtr;
1241 double length;
1242 for (int i = 0; i < np; ++i) {
1243 polymerPtr = &mixture_.polymer(i);
1244 phi = polymerPtr->phi();
1245 mu = polymerPtr->mu();
1246 length = polymerPtr->length();
1247 if (phi > 1.0E-08) {
1248 pressure_ += mu * phi / length;
1249 }
1250 }
1251 }
1252
1253 // Solvent corrections to pressure
1254 if (ns > 0) {
1255 Solvent<D>* solventPtr;
1256 double size;
1257 for (int i = 0; i < ns; ++i) {
1258 solventPtr = &mixture_.solvent(i);
1259 phi = solventPtr->phi();
1260 mu = solventPtr->mu();
1261 size = solventPtr->size();
1262 if (phi > 1.0E-08) {
1263 pressure_ += mu * phi / size;
1264 }
1265 }
1266 }
1267
1268 hasFreeEnergy_ = true;
1269 }
1270
1271 // Output Operations
1272
1273 /*
1274 * Write timer values to output stream (computational cost).
1275 */
1276 template <int D>
1277 void System<D>::writeTimers(std::ostream& out)
1278 {
1279 if (iteratorPtr_) {
1280 iterator().outputTimers(Log::file());
1281 iterator().outputTimers(out);
1282 }
1283 if (hasSimulator()){
1284 simulator().outputTimers(Log::file());
1285 simulator().outputTimers(out);
1286 }
1287 }
1288
1289 /*
1290 * Clear state of all timers.
1291 */
1292 template <int D>
1294 {
1295 if (iteratorPtr_) {
1296 iterator().clearTimers();
1297 }
1298 if (hasSimulator()){
1299 simulator().clearTimers();
1300 }
1301 }
1302
1303 /*
1304 * Write parameter file for SCFT, omitting any sweep block.
1305 */
1306 template <int D>
1307 void System<D>::writeParamNoSweep(std::ostream& out) const
1308 {
1309 out << "System{" << std::endl;
1310 mixture_.writeParam(out);
1311 interaction().writeParam(out);
1312 domain_.writeParam(out);
1313 if (iteratorPtr_) {
1314 iterator().writeParam(out);
1315 }
1316 out << "}" << std::endl;
1317 }
1318
1319 /*
1320 * Write thermodynamic properties to file.
1321 */
1322 template <int D>
1323 void System<D>::writeThermo(std::ostream& out)
1324 {
1325 if (!hasFreeEnergy_) {
1326 computeFreeEnergy();
1327 }
1328
1329 out << std::endl;
1330 out << "fHelmholtz " << Dbl(fHelmholtz(), 18, 11) << std::endl;
1331 out << "pressure " << Dbl(pressure(), 18, 11) << std::endl;
1332 out << std::endl;
1333 out << "fIdeal " << Dbl(fIdeal_, 18, 11) << std::endl;
1334 out << "fInter " << Dbl(fInter_, 18, 11) << std::endl;
1335 if (hasExternalFields()) {
1336 out << "fExt " << Dbl(fExt_, 18, 11) << std::endl;
1337 }
1338 out << std::endl;
1339
1340 int np = mixture().nPolymer();
1341 int ns = mixture().nSolvent();
1342
1343 if (np > 0) {
1344 out << "polymers:" << std::endl;
1345 out << " "
1346 << " phi "
1347 << " mu "
1348 << std::endl;
1349 for (int i = 0; i < np; ++i) {
1350 out << Int(i, 5)
1351 << " " << Dbl(mixture_.polymer(i).phi(),18, 11)
1352 << " " << Dbl(mixture_.polymer(i).mu(), 18, 11)
1353 << std::endl;
1354 }
1355 out << std::endl;
1356 }
1357
1358 if (ns > 0) {
1359 out << "solvents:" << std::endl;
1360 out << " "
1361 << " phi "
1362 << " mu "
1363 << std::endl;
1364 for (int i = 0; i < ns; ++i) {
1365 out << Int(i, 5)
1366 << " " << Dbl(mixture_.solvent(i).phi(),18, 11)
1367 << " " << Dbl(mixture_.solvent(i).mu(), 18, 11)
1368 << std::endl;
1369 }
1370 out << std::endl;
1371 }
1372
1373 out << "cellParams:" << std::endl;
1374 for (int i = 0; i < domain().unitCell().nParameter(); ++i) {
1375 out << Int(i, 5)
1376 << " "
1377 << Dbl(domain().unitCell().parameter(i), 18, 11)
1378 << std::endl;
1379 }
1380 out << std::endl;
1381 }
1382
1383 /*
1384 * Write stress properties to file.
1385 */
1386 template <int D>
1387 void System<D>::writeStress(std::ostream& out)
1388 {
1389 out << "stress:" << std::endl;
1390 for (int i = 0; i < domain().unitCell().nParameter(); ++i) {
1391 out << Int(i, 5)
1392 << " "
1393 << Dbl(mixture_.stress(i), 18, 11)
1394 << std::endl;
1395 }
1396 out << std::endl;
1397 }
1398
1399 /*
1400 * Write w-fields in symmetry-adapted basis format.
1401 */
1402 template <int D>
1403 void System<D>::writeWBasis(std::string const & filename) const
1404 {
1405 UTIL_CHECK(domain_.basis().isInitialized());
1406 UTIL_CHECK(isAllocatedBasis_);
1407 UTIL_CHECK(w_.hasData());
1408 UTIL_CHECK(w_.isSymmetric());
1409 domain_.fieldIo().writeFieldsBasis(filename, w_.basis(),
1410 domain().unitCell());
1411 }
1412
1413 /*
1414 * Write w-fields in real space grid file format.
1415 */
1416 template <int D>
1417 void System<D>::writeWRGrid(const std::string & filename) const
1418 {
1419 UTIL_CHECK(isAllocatedGrid_);
1420 UTIL_CHECK(w_.hasData());
1421 domain_.fieldIo().writeFieldsRGrid(filename, w_.rgrid(),
1422 domain().unitCell(),
1423 w_.isSymmetric());
1424 }
1425
1426 /*
1427 * Write concentration fields in symmetry-adapted basis format.
1428 */
1429 template <int D>
1430 void System<D>::writeCBasis(const std::string & filename) const
1431 {
1432 UTIL_CHECK(domain_.basis().isInitialized());
1433 UTIL_CHECK(isAllocatedBasis_);
1434 UTIL_CHECK(hasCFields_);
1435 UTIL_CHECK(w_.isSymmetric());
1436 domain_.fieldIo().writeFieldsBasis(filename, c_.basis(),
1437 domain().unitCell());
1438 }
1439
1440 /*
1441 * Write concentration fields in real space (r-grid) format.
1442 */
1443 template <int D>
1444 void System<D>::writeCRGrid(const std::string & filename) const
1445 {
1446 UTIL_CHECK(isAllocatedGrid_);
1447 UTIL_CHECK(hasCFields_);
1448 domain_.fieldIo().writeFieldsRGrid(filename, c_.rgrid(),
1449 domain().unitCell(),
1450 w_.isSymmetric());
1451 }
1452
1453 /*
1454 * Write all concentration fields in real space (r-grid) format, for
1455 * each block and solvent species, rather than for each monomer type.
1456 */
1457 template <int D>
1458 void System<D>::writeBlockCRGrid(const std::string & filename) const
1459 {
1460 UTIL_CHECK(isAllocatedGrid_);
1461 UTIL_CHECK(hasCFields_);
1462
1463 // Construct array to hold block and solvent c-field data
1464 DArray< RField<D> > blockCFields;
1465
1466 // Get c-field data from the Mixture
1467 // Note: blockCFields container is allocated within this function
1468 mixture_.createBlockCRGrid(blockCFields);
1469
1470 // Write block and solvent c-field data to file
1471 domain().fieldIo().writeFieldsRGrid(filename, blockCFields,
1472 domain().unitCell(),
1473 w_.isSymmetric());
1474 }
1475
1476 /*
1477 * Write the last time slice of the propagator in r-grid format.
1478 */
1479 template <int D>
1480 void System<D>::writeQSlice(const std::string & filename,
1481 int polymerId, int blockId,
1482 int directionId, int segmentId)
1483 const
1484 {
1485 UTIL_CHECK(polymerId >= 0);
1486 UTIL_CHECK(polymerId < mixture().nPolymer());
1487 Polymer<D> const& polymer = mixture_.polymer(polymerId);
1488 UTIL_CHECK(blockId >= 0);
1489 UTIL_CHECK(blockId < polymer.nBlock());
1490 UTIL_CHECK(directionId >= 0);
1491 UTIL_CHECK(directionId <= 1);
1492 Propagator<D> const &
1493 propagator = polymer.propagator(blockId, directionId);
1494 RField<D> const& field = propagator.q(segmentId);
1495 domain().fieldIo().writeFieldRGrid(filename, field,
1496 domain().unitCell(),
1497 w_.isSymmetric());
1498 }
1499
1500 /*
1501 * Write the last (tail) slice of the propagator in r-grid format.
1502 */
1503 template <int D>
1504 void System<D>::writeQTail(const std::string & filename,
1505 int polymerId, int blockId, int directionId)
1506 const
1507 {
1508 UTIL_CHECK(polymerId >= 0);
1509 UTIL_CHECK(polymerId < mixture().nPolymer());
1510 Polymer<D> const& polymer = mixture_.polymer(polymerId);
1511 UTIL_CHECK(blockId >= 0);
1512 UTIL_CHECK(blockId < polymer.nBlock());
1513 UTIL_CHECK(directionId >= 0);
1514 UTIL_CHECK(directionId <= 1);
1515 RField<D> const&
1516 field = polymer.propagator(blockId, directionId).tail();
1517 domain().fieldIo().writeFieldRGrid(filename, field,
1518 domain().unitCell(),
1519 w_.isSymmetric());
1520 }
1521
1522 /*
1523 * Write the propagator for a block and direction.
1524 */
1525 template <int D>
1526 void System<D>::writeQ(const std::string & filename,
1527 int polymerId, int blockId, int directionId)
1528 const
1529 {
1530 UTIL_CHECK(polymerId >= 0);
1531 UTIL_CHECK(polymerId < mixture().nPolymer());
1532 Polymer<D> const& polymer = mixture_.polymer(polymerId);
1533 UTIL_CHECK(blockId >= 0);
1534 UTIL_CHECK(blockId < polymer.nBlock());
1535 UTIL_CHECK(directionId >= 0);
1536 UTIL_CHECK(directionId <= 1);
1537 Propagator<D> const& propagator
1538 = polymer.propagator(blockId, directionId);
1539 int ns = propagator.ns();
1540
1541 // Open file
1542 std::ofstream file;
1543 fileMaster_.openOutputFile(filename, file);
1544
1545 // Write header
1546 domain().fieldIo().writeFieldHeader(file, 1, domain().unitCell(),
1547 w_.isSymmetric());
1548 file << "mesh " << std::endl
1549 << " " << domain().mesh().dimensions() << std::endl
1550 << "nslice" << std::endl
1551 << " " << ns << std::endl;
1552
1553 // Write data
1554 bool hasHeader = false;
1555 for (int i = 0; i < ns; ++i) {
1556 file << "slice " << i << std::endl;
1557 domain().fieldIo().writeFieldRGrid(file, propagator.q(i),
1558 domain().unitCell(),
1559 hasHeader);
1560 }
1561 file.close();
1562 }
1563
1564 /*
1565 * Write propagators for all blocks of all polymers to files.
1566 */
1567 template <int D>
1568 void System<D>::writeQAll(std::string const & basename)
1569 {
1570 std::string filename;
1571 int np, nb, ip, ib, id;
1572 np = mixture().nPolymer();
1573 for (ip = 0; ip < np; ++ip) {
1574 nb = mixture_.polymer(ip).nBlock();
1575 for (ib = 0; ib < nb; ++ib) {
1576 for (id = 0; id < 2; ++id) {
1577 filename = basename;
1578 filename += "_";
1579 filename += toString(ip);
1580 filename += "_";
1581 filename += toString(ib);
1582 filename += "_";
1583 filename += toString(id);
1584 filename += ".rq";
1585 writeQ(filename, ip, ib, id);
1586 }
1587 }
1588 }
1589 }
1590
1591 /*
1592 * Write description of symmetry-adapted stars and basis to file.
1593 */
1594 template <int D>
1595 void System<D>::writeStars(std::string const & filename) const
1596 {
1597 UTIL_CHECK(domain_.basis().isInitialized());
1598 std::ofstream file;
1599 fileMaster_.openOutputFile(filename, file);
1600 bool isSymmetric = true;
1601 domain().fieldIo().writeFieldHeader(file, mixture().nMonomer(),
1602 domain().unitCell(),
1603 isSymmetric);
1604 domain_.basis().outputStars(file);
1605 file.close();
1606 }
1607
1608 /*
1609 * Write a list of waves and associated stars to file.
1610 */
1611 template <int D>
1612 void System<D>::writeWaves(std::string const & filename) const
1613 {
1614 UTIL_CHECK(domain_.hasGroup());
1615 UTIL_CHECK(domain_.basis().isInitialized());
1616 std::ofstream file;
1617 fileMaster_.openOutputFile(filename, file);
1618 bool isSymmetric = true;
1619 domain().fieldIo().writeFieldHeader(file, mixture().nMonomer(),
1620 domain().unitCell(),
1621 isSymmetric);
1622 domain_.basis().outputWaves(file);
1623 file.close();
1624 }
1625
1626 /*
1627 * Write all elements of the space group to a file.
1628 */
1629 template <int D>
1630 void System<D>::writeGroup(std::string const & filename) const
1631 {
1632 UTIL_CHECK(domain_.hasGroup());
1633 Pscf::Prdc::writeGroup(filename, domain_.group());
1634 }
1635
1636 // Field format conversion functions
1637
1638 /*
1639 * Convert fields from symmetry-adapted basis to real-space grid format.
1640 */
1641 template <int D>
1642 void System<D>::basisToRGrid(std::string const & inFileName,
1643 std::string const & outFileName)
1644 {
1645 UTIL_CHECK(domain_.hasGroup());
1646
1647 // If basis fields are not allocated, peek at field file header to
1648 // get unit cell parameters, initialize basis and allocate fields.
1649 if (!isAllocatedBasis_) {
1650 readFieldHeader(inFileName);
1651 allocateFieldsBasis();
1652 }
1653
1654 // Read, convert, and write fields
1655 UnitCell<D> tmpUnitCell;
1656 FieldIo<D> const & fieldIo = domain().fieldIo();
1657 fieldIo.readFieldsBasis(inFileName, tmpFieldsBasis_, tmpUnitCell);
1658 fieldIo.convertBasisToRGrid(tmpFieldsBasis_, tmpFieldsRGrid_);
1659 fieldIo.writeFieldsRGrid(outFileName, tmpFieldsRGrid_,
1660 tmpUnitCell);
1661 }
1662
1663 /*
1664 * Convert fields from real-space grid to symmetry-adapted basis format.
1665 */
1666 template <int D>
1667 void System<D>::rGridToBasis(std::string const & inFileName,
1668 std::string const & outFileName)
1669 {
1670 UTIL_CHECK(domain_.hasGroup());
1671
1672 // If basis fields are not allocated, peek at field file header to
1673 // get unit cell parameters, initialize basis and allocate fields.
1674 if (!isAllocatedBasis_) {
1675 readFieldHeader(inFileName);
1676 allocateFieldsBasis();
1677 }
1678
1679 // Read, convert and write fields
1680 UnitCell<D> tmpUnitCell;
1681 FieldIo<D> const & fieldIo = domain().fieldIo();
1682 fieldIo.readFieldsRGrid(inFileName, tmpFieldsRGrid_, tmpUnitCell);
1683 fieldIo.convertRGridToBasis(tmpFieldsRGrid_, tmpFieldsBasis_);
1684 fieldIo.writeFieldsBasis(outFileName, tmpFieldsBasis_,
1685 tmpUnitCell);
1686 }
1687
1688 /*
1689 * Convert fields from Fourier (k-grid) to real-space (r-grid) format.
1690 */
1691 template <int D>
1692 void System<D>::kGridToRGrid(std::string const & inFileName,
1693 std::string const & outFileName)
1694 {
1695 // If basis fields are not allocated, peek at field file header to
1696 // get unit cell parameters, initialize basis and allocate fields.
1697 if (domain_.hasGroup() && !isAllocatedBasis_) {
1698 readFieldHeader(inFileName);
1699 allocateFieldsBasis();
1700 }
1701
1702 // Read, convert and write fields
1703 UnitCell<D> tmpUnitCell;
1704 FieldIo<D> const & fieldIo = domain().fieldIo();
1705 fieldIo.readFieldsKGrid(inFileName, tmpFieldsKGrid_, tmpUnitCell);
1706 for (int i = 0; i < mixture().nMonomer(); ++i) {
1707 domain().fft().inverseTransformUnsafe(tmpFieldsKGrid_[i],
1708 tmpFieldsRGrid_[i]);
1709 }
1710 fieldIo.writeFieldsRGrid(outFileName, tmpFieldsRGrid_,
1711 tmpUnitCell);
1712 }
1713
1714 /*
1715 * Convert fields from real-space (r-grid) to Fourier (k-grid) format.
1716 */
1717 template <int D>
1718 void System<D>::rGridToKGrid(std::string const & inFileName,
1719 std::string const & outFileName)
1720 {
1721 // If basis fields are not allocated, peek at field file header to
1722 // get unit cell parameters, initialize basis and allocate fields.
1723 if (domain_.hasGroup() && !isAllocatedBasis_) {
1724 readFieldHeader(inFileName);
1725 allocateFieldsBasis();
1726 }
1727
1728 // Read, convert and write fields
1729 UnitCell<D> tmpUnitCell;
1730 FieldIo<D> const & fieldIo = domain().fieldIo();
1731 fieldIo.readFieldsRGrid(inFileName, tmpFieldsRGrid_,
1732 tmpUnitCell);
1733 for (int i = 0; i < mixture().nMonomer(); ++i) {
1734 domain().fft().forwardTransform(tmpFieldsRGrid_[i],
1735 tmpFieldsKGrid_[i]);
1736 }
1737 fieldIo.writeFieldsKGrid(outFileName, tmpFieldsKGrid_,
1738 tmpUnitCell);
1739 }
1740
1741 /*
1742 * Convert fields from Fourier (k-grid) to symmetry-adapted basis format.
1743 */
1744 template <int D>
1745 void System<D>::kGridToBasis(std::string const & inFileName,
1746 std::string const & outFileName)
1747 {
1748 UTIL_CHECK(domain_.hasGroup());
1749
1750 // If basis fields are not allocated, peek at field file header to
1751 // get unit cell parameters, initialize basis and allocate fields.
1752 if (!isAllocatedBasis_) {
1753 readFieldHeader(inFileName);
1754 allocateFieldsBasis();
1755 }
1756
1757 // Read, convert and write fields
1758 UnitCell<D> tmpUnitCell;
1759 domain_.fieldIo().readFieldsKGrid(inFileName, tmpFieldsKGrid_,
1760 tmpUnitCell);
1761 domain_.fieldIo().convertKGridToBasis(tmpFieldsKGrid_,
1762 tmpFieldsBasis_);
1763 domain_.fieldIo().writeFieldsBasis(outFileName,
1764 tmpFieldsBasis_, tmpUnitCell);
1765 }
1766
1767 /*
1768 * Convert fields from symmetry-adapted basis to Fourier (k-grid) format.
1769 */
1770 template <int D>
1771 void System<D>::basisToKGrid(std::string const & inFileName,
1772 std::string const & outFileName)
1773 {
1774 UTIL_CHECK(domain_.hasGroup());
1775
1776 // If basis fields are not allocated, peek at field file header to
1777 // get unit cell parameters, initialize basis and allocate fields.
1778 if (!isAllocatedBasis_) {
1779 readFieldHeader(inFileName);
1780 allocateFieldsBasis();
1781 }
1782
1783 // Read, convert and write fields
1784 UnitCell<D> tmpUnitCell;
1785 domain_.fieldIo().readFieldsBasis(inFileName,
1786 tmpFieldsBasis_, tmpUnitCell);
1787 domain_.fieldIo().convertBasisToKGrid(tmpFieldsBasis_,
1788 tmpFieldsKGrid_);
1789 domain_.fieldIo().writeFieldsKGrid(outFileName,
1790 tmpFieldsKGrid_, tmpUnitCell);
1791 }
1792
1793 /*
1794 * Check if r-grid fields have declared space group symmetry.
1795 */
1796 template <int D>
1797 bool
1798 System<D>::checkRGridFieldSymmetry(std::string const & inFileName,
1799 double epsilon)
1800 {
1801 UTIL_CHECK(domain_.hasGroup());
1802
1803 // If basis fields are not allocated, peek at field file header to
1804 // get unit cell parameters, initialize basis and allocate fields.
1805 if (!isAllocatedBasis_) {
1806 readFieldHeader(inFileName);
1807 allocateFieldsBasis();
1808 }
1809
1810 // Read fields
1811 UnitCell<D> tmpUnitCell;
1812 domain_.fieldIo().readFieldsRGrid(inFileName,
1813 tmpFieldsRGrid_, tmpUnitCell);
1814
1815 // Check symmetry for all fields
1816 for (int i = 0; i < mixture().nMonomer(); ++i) {
1817 bool symmetric;
1818 symmetric = domain_.fieldIo().hasSymmetry(tmpFieldsRGrid_[i],
1819 epsilon);
1820 if (!symmetric) {
1821 return false;
1822 }
1823 }
1824 return true;
1825
1826 }
1827
1828 /*
1829 * Rescale a field by a constant, write rescaled field to a file.
1830 */
1831 template <int D>
1832 void System<D>::scaleFieldsBasis(std::string const & inFileName,
1833 std::string const & outFileName,
1834 double factor)
1835 {
1836 // If basis fields are not allocated, peek at field file header to
1837 // get unit cell parameters, initialize basis and allocate fields.
1838 if (!isAllocatedBasis_) {
1839 readFieldHeader(inFileName);
1840 allocateFieldsBasis();
1841 }
1842
1843 UnitCell<D> tmpUnitCell;
1844 FieldIo<D> const & fieldIo = domain().fieldIo();
1845 fieldIo.readFieldsBasis(inFileName, tmpFieldsBasis_, tmpUnitCell);
1846 fieldIo.scaleFieldsBasis(tmpFieldsBasis_, factor);
1847 fieldIo.writeFieldsBasis(outFileName, tmpFieldsBasis_,
1848 tmpUnitCell);
1849 }
1850
1851 /*
1852 * Rescale a field by a constant, write rescaled field to a file.
1853 */
1854 template <int D>
1855 void System<D>::scaleFieldsRGrid(std::string const & inFileName,
1856 std::string const & outFileName,
1857 double factor) const
1858 {
1859 UnitCell<D> tmpUnitCell;
1860 FieldIo<D> const & fieldIo = domain().fieldIo();
1861 fieldIo.readFieldsRGrid(inFileName, tmpFieldsRGrid_,
1862 tmpUnitCell);
1863 fieldIo.scaleFieldsRGrid(tmpFieldsRGrid_, factor);
1864 fieldIo.writeFieldsRGrid(outFileName, tmpFieldsRGrid_,
1865 tmpUnitCell);
1866 }
1867
1868 /*
1869 * Compare two fields in basis format, write report to Log file.
1870 */
1871 template <int D>
1873 const DArray< DArray<double> > field2)
1874 {
1875 UTIL_CHECK(domain_.hasGroup());
1876
1877 BFieldComparison comparison(1);
1878 comparison.compare(field1,field2);
1879
1880 Log::file() << "\n Basis expansion field comparison results"
1881 << std::endl;
1882 Log::file() << " Maximum Absolute Difference: "
1883 << comparison.maxDiff() << std::endl;
1884 Log::file() << " Root-Mean-Square Difference: "
1885 << comparison.rmsDiff() << "\n" << std::endl;
1886 }
1887
1888 /*
1889 * Compare two fields in r-grid format, output report to Log file.
1890 */
1891 template <int D>
1893 const DArray< RField<D> > field2)
1894 {
1895 RFieldComparison<D> comparison;
1896 comparison.compare(field1, field2);
1897
1898 Log::file() << "\n Real-space field comparison results"
1899 << std::endl;
1900 Log::file() << " Maximum Absolute Difference: "
1901 << comparison.maxDiff() << std::endl;
1902 Log::file() << " Root-Mean-Square Difference: "
1903 << comparison.rmsDiff() << "\n" << std::endl;
1904 }
1905
1906 // Private member functions
1907
1908 /*
1909 * Allocate memory for fields in grid format.
1910 */
1911 template <int D>
1913 {
1914 // Preconditions
1915 UTIL_CHECK(hasMixture_);
1916 int nMonomer = mixture().nMonomer();
1917 UTIL_CHECK(nMonomer > 0);
1918 UTIL_CHECK(domain().mesh().size() > 0);
1919 UTIL_CHECK(!isAllocatedGrid_);
1920
1921 // Alias for mesh dimensions
1922 IntVec<D> const & dimensions = domain().mesh().dimensions();
1923
1924 // Allocate w (chemical potential) fields
1925 w_.setNMonomer(nMonomer);
1926 w_.allocateRGrid(dimensions);
1927
1928 // Allocate c (monomer concentration) fields
1929 c_.setNMonomer(nMonomer);
1930 c_.allocateRGrid(dimensions);
1931
1932 h_.setNMonomer(nMonomer);
1933
1934 // Allocate work space field arrays
1935 tmpFieldsRGrid_.allocate(nMonomer);
1936 tmpFieldsKGrid_.allocate(nMonomer);
1937 for (int i = 0; i < nMonomer; ++i) {
1938 tmpFieldsRGrid_[i].allocate(dimensions);
1939 tmpFieldsKGrid_[i].allocate(dimensions);
1940 }
1941
1942 isAllocatedGrid_ = true;
1943 }
1944
1945 /*
1946 * Allocate memory for fields in basis format.
1947 */
1948 template <int D>
1949 void System<D>::allocateFieldsBasis()
1950 {
1951 // Preconditions and constants
1952 UTIL_CHECK(hasMixture_);
1953 const int nMonomer = mixture().nMonomer();
1954 UTIL_CHECK(nMonomer > 0);
1955 UTIL_CHECK(isAllocatedGrid_);
1956 UTIL_CHECK(!isAllocatedBasis_);
1957 UTIL_CHECK(domain_.hasGroup());
1958 UTIL_CHECK(domain_.basis().isInitialized());
1959 const int nBasis = domain_.basis().nBasis();
1960 UTIL_CHECK(nBasis > 0);
1961
1962 w_.allocateBasis(nBasis);
1963 c_.allocateBasis(nBasis);
1964
1965 // Temporary work space
1966 tmpFieldsBasis_.allocate(nMonomer);
1967 for (int i = 0; i < nMonomer; ++i) {
1968 tmpFieldsBasis_[i].allocate(nBasis);
1969 }
1970
1971 isAllocatedBasis_ = true;
1972 }
1973
1974 /*
1975 * Peek at field file header, initialize unit cell parameters and basis.
1976 */
1977 template <int D>
1978 void System<D>::readFieldHeader(std::string filename)
1979 {
1980 UTIL_CHECK(hasMixture_);
1981 UTIL_CHECK(mixture().nMonomer() > 0);
1982
1983 // Open field file
1984 std::ifstream file;
1985 fileMaster_.openInputFile(filename, file);
1986
1987 // Read field file header, and initialize basis if needed
1988 int nMonomer;
1989 bool isSymmetric;
1990 domain_.fieldIo().readFieldHeader(file, nMonomer,
1991 domain_.unitCell(), isSymmetric);
1992 // Function FieldIo::readFieldHeader initializes a basis if
1993 // hasGroup is true and the header contains a space group
1994
1995 // Close field file
1996 file.close();
1997
1998 // Postconditions
1999 UTIL_CHECK(mixture().nMonomer() == nMonomer);
2000 UTIL_CHECK(domain().unitCell().nParameter() > 0);
2001 UTIL_CHECK(domain().unitCell().lattice() != UnitCell<D>::Null);
2002 UTIL_CHECK(domain().unitCell().isInitialized());
2003 if (domain_.hasGroup() && isSymmetric) {
2004 UTIL_CHECK(domain_.basis().isInitialized());
2005 UTIL_CHECK(domain_.basis().nBasis() > 0);
2006 }
2007
2008 }
2009
2010 /*
2011 * Read a filename string and echo to log file (used in readCommands).
2012 */
2013 template <int D>
2014 void System<D>::readEcho(std::istream& in, std::string& string) const
2015 {
2016 in >> string;
2017 if (in.fail()) {
2018 UTIL_THROW("Unable to read string parameter.");
2019 }
2020 Log::file() << " " << Str(string, 20) << std::endl;
2021 }
2022
2023 /*
2024 * Read floating point number, echo to log file (used in readCommands).
2025 */
2026 template <int D>
2027 void System<D>::readEcho(std::istream& in, double& value) const
2028 {
2029 in >> value;
2030 if (in.fail()) {
2031 UTIL_THROW("Unable to read floating point parameter.");
2032 }
2033 Log::file() << " " << Dbl(value, 20) << std::endl;
2034 }
2035
2036 /*
2037 * Initialize the homogeneous_ member object, which describes
2038 * thermodynamics of a homogeneous reference system.
2039 */
2040 template <int D>
2041 void System<D>::initHomogeneous()
2042 {
2043 // Set number of molecular species and monomers
2044 UTIL_CHECK(hasMixture_);
2045 int nm = mixture().nMonomer();
2046 int np = mixture().nPolymer();
2047 int ns = mixture().nSolvent();
2048 UTIL_CHECK(homogeneous_.nMolecule() == np + ns);
2049 UTIL_CHECK(homogeneous_.nMonomer() == nm);
2050
2051 int i; // molecule index
2052 int j; // monomer index
2053
2054 // Loop over polymer molecule species
2055 if (np > 0) {
2056
2057 // Allocate array of clump sizes
2059 s.allocate(nm);
2060
2061 int k; // block or clump index
2062 int nb; // number of blocks
2063 int nc; // number of clumps
2064
2065 // Loop over polymer species
2066 for (i = 0; i < np; ++i) {
2067
2068 // Initial array of clump sizes for this polymer
2069 for (j = 0; j < nm; ++j) {
2070 s[j] = 0.0;
2071 }
2072
2073 // Compute clump sizes for all monomer types.
2074 nb = mixture_.polymer(i).nBlock();
2075 for (k = 0; k < nb; ++k) {
2076 Block<D>& block = mixture_.polymer(i).block(k);
2077 j = block.monomerId();
2078 s[j] += block.length();
2079 }
2080
2081 // Count the number of clumps of nonzero size
2082 nc = 0;
2083 for (j = 0; j < nm; ++j) {
2084 if (s[j] > 1.0E-8) {
2085 ++nc;
2086 }
2087 }
2088 homogeneous_.molecule(i).setNClump(nc);
2089
2090 // Set clump properties for this Homogeneous::Molecule
2091 k = 0; // Clump index
2092 for (j = 0; j < nm; ++j) {
2093 if (s[j] > 1.0E-8) {
2094 homogeneous_.molecule(i).clump(k).setMonomerId(j);
2095 homogeneous_.molecule(i).clump(k).setSize(s[j]);
2096 ++k;
2097 }
2098 }
2099 homogeneous_.molecule(i).computeSize();
2100
2101 }
2102 }
2103
2104 // Add solvent contributions
2105 if (ns > 0) {
2106 double size;
2107 int monomerId;
2108 for (int is = 0; is < ns; ++is) {
2109 i = is + np;
2110 monomerId = mixture_.solvent(is).monomerId();
2111 size = mixture_.solvent(is).size();
2112 homogeneous_.molecule(i).setNClump(1);
2113 homogeneous_.molecule(i).clump(0).setMonomerId(monomerId);
2114 homogeneous_.molecule(i).clump(0).setSize(size);
2115 homogeneous_.molecule(i).computeSize();
2116 }
2117 }
2118 }
2119
2120} // namespace Rpc
2121} // namespace Pscf
2122#endif
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 convertBasisToRGrid(DArray< double > const &in, RFRT &out) const
Convert a field from symmetrized basis to spatial grid (r-grid).
void convertRGridToBasis(RFRT const &in, DArray< double > &out, bool checkSymmetry=true, double epsilon=1.0e-8) const
Convert a field from spatial grid (r-grid) to symmetrized basis.
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 writeFieldsKGrid(std::ostream &out, DArray< RFieldDft< D > > const &fields, UnitCell< D > const &unitCell, bool isSymmetric=true) const override
Write array of RFieldDft objects (k-space fields) to file.
void readFieldsRGrid(std::istream &in, DArray< RField< D > > &fields, UnitCell< D > &unitCell) const override
Read array of RField objects (r-grid fields) from a stream.
void readFieldsKGrid(std::istream &in, DArray< RFieldDft< D > > &fields, UnitCell< D > &unitCell) const override
Read array of RFieldDft objects (k-space fields) from a stream.
void writeFieldsRGrid(std::ostream &out, DArray< RField< D > > const &fields, UnitCell< D > const &unitCell, bool writeHeader=true, bool isSymmetric=true, bool writeMeshSize=true) const override
Write array of RField objects (fields on r-space grid) to a stream.
Factory for subclasses of Iterator.
Definition rpc/System.h:35
Descriptor and solver for one polymer species.
MDE solver for one direction of one block.
const QField & q(int i) const
Return q-field at specified step.
int ns() const
Number of values of s (or slices), including head and tail.
const QField & tail() const
Return q-field at the end of the block.
Factory for subclasses of Simulator.
Definition rpc/System.h:39
Solver and descriptor for a solvent species.
Default Factory for subclasses of Sweep.
Definition rpc/System.h:37
Main class for SCFT or PS-FTS simulation of one system.
Definition rpc/System.h:100
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.
int iterate(bool isContinuation=false)
Iteratively solve a SCFT problem.
void writeGroup(std::string const &filename) const
Output all elements of the space group.
void scaleFieldsBasis(const std::string &inFileName, const std::string &outFileName, double factor)
Multiply all components of an array of basis fields by a scalar.
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 readWBasis(const std::string &filename)
Read chemical potential fields in symmetry adapted basis format.
void writeWBasis(const std::string &filename) const
Write chemical potential fields in symmetrized basis 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 writeWRGrid(const std::string &filename) const
Write chemical potential fields in real space grid (r-grid) format.
void setOptions(int argc, char **argv)
Process command line options.
void readCommands()
Read and process commands from the default command file.
void simulate(int nStep)
Perform a field theoretic simulation (PS-FTS).
void writeQAll(std::string const &basename)
Write all propagators of all blocks, each to a separate file.
bool checkRGridFieldSymmetry(const std::string &inFileName, double epsilon=1.0E-8)
Check if r-grid fields have the declared space group symmetry.
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 writeTimers(std::ostream &out)
Write timer information to an output stream.
void compare(const DArray< DArray< double > > field1, const DArray< DArray< double > > field2)
Compare arrays of fields in basis format, output a report.
void writeCBasis(const std::string &filename) const
Write concentration fields in symmetrized basis format.
void writeThermo(std::ostream &out)
Write thermodynamic properties to a file.
void setUnitCell(UnitCell< D > const &unitCell)
Set parameters of the associated unit cell.
~System()
Destructor.
void rGridToKGrid(const std::string &inFileName, const std::string &outFileName)
Convert fields from real-space (r-grid) to Fourier (k-grid) format.
void basisToKGrid(const std::string &inFileName, const std::string &outFileName)
Convert fields from symmetrized basis to Fourier (k-grid) format.
void computeFreeEnergy()
Compute free energy density and pressure for current fields.
void setWBasis(DArray< DArray< double > > const &fields)
Set chemical potential fields, in symmetry-adapted basis format.
System()
Constructor.
void clearTimers()
Clear timers.
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 writeQ(std::string const &filename, int polymerId, int blockId, int directionId) const
Write one propagator for one block, in r-grid format.
void writeBlockCRGrid(const std::string &filename) const
Write c-fields for all blocks and solvents in r-grid format.
virtual void readParameters(std::istream &in)
Read body of parameter block (without opening and closing lines).
void writeParamNoSweep(std::ostream &out) const
Write parameter file to an ostream, omitting any sweep block.
void sweep()
Sweep in parameter space, solving an SCF problem at each point.
void kGridToBasis(const std::string &inFileName, const std::string &outFileName)
Convert fields from Fourier (k-grid) to symmetrized basis format.
void readWRGrid(const std::string &filename)
Read chemical potential fields in real space grid (r-grid) format.
void kGridToRGrid(const std::string &inFileName, const std::string &outFileName)
Convert fields from Fourier (k-grid) to real-space (r-grid) format.
void estimateWfromC(const std::string &filename)
Construct trial w-fields from c-fields in symmetry-adapted form.
void writeCRGrid(const std::string &filename) const
Write concentration fields in real space grid (r-grid) format.
Domain< D > const & domain() const
Get Domain by const reference.
void writeQSlice(std::string const &filename, int polymerId, int blockId, int directionId, int segmentId) const
Write slice of a propagator at fixed s in r-grid format.
void basisToRGrid(const std::string &inFileName, const std::string &outFileName)
Convert a field from symmetrized basis format to r-grid format.
void writeStars(std::string const &filename) const
Output information about stars and symmetrized basis functions.
void readParam()
Read input parameters from default param file.
void writeStress(std::ostream &out)
Write stress properties to a file.
void rGridToBasis(const std::string &inFileName, const std::string &outFileName)
Convert a field from real-space grid to symmetrized basis format.
void writeWaves(std::string const &filename) const
Output information about waves.
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
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.
Fields and FFTs for periodic boundary conditions (CPU)
Definition CField.cpp: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.