Skip to content

Commit

Permalink
merged development into this branch -- resolved conclicts from decima…
Browse files Browse the repository at this point in the history
…tion
  • Loading branch information
saurabh-s-sawant committed Dec 13, 2024
2 parents 906fc0c + a8a690c commit a0b0dda
Show file tree
Hide file tree
Showing 14 changed files with 835 additions and 232 deletions.
1 change: 1 addition & 0 deletions Source/Code_Definitions.H
Original file line number Diff line number Diff line change
Expand Up @@ -6,6 +6,7 @@
#define VFRAC_THREASHOLD 1e-5

#define NUM_MODES 1
#define BLOCK_SIZE 1
#define NUM_CONTACTS 2
#define NUM_ENERGY_PTS_REAL 10

Expand Down
7 changes: 4 additions & 3 deletions Source/Solver/Transport/NEGF/CNT.H
Original file line number Diff line number Diff line change
Expand Up @@ -10,9 +10,9 @@

#include "NEGF_Common.H"

class c_CNT : public c_NEGF_Common<ComplexType[NUM_MODES]>
class c_CNT : public c_NEGF_Common<ComplexType[BLOCK_SIZE]>
{
using BlkType = ComplexType[NUM_MODES];
using BlkType = ComplexType[BLOCK_SIZE];

amrex::Real acc = 1.42e-10;
amrex::Real gamma = 2.5;
Expand Down Expand Up @@ -75,7 +75,8 @@ class c_CNT : public c_NEGF_Common<ComplexType[NUM_MODES]>
virtual int get_offDiag_repeatBlkSize() final { return 2; }

virtual AMREX_GPU_HOST_DEVICE void Compute_SurfaceGreensFunction(
MatrixBlock<BlkType> &gr, const ComplexType EmU) final;
MatrixBlock<BlkType> &gr, const ComplexType E,
const ComplexType U) final;

public:
static int get_1D_site_id(int par_id_local_to_NS)
Expand Down
59 changes: 35 additions & 24 deletions Source/Solver/Transport/NEGF/CNT.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -66,8 +66,8 @@ void c_CNT::Set_MaterialSpecificParameters()

void c_CNT::Set_BlockDegeneracyVector(amrex::Vector<int> &vec)
{
vec.resize(NUM_MODES);
for (int m = 0; m < NUM_MODES; ++m)
vec.resize(BLOCK_SIZE);
for (int m = 0; m < BLOCK_SIZE; ++m)
{
vec[m] = mode_degen_vec[m];
}
Expand Down Expand Up @@ -293,7 +293,8 @@ ComplexType c_CNT::get_beta(int J)

void c_CNT::Define_MPI_BlkType()
{
MPI_Type_vector(1, NUM_MODES, NUM_MODES, MPI_DOUBLE_COMPLEX, &MPI_BlkType);
MPI_Type_vector(1, BLOCK_SIZE, BLOCK_SIZE, MPI_DOUBLE_COMPLEX,
&MPI_BlkType);
MPI_Type_commit(&MPI_BlkType);
}

Expand All @@ -310,7 +311,7 @@ void c_CNT::Construct_Hamiltonian()
h_minusHa(i) = 0.;
}

for (int j = 0; j < NUM_MODES; ++j)
for (int j = 0; j < BLOCK_SIZE; ++j)
{
int J = mode_vec[j];
beta.block[j] = get_beta(J);
Expand Down Expand Up @@ -352,32 +353,42 @@ void c_CNT::Define_ContactInfo()

AMREX_GPU_HOST_DEVICE
void c_CNT::Compute_SurfaceGreensFunction(MatrixBlock<BlkType> &gr,
const ComplexType EmU)
const ComplexType E,
const ComplexType U)
{
auto EmU_sq = pow(EmU, 2.);
auto gamma_sq = pow(gamma, 2.);

for (int i = 0; i < NUM_MODES; ++i)
if (use_decimation)
{
c_NEGF_Common<BlkType>::DecimationTechnique(gr, E - U);
}
else
{
auto Factor = EmU_sq + gamma_sq - pow(beta.block[i], 2);
ComplexType EmU = E - U;
auto EmU_sq = pow(EmU, 2.);
auto gamma_sq = pow(gamma, 2.);

auto Sqrt = sqrt(pow(Factor, 2) - 4. * EmU_sq * gamma_sq);
for (int i = 0; i < BLOCK_SIZE; ++i)
{
auto Factor = EmU_sq + gamma_sq - pow(beta.block[i], 2);

auto Sqrt = sqrt(pow(Factor, 2) - 4. * EmU_sq * gamma_sq);

auto Denom = 2. * gamma_sq * EmU;
auto Denom = 2. * gamma_sq * EmU;

auto val1 = (Factor + Sqrt) / Denom;
auto val2 = (Factor - Sqrt) / Denom;
auto val1 = (Factor + Sqrt) / Denom;
auto val2 = (Factor - Sqrt) / Denom;

if (val1.imag() < 0.)
gr.block[i] = val1;
else if (val2.imag() < 0.)
gr.block[i] = val2;
if (val1.imag() < 0.)
gr.block[i] = val1;
else if (val2.imag() < 0.)
gr.block[i] = val2;

// amrex::Print() << "EmU: " << EmU << "\n";
// amrex::Print() << "Factor: " << Factor << "\n";
// amrex::Print() << "Sqrt: " << Sqrt << "\n";
// amrex::Print() << "Denom: " << Denom << "\n";
// amrex::Print() << "Numerator: " << Factor+Sqrt << "\n";
// amrex::Print() << "Value: " << (Factor+Sqrt)/Denom << "\n";
// amrex::Print() << "EmU: " << EmU << "\n";
// amrex::Print() << "Factor: " << Factor << "\n";
// amrex::Print() << "Sqrt: " << Sqrt << "\n";
// amrex::Print() << "Denom: " << Denom << "\n";
// amrex::Print() << "Numerator: " << Factor+Sqrt << "\n";
// amrex::Print() << "Value: " << (Factor+Sqrt)/Denom << "\n";
}
// amrex::Print() << "Using quadratic, gr: " << gr << "\n";
}
}
4 changes: 2 additions & 2 deletions Source/Solver/Transport/NEGF/Graphene.H
Original file line number Diff line number Diff line change
Expand Up @@ -10,9 +10,9 @@

