Class wavefunc_class
In: MeasureClass/wavefunc_class.F90
comlib lattice_class wavefunc_param_class error_class print_status_class wavefunc_class dot/f_79.png

This defines wavefunction for quark propagotors Several wavefunction with Coulomb-gauge fixing are defined.

Methods

EXP_WAV   LOCAL_WAV   MAXPOL_WAV   POLY_WAV   USER_WAV   WALL_WAV   delete   delete   get   get   imom   ispace   new   new   print   print   read   snk_wavefunc_obj   src_wavefunc_obj   wavefunc_param_obj  

Included Modules

comlib lattice_class wavefunc_param_class error_class print_status_class

Public Instance methods

EXP_WAV
Constant :
EXP_WAV =2 :integer, parameter

Original external subprogram is wavefunc_param_class#EXP_WAV

LOCAL_WAV
Constant :
LOCAL_WAV =1 :integer, parameter

Original external subprogram is wavefunc_param_class#LOCAL_WAV

MAXPOL_WAV
Constant :
MAXPOL_WAV =50 :integer, parameter

Original external subprogram is wavefunc_param_class#MAXPOL_WAV

POLY_WAV
Constant :
POLY_WAV =4 :integer, parameter

Original external subprogram is wavefunc_param_class#POLY_WAV

USER_WAV
Constant :
USER_WAV =5 :integer, parameter

Original external subprogram is wavefunc_param_class#USER_WAV

WALL_WAV
Constant :
WALL_WAV =3 :integer, parameter

Original external subprogram is wavefunc_param_class#WALL_WAV

Subroutine :
this :type(snk_wavefunc_obj), intent(inout)

[Source]

subroutine delete_snk_wavefunc(this)
  implicit none
  type(snk_wavefunc_obj), intent(inout) :: this

  call delete(this%param)
  if (allocated(this%wavefunc)) deallocate(this%wavefunc)

  return
end subroutine
Subroutine :
this :type(src_wavefunc_obj), intent(inout)

[Source]

subroutine delete_src_wavefunc(this)
  implicit none
  type(src_wavefunc_obj), intent(inout) :: this

  call delete(this%param)

  if (allocated(this%wavefunc)) deallocate(this%wavefunc)

  return
end subroutine
Subroutine :
this :type(snk_wavefunc_obj), intent(inout)

[Source]

subroutine get_snk_wavefunc(this)
  use print_status_class
  implicit none
  type(snk_wavefunc_obj), intent(inout) :: this
  type(src_wavefunc_obj), allocatable :: twf
  integer :: itt,isp,it,it0
  integer :: ix,iy,iz
  integer :: itx,ity,itz
  character(len=CHARLEN) :: str

  if (allocated(this%wavefunc)) deallocate(this%wavefunc)
  allocate(this%wavefunc(NSPACE))

  allocate(twf)
  call new(twf)
  twf%param = this%param
  twf%itx0 = 1
  twf%ity0 = 1
  twf%itz0 = 1
  twf%itt0 = 1

  call get(twf)
  do it=1,NT
    itt = it
    if (itt == twf%itt0) then
!$OMP PARALLEL DO PRIVATE(isp)
      do isp=1,NSPACE
        this%wavefunc(isp) = twf%wavefunc(1,1,isp,it)
      enddo
    endif
  enddo

#ifdef _DDD_
!============
do ix=1,NX
do iy=1,NY
do iz=1,NZ
itx = ix + ipsite(1)*NX
ity = iy + ipsite(2)*NY
itz = iz + ipsite(3)*NZ
isp=ispace(ix,iy,iz)
write(str,'(3I3," :",I9,I2,2E24.15,I11," SWF")') itx,ity,itz,isp,1,this%wavefunc(isp), itx-1+(ity-1)*NTX+(itz-1)*NTX*NTY
call print_status(str)
enddo
enddo
enddo
!============
#endif


  call delete(twf)
  deallocate(twf)
  return
end subroutine
Subroutine :
this :type(src_wavefunc_obj), intent(inout)

[Source]

