corespec_eval.f90 25.5 KB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306
!--------------------------------------------------------------------------------
! Copyright (c) 2017 Peter Grünberg Institut, Forschungszentrum Jülich, Germany
! This file is part of FLEUR and available as free software under the conditions
! of the MIT license as expressed in the LICENSE file in more detail.
!--------------------------------------------------------------------------------

MODULE m_corespec_eval

  USE m_corespec
  USE m_types
  USE m_constants

  IMPLICIT NONE

  CONTAINS

!===============================================================================
!
!  S U B R O U T I N E   C O R E S P E C _ G A U N T
!
!-------------------------------------------------------------------------------
!
  SUBROUTINE corespec_gaunt()

!    use factorials
    
    use m_clebsch

    implicit none

!    real :: threejsymbol

    logical :: cmsum,clevn,ctiq1,ctiq2,ctiq3
    real :: twol1p1,twola1p1,twolip1

    smeno = "corespec_gaunt"

    write(*,'(/,a)') trim(smeno)//ssep

!    call init_factorials(6*(lmaxd+1)+1)

    ln = min(0,minval(csv%lc)-1)
    lx = max(csi%lx,maxval(csv%lc)+1)

    lan = 0
    lax = csi%lx+maxval(csv%lc)+1

    lin = minval(csv%lc)-1
    lix = maxval(csv%lc)+1

!!$    print*,"ln,lx,lan,lax,lin,lix"
!!$    print*,ln,lx,lan,lax,lin,lix

    if(.not.allocated(csv%gaunt)) &
         &allocate(csv%gaunt(ln:lx,-lx:lx,lan:lax,-lax:lax,lin:lix,-lix:lix))
    csv%gaunt = 0.d0

! m<=l condition fulfilled by looping m within l value interval {-l,...,+l}

    csv%gaunt = 0.d0
    do l1 = ln,lx
      do m1 = -l1,l1
        do la1 = lan,lax
          do mu1 = -la1,la1
            do li = lin,lix
              do mi = -li,li
                cmsum = (m1+mu1-mi).eq.0  ! sum of m q-nos. = 0
                clevn = mod((l1+la1+li),2).eq.0  ! sum of l q-nos. is even
                ctiq1 = (la1+li-l1).ge.0  ! triangle inequality 1
                ctiq2 = (l1+li-la1).ge.0  ! triangle inequality 2
                ctiq3 = (l1+la1-li).ge.0  ! triangle inequality 3
                twol1p1 = dble(2*l1+1)
                twola1p1 = dble(2*la1+1)
                twolip1 = dble(2*li+1)
                if(cmsum.and.clevn.and.ctiq1.and.ctiq2.and.ctiq3) then
                  csv%gaunt(l1,m1,la1,mu1,li,mi) = &
!                       &threejsymbol((l1),(la1),0,0,(li),0)*&
!                       &threejsymbol((l1),(la1),(m1),(mu1),(li),-(mi)))&
                       &clebsch(real(l1),real(la1),0.0,0.0,real(li),0.0)*&
                       &clebsch(real(l1),real(la1),real(m1),real(mu1),real(li),real(mi))*&
                       &sqrt(twol1p1*twola1p1/(4.d0*pi_const*twolip1))*&
                       &(-1)**(mi)
                  if(csv%gaunt(l1,m1,la1,mu1,li,mi).ne.0.d0) &
                       &write(53,'(6i5,f12.6)') l1,m1,la1,mu1,li,-mi,csv%gaunt(l1,m1,la1,mu1,li,mi)
!!$                  if(abs(csv%gaunt(l1,m1,la1,mu1,li,mi)).lt.1.d-6) &
!!$                       &write(*,'(6i5,f24.20)') l1,m1,la1,mu1,li,-mi,csv%gaunt(l1,m1,la1,mu1,li,mi)
                endif
              enddo
            enddo
          enddo
        enddo
      enddo
    enddo

    if(csi%verb.eq.1) write(*,*) ""

  end subroutine corespec_gaunt
