Skip to content

Commit 99e4a68

Browse files
committed
Add Skala specific rountines to Fortran API
1 parent 33b4ecb commit 99e4a68

File tree

1 file changed

+61
-0
lines changed

1 file changed

+61
-0
lines changed

src/xc_integrator.f90

Lines changed: 61 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -32,6 +32,7 @@ module gauxc_integrator
3232
& gauxc_integrator_eval_exc_gks, &
3333
& gauxc_integrator_eval_exc_vxc_rks, &
3434
& gauxc_integrator_eval_exc_vxc_uks, &
35+
& gauxc_integrator_eval_exc_vxc_onedft_uks, &
3536
& gauxc_integrator_eval_exc_vxc_gks
3637

3738
public :: &
@@ -282,6 +283,37 @@ subroutine gauxc_integrator_eval_exc_vxc_gks(status, integrator, &
282283
end subroutine gauxc_integrator_eval_exc_vxc_gks
283284
end interface gauxc_eval_exc_vxc
284285

286+
interface
287+
!> @brief Evaluate the exchange-correlation energy and potential for UKS.
288+
subroutine gauxc_integrator_eval_exc_vxc_onedft_uks_c(status, integrator, &
289+
density_matrix_s, density_matrix_z, model, exc, vxc_matrix_s, vxc_matrix_z) &
290+
bind(c, name="gauxc_integrator_eval_exc_vxc_onedft_uks")
291+
import :: gauxc_status_type, gauxc_integrator_type, gauxc_matrix_type, c_double, &
292+
& c_char
293+
implicit none
294+
!> @param status Status object to capture any errors.
295+
type(gauxc_status_type), intent(out) :: status
296+
!> @param integrator Handle to the XCIntegrator.
297+
type(gauxc_integrator_type), value :: integrator
298+
!> @param density_matrix_s Density matrix container for total density.
299+
type(gauxc_matrix_type), value :: density_matrix_s
300+
!> @param density_matrix_z Density matrix container for spin density.
301+
type(gauxc_matrix_type), value :: density_matrix_z
302+
!> @param String specifying the OneDFT model to use.
303+
character(kind=c_char), intent(in) :: model(*)
304+
!> @param exc Pointer to store the exchange-correlation energy.
305+
real(c_double), intent(out) :: exc
306+
!> @param vxc_matrix_s Matrix container for total density potential.
307+
type(gauxc_matrix_type), intent(out) :: vxc_matrix_s
308+
!> @param vxc_matrix_z Matrix container for spin density potential.
309+
type(gauxc_matrix_type), intent(out) :: vxc_matrix_z
310+
end subroutine gauxc_integrator_eval_exc_vxc_onedft_uks_c
311+
end interface
312+
313+
interface gauxc_eval_exc_vxc
314+
module procedure gauxc_integrator_eval_exc_vxc_onedft_uks
315+
end interface gauxc_eval_exc_vxc
316+
285317
contains
286318

287319
!> @brief Create a new XCIntegratorFactory instance.
@@ -321,4 +353,33 @@ function gauxc_integrator_factory_new(status, execution_space, &
321353
c_integrator_input_type, c_integrator_kernel_name, &
322354
c_local_work_kernel_name, c_reduction_kernel_name)
323355
end function gauxc_integrator_factory_new
356+
357+
!> @brief Evaluate the exchange-correlation energy and potential for UKS.
358+
subroutine gauxc_integrator_eval_exc_vxc_onedft_uks(status, integrator, &
359+
density_matrix_s, density_matrix_z, model, exc, vxc_matrix_s, vxc_matrix_z)
360+
!> @param status Status object to capture any errors.
361+
type(gauxc_status_type), intent(out) :: status
362+
!> @param integrator Handle to the XCIntegrator.
363+
type(gauxc_integrator_type), value :: integrator
364+
!> @param density_matrix_s Density matrix container for total density.
365+
type(gauxc_matrix_type), value :: density_matrix_s
366+
!> @param density_matrix_z Density matrix container for spin density.
367+
type(gauxc_matrix_type), value :: density_matrix_z
368+
!> @param String specifying the OneDFT model to use.
369+
character(kind=c_char, len=*), intent(in) :: model
370+
!> @param exc Pointer to store the exchange-correlation energy.
371+
real(c_double), intent(inout) :: exc
372+
!> @param vxc_matrix_s Matrix container for total density potential.
373+
type(gauxc_matrix_type), intent(inout) :: vxc_matrix_s
374+
!> @param vxc_matrix_z Matrix container for spin density potential.
375+
type(gauxc_matrix_type), intent(inout) :: vxc_matrix_z
376+
377+
character(kind=c_char), allocatable :: c_model(:)
378+
379+
c_model = transfer(model//c_null_char, [character(kind=c_char) ::], &
380+
& len(model)+1)
381+
382+
call gauxc_integrator_eval_exc_vxc_onedft_uks_c(status, integrator, &
383+
density_matrix_s, density_matrix_z, c_model, exc, vxc_matrix_s, vxc_matrix_z)
384+
end subroutine gauxc_integrator_eval_exc_vxc_onedft_uks
324385
end module gauxc_integrator

0 commit comments

Comments
 (0)