subroutine get_src_wavefunc(this)
  use print_status_class
  implicit none
  type(src_wavefunc_obj), intent(inout) :: this
  integer  :: it,iz,iy,ix,ic,is,jc,js,isp
  integer  :: itt,itz,ity,itx
  integer  :: rx,ry,rz,ipol
  integer  :: ipx,ipy,ipz,ipt
  integer  :: ix0,iy0,iz0,it0
  real(DP) :: rsrc,vsrc,rmax
  character(len=CHARLEN) :: str

  if ( this%param%type == USER_WAV) return

  if (allocated(this%wavefunc)) deallocate(this%wavefunc)
  allocate(this%wavefunc(COL,SPIN,NSPACE,NT))
!$OMP PARALLEL WORKSHARE
  this%wavefunc(:,:,:,:) = Z0
!$OMP END PARALLEL WORKSHARE

  select case (this%param%type)
  case (LOCAL_WAV)

    ipt = (this%itt0-1)/NT
    ipz = (this%itz0-1)/NZ
    ipy = (this%ity0-1)/NY
    ipx = (this%itx0-1)/NX
    if ( (ipt == 0)         .and. (ipz == ipsite(3)) .and. (ipy == ipsite(2)) .and. (ipx == ipsite(1)) ) then
      it0 = mod(this%itt0-1,NT)+1
      iz0 = mod(this%itz0-1,NZ)+1
      iy0 = mod(this%ity0-1,NY)+1
      ix0 = mod(this%itx0-1,NX)+1
      isp=ispace(ix0,iy0,iz0)
      do is=1,SPIN
      do ic=1,COL
        this%wavefunc(ic,is,isp,it0) = Z1
      enddo
      enddo
    endif

#ifdef _FULL_
    do it=1,NT
    do ix=1,NX
    do iy=1,NY
    do iz=1,NZ
      itt=it
      itz=iz+ipsite(3)*NZ
      ity=iy+ipsite(2)*NY
      itx=ix+ipsite(1)*NX
      if ( (itz==this%itz0) .and. (ity==this%ity0) .and. (itx==this%itx0) .and. (itt==this%itt0) ) then
        isp=ispace(ix,iy,iz)
        do is=1,SPIN
        do ic=1,COL
          this%wavefunc(ic,is,isp,it) = Z1
        enddo
        enddo
      endif
    enddo
    enddo
    enddo
    enddo
#endif

  case (EXP_WAV)

    rmax = MIN(NTX,NTZ,NTY)
    rmax = (rmax/2.0_DP - 1.0_DP) + 0.5_DP
    
    ipt = (this%itt0-1)/NT
    it0 = mod(this%itt0-1,NT)+1

#ifdef _DDD_
!============
write(str,'("EXP_WAVE:",4I3," :",3I2," :",2E12.4)') this%itx0,this%ity0,this%itz0,this%itt0, ipsite(1),ipsite(2),ipsite(3), this%param%Asmear,this%param%Bsmear
call print_status(str)
!============
#endif

    do it=1,NT
      itt=it
      if (itt == this%itt0) then

!$OMP PARALLEL DO PRIVATE(ix,iy,iz,itx,ity,itz,isp,rx,ry,rz,rsrc,vsrc,is,ic)
        do ix=1,NX
        do iy=1,NY
        do iz=1,NZ
          itz=iz+ipsite(3)*NZ
          ity=iy+ipsite(2)*NY
          itx=ix+ipsite(1)*NX

          isp=ispace(ix,iy,iz)

          !
          ! 2-point distance with periodic boundary condition 
          !
          rx=MIN(IABS(itx-this%itx0),IABS(itx-this%itx0+NTX),IABS(itx-this%itx0-NTX))
          ry=MIN(IABS(ity-this%ity0),IABS(ity-this%ity0+NTY),IABS(ity-this%ity0-NTY))
          rz=MIN(IABS(itz-this%itz0),IABS(itz-this%itz0+NTZ),IABS(itz-this%itz0-NTZ))
          rsrc = sqrt(real(rx**2+ry**2+rz**2,kind=KIND(rsrc)))
 
          if ( rsrc < rmax) then
            vsrc = this%param%Asmear*exp(-rsrc*this%param%Bsmear)
            do is=1,SPIN
            do ic=1,COL
              this%wavefunc(ic,is,isp,it)=cmplx(vsrc,0.0_DP,kind=KIND(vsrc))
            enddo
            enddo
          endif

          if ( (rx==0) .and. (ry==0) .and. (rz==0) ) then
            do is=1,SPIN
            do ic=1,COL
              this%wavefunc(ic,is,isp,it) = Z1
            enddo
            enddo
          endif

