Octopus
xc.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
21module xc_oct_m
22 use comm_oct_m
23 use debug_oct_m
28 use global_oct_m
29 use grid_oct_m
30 use io_oct_m
32 use iso_c_binding
33 use, intrinsic :: iso_fortran_env
38 use math_oct_m
39 use mesh_oct_m
42 use mpi_oct_m
44 use parser_oct_m
47 use space_oct_m
51 use xc_f03_lib_m
52 use xc_fbe_oct_m
55
56 implicit none
57
58 private
59 public :: &
60 xc_t, &
61 xc_init, &
62 xc_end, &
64 xc_get_vxc, &
66 xc_get_fxc, &
67 xc_get_kxc, &
75
76 ! A Structure that contains the quantities needed to compute the functionals
78 real(real64), pointer :: rho(:,:) ! A pointer to the full density
79
80 real(real64), allocatable :: dens(:,:) ! Density (in the local frame of the magnetization)
81 real(real64), allocatable :: gdens(:,:,:) ! Gradient of the density
82 real(real64), allocatable :: ldens(:,:) ! Laplacian of the density
83 real(real64), allocatable :: tau(:,:) ! Kinetic energy density
85
86
87 type xc_t
88 private
89 integer, public :: family
90 integer, public :: flags
91 integer, public :: kernel_family
92 type(xc_functional_t), public :: functional(2,2)
94
95 type(xc_functional_t), public :: kernel(2,2)
96 real(real64), public :: kernel_lrc_alpha
97
98 real(real64), public :: cam_omega
99 real(real64), public :: cam_alpha
100 real(real64), public :: cam_beta
101 real(real64), public :: cam_ext(3)
102
103 logical :: use_gi_ked
104
105 integer :: xc_density_correction
106 logical :: xcd_optimize_cutoff
107 real(real64) :: xcd_ncutoff
108 logical :: xcd_minimum
109 logical :: xcd_normalize
110 logical :: parallel
111
112 type(internal_quantities_t) :: quantities
113 end type xc_t
115 real(real64), parameter :: tiny = 1.0e-12_real64
116
117 integer, parameter :: &
118 LR_NONE = 0, &
119 lr_x = 1
120
121contains
122
123 ! ---------------------------------------------------------
124 subroutine xc_write_info(xcs, iunit, namespace)
125 type(xc_t), intent(in) :: xcs
126 integer, optional, intent(in) :: iunit
127 type(namespace_t), optional, intent(in) :: namespace
128
129 integer :: ifunc
130
131 push_sub(xc_write_info)
132
133 write(message(1), '(a)') "Exchange-correlation:"
134 call messages_info(1, iunit=iunit, namespace=namespace)
135
136 do ifunc = func_x, func_c
137 call xc_functional_write_info(xcs%functional(ifunc, 1), iunit, namespace)
138 end do
139
140 if (abs(xcs%cam_alpha + xcs%cam_beta) > m_epsilon) then
141 write(message(1), '(1x)')
142 write(message(2), '(a,f8.5)') "Exact exchange mixing = ", xcs%cam_alpha
143 write(message(3), '(a,f8.5)') "Exact exchange for short-range beta = ", xcs%cam_beta
144 write(message(4), '(a,f8.5)') "Exact exchange range-separate omega = ", xcs%cam_omega
145 call messages_info(4, iunit=iunit, namespace=namespace)
146 end if
147
148
149 pop_sub(xc_write_info)
150 end subroutine xc_write_info
151
152
153 ! ---------------------------------------------------------
154 subroutine xc_init(xcs, namespace, ndim, periodic_dim, nel, x_id, c_id, xk_id, ck_id, hartree_fock, ispin)
155 type(xc_t), intent(out) :: xcs
156 type(namespace_t), intent(in) :: namespace
157 integer, intent(in) :: ndim
158 integer, intent(in) :: periodic_dim
159 real(real64), intent(in) :: nel
160 integer, intent(in) :: x_id
161 integer, intent(in) :: c_id
162 integer, intent(in) :: xk_id
163 integer, intent(in) :: ck_id
164 logical, intent(in) :: hartree_fock
165 integer, intent(in) :: ispin
166
167 integer :: isp, xc_major, xc_minor, xc_micro
168 logical :: ll
169 type(block_t) :: blk
171 push_sub(xc_init)
172
173 call xc_f03_version(xc_major, xc_minor, xc_micro)
175 xcs%family = 0
176 xcs%flags = 0
177 xcs%kernel_family = 0
178
179 call parse()
181 !we also need XC functionals that do not depend on the current
182 !get both spin-polarized and unpolarized
184 ! TODO: check that nel should not be spin polarized here
185 do isp = 1, 2
186
187 call xc_functional_init(xcs%functional(func_x, isp), namespace, x_id, ndim, nel, isp)
188 call xc_functional_init(xcs%functional(func_c, isp), namespace, c_id, ndim, nel, isp)
190 call xc_functional_init(xcs%kernel(func_x, isp), namespace, xk_id, ndim, nel, isp)
191 call xc_functional_init(xcs%kernel(func_c, isp), namespace, ck_id, ndim, nel, isp)
193 end do
195 xcs%family = ior(xcs%family, xcs%functional(func_x,1)%family)
196 xcs%family = ior(xcs%family, xcs%functional(func_c,1)%family)
197
198 xcs%flags = ior(xcs%flags, xcs%functional(func_x,1)%flags)
199 xcs%flags = ior(xcs%flags, xcs%functional(func_c,1)%flags)
201 xcs%kernel_family = ior(xcs%kernel_family, xcs%kernel(func_x,1)%family)
202 xcs%kernel_family = ior(xcs%kernel_family, xcs%kernel(func_c,1)%family)
204
205 if (xc_is_not_size_consistent(xcs, namespace) .and. periodic_dim > 0) then
206 message(1) = "Cannot perform a periodic calculation with a functional"
207 message(2) = "that depends on the number of electrons."
208 call messages_fatal(2, namespace=namespace)
209 end if
211 ! Take care of hybrid functionals (they appear in the correlation functional)
212 xcs%cam_omega = m_zero
213 xcs%cam_alpha = m_zero
214 xcs%cam_beta = m_zero
215
216 ll = (hartree_fock) &
217 .or.(xcs%functional(func_x,1)%id == xc_oep_x) &
218 .or. family_is_hybrid(xcs)
219 if (ll) then
220 if ((xcs%functional(func_x,1)%id /= 0).and.(xcs%functional(func_x,1)%id /= xc_oep_x)) then
221 message(1) = "You cannot use an exchange functional when performing"
222 message(2) = "a Hartree-Fock calculation or using a hybrid functional."
223 call messages_fatal(2, namespace=namespace)
224 end if
225
226 if (periodic_dim == ndim) then
227 call messages_experimental("Fock operator (Hartree-Fock, OEP, hybrids) in fully periodic systems", namespace=namespace)
228 end if
229
230 ! get the mixing coefficient for hybrids
231 if (family_is_hybrid(xcs)) then
232 if( any(abs(xcs%cam_ext) > m_epsilon) ) call set_hybrid_params(xcs,namespace)
233
234 call xc_f03_hyb_cam_coef(xcs%functional(func_c,1)%conf, xcs%cam_omega, &
235 xcs%cam_alpha, xcs%cam_beta)
236 call xc_f03_hyb_cam_coef(xcs%functional(func_c,2)%conf, xcs%cam_omega, &
237 xcs%cam_alpha, xcs%cam_beta)
238 else
239 ! we are doing Hartree-Fock plus possibly a correlation functional
240 xcs%cam_omega = m_zero
241 xcs%cam_alpha = m_one
242 xcs%cam_beta = m_zero
243 end if
244
245 ! reset certain variables
246 xcs%functional(func_x,1)%family = xc_family_oep
247 xcs%functional(func_x,1)%id = xc_oep_x
248 xcs%functional(func_x,2)%family = xc_family_oep
249 xcs%functional(func_x,2)%id = xc_oep_x
250 if (.not. hartree_fock) then
251 xcs%family = ior(xcs%family, xc_family_oep)
252 end if
253 end if
254
255 if (in_family(xcs%family, [xc_family_lca])) then
256 call messages_not_implemented("LCA current functionals", namespace) ! not even in libxc!
257 end if
258
259 call messages_obsolete_variable(namespace, 'MGGAimplementation')
260 call messages_obsolete_variable(namespace, 'CurrentInTau', 'XCUseGaugeIndependentKED')
261
262 if (xcs%functional(func_x, 1)%id == xc_mgga_x_tb09 .and. periodic_dim /= 3) then
263 message(1) = "mgga_x_tb09 functional can only be used for 3D periodic systems"
264 call messages_fatal(1, namespace=namespace)
265 end if
266
267 if (family_is_mgga(xcs%family) .or. family_is_nc_mgga(xcs%family)) then
268 !%Variable XCUseGaugeIndependentKED
269 !%Type logical
270 !%Default yes
271 !%Section Hamiltonian::XC
272 !%Description
273 !% If true, when evaluating the XC functional, a term including the (paramagnetic or total) current
274 !% is added to the kinetic-energy density such as to make it gauge-independent.
275 !% Applies only to meta-GGA (and hybrid meta-GGA) functionals.
276 !%End
277 call parse_variable(namespace, 'XCUseGaugeIndependentKED', .true., xcs%use_gi_ked)
278 end if
279
280 pop_sub(xc_init)
281
282 contains
283
284 subroutine parse()
285
286 push_sub(xc_init.parse)
287
288 ! the values of x_id, c_id, xk_id, and c_id are read outside the routine
289
290 !%Variable XCKernelLRCAlpha
291 !%Type float
292 !%Default 0.0
293 !%Section Hamiltonian::XC
294 !%Description
295 !% Set to a non-zero value to add a long-range correction for solids to the kernel.
296 !% This is the <math>\alpha</math> parameter defined in S. Botti <i>et al.</i>, <i>Phys. Rev. B</i>
297 !% 69, 155112 (2004). The <math>\Gamma = \Gamma` = 0</math> term <math>-\alpha/q^2</math> is taken
298 !% into account by introducing an additional pole to the polarizability (see R. Stubner
299 !% <i>et al.</i>, <i>Phys. Rev. B</i> 70, 245119 (2004)). The rest of the terms are included by
300 !% multiplying the Hartree term by <math>1 - \alpha / 4 \pi</math>. The use of non-zero
301 !% <math>\alpha</math> in combination with <tt>HamiltonianVariation</tt> = <tt>V_ext_only</tt>
302 !% corresponds to account of only the <math>\Gamma = \Gamma` = 0</math> term.
303 !% Applicable only to isotropic systems. (Experimental)
304 !%End
305
306 call parse_variable(namespace, 'XCKernelLRCAlpha', m_zero, xcs%kernel_lrc_alpha)
307 if (abs(xcs%kernel_lrc_alpha) > m_epsilon) then
308 call messages_experimental("Long-range correction to kernel", namespace=namespace)
309 end if
310
311 !%Variable XCDensityCorrection
312 !%Type integer
313 !%Default none
314 !%Section Hamiltonian::XC::DensityCorrection
315 !%Description
316 !% This variable controls the long-range correction of the XC
317 !% potential using the <a href=http://arxiv.org/abs/1107.4339>XC density representation</a>.
318 !%Option none 0
319 !% No correction is applied.
320 !%Option long_range_x 1
321 !% The correction is applied to the exchange potential.
322 !%End
323 call parse_variable(namespace, 'XCDensityCorrection', lr_none, xcs%xc_density_correction)
324
325 if (xcs%xc_density_correction /= lr_none) then
326 call messages_experimental('XC density correction', namespace=namespace)
327
328 if(ispin /= unpolarized) then
329 call messages_not_implemented('XCDensityCorrection with SpinComponents /= unpolarized', namespace)
330 end if
331
332 !%Variable XCDensityCorrectionOptimize
333 !%Type logical
334 !%Default true
335 !%Section Hamiltonian::XC::DensityCorrection
336 !%Description
337 !% When enabled, the density cutoff will be
338 !% optimized to replicate the boundary conditions of the exact
339 !% XC potential. If the variable is set to no, the value of
340 !% the cutoff must be given by the <tt>XCDensityCorrectionCutoff</tt>
341 !% variable.
342 !%End
343 call parse_variable(namespace, 'XCDensityCorrectionOptimize', .true., xcs%xcd_optimize_cutoff)
344
345 !%Variable XCDensityCorrectionCutoff
346 !%Type float
347 !%Default 0.0
348 !%Section Hamiltonian::XC::DensityCorrection
349 !%Description
350 !% The value of the cutoff applied to the XC density.
351 !%End
352 call parse_variable(namespace, 'XCDensityCorrectionCutoff', m_zero, xcs%xcd_ncutoff)
353
354 !%Variable XCDensityCorrectionMinimum
355 !%Type logical
356 !%Default true
357 !%Section Hamiltonian::XC::DensityCorrection
358 !%Description
359 !% When enabled, the cutoff optimization will
360 !% return the first minimum of the <math>q_{xc}</math> function if it does
361 !% not find a value of -1 (<a href=http://arxiv.org/abs/1107.4339>details</a>).
362 !% This is required for atoms or small
363 !% molecules, but may cause numerical problems.
364 !%End
365 call parse_variable(namespace, 'XCDensityCorrectionMinimum', .true., xcs%xcd_minimum)
366
367 !%Variable XCDensityCorrectionNormalize
368 !%Type logical
369 !%Default true
370 !%Section Hamiltonian::XC::DensityCorrection
371 !%Description
372 !% When enabled, the correction will be
373 !% normalized to reproduce the exact boundary conditions of
374 !% the XC potential.
375 !%End
376 call parse_variable(namespace, 'XCDensityCorrectionNormalize', .true., xcs%xcd_normalize)
378 end if
379
380 !%Variable ParallelXC
381 !%Type logical
382 !%Default true
383 !%Section Execution::Parallelization
384 !%Description
385 !% When enabled, additional parallelization
386 !% will be used for the calculation of the XC functional.
387 !%End
388 call messages_obsolete_variable(namespace, 'XCParallel', 'ParallelXC')
389 call parse_variable(namespace, 'ParallelXC', .true., xcs%parallel)
390
391
392 !%Variable HybridCAMParameters
393 !%Type block
394 !%Section Hamiltonian::XC
395 !%Description
396 !% This variable specifies the <math>\alpha, \beta, \omega<math> for CAM-type
397 !% hybrid functionals. Defaults are zero.
398 !%End
399 xcs%cam_ext = m_zero
400
401 if(parse_block(namespace, 'HybridCamParameters', blk) == 0) then
402 do isp = 1, size(xcs%cam_ext)
403 call parse_block_float(blk, 0, isp - 1, xcs%cam_ext(isp))
404 end do
405 call parse_block_end(blk)
406 end if
407
408 if( any(abs(xcs%cam_ext) > m_epsilon) ) then
409 write(message(1), '(1x)')
410 write(message(2), '(a,f8.5)') "Info: Setting external cam_alpha = ", xcs%cam_ext(1)
411 write(message(3), '(a,f8.5)') "Info: Setting external cam_beta = ", xcs%cam_ext(2)
412 write(message(4), '(a,f8.5)') "Info: Setting external cam_omega = ", xcs%cam_ext(3)
413 call messages_info(4, namespace=namespace)
414 end if
415
416 pop_sub(xc_init.parse)
417 end subroutine parse
418
419 end subroutine xc_init
420
421
422 ! ---------------------------------------------------------
423 subroutine xc_end(xcs)
424 type(xc_t), intent(inout) :: xcs
425
426 integer :: isp
427
428 push_sub(xc_end)
429
430 do isp = 1, 2
431 call xc_functional_end(xcs%functional(func_x, isp))
432 call xc_functional_end(xcs%functional(func_c, isp))
433 call xc_functional_end(xcs%kernel(func_x, isp))
434 call xc_functional_end(xcs%kernel(func_c, isp))
435 end do
436 xcs%family = 0
437 xcs%flags = 0
438
439 pop_sub(xc_end)
440 end subroutine xc_end
441
442 ! ---------------------------------------------------------
447 logical pure function xc_is_orbital_dependent(xcs)
448 type(xc_t), intent(in) :: xcs
449
450 xc_is_orbital_dependent = family_is_hybrid(xcs) .or. &
451 in_family(xcs%functional(func_x,1)%family, [xc_family_oep]) .or. &
452 in_family(xcs%family, [xc_family_mgga, xc_family_nc_mgga])
453
454 end function xc_is_orbital_dependent
455
456 ! ---------------------------------------------------------
458 pure logical function family_is_gga(family, only_collinear)
459 integer, intent(in) :: family
460 logical, optional, intent(in) :: only_collinear
461
462 if(optional_default(only_collinear, .false.)) then
463 family_is_gga = in_family(family, [xc_family_gga, xc_family_hyb_gga, &
464 xc_family_mgga, xc_family_hyb_mgga, xc_family_libvdwxc])
465 else
466 family_is_gga = in_family(family, [xc_family_gga, xc_family_hyb_gga, &
467 xc_family_mgga, xc_family_hyb_mgga, xc_family_libvdwxc, xc_family_nc_mgga])
468 end if
469 end function family_is_gga
470
471 !----------------------------------------------------------------------
476 pure logical function family_is_supported(family)
477 integer, intent(in) :: family
478
479 family_is_supported = in_family(family, [xc_family_lda, xc_family_hyb_lda, xc_family_gga, xc_family_hyb_gga, &
480 xc_family_mgga, xc_family_hyb_mgga, xc_family_libvdwxc])
481 end function family_is_supported
482
483 ! ---------------------------------------------------------
485 pure logical function family_is_mgga(family, only_collinear)
486 integer, optional, intent(in) :: family
487 logical, optional, intent(in) :: only_collinear
488
489 if(optional_default(only_collinear, .false.)) then
490 family_is_mgga = in_family(family, [xc_family_mgga, xc_family_hyb_mgga])
491 else
492 family_is_mgga = in_family(family, [xc_family_mgga, xc_family_hyb_mgga, xc_family_nc_mgga])
493 end if
494 end function family_is_mgga
495
496 ! ---------------------------------------------------------
498 logical pure function family_is_mgga_with_exc(xcs)
499 type(xc_t), intent(in) :: xcs
500
501 integer :: ixc
502
504 do ixc = 1, 2
505 if (in_family(xcs%functional(ixc, 1)%family, [xc_family_mgga, xc_family_hyb_mgga, xc_family_nc_mgga]) &
506 .and. xc_functional_is_energy_functional(xcs%functional(ixc, 1))) then
508 end if
509 end do
510 end function family_is_mgga_with_exc
511
513 logical pure function family_is_hybrid(xcs)
514 type(xc_t), intent(in) :: xcs
515
516 integer :: ixc
517
518 family_is_hybrid = .false.
519 do ixc = 1, 2
520 if (in_family(xcs%functional(ixc, 1)%family, [xc_family_hyb_lda, xc_family_hyb_gga, xc_family_hyb_mgga])) then
522 end if
523 end do
524 end function family_is_hybrid
525
526 pure logical function in_family(family, xc_families)
527 integer, intent(in) :: family
528 integer, intent(in) :: xc_families(:)
529
530 in_family = bitand(family, sum(xc_families)) /= 0
531 end function in_family
532
533 ! ---------------------------------------------------------
535 subroutine copy_global_to_local(global, local, n_block, nspin, ip)
536 real(real64), intent(in) :: global(:,:)
537 real(real64), intent(out) :: local(:,:)
538 integer, intent(in) :: n_block
539 integer, intent(in) :: nspin
540 integer, intent(in) :: ip
541
542 integer :: ib, is
543
544 push_sub(copy_global_to_local)
545
546 do is = 1, nspin
547 !$omp parallel do
548 do ib = 1, n_block
549 local(is, ib) = global(ib + ip - 1, is)
550 end do
551 end do
552
553 pop_sub(copy_global_to_local)
554 end subroutine copy_global_to_local
555
556 ! ---------------------------------------------------------
557 subroutine copy_local_to_global(local, global, n_block, spin_channels, ip)
558 real(real64), intent(in) :: local(:,:)
559 real(real64), intent(inout) :: global(:,:)
560 integer, intent(in) :: n_block
561 integer, intent(in) :: spin_channels
562 integer, intent(in) :: ip
563
564 integer :: ib, is
565
567
568 do is = 1, spin_channels
569 !$omp parallel do
570 do ib = 1, n_block
571 global(ib + ip - 1, is) = global(ib + ip - 1, is) + local(is, ib)
572 end do
573 end do
574
576 end subroutine copy_local_to_global
577
578 ! ---------------------------------------------------------
580 subroutine set_hybrid_params(xcs,namespace)
581 type(namespace_t), intent(in) :: namespace
582 type(xc_t), intent(inout) :: xcs
583
584 real(real64) :: parameters(3)
585
586 push_sub(set_hybrid_params)
587
588 parameters = xcs%cam_ext
589 xcs%cam_alpha = parameters(1)
590
591 select case(xcs%functional(func_c,1)%id)
592 case(xc_hyb_gga_xc_pbeh, xc_hyb_lda_xc_lda0) ! original PBE0/LDA0 in libxc
593 ! If the value is negative, we use the default value of PBE0/LDA0
594 if(parameters(1) < m_zero) parameters(1) = 0.25_real64
595
596 call xc_f03_func_set_ext_params(xcs%functional(func_c,1)%conf, parameters)
597 call xc_f03_func_set_ext_params(xcs%functional(func_c,2)%conf, parameters)
598
599 write(message(1), '(a,f6.3,a)') 'Info: Setting mixing parameter (' , parameters(1) ,').'
600 call messages_info(1)
601
602 case(xc_hyb_gga_xc_cam_pbeh, xc_hyb_lda_xc_cam_lda0)
603 xcs%cam_beta = parameters(2)
604 xcs%cam_omega = parameters(3)
605 ! check parameters
606 call xc_f03_hyb_cam_coef(xcs%functional(func_c,1)%conf, xcs%cam_omega, &
607 xcs%cam_alpha, xcs%cam_beta)
608 call xc_f03_hyb_cam_coef(xcs%functional(func_c,2)%conf, xcs%cam_omega, &
609 xcs%cam_alpha, xcs%cam_beta)
610 write(message(1), '(a,f6.3,a)') 'Info: Setting alpha parameter (' , xcs%cam_alpha ,').'
611 write(message(2), '(a,f6.3,a)') 'Info: Setting beta parameter (' , xcs%cam_beta ,').'
612 write(message(3), '(a,f6.3,a)') 'Info: Setting omega parameter (' , xcs%cam_omega ,').'
613 call messages_info(3)
614
615 case default
616 assert(.false.)
617 end select
618
620 end subroutine set_hybrid_params
621
622 ! ---------------------------------------------------------
624 logical function xc_is_not_size_consistent(xcs, namespace)
625 type(xc_t), intent(in) :: xcs
626 type(namespace_t), intent(in) :: namespace
627
628 xc_is_not_size_consistent = xc_functional_is_not_size_consistent(xcs%functional(func_x,1), namespace) &
629 .or. xc_functional_is_not_size_consistent(xcs%functional(func_c,1), namespace)
630 end function xc_is_not_size_consistent
631
632 ! ---------------------------------------------------------
634 logical pure function xc_is_energy_functional(xcs)
635 type(xc_t), intent(in) :: xcs
636 xc_is_energy_functional = xc_functional_is_energy_functional(xcs%functional(func_x,1)) &
637 .or. xc_functional_is_energy_functional(xcs%functional(func_c,1))
638 end function xc_is_energy_functional
639
640#include "xc_vxc_inc.F90"
641#include "xc_fxc_inc.F90"
642#include "xc_kxc_inc.F90"
643
644#include "xc_vxc_nc_inc.F90"
645
646end module xc_oct_m
647
648
649!! Local Variables:
650!! mode: f90
651!! coding: utf-8
652!! End:
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
integer, parameter, public unpolarized
Parameters...
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_one
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
Definition: io.F90:114
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
This module handles spin dimensions of the states and the k-point distribution.
This module defines the unit system, used for input and output.
subroutine, public xc_functional_write_info(functl, iunit, namespace)
Write functional information.
subroutine, public xc_functional_init(functl, namespace, id, ndim, nel, spin_channels)
integer, parameter, public xc_family_nc_mgga
integer, parameter, public xc_oep_x
Exact exchange.
subroutine, public xc_functional_end(functl)
integer, parameter, public func_c
integer, parameter, public func_x
Definition: xc.F90:114
subroutine, public xc_get_fxc(xcs, mesh, namespace, rho, ispin, fxc, zfxc)
Definition: xc.F90:1974
subroutine, public xc_write_info(xcs, iunit, namespace)
Definition: xc.F90:218
subroutine xc_compute_vxc(der, xcs, st, psolver, namespace, space, quantities, ispin, vxc, ex, ec, deltaxc, vtau, ex_density, ec_density, stress_xc)
Definition: xc.F90:979
subroutine, public xc_init(xcs, namespace, ndim, periodic_dim, nel, x_id, c_id, xk_id, ck_id, hartree_fock, ispin)
Definition: xc.F90:248
subroutine, public xc_get_kxc(xcs, mesh, namespace, rho, ispin, kxc)
Definition: xc.F90:2227
pure logical function family_is_supported(family)
Is the xc family internally supported by Octopus.
Definition: xc.F90:570
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:579
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Definition: xc.F90:592
integer, parameter lr_x
Definition: xc.F90:210
subroutine, public xc_end(xcs)
Definition: xc.F90:517
logical pure function, public xc_is_energy_functional(xcs)
Is one of the x or c functional is not an energy functional.
Definition: xc.F90:728
subroutine, public xc_get_vxc(gr, xcs, st, kpoints, psolver, namespace, space, rho, ispin, rcell_volume, vxc, ex, ec, deltaxc, vtau, ex_density, ec_density, stress_xc, force_orbitalfree)
Definition: xc.F90:754
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
Definition: xc.F90:607
logical function, public xc_is_not_size_consistent(xcs, namespace)
Is one of the x or c functional is not size consistent.
Definition: xc.F90:718
subroutine copy_local_to_global(local, global, n_block, spin_channels, ip)
Definition: xc.F90:651
subroutine, public xc_get_nc_vxc(gr, xcs, st, kpoints, space, namespace, rho, vxc, ex, ec, vtau, ex_density, ec_density)
This routines is similar to xc_get_vxc but for noncollinear functionals, which are not implemented in...
Definition: xc.F90:2413
pure logical function family_is_nc_mgga(family)
Returns true is the functional is a noncollinear functional.
Definition: xc.F90:2986
pure logical function family_is_gga(family, only_collinear)
Is the xc function part of the GGA family.
Definition: xc.F90:552
subroutine set_hybrid_params(xcs, namespace)
Sets external parameters for some hybrid functionals.
Definition: xc.F90:674
pure logical function, public in_family(family, xc_families)
Definition: xc.F90:620
subroutine copy_global_to_local(global, local, n_block, nspin, ip)
make a local copy with the correct memory order for libxc
Definition: xc.F90:629
logical pure function, public xc_is_orbital_dependent(xcs)
Is the xc family orbital dependent.
Definition: xc.F90:541
int true(void)
subroutine parse()
Definition: xc.F90:378