Line data Source code
1 : !--------------------------------------------------------------------------------------------------!
2 : ! CP2K: A general program to perform molecular dynamics simulations !
3 : ! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4 : ! !
5 : ! SPDX-License-Identifier: GPL-2.0-or-later !
6 : !--------------------------------------------------------------------------------------------------!
7 :
8 : ! **************************************************************************************************
9 : !> \brief Print basis sets in CP2K format
10 : !> \par History
11 : !> \author JGH (12.2017)
12 : ! **************************************************************************************************
13 : MODULE basis_set_output
14 : USE basis_set_types, ONLY: get_gto_basis_set,&
15 : gto_basis_set_type
16 : USE cp2k_info, ONLY: compile_revision,&
17 : cp2k_version,&
18 : r_datx,&
19 : r_host_name,&
20 : r_user_name
21 : USE cp_files, ONLY: close_file,&
22 : open_file
23 : USE cp_log_handling, ONLY: cp_get_default_logger,&
24 : cp_logger_get_default_io_unit,&
25 : cp_logger_type
26 : USE input_section_types, ONLY: section_vals_type,&
27 : section_vals_val_get
28 : USE kinds, ONLY: default_string_length,&
29 : dp
30 : USE qs_environment_types, ONLY: get_qs_env,&
31 : qs_environment_type
32 : USE qs_kind_types, ONLY: get_qs_kind,&
33 : qs_kind_type
34 : #include "./base/base_uses.f90"
35 :
36 : IMPLICIT NONE
37 : PRIVATE
38 :
39 : ! Global parameters
40 : CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'basis_set_output'
41 : PUBLIC :: print_basis_set_file
42 :
43 : ! **************************************************************************************************
44 :
45 : CONTAINS
46 :
47 : ! **************************************************************************************************
48 : !> \brief ...
49 : !> \param qs_env ...
50 : !> \param base_section ...
51 : ! **************************************************************************************************
52 2 : SUBROUTINE print_basis_set_file(qs_env, base_section)
53 :
54 : TYPE(qs_environment_type), POINTER :: qs_env
55 : TYPE(section_vals_type), POINTER :: base_section
56 :
57 : CHARACTER(LEN=2) :: element_symbol
58 : CHARACTER(LEN=default_string_length) :: bname, filename
59 : INTEGER :: ikind, iunit, nkind, ounit
60 : INTEGER, SAVE :: ncalls = 0
61 : TYPE(cp_logger_type), POINTER :: logger
62 : TYPE(gto_basis_set_type), POINTER :: aux_fit_basis, lri_aux_basis, orb_basis, &
63 : p_lri_aux_basis, ri_aux_basis, ri_hfx_basis, ri_hxc_basis, ri_xas_basis
64 2 : TYPE(qs_kind_type), DIMENSION(:), POINTER :: qs_kind_set
65 : TYPE(qs_kind_type), POINTER :: qs_kind
66 :
67 0 : IF (ncalls > 0) RETURN
68 2 : ncalls = ncalls + 1
69 :
70 2 : logger => cp_get_default_logger()
71 2 : ounit = cp_logger_get_default_io_unit(logger)
72 :
73 2 : CALL section_vals_val_get(base_section, "FILENAME", c_val=filename)
74 :
75 2 : IF (ounit > 0) THEN
76 1 : WRITE (UNIT=ounit, FMT='(/,(T2,A))') REPEAT("-", 79)
77 1 : WRITE (UNIT=ounit, FMT='((T2,A,A))') "Print Basis Set File: ", TRIM(filename)
78 1 : WRITE (UNIT=ounit, FMT='((T2,A))') REPEAT("-", 79)
79 1 : CALL open_file(filename, unit_number=iunit, file_status="UNKNOWN", file_action="WRITE")
80 : WRITE (UNIT=iunit, FMT="(A8,T11,A)") &
81 1 : "# TITLE ", "Basis set file created by "//TRIM(cp2k_version)//" (revision "//TRIM(compile_revision)//")", &
82 2 : "# AUTHOR", TRIM(r_user_name)//"@"//TRIM(r_host_name)//" "//r_datx(1:19)
83 :
84 : END IF
85 :
86 2 : CALL get_qs_env(qs_env, qs_kind_set=qs_kind_set, nkind=nkind)
87 6 : DO ikind = 1, nkind
88 4 : qs_kind => qs_kind_set(ikind)
89 4 : CALL get_qs_kind(qs_kind, element_symbol=element_symbol)
90 4 : NULLIFY (orb_basis, ri_aux_basis, lri_aux_basis, p_lri_aux_basis, aux_fit_basis)
91 4 : CALL get_qs_kind(qs_kind, basis_set=orb_basis, basis_type="ORB")
92 4 : CALL get_qs_kind(qs_kind, basis_set=ri_aux_basis, basis_type="RI_AUX")
93 4 : CALL get_qs_kind(qs_kind, basis_set=ri_hxc_basis, basis_type="RI_HXC")
94 4 : CALL get_qs_kind(qs_kind, basis_set=ri_hfx_basis, basis_type="RI_HFX")
95 4 : CALL get_qs_kind(qs_kind, basis_set=lri_aux_basis, basis_type="LRI_AUX")
96 4 : CALL get_qs_kind(qs_kind, basis_set=p_lri_aux_basis, basis_type="P_LRI_AUX")
97 4 : CALL get_qs_kind(qs_kind, basis_set=aux_fit_basis, basis_type="AUX_FIT")
98 4 : CALL get_qs_kind(qs_kind, basis_set=ri_xas_basis, basis_type="RI_XAS")
99 6 : IF (ounit > 0) THEN
100 2 : IF (ASSOCIATED(orb_basis)) THEN
101 2 : bname = "local_orbital"
102 2 : CALL basis_out(orb_basis, element_symbol, bname, iunit)
103 : END IF
104 2 : IF (ASSOCIATED(ri_aux_basis)) THEN
105 0 : bname = "local_ri_aux"
106 0 : CALL basis_out(ri_aux_basis, element_symbol, bname, iunit)
107 : END IF
108 2 : IF (ASSOCIATED(ri_hxc_basis)) THEN
109 0 : bname = "local_ri_hxc"
110 0 : CALL basis_out(ri_hxc_basis, element_symbol, bname, iunit)
111 : END IF
112 2 : IF (ASSOCIATED(lri_aux_basis)) THEN
113 2 : bname = "local_lri_aux"
114 2 : CALL basis_out(lri_aux_basis, element_symbol, bname, iunit)
115 : END IF
116 2 : IF (ASSOCIATED(p_lri_aux_basis)) THEN
117 2 : bname = "local_p_lri_aux"
118 2 : CALL basis_out(p_lri_aux_basis, element_symbol, bname, iunit)
119 : END IF
120 2 : IF (ASSOCIATED(aux_fit_basis)) THEN
121 0 : bname = "local_aux_fit"
122 0 : CALL basis_out(aux_fit_basis, element_symbol, bname, iunit)
123 : END IF
124 2 : IF (ASSOCIATED(ri_xas_basis)) THEN
125 0 : bname = "local_ri_xas"
126 0 : CALL basis_out(ri_xas_basis, element_symbol, bname, iunit)
127 : END IF
128 2 : IF (ASSOCIATED(ri_hfx_basis)) THEN
129 0 : bname = "local_ri_hfx"
130 0 : CALL basis_out(ri_hfx_basis, element_symbol, bname, iunit)
131 : END IF
132 : END IF
133 : END DO
134 :
135 2 : IF (ounit > 0) THEN
136 1 : CALL close_file(iunit)
137 : END IF
138 :
139 2 : END SUBROUTINE print_basis_set_file
140 :
141 : ! **************************************************************************************************
142 :
143 : ! **************************************************************************************************
144 : !> \brief ...
145 : !> \param basis ...
146 : !> \param element_symbol ...
147 : !> \param bname ...
148 : !> \param iunit ...
149 : ! **************************************************************************************************
150 12 : SUBROUTINE basis_out(basis, element_symbol, bname, iunit)
151 : TYPE(gto_basis_set_type), POINTER :: basis
152 : CHARACTER(LEN=*), INTENT(IN) :: element_symbol, bname
153 : INTEGER, INTENT(IN) :: iunit
154 :
155 : INTEGER :: ipgf, iset, ishell, ll, nset
156 : INTEGER, DIMENSION(0:9) :: lset
157 6 : INTEGER, DIMENSION(:), POINTER :: lmax, lmin, npgf, nshell
158 6 : INTEGER, DIMENSION(:, :), POINTER :: l, n
159 6 : REAL(KIND=dp), DIMENSION(:, :), POINTER :: zet
160 6 : REAL(KIND=dp), DIMENSION(:, :, :), POINTER :: gcc
161 :
162 6 : WRITE (iunit, "(A1)") "#"
163 6 : WRITE (iunit, "(A2,T5,A)") element_symbol, ADJUSTL(TRIM(bname))
164 :
165 : CALL get_gto_basis_set(basis, nset=nset, npgf=npgf, lmax=lmax, lmin=lmin, &
166 : nshell=nshell, n=n, l=l, &
167 6 : gcc=gcc, zet=zet)
168 :
169 6 : WRITE (iunit, "(I5)") nset
170 58 : DO iset = 1, nset
171 52 : lset = 0
172 246 : DO ishell = 1, nshell(iset)
173 194 : ll = l(ishell, iset)
174 246 : lset(ll) = lset(ll) + 1
175 : END DO
176 52 : WRITE (iunit, "(I5,2I3,I5,2X,10(I3))") n(1, iset), lmin(iset), lmax(iset), npgf(iset), &
177 104 : (lset(ll), ll=lmin(iset), lmax(iset))
178 122 : DO ipgf = 1, npgf(iset)
179 116 : WRITE (iunit, "(F20.10,50(F15.10))") zet(ipgf, iset), (gcc(ipgf, ishell, iset), ishell=1, nshell(iset))
180 : END DO
181 : END DO
182 :
183 6 : END SUBROUTINE basis_out
184 :
185 : ! **************************************************************************************************
186 :
187 : END MODULE basis_set_output
|