#ifdef _DDD_
!================
do is=1,SPIN
do ic=1,COL
write(str,'(3I3," :",I9,I2,2E24.15,I11,2I2" WWW")') itx,ity,itz,isp,it,this%wavefunc(ic,is,isp,it), itx-1+(ity-1)*NTX+(itz-1)*NTX*NTY,ic,is
call print_status(str)
enddo
enddo
!================
#endif

        enddo
        enddo
        enddo
      endif
    enddo

  case (WALL_WAV)

    vsrc = 1.0_DP/real(NTX*NTY*NTZ,kind=KIND(vsrc))

    do it=1,NT
      itt=it
      if (itt==this%itt0) then
!$OMP PARALLEL DO PRIVATE(ix,iy,iz,itx,ity,itz,isp,is,ic)
        do ix=1,NX
        do iy=1,NY
        do iz=1,NZ
          itz=iz+ipsite(3)*NZ
          ity=iy+ipsite(2)*NY
          itx=ix+ipsite(1)*NX

          isp=ispace(ix,iy,iz)

          do is=1,SPIN
          do ic=1,COL
            this%wavefunc(ic,is,isp,it) = cmplx(vsrc,0.0_DP,kind=KIND(vsrc))
          enddo
          enddo

        enddo
        enddo
        enddo
      endif
    enddo
  
  case (POLY_WAV)

    rmax = MIN(NTX,NTZ,NTY)
    rmax = (rmax/2.0_DP - 1.0_DP) + 0.5_DP
    
    do it=1,NT
      itt=it
      if (itt==this%itt0) then
!$OMP PARALLEL DO PRIVATE(ix,iy,iz,itx,ity,itz,isp,rx,ry,rz,rsrc,vsrc,ipol,is,ic)
        do ix=1,NX
        do iy=1,NY
        do iz=1,NZ
          itz=iz+ipsite(3)*NZ
          ity=iy+ipsite(2)*NY
          itx=ix+ipsite(1)*NX

          isp=ispace(ix,iy,iz)

          !
          ! 2-point distance with periodic boundary condition 
          !
          rx=MIN(IABS(itx-this%itx0),IABS(itx-this%itx0+NTX),IABS(itx-this%itx0-NTX))
          ry=MIN(IABS(ity-this%ity0),IABS(ity-this%ity0+NTY),IABS(ity-this%ity0-NTY))
          rz=MIN(IABS(itz-this%itz0),IABS(itz-this%itz0+NTZ),IABS(itz-this%itz0-NTZ))
          rsrc = sqrt(real(rx**2+ry**2+rz**2,kind=KIND(rsrc)))
 
          vsrc = 1.0_DP
          do ipol = 1,this%param%Npolsf
            vsrc = vsrc - this%param%Psmear(ipol)*rsrc**ipol
          enddo

          do is=1,SPIN
          do ic=1,COL
            this%wavefunc(ic,is,isp,it)=cmplx(vsrc,0.0_DP,kind=KIND(vsrc))
          enddo
          enddo

        enddo
        enddo
        enddo
      endif
    enddo

  end select

#ifdef _DDD_
!==============
vsrc = 0.0_DP
do ix=1,NX
do iy=1,NY
do iz=1,NZ
  isp=ispace(ix,iy,iz)
  vsrc = vsrc + REAL(this%wavefunc(1,1,isp,this%itt0),kind=DP)**2
enddo
enddo
enddo
call comlib_sumcast(vsrc)
if (0==nodeid) then
  write(*,'("|WF|^2=",E24.15)')vsrc
endif
!==============
#endif

  return
end subroutine
Function :
ipp :integer
ipx :integer
ipy :integer
ipz :integer

periodic boundary condition

[Source]

function imom(ipx,ipy,ipz) result(ipp)
!
! periodic boundary condition
!
  implicit none
  integer :: ipx,ipy,ipz,ipp
  ipp = 1 + mod(ipz-1+NZ,NZ) + mod(ipy-1+NY,NY)*NZ + mod(ipx-1+NX,NX)*NZ*NY
  return
end function
Function :
isp :integer
ix :integer
iy :integer
iz :integer

[Source]

function ispace(ix,iy,iz) result(isp)
  implicit none
  integer :: ix,iy,iz,isp
  isp=iz + (iy-1)*NZ + (ix-1)*NZ*NY
  return
