oILAB
Loading...
Searching...
No Matches
LatticeFunctionImplementation.h
Go to the documentation of this file.
1//
2// Created by Nikhil Chandra Admal on 5/31/24.
3//
4#ifndef OILAB_LATTICEFUNCTIONIMPLEMENTATION_H
5#define OILAB_LATTICEFUNCTIONIMPLEMENTATION_H
6#include <FFT.h>
7#include <iostream>
8
9namespace gbLAB
10{
11 template<typename Scalar, int dim>
12 LatticeFunction<Scalar,dim>::LatticeFunction(const Eigen::array<Eigen::Index,dim>& n,
13 const Eigen::Matrix<double,Eigen::Dynamic,dim>& _basisVectors) :
14 values(n), basisVectors(_basisVectors)
15 {
16 values.setZero();
17 }
18
19 template<typename Scalar, int dim>
20 template<typename T, typename, typename, int dm, typename>
21 LatticeFunction<Scalar,dim>::LatticeFunction(const Eigen::array<Eigen::Index, dim> &n,
22 const Eigen::Matrix<double, Eigen::Dynamic, dim> &_basisVectors,
23 const Function<T,Scalar>& fun) :
24 values(n), basisVectors(_basisVectors) {
25 for (int i = 0; i < n[0]; i++) {
26 int in= i > n[0]/2 ? i-n[0] : i;
27 values(i) = fun(in * basisVectors.col(0));
28 }
29 }
30
31 template<typename Scalar, int dim>
32 template<typename T, typename, int dm, typename>
33 LatticeFunction<Scalar,dim>::LatticeFunction(const Eigen::array<Eigen::Index, dim> &n,
34 const Eigen::Matrix<double,Eigen::Dynamic,dim>& _basisVectors,
35 const Function<T,Scalar>& fun) :
36 values(n), basisVectors(_basisVectors) {
37 for (int i = 0; i < n[0]; i++) {
38 for (int j = 0; j < n[1]; j++) {
39 //values(i, j) = fun(i * basisVectors.col(0) + j * basisVectors.col(1));
40 int in= i > n[0]/2 ? i-n[0] : i;
41 int jn= j > n[1]/2 ? j-n[1] : j;
42 values(i,j) = fun(in * basisVectors.col(0) + jn * basisVectors.col(1));
43 }
44 }
45 }
46
47 template<typename Scalar, int dim>
48 template<typename T, int dm, typename>
49 LatticeFunction<Scalar,dim>::LatticeFunction(const Eigen::array<Eigen::Index, dim> &n,
50 const Eigen::Matrix<double, Eigen::Dynamic, dim> &_basisVectors,
51 const Function<T, Scalar> &fun) :
52 values(n), basisVectors(_basisVectors) {
53 for (int i = 0; i < n[0]; i++) {
54 for (int j = 0; j < n[1]; j++) {
55 for (int k = 0; k < n[2]; k++) {
56 //values(i, j, k) = fun(i * basisVectors.col(0) + j * basisVectors.col(1) + k * basisVectors.col(2));
57 int in= i>n[0]/2 ? i-n[0] : i;
58 int jn= j>n[1]/2 ? j-n[1] : j;
59 int kn= k>n[2]/2 ? k-n[2] : k;
60 values(i, j, k) = fun(in * basisVectors.col(0) +
61 jn * basisVectors.col(1) +
62 kn * basisVectors.col(2));
63 }
64 }
65 }
66 }
67
68 template<typename Scalar, int dim>
69 std::complex<double> LatticeFunction<Scalar,dim>::dot(const LatticeFunction<std::complex<double>,dim>& other) const
70 {
71 Eigen::Tensor<std::complex<double>,0> sum((this->values * other.values.conjugate()).sum());
72 Eigen::Matrix<double,dim,dim> gramMatrix;
73 for(int i=0; i<dim; ++i)
74 for(int j=0; j<dim; ++j)
75 gramMatrix(i,j)= basisVectors.col(i).dot(basisVectors.col(j));
76 return sum(0) * sqrt(gramMatrix.determinant());
77
78 }
79
80 template<typename Scalar, int dim>
82 {
83 Eigen::Matrix<double,Eigen::Dynamic,dim> unitCell(basisVectors.transpose().completeOrthogonalDecomposition().pseudoInverse());
84 PeriodicFunction<dcomplex,dim> pf(values.dimensions(),unitCell);
85 FFT::ifft(values.template cast<dcomplex>(),pf.values);
86 //return pf;
87 // Calculate the area spanned by the unit cell vectors
88 Eigen::Matrix<double,dim,dim> unitCellGramMatrix;
89 for(int i=0; i<dim; ++i)
90 for(int j=0; j<dim; ++j)
91 unitCellGramMatrix(i,j)= unitCell.col(i).dot(unitCell.col(j));
92
93 Eigen::array<Eigen::Index,dim> n= this->values.dimensions();
94 int prod= std::accumulate(std::begin(n),
95 std::begin(n) + dim,
96 1,
97 std::multiplies<>{});
98 pf.values= pf.values * (dcomplex)(prod/sqrt(unitCellGramMatrix.determinant()));
99 return pf;
100 }
101
102 template<typename Scalar, int dim>
104 {
105 // assert that the two lattice functions have the same domain
106 const auto n= lf1.values.dimensions();
108 output.values= lf1.values * lf2.values;
109 return output;
110 }
111}
112#endif // OILAB_LATTICEFUNCTIONIMPLEMENTATION_H
static void ifft(const Eigen::Tensor< dcomplex, 3 > &in, Eigen::Tensor< dcomplex, 3 > &out)
Definition FFT.h:54
Eigen::Tensor< Scalar, dim > values
PeriodicFunction< dcomplex, dim > ifft() const
const Eigen::Matrix< double, Eigen::Dynamic, dim > basisVectors
std::complex< double > dcomplex
LatticeFunction(const Eigen::array< Eigen::Index, dim > &n, const Eigen::Matrix< double, Eigen::Dynamic, dim > &_basisVectors)
std::complex< double > dot(const LatticeFunction< std::complex< double >, dim > &other) const
Eigen::Tensor< Scalar, dim > values
LatticeVector< dim > operator*(const typename LatticeVector< dim >::IntScalarType &scalar, const LatticeVector< dim > &L)