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
23  implicit none
24
25  private
26
27  public ::                              &
28    pseudo_t,                            &
29    pseudo_detect_format,                &
30    pseudo_init,                         &
31    pseudo_end,                          &
32    pseudo_type,                         &
33    pseudo_format,                       &
34    pseudo_valence_charge,               &
35    pseudo_mesh_spacing,                 &
36    pseudo_mesh_size,                    &
37    pseudo_mass,                         &
38    pseudo_lmax,                         &
39    pseudo_llocal,                       &
40    pseudo_nchannels,                    &
41    pseudo_nprojectors,                  &
42    pseudo_grid,                         &
43    pseudo_grid_weights,                 &
44    pseudo_local_potential,              &
45    pseudo_projector,                    &
46    pseudo_has_radial_function,          &
47    pseudo_radial_function,              &
48    pseudo_radial_potential,             &
49    pseudo_has_nlcc,                     &
50    pseudo_nlcc_density,                 &
51    pseudo_dij,                          &
52    pseudo_has_density,                  &
53    pseudo_density,                      &
54    pseudo_nwavefunctions,               &
55    pseudo_wavefunction,                 &
56    pseudo_exchange,                     &
57    pseudo_correlation,                  &
58    pseudo_has_total_angular_momentum,   &
59    pseudo_projector_2j,                 &
60    pseudo_wavefunction_2j
61
62  !the following sets of values have to match with those on base.hpp
63  integer, parameter, public ::               &
64    PSEUDO_TYPE_ULTRASOFT         = 30,       &
65    PSEUDO_TYPE_SEMILOCAL         = 31,       &
66    PSEUDO_TYPE_KLEINMAN_BYLANDER = 32,       &
67    PSEUDO_TYPE_PAW               = 33
68
69  integer, parameter, public ::                       &
70    PSEUDO_STATUS_SUCCESS                      = 0,   &
71    PSEUDO_STATUS_FILE_NOT_FOUND               = 455, &
72    PSEUDO_STATUS_FORMAT_NOT_SUPPORTED         = 456, &
73    PSEUDO_STATUS_UNKNOWN_FORMAT               = 457, &
74    PSEUDO_STATUS_UNSUPPORTED_TYPE_ULTRASOFT   = 458, &
75    PSEUDO_STATUS_UNSUPPORTED_TYPE_PAW         = 459, &
76    PSEUDO_STATUS_UNSUPPORTED_TYPE             = 460
77
78  integer, parameter, public ::                       &
79    PSEUDO_FORMAT_FILE_NOT_FOUND             = 773,   &
80    PSEUDO_FORMAT_UNKNOWN                    = 774,   &
81    PSEUDO_FORMAT_UPF1                       = 775,   &
82    PSEUDO_FORMAT_UPF2                       = 776,   &
83    PSEUDO_FORMAT_QSO                        = 777,   &
84    PSEUDO_FORMAT_PSML                       = 778,   &
85    PSEUDO_FORMAT_PSF                        = 779,   &
86    PSEUDO_FORMAT_CPI                        = 780,   &
87    PSEUDO_FORMAT_FHI                        = 781,   &
88    PSEUDO_FORMAT_HGH                        = 782,   &
89    PSEUDO_FORMAT_PSP8                       = 783
90
91  ! we only define these values here, the specific functionals are
92  ! obtained from libxc
93  integer, parameter, public ::                       &
94    PSEUDO_EXCHANGE_UNKNOWN                  = -2,    &
95    PSEUDO_EXCHANGE_ANY                      = -1
96
97  integer, parameter, public ::                       &
98    PSEUDO_CORRELATION_UNKNOWN               = -2,    &
99    PSEUDO_CORRELATION_ANY                   = -1
100
101  type pseudo_t
102    private
103    integer(8) :: dummy
104  end type pseudo_t
105
106  interface
107
108    ! -------------------------------------------------
109
110    integer function pseudo_detect_format(filename)
111      import :: pseudo_t
112      implicit none
113
114      character(len=*), intent(in)    :: filename
115    end function pseudo_detect_format
116
117    ! -------------------------------------------------
118
119    subroutine pseudo_init(pseudo, filename, fmt, ierr)
120      import :: pseudo_t
121      implicit none
122
123      type(pseudo_t),   intent(out)   :: pseudo
124      character(len=*), intent(in)    :: filename
125      integer,          intent(in)    :: fmt
126      integer,          intent(out)   :: ierr
127    end subroutine pseudo_init
128
129    ! -------------------------------------------------
130
131    subroutine pseudo_end(pseudo)
132      import :: pseudo_t
133      implicit none
134
135      type(pseudo_t),   intent(inout) :: pseudo
136    end subroutine pseudo_end
137
138    ! -------------------------------------------------
139
140    integer function pseudo_type(pseudo)
141      import :: pseudo_t
142      implicit none
143
144      type(pseudo_t),   intent(in)    :: pseudo
145    end function pseudo_type
146
147    ! -------------------------------------------------
148
149    integer function pseudo_format(pseudo)
150      import :: pseudo_t
151      implicit none
152
153      type(pseudo_t),   intent(in)    :: pseudo
154    end function pseudo_format
155
156    ! -------------------------------------------------
157
158    FLOAT function pseudo_valence_charge(pseudo)
159      import :: pseudo_t
160      implicit none
161
162      type(pseudo_t),   intent(in)    :: pseudo
163    end function pseudo_valence_charge
164
165    ! -------------------------------------------------
166
167    FLOAT function pseudo_mesh_spacing(pseudo)
168      import :: pseudo_t
169      implicit none
170
171      type(pseudo_t),   intent(in)    :: pseudo
172    end function pseudo_mesh_spacing
173
174    ! -------------------------------------------------
175
176    integer function pseudo_mesh_size(pseudo)
177      import :: pseudo_t
178      implicit none
179
180      type(pseudo_t),   intent(in)    :: pseudo
181    end function pseudo_mesh_size
182
183    ! -------------------------------------------------
184
185    FLOAT function pseudo_mass(pseudo)
186      import :: pseudo_t
187      implicit none
188
189      type(pseudo_t),   intent(in)    :: pseudo
190    end function pseudo_mass
191
192    ! -------------------------------------------------
193
194    integer function pseudo_lmax(pseudo)
195      import :: pseudo_t
196      implicit none
197
198      type(pseudo_t),   intent(in)    :: pseudo
199    end function pseudo_lmax
200
201    ! -------------------------------------------------
202
203    integer function pseudo_llocal(pseudo)
204      import :: pseudo_t
205      implicit none
206
207      type(pseudo_t),   intent(in)    :: pseudo
208    end function pseudo_llocal
209
210    ! -------------------------------------------------
211
212    integer function pseudo_nchannels(pseudo)
213      import :: pseudo_t
214      implicit none
215
216      type(pseudo_t),   intent(in)    :: pseudo
217    end function pseudo_nchannels
218
219    ! -------------------------------------------------
220
221    integer function pseudo_nprojectors(pseudo)
222      import :: pseudo_t
223      implicit none
224
225      type(pseudo_t),   intent(in)    :: pseudo
226    end function pseudo_nprojectors
227
228    ! -------------------------------------------------
229
230    subroutine pseudo_grid(pseudo, grid)
231      import :: pseudo_t
232      implicit none
233
234      type(pseudo_t),   intent(in)    :: pseudo
235      FLOAT,            intent(out)   :: grid
236    end subroutine pseudo_grid
237
238    ! -------------------------------------------------
239
240    subroutine pseudo_grid_weights(pseudo, weight)
241      import :: pseudo_t
242      implicit none
243
244      type(pseudo_t),   intent(in)    :: pseudo
245      FLOAT,            intent(out)   :: weight
246    end subroutine pseudo_grid_weights
247
248    ! -------------------------------------------------
249
250    subroutine pseudo_local_potential(pseudo, local_potential)
251      import :: pseudo_t
252      implicit none
253
254      type(pseudo_t),   intent(in)    :: pseudo
255      FLOAT,            intent(out)   :: local_potential
256    end subroutine pseudo_local_potential
257
258    ! -------------------------------------------------
259
260    subroutine pseudo_projector(pseudo, l, ic, projector)
261      import :: pseudo_t
262      implicit none
263
264      type(pseudo_t),   intent(in)    :: pseudo
265      integer,          intent(in)    :: l
266      integer,          intent(in)    :: ic
267      FLOAT,            intent(out)   :: projector
268    end subroutine pseudo_projector
269
270    ! -------------------------------------------------
271
272    FLOAT function pseudo_dij(pseudo, l, ic, jc)
273      import :: pseudo_t
274      implicit none
275
276      type(pseudo_t),   intent(in)    :: pseudo
277      integer,          intent(in)    :: l
278      integer,          intent(in)    :: ic
279      integer,          intent(in)    :: jc
280    end function pseudo_dij
281
282    ! -------------------------------------------------
283
284    subroutine pseudo_radial_potential(pseudo, l, radial_potential)
285      import :: pseudo_t
286      implicit none
287
288      type(pseudo_t),   intent(in)    :: pseudo
289      integer,          intent(in)    :: l
290      FLOAT,            intent(inout) :: radial_potential
291    end subroutine pseudo_radial_potential
292
293    ! -------------------------------------------------
294
295    subroutine pseudo_radial_function(pseudo, l, radial_function)
296      import :: pseudo_t
297      implicit none
298
299      type(pseudo_t),   intent(in)    :: pseudo
300      integer,          intent(in)    :: l
301      FLOAT,            intent(out)   :: radial_function
302    end subroutine pseudo_radial_function
303
304    ! -------------------------------------------------
305
306    subroutine pseudo_nlcc_density(pseudo, nlcc_density)
307      import :: pseudo_t
308      implicit none
309
310      type(pseudo_t),   intent(in)    :: pseudo
311      FLOAT,            intent(out)   :: nlcc_density
312    end subroutine pseudo_nlcc_density
313
314    ! -------------------------------------------------
315
316    subroutine pseudo_density(pseudo, density)
317      import :: pseudo_t
318      implicit none
319
320      type(pseudo_t),   intent(in)    :: pseudo
321      FLOAT,            intent(out)   :: density
322    end subroutine pseudo_density
323
324    ! -------------------------------------------------
325
326    integer function pseudo_nwavefunctions(pseudo)
327      import :: pseudo_t
328      implicit none
329
330      type(pseudo_t),   intent(in)    :: pseudo
331    end function pseudo_nwavefunctions
332
333    ! -------------------------------------------------
334
335    subroutine pseudo_wavefunction(pseudo, index, n, l, occ, wf)
336      import :: pseudo_t
337      implicit none
338
339      type(pseudo_t),   intent(in)    :: pseudo
340      integer,          intent(in)    :: index
341      integer,          intent(out)   :: n
342      integer,          intent(out)   :: l
343      FLOAT,            intent(out)   :: occ
344      FLOAT,            intent(out)   :: wf
345    end subroutine pseudo_wavefunction
346
347    ! -------------------------------------------------
348
349    integer function pseudo_exchange(pseudo)
350      import :: pseudo_t
351      implicit none
352
353      type(pseudo_t),   intent(in)    :: pseudo
354    end function pseudo_exchange
355
356    ! -------------------------------------------------
357
358    integer function pseudo_correlation(pseudo)
359      import :: pseudo_t
360      implicit none
361
362      type(pseudo_t),   intent(in)    :: pseudo
363    end function pseudo_correlation
364
365    ! -------------------------------------------------
366
367    integer function pseudo_projector_2j(pseudo, l, ic)
368      import :: pseudo_t
369      implicit none
370
371      type(pseudo_t),   intent(in)    :: pseudo
372      integer,          intent(in)    :: l
373      integer,          intent(in)    :: ic
374    end function pseudo_projector_2j
375
376    ! -------------------------------------------------
377
378    integer function pseudo_wavefunction_2j(pseudo, ii)
379      import :: pseudo_t
380      implicit none
381
382      type(pseudo_t),   intent(in)    :: pseudo
383      integer,          intent(in)    :: ii
384    end function pseudo_wavefunction_2j
385
386  end interface
387
388contains
389
390  ! -------------------------------------------------
391
392  logical function pseudo_has_nlcc(pseudo)
393    type(pseudo_t),   intent(in)      :: pseudo
394
395    interface
396      integer function pseudo_has_nlcc_low(pseudo)
397        import :: pseudo_t
398        implicit none
399
400        type(pseudo_t),   intent(in)      :: pseudo
401      end function pseudo_has_nlcc_low
402    end interface
403
404    pseudo_has_nlcc = (pseudo_has_nlcc_low(pseudo) /= 0)
405
406  end function pseudo_has_nlcc
407
408  ! -------------------------------------------------
409
410  logical function pseudo_has_density(pseudo)
411    type(pseudo_t),   intent(in)      :: pseudo
412
413    interface
414      integer function pseudo_has_density_low(pseudo)
415        import :: pseudo_t
416        implicit none
417
418        type(pseudo_t),   intent(in)      :: pseudo
419      end function pseudo_has_density_low
420    end interface
421
422    pseudo_has_density = (pseudo_has_density_low(pseudo) /= 0)
423
424  end function pseudo_has_density
425    ! -------------------------------------------------
426
427  logical function pseudo_has_total_angular_momentum(pseudo)
428    type(pseudo_t),   intent(in)      :: pseudo
429
430    interface
431      integer function pseudo_has_total_angular_momentum_low(pseudo)
432        import :: pseudo_t
433        implicit none
434
435        type(pseudo_t),   intent(in)      :: pseudo
436      end function pseudo_has_total_angular_momentum_low
437    end interface
438
439    pseudo_has_total_angular_momentum = (pseudo_has_total_angular_momentum_low(pseudo) /= 0)
440
441  end function pseudo_has_total_angular_momentum
442
443  ! -------------------------------------------------
444
445  logical function pseudo_has_radial_function(pseudo, l)
446    type(pseudo_t),   intent(in)      :: pseudo
447    integer,          intent(in)      :: l
448
449    interface
450      integer function pseudo_has_radial_function_low(pseudo, l)
451        import :: pseudo_t
452        implicit none
453
454        type(pseudo_t),   intent(in)      :: pseudo
455        integer,          intent(in)      :: l
456      end function pseudo_has_radial_function_low
457    end interface
458
459    pseudo_has_radial_function = (pseudo_has_radial_function_low(pseudo, l) /= 0)
460
461  end function pseudo_has_radial_function
462
463end module pseudo_oct_m
464
465!! Local Variables:
466!! mode: f90
467!! coding: utf-8
468!! End:
469