!
!===============================================================================
!===============================================================================
!
!  S U B R O U T I N E   C O R E S P E C _ R M E
!
!-------------------------------------------------------------------------------
!
  subroutine corespec_rme(atoms,input,itype,nstd,&
                          jspins,jspin,efermi,&
                          msh,vr,f,g)

    USE m_constants, ONLY : c_light
    USE m_setcor
    USE m_differ
    USE m_intgr, ONLY : intgr3
    USE m_dr2fdr
    USE m_sphbes
    USE m_intgr, ONLY : intgr3

    implicit none

    TYPE(t_atoms),INTENT(IN)   :: atoms
    TYPE(t_input),INTENT(IN)   :: input

    integer, intent(in) :: itype  ! call in ntype loop with itype = n
    integer, intent(in) :: nstd
    integer, intent(in) :: jspins,jspin
    real, intent(in) :: efermi
    integer, intent(in) :: msh
    real, intent    (in) :: vr(atoms%jmtd,atoms%ntype,jspins)
    real, intent (in) :: f(atoms%jmtd,2,0:atoms%lmaxd,jspin:jspin)
    real, intent (in) :: g(atoms%jmtd,2,0:atoms%lmaxd,jspin:jspin)

    integer :: nr,lx,lax,lin,lix,nqv,nen,nex

    integer :: ir,id,iljc,ic,il,ila,iqv,ie,ierr
    integer :: nst,kappa(nstd),nprnc(nstd)
    real :: nc,nlc,njc
    real :: c,bmu,t2,weight,e,d,rn,res,qr
    real :: vrd(msh),occ(nstd,jspins),a(msh),b(msh)
    real :: resd

    real, allocatable :: fpd(:)
    real, allocatable :: fp(:),fc(:),fsb(:)
    real :: sum1,sum2,sum3,sum1d,sum2d

    smeno = "corespec_rme"

    if(itype.ne.csi%type) return

    write(*,'(/,a)') trim(smeno)//ssep

    c = c_light(1.0)

    nr = atoms%jri(itype)

    allocate(fp(nr),fpd(nr),fc(nr))

  ! CORE functions
  ! csv%fc(ir,:,:,:) : ir = 1:nr
  ! csv%fc(:,id,:,:) : id = 1 { r*fc(r) } or 2 { r*[dfc(r)/dr] }
  ! csv%fc(:,:,iljc,:) : iljc = 1:csv%nljc
  ! csv%fc(:,:,:,ic) : ic = 1 { large component } or 2 { small component }
    if(.not.allocated(csv%fc)) allocate(csv%fc(nr,2,csv%nljc,2))
    csv%fc = 0.d0

  ! core setup
    bmu = 0.0

    CALL setcor(itype,jspins,atoms,input,bmu,nst,kappa,nprnc,occ)

  ! extend core potential
    vrd(1:nr) = vr(1:nr,itype,jspin)
    t2 = vrd(nr)/(nr-msh)
    do ir = nr+1,msh
      vrd(ir) = vrd(nr)+t2*(ir-nr)
    enddo

  ! calculate core radial functions
    nc = real(csv%nc)
    do iljc = 1,csv%nljc
      njc = real(edgej(csi%edgeidx(iljc)))/2.0
      nlc = real(edgel(csi%edgeidx(iljc)))
      weight = 2*njc+1.0
      csv%eedge(iljc) = -2*(atoms%zatom(itype)/(nc+nlc))**2
      d = exp(atoms%dx(itype))
      rn = atoms%rmsh(1,itype)*(d**(msh-1))

      CALL differ(nc,nlc,njc,c,atoms%zatom(itype),atoms%dx(itype),&
                  atoms%rmsh(1,itype),rn,d,msh,vrd,&
                  e,&
                  a,b,ierr)

      csv%eedge(iljc)=dble(e)
      csv%fc(:,1,iljc,1) = a(1:nr)  ! large component
      csv%fc(:,1,iljc,2) = b(1:nr)  ! small component
      do ic = 1,2
        fp(:) = real(csv%fc(:,1,iljc,ic)*atoms%rmsh(1:nr,itype))
        CALL dr2fdr(fp,atoms%rmsh(1,itype),nr,fc)
        csv%fc(:,2,iljc,ic)=dble(fc(:)/atoms%rmsh(1:nr,itype))

        if(ic.eq.1) then
        do ir=1,nr
        write(90,'(2i5,16e12.4)') iljc,ir,atoms%rmsh(ir,itype),csv%fc(ir,1,iljc,ic),csv%fc(ir,2,iljc,ic)
        enddo
        write(90,*) ''
        write(90,*) ''
        endif

      enddo
      
      fp = csv%fc(:,1,iljc,1)**2
      CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,sum1)
      fp = csv%fc(:,1,iljc,2)**2
      CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,sum2)
      write(*,'(a,i5,3f8.4)') "ui",0,sum1,sum2,sum1+sum2

      fp = csv%fc(:,2,iljc,1)**2
      CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,sum1)
      fp = csv%fc(:,2,iljc,2)**2
      CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,sum2)
      write(*,'(a,i5,3f8.4)') "ui",0,sum1,sum2,sum1+sum2

      write(60,*) ""
      csv%occ(iljc) = dble(occ((nc-1)**2+csi%edgeidx(iljc),jspin))
      write(*,"(a,2(a,i2),a,f3.1,2(a,i2),a,f16.8,a)") trim(smeno)//ssep,&
           &"core state: iljc = ",iljc,&
           &", nc = ",nint(nc),&
           &", njc = ",njc,&
           &", nlc = ",nint(nlc),&
           &", occ. csv%occ = ",nint(csv%occ(iljc)),&
           &", energy csv%eedge(iljc) = ",csv%eedge(iljc)," Ha found"
      if(efermi-csv%eedge(iljc).lt.ecoredeep) then
        write(*,csmsgsfs)  trim(smeno),&
             &"core state energy found not very deep: ",&
             &"efermi-csv%eedge(iljc) = ",&
             &(efermi-csv%eedge(iljc))*hartree_to_ev_const,"eV ; are you sure ? "//csmsgwar
      endif
    enddo

    CALL corespec_eloss_qv(efermi)  ! set-up csv%eloss and csv%qv arrays

    lx = csi%lx  ! lmax for l index
    lax = lx+maxval(csv%lc)+1  ! lmax for la index
    lin = minval(csv%lc)  ! minimum lc q-no.
    lix = maxval(csv%lc)  ! maximum lc q-no.
    nqv = csv%nqv
    nen = csv%nen
    nex = csv%nex

    allocate(fsb(0:lax))

  ! VALENCE functions
  ! csv%fv(ir,:,:,:) : ir = 1:nr
  ! csv%fv(:,il,:,:) : il = 0:csi%lx
  ! csv%fv(:,:,id,:) : id = 1 { a.u } or 2 { b.u' }
  ! csv%fv(:,:,:,ic) : ic = 1 { large component } or 2 { small component }
    if(.not.allocated(csv%fv)) allocate(csv%fv(nr,0:lx,2,2))
    csv%fv = 0.d0

    do ic = 1,2
      do il = 0,lx
        csv%fv(:,il,1,ic) = f(1:nr,ic,il,jspin)
        csv%fv(:,il,2,ic) = g(1:nr,ic,il,jspin)

        if(ic.eq.1) then
        do ir=1,nr
        write(70,'(3i5,16e12.4)') ic,il,ir,atoms%rmsh(ir,itype),csv%fv(ir,il,1,ic),csv%fv(ir,il,2,ic)
        enddo
        write(70,*) ''
        write(70,*) ''
        endif

       fp(:) = csv%fv(:,il,1,ic)**2
        CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,sum1)
        fp(:) = csv%fv(:,il,2,ic)**2
        CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,sum2)
        fp(:) = csv%fv(:,il,1,ic)*csv%fv(:,il,2,ic)
        CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,sum3)
        write(*,'(a,i5,3f8.4)') "u ",il,sum1,sum2,sum3

      enddo
    enddo

  ! BESSEL functions
  ! csv%fb(ir,:,:,:,:) : ir = 1:nr
  ! csv%fb(:,il,:,:,:) : il = 0:lax
  ! csv%fb(:,,:,iljc,:,:) : iljc = 1:csv%nljc
  ! csv%fb(:,:,:,iqv,:) : iqv = 1:nqv
  ! csv%fb(:,:,:,:,ie) : ie = nen:nex
    if(.not.allocated(csv%fb)) allocate(csv%fb(nr,0:lax,csv%nljc,nqv,nen:nex))
    csv%fb = 0.d0

    do ie = nen,nex
      do iqv = 1,nqv
        do iljc = 1,csv%nljc
          do ir = 1,nr
            fsb=0.d0
            qr = real(csv%qv(0,iljc,iqv,ie)*atoms%rmsh(ir,itype))
            CALL sphbes(lax,qr,fsb)
            csv%fb(ir,:,iljc,iqv,ie) = dble(fsb)
