StatMech
Loading...
Searching...
No Matches
Ensemble.cuh File Reference

Classes

class  SubSpaceManager
 
class  ScalableOperator< MatrixType >
 
class  HamiltonianEnsemble
 
class  ObservableEnsemble
 

Typedefs

using RealScalar = Real_t
 
using Scalar = Complex_t< RealScalar >
 

Functions

template<class Derived1 , class Derived2 , class TotalSpace , typename Scalar >
__global__ void construct_globalOp_kernel (Eigen::MatrixBase< Derived1 > *dTotMat, Eigen::MatrixBase< Derived2 > const *dLocMat, SubSpace< TotalSpace, Scalar > const *subSpace)
 
template<class Derived1 , class Derived2 , class TotalSpace >
__global__ void ObservableEnsemble_sample_kernel (Eigen::MatrixBase< Derived1 > *resPtr, Eigen::MatrixBase< Derived2 > const *opTotPtr, SubSpace< TotalSpace, typename Derived1::Scalar > const *subSpacePtr)
 

Variables

int g_dimLoc = 2
 
int g_momentum = 0
 

Typedef Documentation

◆ RealScalar

using RealScalar = Real_t

◆ Scalar

Function Documentation

◆ construct_globalOp_kernel()

template<class Derived1 , class Derived2 , class TotalSpace , typename Scalar >
__global__ void construct_globalOp_kernel ( Eigen::MatrixBase< Derived1 > *  dTotMat,
Eigen::MatrixBase< Derived2 > const *  dLocMat,
SubSpace< TotalSpace, Scalar > const *  subSpace 
)
66 {
67 int const idx = blockIdx.x * blockDim.x + threadIdx.x;
68 int const idy = blockIdx.y * blockDim.y + threadIdx.y;
69 if(idx >= dTotMat->rows() || idy >= dTotMat->cols()) return;
70 if(idx < idy) return;
71
72 int const dimLoc = dLocMat->cols();
73 int const dimTot = subSpace->basis().outerIndexPtr()[dTotMat->rows()];
74
75 (*dTotMat)(idx, idy) = {0, 0};
76 int innerX, innerY, locIdX, locIdY;
77 for(auto valuePosX = subSpace->basis().outerIndexPtr()[idx];
78 valuePosX != subSpace->basis().outerIndexPtr()[idx + 1]; ++valuePosX) {
79 innerX = subSpace->basis().innerIndexPtr()[valuePosX];
80 auto conjCoeffX = conj(subSpace->basis().valuePtr()[valuePosX]);
81
82 for(auto valuePosY = subSpace->basis().outerIndexPtr()[idy];
83 valuePosY != subSpace->basis().outerIndexPtr()[idy + 1]; ++valuePosY) {
84 innerY = subSpace->basis().innerIndexPtr()[valuePosY];
85 auto CoeffY = subSpace->basis().valuePtr()[valuePosY];
86
87 if(innerX / dimLoc != innerY / dimLoc) continue;
88 locIdX = innerX % dimLoc;
89 locIdY = innerY % dimLoc;
90 (*dTotMat)(idx, idy) += conjCoeffX * (*dLocMat)(locIdX, locIdY) * CoeffY;
91 }
92 }
93
94 if(idx == idy) { (*dTotMat)(idx, idy) = real((*dTotMat)(idx, idy)); }
95 else { (*dTotMat)(idy, idx) = conj((*dTotMat)(idx, idy)); }
96 // printf("LDT=%d, id=(%d,%d), val=%+f%+f*i\t %s is yet to be tested.\n", LDT, idx, idy,
97 // real(dTotMat[idx + LDT * idy]), imag(dTotMat[idx + LDT * idy]), __func__);
98 // printf("\t id=(%d,%d)\n", idx, idy);
99 return;
100}
__host__ __device__ int * innerIndexPtr() const
Definition MatrixUtils.cuh:421
__host__ __device__ int * outerIndexPtr() const
Definition MatrixUtils.cuh:420
__host__ __device__ Scalar_t * valuePtr() const
Definition MatrixUtils.cuh:422
__host__ __device__ Matrix_t & basis()
Definition HilbertSpace.hpp:668

◆ ObservableEnsemble_sample_kernel()

template<class Derived1 , class Derived2 , class TotalSpace >
__global__ void ObservableEnsemble_sample_kernel ( Eigen::MatrixBase< Derived1 > *  resPtr,
Eigen::MatrixBase< Derived2 > const *  opTotPtr,
SubSpace< TotalSpace, typename Derived1::Scalar > const *  subSpacePtr 
)
327 {
328 int const idx = blockIdx.x * blockDim.x + threadIdx.x;
329 int const idy = blockIdx.y * blockDim.y + threadIdx.y;
330 auto& basis = subSpacePtr->basis();
331 if(idx >= resPtr->rows() || idy >= resPtr->cols()) return;
332 if(idx > idy) return;
333
334 (*resPtr)(idx, idy) = {0, 0};
335 int innerX, innerY;
336 for(auto valuePosX = basis.outerIndexPtr()[idx]; valuePosX != basis.outerIndexPtr()[idx + 1];
337 ++valuePosX) {
338 innerX = basis.innerIndexPtr()[valuePosX];
339 auto conjCoeffX = conj(basis.valuePtr()[valuePosX]);
340
341 for(auto valuePosY = basis.outerIndexPtr()[idy];
342 valuePosY != basis.outerIndexPtr()[idy + 1]; ++valuePosY) {
343 innerY = basis.innerIndexPtr()[valuePosY];
344 auto CoeffY = basis.valuePtr()[valuePosY];
345
346 (*resPtr)(idx, idy) += conjCoeffX * (*opTotPtr)(innerX, innerY) * CoeffY;
347 }
348 }
349
350 if(idx == idy) (*resPtr)(idy, idx) = real((*resPtr)(idy, idx));
351
352 (*resPtr)(idy, idx) = conj((*resPtr)(idx, idy));
353}

Variable Documentation

◆ g_dimLoc

int g_dimLoc = 2

◆ g_momentum

int g_momentum = 0