Octopus
pseudo.F90
Go to the documentation of this file.
1!! Copyright (C) 2018 X. Andrade
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the Lesser GNU General Public License as published by
5!! the Free Software Foundation; either version 3, 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 pseudo_oct_m
22 use, intrinsic :: iso_fortran_env
23
24 implicit none
25
26 private
27
28 public :: &
29 pseudo_t, &
32 pseudo_end, &
53 pseudo_dij, &
63
64 !the following sets of values have to match with those on base.hpp
65 integer, parameter, public :: &
66 PSEUDO_TYPE_ULTRASOFT = 30, &
70
71 integer, parameter, public :: &
72 PSEUDO_STATUS_SUCCESS = 0, &
79
80 integer, parameter, public :: &
81 PSEUDO_FORMAT_FILE_NOT_FOUND = 773, &
83 pseudo_format_upf1 = 775, &
84 pseudo_format_upf2 = 776, &
85 pseudo_format_qso = 777, &
86 pseudo_format_psml = 778, &
87 pseudo_format_psf = 779, &
88 pseudo_format_cpi = 780, &
89 pseudo_format_fhi = 781, &
90 pseudo_format_hgh = 782, &
92
93 ! we only define these values here, the specific functionals are
94 ! obtained from libxc
95 integer, parameter, public :: &
96 PSEUDO_EXCHANGE_UNKNOWN = -2, &
98
99 integer, parameter, public :: &
100 PSEUDO_CORRELATION_UNKNOWN = -2, &
102
103 type pseudo_t
104 private
105 integer(int64) :: dummy
106 end type pseudo_t
107
108 interface
109
110 ! -------------------------------------------------
111
112 integer function pseudo_detect_format(filename)
113 import :: pseudo_t
114 implicit none
115
116 character(len=*), intent(in) :: filename
117 end function pseudo_detect_format
118
119 ! -------------------------------------------------
120
121 subroutine pseudo_init(pseudo, filename, fmt, ierr)
122 import :: pseudo_t
123 implicit none
124
125 type(pseudo_t), intent(out) :: pseudo
126 character(len=*), intent(in) :: filename
127 integer, intent(in) :: fmt
128 integer, intent(out) :: ierr
129 end subroutine pseudo_init
130
131 ! -------------------------------------------------
132
133 subroutine pseudo_end(pseudo)
134 import :: pseudo_t
135 implicit none
136
137 type(pseudo_t), intent(inout) :: pseudo
138 end subroutine pseudo_end
139
140 ! -------------------------------------------------
141
142 integer function pseudo_type(pseudo)
143 import :: pseudo_t
144 implicit none
145
146 type(pseudo_t), intent(in) :: pseudo
147 end function pseudo_type
148
149 ! -------------------------------------------------
150
151 integer function pseudo_format(pseudo)
152 import :: pseudo_t
153 implicit none
154
155 type(pseudo_t), intent(in) :: pseudo
156 end function pseudo_format
157
158 ! -------------------------------------------------
159
160 real(real64) function pseudo_valence_charge(pseudo)
161 import :: pseudo_t
162 import :: real64
163 implicit none
165 type(pseudo_t), intent(in) :: pseudo
166 end function pseudo_valence_charge
167
168 ! -------------------------------------------------
169
170 real(real64) function pseudo_mesh_spacing(pseudo)
171 import :: pseudo_t
172 import :: real64
173 implicit none
174
175 type(pseudo_t), intent(in) :: pseudo
176 end function pseudo_mesh_spacing
177
178 ! -------------------------------------------------
179
180 integer function pseudo_mesh_size(pseudo)
181 import :: pseudo_t
182 implicit none
183
184 type(pseudo_t), intent(in) :: pseudo
185 end function pseudo_mesh_size
186
187 ! -------------------------------------------------
189 real(real64) function pseudo_mass(pseudo)
190 import :: pseudo_t
191 import :: real64
192 implicit none
193
194 type(pseudo_t), intent(in) :: pseudo
195 end function pseudo_mass
197 ! -------------------------------------------------
199 integer function pseudo_lmax(pseudo)
200 import :: pseudo_t
201 implicit none
202
203 type(pseudo_t), intent(in) :: pseudo
204 end function pseudo_lmax
206 ! -------------------------------------------------
207
208 integer function pseudo_llocal(pseudo)
209 import :: pseudo_t
210 implicit none
211
212 type(pseudo_t), intent(in) :: pseudo
213 end function pseudo_llocal
215 ! -------------------------------------------------
216
217 integer function pseudo_nchannels(pseudo)
218 import :: pseudo_t
219 implicit none
220
221 type(pseudo_t), intent(in) :: pseudo
222 end function pseudo_nchannels
223
224 ! -------------------------------------------------
225
226 integer function pseudo_nprojectors(pseudo)
227 import :: pseudo_t
228 implicit none
229
230 type(pseudo_t), intent(in) :: pseudo
231 end function pseudo_nprojectors
232
233 ! -------------------------------------------------
234
235 integer function pseudo_nprojectors_per_l(pseudo, l)
236 import :: pseudo_t
237 implicit none
238
239 type(pseudo_t), intent(in) :: pseudo
240 integer, intent(in) :: l
241 end function pseudo_nprojectors_per_l
242
243 ! -------------------------------------------------
245 subroutine pseudo_grid(pseudo, grid)
246 import :: pseudo_t
247 import :: real64
248 implicit none
249
250 type(pseudo_t), intent(in) :: pseudo
251 real(real64), intent(out) :: grid
252 end subroutine pseudo_grid
254 ! -------------------------------------------------
255
256 subroutine pseudo_grid_weights(pseudo, weight)
257 import :: pseudo_t
258 import :: real64
259 implicit none
260
261 type(pseudo_t), intent(in) :: pseudo
262 real(real64), intent(out) :: weight
263 end subroutine pseudo_grid_weights
264
265 ! -------------------------------------------------
266
267 subroutine pseudo_local_potential(pseudo, local_potential)
268 import :: pseudo_t
269 import :: real64
270 implicit none
271
272 type(pseudo_t), intent(in) :: pseudo
273 real(real64), intent(out) :: local_potential
274 end subroutine pseudo_local_potential
275
276 ! -------------------------------------------------
277
278 subroutine pseudo_projector(pseudo, l, ic, projector)
279 import :: pseudo_t
280 import :: real64
281 implicit none
283 type(pseudo_t), intent(in) :: pseudo
284 integer, intent(in) :: l
285 integer, intent(in) :: ic
286 real(real64), intent(out) :: projector
287 end subroutine pseudo_projector
288
289 ! -------------------------------------------------
290
291 real(real64) function pseudo_dij(pseudo, l, ic, jc)
292 import :: pseudo_t
293 import :: real64
294 implicit none
295
296 type(pseudo_t), intent(in) :: pseudo
297 integer, intent(in) :: l
298 integer, intent(in) :: ic
299 integer, intent(in) :: jc
300 end function pseudo_dij
302 ! -------------------------------------------------
303
304 subroutine pseudo_radial_potential(pseudo, l, radial_potential)
305 import :: pseudo_t
306 import :: real64
307 implicit none
308
309 type(pseudo_t), intent(in) :: pseudo
310 integer, intent(in) :: l
311 real(real64), intent(inout) :: radial_potential
312 end subroutine pseudo_radial_potential
313
314 ! -------------------------------------------------
315
316 subroutine pseudo_radial_function(pseudo, l, radial_function)
317 import :: pseudo_t
318 import :: real64
319 implicit none
320
321 type(pseudo_t), intent(in) :: pseudo
322 integer, intent(in) :: l
323 real(real64), intent(out) :: radial_function
324 end subroutine pseudo_radial_function
325
326 ! -------------------------------------------------
327
328 subroutine pseudo_nlcc_density(pseudo, nlcc_density)
329 import :: pseudo_t
330 import :: real64
331 implicit none
332
333 type(pseudo_t), intent(in) :: pseudo
334 real(real64), intent(out) :: nlcc_density
335 end subroutine pseudo_nlcc_density
336
337 ! -------------------------------------------------
339 subroutine pseudo_density(pseudo, density)
340 import :: pseudo_t
341 import :: real64
342 implicit none
343
344 type(pseudo_t), intent(in) :: pseudo
345 real(real64), intent(out) :: density
346 end subroutine pseudo_density
347
348 ! -------------------------------------------------
350 integer function pseudo_nwavefunctions(pseudo)
351 import :: pseudo_t
352 implicit none
353
354 type(pseudo_t), intent(in) :: pseudo
355 end function pseudo_nwavefunctions
356
357 ! -------------------------------------------------
358
359 subroutine pseudo_wavefunction(pseudo, index, n, l, occ, wf)
360 import :: pseudo_t
361 import :: real64
362 implicit none
363
364 type(pseudo_t), intent(in) :: pseudo
365 integer, intent(in) :: index
366 integer, intent(out) :: n
367 integer, intent(out) :: l
368 real(real64), intent(out) :: occ
369 real(real64), intent(out) :: wf
370 end subroutine pseudo_wavefunction
372 ! -------------------------------------------------
373
374 integer function pseudo_exchange(pseudo)
375 import :: pseudo_t
376 implicit none
377
378 type(pseudo_t), intent(in) :: pseudo
379 end function pseudo_exchange
380
381 ! -------------------------------------------------
382
383 integer function pseudo_correlation(pseudo)
384 import :: pseudo_t
385 implicit none
386
387 type(pseudo_t), intent(in) :: pseudo
388 end function pseudo_correlation
389
390 ! -------------------------------------------------
391
392 integer function pseudo_projector_2j(pseudo, l, ic)
393 import :: pseudo_t
394 implicit none
395
396 type(pseudo_t), intent(in) :: pseudo
397 integer, intent(in) :: l
398 integer, intent(in) :: ic
399 end function pseudo_projector_2j
400
401 ! -------------------------------------------------
402
403 integer function pseudo_wavefunction_2j(pseudo, ii)
404 import :: pseudo_t
405 implicit none
406
407 type(pseudo_t), intent(in) :: pseudo
408 integer, intent(in) :: ii
410
411 end interface
412
413contains
414
415 ! -------------------------------------------------
416
417 logical function pseudo_has_nlcc(pseudo)
418 type(pseudo_t), intent(in) :: pseudo
419
420 interface
421 integer function pseudo_has_nlcc_low(pseudo)
422 import :: pseudo_t
423 implicit none
424
425 type(pseudo_t), intent(in) :: pseudo
426 end function pseudo_has_nlcc_low
427 end interface
428
429 pseudo_has_nlcc = (pseudo_has_nlcc_low(pseudo) /= 0)
430
431 end function pseudo_has_nlcc
433 ! -------------------------------------------------
434
435 logical function pseudo_has_density(pseudo)
436 type(pseudo_t), intent(in) :: pseudo
437
438 interface
439 integer function pseudo_has_density_low(pseudo)
440 import :: pseudo_t
441 implicit none
442
443 type(pseudo_t), intent(in) :: pseudo
444 end function pseudo_has_density_low
445 end interface
446
447 pseudo_has_density = (pseudo_has_density_low(pseudo) /= 0)
448
449 end function pseudo_has_density
450
451
452 ! -------------------------------------------------
453
454 logical pure function pseudo_has_total_angular_momentum(pseudo)
455 type(pseudo_t), intent(in) :: pseudo
456
457 interface
458 integer pure function pseudo_has_total_angular_momentum_low(pseudo)
459 import :: pseudo_t
460 implicit none
461
462 type(pseudo_t), intent(in) :: pseudo
463 end function pseudo_has_total_angular_momentum_low
464 end interface
465
466 pseudo_has_total_angular_momentum = (pseudo_has_total_angular_momentum_low(pseudo) /= 0)
469
470 ! -------------------------------------------------
471
472 logical function pseudo_has_radial_function(pseudo, l)
473 type(pseudo_t), intent(in) :: pseudo
474 integer, intent(in) :: l
475
476 interface
477 integer function pseudo_has_radial_function_low(pseudo, l)
478 import :: pseudo_t
479 implicit none
480
481 type(pseudo_t), intent(in) :: pseudo
482 integer, intent(in) :: l
483 end function pseudo_has_radial_function_low
484 end interface
486 pseudo_has_radial_function = (pseudo_has_radial_function_low(pseudo, l) /= 0)
487
488 end function pseudo_has_radial_function
489
490end module pseudo_oct_m
491
492!! Local Variables:
493!! mode: f90
494!! coding: utf-8
495!! End:
integer, parameter, public pseudo_type_kleinman_bylander
Definition: pseudo.F90:158
logical function, public pseudo_has_radial_function(pseudo, l)
Definition: pseudo.F90:566
integer, parameter, public pseudo_format_qso
Definition: pseudo.F90:173
logical function, public pseudo_has_density(pseudo)
Definition: pseudo.F90:529
integer, parameter, public pseudo_format_hgh
Definition: pseudo.F90:173
integer, parameter, public pseudo_format_upf2
Definition: pseudo.F90:173
integer, parameter, public pseudo_type_paw
Definition: pseudo.F90:158
logical function, public pseudo_has_nlcc(pseudo)
Definition: pseudo.F90:511
integer, parameter, public pseudo_format_cpi
Definition: pseudo.F90:173
integer, parameter, public pseudo_status_unsupported_type
Definition: pseudo.F90:164
integer, parameter, public pseudo_status_file_not_found
Definition: pseudo.F90:164
integer, parameter, public pseudo_type_semilocal
Definition: pseudo.F90:158
integer, parameter, public pseudo_format_psp8
Definition: pseudo.F90:173
integer, parameter, public pseudo_status_unsupported_type_ultrasoft
Definition: pseudo.F90:164
logical pure function, public pseudo_has_total_angular_momentum(pseudo)
Definition: pseudo.F90:548
integer, parameter, public pseudo_format_psf
Definition: pseudo.F90:173
integer, parameter, public pseudo_status_unknown_format
Definition: pseudo.F90:164
integer, parameter, public pseudo_status_unsupported_type_paw
Definition: pseudo.F90:164
integer, parameter, public pseudo_format_fhi
Definition: pseudo.F90:173
integer, parameter, public pseudo_correlation_any
Definition: pseudo.F90:192
integer, parameter, public pseudo_format_upf1
Definition: pseudo.F90:173
integer, parameter, public pseudo_format_unknown
Definition: pseudo.F90:173
integer, parameter, public pseudo_status_format_not_supported
Definition: pseudo.F90:164
integer, parameter, public pseudo_exchange_any
Definition: pseudo.F90:188
integer, parameter, public pseudo_format_psml
Definition: pseudo.F90:173