!            write(70,'(4i5,16e12.4)') ie,iqv,iljc,ir,atoms%rmsh(ir,itype),fsb
          enddo
!          write(70,*) ''
        enddo
      enddo
    enddo

307 308 309 310 311 312 313 314
    if(.NOT.ALLOCATED(csv%rmeA)) THEN
       ALLOCATE(csv%rmeA(2,0:lx,0:lax,csv%nljc,2,nqv,nen:nex))
       ALLOCATE(csv%rmeB(2,0:lx,0:lax,csv%nljc,2,nqv,nen:nex))
       ALLOCATE(csv%rmeC(2,0:lx,0:lax,csv%nljc,2,nqv,nen:nex))
    END IF
    csv%rmeA = 0.0
    csv%rmeB = 0.0
    csv%rmeC = 0.0
315 316 317 318 319 320 321 322 323 324 325 326

    do ie = nen,nex
      do iqv = 1,nqv
        do ic = 1,2
          do iljc = 1,csv%nljc
            do ila = 0,lax
              do il = 0,lx
                do id = 1,2
                  fp(:)=csv%fc(1:nr,1,iljc,ic)*&
                       &csv%fv(1:nr,il,id,ic)*&
                       &csv%fb(1:nr,ila,iljc,iqv,ie)
                  CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,res)