end function
Subroutine :
this :type(snk_wavefunc_obj), intent(inout)

[Source]

subroutine new_snk_wavefunc(this)
  implicit none
  type(snk_wavefunc_obj), intent(inout) :: this

  call new(this%param)
  if (allocated(this%wavefunc)) deallocate(this%wavefunc)

  return
end subroutine
Subroutine :
this :type(src_wavefunc_obj), intent(inout)

[Source]

subroutine new_src_wavefunc(this)
  implicit none
  type(src_wavefunc_obj), intent(inout) :: this

  call new(this%param)

  this%itx0 = 1
  this%ity0 = 1
  this%itz0 = 1
  this%itt0 = 1

  if (allocated(this%wavefunc)) deallocate(this%wavefunc)

  return
end subroutine
Subroutine :
this :type(snk_wavefunc_obj), intent(inout)

[Source]

subroutine print_snk_wavefunc(this)
  implicit none
  type(snk_wavefunc_obj), intent(inout) :: this

  if (nodeid==0) then
    write(*,'("==== Sink Parameters ====")')
  endif
  call print(this%param)

  return
end subroutine
Subroutine :
this :type(src_wavefunc_obj), intent(inout)

[Source]

subroutine print_src_wavefunc(this)
  implicit none
  type(src_wavefunc_obj), intent(inout) :: this

  if (nodeid==0) then
    write(*,'("==== Source Parameters ====")')
    write(*,'(" Source Center (xyzt):",4I3)') this%itx0,this%ity0,this%itz0,this%itt0 
  endif
  call print(this%param)

  return
end subroutine
Subroutine :
this :type(src_wavefunc_obj), intent(inout)
iout :integer, intent(in)

[Source]

subroutine read_src_wavefunc(this,iout)
  use error_class
  implicit none
  type(src_wavefunc_obj), intent(inout) :: this
  integer,                intent(in)    :: iout

  !
  ! read source center,time
  !
  if (nodeid==0) then
    read(iout,*) this%itx0,this%ity0,this%itz0,this%itt0
  endif
  if (NPU>1) then
    call comlib_bcast(this%itx0,0)
    call comlib_bcast(this%ity0,0)
    call comlib_bcast(this%itz0,0)
    call comlib_bcast(this%itt0,0)
  endif

  if (   (this%itx0 < 1) .or. (this%itx0 > NTX) .or. (this%ity0 < 1) .or. (this%ity0 > NTY) .or. (this%itz0 < 1) .or. (this%itz0 > NTZ) .or. (this%itt0 < 1) .or. (this%itt0 > NTT)) then
    call error_message(" Source center (x,y,z,t) error.")
    call    error_stop(" x<=[1,NTX], y<=[1,NTY], z<=[1,NTZ] t<=[1,NTT]")
  endif

  !
  ! read source wave function parameters
  !
  call read(this%param,iout)

  return
end subroutine
snk_wavefunc_obj
Derived Type :
param :type(wavefunc_param_obj)
wavefunc(:) :complex(DP), allocatable

quark sink wavefunction

src_wavefunc_obj
Derived Type :
param :type(wavefunc_param_obj)
wavefunc(:,:,:,:) :complex(DP), allocatable
: wf(col,spin,zxy,t) wave function
itx0 :integer
: source center (1..NTX) (1..NTY) (1..NTZ) (1..NTT)
ity0 :integer
: source center (1..NTX) (1..NTY) (1..NTZ) (1..NTT)
itz0 :integer
: source center (1..NTX) (1..NTY) (1..NTZ) (1..NTT)
itt0 :integer
: source center (1..NTX) (1..NTY) (1..NTZ) (1..NTT)

quark source wavefunction

wavefunc_param_obj
Derived Type :
type :integer
: 1 for local, 2 for exponential, 3 for wall, 4 for polynomial, 5 user supplied
Npolsf :integer
: number of polynomial coeffciient
idummy(2) :integer
Asmear :real(DP)
: exponetial parameter, psi(r)=A*exp(-r*B)
Bsmear :real(DP)
: exponetial parameter, psi(r)=A*exp(-r*B)
Psmear(MAXPOL_WAV) :real(DP)
: polynomial coeffcient, psi(r)= 1 - sum_(i=1,Npol) P(i)*r^i

quark source function parameters

Original external subprogram is wavefunc_param_class#wavefunc_param_obj