VYukawaFilm.f90 42.5 KB
Newer Older
1 2
module m_VYukawaFilm

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
  ! Computation of the film-case Yukawa potential for the preconditioning of the
  ! residual charge density in 5 steps:
  ! 1. pseudo-charge density generation
  ! 2. vacuum potential generation
  ! 3. interstitial potential generation
  ! 4. muffin-tin potential generation
  ! 5. modification for charge neutrality

  ! The Yukawa potential is the solution to the modified Helmholtz equation
  ! ( Delta - lambda^2 ) V_lambda = -4 pi ( rho_out - rho_in )
  ! subject to some conditions. 
  ! The general scheme (steps 1 to 4) is the same as for the Poisson equation --
  ! we use Green function methods for the z-dependent vacuum and interstitial
  ! potentials as well as for the muffin-tin potential and apply Weinert's
  ! method.
  ! You can choose between two variants:
  ! 1. variant: 
  ! zero Dirichlet boundary conditions at +/- infinity; 
  ! multiplication with a decaying exponential in vacuum;
  ! modification in the film for charge neutrality (step 5)
  ! 2. variant:
  ! zero Dirichlet boundary conditions near the film boundary in vacuum (D/2+2R);
  ! modification in the film for charge neutrality (step 5) 
  ! In both cases charge neutrality is broken.
  ! To restore charge neutrality, we need the integral over the potential to be
  ! zero.
  ! In step 5 we therefore solve the modified Helmholtz equation again with
  ! constant right-hand side, for an additive correction to the potential.
  ! The constant is chosen such that the integral over the final potential is
  ! zero.

  contains 



  subroutine VYukawaFilm( stars, vacuum, cell, sym, input, mpi, atoms, sphhar, oneD, noco, den, &
39 40 41 42 43
                          VYukawa )

    use m_constants
    use m_types
    use m_psqpw
44
    use m_vmts
45 46 47 48 49 50 51 52 53 54 55
    implicit none

    type(t_stars),      intent(in)    :: stars
    type(t_vacuum),     intent(in)    :: vacuum
    type(t_cell),       intent(in)    :: cell
    type(t_sym),        intent(in)    :: sym
    type(t_input),      intent(in)    :: input
    type(t_mpi),        intent(in)    :: mpi
    type(t_atoms),      intent(in)    :: atoms 
    type(t_sphhar),     intent(in)    :: sphhar
    type(t_oneD),       intent(in)    :: oneD
56 57
    type(t_noco),       intent(in)    :: noco
    type(t_potden),     intent(inout) :: den
58 59 60 61

    type(t_potden),     intent(inout) :: VYukawa

    complex                           :: psq(stars%ng3)
62 63 64
    complex                           :: alphm(stars%ng2,2)
    real                              :: dh_prec
    real                              :: coshdh(stars%ng2)
65

66
 
67 68
    ! PSEUDO-CHARGE DENSITY

69
    call psqpw( mpi, atoms, sphhar, stars, vacuum, cell, input, sym, oneD, &
70 71
                den%pw(:,1), den%mt(:,:,:,1), den%vacz(:,:,1), .false., VYukawa%potdenType, &
                psq )
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
    ChooseVariant: if ( .true. ) then

      ! VACUUM POTENTIAL

      call VYukawaFilmVacuumVariant1( &
              stars, vacuum, cell, sym, input, atoms%rmt(1), &
              psq, den%vacxy(:,:,:,1), den%vacz(:,:,1), &
              VYukawa%vacxy, VYukawa%vacz, alphm )


      ! INTERSTITIAL POTENTIAL

      call VYukawaFilmInterstitialVariant1( &
              stars, vacuum, cell, sym, input, &
              psq, VYukawa%vacxy, VYukawa%vacz, alphm, &
              VYukawa%pw(:,1) )

    else ChooseVariant

      ! VACUUM POTENTIAL

      call VYukawaFilmVacuumVariant2( &
              stars, vacuum, cell, sym, input, 2*atoms%rmt(1), &
              psq, den%vacxy(:,:,:,1), den%vacz(:,:,1), &
              VYukawa%vacxy, VYukawa%vacz, alphm, dh_prec, coshdh )


      ! INTERSTITIAL POTENTIAL

      call VYukawaFilmInterstitialVariant2( &
              stars, vacuum, cell, sym, input, &
              psq, VYukawa%vacxy, VYukawa%vacz, alphm, dh_prec, coshdh, &
              VYukawa%pw(:,1) )
107

108
    end if ChooseVariant
109 110 111 112


    ! MUFFIN-TIN POTENTIAL

113 114 115
    call Vmts( input, mpi, stars, sphhar, atoms, sym, cell, oneD, &
               VYukawa%pw(:,1), den%mt(:,0:,:,1), VYukawa%potdenType, &
               VYukawa%mt(:,0:,:,1) )
116

117 118
 
    ! MODIFICATION FOR CHARGE NEUTRALITY
119

120 121
    call VYukawaModify( stars, vacuum, cell, sym, input, mpi, atoms, sphhar, oneD, noco, &
                        den, &
122
                        VYukawa )
123 124


125
  end subroutine VYukawaFilm
126 127


128