327
                  csv%rmeA(id,il,ila,iljc,ic,iqv,ie)=dble(res)
328 329
                  fp(:)=fp(:)/atoms%rmsh(1:nr,itype)
                  CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,res)
330
                  csv%rmeC(id,il,ila,iljc,ic,iqv,ie)=dble(res)
331 332 333 334
                  fp(:)=csv%fc(1:nr,2,iljc,ic)*&
                       &csv%fv(1:nr,il,id,ic)*&
                       &csv%fb(1:nr,ila,iljc,iqv,ie)!/atoms%rmsh(1:nr,itype)
                  CALL intgr3(fp,atoms%rmsh(1,itype),atoms%dx(itype),nr,res)
335
                  csv%rmeB(id,il,ila,iljc,ic,iqv,ie)=dble(res)
336
                  write(41,'(7(a,i5),3f12.6)') 'ie=',ie,' iqv=',iqv,' ic=',ic,&
337 338 339 340
                        ' iljc=',iljc,' id=',id,' ila=',ila,' il=',il,&
                        csv%rmeA(id,il,ila,iljc,ic,iqv,ie),&
                        csv%rmeB(id,il,ila,iljc,ic,iqv,ie),&
                        csv%rmeC(id,il,ila,iljc,ic,iqv,ie)
341 342 343 344 345 346 347 348
                enddo  ! id
              enddo  ! il
            enddo  ! ila
          enddo  ! iljc
        enddo  ! ic
      enddo  ! iqv
    enddo  ! ie

349
    print*,size(3*csv%rmeA)
350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578
 
    deallocate(fsb,fc,fpd,fp)

    if(csi%verb.eq.1) write(*,*) ""

  end subroutine corespec_rme
!
!===============================================================================
!===============================================================================
!
!  S U B R O U T I N E   C O R E S P E C _ D O S
!
!-------------------------------------------------------------------------------
!
  subroutine corespec_dos(atoms,usdus,ispin,lmd,nkpt,ikpt,&
                          neigd,noccbd,efermi,sig_dos,&
                          eig,we,acof,bcof,&
                          ccof)

    IMPLICIT NONE

    TYPE (t_atoms), INTENT(IN)   :: atoms
    TYPE (t_usdus), INTENT(IN)   :: usdus

!     .. Scalar Arguments ..
    integer, intent(in) :: ispin,lmd,nkpt,ikpt
    integer, intent(in) :: neigd,noccbd
    real, intent(in) :: efermi,sig_dos
