(git:d18deda)
Loading...
Searching...
No Matches
qs_basis_rotation_methods.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
12 USE cell_types, ONLY: cell_type
15 USE kinds, ONLY: dp
16 USE kpoint_types, ONLY: kpoint_sym_type,&
18 USE orbital_pointers, ONLY: nso
24 USE qs_kind_types, ONLY: get_qs_kind,&
27#include "./base/base_uses.f90"
28
29 IMPLICIT NONE
30
31 PRIVATE
32
33 ! Global parameters (only in this module)
34
35 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qs_basis_rotation_methods'
36
37 ! Public subroutines
38
39 PUBLIC :: qs_basis_rotation
40
41CONTAINS
42
43! **************************************************************************************************
44!> \brief Construct basis set rotation matrices
45!> \param qs_env ...
46!> \param kpoints ...
47! **************************************************************************************************
48 SUBROUTINE qs_basis_rotation(qs_env, kpoints)
49
50 TYPE(qs_environment_type), POINTER :: qs_env
51 TYPE(kpoint_type), POINTER :: kpoints
52
53 INTEGER :: ik, ikind, ir, ira, irot, jr, lval, &
54 nkind, nrot
55 REAL(kind=dp), DIMENSION(3, 3) :: rotmat
56 TYPE(cell_type), POINTER :: cell
57 TYPE(dft_control_type), POINTER :: dft_control
58 TYPE(gto_basis_set_type), POINTER :: orb_basis
59 TYPE(kpoint_sym_type), POINTER :: kpsym
60 TYPE(orbrotmat_type), DIMENSION(:), POINTER :: orbrot
61 TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
62
63 cpassert(ASSOCIATED(qs_env))
64 cpassert(ASSOCIATED(kpoints))
65 IF (ASSOCIATED(kpoints%kind_rotmat)) THEN
66 CALL get_qs_env(qs_env, cell=cell)
67 CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set)
68 CALL get_qs_kind_set(qs_kind_set, maxlgto=lval)
69 nrot = SIZE(kpoints%kind_rotmat, 1)
70 nkind = SIZE(kpoints%kind_rotmat, 2)
71 ! remove possible old rotation matrices
72 DO irot = 1, nrot
73 DO ikind = 1, nkind
74 IF (ASSOCIATED(kpoints%kind_rotmat(irot, ikind)%rmat)) THEN
75 DEALLOCATE (kpoints%kind_rotmat(irot, ikind)%rmat)
76 END IF
77 END DO
78 END DO
79 ! check all rotations needed
80 NULLIFY (orbrot)
81 CALL get_qs_env(qs_env, dft_control=dft_control)
82 DO ik = 1, kpoints%nkp
83 kpsym => kpoints%kp_sym(ik)%kpoint_sym
84 IF (kpsym%apply_symmetry) THEN
85 DO irot = 1, SIZE(kpsym%rotp)
86 ir = kpsym%rotp(irot)
87 ira = 0
88 DO jr = 1, SIZE(kpoints%ibrot)
89 IF (ir == kpoints%ibrot(jr)) ira = jr
90 END DO
91 IF (ira > 0) THEN
92 IF (.NOT. ASSOCIATED(kpoints%kind_rotmat(ira, 1)%rmat)) THEN
93 rotmat(1:3, 1:3) = matmul(cell%h_inv, &
94 matmul(kpsym%rot(:, :, irot), cell%hmat))
95 CALL calculate_rotmat(orbrot, rotmat, lval)
96 IF (dft_control%qs_control%method_id == do_method_dftb) THEN
97 cpabort("ROTMAT")
98 ELSE
99 DO ikind = 1, nkind
100 CALL get_qs_kind(qs_kind_set(ikind), basis_set=orb_basis)
101 NULLIFY (kpoints%kind_rotmat(ira, ikind)%rmat)
102 CALL set_rotmat_basis(kpoints%kind_rotmat(ira, ikind)%rmat, orbrot, orb_basis)
103 END DO
104 END IF
105 END IF
106 END IF
107 END DO
108 END IF
109 END DO
110 CALL release_rotmat(orbrot)
111 END IF
112
113 END SUBROUTINE qs_basis_rotation
114
115! **************************************************************************************************
116!> \brief ...
117!> \param rmat ...
118!> \param orbrot ...
119!> \param basis ...
120! **************************************************************************************************
121 SUBROUTINE set_rotmat_basis(rmat, orbrot, basis)
122 REAL(kind=dp), DIMENSION(:, :), POINTER :: rmat
123 TYPE(orbrotmat_type), DIMENSION(:), POINTER :: orbrot
124 TYPE(gto_basis_set_type), POINTER :: basis
125
126 INTEGER :: fs1, fs2, iset, ishell, l, nset, nsgf
127 INTEGER, DIMENSION(:), POINTER :: nshell
128 INTEGER, DIMENSION(:, :), POINTER :: first_sgf, lshell
129
130 CALL get_gto_basis_set(gto_basis_set=basis, nsgf=nsgf)
131 ALLOCATE (rmat(nsgf, nsgf))
132 rmat = 0.0_dp
133
134 CALL get_gto_basis_set(gto_basis_set=basis, nset=nset, nshell=nshell, l=lshell, &
135 first_sgf=first_sgf)
136 DO iset = 1, nset
137 DO ishell = 1, nshell(iset)
138 l = lshell(ishell, iset)
139 fs1 = first_sgf(ishell, iset)
140 fs2 = fs1 + nso(l) - 1
141 rmat(fs1:fs2, fs1:fs2) = orbrot(l)%mat(1:nso(l), 1:nso(l))
142 END DO
143 END DO
144
145 END SUBROUTINE set_rotmat_basis
146
subroutine, public get_gto_basis_set(gto_basis_set, name, aliases, norm_type, kind_radius, ncgf, nset, nsgf, cgf_symbol, sgf_symbol, norm_cgf, set_radius, lmax, lmin, lx, ly, lz, m, ncgf_set, npgf, nsgf_set, nshell, cphi, pgf_radius, sphi, scon, zet, first_cgf, first_sgf, l, last_cgf, last_sgf, n, gcc, maxco, maxl, maxpgf, maxsgf_set, maxshell, maxso, nco_sum, npgf_sum, nshell_sum, maxder, short_kind_radius, npgf_seg_sum)
...
Handles all functions related to the CELL.
Definition cell_types.F:15
Defines control structures, which contain the parameters and the settings for the DFT-based calculati...
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_method_dftb
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
Types and basic routines needed for a kpoint calculation.
K-points and crystal symmetry routines based on.
Definition kpsym.F:28
Provides Cartesian and spherical orbital pointers and indices.
integer, dimension(:), allocatable, public nso
Calculation of the spherical harmonics and the corresponding orbital transformation matrices.
subroutine, public calculate_rotmat(orbrotmat, rotmat, lval)
Calculate rotation matrices for spherical harmonics up to value l Joseph Ivanic and Klaus Ruedenberg ...
subroutine, public release_rotmat(orbrotmat)
Release rotation matrices.
subroutine, public qs_basis_rotation(qs_env, kpoints)
Construct basis set rotation matrices.
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.
subroutine, public get_qs_kind_set(qs_kind_set, all_potential_present, tnadd_potential_present, gth_potential_present, sgp_potential_present, paw_atom_present, dft_plus_u_atom_present, maxcgf, maxsgf, maxco, maxco_proj, maxgtops, maxlgto, maxlprj, maxnset, maxsgf_set, ncgf, npgf, nset, nsgf, nshell, maxpol, maxlppl, maxlppnl, maxppnl, nelectron, maxder, max_ngrid_rad, max_sph_harm, maxg_iso_not0, lmax_rho0, basis_rcut, basis_type, total_zeff_corr, npgf_seg)
Get attributes of an atomic kind set.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
Keeps symmetry information about a specific k-point.
Contains information about kpoints.
Provides all information about a quickstep kind.