Hamiltonians::J1J2XXZModel_SquareLattice Class Reference

Detailed Description

Implements the Hamiltonian for the J1-J2 XY model on the square lattice.

Options Database
Command line argument Description
-Lx [int] lattice dimension in the longitudinal direction (growing) [def: 4]
-Ly [int] lattice dimension in the transverse direction [def: 4]
-J1 [float] coupling constant for the nearest neighbor interaction
-Jz1 [float] anisotropy in the z-direction for the nearest neighbor interaction
-J2 [float] coupling constant for the next-nearest neighbor interaction
-Jz2 [float] anisotropy in the z-direction for the next-nearest neighbor interaction
-heisenberg [float] Overrides other coupling constants and reverts to an XXZ Heisenberg model where the given parameter is the coupling constant along the z-direction
-BCopen Either for open boundary conditions or for periodic (toroidal) boundary conditions. Only at most one of them can be true. If none of them are set, the default cylindrical boundary conditions are used.
-BCperiodic
See also
Options regarding the DMRG runtime can be found at the DMRGBlockContainer documentation.

Definition at line 77 of file Hamiltonians.hpp.

#include <Hamiltonians.hpp>

Public Member Functions

 J1J2XXZModel_SquareLattice ()
 Constructor. More...
 
PetscErrorCode SetFromOptions ()
 Set attributes from command line arguments. More...
 
PetscErrorCode SaveAsOptions (const std::string &filename)
 Dumps options keys and values to file. More...
 
PetscInt NumSites () const
 Returns the number of sites in the square lattice. More...
 
PetscInt NumEnvSites () const
 Returns the number of sites in a single cluster/column of the environment block as suggested by Liang and Pang, 1994, for the square lattice. More...
 
std::vector< TermH (const PetscInt &nsites)
 Returns the object used to construct the Hamiltonian using sites counted from the left side of the lattice. More...
 
void PrintOut () const
 Prints out some information to stdout. More...
 
void SaveOut (FILE *fp) const
 Writes out some JSON information to stdout. More...
 
PetscInt Lx () const
 
PetscInt Ly () const
 
PetscInt To1D (const PetscInt ix, const PetscInt jy) const
 
PetscErrorCode To2D (const PetscInt idx, PetscInt &ix, PetscInt &jy) const
 
std::vector< std::vector< PetscInt > > NeighborPairs (const PetscInt d=1) const
 Returns all pairs of nearest-neighbors in the full square lattice. More...
 

Private Member Functions

PetscErrorCode PrintInfo () const
 Prints some information related to the parameters of the Hamiltonian. More...
 
std::vector< PetscInt > GetNearestNeighbors (const PetscInt &ix, const PetscInt &jy, const PetscInt &nsites_in) const
 Given the 2d lattice coordinates, this function retrieves the nearest neighbors' 1D coordinates. More...
 
std::vector< PetscInt > GetNextNearestNeighbors (const PetscInt &ix, const PetscInt &jy, const PetscInt &nsites_in) const
 Given the 2d lattice coordinates, this function retrieves the next-nearest neighbors' 1D coordinates. More...
 

Private Attributes

PetscBool heisenberg = PETSC_FALSE
 Whether to use the Heisenberg model. More...
 
PetscBool set_from_options = PETSC_FALSE
 Set from options. More...
 
PetscScalar _Jz1 = 0.0
 Transverse anisotropy for nearest-neighbor interactions. More...
 
PetscScalar _Jz2 = 0.0
 Transverse anisotropy for next-nearest-neighbor interactions. More...
 
PetscScalar _J1 = 1.0
 Coupling strength for nearest-neighbor interactions. More...
 
PetscScalar _J2 = 1.0
 Coupling strength for next-nearest-neighbor interactions. More...
 
PetscInt _Lx = 4
 Length along (growing) longitudinal direction. More...
 
PetscInt _Ly = 4
 Length along (fixed) transverse direction. More...
 
PetscBool verbose = PETSC_FALSE
 Flag that tells whether to print some information. More...
 
BC_t _BCx = OpenBC
 Boundary conditions for the longitudinal direction (defaults to cylindrical BC) More...
 
BC_t _BCy = PeriodicBC
 Boundary conditions for the transverse direction (defaults to cylindrical BC) More...
 
std::vector< TermH_full
 Contains the Hamiltonian involving the entire lattice. More...
 
PetscBool H_full_filled = PETSC_FALSE
 Tells whether H_full has been constructed. More...
 

The documentation for this class was generated from the following files:
This site was generated by Sphinx using Doxygen with a customized theme from doxygen-bootstrapped.