!     .. Array Arguments ..
    real, intent (in) :: eig(neigd),we(noccbd)
    complex, intent (in) :: acof(noccbd,0:lmd,atoms%nat)
    complex, intent (in) :: bcof(noccbd,0:lmd,atoms%nat)
    complex, intent (in) :: ccof(-atoms%llod:atoms%llod,noccbd,atoms%nlod,atoms%nat)

! local variables
    integer :: lx,lmx,nen,nex
    integer :: iatom,iband,l1,m1,l2,m2,lm1,lm2,ie!,ljc,iqv
    real :: sigma,eigos(noccbd)
    real :: sum11,sum22

    smeno = "corespec_dos"

    lx = csi%lx
    lmx = lx*(lx+2)
    nen = csv%nen
    nex = csv%nex
    iatom = atoms%neq(csi%type)
    sigma = sqrt(2.d0)*sig_dos*hartree_to_ev_const
    sigma = sig_dos*hartree_to_ev_const
    eigos(1:noccbd) = (eig(1:noccbd)-efermi)*hartree_to_ev_const/dble(sigma)

    if(ikpt.eq.1) then
      write(*,'(/,a)') trim(smeno)//ssep
      if(.not.allocated(csv%dose)) allocate(csv%dose(2,2,0:lmx,0:lmx,0:nex))
      if(.not.allocated(csv%dosb)) allocate(csv%dosb(2,2,0:lmx,0:lmx,noccbd))
      if(.not.allocated(csv%eos)) then
        allocate(csv%eos(0:nex))
        csv%eos(:) = csv%egrid(:)/dble(sigma)
      endif
      csv%dose = 0.d0
    endif
    csv%dosb = 0.d0

    do iband = 1,noccbd
      do l1 = 0,lx
        do m1 = -l1,l1
          lm1 = l1*(l1+1)+m1
!!$          do l2 = 0,lx
!!$            do m2 = -l2,l2
!!$              lm2 = l2*(l2+1)+m2
          csv%dosb(1,1,lm1,lm1,iband) = dble(acof(iband,lm1,iatom)*&
               &conjg(acof(iband,lm1,iatom)))*we(1)
          csv%dosb(1,2,lm1,lm1,iband) = dble(acof(iband,lm1,iatom)*&
               &conjg(bcof(iband,lm1,iatom)))
          csv%dosb(2,1,lm1,lm1,iband) = dble(bcof(iband,lm1,iatom)*&
               &conjg(acof(iband,lm1,iatom)))
          csv%dosb(2,2,lm1,lm1,iband) = dble(bcof(iband,lm1,iatom)*&
               &conjg(bcof(iband,lm1,iatom)))*we(1)*usdus%ddn(l1,csi%type,ispin)
!!!!! this has to be checked: is >> ddn << factor necessary !!!!!
!!$        enddo
!!$        enddo
        enddo
      enddo
      if(eigos(iband)+3.d0*sigma.ge.csv%eos(0).and.&
           &eigos(iband)-3.d0*sigma.le.csv%eos(nex)) then
        do ie = 0,nex
          csv%dose(:,:,:,:,ie) = csv%dose(:,:,:,:,ie)+&
               &csv%dosb(:,:,:,:,iband)*exp(-(eigos(iband)-csv%eos(ie))**2)
        enddo
      endif
    enddo

    if(ikpt.eq.nkpt) then
      csv%dose = csv%dose/(sqrt(pi_const)*sigma)
      do ie=0,nex
        write(36,*) csv%egrid(ie),sum(csv%dose(1,1,:,:,ie)+csv%dose(2,2,:,:,ie))
      enddo
      write(36,*) ""
      write(*,'(10i8)') atoms%llod,noccbd,atoms%nlod,atoms%nat,neigd,atoms%ntype,atoms%lmaxd
      write(*,'(10i8)') lmd,atoms%ntype

      if(csi%verb.eq.1) write(*,*) ""
    endif

  end subroutine corespec_dos
