PSCF v1.2
Basis.tpp
1#ifndef PRDC_BASIS_TPP
2#define PRDC_BASIS_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 "Basis.h"
12#include "TWave.h"
13#include "groupFile.h"
14#include <prdc/crystal/UnitCell.h>
15#include <prdc/crystal/SpaceGroup.h>
16#include <prdc/crystal/shiftToMinimum.h>
17#include <pscf/mesh/Mesh.h>
18#include <pscf/mesh/MeshIterator.h>
19#include <algorithm>
20#include <vector>
21#include <set>
22#include <fstream>
23
24namespace Pscf {
25namespace Prdc {
26
27 /*
28 * Constructor.
29 */
30 template <int D>
32 : waves_(),
33 stars_(),
34 waveIds_(),
35 starIds_(),
36 nWave_(0),
37 nBasisWave_(0),
38 nStar_(0),
39 nBasis_(0),
40 unitCellPtr_(0),
41 meshPtr_(0),
42 isInitialized_(false)
43 {}
44
45 /*
46 * Destructor.
47 */
48 template <int D>
51
52 /*
53 * Construct basis for pseudo-spectral scft.
54 */
55 template <int D>
56 void Basis<D>::makeBasis(Mesh<D> const & mesh,
57 UnitCell<D> const & unitCell,
58 std::string groupName)
59 {
60 SpaceGroup<D> group;
61 readGroup(groupName, group);
62 makeBasis(mesh, unitCell, group);
63 }
64
65 /*
66 * Construct a symmetry-adapted basis for pseudo-spectral scft.
67 */
68 template <int D>
69 void Basis<D>::makeBasis(Mesh<D> const & mesh,
70 UnitCell<D> const & unitCell,
71 SpaceGroup<D> const & group)
72 {
73 // Precondition: Check compatibility of mesh with space group
74 group.checkMeshDimensions(mesh.dimensions());
75
76 // Save pointers to mesh and unit cell
77 meshPtr_ = &mesh;
78 unitCellPtr_ = &unitCell;
79
80 // Allocate arrays
81 nWave_ = mesh.size();
82 waves_.allocate(nWave_);
83 waveIds_.allocate(nWave_);
84
85 // Make sorted array of waves
86 makeWaves();
87
88 // Identify stars of waves that are related by symmetry
89 makeStars(group);
90
91 // Apply validity test suite
92 bool valid = isValid();
93 if (!valid) {
94 UTIL_THROW("Basis failed validity check");
95 }
96
97 isInitialized_ = true;
98 }
99
100 /*
101 * Construct ordered list of waves.
102 *
103 * On exit:
104 * - Array waves_ contains list of waves ordered by sqNorm.
105 * - Each wave has indicesDft, indicesBz and sqNorm set.
106 * - Array stars_ is still empty.
107 */
108 template <int D>
110 {
111 IntVec<D> meshDimensions = mesh().dimensions();
112 std::vector< TWave<D> > twaves;
113 twaves.reserve(nWave_);
114
115 // Loop over dft mesh to generate all waves, add to twaves
116 TWave<D> w;
117 IntVec<D> v;
118 MeshIterator<D> itr(mesh().dimensions());
119 for (itr.begin(); !itr.atEnd(); ++itr) {
120 w.indicesDft = itr.position();
121 v = shiftToMinimum(w.indicesDft, meshDimensions, *unitCellPtr_);
122 w.indicesBz = v;
123 w.sqNorm = unitCell().ksq(v);
124 twaves.push_back(w);
125 }
126
127 // Sort temporary array twaves
128 TWaveNormComp<D> comp;
129 std::sort(twaves.begin(), twaves.end(), comp);
130
131 // Copy temporary array twaves into member variable waves_
132 for (int i = 0; i < nWave_; ++i) {
133 waves_[i].sqNorm = twaves[i].sqNorm;
134 waves_[i].indicesDft = twaves[i].indicesDft;
135 waves_[i].indicesBz = twaves[i].indicesBz;
136 }
137
138 }
139
140 /*
141 * Complete construction of a basis by grouping presorted waves into
142 * stars and completing initialization of all Wave and Star objects.
143 */
144 template <int D>
145 void Basis<D>::makeStars(SpaceGroup<D> const & group)
146 {
147
148 /*
149 * Conceptual definitions:
150 *
151 * - A "list" is a set of wavevectors of equal magnitude.
152 * - A "star" is a set of wavevectors that are related by symmetry.
153 *
154 * Each list may contain one or more complete stars. Complete lists
155 * are identified as an intermediate step in identification of stars.
156 *
157 * During initial processing, wavevectors are temporarily stored in
158 * TWave<D> objects. The following local containers of TWave<D>
159 * objects are used:
160 *
161 * list - a std::set of waves of equal norm (a "list")
162 * star - a std::set of symmetry-related waves (a "star")
163 * tempStar - a sorted star, sorted by descending indicesBz
164 * tempList - a sorted list, with contiguous sorted stars
165 *
166 * Reasons for choice of some C++ standard lib container types:
167 * std::set for "list", "star" prevents duplicates, allows removal
168 * std::vector for tempStar allows use of std::sort
169 */
170
171 // Local TWave<D> containers and associated iterators
172 std::set< TWave<D>, TWaveDftComp<D> > list;
173 std::set< TWave<D>, TWaveDftComp<D> > star;
174 std::vector< TWave<D> > tempStar;
175 GArray< TWave<D> > tempList;
176 typename std::set< TWave<D>, TWaveDftComp<D> >::iterator rootItr;
177 typename std::set< TWave<D>, TWaveDftComp<D> >::iterator setItr;
178
179 // Local variables
180 TWave<D> wave;
181 Basis<D>::Star newStar;
182 std::complex<double> coeff;
183 double Gsq;
184 double Gsq_max;
185 double phase_diff;
186 const double twoPi = 2.0*Constants::Pi;
187 const double epsilon = 1.0E-8;
188 IntVec<D> meshDimensions = mesh().dimensions();
189 IntVec<D> rootVecBz; // BZ indices for root of this star
190 IntVec<D> rootVecDft; // DFT indices for root of this star
191 IntVec<D> vec; // Indices of temporary wavevector
192 IntVec<D> nVec; // Indices of inverse of a wavevector
193 int listBegin = 0; // id of first wave in this list
194 int listEnd = 0; // (id of last wave in this list) + 1
195 int listSize; // listEnd - listBegin
196 int starBegin = 0; // id of first wave in this star
197 int i, j, k;
198 bool cancel;
199
200 /*
201 * Overview of algorithm:
202 *
203 * Precondition: Wavevectors in the array waves_ are sorted in
204 * nondecreasing order by wavevector norm.
205 *
206 * Loop over index i of array waves_ {
207 *
208 * Search for end of a "list" (i.e., contiguous block of waves
209 * of equal magnitude) by identifying changes in magnitude.
210 * The resulting list has indices [listBegin,listEnd-1].
211 * Set newList true.
212 *
213 * // Each list may contain one or more stars.
214 * // Process the newly identified list to identify stars
215 * If (newList) {
216 *
217 * Copy all waves in the range into std::set list
218 *
219 * Set rootItr to the first wave in the list
220 *
221 * // Loop over stars within the list
222 * while (list.size() > 0) {
223 *
224 * // To generate a star from a root wave rootItr,
225 * // loop over symmetry operations of space group.
226 * For each symmetry operation group[j] {
227 * Compute vec = (rootItr->indicesBz)*group[j]
228 * Set phase = rootItr->indicesBz .dot. (group[j].t)
229 * Check for cancellation of the star, set "cancel" flag
230 * Add wave to std::set<TWave> star if not added before
231 * // Here, use of a std::set simplifies test of uniqueness
232 * }
233 *
234 * Copy all waves from star to std::vector<TWave> tempStar
235 * Sort tempStar by indicesBz, in descending order
236 * // Here, use of a std::vector for tempStar allows sorting
237 *
238 * // Add waves in star to tempList and remove from list
239 * For each wave in tempStar {
240 * Append the wave to GArray<TWave> tempList
241 * Erase the wave from std::set<TWave> list
242 * // Here, use of a std::set for list simplifies erasure
243 * }
244 *
245 * Initialize a Star object named newStar
246 * Assign values to members beginId, endId, size, cancel
247 *
248 * // Assign values of newStar.invertFlag, rootItr, nextInvert
249 * if (nextInvert == -1) {
250 * // This is the second star in pair
251 * newStar.invertFlag = -1;
252 * nextInvert = 1;
253 * Set rootItr to the first wave in remaining list
254 * } else {
255 * Search for inverse of rootItr in this star
256 * if inverse is in this star {
257 * newStar.invertFlag = 0
258 * nextInvert = 1;
259 * Set rootItr to the first wave in remaining list
260 * } else
261 * Search for inverse of rootItr in this remaining list
262 * if the inverse is in the remaining list {
263 * newStar.invertFlag = 1
264 * nextInvert = -1;
265 * set rootItr to inverse of current root
266 * }
267 * }
268 *
269 * Append newStar object to GArray<Star> stars_
270 *
271 * } // end loop over stars in a single list
272 *
273 * // At this point, tempList contains the contents of the
274 * // waves_ array occupying the range [beginId, endId-1],
275 * // grouped by stars, with waves within each star sorted
276 * // by indexBz.
277 *
278 * // Overwrite the block of array waves_ with indices in the
279 * // range [beginId, endId-1] with the contents of tempList.
280 * For each wave in tempList {
281 * Copy a TWave in tempList to a Basis:Wave in waves_
282 * Assign a complex coefficient of unit norm to the Wave
283 * }
284 *
285 * // At this point, coefficients of waves have unit magnitude
286 * // and correct relative phases within each star, but not the
287 * // final absolute phases or magnitude.
288 *
289 * } // finish processing of one list (waves of equal norm)
290 *
291 * } // End initial processing of all waves and stars
292 *
293 * // Set phases of wave coefficients
294 * For each star in array stars_ {
295 * if star is closed under inversion (star.invertFlag == 0) {
296 * Find the inverse of every wave in this star
297 * if star is cancelled {
298 * Set coefficients of all waves to zero
299 * } else {
300 * Set the root to the first wave in the star
301 * For each wave in star:
302 * Divide coeff by the root coefficient
303 * }
304 * if (coeffs of root & inverse are not complex conjugates){
305 * Divide all coeffs by a common phasor chosen to obtain
306 * complex conjugate coefficients for root and inverse
307 * }
308 * }
309 * } else
310 * if (star.invertFlag == 1) {
311 * Find the inverse of every wave in this star and next star
312 * Set root of this star to the 1st wave in this star
313 * Set partner to the inverse of the root of this star
314 * If this star is cancelled {
315 * Set coefficients in this star and next to zero
316 * } else {
317 * For each wave in this star:
318 * Divide coeff by the root coefficient
319 * }
320 * For each wave in the next star:
321 * Divide coeff by the partner coefficient
322 * }
323 * }
324 * }
325 * // Note: If star.invertFlag = -1, do nothing because properties
326 * // of this star were all set when processing its partner.
327 * }
328 *
329 * // For all waves, normalize coefficients and set starId
330 * For each star in array stars_ {
331 * For each wave in this star {
332 * Set Wave.starId
333 * Divide coefficient by sqrt(double(star.size))
334 * }
335 * }
336 *
337 * // For all waves, set implicit member and add to look up table
338 * For each wave in array waves_ {
339 * Set Wave::implicit attribute
340 * Set waveIds_[rank] = i
341 * }
342 */
343
344 // Loop over all waves (initial processing of waves)
345 nBasis_ = 0;
346 nBasisWave_ = 0;
347 Gsq_max = waves_[0].sqNorm;
348 for (i = 1; i <= nWave_; ++i) {
349
350 // Determine if this wave begins a new list
351 bool newList = false;
352 if (i == nWave_) {
353 listEnd = i;
354 listSize = listEnd - listBegin;
355 newList = true;
356 } else {
357 Gsq = waves_[i].sqNorm;
358 if (Gsq > Gsq_max + epsilon) {
359 Gsq_max = Gsq;
360 listEnd = i;
361 listSize = listEnd - listBegin;
362 newList = true;
363 }
364 }
365
366 // Process completed list of wavectors of equal norm
367 if (newList) {
368
369 // Copy waves of equal norm into std::set "list"
370 list.clear();
371 tempList.clear();
372 for (j = listBegin; j < listEnd; ++j) {
373 wave.indicesDft = waves_[j].indicesDft;
374 wave.indicesBz = waves_[j].indicesBz;
375 wave.sqNorm = waves_[j].sqNorm;
376 if (j > listBegin) {
377 UTIL_CHECK( std::abs(wave.sqNorm-waves_[j].sqNorm)
378 < 2.0*epsilon );
379 }
380 list.insert(wave);
381 }
382
383 // On entry to each iteration of the loop over stars,
384 // rootItr and nextInvert are known. The iterator rootItr
385 // points to the wave in the remaining list that will be
386 // used as the root of the next star. The flag nextInvert
387 // is equal to -1 iff the previous star was the first of
388 // a pair that are open under inversion, and is equal
389 // to + 1 otherwise.
390
391 // Initial values for first star in this list
392 rootItr = list.begin();
393 int nextInvert = 1;
394
395 // Loop over stars with a list of waves of equal norm,
396 // removing each star from the list as it is identified.
397 // The root of the next star must have been chosen on
398 // entry to each iteration of this loop.
399
400 while (list.size() > 0) {
401
402 rootVecBz = rootItr->indicesBz;
403 rootVecDft = rootItr->indicesDft;
404 Gsq = rootItr->sqNorm;
405 cancel = false;
406 star.clear();
407
408 // Construct a star from root vector, by applying every
409 // symmetry operation in the group to the root wavevector.
410 for (j = 0; j < group.size(); ++j) {
411
412 // Apply symmetry (i.e., multiply by rotation matrix)
413 // vec = rotated wavevector.
414 vec = rootVecBz*group[j];
415
416 // Check that rotated vector has same norm as root.
417 UTIL_CHECK(std::abs(Gsq - unitCell().ksq(vec)) < epsilon);
418
419 // Initialize TWave object associated with rotated wave
420 wave.sqNorm = Gsq;
421 wave.indicesBz = shiftToMinimum(vec, meshDimensions,
422 *unitCellPtr_);
423 wave.indicesDft = vec;
424 mesh().shift(wave.indicesDft);
425
426 // Compute phase for coeff. of wave in basis function.
427 // Convention -pi < phase <= pi.
428 wave.phase = 0.0;
429 for (k = 0; k < D; ++k) {
430 wave.phase += rootVecBz[k]*(group[j].t(k));
431 }
432 while (wave.phase > 0.5) {
433 wave.phase -= 1.0;
434 }
435 while (wave.phase <= -0.5) {
436 wave.phase += 1.0;
437 }
438 wave.phase *= twoPi;
439
440 // Check for cancellation of star: The star is
441 // cancelled if application of any symmetry operation
442 // in the group to the root vector yields a rotated
443 // vector equivalent to the root vector but with a
444 // nonzero phase, creating a contradiction.
445
446 if (wave.indicesDft == rootVecDft) {
447 if (std::abs(wave.phase) > 1.0E-6) {
448 cancel = true;
449 }
450 }
451
452 // Search for an equivalent wave already in the star.
453 // Note: Equivalent waves have equal DFT indices.
454 setItr = star.find(wave);
455
456 if (setItr == star.end()) {
457
458 // If no equivalent wave is found in the star,
459 // then add this wave to the star
460 star.insert(wave);
461
462 } else {
463
464 // If an equivalent wave is found, check if the
465 // phases are equivalent. If not, the star is
466 // cancelled.
467
468 phase_diff = setItr->phase - wave.phase;
469 while (phase_diff > 0.5) {
470 phase_diff -= 1.0;
471 }
472 while (phase_diff <= -0.5) {
473 phase_diff += 1.0;
474 }
475 if (std::abs(phase_diff) > 1.0E-6) {
476 cancel = true;
477 }
478
479 }
480
481 }
482
483 // Copy all waves from std::set star to std::vector tempStar
484 tempStar.clear();
485 setItr = star.begin();
486 for ( ; setItr != star.end(); ++setItr) {
487 tempStar.push_back(*setItr);
488 }
489
490 // Sort tempStar, in descending order by indicesBz.
491 TWaveBzComp<D> waveBzComp;
492 std::sort(tempStar.begin(), tempStar.end(), waveBzComp);
493
494 // Append contents of tempStar to tempList, erase from list
495 int tempStarSize = tempStar.size();
496 for (j = 0; j < tempStarSize; ++j) {
497 list.erase(tempStar[j]);
498 tempList.append(tempStar[j]);
499 }
500 UTIL_CHECK((int)(tempList.size()+list.size()) == listSize);
501
502 // If this star is not cancelled, increment the number of
503 // basis functions (nBasis_) & waves in basis (nBasisWave_)
504 if (!cancel) {
505 ++nBasis_;
506 nBasisWave_ += star.size();
507 }
508
509 // Initialize a Star object
510 // newStar.eigen = Gsq;
511 newStar.beginId = starBegin;
512 newStar.endId = newStar.beginId + star.size();
513 newStar.size = star.size();
514 newStar.cancel = cancel;
515 // Note: newStar.starInvert is not yet known
516
517 // Determine invertFlag, rootItr and nextInvert
518 if (nextInvert == -1) {
519
520 // If this star is 2nd of a pair related by inversion,
521 // set root for next star to 1st wave of remaining list.
522
523 newStar.invertFlag = -1;
524 rootItr = list.begin();
525 nextInvert = 1;
526
527 } else {
528
529 // If this star is not the 2nd of a pair of partners,
530 // then determine if it is closed under inversion.
531
532 // Compute inverse nVec of root vector in FBZ
533 nVec.negate(rootVecBz);
534
535 // Shift inverse nVec to the DFT mesh
536 (*meshPtr_).shift(nVec);
537
538 // Search for inverse of root vector within this star
539 bool inverseFound = false;
540 setItr = star.begin();
541 for ( ; setItr != star.end(); ++setItr) {
542 if (nVec == setItr->indicesDft) {
543 inverseFound = true;
544 break;
545 }
546 }
547
548 if (inverseFound) {
549
550 // If this star is closed under inversion, the root
551 // of next star is the 1st vector of remaining list.
552
553 newStar.invertFlag = 0;
554 rootItr = list.begin();
555 nextInvert = 1;
556
557 } else {
558
559 // This star is open under inversion, and is the
560 // first star of a pair related by inversion
561
562 newStar.invertFlag = 1;
563 nextInvert = -1;
564
565 // Find inverse of the root of this star in the
566 // remaining list, and use this inverse as the
567 // root of the next star.
568
569 setItr = list.begin();
570 for ( ; setItr != list.end(); ++setItr) {
571 if (nVec == setItr->indicesDft) {
572 inverseFound = true;
573 rootItr = setItr;
574 break;
575 }
576 }
577 // If inverseFound, then rootVecDft = nVec
578
579 // Failure to find the inverse here is an error:
580 // It must be either in this star or remaining list
581
582 if (!inverseFound) {
583 std::cout << "Inverse not found for: " << "\n";
584 std::cout << " vec (ft):"
585 << rootVecDft <<"\n";
586 std::cout << " vec (bz):"
587 << rootVecBz <<"\n";
588 std::cout << "-vec (dft):" << nVec << "\n";
589 UTIL_CHECK(inverseFound);
590 }
591
592 }
593
594 }
595
596 stars_.append(newStar);
597 starBegin = newStar.endId;
598
599 }
600 // End loop over stars within a list.
601
602 UTIL_CHECK(list.size() == 0);
603 UTIL_CHECK(tempList.size() == listEnd - listBegin);
604
605 // Copy tempList into corresponding section of waves_,
606 // overwriting the section of waves_ used to create the list.
607 // Compute a complex coefficient of unit norm for each wave.
608 for (j = 0; j < tempList.size(); ++j) {
609 k = j + listBegin;
610 waves_[k].indicesDft = tempList[j].indicesDft;
611 waves_[k].indicesBz = tempList[j].indicesBz;
612 waves_[k].sqNorm = tempList[j].sqNorm;
613 coeff = std::complex<double>(0.0, tempList[j].phase);
614 coeff = exp(coeff);
615 if (std::abs(imag(coeff)) < 1.0E-6) {
616 coeff = std::complex<double>(real(coeff), 0.0);
617 }
618 if (std::abs(real(coeff)) < 1.0E-6) {
619 coeff = std::complex<double>(0.0, imag(coeff));
620 }
621 waves_[k].coeff = coeff;
622 }
623
624 // Processing of list is now complete.
625 // Here, waves_[k].coeff has unit absolute magnitude, and
626 // correct relative phases for waves within a star, but
627 // the coeff may not be unity for the first or last wave in
628 // the star.
629
630 listBegin = listEnd;
631 }
632 // Finished processing a list of waves of equal norm
633
634 } // End loop over all waves
635 nStar_ = stars_.size();
636 // Complete initial processing of all lists and stars
637
638 /*
639 * Conventions for phases of wave coefficients (imposed below):
640 * - Coefficients of the root of each star and its inverse must
641 * be complex conjugates.
642 * - In a closed star (starInvert = 0), the coefficient of the
643 * root must have a non-negative real part. If the root
644 * coefficient is pure imaginary, the imaginary part must
645 * be negative.
646 * - In a pair of open stars that are related by inversion
647 * symmetry, the coefficients of the first wave of the first
648 * star and its inverse must have real coefficients.
649 */
650
651 // Final processing of phases of of waves in stars:
652 std::complex<double> rootCoeff; // Coefficient of root wave
653 std::complex<double> partCoeff; // Coefficient of partner of root
654 std::complex<double> d;
655 int rootId, partId;
656 for (i = 0; i < nStar_; ++i) {
657
658 // Treat open and closed stars differently
659
660 if (stars_[i].invertFlag == 0) {
661
662 // First, assign inverseId for each wave in star
663 for (j = stars_[i].beginId; j < stars_[i].endId; ++j) {
664 if (waves_[j].inverseId < 0) { // if inverseId is unassigned
665
666 // Compute nVec = inverse of root, shifted to DFT mesh
667 nVec.negate(waves_[j].indicesBz);
668 (*meshPtr_).shift(nVec);
669
670 // Find inverse
671
672 // Check in the position that the inverse is
673 // expected to be located for a typical star
674 k = stars_[i].endId - 1 - (j - stars_[i].beginId);
675 if (nVec == waves_[k].indicesDft) {
676 waves_[j].inverseId = k;
677 waves_[k].inverseId = j;
678 } else {
679 // Inverse not in expected position, search full star
680 // (this usually occurs for stars on the edge of the
681 // Brillouin zone)
682 for (k = j; k < stars_[i].endId; ++k) {
683 if (nVec == waves_[k].indicesDft) {
684 waves_[j].inverseId = k;
685 waves_[k].inverseId = j;
686 break;
687 }
688 }
689 }
690
691 // For invertFlag == 0, failure to find nVec in this
692 // star is a fatal error
693 if (waves_[j].inverseId < 0) {
694 std::cout << "\n";
695 std::cout << "Inverse not found in closed star"
696 << std::endl;
697 std::cout << "G = " << waves_[j].indicesBz
698 << ", coeff = " << waves_[j].coeff
699 << std::endl;
700 std::cout << "All waves in star " << i
701 << std::endl;
702 for (k=stars_[i].beginId; k < stars_[i].endId; ++k) {
703 std::cout << waves_[k].indicesBz << " "
704 << waves_[k].coeff << std::endl;
705 }
706 UTIL_CHECK(waves_[j].inverseId >= 0);
707 }
708
709 }
710 }
711
712 // Identify root of this star (star i)
713 // Set the root to be the first wave in the star
714
715 rootId = stars_[i].beginId;
716 stars_[i].waveBz = waves_[rootId].indicesBz;
717
718 if (stars_[i].cancel) {
719
720 // If the star is cancelled, set all coefficients to zero
721 std::complex<double> czero(0.0, 0.0);
722 for (j = stars_[i].beginId; j < stars_[i].endId; ++j) {
723 waves_[j].coeff = czero;
724 }
725
726 } else { // if not cancelled
727
728 // Set partId to index of the inverse of the root
729 partId = waves_[rootId].inverseId;
730
731 // Divide all coefficients by the root coefficient
732 rootCoeff = waves_[rootId].coeff;
733 for (j = stars_[i].beginId; j < stars_[i].endId; ++j) {
734 waves_[j].coeff /= rootCoeff;
735 }
736 rootCoeff = waves_[rootId].coeff;
737 UTIL_CHECK(std::abs(real(rootCoeff) - 1.0) < 1.0E-9);
738 UTIL_CHECK(std::abs(imag(rootCoeff)) < 1.0E-9);
739
740 // Require coefficients of root and inverse are conjugates
741 if (partId != rootId) {
742
743 partCoeff = waves_[partId].coeff;
744 UTIL_CHECK(std::abs(std::abs(partCoeff) - 1.0) < 1.0E-9);
745 if (std::abs(partCoeff - rootCoeff) > 1.0E-6) {
746 d = sqrt(partCoeff);
747 if (real(d) < -1.0E-4) {
748 d = -d;
749 } else
750 if (std::abs(real(d)) <= 1.0E-4) {
751 if (imag(d) < 0.0) {
752 d = -d;
753 }
754 }
755 for (j=stars_[i].beginId; j < stars_[i].endId; ++j){
756 waves_[j].coeff /= d;
757 }
758 }
759
760 }
761
762 } // end if (cancel) ... else ...
763
764 } // end if (stars_[i].invertFlag == 0)
765 else
766 if (stars_[i].invertFlag == 1) {
767
768 // Process a pair of open stars related by inversion.
769
770 // Preconditions:
771 UTIL_CHECK(stars_[i].size == stars_[i+1].size);
772 UTIL_CHECK(stars_[i].cancel == stars_[i+1].cancel);
773 // UTIL_CHECK(stars_[i+1].invertFlag == -1);
774
775 // First, assign inverseId for each wave in pair of stars
776 for (j = stars_[i].beginId; j < stars_[i].endId; ++j) {
777 if (waves_[j].inverseId < 0) { // if inverseId is unassigned
778
779 // Compute nVec = inverse of root, shifted to DFT mesh
780 nVec.negate(waves_[j].indicesBz);
781 (*meshPtr_).shift(nVec);
782
783 // Find inverse
784
785 // Check in the position that the inverse is expected
786 // to be located for a typical pair of stars
787 k = stars_[i+1].endId - 1 - (j - stars_[i].beginId);
788 if (nVec == waves_[k].indicesDft) {
789 waves_[j].inverseId = k;
790 waves_[k].inverseId = j;
791 } else {
792 // Inverse not in expected position, search full star
793 // (this usually occurs for stars on the edge of the
794 // Brillouin zone)
795 k = stars_[i+1].beginId;
796 for ( ; k < stars_[i+1].endId; ++k) {
797 if (nVec == waves_[k].indicesDft) {
798 waves_[j].inverseId = k;
799 waves_[k].inverseId = j;
800 break;
801 }
802 }
803 }
804
805 // For invertFlag = 1, failure to find nVec in the
806 // next star is a fatal error
807 UTIL_CHECK(waves_[j].inverseId >= 0);
808
809 }
810 }
811
812 // Check that inverseId was assigned for all waves in the
813 // next star
814 for (j = stars_[i+1].beginId; j < stars_[i+1].endId; ++j) {
815 UTIL_CHECK(waves_[j].inverseId >= 0);
816 }
817
818 // Identify root of this star (star i)
819 // Set the root to be the first wave in the star
820 rootId = stars_[i].beginId;
821 stars_[i].waveBz = waves_[rootId].indicesBz;
822
823 // Identify root of the next star (star i+1)
824 // Set the root to be the inverse of the root of star i
825 partId = waves_[rootId].inverseId;
826 stars_[i+1].waveBz = waves_[partId].indicesBz;
827
828 if (stars_[i].cancel) {
829
830 std::complex<double> czero(0.0, 0.0);
831 for (j = stars_[i].beginId; j < stars_[i].endId; ++j) {
832 waves_[j].coeff = czero;
833 }
834 for (j = stars_[i+1].beginId; j < stars_[i+1].endId; ++j) {
835 waves_[j].coeff = czero;
836 }
837
838 } else { // if star is not cancelled
839
840 // Divide all coefficients in this star by root coeff
841 rootCoeff = waves_[rootId].coeff;
842 for (j = stars_[i].beginId; j < stars_[i].endId; ++j) {
843 waves_[j].coeff /= rootCoeff;
844 }
845
846 // Divide coefficients in next star by a partner coeff.
847 partCoeff = waves_[partId].coeff;
848 for (j = stars_[i+1].beginId; j < stars_[i+1].endId; ++j) {
849 waves_[j].coeff /= partCoeff;
850 }
851
852 } // end if (cancel) ... else ...
853
854 } // end if (invertFlag==0) ... else if (invertFlag==1) ...
855
856 // Note: If invertFlag == -1, do nothing and continue
857 // Related stars with invertFlag == 1 and -1 are treated together
858
859 } // end loop over stars
860
861 // For all waves, normalize coefficients and set starId
862 for (i = 0; i < nStar_; ++i) {
863 double snorm = 1.0/sqrt(double(stars_[i].size));
864 for (j = stars_[i].beginId; j < stars_[i].endId; ++j) {
865 waves_[j].coeff *= snorm;
866 waves_[j].starId = i;
867 }
868 }
869
870 // Set tiny real and imaginary parts to zero (due to round-off)
871 for (i = 0; i < nWave_; ++i) {
872 if (std::abs(real(waves_[i].coeff)) < 1.0E-8) {
873 waves_[i].coeff
874 = std::complex<double>(0.0, imag(waves_[i].coeff));
875 }
876 if (std::abs(imag(waves_[i].coeff)) < 1.0E-8) {
877 waves_[i].coeff
878 = std::complex<double>(real(waves_[i].coeff), 0.0);
879 }
880 }
881
882 // For each wave, set implicit attribute and add to look-up table
883 for (i = 0; i < nWave_; ++i) {
884 vec = waves_[i].indicesDft;
885
886 // Validity check - check that vec is in dft mesh
887 for (j = 0; j < D; ++j) {
888 UTIL_CHECK(vec[j] >= 0);
889 UTIL_CHECK(vec[j] < meshDimensions[j]);
890 }
891
892 // Set implicit attribute
893 if ((vec[D-1] + 1) > (meshDimensions[D-1]/2 + 1)) {
894 waves_[i].implicit = true;
895 } else {
896 waves_[i].implicit = false;
897 }
898
899 // Look up table for waves
900 waveIds_[mesh().rank(vec)] = i;
901 }
902
903 // Set Star::starId and Star::basisId, add to look-up table starIds_.
904 starIds_.allocate(nBasis_);
905 j = 0;
906 for (i = 0; i < nStar_; ++i) {
907 stars_[i].starId = i;
908 if (stars_[i].cancel) {
909 stars_[i].basisId = -1;
910 } else {
911 stars_[i].basisId = j;
912 UTIL_CHECK(j < nBasis_);
913 starIds_[j] = i;
914 ++j;
915 }
916 }
917 UTIL_CHECK(j == nBasis_);
918
919 }
920
921 // Return value of nBasis
922 template <int D>
924 { return nBasis_; }
925
926
927 template <int D>
928 void Basis<D>::outputWaves(std::ostream& out, bool outputAll) const
929 {
930 out << "N_wave" << std::endl;
931 if (outputAll) {
932 out << " " << nWave_ << std::endl;
933 } else {
934 out << " " << nBasisWave_ << std::endl;
935 }
936 int i, j, k, starId;
937 k = 0;
938 for (i = 0; i < nWave_; ++i) {
939 starId = waves_[i].starId;
940 if (outputAll || (!stars_[starId].cancel)) {
941 out << Int(k, 8);
942 out << Int(i, 8);
943 for (j = 0; j < D; ++j) {
944 out << Int(waves_[i].indicesBz[j], 5);
945 }
946 out << Int(waves_[i].starId, 6);
947 out << " " << Dbl(waves_[i].coeff.real(), 15);
948 out << " " << Dbl(waves_[i].coeff.imag(), 15);
949 out << std::endl;
950 k++;
951 }
952 }
953 }
954
955
956 template <int D>
957 void Basis<D>::outputStars(std::ostream& out, bool outputAll) const
958 {
959 // Output number of stars in appropriate format
960 if (outputAll) {
961 out << "N_star" << std::endl
962 << " " << nStar_ << std::endl;
963 } else {
964 out << "N_basis" << std::endl
965 << " " << nBasis_ << std::endl;
966 }
967
968 // Loop over stars
969 int i, j;
970 for (i = 0; i < nStar_; ++i) {
971 if (outputAll || (!stars_[i].cancel)) {
972 out << Int(stars_[i].basisId, 6); // basisId
973 out << Int(i, 6); // starId
974 out << Int(stars_[i].size, 5)
975 << Int(stars_[i].beginId, 8)
976 << Int(stars_[i].endId, 8)
977 << Int(stars_[i].invertFlag, 4);
978 if (outputAll) {
979 out << Int(stars_[i].cancel, 4);
980 }
981 for (j = 0; j < D; ++j) {
982 out << Int(stars_[i].waveBz[j], 6);
983 }
984 out << std::endl;
985 }
986 }
987 }
988
989 template <int D>
990 bool Basis<D>::isValid() const
991 {
992 double Gsq;
993 IntVec<D> v;
994 int is, ib, iw, iwp, j;
995
996 // Check total number of waves == # of grid points
997 if (nWave_ != mesh().size()) {
998 std::cout << "nWave != size of mesh" << std::endl;
999 return false;
1000 }
1001
1002 // Loop over dft mesh to check consistency of waveIds_ and waves_
1003 MeshIterator<D> itr(mesh().dimensions());
1004 for (itr.begin(); !itr.atEnd(); ++itr) {
1005 v = itr.position();
1006 iw = waveId(v);
1007 if (wave(iw).indicesDft != v) {
1008 std::cout << "Inconsistent waveId and Wave::indicesDft"
1009 << std::endl;
1010 return false;
1011 }
1012 }
1013
1014 // Loop over elements of waves_, check consistency of wave data.
1015 for (iw = 0; iw < nWave_; ++iw) {
1016
1017 // Check sqNorm
1018 v = waves_[iw].indicesBz;
1019 Gsq = unitCell().ksq(v);
1020 if (std::abs(Gsq - waves_[iw].sqNorm) > 1.0E-8) {
1021 std::cout << "\n";
1022 std::cout << "Incorrect sqNorm:" << "\n"
1023 << "wave.indicesBz = " << "\n"
1024 << "wave.sqNorm = " << waves_[iw].sqNorm << "\n"
1025 << "|v|^{2} = " << Gsq << "\n";
1026 return false;
1027 }
1028
1029 // Check that wave indicesBz is an image of indicesDft
1030 mesh().shift(v);
1031 if (v != waves_[iw].indicesDft) {
1032 std::cout << "\n";
1033 std::cout << "shift(indicesBz) != indicesDft" << std::endl;
1034 return false;
1035 }
1036
1037 // Compare Wave::starId to Star::beginId and Star::endId
1038 is = waves_[iw].starId;
1039 if (iw < stars_[is].beginId) {
1040 std::cout << "\n";
1041 std::cout << "Wave::starId < Star::beginId" << std::endl;
1042 return false;
1043 }
1044 if (iw >= stars_[is].endId) {
1045 std::cout << "\n";
1046 std::cout << "Wave::starId >= Star::endId" << std::endl;
1047 return false;
1048 }
1049
1050 // Check that inverseId has been assigned
1051 if (waves_[iw].inverseId < 0) {
1052 std::cout << "\n";
1053 std::cout << "Wave::inverseId not assigned\n";
1054 std::cout << "G = " << waves_[iw].indicesBz << std::endl;
1055 return false;
1056 }
1057
1058 // Check that inverseId points to the correct wave
1059 v.negate(waves_[iw].indicesBz);
1060 mesh().shift(v);
1061 iwp = waves_[iw].inverseId;
1062 if (waves_[iwp].indicesDft != v) {
1063 std::cout << "\n";
1064 std::cout << "Wave::inverseId is not inverse" << std::endl;
1065 std::cout << "G = " << waves_[iw].indicesBz << std::endl;
1066 std::cout << "-G (from inverseId) = "
1067 << waves_[iwp].indicesBz << std::endl;
1068 return false;
1069 }
1070
1071 // Check that inverseId of the inverse wave is correct
1072 if (waves_[iwp].inverseId != iw) {
1073 std::cout << "\n";
1074 std::cout << "Wave::inverseId values do not agree\n";
1075 std::cout << "+G = " << waves_[iw].indicesBz << std::endl;
1076 std::cout << "-G = " << waves_[iwp].indicesBz << std::endl;
1077 return false;
1078 }
1079
1080 // Check that either this wave or its inverse is explicit
1081 if (waves_[iw].implicit == true && waves_[iwp].implicit == true)
1082 {
1083 std::cout << "\n";
1084 std::cout << "Wave and its inverse are both implicit";
1085 std::cout << "+G = " << waves_[iw].indicesBz << std::endl;
1086 std::cout << "-G = " << waves_[iwp].indicesBz << std::endl;
1087 return false;
1088 }
1089 }
1090
1091 // Loop over all stars (elements of stars_ array)
1092 int nWave = 0;
1093 for (is = 0; is < nStar_; ++is) {
1094
1095 // Check star size
1096 nWave += stars_[is].size;
1097 if (stars_[is].size != stars_[is].endId - stars_[is].beginId) {
1098 std::cout << "\n";
1099 std::cout << "Inconsistent Star::size:" << std::endl;
1100 std::cout << "Star id " << is << std::endl;
1101 std::cout << "star size " << stars_[is].size << std::endl;
1102 std::cout << "Star begin " << stars_[is].beginId << std::endl;
1103 std::cout << "Star end " << stars_[is].endId << std::endl;
1104 return false;
1105 }
1106 if (is > 0) {
1107 if (stars_[is].beginId != stars_[is-1].endId) {
1108 std::cout << "\n";
1109 std::cout << "Star ranges not consecutive:" << std::endl;
1110 std::cout << "Star id " << is << std::endl;
1111 std::cout << "stars_[" << is << "]" << ".beginId = "
1112 << stars_[is].beginId << std::endl;
1113 std::cout << "stars_[" << is-1 << "]" << ".endId = "
1114 << stars_[is-1].endId << std::endl;
1115 return false;
1116 }
1117 }
1118
1119 // Check waveBz indices of star
1120 if (stars_[is].invertFlag == -1) {
1121 v.negate(stars_[is-1].waveBz);
1122 v = shiftToMinimum(v, mesh().dimensions(), *unitCellPtr_);
1123 if (stars_[is].waveBz != v) {
1124 std::cout << "\n";
1125 std::cout << "waveBz of star is not inverse of waveBz "
1126 << "of previous star" << std::endl;
1127 std::cout << "star id " << is << std::endl;
1128 std::cout << "waveBz " << stars_[is].waveBz << std::endl;
1129 std::cout << "waveBz (previous star) "
1130 << stars_[is-1].waveBz << std::endl;
1131 return false;
1132 }
1133 } else {
1134 v = waves_[stars_[is].beginId].indicesBz;
1135 if (stars_[is].waveBz != v) {
1136 std::cout << "\n";
1137 std::cout << "waveBz of star != first wave of star"
1138 << std::endl;
1139 std::cout << "star id " << is << std::endl;
1140 std::cout << "waveBz " << stars_[is].waveBz
1141 << std::endl;
1142 std::cout << "first wave " << v << std::endl;
1143 return false;
1144 }
1145 }
1146
1147 // Check star ids of waves in star
1148 for (iw = stars_[is].beginId; iw < stars_[is].endId; ++iw) {
1149 if (waves_[iw].starId != is) {
1150 std::cout << "\n";
1151 std::cout << "Inconsistent Wave::starId :" << std::endl;
1152 std::cout << "star id " << is << std::endl;
1153 std::cout << "star beginId " << stars_[is].beginId << "\n";
1154 std::cout << "star endId " << stars_[is].endId << "\n";
1155 std::cout << "wave id " << iw << "\n";
1156 std::cout << "wave starId " << waves_[iw].starId << "\n";
1157 return false;
1158 }
1159 }
1160
1161 // Check Star::starId is equal to array index
1162 if (stars_[is].starId != is) {
1163 std::cout << "\n";
1164 std::cout << "stars_[is].starId != is for "
1165 << "is = " << is << "\n";
1166 return false;
1167 }
1168
1169 // Check Star::basisId and starIds_ look up table
1170 ib = stars_[is].basisId;
1171 if (stars_[is].cancel) {
1172 if (ib != -1) {
1173 std::cout << "\n";
1174 std::cout << "basisId != -1 for cancelled star\n";
1175 std::cout << "star id = " << is << "\n";
1176 return false;
1177 }
1178 } else {
1179 if (starIds_[ib] != is) {
1180 std::cout << "\n";
1181 std::cout << "starIds_[stars_[is].basisId] != is for: \n";
1182 std::cout << "is = " << is << "\n";
1183 std::cout << "ib = stars_[is].basisId = " << ib << "\n";
1184 std::cout << "starIds_[ib] = " << starIds_[ib]
1185 << "\n";
1186 return false;
1187 }
1188 }
1189
1190 // Check ordering of waves in star
1191 for (iw = stars_[is].beginId + 1; iw < stars_[is].endId; ++iw) {
1192 if (waves_[iw].indicesBz > waves_[iw-1].indicesBz) {
1193 std::cout << "\n";
1194 std::cout << "Failure of ordering by indicesB within star"
1195 << std::endl;
1196 return false;
1197 }
1198 if (waves_[iw].indicesBz == waves_[iw-1].indicesBz) {
1199 std::cout << "\n";
1200 std::cout << "Equal values of indicesBz within star"
1201 << std::endl;
1202 return false;
1203 }
1204 }
1205
1206 // Check that all coefficients are zero if star is cancelled
1207 if (stars_[is].cancel) {
1208 for (iw = stars_[is].beginId + 1; iw < stars_[is].endId; ++iw)
1209 {
1210 if (std::abs(waves_[iw].coeff) > 1.0E-8) {
1211 std::cout << "\n";
1212 std::cout << "Nonzero coefficient in a cancelled star"
1213 << "\n";
1214 std::cout << "G = " << waves_[iw].indicesBz
1215 << " coeff = " << waves_[iw].coeff
1216 << "\n";
1217 return false;
1218 }
1219 }
1220 }
1221
1222 } // End do loop over all stars
1223
1224 // Check that all waves in mesh are accounted for in stars
1225 if (stars_[nStar_-1].endId != mesh().size()) {
1226 std::cout << "\n";
1227 std::cout << "Star endId of last star != mesh size" << std::endl;
1228 return false;
1229 }
1230 if (nWave != mesh().size()) {
1231 std::cout << "\n";
1232 std::cout << "Sum of star sizes != mesh size" << std::endl;
1233 return false;
1234 }
1235
1236 // Loop over closed stars and related pairs of stars.
1237 // Test closure under inversion and conjugacy of coefficients.
1238 std::complex<double> cdel;
1239 bool cancel;
1240 is = 0;
1241 while (is < nStar_) {
1242 cancel = stars_[is].cancel;
1243
1244 if (stars_[is].invertFlag == 0) {
1245
1246 // Test that star is closed under inversion and real
1247 int begin = stars_[is].beginId;
1248 int end = stars_[is].endId;
1249 for (iw = begin; iw < end; ++iw) {
1250 iwp = waves_[iw].inverseId;
1251 if (waves_[iwp].starId != is) {
1252 std::cout << "\n";
1253 std::cout << "Inverse not found in closed star"
1254 << std::endl;
1255 std::cout << "G = " << waves_[iw].indicesBz
1256 << "coeff = " << waves_[iw].coeff
1257 << std::endl;
1258 std::cout << "All waves in star " << is << "\n";
1259 for (j=begin; j < end; ++j) {
1260 std::cout << waves_[j].indicesBz << " "
1261 << waves_[j].coeff << "\n";
1262 }
1263 return false;
1264 }
1265 if (!cancel) {
1266 cdel = std::conj(waves_[iwp].coeff);
1267 cdel -= waves_[iw].coeff;
1268 if (std::abs(cdel) > 1.0E-8) {
1269 std::cout << "\n";
1270 std::cout << "Function for closed star is not real:"
1271 << "\n";
1272 std::cout << "+G = " << waves_[iw].indicesBz
1273 << " coeff = " << waves_[iw].coeff
1274 << "\n";
1275 std::cout << "-G = " << waves_[iwp].indicesBz
1276 << " coeff = " << waves_[iwp].coeff
1277 << "\n";
1278 std::cout << "Coefficients are not conjugates."
1279 << "\n";
1280 std::cout << "All waves in star " << is << "\n";
1281 for (j=begin; j < end; ++j) {
1282 std::cout << waves_[j].indicesBz << " "
1283 << waves_[j].coeff << "\n";
1284 }
1285 return false;
1286 }
1287 }
1288 }
1289
1290 // Finished processing a closed star, increment counter is
1291 ++is;
1292
1293 } else {
1294
1295 // Test pairs of open stars
1296
1297 if (stars_[is].invertFlag != 1) {
1298 std::cout << "\n";
1299 std::cout << "Expected invertFlag == 1" << std::endl;
1300 return false;
1301 }
1302 if (stars_[is+1].invertFlag != -1) {
1303 std::cout << "\n";
1304 std::cout << "Expected invertFlag == -1" << std::endl;
1305 return false;
1306 }
1307 if (stars_[is+1].size != stars_[is].size) {
1308 std::cout << "\n";
1309 std::cout << "Partner stars of different size" << std::endl;
1310 return false;
1311 }
1312 if (stars_[is+1].cancel != stars_[is].cancel) {
1313 std::cout << "\n";
1314 std::cout << "Partners stars with different cancel flags"
1315 << std::endl;
1316 return false;
1317 }
1318
1319 // Begin and end wave ids for the first and second stars
1320 int begin1 = stars_[is].beginId;
1321 int end1 = stars_[is].endId;
1322 int begin2 = stars_[is+1].beginId;
1323 int end2 = stars_[is+1].endId;
1324
1325 // Check that inverse is in next star and check for
1326 // conjugate coefficients
1327
1328 // Loop over waves in first star
1329 for (iw = begin1; iw < end1; ++iw) {
1330 iwp = waves_[iw].inverseId;
1331 if (waves_[iwp].starId != is + 1) {
1332 std::cout << "\n";
1333 std::cout << "Inverse not found for G in open star"
1334 << std::endl;
1335 std::cout << "First star id = " << is << std::endl;
1336 std::cout << "+G = " << waves_[iw].indicesBz
1337 << "coeff = " << waves_[iw].coeff
1338 << std::endl;
1339 std::cout << "Waves in star " << is
1340 << " (starInvert ==1):" << "\n";
1341 for (j = begin1; j < end1; ++j) {
1342 std::cout << waves_[j].indicesBz << " "
1343 << waves_[j].coeff << "\n";
1344 }
1345 std::cout << "Waves in star " << is+1
1346 << " (starInvert == -1):" << "\n";
1347 for (j=begin2; j < end2; ++j) {
1348 std::cout << waves_[j].indicesBz << " "
1349 << waves_[j].coeff << "\n";
1350 }
1351 return false;
1352 }
1353 if (!cancel) {
1354 cdel = std::conj(waves_[iwp].coeff);
1355 cdel -= waves_[iw].coeff;
1356 if (std::abs(cdel) > 1.0E-8) {
1357 std::cout << "\n";
1358 std::cout << "Error of coefficients in open stars:"
1359 << "\n";
1360 std::cout << "First star id = " << is << std::endl;
1361 std::cout << "+G = " << waves_[iw].indicesBz
1362 << " coeff = " << waves_[iw].coeff
1363 << "\n";
1364 std::cout << "-G = " << waves_[iwp].indicesBz
1365 << " coeff = " << waves_[iwp].coeff
1366 << "\n";
1367 std::cout << "Coefficients are not conjugates."
1368 << "\n";
1369 std::cout << "Waves in star " << is
1370 << " (starInvert ==1):" << "\n";
1371 for (j = begin1; j < end1; ++j) {
1372 std::cout << waves_[j].indicesBz << " "
1373 << waves_[j].coeff << "\n";
1374 }
1375 std::cout << "Waves in star " << is+1
1376 << " (starInvert == -1):" << "\n";
1377 for (j=begin2; j < end2; ++j) {
1378 std::cout << waves_[j].indicesBz << " "
1379 << waves_[j].coeff << "\n";
1380 }
1381 return false;
1382 }
1383 }
1384 }
1385
1386 // Finished processing a pair, increment star counter by 2
1387 is += 2;
1388
1389 } // end if (stars_[is].invertFlag == 0) ... else ...
1390
1391 } // end while (is < nStar_) loop over stars
1392
1393 // Loop over basis functions
1394 for (ib = 0; ib < nBasis_; ++ib) {
1395 is = starIds_[ib];
1396 if (stars_[is].cancel) {
1397 std::cout << "\n";
1398 std::cout << "Star referred to by starIds_ is cancelled\n";
1399 return false;
1400 }
1401 if (stars_[is].basisId != ib) {
1402 std::cout << "\n";
1403 std::cout << "Error: stars_[starIds_[ib]].basisId != ib\n";
1404 std::cout << "Basis function index ib = " << ib << "\n";
1405 std::cout << "is = starIds_[ib] = " << is << "\n";
1406 std::cout << "stars_[is].basisId = "
1407 << stars_[is].basisId << "\n";
1408 return false;
1409 }
1410 }
1411
1412 // The end of this function is reached iff all tests passed.
1413 return true;
1414 }
1415
1416}
1417}
1418#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Iterator over points in a Mesh<D>.
void begin()
Set iterator to the first point in the mesh.
bool atEnd() const
Is this the end (i.e., one past the last point)?
IntVec< D > position() const
Get current position in the grid, as integer vector.
Description of a regular grid of points in a periodic domain.
IntVec< D > dimensions() const
Get an IntVec<D> of the grid dimensions.
Definition Mesh.h:217
int size() const
Get total number of grid points.
Definition Mesh.h:229
Symmetry-adapted Fourier basis for pseudo-spectral scft.
Definition fieldIoUtil.h:21
int nBasis() const
Total number of nonzero symmetry-adapted basis functions.
Definition Basis.tpp:923
Basis()
Default constructor.
Definition Basis.tpp:31
~Basis()
Destructor.
Definition Basis.tpp:49
void makeBasis(Mesh< D > const &mesh, UnitCell< D > const &unitCell, SpaceGroup< D > const &group)
Construct basis for a specific mesh and space group.
Definition Basis.tpp:69
void outputWaves(std::ostream &out, bool outputAll=false) const
Print a list of all waves to an output stream.
Definition Basis.tpp:928
bool isValid() const
Returns true if valid, false otherwise.
Definition Basis.tpp:990
void outputStars(std::ostream &out, bool outputAll=false) const
Print a list of all stars to an output stream.
Definition Basis.tpp:957
Crystallographic space group.
Definition FieldIoReal.h:23
void checkMeshDimensions(IntVec< D > const &dimensions) const
Check if input mesh dimensions are compatible with space group.
Base template for UnitCell<D> classes, D=1, 2 or 3.
Definition rpg/System.h:34
Vec< D, T > & negate(const Vec< D, T > &v)
Return negative of vector v.
Definition Vec.h:520
static const double Pi
Trigonometric constant Pi.
Definition Constants.h:35
Wrapper for a double precision number, for formatted ostream output.
Definition Dbl.h:40
An automatically growable array, analogous to a std::vector.
Definition GArray.h:34
int size() const
Return logical size of this array (i.e., current number of elements).
Definition GArray.h:455
void clear()
Reset to empty state.
Definition GArray.h:299
void append(Data const &data)
Append an element to the end of the sequence.
Definition GArray.h:306
Wrapper for an int, for formatted ostream output.
Definition Int.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
RT imag(CT const &a)
Return the imaginary part of a complex number.
RT real(CT const &a)
Return the real part of a complex number.
PSCF package top-level namespace.
Definition param_pc.dox:1
Simple wave struct for use within Basis construction.
Definition TWave.h:24