61#include "./base/base_uses.f90"
66 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'hirshfeld_methods'
86 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(inout) :: charges
90 INTEGER,
INTENT(IN) :: unit_nr
92 CHARACTER(len=2) :: element_symbol
93 INTEGER :: iatom, ikind, natom, nspin
94 REAL(kind=
dp) :: refc, tc1, zeff
96 natom =
SIZE(charges, 1)
97 nspin =
SIZE(charges, 2)
98 WRITE (unit_nr,
'(/,T2,A)')
'!-----------------------------------------------------------------------------!'
99 WRITE (unit=unit_nr, fmt=
"(T28,A)")
"Hirshfeld Charges"
101 WRITE (unit=unit_nr, fmt=
"(/,T3,A,A)") &
102 "#Atom Element Kind ",
" Ref Charge Population Net charge"
104 WRITE (unit=unit_nr, fmt=
"(/,T3,A,A)") &
105 "#Atom Element Kind ",
" Ref Charge Population Spin moment Net charge"
110 element_symbol=element_symbol, kind_number=ikind)
111 refc = hirshfeld_env%charges(iatom)
114 WRITE (unit=unit_nr, fmt=
"(i7,T15,A2,T20,i3,T27,F8.3,T42,F8.3,T72,F8.3)") &
115 iatom, element_symbol, ikind, refc, charges(iatom, 1), zeff - charges(iatom, 1)
117 WRITE (unit=unit_nr, fmt=
"(i7,T15,A2,T20,i3,T27,F8.3,T36,2F8.3,T61,F8.3,T72,F8.3)") &
118 iatom, element_symbol, ikind, refc, charges(iatom, 1), charges(iatom, 2), &
119 charges(iatom, 1) - charges(iatom, 2), zeff - sum(charges(iatom, :))
121 tc1 = tc1 + (zeff - sum(charges(iatom, :)))
123 WRITE (unit=unit_nr, fmt=
"(/,T3,A,T72,F8.3)")
"Total Charge ", tc1
124 WRITE (unit_nr,
'(T2,A)')
'!-----------------------------------------------------------------------------!'
136 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(inout) :: charges
138 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
141 CHARACTER(LEN=default_string_length) :: description
142 INTEGER :: iatom, ikind, natom
143 REAL(kind=
dp) :: zeff
144 REAL(kind=
dp),
DIMENSION(:),
POINTER :: charges_save
150 natom =
SIZE(charges, 1)
151 ALLOCATE (charges_save(natom))
157 charges_save(iatom) = zeff - sum(charges(iatom, :))
161 description =
"[HIRSHFELD-CHARGES]"
163 CALL put_results(results=results, description=description, &
166 DEALLOCATE (charges_save)
180 TYPE(
qs_kind_type),
DIMENSION(:),
POINTER :: qs_kind_set
182 REAL(kind=
dp),
OPTIONAL :: radius
183 REAL(kind=
dp),
DIMENSION(:),
OPTIONAL,
POINTER :: radii_list
185 INTEGER,
PARAMETER :: ngto = 8
187 CHARACTER(len=2) :: esym
188 INTEGER :: ikind, nkind
190 REAL(kind=
dp) :: al, rco, zeff
191 REAL(kind=
dp),
DIMENSION(ngto, 2) :: ppdens
195 cpassert(
ASSOCIATED(hirshfeld_env))
197 nkind =
SIZE(qs_kind_set)
198 ALLOCATE (hirshfeld_env%kind_shape_fn(nkind))
200 SELECT CASE (hirshfeld_env%shape_function_type)
203 hirshfeld_env%kind_shape_fn(ikind)%numexp = 1
204 ALLOCATE (hirshfeld_env%kind_shape_fn(ikind)%zet(1))
205 ALLOCATE (hirshfeld_env%kind_shape_fn(ikind)%coef(1))
206 CALL get_qs_kind(qs_kind_set(ikind), element_symbol=esym)
208 SELECT CASE (hirshfeld_env%radius_type)
211 rco = max(rco, 1.0_dp)
213 cpassert(
PRESENT(radii_list))
214 cpassert(
ASSOCIATED(radii_list))
215 cpassert(
SIZE(radii_list) == nkind)
217 rco = radii_list(ikind)
220 IF (.NOT. found)
THEN
221 rco = max(rco, 1.0_dp)
223 IF (hirshfeld_env%use_bohr) &
228 IF (.NOT. found)
THEN
229 rco = max(rco, 1.0_dp)
231 IF (hirshfeld_env%use_bohr) &
235 cpassert(
PRESENT(radius))
239 hirshfeld_env%kind_shape_fn(ikind)%zet(1) = al
240 hirshfeld_env%kind_shape_fn(ikind)%coef(1) = (al/
pi)**1.5_dp
245 atomic_kind => atomic_kind_set(ikind)
246 qs_kind => qs_kind_set(ikind)
249 hirshfeld_env%kind_shape_fn(ikind)%numexp = ngto
250 ALLOCATE (hirshfeld_env%kind_shape_fn(ikind)%zet(ngto))
251 ALLOCATE (hirshfeld_env%kind_shape_fn(ikind)%coef(ngto))
252 hirshfeld_env%kind_shape_fn(ikind)%zet(:) = ppdens(:, 1)
254 hirshfeld_env%kind_shape_fn(ikind)%coef(:) = ppdens(:, 2)/zeff
258 cpabort(
"Unknown shape function")
272 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(inout) :: charges
275 LOGICAL :: rho_r_valid
276 REAL(kind=
dp) :: tnfun
285 CALL calculate_hirshfeld_normalization(qs_env, hirshfeld_env)
288 tnfun = abs(tnfun - sum(hirshfeld_env%charges))
290 CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, rho=rho)
291 CALL qs_rho_get(rho, rho_r=rho_r, rho_r_valid=rho_r_valid)
292 CALL pw_env_get(pw_env=pw_env, auxbas_pw_pool=auxbas_pw_pool)
293 CALL auxbas_pw_pool%create_pw(rhonorm)
295 DO is = 1,
SIZE(rho_r)
296 IF (rho_r_valid)
THEN
297 CALL hfun_scale(rhonorm%array, rho_r(is)%array, &
298 hirshfeld_env%fnorm%array)
300 cpabort(
"We need rho in real space")
302 CALL hirshfeld_integration(qs_env, hirshfeld_env, rhonorm, charges(:, is))
303 charges(:, is) = charges(:, is)*hirshfeld_env%charges(:)
305 CALL auxbas_pw_pool%give_back_pw(rhonorm)
314 SUBROUTINE hfun_scale(fout, fun1, fun2)
315 REAL(kind=
dp),
DIMENSION(:, :, :),
INTENT(OUT) :: fout
316 REAL(kind=
dp),
DIMENSION(:, :, :),
INTENT(IN) :: fun1, fun2
318 REAL(kind=
dp),
PARAMETER :: small = 1.0e-12_dp
320 INTEGER :: i1, i2, i3, n1, n2, n3
325 cpassert(n1 ==
SIZE(fun1, 1))
326 cpassert(n2 ==
SIZE(fun1, 2))
327 cpassert(n3 ==
SIZE(fun1, 3))
328 cpassert(n1 ==
SIZE(fun2, 1))
329 cpassert(n2 ==
SIZE(fun2, 2))
330 cpassert(n3 ==
SIZE(fun2, 3))
335 IF (fun2(i1, i2, i3) > small)
THEN
336 fout(i1, i2, i3) = fun1(i1, i2, i3)/fun2(i1, i2, i3)
338 fout(i1, i2, i3) = 0.0_dp
344 END SUBROUTINE hfun_scale
356 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(inout) :: charges
357 INTEGER,
INTENT(IN) :: ounit
359 INTEGER,
PARAMETER :: maxloop = 100
360 REAL(kind=
dp),
PARAMETER :: maxres = 1.0e-2_dp
362 CHARACTER(len=3) :: yesno
363 INTEGER :: iat, iloop, is, natom
364 LOGICAL :: rho_r_valid
365 REAL(kind=
dp) :: res, tnfun
374 natom =
SIZE(charges, 1)
376 IF (ounit > 0)
WRITE (ounit,
"(/,T2,A)")
"Hirshfeld charge iterations: Residuals ..."
378 CALL get_qs_env(qs_env=qs_env, pw_env=pw_env, rho=rho)
379 CALL qs_rho_get(rho, rho_r=rho_r, rho_r_valid=rho_r_valid)
380 CALL pw_env_get(pw_env=pw_env, auxbas_pw_pool=auxbas_pw_pool)
381 CALL auxbas_pw_pool%create_pw(rhonorm)
383 DO iloop = 1, maxloop
386 CALL calculate_hirshfeld_normalization(qs_env, hirshfeld_env)
389 tnfun = abs(tnfun - sum(hirshfeld_env%charges))
391 DO is = 1,
SIZE(rho_r)
392 IF (rho_r_valid)
THEN
393 CALL hfun_scale(rhonorm%array, rho_r(is)%array, &
394 hirshfeld_env%fnorm%array)
396 cpabort(
"We need rho in real space")
398 CALL hirshfeld_integration(qs_env, hirshfeld_env, rhonorm, charges(:, is))
399 charges(:, is) = charges(:, is)*hirshfeld_env%charges(:)
404 res = res + (sum(charges(iat, :)) - hirshfeld_env%charges(iat))**2
406 res = sqrt(res/real(natom, kind=
dp))
409 IF (mod(iloop, 10) == 0) yesno =
"YES"
410 WRITE (ounit, fmt=
"(F8.3)", advance=yesno) res
414 hirshfeld_env%charges(iat) = sum(charges(iat, :))
416 IF (res < maxres)
EXIT
420 CALL auxbas_pw_pool%give_back_pw(rhonorm)
429 SUBROUTINE calculate_hirshfeld_normalization(qs_env, hirshfeld_env)
434 CHARACTER(len=*),
PARAMETER :: routinen =
'calculate_hirshfeld_normalization'
436 INTEGER :: atom_a, handle, iatom, iex, ikind, &
437 ithread, j, natom, npme, nthread, &
438 numexp, subpatch_pattern
439 INTEGER,
DIMENSION(:),
POINTER :: atom_list, cores
440 REAL(kind=
dp) :: alpha, coef, eps_rho_rspace, radius
441 REAL(kind=
dp),
DIMENSION(3) :: ra
442 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: pab
453 CALL timeset(routinen, handle)
455 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set, cell=cell, &
456 dft_control=dft_control, particle_set=particle_set, pw_env=pw_env)
457 CALL pw_env_get(pw_env, auxbas_rs_desc=auxbas_rs_desc, auxbas_rs_grid=rs_rho, &
458 auxbas_pw_pool=auxbas_pw_pool)
462 eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
467 DO ikind = 1,
SIZE(atomic_kind_set)
468 numexp = hirshfeld_env%kind_shape_fn(ikind)%numexp
469 IF (numexp <= 0) cycle
470 CALL get_atomic_kind(atomic_kind_set(ikind), natom=natom, atom_list=atom_list)
471 ALLOCATE (cores(natom))
474 alpha = hirshfeld_env%kind_shape_fn(ikind)%zet(iex)
475 coef = hirshfeld_env%kind_shape_fn(ikind)%coef(iex)
479 atom_a = atom_list(iatom)
480 ra(:) =
pbc(particle_set(atom_a)%r, cell)
481 IF (rs_rho%desc%parallel .AND. .NOT. rs_rho%desc%distributed)
THEN
483 IF (
modulo(iatom, rs_rho%desc%group_size) == rs_rho%desc%my_pos)
THEN
494 atom_a = atom_list(iatom)
495 pab(1, 1) = hirshfeld_env%charges(atom_a)*coef
496 ra(:) =
pbc(particle_set(atom_a)%r, cell)
499 ra=ra, rb=ra, rp=ra, zetp=alpha, eps=eps_rho_rspace, &
500 pab=pab, o1=0, o2=0, &
501 prefactor=1.0_dp, cutoff=0.0_dp)
505 (/0.0_dp, 0.0_dp, 0.0_dp/), 1.0_dp, pab, 0, 0, rs_rho, &
507 use_subpatch=.true., subpatch_pattern=subpatch_pattern)
517 IF (
ASSOCIATED(fnorm))
THEN
522 CALL auxbas_pw_pool%create_pw(fnorm)
527 CALL timestop(handle)
529 END SUBROUTINE calculate_hirshfeld_normalization
539 SUBROUTINE hirshfeld_integration(qs_env, hirshfeld_env, rfun, fval, fderiv)
544 REAL(kind=
dp),
DIMENSION(:),
INTENT(inout) :: fval
545 REAL(kind=
dp),
DIMENSION(:, :),
INTENT(inout), &
548 CHARACTER(len=*),
PARAMETER :: routinen =
'hirshfeld_integration'
550 INTEGER :: atom_a, handle, iatom, iex, ikind, &
551 ithread, j, natom, npme, nthread, &
553 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: cores
554 INTEGER,
DIMENSION(:),
POINTER :: atom_list
556 REAL(kind=
dp) :: alpha, coef, dvol, eps_rho_rspace, radius
557 REAL(kind=
dp),
DIMENSION(3) :: force_a, force_b, ra
558 REAL(kind=
dp),
DIMENSION(:, :),
POINTER :: hab, pab
568 CALL timeset(routinen, handle)
570 do_force =
PRESENT(fderiv)
572 dvol = rfun%pw_grid%dvol
574 NULLIFY (pw_env, auxbas_rs_desc)
576 CALL pw_env_get(pw_env=pw_env, auxbas_rs_desc=auxbas_rs_desc, &
580 CALL get_qs_env(qs_env=qs_env, atomic_kind_set=atomic_kind_set, cell=cell, &
581 dft_control=dft_control, particle_set=particle_set)
582 eps_rho_rspace = dft_control%qs_control%eps_rho_rspace
586 ALLOCATE (hab(1, 1), pab(1, 1))
588 DO ikind = 1,
SIZE(atomic_kind_set)
589 numexp = hirshfeld_env%kind_shape_fn(ikind)%numexp
590 IF (numexp <= 0) cycle
591 CALL get_atomic_kind(atomic_kind_set(ikind), natom=natom, atom_list=atom_list)
592 ALLOCATE (cores(natom))
595 alpha = hirshfeld_env%kind_shape_fn(ikind)%zet(iex)
596 coef = hirshfeld_env%kind_shape_fn(ikind)%coef(iex)
600 atom_a = atom_list(iatom)
601 ra(:) =
pbc(particle_set(atom_a)%r, cell)
602 IF (rs_v%desc%parallel .AND. .NOT. rs_v%desc%distributed)
THEN
604 IF (
modulo(iatom, rs_v%desc%group_size) == rs_v%desc%my_pos)
THEN
616 atom_a = atom_list(iatom)
617 ra(:) =
pbc(particle_set(atom_a)%r, cell)
624 ra=ra, rb=ra, rp=ra, &
625 zetp=alpha, eps=eps_rho_rspace, &
626 pab=pab, o1=0, o2=0, &
627 prefactor=1.0_dp, cutoff=1.0_dp)
630 0, 0.0_dp, 0, ra, (/0.0_dp, 0.0_dp, 0.0_dp/), &
631 rs_v, hab, pab=pab, o1=0, o2=0, &
632 radius=radius, calculate_forces=do_force, &
633 force_a=force_a, force_b=force_b, use_virial=.false., &
634 use_subpatch=.true., subpatch_pattern=0)
635 fval(atom_a) = fval(atom_a) + hab(1, 1)*dvol*coef
637 fderiv(:, atom_a) = fderiv(:, atom_a) + force_a(:)*dvol
646 DEALLOCATE (hab, pab)
648 CALL get_qs_env(qs_env=qs_env, para_env=para_env)
649 CALL para_env%sum(fval)
651 CALL timestop(handle)
653 END SUBROUTINE hirshfeld_integration
static GRID_HOST_DEVICE int modulo(int a, int m)
Equivalent of Fortran's MODULO, which always return a positive number. https://gcc....
All kind of helpful little routines.
real(kind=dp) function, public exp_radius_very_extended(la_min, la_max, lb_min, lb_max, pab, o1, o2, ra, rb, rp, zetp, eps, prefactor, cutoff, epsabs)
computes the radius of the Gaussian outside of which it is smaller than eps
calculate the orbitals for a given atomic kind type
subroutine, public calculate_atomic_density(density, atomic_kind, qs_kind, ngto, iunit, optbasis, allelectron, confine)
...
Define the atomic kind types and their sub types.
subroutine, public get_atomic_kind(atomic_kind, fist_potential, element_symbol, name, mass, kind_number, natom, atom_list, rcov, rvdw, z, qeff, apol, cpol, mm_radius, shell, shell_active, damping)
Get attributes of an atomic kind.
Handles all functions related to the CELL.
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
set of type/routines to handle the storage of results in force_envs
subroutine, public cp_results_erase(results, description, nval)
erase a part of result_list
set of type/routines to handle the storage of results in force_envs
real(kind=dp) function, public cp_unit_to_cp2k(value, unit_str, defaults, power)
converts to the internal cp2k units to the given unit
Fortran API for the grid package, which is written in C.
integer, parameter, public grid_func_ab
subroutine, public integrate_pgf_product(la_max, zeta, la_min, lb_max, zetb, lb_min, ra, rab, rsgrid, hab, pab, o1, o2, radius, calculate_forces, force_a, force_b, compute_tau, use_virial, my_virial_a, my_virial_b, hdab, hadb, a_hdab, use_subpatch, subpatch_pattern)
low level function to compute matrix elements of primitive gaussian functions
subroutine, public collocate_pgf_product(la_max, zeta, la_min, lb_max, zetb, lb_min, ra, rab, scale, pab, o1, o2, rsgrid, ga_gb_function, radius, use_subpatch, subpatch_pattern)
low level collocation of primitive gaussian functions
Calculate Hirshfeld charges and related functions.
subroutine, public comp_hirshfeld_charges(qs_env, hirshfeld_env, charges)
...
subroutine, public create_shape_function(hirshfeld_env, qs_kind_set, atomic_kind_set, radius, radii_list)
creates kind specific shape functions for Hirshfeld charges
subroutine, public write_hirshfeld_charges(charges, hirshfeld_env, particle_set, qs_kind_set, unit_nr)
...
subroutine, public comp_hirshfeld_i_charges(qs_env, hirshfeld_env, charges, ounit)
...
subroutine, public save_hirshfeld_charges(charges, particle_set, qs_kind_set, qs_env)
saves the Hirshfeld charges to the results structure
The types needed for the calculation of Hirshfeld charges and related functions.
subroutine, public get_hirshfeld_info(hirshfeld_env, shape_function_type, iterative, ref_charge, fnorm, radius_type, use_bohr)
Get information from a Hirshfeld env.
subroutine, public set_hirshfeld_info(hirshfeld_env, shape_function_type, iterative, ref_charge, fnorm, radius_type, use_bohr)
Set values of a Hirshfeld env.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
Interface to the message passing library MPI.
Define the data structure for the particle information.
Periodic Table related data definitions.
subroutine, public get_ptable_info(symbol, number, amass, ielement, covalent_radius, metallic_radius, vdw_radius, found)
Pass information about the kind given the element symbol.
container for various plainwaves related things
subroutine, public pw_env_get(pw_env, pw_pools, cube_info, gridlevel_info, auxbas_pw_pool, auxbas_grid, auxbas_rs_desc, auxbas_rs_grid, rs_descs, rs_grids, xc_pw_pool, vdw_pw_pool, poisson_env, interp_section)
returns the various attributes of the pw env
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
subroutine, public get_qs_env(qs_env, atomic_kind_set, qs_kind_set, cell, super_cell, cell_ref, use_ref_cell, kpoints, dft_control, mos, sab_orb, sab_all, qmmm, qmmm_periodic, sac_ae, sac_ppl, sac_lri, sap_ppnl, sab_vdw, sab_scp, sap_oce, sab_lrc, sab_se, sab_xtbe, sab_tbe, sab_core, sab_xb, sab_xtb_pp, sab_xtb_nonbond, sab_almo, sab_kp, sab_kp_nosym, particle_set, energy, force, matrix_h, matrix_h_im, matrix_ks, matrix_ks_im, matrix_vxc, run_rtp, rtp, matrix_h_kp, matrix_h_im_kp, matrix_ks_kp, matrix_ks_im_kp, matrix_vxc_kp, kinetic_kp, matrix_s_kp, matrix_w_kp, matrix_s_ri_aux_kp, matrix_s, matrix_s_ri_aux, matrix_w, matrix_p_mp2, matrix_p_mp2_admm, rho, rho_xc, pw_env, ewald_env, ewald_pw, active_space, mpools, input, para_env, blacs_env, scf_control, rel_control, kinetic, qs_charges, vppl, rho_core, rho_nlcc, rho_nlcc_g, ks_env, ks_qmmm_env, wf_history, scf_env, local_particles, local_molecules, distribution_2d, dbcsr_dist, molecule_kind_set, molecule_set, subsys, cp_subsys, oce, local_rho_set, rho_atom_set, task_list, task_list_soft, rho0_atom_set, rho0_mpole, rhoz_set, ecoul_1c, rho0_s_rs, rho0_s_gs, do_kpoints, has_unit_metric, requires_mo_derivs, mo_derivs, mo_loc_history, nkind, natom, nelectron_total, nelectron_spin, efield, neighbor_list_id, linres_control, xas_env, virial, cp_ddapc_env, cp_ddapc_ewald, outer_scf_history, outer_scf_ihistory, x_data, et_coupling, dftb_potential, results, se_taper, se_store_int_env, se_nddo_mpole, se_nonbond_env, admm_env, lri_env, lri_density, exstate_env, ec_env, harris_env, dispersion_env, gcp_env, vee, rho_external, external_vxc, mask, mp2_env, bs_env, kg_env, wanniercentres, atprop, ls_scf_env, do_transport, transport_env, v_hartree_rspace, s_mstruct_changed, rho_changed, potential_changed, forces_up_to_date, mscfg_env, almo_scf_env, gradient_history, variable_history, embed_pot, spin_embed_pot, polar_env, mos_last_converged, eeq, rhs)
Get the QUICKSTEP environment.
Define the quickstep kind type and their sub types.
subroutine, public get_qs_kind(qs_kind, basis_set, basis_type, ncgf, nsgf, all_potential, tnadd_potential, gth_potential, sgp_potential, upf_potential, se_parameter, dftb_parameter, xtb_parameter, dftb3_param, zatom, zeff, elec_conf, mao, lmax_dftb, alpha_core_charge, ccore_charge, core_charge, core_charge_radius, paw_proj_set, paw_atom, hard_radius, hard0_radius, max_rad_local, covalent_radius, vdw_radius, gpw_type_forced, harmonics, max_iso_not0, max_s_harm, grid_atom, ngrid_ang, ngrid_rad, lmax_rho0, dft_plus_u_atom, l_of_dft_plus_u, n_of_dft_plus_u, u_minus_j, u_of_dft_plus_u, j_of_dft_plus_u, alpha_of_dft_plus_u, beta_of_dft_plus_u, j0_of_dft_plus_u, occupation_of_dft_plus_u, dispersion, bs_occupation, magnetization, no_optimize, addel, laddel, naddel, orbitals, max_scf, eps_scf, smear, u_ramping, u_minus_j_target, eps_u_ramping, init_u_ramping_each_scf, reltmat, ghost, floating, name, element_symbol, pao_basis_size, pao_model_file, pao_potentials, pao_descriptors, nelec)
Get attributes of an atomic kind.
superstucture that hold various representations of the density and keeps track of which ones are vali...
subroutine, public qs_rho_get(rho_struct, rho_ao, rho_ao_im, rho_ao_kp, rho_ao_im_kp, rho_r, drho_r, rho_g, drho_g, tau_r, tau_g, rho_r_valid, drho_r_valid, rho_g_valid, drho_g_valid, tau_r_valid, tau_g_valid, tot_rho_r, tot_rho_g, rho_r_sccs, soft_valid, complex_rho_ao)
returns info about the density described by this object. If some representation is not available an e...
subroutine, public transfer_pw2rs(rs, pw)
...
subroutine, public transfer_rs2pw(rs, pw)
...
subroutine, public rs_grid_zero(rs)
Initialize grid to zero.
Provides all information about an atomic kind.
Type defining parameters related to the simulation cell.
contains arbitrary information which need to be stored
quantities needed for a Hirshfeld based partitioning of real space
stores all the informations relevant to an mpi environment
contained for different pw related things
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Provides all information about a quickstep kind.
keeps the density in various representations, keeping track of which ones are valid.