!
!===============================================================================
!===============================================================================
!
!  S U B R O U T I N E   C O R E S P E C _ D D S C S
!
!-------------------------------------------------------------------------------
!
  subroutine corespec_ddscs(jspin,jspins)

    use m_ylm

    implicit none

    integer, intent(in) :: jspin,jspins

    integer :: lx,lmx,lan,lax,nqv,nen,nex,nor

    integer :: ic,ie,iqv,ior,it,iljc,imi,id1,id2,ip1,ip2
    integer :: l1,l2,m1,m2,lm1,lm2
    integer :: la1,la2,mu1,mu2
    integer :: li,mi
    integer :: lamu,lamu1,lamu2

    real :: gamma,beta,rho,qepref,orvec(3)
    real :: ga(0:2,2)
    real :: prd(0:2,0:2)
    complex :: td(2),orfac,ila1la2
    complex, allocatable :: tdy(:,:),orpref(:),ylm(:,:)

    smeno = "corespec_ddscs"

    write(*,'(/,a)') trim(smeno)//ssep

    lx = csi%lx
    lmx = lx*(lx+2)
    lan = 0
    lax = csi%lx+maxval(csv%lc)+1
    nqv = csv%nqv
    nen = csv%nen
    nex = csv%nex

    nor = 1
    orvec = (/1.d0,0.d0,0.d0/)

    if(.not.allocated(csv%ddscs)) then
      allocate(csv%ddscs(2,0:nor,1:csv%nljc,nqv,0:nex))
      csv%ddscs = cmplx(0.d0,0.d0)
    endif
    if(.not.allocated(tdy)) allocate(tdy(0:nor,2))
    if(.not.allocated(orpref)) then
      allocate(orpref(0:nor))
      orpref(0) = 1.d0
      if(nor.gt.0) then
        orpref(1:nor) = (4.d0*pi_const)**2
        if(.not.allocated(ylm)) allocate(ylm(0:lax*(lax+2),nor))
        do ior = 1,nor
          CALL ylm4(lax,orvec,ylm(0,ior))
          do la1 = lan,lax ; do mu1 = -la1,la1
            lamu = la1*(la1+1)+mu1
            write(98,'(3i5,2f12.8)') la1,mu1,lamu,ylm(lamu,ior)
          enddo; enddo
        enddo
      endif
    endif

    ic = 1
    gamma = csv%gamma
    beta = csv%beta

    rho = alpha*beta*sqrt(4.d0*pi_const/3.d0)
    print*,gamma,beta,rho
!    rho = 0.d0

    do ie = nen,nex  ! energy
      do iqv = 1,nqv  ! q-vector
        do iljc = 1,csv%nljc  ! core levels
          li = edgel(csi%edgeidx(iljc))
          qepref = 4.d0*gamma**2*csv%qv1(iljc,iqv,ie)/csv%qv0/&
               &(csv%qv(0,iljc,iqv,ie)**2-(csv%eloss(iljc,ie)*alpha)**2)**2