129 130 131 132
  subroutine VYukawaFilmVacuumVariant1( &
                stars, vacuum, cell, sym, input, rmt, &
                psq, rhtxy, rht, &
                VVxy, VVz, alphm )
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
    ! 1. part: Compute the contribution from the interstitial charge density to the vacuum potential as a function of q_xy and z (analytic expression for integral)
    ! 2. part: Compute the contribution from the vacuum charge density to the vacuum potential as a function of q_xy and z by numerical integration

    use m_ExpSave
    use m_constants
    use m_types
    use m_intgr, only: intgz1Reverse
    use m_qsf
    implicit none

    type(t_stars),  intent(in)  :: stars
    type(t_vacuum), intent(in)  :: vacuum
    type(t_cell),   intent(in)  :: cell
    type(t_sym),    intent(in)  :: sym
    type(t_input),  intent(in)  :: input
    real,           intent(in)  :: rmt
    complex,        intent(in)  :: psq(stars%ng3)
    complex,        intent(in)  :: rhtxy(vacuum%nmzxyd,stars%ng2-1,2)
    real,           intent(in)  :: rht(vacuum%nmzd,2) 

    complex,        intent(out) :: VVxy(vacuum%nmzxyd,2:stars%ng2,2) ! this is the qxy /= 0 part of the vacuum potential
    real,           intent(out) :: VVz(vacuum%nmzd,2)                ! this is the qxy  = 0 part of the vacuum potential
    complex,        intent(out) :: alphm(stars%ng2,2)                ! these are the integrals in upper and lower vacuum, now including the first star---integral for star ig2 is in alphm(ig2,ivac) 

    complex                     :: sum_qz(2,stars%ng2)
    complex                     :: c_ph(-stars%mx3:stars%mx3,stars%ng2)
    complex                     :: signIqz
    real                        :: g_damped(stars%ng2), qz, sign, vcons(stars%ng2) 
    real                        :: exp_m(vacuum%nmzd,stars%ng2), exp_p(vacuum%nmzd,stars%ng2)
    real                        :: expDhg(stars%ng2), expDg(stars%ng2)
    real                        :: z(vacuum%nmzd)
    integer                     :: iz, irec2, irec3, ivac, iqz
    complex                     :: fa(vacuum%nmzxyd,2:stars%ng2), fb(vacuum%nmzxyd,2:stars%ng2)
    complex                     :: alpha(vacuum%nmzxyd,2:stars%ng2,2), beta(vacuum%nmzxyd,2:stars%ng2,2), gamma(vacuum%nmzxyd,2:stars%ng2)
    real                        :: ga(vacuum%nmzd), gb(vacuum%nmzd)
    real                        :: delta(vacuum%nmzd,2), epsilon(vacuum%nmzd,2), zeta(vacuum%nmzd)


    ! DEFINITIONS / ALLOCATIONS / INITIALISATIONS    
    
    do iz = 1, vacuum%nmz
      z(iz) = ( iz - 1 ) * vacuum%delz
    end do
    do irec2 = 1, stars%ng2
      g_damped(irec2) = sqrt( stars%sk2(irec2) ** 2 + input%preconditioning_param ** 2 )
      vcons(irec2) = tpi_const / g_damped(irec2)
      do iz = 1, vacuum%nmz
        exp_m(iz,irec2) = exp_save( - g_damped(irec2) * z(iz) )
        exp_p(iz,irec2) = exp_save(   g_damped(irec2) * z(iz) )
      end do
      expDhg(irec2) = exp_save( - cell%z1 * g_damped(irec2) )
      expDg(irec2)  = exp_save( -2 * cell%z1 * g_damped(irec2) )
    end do
    sum_qz = (0.,0.)
    VVxy = (0.,0.)
    VVz = 0.
   

    ! CONTRIBUTION FROM THE INTERSTITIAL CHARGE DENSITY

    do irec2 = 1, stars%ng2
      do ivac = 1, vacuum%nvac
        sign = 3. - 2. * ivac
        do iqz = -stars%mx3, stars%mx3
          irec3 = stars%ig(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
          ! use only stars within the g_max sphere -> stars outside the sphere have per definition index ig3n = 0
          if ( irec3 /= 0 ) then
            c_ph(iqz,irec2) = stars%rgphs(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
            qz = iqz * cell%bmat(3,3)
            signIqz = sign * ImagUnit * qz 
            sum_qz(ivac,irec2) = sum_qz(ivac,irec2) + c_ph(iqz,irec2) * psq(irec3) * ( exp( signIqz * cell%z1 ) - exp( - signIqz * cell%z1 ) * expDg(irec2) ) / ( signIqz + g_damped(irec2) )
          endif
        enddo
        if( irec2 /= 1 ) then
          VVxy(1:vacuum%nmzxy,irec2,ivac) = vcons(irec2) * sum_qz(ivac,irec2) * exp_m(1:vacuum%nmzxy,irec2)
        else
          VVz(1:vacuum%nmz,ivac)          = vcons(1)     * sum_qz(ivac,1)     * exp_m(1:vacuum%nmz,1)
        end if
      enddo
    enddo


    ! CONTRIBUTION FROM THE VACUUM CHARGE DENSITY

    ! shifting z:
    do irec2 = 1, stars%ng2
      exp_m(1:vacuum%nmz,irec2) = exp_m(1:vacuum%nmz,irec2) * expDhg(irec2)
      exp_p(1:vacuum%nmz,irec2) = exp_p(1:vacuum%nmz,irec2) / expDhg(irec2)
    end do

    ! case irec2 > 1:
    do irec2 = 2, stars%ng2
      do ivac = 1, vacuum%nvac
        ! integrands:
        fa(1:vacuum%nmzxy,irec2) = rhtxy(1:vacuum%nmzxy,irec2-1,ivac) * exp_m(1:vacuum%nmzxy,irec2)
        fb(1:vacuum%nmzxy,irec2) = rhtxy(1:vacuum%nmzxy,irec2-1,ivac) * exp_p(1:vacuum%nmzxy,irec2)
        ! integrals:
        ! alpha(z,q_xy,ivac) = int_z^infty rho(z',q_xy,ivac) exp(-sqrt(q_xy**2+prec_param**2)*z') dz'
        ! beta (z,q_xy,ivac) = int_{D/2}^z rho(z',q_xy,ivac) exp(+sqrt(q_xy**2+prec_param**2)*z') dz'
        ! where for z < 0 the lower vacuum charge density (ivac=2) is defined by rho(q_xy,z,ivac=2) := rho(q_xy,-z,ivac=2)
        call intgz1Reverse( fa(:,irec2), vacuum%delz, vacuum%nmzxy, alpha(:,irec2,ivac), .true. ) 
        call qsfComplex( vacuum%delz, fb(:,irec2), beta(:,irec2,ivac), vacuum%nmzxy, 1 )
        ! alphm(q_xy,ivac) = alpha(D/2,q_xy,ivac) --- these integrals are also needed for the interstitial potential
        alphm(irec2,ivac) = alpha(1,irec2,ivac)
      end do
      if ( vacuum%nvac == 1 ) then
        if ( sym%invs ) then
          alphm(irec2,2) = conjg( alphm(irec2,1) )
        else
          alphm(irec2,2) = alphm(irec2,1)
        end if
      end if
      do ivac = 1, vacuum%nvac
        gamma(1:vacuum%nmzxy,irec2) = exp_m(1:vacuum%nmzxy,irec2) * ( alphm(irec2,mod(ivac,2)+1) + beta(1:vacuum%nmzxy,irec2,ivac) ) &
                                    + exp_p(1:vacuum%nmzxy,irec2) *                               alpha(1:vacuum%nmzxy,irec2,ivac) ! mod(ivac,2)+1 outputs the other ivac value 
        where ( 2. * gamma(:,irec2) == gamma(:,irec2) ) gamma(:,irec2) = cmplx( 0., 0. )
        VVxy(1:vacuum%nmzxy,irec2,ivac) = VVxy(1:vacuum%nmzxy,irec2,ivac) + vcons(irec2) * gamma(1:vacuum%nmzxy,irec2) 
      end do
    end do

    ! case irec2 = 1:
    do ivac = 1, vacuum%nvac
      ga(1:vacuum%nmz) = rht(1:vacuum%nmz,ivac) * exp_m(1:vacuum%nmz,1)
      gb(1:vacuum%nmz) = rht(1:vacuum%nmz,ivac) * exp_p(1:vacuum%nmz,1)
      call intgz1Reverse( ga(:), vacuum%delz, vacuum%nmz, delta(:,ivac), .true. ) ! integrals 
      call qsf( vacuum%delz, gb(:), epsilon(:,ivac), vacuum%nmz, 1 )
      alphm(1,ivac) = delta(1,ivac)
    end do
    if ( vacuum%nvac == 1 ) alphm(1,2) = alphm(1,1)
    do ivac = 1, vacuum%nvac
      zeta(1:vacuum%nmz) = exp_m(1:vacuum%nmz,1) * ( alphm(1,mod(ivac,2)+1) + epsilon(1:vacuum%nmz,ivac) ) &
                         + exp_p(1:vacuum%nmz,1) *                              delta(1:vacuum%nmz,ivac) 
      where ( 2. * zeta == zeta ) zeta = 0.
      VVz(1:vacuum%nmz,ivac) = VVz(1:vacuum%nmz,ivac) + vcons(1) * zeta(1:vacuum%nmz)
    end do

    ! damping in vacuum:
    do ivac = 1, vacuum%nvac
      VVz(1:vacuum%nmz,ivac) = VVz(1:vacuum%nmz,ivac) * exp( -0.1 / rmt * z(1:vacuum%nmz) )
      do irec2 = 2, stars%ng2
        VVxy(1:vacuum%nmzxy,irec2,ivac) = VVxy(1:vacuum%nmzxy,irec2,ivac) * exp( -0.1 / rmt * z(1:vacuum%nmzxy) )
      end do
    end do


  end subroutine VYukawaFilmVacuumVariant1



  subroutine VYukawaFilmInterstitialVariant1( &
                stars, vacuum, cell, sym, input, &
                psq, VVxy, VVz, alphm, &
                VIq )

    ! main parts:
    ! 1. part: Compute the contribution from the interstitial charge density to the interstitial potential as a function of q_xy and z (analytic expression for integral)
    ! 2. part: Add the contribution from the vacuum charge density to the interstitial potential, which had already been computed earlier for the vacuum potential
 
    ! 4. part: Compute the coefficients V^I(q_xy,q_z) from the function V^I(q_xy,z) by a 1D Fourier transform:
    !          V^I(q_xy,z) = sum_{q_z} V^I(q_xy,q_z) * exp( ImagUnit * q_z * z )
    ! In order to be able to match the interstitial and vacuum potentials smoothly at the interface, the Fourier transform is done
    ! in a slightly larger region. -> 3. part
    ! 3. part: Interpolate the vacuum potential in a small region surrounding the slab
297 298 299 300 301 302 303

    use m_ExpSave
    use m_constants
    use m_types
    use m_cfft
    implicit none

304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324
    type(t_stars),  intent(in)  :: stars
    type(t_vacuum), intent(in)  :: vacuum
    type(t_cell),   intent(in)  :: cell
    type(t_sym),    intent(in)  :: sym
    type(t_input),  intent(in)  :: input
    complex,        intent(in)  :: psq(stars%ng3)
    complex,        intent(in)  :: VVxy(vacuum%nmzxyd,2:stars%ng2,2)
    real,           intent(in)  :: VVz(vacuum%nmzd,2)
    complex,        intent(in)  :: alphm(stars%ng2,2)

    complex,        intent(out) :: VIq(stars%ng3)

    real                        :: partitioning, rz, qz, q
    integer                     :: irec2, irec3, iz, jz, ivac, iqz, nfft, nzmax, nzmin, nzdh, nLower, nUpper, jvac
    complex, allocatable        :: VIz(:,:), eta(:,:)
    complex                     :: VIqz(-stars%mx3:stars%mx3,stars%ng2), c_ph(-stars%mx3:stars%mx3,stars%ng2)
    complex                     :: vcons1(stars%ng3)
    real                        :: VIzReal(3*stars%mx3,stars%ng2), VIzImag(3*stars%mx3,stars%ng2)
    real, allocatable           :: exp_m(:,:), exp_p(:,:)
    real, allocatable           :: z(:)
    real                        :: g_damped(stars%ng2), vcons2(stars%ng2), expDhg(stars%ng2)
325 326 327 328
    

    ! DEFINITIONS / ALLOCATIONS / INITIALISATIONS

329 330
    ! grid points z_i:
    nfft = 3 * stars%mx3                                  ! number of grid points for Fourier transform
331
    partitioning = 1. / real( nfft )
332 333 334
    nzmax = nfft / 2                                      ! index of maximal z below D~/2 
    nzmin = nzmax - nfft + 1                              ! index of minimal z above -D~/2
    nzdh = ceiling( cell%z1 / cell%amat(3,3) * nfft ) - 1 ! index of maximal z below D/2
335
    allocate( z(nzmin:nzmax) )
336
    ! definition of z_i:        ! indexing: z_0 = 0; positive indices for z > 0; negative indices for z < 0
337 338
    do iz = nzmin, nzmax
      z(iz) = cell%amat(3,3) * iz * partitioning
339
    end do
340
    ! other variables:
341 342
    allocate( VIz(nzmin:nzmax,stars%ng2), eta(nzmin:nzmax,stars%ng2) )
    allocate( exp_m(nzmin:nzmax,stars%ng2), exp_p(nzmin:nzmax,stars%ng2) )
343 344
    do irec2 = 1, stars%ng2
      g_damped(irec2) = sqrt( stars%sk2(irec2) ** 2 + input%preconditioning_param ** 2 )
345 346 347 348 349 350
      vcons2(irec2) = -1. / ( 2. * g_damped(irec2) )
      do iz = nzmin, nzmax
        exp_m(iz,irec2) = exp_save( - g_damped(irec2) * z(iz) )
        exp_p(iz,irec2) = exp_save(   g_damped(irec2) * z(iz) )
      end do
      expDhg(irec2) = exp_save( - cell%z1 * g_damped(irec2) )
351 352 353 354 355 356 357 358 359
    end do
    do irec3 = 1, stars%ng3
      vcons1(irec3) = fpi_const * psq(irec3) / ( stars%sk3(irec3) ** 2 + input%preconditioning_param ** 2 )
    end do
    VIz = (0.,0.)
    VIq = (0.,0.)


    ! CONTRIBUTION FROM THE INTERSTITIAL CHARGE DENSITY
360
    
361
    ! compute V^I(q_xy,z) as a function of q_xy and z
362
    do irec2 = 1, stars%ng2  
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
      do iz = -nzdh, nzdh
        do iqz = -stars%mx3, stars%mx3
          irec3 = stars%ig(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
          if ( irec3 /= 0 ) then ! use only stars within the g_max sphere
            c_ph(iqz,irec2) = stars%rgphs(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
            qz = iqz * cell%bmat(3,3)
            VIz(iz,irec2) = VIz(iz,irec2) &
                          + vcons1(irec3) * c_ph(iqz,irec2) * &
                            ( exp( ImagUnit * qz * z(iz) ) &
                            + vcons2(irec2) * expDhg(irec2) * &
                              ( ( g_damped(irec2) + ImagUnit * qz ) * exp_p(iz,irec2) * exp(   ImagUnit * qz * cell%z1 ) &
                              + ( g_damped(irec2) - ImagUnit * qz ) * exp_m(iz,irec2) * exp( - ImagUnit * qz * cell%z1 ) ) )
          end if
        enddo
      end do
    ! irec2 loop continues


    ! CONTRIBUTION FROM THE VACUUM CHARGE DENSITY

    ! irec2 loop continues
      eta(-nzdh:nzdh,irec2) = exp_m(-nzdh:nzdh,irec2) * alphm(irec2,2) + exp_p(-nzdh:nzdh,irec2) * alphm(irec2,1)
      where ( 2.0 * eta(:,irec2) == eta(:,irec2) ) eta(:,irec2) = cmplx( 0.0, 0.0 )
      VIz(-nzdh:nzdh,irec2) = VIz(-nzdh:nzdh,irec2) + tpi_const / g_damped(irec2) * eta(-nzdh:nzdh,irec2)
    ! irec2 loop continues
    
    
    ! INTERPOLATION IN VACUUM REGION

    ! use Lagrange polynomials of order 3 to interpolate the vacuum potential outside I
    ! q, q-1 and q-2 (scaled with +/-1 or +/-0.5) are the factors of the Lagrange basis polynomials
    ! irec2 loop continues
      do ivac = 1, 2
        select case( ivac )
          case( 1 )
            nUpper = nzmax; nLower =  nzdh + 1; jvac = 1
          case( 2 )
            nLower = nzmin; nUpper = -nzdh - 1; jvac = 2; if ( sym%invs .or. sym%zrfs ) jvac = 1
        end select
        do iz = nLower, nUpper
          rz = ( abs( z(iz) ) - cell%z1 ) / vacuum%delz + 1.0
          jz = rz      ! index of maximal vacuum grid point below z_i
          q = rz - jz  ! factor in Lagrange basis polynomials
          if ( irec2 == 1 ) then
            VIz(iz,irec2) = 0.5     * ( q - 1. ) * ( q - 2. ) * VVz(jz,  jvac) &
                          -       q *              ( q - 2. ) * VVz(jz+1,jvac) &
                          + 0.5 * q * ( q - 1. )              * VVz(jz+2,jvac)
          else if ( jz + 2 <= vacuum%nmzxy ) then
            VIz(iz,irec2) = 0.5 *     ( q - 1. ) * ( q - 2. ) * VVxy(jz,  irec2,jvac) &
                          -       q              * ( q - 2. ) * VVxy(jz+1,irec2,jvac) &
                          + 0.5 * q * ( q - 1. )              * VVxy(jz+2,irec2,jvac)
            if ( ( sym%invs .and. .not. sym%zrfs ) .and. ivac == 2 ) VIz(iz,irec2) = conjg( VIz(iz,irec2) )
          end if
        end do
      end do
    end do ! irec2


    ! 1D FOURIER TRANSFORM TO FIND THE COEFFICIENTS V^I(q_xy,q_z)

    ! change the indexing for the subroutine cfft, and split real and imaginary parts
    VIzReal(1:nzmax+1,:) =  real( VIz(0:nzmax,:) ); VIzReal(nzmax+2:nfft,:) =  real( VIz(nzmin:-1,:) )
    VIzImag(1:nzmax+1,:) = aimag( VIz(0:nzmax,:) ); VIzImag(nzmax+2:nfft,:) = aimag( VIz(nzmin:-1,:) )
 
    ! V^I(q_xy,z) = sum_{q_z} V^I(q_xy,q_z) * exp( ImagUnit * q_z * z )
    do irec2 = 1, stars%ng2
      call cfft( VIzReal(:,irec2), VIzImag(:,irec2), nfft, nfft, nfft, -1 )
    ! irec2 loop continues

    ! reorder
      VIqz(0,irec2) = cmplx( VIzReal(1,irec2), VIzImag(1,irec2) )
      do iqz = 1, stars%mx3
        VIqz( iqz,irec2) = cmplx( VIzReal(iqz+1,     irec2), VIzImag(iqz+1,     irec2) )
        VIqz(-iqz,irec2) = cmplx( VIzReal(nfft+1-iqz,irec2), VIzImag(nfft+1-iqz,irec2) )
      end do
   
    ! add the computed components to V^I(q_xy,q_z):
      do iqz= -stars%mx3, stars%mx3
        irec3 = stars%ig(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
        if ( irec3 /= 0 ) VIq(irec3) = VIq(irec3) + VIqz(iqz,irec2) * partitioning / ( stars%nstr(irec3) / stars%nstr2(irec2) )
      end do
    end do


  end subroutine VYukawaFilmInterstitialVariant1



  subroutine VYukawaFilmVacuumVariant2( &
                stars, vacuum, cell, sym, input, rmt, &
                psq, rhtxy, rht, &
                VVxy, VVz, alphm, dh_prec, coshdh )

    ! 1. part: Compute the contribution from the interstitial charge density to the vacuum potential as a function of q_xy and z (analytic expression for integral)
    ! 2. part: Compute the contribution from the vacuum charge density to the vacuum potential as a function of q_xy and z by numerical integration

    use m_ExpSave
    use m_constants
    use m_types
    use m_intgr, only: intgz1Reverse
    use m_qsf
    implicit none

    type(t_stars),  intent(in)  :: stars
    type(t_vacuum), intent(in)  :: vacuum
    type(t_cell),   intent(in)  :: cell
    type(t_sym),    intent(in)  :: sym
    type(t_input),  intent(in)  :: input
    real,           intent(in)  :: rmt
    complex,        intent(in)  :: psq(stars%ng3)
    complex,        intent(in)  :: rhtxy(vacuum%nmzxyd,stars%ng2-1,2)
    real,           intent(in)  :: rht(vacuum%nmzd,2) 

    complex,        intent(out) :: VVxy(vacuum%nmzxyd,2:stars%ng2,2) ! this is the qxy /= 0 part of the vacuum potential
    real,           intent(out) :: VVz(vacuum%nmzd,2)                ! this is the qxy  = 0 part of the vacuum potential
    complex,        intent(out) :: alphm(stars%ng2,2)                ! these are the integrals in upper and lower vacuum, now including the first star---integral for star ig2 is in alphm(ig2,ivac) 
    real,           intent(out) :: dh_prec
    real,           intent(out) :: coshdh(stars%ng2)

    integer                                             :: iz, irec2, irec3, ivac, iqz, nzdhprec, sign
    real                                                :: dh, qz, qxy_numerics 
    real, allocatable                                   :: z(:)
    complex, dimension(stars%ng2)                       :: sum_qz
    complex, dimension(stars%ng3)                       :: vcons3
    real, dimension(stars%ng2)                          :: g_damped, vcons2
    complex, dimension(-stars%mx3:stars%mx3,stars%ng2)  :: c_ph, quotq
    real, allocatable                                   :: quotz(:,:), sinhz(:,:)
    real, dimension(-1:1,stars%ng2)                     :: quotvardh
    complex, dimension(-stars%mx3:stars%mx3)            :: expdhqz
    complex, allocatable                                :: fa(:,:), fb(:,:)
    complex, allocatable                                :: alpha(:,:,:), beta(:,:,:), gamma(:,:)
    real, allocatable                                   :: ga(:), gb(:)
    real, allocatable                                   :: delta(:,:), epsilon(:,:), zeta(:)


    ! DEFINITIONS / ALLOCATIONS / INITIALISATIONS    
    
    dh = cell%z1 ! half the thickness of the film
    nzdhprec = ceiling( rmt / vacuum%delz ) ! index of dh_prec, see below
    dh_prec = dh + ( nzdhprec - 1 ) * vacuum%delz ! dh_prec is about dh + R; preconditioning boundary
    qxy_numerics = sqrt( ( 100 / dh_prec ) ** 2 - input%preconditioning_param ** 2 )
    allocate( z(nzdhprec) )
    allocate( quotz(-nzdhprec:nzdhprec,stars%ng2) )
    allocate( sinhz(nzdhprec,stars%ng2) )
    do iz = 1, nzdhprec ! new boundary
      z(iz) = ( iz - 1 ) * vacuum%delz + dh
    end do
    do irec2 = 1, stars%ng2
      g_damped(irec2) = sqrt( stars%sk2(irec2) ** 2 + input%preconditioning_param ** 2 )
      vcons2(irec2) = fpi_const / g_damped(irec2)
      coshdh(irec2) = cosh( g_damped(irec2) * ( dh_prec - dh ) )
514
      if( stars%sk2(irec2) < qxy_numerics ) then ! numerics ok
515 516 517 518 519
        do iz = 1, nzdhprec
          quotz( iz,irec2) = ( cosh( g_damped(irec2) * z(iz) ) / cosh( g_damped(irec2) * dh_prec ) &
                             + sinh( g_damped(irec2) * z(iz) ) / sinh( g_damped(irec2) * dh_prec ) ) / 2
          quotz(-iz,irec2) = ( cosh( g_damped(irec2) * z(iz) ) / cosh( g_damped(irec2) * dh_prec ) &
                             - sinh( g_damped(irec2) * z(iz) ) / sinh( g_damped(irec2) * dh_prec ) ) / 2
520
        end do
521 522 523 524
        quotvardh( 1,irec2) = ( cosh( g_damped(irec2) * dh ) / sinh( g_damped(irec2) * dh_prec ) &
                              + sinh( g_damped(irec2) * dh ) / cosh( g_damped(irec2) * dh_prec ) ) / 2
        quotvardh(-1,irec2) = ( cosh( g_damped(irec2) * dh ) / sinh( g_damped(irec2) * dh_prec ) &
                              - sinh( g_damped(irec2) * dh ) / cosh( g_damped(irec2) * dh_prec ) ) / 2
525
      else ! numerical treatment necessary
526 527 528
        do iz = 1, nzdhprec
          quotz( iz,irec2) = exp_save( g_damped(irec2) * (  z(iz) - dh_prec ) )
          quotz(-iz,irec2) = exp_save( g_damped(irec2) * ( -z(iz) - dh_prec ) )
529
        end do
530 531
        quotvardh( 1,irec2) = quotz( 1,irec2)
        quotvardh(-1,irec2) = quotz(-1,irec2)
532
      end if
533 534 535 536 537 538 539 540 541 542 543 544
      do iz = 1, nzdhprec
        sinhz(iz,irec2) = sinh( g_damped(irec2) * ( dh_prec - z(iz) ) )
      end do
      do iqz = -stars%mx3, stars%mx3
        quotq(iqz,irec2) = ImagUnit * iqz * cell%bmat(3,3) / g_damped(irec2)
      end do
    end do
    sum_qz = (0.,0.)
    VVxy = (0.,0.)
    VVz = 0.
    do irec3 = 1, stars%ng3
      vcons3(irec3) = fpi_const * psq(irec3) / ( stars%sk3(irec3) ** 2 + input%preconditioning_param ** 2 )
545
    end do
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 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 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 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793
    do iqz = -stars%mx3, stars%mx3
      expdhqz(iqz) = exp( ImagUnit * iqz * cell%bmat(3,3) * dh )
    end do


    ! CONTRIBUTION FROM THE INTERSTITIAL CHARGE DENSITY

    do ivac = 1, vacuum%nvac
      sign = 3 - 2 * ivac
      do irec2 = 1, stars%ng2
        do iqz = -stars%mx3, stars%mx3
          irec3 = stars%ig(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
          ! use only stars within the g_max sphere -> stars outside the sphere have per definition index ig3n = 0
          if ( irec3 /= 0 ) then
            c_ph(iqz,irec2) = stars%rgphs(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
            sum_qz(irec2) = sum_qz(irec2) + &
                            c_ph(iqz,irec2) * vcons3(irec3) * &
                            ( ( quotvardh( 1,irec2) - sign * quotq(iqz,irec2) * quotz( 1,irec2) ) * expdhqz( sign*iqz) &
                            - ( quotvardh(-1,irec2) - sign * quotq(iqz,irec2) * quotz(-1,irec2) ) * expdhqz(-sign*iqz) )
          endif
        enddo
        if( irec2 /= 1 ) then
          VVxy(1:nzdhprec,irec2,ivac) = sum_qz(irec2) * sign * sinhz(1:nzdhprec,irec2)
        else
          VVz( 1:nzdhprec,      ivac) = sum_qz(1)     * sign * sinhz(1:nzdhprec,1)
        end if
      enddo
    enddo


    ! CONTRIBUTION FROM THE VACUUM CHARGE DENSITY

    allocate( fa(nzdhprec,2:stars%ng2), fb(nzdhprec,2:stars%ng2) )
    allocate( alpha(nzdhprec,2:stars%ng2,2), beta(nzdhprec,2:stars%ng2,2), gamma(nzdhprec,2:stars%ng2) )
    ! case irec2 > 1:
    do irec2 = 2, stars%ng2
      do ivac = 1, vacuum%nvac
        ! integrands:
        fa(1:nzdhprec,irec2) = rhtxy(1:nzdhprec,irec2-1,ivac) * sinhz(1:nzdhprec,irec2)
        fb(1:nzdhprec,irec2) = rhtxy(1:nzdhprec,irec2-1,ivac) * quotz(1:nzdhprec,irec2)
        ! integrals:
        ! alpha(z,q_xy,ivac) = int_z^infty rho(z',q_xy,ivac) sinhz dz'
        ! beta (z,q_xy,ivac) = int_{D/2}^z rho(z',q_xy,ivac) quotz dz'
        ! where for z < 0 the lower vacuum charge density (ivac=2) is defined by rho(q_xy,z,ivac=2) := rho(q_xy,-z,ivac=2)
        call intgz1Reverse( fa(:,irec2), vacuum%delz, nzdhprec, alpha(:,irec2,ivac), .false. ) 
        call qsfComplex( vacuum%delz, fb(:,irec2), beta(:,irec2,ivac), nzdhprec, 1 )
        ! alphm(q_xy,ivac) = alpha(D/2,q_xy,ivac) --- these integrals are also needed for the interstitial potential
        alphm(irec2,ivac) = alpha(1,irec2,ivac)
      end do
      if ( vacuum%nvac == 1 ) then
        if ( sym%invs ) then
          alphm(irec2,2) = conjg( alphm(irec2,1) )
        else
          alphm(irec2,2) = alphm(irec2,1)
        end if
      end if
      do ivac = 1, vacuum%nvac
        gamma(1:nzdhprec,irec2) = quotz(-1:-nzdhprec:-1,irec2) * alphm(irec2,mod(ivac,2)+1) &
                                + quotz(1:nzdhprec,irec2) * alpha(1:nzdhprec,irec2,ivac) &
                                + sinhz(1:nzdhprec,irec2) * beta(1:nzdhprec,irec2,ivac)
        VVxy(1:nzdhprec,irec2,ivac) = VVxy(1:nzdhprec,irec2,ivac) + vcons2(irec2) * gamma(1:nzdhprec,irec2) 
      end do
    end do


    allocate( ga(nzdhprec), gb(nzdhprec) )
    allocate( delta(nzdhprec,2), epsilon(nzdhprec,2), zeta(nzdhprec) )
    ! case irec2 = 1:
    do ivac = 1, vacuum%nvac
      ga(1:nzdhprec) = rht(1:nzdhprec,ivac) * sinhz(1:nzdhprec,1)
      gb(1:nzdhprec) = rht(1:nzdhprec,ivac) * quotz(1:nzdhprec,1)
      call intgz1Reverse( ga(:), vacuum%delz, nzdhprec, delta(:,ivac), .false. ) 
      call qsf( vacuum%delz, gb(:), epsilon(:,ivac), nzdhprec, 1 )
      alphm(1,ivac) = delta(1,ivac)
    end do
    if ( vacuum%nvac == 1 ) alphm(1,2) = alphm(1,1)
    do ivac = 1, vacuum%nvac
      zeta(1:nzdhprec) = quotz(-1:-nzdhprec:-1,1) * alphm(1,mod(ivac,2)+1) &
                       + quotz(1:nzdhprec,1) * delta(1:nzdhprec,ivac) &
                       + sinhz(1:nzdhprec,1) * epsilon(1:nzdhprec,ivac) 
      VVz(1:nzdhprec,ivac) = VVz(1:nzdhprec,ivac) + vcons2(1) * zeta(1:nzdhprec)
    end do


  end subroutine VYukawaFilmVacuumVariant2



  subroutine VYukawaFilmInterstitialVariant2( &
                stars, vacuum, cell, sym, input, &
                psq, VVxy, VVz, alphm, dh_prec, coshdh, &
                VIq )

    ! main parts:
    ! 1. part: Compute the contribution from the interstitial charge density to the interstitial potential as a function of q_xy and z (analytic expression for integral)
    ! 2. part: Add the contribution from the vacuum charge density to the interstitial potential, which had already been computed earlier for the vacuum potential
 
    ! 4. part: Compute the coefficients V^I(q_xy,q_z) from the function V^I(q_xy,z) by a 1D Fourier transform:
    !          V^I(q_xy,z) = sum_{q_z} V^I(q_xy,q_z) * exp( ImagUnit * q_z * z )
    ! In order to be able to match the interstitial and vacuum potentials smoothly at the interface, the Fourier transform is done
    ! in a slightly larger region. -> 3. part
    ! 3. part: Interpolate the vacuum potential in a small region surrounding the slab

    use m_ExpSave
    use m_constants
    use m_types
    use m_cfft
    implicit none

    type(t_stars),  intent(in)  :: stars
    type(t_vacuum), intent(in)  :: vacuum
    type(t_cell),   intent(in)  :: cell
    type(t_sym),    intent(in)  :: sym
    type(t_input),  intent(in)  :: input
    complex,        intent(in)  :: psq(stars%ng3)
    complex,        intent(in)  :: VVxy(vacuum%nmzxyd,2:stars%ng2,2)
    real,           intent(in)  :: VVz(vacuum%nmzd,2)
    complex,        intent(in)  :: alphm(stars%ng2,2)
    real,           intent(in)  :: dh_prec
    real,           intent(in)  :: coshdh(stars%ng2)

    complex,        intent(out) :: VIq(stars%ng3)

    real                                               :: partitioning, rz, qz, q, qxy_numerics
    integer                                            :: irec2, irec3, iz, jz, ivac, iqz, jvac
    integer                                            :: nfft, nzmax, nzmin, nzdh, nLower, nUpper
    complex, allocatable                               :: VIz(:,:), eta(:,:)
    complex, allocatable                               :: expzqz(:,:)
    complex, dimension(-stars%mx3:stars%mx3,stars%ng2) :: VIqz, c_ph
    complex, dimension(-stars%mx3:stars%mx3,stars%ng2) :: qquot, qquottrigp, qquottrigm
    complex, dimension(stars%ng3)                      :: vcons3
    real, dimension(3*stars%mx3,stars%ng2)             :: VIzReal, VIzImag
    real, allocatable                                  :: quotz(:,:), sinhz(:,:)
    real, allocatable                                  :: z(:)
    real, dimension(stars%ng2)                         :: g_damped, vcons2
    

    ! DEFINITIONS / ALLOCATIONS / INITIALISATIONS

    ! grid points z_i:
    qxy_numerics = sqrt( ( 100 / dh_prec ) ** 2 - input%preconditioning_param ** 2 )
    nfft = 3 * stars%mx3                                  ! number of grid points for Fourier transform
    partitioning = 1. / real( nfft )
    nzmax = nfft / 2                                      ! index of maximal z below D~/2 
    nzmin = nzmax - nfft + 1                              ! index of minimal z above -D~/2
    nzdh = ceiling( cell%z1 / cell%amat(3,3) * nfft ) - 1 ! index of maximal z below D/2
    allocate( z(nzmin:nzmax) )
    ! definition of z_i:        ! indexing: z_0 = 0; positive indices for z > 0; negative indices for z < 0
    do iz = nzmin, nzmax
      z(iz) = cell%amat(3,3) * iz * partitioning
    end do
    ! other variables:
    allocate( VIz(nzmin:nzmax,stars%ng2) )
    allocate( eta(-nzdh:nzdh,stars%ng2) )
    allocate( sinhz(-nzdh:nzdh,stars%ng2) )
    allocate( quotz(-nzdh:nzdh,stars%ng2) )
    allocate( expzqz(-nzdh:nzdh,-stars%mx3:stars%mx3) )
    do irec2 = 1, stars%ng2
      g_damped(irec2) = sqrt( stars%sk2(irec2) ** 2 + input%preconditioning_param ** 2 )
      vcons2(irec2) = fpi_const / g_damped(irec2)
      if( stars%sk2(irec2) < qxy_numerics ) then ! numerics ok
        do iz = -nzdh, nzdh
          quotz( iz,irec2) = ( cosh( g_damped(irec2) * z(iz) ) / cosh( g_damped(irec2) * dh_prec ) &
                             + sinh( g_damped(irec2) * z(iz) ) / sinh( g_damped(irec2) * dh_prec ) ) / 2
        end do
      else ! numerical treatment necessary
        do iz = -nzdh, nzdh
          quotz( iz,irec2) = exp_save( g_damped(irec2) * (  z(iz) - dh_prec ) )
        end do
      end if
      do iz = -nzdh, nzdh
        sinhz(iz,irec2) = sinh( g_damped(irec2) * ( dh_prec - z(iz) ) )
        do iqz = -stars%mx3, stars%mx3
          expzqz(iz,iqz) = exp( ImagUnit * iqz * cell%bmat(3,3) * z(iz) )
        end do
      end do
      do iqz = -stars%mx3, stars%mx3
        qquot(iqz,irec2) = ImagUnit * iqz * cell%bmat(3,3) / g_damped(irec2)
        qquottrigp(iqz,irec2) = coshdh(irec2) + sinhz(nzdh,irec2) * qquot(iqz,irec2)
        qquottrigm(iqz,irec2) = coshdh(irec2) - sinhz(nzdh,irec2) * qquot(iqz,irec2)
      end do
    end do
    do irec3 = 1, stars%ng3
      vcons3(irec3) = fpi_const * psq(irec3) / ( stars%sk3(irec3) ** 2 + input%preconditioning_param ** 2 )
    end do
    VIz = (0.,0.)
    VIq = (0.,0.)


    ! CONTRIBUTION FROM THE INTERSTITIAL CHARGE DENSITY
    
    ! compute V^I(q_xy,z) as a function of q_xy and z
    do irec2 = 1, stars%ng2  
      do iz = -nzdh, nzdh
        do iqz = -stars%mx3, stars%mx3
          irec3 = stars%ig(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
          if ( irec3 /= 0 ) then ! use only stars within the g_max sphere
            c_ph(iqz,irec2) = stars%rgphs(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
            qz = iqz * cell%bmat(3,3)
            VIz(iz,irec2) = VIz(iz,irec2) &
                          + vcons3(irec3) * c_ph(iqz,irec2) * &
                            ( expzqz(iz,iqz) &
                            - qquottrigp(iqz,irec2) * expzqz( nzdh,iqz) * quotz( iz,irec2) &
                            - qquottrigm(iqz,irec2) * expzqz(-nzdh,iqz) * quotz(-iz,irec2) )
          end if
        enddo
      end do
    ! irec2 loop continues


    ! CONTRIBUTION FROM THE VACUUM CHARGE DENSITY

    ! irec2 loop continues
      eta(-nzdh:nzdh,irec2) = quotz(nzdh:-nzdh:-1,irec2) * alphm(irec2,2) &
                            + quotz(-nzdh:nzdh,   irec2) * alphm(irec2,1)
      VIz(-nzdh:nzdh,irec2) = VIz(-nzdh:nzdh,irec2) + vcons2(irec2) * eta(-nzdh:nzdh,irec2)
    ! irec2 loop continues
   
    
    ! INTERPOLATION IN VACUUM REGION

    ! use Lagrange polynomials of order 3 to interpolate the vacuum potential outside I
    ! q, q-1 and q-2 (scaled with +/-1 or +/-0.5) are the factors of the Lagrange basis polynomials
    ! irec2 loop continues
      do ivac = 1, 2
        select case( ivac )
          case( 1 )
            nUpper = nzmax; nLower =  nzdh + 1; jvac = 1
          case( 2 )
            nLower = nzmin; nUpper = -nzdh - 1; jvac = 2; if ( sym%invs .or. sym%zrfs ) jvac = 1
        end select
        do iz = nLower, nUpper
          rz = ( abs( z(iz) ) - cell%z1 ) / vacuum%delz + 1.0
          jz = rz      ! index of maximal vacuum grid point below z_i
          q = rz - jz  ! factor in Lagrange basis polynomials
          if ( irec2 == 1 ) then
            VIz(iz,irec2) = 0.5     * ( q - 1. ) * ( q - 2. ) * VVz(jz,  jvac) &
                          -       q *              ( q - 2. ) * VVz(jz+1,jvac) &
                          + 0.5 * q * ( q - 1. )              * VVz(jz+2,jvac)
          else if ( jz + 2 <= vacuum%nmzxy ) then
            VIz(iz,irec2) = 0.5 *     ( q - 1. ) * ( q - 2. ) * VVxy(jz,  irec2,jvac) &
                          -       q              * ( q - 2. ) * VVxy(jz+1,irec2,jvac) &
                          + 0.5 * q * ( q - 1. )              * VVxy(jz+2,irec2,jvac)
            if ( ( sym%invs .and. .not. sym%zrfs ) .and. ivac == 2 ) VIz(iz,irec2) = conjg( VIz(iz,irec2) )
          end if
        end do
      end do
    end do ! irec2
794

795

796
    ! 1D FOURIER TRANSFORM TO FIND THE COEFFICIENTS V^I(q_xy,q_z)
797

798 799 800 801 802 803 804 805
    ! change the indexing for the subroutine cfft, and split real and imaginary parts
    VIzReal(1:nzmax+1,:) =  real( VIz(0:nzmax,:) ); VIzReal(nzmax+2:nfft,:) =  real( VIz(nzmin:-1,:) )
    VIzImag(1:nzmax+1,:) = aimag( VIz(0:nzmax,:) ); VIzImag(nzmax+2:nfft,:) = aimag( VIz(nzmin:-1,:) )
 
    ! V^I(q_xy,z) = sum_{q_z} V^I(q_xy,q_z) * exp( ImagUnit * q_z * z )
    do irec2 = 1, stars%ng2
      call cfft( VIzReal(:,irec2), VIzImag(:,irec2), nfft, nfft, nfft, -1 )
    ! irec2 loop continues
806

807
    ! reorder
808
      VIqz(0,irec2) = cmplx( VIzReal(1,irec2), VIzImag(1,irec2) )
809 810 811
      do iqz = 1, stars%mx3
        VIqz( iqz,irec2) = cmplx( VIzReal(iqz+1,     irec2), VIzImag(iqz+1,     irec2) )
        VIqz(-iqz,irec2) = cmplx( VIzReal(nfft+1-iqz,irec2), VIzImag(nfft+1-iqz,irec2) )
812
      end do
813 814 815 816 817
   
    ! add the computed components to V^I(q_xy,q_z):
      do iqz= -stars%mx3, stars%mx3
        irec3 = stars%ig(stars%kv2(1,irec2),stars%kv2(2,irec2),iqz)
        if ( irec3 /= 0 ) VIq(irec3) = VIq(irec3) + VIqz(iqz,irec2) * partitioning / ( stars%nstr(irec3) / stars%nstr2(irec2) )
818 819 820
      end do
    end do

821

822
  end subroutine VYukawaFilmInterstitialVariant2
823 824 825



826
  subroutine VYukawaModify( stars, vacuum, cell, sym, input, mpi, atoms, sphhar, oneD, noco, den, &
827 828 829 830 831 832 833
                            VYukawa )

    ! This subroutine adds a potential to the previously computed Yukawa
    ! potential to ensure charge neutrality.
    ! The added potential itself is a solution to the modified Helmholtz
    ! equation with constant right-hand side, where the constant is chosen such
    ! that charge neutrality is obtained.
834 835 836
    ! The charge is distributed only over the film region, and we therefore 
    ! solve the differential equation subject to a boundary condition on the 
    ! film surface, in contrast to the basic Yukawa potential above. 
837 838 839 840 841 842 843 844 845 846 847 848 849 850 851 852 853 854 855 856 857 858 859 860 861 862 863 864 865 866 867 868 869 870 871 872 873 874 875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 893 894 895 896 897 898 899 900 901

    use m_constants
    use m_types
    use m_vmts
    use m_constants
    use m_cdntot
    use m_cfft
    implicit none

    type(t_stars),      intent(in)    :: stars
    type(t_vacuum),     intent(in)    :: vacuum
    type(t_cell),       intent(in)    :: cell
    type(t_sym),        intent(in)    :: sym
    type(t_input),      intent(in)    :: input
    type(t_mpi),        intent(in)    :: mpi
    type(t_atoms),      intent(in)    :: atoms
    type(t_sphhar),     intent(in)    :: sphhar
    type(t_oneD),       intent(in)    :: oneD
    type(t_noco),       intent(in)    :: noco
    type(t_potden),     intent(inout) :: den
    type(t_potden),     intent(inout) :: VYukawa

    integer                           :: n, lh, irec3, iz, iqz, nfft, nzmax, nzmin, nzdh
    real                              :: q0, qhat, qbar, ldh, partitioning, dh
    real                              :: q(input%jspins), qis(input%jspins), qmt(atoms%ntype,input%jspins),qvac(2,input%jspins), qtot, qistot
    complex                           :: psq(stars%ng3)
    type(t_potden)                    :: VYukawaModification
    real, allocatable                 :: z(:)
    real                              :: VIzReal(3*stars%mx3), VIzImag(3*stars%mx3)
    complex, allocatable              :: VIz(:)
    complex                           :: VIqz(-stars%mx3:stars%mx3)


    ! DEFINITIONS / ALLOCATIONS / INITIALISATIONS

    ! constants:
    dh = cell%z1   ! half the width of the film
    ! indexing of grid points z_i:
    nfft = 3 * stars%mx3   ! number of grid points for Fourier transform
    partitioning = 1. / real( nfft )
    nzmax = nfft / 2           ! index of maximal z below D~/2 
    nzmin = nzmax - nfft + 1   ! index of minimal z above -D~/2
    nzdh = ceiling( dh / cell%amat(3,3) * nfft ) - 1   ! index of maximal z below D/2
    allocate( z(nzmin:nzmax) )
    ! definition of grid points z_i:
    ! indexing: z_0 = 0; positive indices for z > 0; negative indices for z < 0
    do iz = nzmin, nzmax
      z(iz) = cell%amat(3,3) * iz * partitioning
    end do


    ! INTEGRATION OF THE PREVIOUSLY COMPUTED YUKAWA POTENTIAL

    ! initialise VYukawaModification with in-going VYukawa and prepare for integration
    call VYukawaModification%init( stars, atoms, sphhar, vacuum, noco, input%jspins, 4 )
    call VYukawaModification%copyPotDen( VYukawa )
    do n = 1, atoms%ntype
      do lh = 0, sphhar%nlhd    
        VYukawaModification%mt(1:atoms%jri(n),lh,n,1) = VYukawaModification%mt(1:atoms%jri(n),lh,n,1) * atoms%rmsh(1:atoms%jri(n),n) ** 2
      end do
    end do

    ! integrate the potential over the film region
    call cdntot_integrate( stars, atoms, sym, vacuum, input, cell, oneD, VYukawaModification, q, qis, qmt, qvac, qtot, qistot  )
    q0 = qtot / cell%area
902
    ldh = input%preconditioning_param * dh
903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 956 957 958 959 960 961 962 963 964 965 966 967
    qhat = ( q0 / ( 2 * dh ) ) / ( sinh(ldh) / ( ldh * cosh( ldh ) ) - 1 )
    qbar = input%preconditioning_param ** 2 / fpi_const *  qhat


    ! SET UP CONSTANT CHARGE DENSITY

    ! instead of den%pw(1,1) = qbar we directly set the pseudo charge density
    den%mt = 0; den%pw = 0; den%vacxy = 0; den%vacz = 0    
    do n = 1, atoms%ntype
      den%mt(1:atoms%jri(n),0,n,1) = sfp_const * qbar * atoms%rmsh(1:atoms%jri(n),n) ** 2
    end do
    psq = cmplx(0.0,0.0); psq(1) = qbar


    ! CALCULATE THE INTERSTITIAL POTENTIAL AS A FUNCTION OF z

    ! initialise and calculate out-going modification potential; reuse VYukawaModification
    VYukawaModification%mt = 0; VYukawaModification%pw = 0; VYukawaModification%vacxy = 0; VYukawaModification%vacz = 0

    allocate( VIz(nzmin:nzmax) )
    VIz = (0.,0.)
    do iz = -nzdh+1, nzdh-1
      VIz(iz) = qhat * ( 1 - cosh( input%preconditioning_param * z(iz) ) / cosh( ldh ) ) 
    end do


    ! 1D FOURIER TRANSFORM TO FIND THE 3D-FOURIER COEFFICIENTS

    ! change the indexing for the subroutine cfft, and split real and imaginary parts
    VIzReal(1:nzmax+1) =  real( VIz(0:nzmax) ); VIzReal(nzmax+2:nfft) = real( VIz(nzmin:-1) )
    VIzImag(1:nzmax+1) = aimag( VIz(0:nzmax) ); VIzImag(nzmax+2:nfft) = aimag( VIz(nzmin:-1) )

    call cfft( VIzReal(:), VIzImag(:), nfft, nfft, nfft, -1 )

    ! reorder
    VIqz = 0
    VIqz(0) = cmplx( VIzReal(1), VIzImag(1) )
    do iqz = 1, stars%mx3
      VIqz( iqz) = cmplx( VIzReal(iqz+1     ), VIzImag(iqz+1     ) )
      VIqz(-iqz) = cmplx( VIzReal(nfft+1-iqz), VIzImag(nfft+1-iqz) )
    end do

    ! add the computed components
    do iqz= -stars%mx3, stars%mx3
      irec3 = stars%ig(stars%kv2(1,1),stars%kv2(2,1),iqz)
      if ( irec3 /= 0 ) VYukawaModification%pw(irec3,1) = VYukawaModification%pw(irec3,1) + VIqz(iqz) * partitioning / ( stars%nstr(irec3) / stars%nstr2(1) )
    end do


    ! MUFFIN-TIN POTENTIAL

    call Vmts( input, mpi, stars, sphhar, atoms, sym, cell, oneD, &
               VYukawaModification%pw(:,1), den%mt(:,0:,:,1), VYukawaModification%potdenType, &
               VYukawaModification%mt(:,0:,:,1) )


    ! APPLYING THE MODIFICATION TO THE YUKAWA POTENTIAL

    call VYukawa%AddPotDen( VYukawa, VYukawaModification )


  end subroutine VYukawaModify



968
end module m_VYukawaFilm