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_cam_oct_m
52 use xc_f03_lib_m
53#ifdef HAVE_LIBXC_FUNCS
54 use xc_f03_funcs_m
55#endif
56 use xc_fbe_oct_m
59 use xc_lrc_oct_m
60
61 implicit none
62
63 private
64 public :: &
65 xc_t, &
67 xc_init, &
68 xc_end, &
81
82 ! A Structure that contains the quantities needed to compute the functionals
84 real(real64), pointer :: rho(:,:) ! A pointer to the full density
85
86 real(real64), allocatable :: dens(:,:) ! Density (in the local frame of the magnetization)
87 real(real64), allocatable :: gdens(:,:,:) ! Gradient of the density
88 real(real64), allocatable :: ldens(:,:) ! Laplacian of the density
89 real(real64), allocatable :: tau(:,:) ! Kinetic energy density
91
92
93 type xc_t
94 private
95 integer, public :: family
96 integer, public :: flags
97 integer, public :: kernel_family
98 type(xc_functional_t), public :: functional(2,2)
100
101 type(xc_functional_t), public :: kernel(2,2)
102 type(xc_cam_t), public :: cam
103 type(xc_lrc_t), public :: lrc
104
105 logical, public :: use_gi_ked
106 integer, public :: xc_density_correction
107 logical, public :: xcd_optimize_cutoff
108 real(real64), public :: xcd_ncutoff
109 logical, public :: xcd_minimum
110 logical, public :: xcd_normalize
111 logical, public :: parallel
112
113 type(internal_quantities_t), public :: quantities
114
115 contains
116 procedure :: compute_exchange => xc_compute_exchange
117
118 end type xc_t
119
121 real(real64), public, parameter :: xc_tiny = 1.0e-12_real64
122
123 integer, public, parameter :: &
124 LR_NONE = 0, &
125 lr_x = 1
126
127contains
128
129 ! ---------------------------------------------------------
130 subroutine xc_write_info(xcs, iunit, namespace)
131 type(xc_t), intent(in) :: xcs
132 integer, optional, intent(in) :: iunit
133 type(namespace_t), optional, intent(in) :: namespace
134
135 integer :: ifunc
136
137 push_sub(xc_write_info)
138
139 write(message(1), '(a)') "Exchange-correlation:"
140 call messages_info(1, iunit=iunit, namespace=namespace)
141
142 do ifunc = func_x, func_c
143 call xc_functional_write_info(xcs%functional(ifunc, 1), iunit, namespace)
144 end do
145
146 if (abs(xcs%cam%alpha + xcs%cam%beta) > m_epsilon) then
147 write(message(1), '(1x)')
148 write(message(2), '(a,f8.5)') "Exact exchange mixing = ", xcs%cam%alpha
149 write(message(3), '(a,f8.5)') "Exact exchange for short-range beta = ", xcs%cam%beta
150 write(message(4), '(a,f8.5)') "Exact exchange range-separate omega = ", xcs%cam%omega
151 call messages_info(4, iunit=iunit, namespace=namespace)
152 end if
153
154
155 pop_sub(xc_write_info)
156 end subroutine xc_write_info
157
158 ! ---------------------------------------------------------
159 subroutine xc_write_fxc_info(xcs, iunit, namespace)
160 type(xc_t), intent(in) :: xcs
161 integer, optional, intent(in) :: iunit
162 type(namespace_t), optional, intent(in) :: namespace
163
164 integer :: ifunc
165
166 push_sub(xc_write_fxc_info)
167
168 write(message(1), '(a)') "Exchange-correlation kernel:"
169 call messages_info(1, iunit=iunit, namespace=namespace)
170
171 do ifunc = func_x, func_c
172 call xc_functional_write_info(xcs%kernel(ifunc, 1), iunit, namespace)
173 end do
174
175 pop_sub(xc_write_fxc_info)
176 end subroutine xc_write_fxc_info
177
179 ! ---------------------------------------------------------
180 subroutine xc_init(xcs, namespace, ndim, periodic_dim, nel, x_id, c_id, xk_id, ck_id, hartree_fock, ispin)
181 type(xc_t), intent(out) :: xcs
182 type(namespace_t), intent(in) :: namespace
183 integer, intent(in) :: ndim
184 integer, intent(in) :: periodic_dim
185 real(real64), intent(in) :: nel
186 integer, intent(in) :: x_id
187 integer, intent(in) :: c_id
188 integer, intent(in) :: xk_id
189 integer, intent(in) :: ck_id
190 logical, intent(in) :: hartree_fock
191 integer, intent(in) :: ispin
193 integer :: isp, xc_major, xc_minor, xc_micro
194 logical :: ll
195 type(block_t) :: blk
196 type(xc_cam_t) :: cam_ext
198 push_sub(xc_init)
199
200 call xc_f03_version(xc_major, xc_minor, xc_micro)
202 xcs%family = 0
203 xcs%flags = 0
204 xcs%kernel_family = 0
206 call parse()
207
208 call xc_lrc_init(xcs%lrc, namespace, ndim, periodic_dim)
209
210 !we also need XC functionals that do not depend on the current
211 !get both spin-polarized and unpolarized
212
213 ! TODO: check that nel should not be spin polarized here
214 do isp = 1, 2
215
216 call xc_functional_init(xcs%functional(func_x, isp), namespace, x_id, ndim, nel, isp)
217 call xc_functional_init(xcs%functional(func_c, isp), namespace, c_id, ndim, nel, isp)
219 call xc_functional_init(xcs%kernel(func_x, isp), namespace, xk_id, ndim, nel, isp)
220 call xc_functional_init(xcs%kernel(func_c, isp), namespace, ck_id, ndim, nel, isp)
221
222 end do
223
224 xcs%family = ior(xcs%family, xcs%functional(func_x,1)%family)
225 xcs%family = ior(xcs%family, xcs%functional(func_c,1)%family)
226
227 xcs%flags = ior(xcs%flags, xcs%functional(func_x,1)%flags)
228 xcs%flags = ior(xcs%flags, xcs%functional(func_c,1)%flags)
229
230 xcs%kernel_family = ior(xcs%kernel_family, xcs%kernel(func_x,1)%family)
231 xcs%kernel_family = ior(xcs%kernel_family, xcs%kernel(func_c,1)%family)
232
233
234 if (xc_is_not_size_consistent(xcs, namespace) .and. periodic_dim > 0) then
235 message(1) = "Cannot perform a periodic calculation with a functional"
236 message(2) = "that depends on the number of electrons."
237 call messages_fatal(2, namespace=namespace)
238 end if
239
240 ! Take care of hybrid functionals (they appear in the correlation functional)
241 xcs%cam = cam_null
242
243 ll = (hartree_fock) &
244 .or.(xcs%functional(func_x,1)%id == xc_oep_x) &
245 .or. family_is_hybrid(xcs)
246 if (ll) then
247 if ((xcs%functional(func_x,1)%id /= 0).and.(xcs%functional(func_x,1)%id /= xc_oep_x)) then
248 message(1) = "You cannot use an exchange functional when performing"
249 message(2) = "a Hartree-Fock calculation or using a hybrid functional."
250 call messages_fatal(2, namespace=namespace)
251 end if
252
253 if (periodic_dim == ndim) then
254 call messages_experimental("Fock operator (Hartree-Fock, OEP, hybrids) in fully periodic systems", namespace=namespace)
255 end if
256
257 ! get the mixing coefficient for hybrids
258 if (family_is_hybrid(xcs)) then
259 if( .not. cam_ext%is_null() ) call set_hybrid_params(xcs, namespace, cam_ext)
260 call xc_f03_hyb_cam_coef(xcs%functional(func_c,1)%conf, xcs%cam%omega, &
261 xcs%cam%alpha, xcs%cam%beta)
262 call xc_f03_hyb_cam_coef(xcs%functional(func_c,2)%conf, xcs%cam%omega, &
263 xcs%cam%alpha, xcs%cam%beta)
264 else
265 ! we are doing Hartree-Fock plus possibly a correlation functional
266 xcs%cam = cam_exact_exchange
267 end if
268
269 ! reset certain variables
270 xcs%functional(func_x,1)%family = xc_family_oep
271 xcs%functional(func_x,1)%id = xc_oep_x
272 xcs%functional(func_x,2)%family = xc_family_oep
273 xcs%functional(func_x,2)%id = xc_oep_x
274 if (.not. hartree_fock) then
275 xcs%family = ior(xcs%family, xc_family_oep)
276 end if
277 end if
278
279 if (in_family(xcs%family, [xc_family_lca])) then
280 call messages_not_implemented("LCA current functionals", namespace) ! not even in libxc!
281 end if
282
283 call messages_obsolete_variable(namespace, 'MGGAimplementation')
284 call messages_obsolete_variable(namespace, 'CurrentInTau', 'XCUseGaugeIndependentKED')
285
286 if (xcs%functional(func_x, 1)%id == xc_mgga_x_tb09 .and. periodic_dim /= 3) then
287 message(1) = "mgga_x_tb09 functional can only be used for 3D periodic systems"
288 call messages_fatal(1, namespace=namespace)
289 end if
290
291 if (family_is_mgga(xcs%family) .or. family_is_nc_mgga(xcs%family)) then
292 !%Variable XCUseGaugeIndependentKED
293 !%Type logical
294 !%Default yes
295 !%Section Hamiltonian::XC
296 !%Description
297 !% If true, when evaluating the XC functional, a term including the (paramagnetic or total) current
298 !% is added to the kinetic-energy density such as to make it gauge-independent.
299 !% Applies only to meta-GGA (and hybrid meta-GGA) functionals.
300 !%End
301 call parse_variable(namespace, 'XCUseGaugeIndependentKED', .true., xcs%use_gi_ked)
302 end if
303
304 pop_sub(xc_init)
305
306 contains
307
308 subroutine parse()
309
310 push_sub(xc_init.parse)
311
312 ! the values of x_id, c_id, xk_id, and c_id are read outside the routine
313
314 !%Variable XCDensityCorrection
315 !%Type integer
316 !%Default none
317 !%Section Hamiltonian::XC::DensityCorrection
318 !%Description
319 !% This variable controls the long-range correction of the XC
320 !% potential using the <a href=http://arxiv.org/abs/1107.4339>XC density representation</a>.
321 !%Option none 0
322 !% No correction is applied.
323 !%Option long_range_x 1
324 !% The correction is applied to the exchange potential.
325 !%End
326 call parse_variable(namespace, 'XCDensityCorrection', lr_none, xcs%xc_density_correction)
327
328 if (xcs%xc_density_correction /= lr_none) then
329 call messages_experimental('XC density correction', namespace=namespace)
330
331 if(ispin /= unpolarized) then
332 call messages_not_implemented('XCDensityCorrection with SpinComponents /= unpolarized', namespace)
333 end if
334
335 !%Variable XCDensityCorrectionOptimize
336 !%Type logical
337 !%Default true
338 !%Section Hamiltonian::XC::DensityCorrection
339 !%Description
340 !% When enabled, the density cutoff will be
341 !% optimized to replicate the boundary conditions of the exact
342 !% XC potential. If the variable is set to no, the value of
343 !% the cutoff must be given by the <tt>XCDensityCorrectionCutoff</tt>
344 !% variable.
345 !%End
346 call parse_variable(namespace, 'XCDensityCorrectionOptimize', .true., xcs%xcd_optimize_cutoff)
347
348 !%Variable XCDensityCorrectionCutoff
349 !%Type float
350 !%Default 0.0
351 !%Section Hamiltonian::XC::DensityCorrection
352 !%Description
353 !% The value of the cutoff applied to the XC density.
354 !%End
355 call parse_variable(namespace, 'XCDensityCorrectionCutoff', m_zero, xcs%xcd_ncutoff)
356
357 !%Variable XCDensityCorrectionMinimum
358 !%Type logical
359 !%Default true
360 !%Section Hamiltonian::XC::DensityCorrection
361 !%Description
362 !% When enabled, the cutoff optimization will
363 !% return the first minimum of the <math>q_{xc}</math> function if it does
364 !% not find a value of -1 (<a href=http://arxiv.org/abs/1107.4339>details</a>).
365 !% This is required for atoms or small
366 !% molecules, but may cause numerical problems.
367 !%End
368 call parse_variable(namespace, 'XCDensityCorrectionMinimum', .true., xcs%xcd_minimum)
369
370 !%Variable XCDensityCorrectionNormalize
371 !%Type logical
372 !%Default true
373 !%Section Hamiltonian::XC::DensityCorrection
374 !%Description
375 !% When enabled, the correction will be
376 !% normalized to reproduce the exact boundary conditions of
377 !% the XC potential.
378 !%End
379 call parse_variable(namespace, 'XCDensityCorrectionNormalize', .true., xcs%xcd_normalize)
380
381 end if
382
383 !%Variable ParallelXC
384 !%Type logical
385 !%Default true
386 !%Section Execution::Parallelization
387 !%Description
388 !% When enabled, additional parallelization
389 !% will be used for the calculation of the XC functional.
390 !%End
391 call messages_obsolete_variable(namespace, 'XCParallel', 'ParallelXC')
392 call parse_variable(namespace, 'ParallelXC', .true., xcs%parallel)
393
394
395 !%Variable HybridCAMParameters
396 !%Type block
397 !%Section Hamiltonian::XC
398 !%Description
399 !% This variable specifies the <math>\alpha, \beta, \omega</math> for CAM-type
400 !% hybrid functionals. Defaults are zero.
401 !%End
402 cam_ext = cam_null
404 if(parse_block(namespace, 'HybridCamParameters', blk) == 0) then
405 call parse_block_float(blk, 0, 0, cam_ext%alpha)
406 call parse_block_float(blk, 0, 1, cam_ext%beta)
407 call parse_block_float(blk, 0, 2, cam_ext%omega)
408 call parse_block_end(blk)
409 end if
410
411 if(.not. cam_ext%is_null()) then
412 call cam_ext%print(namespace, msg="Info: Setting external CAM parameters")
413 endif
414
415 pop_sub(xc_init.parse)
416 end subroutine parse
417
418 end subroutine xc_init
419
420
421 ! ---------------------------------------------------------
422 subroutine xc_end(xcs)
423 type(xc_t), intent(inout) :: xcs
424
425 integer :: isp
426
427 push_sub(xc_end)
428
429 do isp = 1, 2
430 call xc_functional_end(xcs%functional(func_x, isp))
431 call xc_functional_end(xcs%functional(func_c, isp))
432 call xc_functional_end(xcs%kernel(func_x, isp))
433 call xc_functional_end(xcs%kernel(func_c, isp))
434 end do
435 xcs%family = 0
436 xcs%flags = 0
437
438 pop_sub(xc_end)
439 end subroutine xc_end
440
441 ! ---------------------------------------------------------
446 logical pure function xc_is_orbital_dependent(xcs)
447 type(xc_t), intent(in) :: xcs
448
449 xc_is_orbital_dependent = family_is_hybrid(xcs) .or. &
450 in_family(xcs%functional(func_x,1)%family, [xc_family_oep]) .or. &
451 in_family(xcs%family, [xc_family_mgga, xc_family_nc_mgga])
452
453 end function xc_is_orbital_dependent
454
455 ! ---------------------------------------------------------
457 pure logical function family_is_gga(family, only_collinear)
458 integer, intent(in) :: family
459 logical, optional, intent(in) :: only_collinear
460
461 if(optional_default(only_collinear, .false.)) then
462 family_is_gga = in_family(family, [xc_family_gga, xc_family_hyb_gga, &
463 xc_family_mgga, xc_family_hyb_mgga, xc_family_libvdwxc])
464 else
465 family_is_gga = in_family(family, [xc_family_gga, xc_family_hyb_gga, &
466 xc_family_mgga, xc_family_hyb_mgga, xc_family_libvdwxc, xc_family_nc_mgga])
467 end if
468 end function family_is_gga
469
470 !----------------------------------------------------------------------
475 pure logical function family_is_supported(family)
476 integer, intent(in) :: family
477
478 family_is_supported = in_family(family, [xc_family_lda, xc_family_hyb_lda, xc_family_gga, xc_family_hyb_gga, &
479 xc_family_mgga, xc_family_hyb_mgga, xc_family_libvdwxc])
480 end function family_is_supported
481
482 ! ---------------------------------------------------------
484 pure logical function family_is_mgga(family, only_collinear)
485 integer, optional, intent(in) :: family
486 logical, optional, intent(in) :: only_collinear
487
488 if(optional_default(only_collinear, .false.)) then
489 family_is_mgga = in_family(family, [xc_family_mgga, xc_family_hyb_mgga])
490 else
491 family_is_mgga = in_family(family, [xc_family_mgga, xc_family_hyb_mgga, xc_family_nc_mgga])
492 end if
493 end function family_is_mgga
494
495 pure logical function family_is_nc_mgga(family)
496 integer, intent(in) :: family
497
498 family_is_nc_mgga = bitand(family, xc_family_nc_mgga) /= 0
499 end function family_is_nc_mgga
500
501 ! ---------------------------------------------------------
503 logical pure function family_is_mgga_with_exc(xcs)
504 type(xc_t), intent(in) :: xcs
505
506 integer :: ixc
507
509 do ixc = 1, 2
510 if (in_family(xcs%functional(ixc, 1)%family, [xc_family_mgga, xc_family_hyb_mgga, xc_family_nc_mgga]) &
511 .and. xc_functional_is_energy_functional(xcs%functional(ixc, 1))) then
513 end if
514 end do
515 end function family_is_mgga_with_exc
516
518 logical pure function family_is_hybrid(xcs)
519 type(xc_t), intent(in) :: xcs
520
521 integer :: ixc
522
523 family_is_hybrid = .false.
524 do ixc = 1, 2
525 if (in_family(xcs%functional(ixc, 1)%family, [xc_family_hyb_lda, xc_family_hyb_gga, xc_family_hyb_mgga])) then
527 end if
528 end do
529 end function family_is_hybrid
530
531 pure logical function in_family(family, xc_families)
532 integer, intent(in) :: family
533 integer, intent(in) :: xc_families(:)
534
535 in_family = bitand(family, sum(xc_families)) /= 0
536 end function in_family
537
539 logical function xc_compute_exchange(xc, theory_level)
540 class(xc_t), intent(in) :: xc
541 integer, intent(in) :: theory_level
542
543 integer, parameter :: exchange_theory_level(3) = [hartree, hartree_fock, rdmft]
544
545 push_sub(xc_compute_exchange)
546
547 xc_compute_exchange = any(exchange_theory_level == theory_level) .or. &
548 (theory_level == generalized_kohn_sham_dft .and. family_is_hybrid(xc)) .or. &
549 (xc%functional(func_x, 1)%id == xc_oep_x_slater) .or. &
550 (bitand(xc%family, xc_family_oep) /= 0)
551
553
554 end function xc_compute_exchange
555
556 ! ---------------------------------------------------------
558 subroutine set_hybrid_params(xcs, namespace, cam_ext)
559 type(xc_t), intent(inout) :: xcs
560 type(namespace_t), intent(in) :: namespace
561 type(xc_cam_t), intent(in) :: cam_ext
562
563 real(real64), parameter :: default_alpha_pbe0 = 0.25_real64
564 real(real64) :: parameters(3)
565
566 push_sub(set_hybrid_params)
567
568 ! LibXC expects an array of reals, ordered [alpha, beta, omega]
569 parameters = cam_ext%as_array()
570 xcs%cam%alpha = parameters(1)
571
572 select case(xcs%functional(func_c, 1)%id)
573
574 case(xc_hyb_gga_xc_pbeh, xc_hyb_lda_xc_lda0) ! original PBE0/LDA0 in libxc
575 if(parameters(1) < m_zero) parameters(1) = default_alpha_pbe0
576 call xc_f03_func_set_ext_params(xcs%functional(func_c, 1)%conf, parameters)
577 call xc_f03_func_set_ext_params(xcs%functional(func_c, 2)%conf, parameters)
578 write(message(1), '(a,f6.3,a)') 'Info: Setting mixing parameter (' , parameters(1) ,').'
579 call messages_info(1)
580
581 case(xc_hyb_gga_xc_cam_pbeh, xc_hyb_lda_xc_cam_lda0)
582 xcs%cam%beta = parameters(2)
583 xcs%cam%omega = parameters(3)
584 ! check parameters
585 call xc_f03_hyb_cam_coef(xcs%functional(func_c,1)%conf, xcs%cam%omega, &
586 xcs%cam%alpha, xcs%cam%beta)
587 call xc_f03_hyb_cam_coef(xcs%functional(func_c,2)%conf, xcs%cam%omega, &
588 xcs%cam%alpha, xcs%cam%beta)
589 call xcs%cam%print(namespace, msg="Setting CAM parameters:")
591 case default
592 assert(.false.)
593
594 end select
595
596 pop_sub(set_hybrid_params)
597 end subroutine set_hybrid_params
599 ! ---------------------------------------------------------
601 logical function xc_is_not_size_consistent(xcs, namespace)
602 type(xc_t), intent(in) :: xcs
603 type(namespace_t), intent(in) :: namespace
604
605 xc_is_not_size_consistent = xc_functional_is_not_size_consistent(xcs%functional(func_x,1), namespace) &
606 .or. xc_functional_is_not_size_consistent(xcs%functional(func_c,1), namespace)
607 end function xc_is_not_size_consistent
608
609 ! ---------------------------------------------------------
611 logical pure function xc_is_energy_functional(xcs)
612 type(xc_t), intent(in) :: xcs
613 xc_is_energy_functional = xc_functional_is_energy_functional(xcs%functional(func_x,1)) &
614 .or. xc_functional_is_energy_functional(xcs%functional(func_c,1))
615 end function xc_is_energy_functional
616
617end module xc_oct_m
618
619
620!! Local Variables:
621!! mode: f90
622!! coding: utf-8
623!! 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:191
integer, parameter, public hartree_fock
Definition: global.F90:241
real(real64), parameter, public m_epsilon
Definition: global.F90:207
This module implements the underlying real-space grid.
Definition: grid.F90:119
Definition: io.F90:116
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1067
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:999
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:409
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
This module handles spin dimensions of the states and the k-point distribution.
This module defines the unit system, used for input and output.
type(xc_cam_t), parameter, public cam_null
All CAM parameters set to zero.
Definition: xc_cam.F90:152
type(xc_cam_t), parameter, public cam_exact_exchange
Use only Hartree Fock exact exchange.
Definition: xc_cam.F90:155
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
subroutine, public xc_lrc_init(this, namespace, dim, periodic_dim)
Definition: xc_lrc.F90:146
Definition: xc.F90:116
subroutine, public xc_write_info(xcs, iunit, namespace)
Definition: xc.F90:226
subroutine, public xc_init(xcs, namespace, ndim, periodic_dim, nel, x_id, c_id, xk_id, ck_id, hartree_fock, ispin)
Definition: xc.F90:276
pure logical function, public family_is_mgga(family, only_collinear)
Is the xc function part of the mGGA family.
Definition: xc.F90:580
subroutine set_hybrid_params(xcs, namespace, cam_ext)
Sets external parameters for some hybrid functionals.
Definition: xc.F90:654
pure logical function, public family_is_gga(family, only_collinear)
Is the xc function part of the GGA family.
Definition: xc.F90:553
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:599
subroutine, public xc_end(xcs)
Definition: xc.F90:518
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:707
logical pure function, public family_is_hybrid(xcs)
Returns true if the functional is an hybrid functional.
Definition: xc.F90:614
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:697
integer, parameter, public lr_x
Definition: xc.F90:218
pure logical function, public family_is_nc_mgga(family)
Definition: xc.F90:591
logical function xc_compute_exchange(xc, theory_level)
Theory levels and functionals for which exacy exchange is required.
Definition: xc.F90:635
pure logical function, public family_is_supported(family)
Is the xc family internally supported by Octopus.
Definition: xc.F90:571
pure logical function, public in_family(family, xc_families)
Definition: xc.F90:627
subroutine, public xc_write_fxc_info(xcs, iunit, namespace)
Definition: xc.F90:255
logical pure function, public xc_is_orbital_dependent(xcs)
Is the xc family orbital dependent.
Definition: xc.F90:542
Coulomb-attenuating method parameters, used in the partitioning of the Coulomb potential into a short...
Definition: xc_cam.F90:141
int true(void)
subroutine parse()
Definition: xc.F90:404