!!$          write(*,'(2i5,3f20.4)') ie,iljc,csv%qv(0,iljc,iqv,ie),csv%eloss(iljc,ie)*alpha,qepref
          tdy = cmplx(0.d0,0.d0)

          do imi = 1,nint(csv%occ(iljc)*jspins/2)
            mi = sign(jspin)*(edgej(csi%edgeidx(iljc))-4*(imi-1)-1)/2
            print*,jspin,ie,iljc,li,mi
            write(39,*) jspin,ie,iljc,li,mi

            do l1 = 0,lx ; do m1 = -l1,l1
              lm1 = l1*(l1+1)+m1
              do la1 = lan,lax ; do mu1 = -la1,la1
                lamu1 = la1*(la1+1)+mu1

                ga(0,1) = csv%gaunt(l1,-m1,la1,mu1,li,mi)
                ga(1,1) = csv%gaunt(li+1,-mi,li,mi,1,0)*&
                     &csv%gaunt(l1,-m1,la1,mu1,li+1,mi)+&
                     &csv%gaunt(li-1,-mi,li,mi,1,0)*&
                     &csv%gaunt(l1,-m1,la1,mu1,li-1,mi)
                ga(2,1) = csv%gaunt(li+1,-mi,li,mi+1,1,-1)*&
                     &csv%gaunt(l1,-m1,la1,mu1,li+1,mi)+&
                     &csv%gaunt(li-1,-mi,li,mi+1,1,-1)*&
                     &csv%gaunt(l1,-m1,la1,mu1,li-1,mi)*&
                     &sqrt(dble(2*(li-mi)*(li+mi+1)))+mi*ga(1,1)

                do l2 = 0,lx ; do m2 = -l2,l2
                  lm2 = l2*(l2+1)+m2
                  do la2 = lan,lax ; do mu2 = -la2,la2
                    lamu2 = la2*(la2+1)+mu2
                    
                    if(l1.eq.l2.and.m1.eq.m2) then

                    ga(0,2) = csv%gaunt(l2,-m2,la2,mu2,li,mi)
                    ga(1,2) = csv%gaunt(li+1,-mi,li,mi,1,0)*&
                             &csv%gaunt(l2,-m2,la2,mu2,li+1,mi)+&
                             &csv%gaunt(li-1,-mi,li,mi,1,0)*&
                             &csv%gaunt(l2,-m2,la2,mu2,li-1,mi)
                    ga(2,2) = csv%gaunt(li+1,-mi,li,mi+1,1,-1)*&
                             &csv%gaunt(l2,-m2,la2,mu2,li+1,mi)+&
                             &csv%gaunt(li-1,-mi,li,mi+1,1,-1)*&
                             &csv%gaunt(l2,-m2,la2,mu2,li-1,mi)*&
                             &sqrt(dble(2*(li-mi)*(li+mi+1)))+mi*ga(1,2)

                    prd = 0.d0