#include "NEGF_Common.H"

class c_Graphene : public c_NEGF_Common<ComplexType[NUM_MODES][NUM_MODES]>
class c_Graphene : public c_NEGF_Common<ComplexType[BLOCK_SIZE][BLOCK_SIZE]>
{
using BlkType = ComplexType[NUM_MODES][NUM_MODES];
using BlkType = ComplexType[BLOCK_SIZE][BLOCK_SIZE];

static amrex::Array<int, 2> type_id;

Expand Down
45 changes: 44 additions & 1 deletion Source/Solver/Transport/NEGF/Matrix_Block.H
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,44 @@ struct MatrixBlock
public:
T block;

// Non-const .data() for 1D array
template <typename U = T>
auto data() ->
typename std::enable_if<std::rank<U>::value == 1,
typename std::remove_extent<U>::type *>::type
{
return &block[0];
}

// Const .data() for 1D array
template <typename U = T>
auto const_data() const -> typename std::enable_if<
std::rank<U>::value == 1,
const typename std::remove_extent<U>::type *>::type
{
return &block[0];
}

// Non-const .data() for 2D array
template <typename U = T>
auto data() -> typename std::enable_if<
std::rank<U>::value == 2,
typename std::remove_all_extents<U>::type *>::type
{
return reinterpret_cast<typename std::remove_all_extents<U>::type *>(
&block[0][0]);
}

// Const .data() for 2D array
template <typename U = T>
auto const_data() const -> typename std::enable_if<
std::rank<U>::value == 2,
const typename std::remove_all_extents<U>::type *>::type
{
return reinterpret_cast<
const typename std::remove_all_extents<U>::type *>(&block[0][0]);
}

AMREX_GPU_HOST_DEVICE MatrixBlock<T> operator=(const amrex::Real c);
AMREX_GPU_HOST_DEVICE MatrixBlock<T> operator=(const ComplexType c);

Expand Down Expand Up @@ -63,7 +101,11 @@ struct MatrixBlock
AMREX_GPU_HOST_DEVICE MatrixBlock<T> Conj() const; /*Conjugate*/
AMREX_GPU_HOST_DEVICE MatrixBlock<T> Tran() const; /*Transpose*/
AMREX_GPU_HOST_DEVICE MatrixBlock<T> Dagger() const; /*Conjugate-Transpose*/
AMREX_GPU_HOST_DEVICE ComplexType DiagSum() const; /*Trace*/
AMREX_GPU_HOST_DEVICE MatrixBlock<T> Inverse() const; /*Inverse*/
AMREX_GPU_HOST_DEVICE ComplexType DiagSum() const; /*Trace*/
AMREX_GPU_HOST_DEVICE ComplexType FrobeniusNorm() const;
AMREX_GPU_HOST_DEVICE void SetDiag(
const ComplexType c); /*set diagonal element*/

template <typename U>
AMREX_GPU_HOST_DEVICE ComplexType
Expand All @@ -75,6 +117,7 @@ struct MatrixBlock
// MatrixBlock<T> AtomicAdd(MatrixBlock<T>& B) const; /*Atomic add for the
// block*/ template<typename U, std::size_t N> MatrixBlock<amrex::Real>
// Real() const;
//
};

#endif
Loading

0 comments on commit a0b0dda

Please sign in to comment.