PSCF v1.1
SpaceSymmetry.cpp
1/*
2* PSCF - Polymer Self-Consistent Field Theory
3*
4* Copyright 2016 - 2022, The Regents of the University of Minnesota
5* Distributed under the terms of the GNU General Public License.
6*/
7
8#include "SpaceSymmetry.tpp"
9
10namespace Pscf
11{
12
13 using namespace Util;
14
15 // Explicit specializations D = 1
16
17 template<>
19 { return R_(0, 0); }
20
21 template<>
23 {
24 // Compute and check the determinant
25 int det = R_(0,0);
26 if (1 != det*det) {
27 std::cout << "Determinant (1D symmetry) =" << det << std::endl;
28 UTIL_THROW("Invalid SpaceSymmetry<1>: |det| != 1");
29 }
30
31 // Compute the inverse matrix
33 A(0,0) = R_(0,0);
34
35 return A;
36 }
37
38 // Explicit specializations D = 2
39
40 template<>
42 { return (R_(0,0)*R(1,1) - R(0,1)*R(1,0)); }
43
44 template<>
46 {
47 // Compute and check the determinant
48 int det = determinant();
49 if (1 != det*det) {
50 UTIL_THROW("Invalid SpaceSymmetry<2>: |det| != 1");
51 }
52
53 // Compute the inverse matrix
55 A(0,0) = R_(1,1)/det;
56 A(1,1) = R_(0,0)/det;
57 A(0,1) = -R_(0,1)/det;
58 A(1,0) = -R_(1,0)/det;
59
60 return A;
61 }
62
63 // Explicit specializations D = 3
64
65 template<>
67 {
68 int det = 0;
69 det += R_(0,0)*(R(1,1)*R(2,2) - R(1,2)*R(2,1));
70 det += R_(0,1)*(R(1,2)*R(2,0) - R(1,0)*R(2,2));
71 det += R_(0,2)*(R(1,0)*R(2,1) - R(1,1)*R(2,0));
72 return det;
73 }
74
75 template<>
77 {
78 // Compute and check the determinant
79 int det = determinant();
80 if (1 != det*det) {
81 UTIL_THROW("Invalid SpaceSymmetry<3>: |det| != 1");
82 }
83
84 // Compute the inverse matrix
86 A(0,0) = (R_(1,1)*R(2,2) - R(1,2)*R(2,1))/det;
87 A(1,0) = (R_(1,2)*R(2,0) - R(1,0)*R(2,2))/det;
88 A(2,0) = (R_(1,0)*R(2,1) - R(1,1)*R(2,0))/det;
89
90 A(0,1) = (R_(2,1)*R(0,2) - R(2,2)*R(0,1))/det;
91 A(1,1) = (R_(2,2)*R(0,0) - R(2,0)*R(0,2))/det;
92 A(2,1) = (R_(2,0)*R(0,1) - R(2,1)*R(0,0))/det;
93
94 A(0,2) = (R_(0,1)*R(1,2) - R(0,2)*R(1,1))/det;
95 A(1,2) = (R_(0,2)*R(1,0) - R(0,0)*R(1,2))/det;
96 A(2,2) = (R_(0,0)*R(1,1) - R(0,1)*R(1,0))/det;
97
98 return A;
99 }
100
101 // Explicit instantiation of required class instances
102
103 template class SpaceSymmetry<1>;
104 template class SpaceSymmetry<2>;
105 template class SpaceSymmetry<3>;
106
107}
int determinant() const
Compute and return the determinant of the rotation matrix.
SpaceSymmetry< D >::Rotation inverseRotation() const
Compute and return the inverse of the rotation matrix.
FMatrix< int, D, D > Rotation
Typedef for matrix used to represent point group operation.
#define UTIL_THROW(msg)
Macro for throwing an Exception, reporting function, file and line number.
Definition: global.h:51
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1