579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600
                    do id1 = 1,2 ; 
                      do id2 = 1,2
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeA(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeA(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeA(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeB(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeA(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeC(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeB(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeA(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeB(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeB(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeB(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeC(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeC(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeA(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeC(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeB(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                        prd(ip1,ip2) = prd(ip1,ip2)+ &
                           csv%rmeC(id1,l1,la1,iljc,ic,iqv,ie)*csv%rmeC(id1,l1,la1,iljc,ic,iqv,ie)*csv%dose(id1,id2,lm1,lm2,ie)
                      enddo 
                    enddo
601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754

                    td(1) = prd(0,0)*ga(0,1)*ga(0,2)
                    td(2) = cone*rho**2*(&
                           &prd(1,1)*ga(1,1)*ga(1,2)&
                          &+prd(2,2)*ga(2,1)*ga(2,2)&
                          &-prd(1,2)*ga(1,1)*ga(2,2)&
                          &-prd(2,1)*ga(2,1)*ga(1,2))&
                          &+cimu*rho*(&
                          &-prd(0,1)*ga(0,1)*ga(1,2)&
                          &+prd(0,2)*ga(0,1)*ga(2,2)&
                          &+prd(1,0)*ga(1,1)*ga(0,2)&
                          &-prd(2,0)*ga(2,1)*ga(0,2))
                    
                    ila1la2 = cimu**(la1-la2)
                    
                    if(abs(real(td(1))).gt.0.d0.or.abs(real(td(2))).gt.0.d0.or.abs(imag(td(2))).gt.0.d0) then
                    write(39,'(2f4.0,i2,6i4,a,6i4,a,6f7.3,a,4f10.6)') ila1la2,la1-la2,l1,-m1,la1,mu1,li,mi,'  ',l2,-m2,la2,mu2,li,mi,'  ',ga(0,1),ga(0,2),ga(1,1),ga(1,2),ga(2,1),ga(2,2),'  ',1000000*td
                    endif

                    do ior = 0,nor  ! orientation
                      if(ior.eq.0) then
                        orfac = cone
                      else
                        orfac = ylm(lamu1,ior)*conjg(ylm(lamu2,ior))
                      endif
                      
                      tdy(ior,1:2) = tdy(ior,1:2)+td(1:2)*orfac*ila1la2

                    enddo  ! ior

                  endif

                  enddo; enddo
                enddo; enddo

              enddo; enddo
            enddo; enddo

          enddo  ! mi

          do it = 1,2
            do ior = 0,nor
              csv%ddscs(it,ior,iljc,iqv,ie) = csv%ddscs(it,ior,iljc,iqv,ie)+&
                   &qepref*orpref(ior)*tdy(ior,it)
            enddo
          enddo

        enddo  ! iljc
      enddo  ! iqv
    enddo  ! ie

    if(jspin.eq.1) then
      do ior = 0,nor
      do iljc = 1,csv%nljc
        do ie = nen,nex
          write(37,'(2i5,f8.3,4es16.4)') ior,iljc,csv%eloss(iljc,ie)*hartree_to_ev_const,csv%ddscs(1,ior,iljc,1,ie),csv%ddscs(2,ior,iljc,1,ie)
        enddo
        write(37,*) ""
      enddo
      write(37,*) ""
      enddo
    endif
    if(jspin.eq.2) then
      do ior = 0,nor
      do iljc = 1,csv%nljc
        do ie = nen,nex
          write(38,'(2i5,f8.3,4es16.4)') ior,iljc,csv%eloss(iljc,ie)*hartree_to_ev_const,csv%ddscs(1,ior,iljc,1,ie),csv%ddscs(2,ior,iljc,1,ie)
        enddo
        write(38,*) ""
      enddo
      write(38,*) ""
      enddo
    endif

    if(csi%verb.eq.1) write(*,*) ""

  end subroutine corespec_ddscs
!
!===============================================================================
!===============================================================================
!
!  S U B R O U T I N E   C O R E S P E C _ E L O S S _ Q V
!
!-------------------------------------------------------------------------------
!
  subroutine corespec_eloss_qv(efermi)

    implicit none

    real, intent(in) :: efermi

    integer :: ie,iljc,iqv
    real :: eout,relfac

    smeno = "corespec_eloss_qv"

    write(*,'(/,a)') trim(smeno)//ssep

    csv%nqv = 1

    if(.not.allocated(csv%eloss)) &
         &allocate(csv%eloss(csv%nljc,csv%nen:csv%nex))
    if(.not.allocated(csv%qv1)) &
         &allocate(csv%qv1(csv%nljc,csv%nqv,csv%nen:csv%nex))
    do ie = csv%nen,csv%nex
      do iljc = 1,csv%nljc
        csv%eloss(iljc,ie) = csv%egrid(ie)/hartree_to_ev_const+dble(efermi)-csv%eedge(iljc)
!!$        print*,iljc,ie,csv%egrid(ie),csv%eloss(iljc,ie)
      enddo
    enddo

    csv%qv0 = e2q(csi%ek0/hartree_to_ev_const)
    relfac = (mec2)**2/(csi%ek0+mec2)**2
!!$    print*,csi%ek0,csv%qv0

    if(.not.allocated(csv%qv)) &
         &allocate(csv%qv(0:3,csv%nljc,csv%nqv,csv%nen:csv%nex))
    csv%qv=0.d0
    do ie = csv%nen,csv%nex
      do iqv = 1,csv%nqv
        do iljc = 1,csv%nljc
          eout = csi%ek0/hartree_to_ev_const-csv%eloss(iljc,ie)
          csv%qv1(iljc,iqv,ie) = e2q(eout)
          csv%qv(3,iljc,iqv,ie) = (csv%qv1(iljc,iqv,ie)-csv%qv0)*relfac
          csv%qv(0,iljc,iqv,ie) = sqrt(&
               &dot_product(csv%qv(1:3,iljc,iqv,ie),csv%qv(1:3,iljc,iqv,ie)))
!!$          write(*,'(3i5,2f16.2,6f16.6)') ie,iqv,iljc,csi%ek0,eout*hartree_to_ev_const,csv%qv1(iljc,iqv,ie),csv%eloss(iljc,ie),csv%qv(:,iljc,iqv,ie)
        enddo
      enddo
    enddo

    if(csi%verb.eq.1) write(*,*) ""

  end subroutine corespec_eloss_qv
!
!===============================================================================
!===============================================================================
!  F U N C T I O N   E 2 Q
!-------------------------------------------------------------------------------
!
  real function e2q(e)

    use m_corespec, only : mec2,alpha
    implicit none
    real, intent(in) :: e

    e2q=sqrt(e**2+2.d0*e*mec2/hartree_to_ev_const)*alpha

  end function e2q
!
!===============================================================================


end module m_corespec_eval