Octopus
species_factory.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2023-2024 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19#include "global.h"
20
23 use debug_oct_m
25 use global_oct_m
26 use iihash_oct_m
27 use io_oct_m
31 use parser_oct_m
36 use string_oct_m
37 use unit_oct_m
39
40 implicit none
41
42 private
43 public :: &
47
48 type :: species_factory_t
49 logical :: initialized = .false.
50 integer :: default_pseudopotential_set_id
51 type(pseudo_set_t) :: default_pseudopotential_set
52 integer :: default_allelectron_type
53 real(real64) :: default_sigma
54 real(real64) :: default_anc_a
55 contains
56 procedure :: init => species_factory_init
57 procedure :: end => species_factory_end
58 procedure :: create_from_input => species_factory_create_from_input
59 procedure :: create_from_block => read_from_block
60 end type species_factory_t
61
62contains
63
64 ! ---------------------------------------------------------
65 subroutine species_factory_init(factory, namespace)
66 class(species_factory_t), intent(inout) :: factory
67 type(namespace_t), intent(in) :: namespace
68
69 integer :: ierr, default_val
70
71 if (factory%initialized) return
72
73 push_sub(species_factory_init)
74
75 factory%initialized = .true.
76
77 call share_directory_set(conf%share)
78
79 !%Variable AllElectronType
80 !%Type integer
81 !%Default no
82 !%Section System::Species
83 !%Description
84 !% Selects the type of all-electron species that applies by default to all
85 !% atoms. This is not compatible with <tt>PseudopotentialSet</tt>, but it is
86 !% compatible with the <tt>Species</tt> block.
87 !%
88 !%Option no 0
89 !% Do not specify any default all-electron type of species. All species must be
90 !% specified in the Species block.
91 !%Option full_delta 1
92 !% All atoms are supposed to be by default of type <tt>species_full_delta</tt>.
93 !%Option full_gaussian 2
94 !% All atoms are supposed to be by default of type <tt>species_full_gaussian</tt>.
95 !%Option full_anc 3
96 !% All atoms are supposed to be by default of type <tt>species_full_anc</tt>.
97 !%End
98 call parse_variable(namespace, 'AllElectronType', option__allelectrontype__no, factory%default_allelectron_type)
99 call messages_print_var_option('AllElectronType', factory%default_allelectron_type, namespace=namespace)
100
101 !%Variable AllElectronSigma
102 !%Type integer
103 !%Default 0.6
104 !%Section System::Species
105 !%Description
106 !% Default value for the parameter <tt>gaussian_width</tt>. This is useful
107 !% for specifying multiple atoms without specifying the species block. The
108 !% default value is taken from the recommendation in
109 !% <i>Phys. Rev. B</i> <b>55</b>, 10289 (1997).
110 !%
111 !%End
112 call parse_variable(namespace, 'AllElectronSigma', 0.6_real64, factory%default_sigma)
113
114 !%Variable AllElectronANCParam
115 !%Type integer
116 !%Default 4
117 !%Section System::Species
118 !%Description
119 !% Default values for the parameter <tt>anc_a</tt>. This is usefull
120 !% for specifying multiple atoms without specifying the species block.
121 !%
122 !%End
123 call parse_variable(namespace, 'AllElectronANCParam', 4.0_real64, factory%default_anc_a)
124
125
126 !%Variable PseudopotentialSet
127 !%Type integer
128 !%Default standard
129 !%Section System::Species
130 !%Description
131 !% Selects the set of pseudopotentials used by default for species
132 !% not defined in the <tt>Species</tt> block.
133 !%
134 !% These sets of pseudopotentials come from different
135 !% sources. Octopus developers have not validated them. We include
136 !% them with the code for convenience of the users, but you are
137 !% expected to check the quality and suitability of the
138 !% pseudopotential for your application.
139 !%
140 !%Option none 0
141 !% Do not load any pseudopotential by default. All species must be
142 !% specified in the Species block.
143 !%Option standard 1
144 !% The standard set of Octopus that provides LDA pseudopotentials
145 !% in the PSF format for some elements: H, Li, C, N, O, Na, Si, S, Ti, Se, Cd.
146 !%Option sg15 2
147 !% The set of Optimized Norm-Conserving Vanderbilt
148 !% PBE pseudopotentials. Ref: M. Schlipf and F. Gygi, <i>Comp. Phys. Commun.</i> <b>196</b>, 36 (2015).
149 !% This set provides pseudopotentials for elements up to Z = 83
150 !% (Bi), excluding Lanthanides.
151 !% Current version of the set is 1.2.
152 !%Option hgh_lda 3
153 !% The set of Hartwigsen-Goedecker-Hutter LDA pseudopotentials for elements from H to Rn.
154 !% Ref: C. Hartwigsen, S. Goedecker, and J. Hutter, <i>Phys. Rev. B</i> <b>58</b>, 3641 (1998).
155 !%Option hgh_lda_sc 31
156 !% The semicore set of Hartwigsen-Goedecker-Hutter LDA pseudopotentials.
157 !% Ref: C. Hartwigsen, S. Goedecker, and J. Hutter, <i>Phys. Rev. B</i> <b>58</b>, 3641 (1998).
158 !%Option hscv_lda 4
159 !% The set of Hamann-Schlueter-Chiang-Vanderbilt (HSCV) potentials
160 !% for LDA exchange and correlation downloaded from http://fpmd.ucdavis.edu/potentials/index.htm.
161 !% These pseudopotentials were originally intended for the QBox
162 !% code. They were generated using the method of Hamann, Schluter
163 !% and Chiang. Ref: D. Vanderbilt, <i>Phys. Rev. B</i> <b>32</b>, 8412 (1985).
164 !% Warning from the original site: The potentials provided in this
165 !% site are distributed without warranty. In most cases,
166 !% potentials were not tested. Potentials should be thoroughly
167 !% tested before being used in simulations.
168 !%Option hscv_pbe 5
169 !% PBE version of the HSCV pseudopotentials. Check the
170 !% documentation of the option <tt>hscv_lda</tt> for details and warnings.
171 !%Option pseudodojo_pbe 100
172 !% PBE version of the pseudopotentials of http://pseudo-dojo.org. Version 0.4.
173 !%Option pseudodojo_lda 103
174 !% LDA pseudopotentials of http://pseudo-dojo.org. Version 0.4.
175 !%Option pseudodojo_pbesol 105
176 !% PBEsol version of the pseudopotentials of http://pseudo-dojo.org. Version 0.3.
177 !%End
178 default_val = option__pseudopotentialset__standard
179 if(factory%default_allelectron_type /= option__allelectrontype__no) default_val = option__pseudopotentialset__none
180 call parse_variable(namespace, 'PseudopotentialSet', default_val, factory%default_pseudopotential_set_id)
181 call messages_print_var_option('PseudopotentialSet', factory%default_pseudopotential_set_id, namespace=namespace)
182 if (factory%default_pseudopotential_set_id /= option__pseudopotentialset__none) then
183 call pseudo_set_init(factory%default_pseudopotential_set, get_set_directory(factory%default_pseudopotential_set_id), ierr)
184 else
185 call pseudo_set_nullify(factory%default_pseudopotential_set)
186 end if
187
188 if (factory%default_pseudopotential_set_id /= option__pseudopotentialset__none .and. factory%default_allelectron_type /= option__allelectrontype__no) then
189 message(1) = "PseudopotentialSet /= none cannot be used with AllElectronType /= no."
190 call messages_fatal(1, namespace=namespace)
191 end if
192
193 pop_sub(species_factory_init)
194 end subroutine species_factory_init
195
196 ! ---------------------------------------------------------
197
198 subroutine species_factory_end(factory)
199 class(species_factory_t), intent(inout) :: factory
200
201 push_sub(species_factory_end)
202
203 if (factory%initialized) then
204 call pseudo_set_end(factory%default_pseudopotential_set)
205 factory%initialized = .false.
206 end if
207
208 pop_sub(species_factory_end)
209 end subroutine species_factory_end
210
211
212 ! ---------------------------------------------------------
217 ! ---------------------------------------------------------
218 function species_factory_create_from_input(factory, namespace, label, index) result(spec)
219 class(species_factory_t), intent(in) :: factory
220 type(namespace_t), intent(in) :: namespace
221 character(len=*), intent(in) :: label
222 integer, intent(in) :: index
223 class(species_t), pointer :: spec
224
225 character(len=LABEL_LEN) :: lab
226 integer :: ib, row, n_spec_block, read_data
227 type(block_t) :: blk
228 type(element_t) :: element
229
231
232 read_data = 0
233
234 !%Variable Species
235 !%Type block
236 !%Section System::Species
237 !%Description
238 !% A species is by definition either an "ion" (nucleus + core electrons) described
239 !% through a pseudopotential, or a model potential.
240 !%
241 !% Note that some sets of pseudopotentials are distributed with
242 !% the code. To use these pseudopotentials, you do not need to define them
243 !% explicitly in the <tt>Species</tt> block, as default parameters
244 !% are provided.
245 !% You can select the set for default pseudopotentials using the
246 !% <tt>PseudopotentialSet</tt> variable.
247 !%
248 !% Supported norm-conserving pseudopotential formats are
249 !% detected by the file extension: UPF (<tt>.upf</tt>), PSF (SIESTA, <tt>.psf</tt>), FHI (ABINIT 6, <tt>.fhi</tt>),
250 !% CPI (Fritz-Haber, <tt>.cpi</tt>), QSO (quantum-simulation.org, for Qbox, <tt>.xml</tt>),
251 !% HGH (Hartwigsen-Goedecker-Hutter, <tt>.hgh</tt>).
252 !% PSPIO format can also be used via <tt>species_pspio</tt> if that library is linked.
253 !% Note: pseudopotentials may only be used in 3D.
254 !%
255 !% The format of this block is the following: The first field is a
256 !% string that defines the name of the species. The second field
257 !% defines the type of species (the valid options are detailed
258 !% below).
259 !%
260 !% Then a list of parameters follows. The parameters are specified
261 !% by a first field with the parameter name and the field that
262 !% follows with the value of the parameter. Some parameters are
263 !% specific to a certain species while others are accepted by all
264 !% species. These are <tt>mass</tt>, <tt>max_spacing</tt>, and <tt>min_radius</tt>.
265 !%
266 !% These are examples of possible species:
267 !%
268 !% <tt>%Species
269 !% <br>&nbsp;&nbsp;'O' | species_pseudo | file | 'O.psf' | lmax | 1 | lloc | 1
270 !% <br>&nbsp;&nbsp;'H' | species_pseudo | file | '../H.hgh'
271 !% <br>&nbsp;&nbsp;'Xe' | species_pseudo | set | pseudojo_pbe_stringent
272 !% <br>&nbsp;&nbsp;'C' | species_pseudo | file | "carbon.xml"
273 !% <br>&nbsp;&nbsp;'jlm' | species_jellium | jellium_radius | 5.0
274 !% <br>&nbsp;&nbsp;'rho' | species_charge_density | density_formula | "exp(-r/a)" | mass | 17.0 | valence | 6
275 !% <br>&nbsp;&nbsp;'udf' | species_user_defined | potential_formula | "1/2*r^2" | valence | 8
276 !% <br>&nbsp;&nbsp;'He_all' | species_full_delta
277 !% <br>&nbsp;&nbsp;'H_all' | species_full_gaussian | gaussian_width | 0.2
278 !% <br>&nbsp;&nbsp;'Li1D' | species_soft_coulomb | softening | 1.5 | valence | 3
279 !% <br>&nbsp;&nbsp;'H_all' | species_full_anc | anc_a | 4
280 !% <br>%</tt>
281 !%Option species_pseudo -7
282 !% The species is a pseudopotential. How to get the
283 !% pseudopotential can be specified by the <tt>file</tt> or
284 !% the <tt>set</tt> parameters. If both are missing, the
285 !% pseudopotential will be taken from the <tt>PseudopotentialSet</tt>
286 !% specified for the run, this is useful if you want to change
287 !% some parameters of the pseudo, like the <tt>mass</tt>.
288 !%
289 !% The optional parameters for this type of species are
290 !% <tt>lmax</tt>, that defines the maximum angular momentum
291 !% component to be used, and <tt>lloc</tt>, that defines the
292 !% angular momentum to be considered as local. When these
293 !% parameters are not set, the value for lmax is the maximum
294 !% angular component from the pseudopotential file. The default
295 !% value for <tt>lloc</tt> is taken from the pseudopotential if
296 !% available, if not, it is set to 0. Note that, depending on the
297 !% type of pseudopotential, it might not be possible to select
298 !% <tt>lmax</tt> and <tt>lloc</tt>, if that is the case the
299 !% parameters will be ignored.
300 !%
301 !%Option species_pspio -110
302 !% (experimental) Alternative method to read pseudopotentials
303 !% using the PSPIO library. This species uses the same parameters
304 !% as <tt>species_pseudo</tt>.
305 !%Option species_user_defined -123
306 !% Species with user-defined potential. The potential for the
307 !% species is defined by the formula given by the <tt>potential_formula</tt>
308 !% parameter.
309 !% The
310 !% <tt>valence</tt> parameter determines the number of electrons
311 !% associated with the species. By default, a valence of 0 is assumed.
312 !%Option species_charge_density -125
313 !% The potential for this species is created from the distribution
314 !% of charge given by the <tt>density_formula</tt> parameter.
315 !% The
316 !% <tt>valence</tt> parameter determines the number of electrons
317 !% associated with the species. By default, a valence of 0 is assumed.
318 !%Option species_point -3
319 !%Option species_jellium -3
320 !% Jellium sphere.
321 !% The charge associated with this species must be given by the <tt>valence</tt> parameter.
322 !%Option species_jellium_slab -4
323 !% A slab of jellium that extends across the simulation box in the
324 !% <i>xy</i>-plane. The dimension along the <i>z</i> direction is
325 !% determined by the required parameter <tt>thickness</tt>.
326 !% The charge associated with this species must be given by the <tt>valence</tt> parameter.
327 !%Option species_full_delta -127
328 !% Full atomic potential represented by a delta charge
329 !% distribution. The atom will be displaced to the nearest grid
330 !% point. The atomic number is determined from the name of the species.
331 !%Option species_full_anc -130
332 !% Analytical norm-conserving regulized Coulomb potential from
333 !% [Gygi J. Chem. Theory Comput. 2023, 19, 1300−1309].
334 !%Option species_full_gaussian -124
335 !% A full-potential atom is defined by a Gaussian accumulation of
336 !% positive charge (distorted if curvilinear coordinates are
337 !% used), in the form:
338 !%
339 !% <math>q(r) = z \beta \exp[ - (\vec{r}-\vec{r_0})^2 / (\sqrt{2} \delta \sigma) ] </math>
340 !%
341 !% <math>\beta</math> is chosen in order to maintain proper
342 !% normalization (the integral of <math>q</math> should sum up to
343 !% <math>z</math>). <math>\delta</math> is the grid spacing (the
344 !% grid spacing in the first dimension, to be precise).
345 !% <math>\vec{r_0}</math> is calculated in such a way that the the
346 !% first moment of <math>q(r)/z</math> is equal to the atomic
347 !% position. For a precise description, see N. A. Modine,
348 !% <i>Phys. Rev. B</i> <b>55</b>, 10289 (1997). The width of the
349 !% Gaussian is set by parameter <tt>gaussian_width</tt>. The
350 !% atomic number is determined from the name of the species.
351 !%Option species_from_file -126
352 !% The potential is read from a file. Accepted file formats, detected by extension: obf, ncdf and csv.
353 !% The
354 !% <tt>valence</tt> parameter determines the number of electrons
355 !% associated with the species. By default, a valence of 0 is assumed.
356 !%Option species_soft_coulomb -128
357 !% The potential is a soft-Coulomb function, <i>i.e.</i> a function in the form:
358 !%
359 !% <math>v(r) = - z_{val} / \sqrt{a^2 + r^2}</math>
360 !%
361 !% The value of <i>a</i> should be given by the mandatory <tt>softening</tt> parameter.
362 !% The charge associated with this species must be given by the <tt>valence</tt> parameter.
363 !%Option species_jellium_charge_density -129
364 !% The parameter is the name of a volume block specifying the shape of the jellium.
365 !%Option lmax -10003
366 !% The maximum angular-momentum channel that will be used for the pseudopotential.
367 !%Option lloc -10004
368 !% The angular-momentum channel of the pseudopotential to be considered local.
369 !%Option mass -10005
370 !% The mass of the species in atomic mass units, <i>i.e.</i> the mass of a proton is
371 !% roughly one. It is set automatically for pseudopotentials from the
372 !% <a href=http://www.nist.gov/pml/data/comp.cfm>NIST values</a>.
373 !% For other species, the default is 1.0.
374 !%Option valence -10006
375 !% The number of electrons of the species. It is set automatically from the name of the species.
376 !% if it correspond to the name in the periodic table. If not specified and if the name
377 !% does not match an atom name, a value of 0 is assumed.
378 !%Option jellium_radius -10007
379 !% The radius of the sphere for <tt>species_jellium</tt>. If this value is not specified,
380 !% the default of 0.5 bohr is used.
381 !%Option set -10017
382 !% For a <tt>species_pseudo</tt>, get the pseudopotential from a
383 !% particular set. This flag must be followed with one of the
384 !% valid values for the variable <tt>PseudopotentialSet</tt>.
385 !%Option gaussian_width -10008
386 !% The width of the Gaussian (in units of spacing) used to represent
387 !% the nuclear charge for <tt>species_full_gaussian</tt>. If not present,
388 !% the default is 0.6.
389 !%Option softening -10009
390 !% The softening parameter <i>a</i> for <tt>species_soft_coulomb</tt> in units of length.
391 !%Option file -10010
392 !% The path for the file that describes the species.
393 !%Option db_file -10011
394 !% Obsolete. Use the <tt>set</tt> option of the <tt>PseudopotentialSet</tt> variable instead.
395 !%Option potential_formula -10012
396 !% Mathematical expression that defines the potential for <tt>species_user_defined</tt>. You can use
397 !% any of the <i>x</i>, <i>y</i>, <i>z</i> or <i>r</i> variables.
398 !%Option density_formula -10013
399 !% Mathematical expression that defines the charge density for <tt>species_charge_density</tt>. You can use
400 !% any of the <i>x</i>, <i>y</i>, <i>z</i> or <i>r</i> variables.
401 !%Option thickness -10014
402 !% The thickness of the slab for species_jellium_slab. Must be positive.
403 !%Option vdw_radius -10015
404 !% The van der Waals radius that will be used for this species.
405 !%Option volume -10016
406 !% Name of a volume block
407 !%Option hubbard_l -10018
408 !% The angular-momentum for which the effective U will be applied.
409 !%Option hubbard_u -10019
410 !% The effective U that will be used for the DFT+U calculations.
411 !%Option hubbard_j -10020
412 !% The value of j (hubbard_l-1/2 or hubbard_l+1/2) on which the effective U is applied.
413 !%Option hubbard_alpha -10021
414 !% The strength of the potential constraining the occupations of the localized subspace
415 !% as defined in PRB 71, 035105 (2005)
416 !%Option anc_a -10022
417 !% The value of the parameter a of the ANC potential, as defined in [Gygi, JCTC 2023, 19, 1300−1309].
418 !% This parameter has the unit of inverse length and determines the range of regularization.
419 !%End
420
421 call messages_obsolete_variable(namespace, 'SpecieAllElectronSigma', 'Species')
422 call messages_obsolete_variable(namespace, 'SpeciesAllElectronSigma', 'Species')
423
424 ! First, find out if there is a Species block.
425 n_spec_block = 0
426 if (parse_block(namespace, 'Species', blk) == 0) then
427 n_spec_block = parse_block_n(blk)
428 end if
429
430 ! Find out if the sought species is in the block
431 row = -1
432 block: do ib = 1, n_spec_block
433 call parse_block_string(blk, ib-1, 0, lab)
434 if (trim(lab) == trim(label)) then
435 row = ib - 1
436 exit block
437 end if
438 end do block
439
440 ! Read whatever may be read from the block
441 if (row >= 0) then
442 spec => factory%create_from_block(namespace, blk, row, label, index, read_data)
443 call parse_block_end(blk)
444
445 assert(read_data > 0)
446
448 return
449 end if
450
451 ! We get here if there is a Species block but it does not contain
452 ! the species we are looking for.
453 if (n_spec_block > 0) call parse_block_end(blk)
454
455 ! Initialize all electron species (except soft-Coulomb) from specified allelectron type
456 if(factory%default_allelectron_type /= 0) then
457 select case(factory%default_allelectron_type)
458 case(option__allelectrontype__full_delta)
459 spec => full_delta_t(label, index, factory%default_sigma)
460 case(option__allelectrontype__full_gaussian)
461 spec => full_gaussian_t(label, index, factory%default_sigma)
462 case(option__allelectrontype__full_anc)
463 spec => full_anc_t(label, index, factory%default_anc_a)
464 case default
465 assert(.false.)
466 end select
467
468 ! get the mass, vdw radius and atomic number for this element
469 call element_init(element, get_symbol(spec%get_label()))
470
471 assert(element_valid(element))
472
473 call spec%set_z(real(element_atomic_number(element), real64))
474 call spec%set_zval(spec%get_z())
475 call spec%set_mass(element_mass(element))
476 call spec%set_vdw_radius(element_vdw_radius(element))
477
478 call element_end(element)
479
480 else ! Pseudopotential from specified set or the default one
481 spec => pseudopotential_t(label, index)
482 select type(spec)
483 type is(pseudopotential_t)
484 call read_from_set(spec, factory%default_pseudopotential_set_id, factory%default_pseudopotential_set, read_data)
485
486 if (read_data == 0) then
487 call messages_write( 'Species '//trim(spec%get_label())//' not found in default pseudopotential set.', new_line=.true. )
488 call messages_write('( '//trim(get_set_directory(factory%default_pseudopotential_set_id))//' )')
489 call messages_fatal(namespace=namespace)
490 end if
491 end select
492 end if
493
496
497!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
498! Private procedures
499
500 ! ---------------------------------------------------------
502 function read_from_block(factory, namespace, blk, row, label, index, read_data) result(spec)
503 class(species_factory_t), intent(in) :: factory
504 type(namespace_t), intent(in) :: namespace
505 type(block_t), intent(in) :: blk
506 integer, intent(in) :: row
507 character(len=*), intent(in) :: label
508 integer, intent(in) :: index
509 integer, intent(out):: read_data
510 class(species_t), pointer :: spec
511
512 integer :: ncols, icol, flag, set_read_data, ierr, type
513 type(element_t) :: element
514 type(iihash_t) :: read_parameters
515 integer :: user_lmax, user_llocal, hubbard_l, pseudopotential_set_id
516 real(real64) :: hubbard_u, hubbard_j, hubbard_alpha, mass, z_val, jradius, jthick, vdw_radius, aa
517 real(real64) :: sigma, softening
518 character(len=MAX_PATH_LEN) :: filename
519
520 integer, parameter :: &
521 species_jellium = 3, & !< jellium sphere.
524 species_pseudo = 7, &
525 species_pspio = 110, &
526 species_usdef = 123, &
527 species_full_gaussian = 124, &
528 species_charge_density = 125, &
529 species_from_file = 126, &
530 species_full_delta = 127, &
531 species_soft_coulomb = 128, &
532 species_full_anc = 130
533
534 push_sub(read_from_block)
535
536 ncols = parse_block_cols(blk, row)
537 read_data = 0
538
539 call parse_block_integer(blk, row, 1, type)
540
541 ! To detect the old species block format, options are represented
542 ! as negative values. If we get a non-negative value we know we
543 ! are reading a mass.
544 if (type >= 0) then
545 call messages_write('Found a species with the old format. Please update', new_line = .true.)
546 call messages_write('the Species block to the new format, where the second', new_line = .true.)
547 call messages_write('column indicates the type of the species.')
548 call messages_fatal(namespace=namespace)
549 end if
550
551 ! now we convert back to positive
552 type = -type
553
554 read_data = 2
555
556
557 select case (type)
558
559 case (species_soft_coulomb)
560 spec => soft_coulomb_t(label, index)
561
562 case (species_usdef) ! user-defined
563 spec => species_user_defined_t(label, index)
564
565 case (species_from_file)
566 spec => species_from_file_t(label, index)
567
568 case (species_jellium)
569 spec => jellium_sphere_t(label, index)
570
572 spec => jellium_slab_t(label, index)
573
574 case (species_full_delta)
575 spec => full_delta_t(label, index, factory%default_sigma)
576
577 case (species_full_gaussian)
578 spec => full_gaussian_t(label, index, factory%default_sigma)
579
580 case (species_full_anc)
581 spec => full_anc_t(label, index, factory%default_anc_a)
582
584 spec => species_charge_density_t(label, index)
585
587 spec => jellium_charge_t(label, index)
588
589 case (species_pseudo)
590 spec => pseudopotential_t(label, index)
591
592 case (species_pspio) ! a pseudopotential file to be handled by the pspio library
593 spec => pseudopotential_t(label, index)
594
595 case default
596 call messages_input_error(namespace, 'Species', "Unknown type for species '"//trim(spec%get_label())//"'", row=row, column=1)
597 end select
598
599 call spec%set_mass(-m_one)
600 call spec%set_vdw_radius(-m_one)
601 call spec%set_zval(-m_one)
602
603 call iihash_init(read_parameters)
604
605 icol = read_data
606 do
607 if (icol >= ncols) exit
608
609 call parse_block_integer(blk, row, icol, flag)
610
611 select case (flag)
612
613 case (option__species__lmax)
614 call check_duplication(option__species__lmax)
615 call parse_block_integer(blk, row, icol + 1, user_lmax)
616
617 select type(spec)
618 class is(pseudopotential_t)
619 call spec%set_user_lmax(user_lmax)
620 class default
621 call messages_input_error(namespace, 'Species', &
622 "The 'lmax' parameter in species "//trim(spec%get_label())//" can only be used with pseudopotential species", &
623 row=row, column=icol+1)
624 end select
625
626 if (user_lmax < 0) then
627 call messages_input_error(namespace, 'Species', &
628 "The 'lmax' parameter in species "//trim(spec%get_label())//" cannot be negative", &
629 row=row, column=icol+1)
630 end if
631
632 case (option__species__lloc)
633 call check_duplication(option__species__lloc)
634 call parse_block_integer(blk, row, icol + 1, user_llocal)
635
636 select type(spec)
637 class is(pseudopotential_t)
638 call spec%set_user_lloc(user_llocal)
639 class default
640 call messages_input_error(namespace, 'Species', &
641 "The 'lloc' parameter in species "//trim(spec%get_label())//" can only be used with pseudopotential species", &
642 row=row, column=icol+1)
643 end select
644
645 if (user_llocal < 0) then
646 call messages_input_error(namespace, 'Species', &
647 "The 'lloc' parameter in species "//trim(spec%get_label())//" cannot be negative", row=row, column=icol+1)
648 end if
649
650 case (option__species__hubbard_l)
651 call check_duplication(option__species__hubbard_l)
652 call parse_block_integer(blk, row, icol + 1, hubbard_l)
653
654 select type(spec)
655 class is(pseudopotential_t)
656 call spec%set_hubbard_l(hubbard_l)
657 class default
658 call messages_input_error(namespace, 'Species', &
659 "The 'hubbard_l' parameter in species "//trim(spec%get_label())//" can only be used with pseudopotential species", &
660 row=row, column=icol+1)
661 end select
662
663 if (hubbard_l < 0) then
664 call messages_input_error(namespace, 'Species', &
665 "The 'hubbard_l' parameter in species "//trim(spec%get_label())//" cannot be negative", row=row, column=icol+1)
666 end if
667
668 case (option__species__hubbard_u)
669 call check_duplication(option__species__hubbard_u)
670 call parse_block_float(blk, row, icol + 1, hubbard_u, unit = units_inp%energy)
671 call spec%set_hubbard_u(hubbard_u)
672
673 case (option__species__hubbard_alpha)
674 call check_duplication(option__species__hubbard_alpha)
675 call parse_block_float(blk, row, icol + 1, hubbard_alpha, unit = units_inp%energy)
676 call spec%set_hubbard_u(hubbard_alpha)
677
678 case (option__species__hubbard_j)
679 call check_duplication(option__species__hubbard_j)
680 call parse_block_float(blk, row, icol + 1, hubbard_j)
681 call spec%set_hubbard_u(hubbard_alpha)
682
683 if (abs(abs(spec%get_hubbard_j()-spec%get_hubbard_l())-m_half) <= m_epsilon) then
684 call messages_input_error(namespace, 'Species', "The 'hubbard_j' parameter in species "// &
685 trim(spec%get_label())//" can only be hubbard_l +/- 1/2", row=row, column=icol+1)
686 end if
687
688 case (option__species__mass)
689 call check_duplication(option__species__mass)
690 call parse_block_float(blk, row, icol + 1, mass, unit = units_inp%mass)
691 call spec%set_mass(mass)
692
693 case (option__species__valence)
694 call check_duplication(option__species__valence)
695 call parse_block_float(blk, row, icol + 1, z_val)
696 call spec%set_zval(z_val)
697 call spec%set_z(z_val)
698
699 case (option__species__jellium_radius)
700 call check_duplication(option__species__jellium_radius)
701 call parse_block_float(blk, row, icol + 1, jradius)
702 if (jradius <= m_zero) call messages_input_error(namespace, 'Species', 'jellium_radius must be positive', &
703 row=row, column=icol+1)
704
705 select type(spec)
706 type is(jellium_sphere_t)
707 call spec%set_radius(jradius)
708 class default
709 call messages_input_error(namespace, 'Species', 'jellium_radius can only be used with species_jellium', &
710 row=row, column=icol+1)
711 end select
712
713 case (option__species__gaussian_width)
714 call check_duplication(option__species__gaussian_width)
715 call parse_block_float(blk, row, icol + 1, sigma)
716 if (sigma <= m_zero) call messages_input_error(namespace, 'Species', 'gaussian_width must be positive', &
717 row=row, column=icol+1)
718 select type(spec)
719 type is(full_gaussian_t)
720 call spec%set_sigma(sigma)
721 class default
722 call messages_input_error(namespace, 'Species', 'gaussian_width can only be used with species_full_gaussian', &
723 row=row, column=icol+1)
724 end select
725
726 case (option__species__anc_a)
727 call check_duplication(option__species__anc_a)
728 call parse_block_float(blk, row, icol + 1, aa)
729 if (aa <= m_zero) call messages_input_error(namespace, 'Species', 'anc_a must be positive', &
730 row=row, column=icol+1)
731
732 select type(spec)
733 type is(full_anc_t)
734 call spec%set_a(aa)
735 class default
736 call messages_input_error(namespace, 'Species', 'anc_a can only be used with species_full_anc', &
737 row=row, column=icol+1)
738 end select
739
740 case (option__species__softening)
741 call check_duplication(option__species__softening)
742 call parse_block_float(blk, row, icol + 1, softening)
743 softening = softening**2
744
745 select type(spec)
746 type is(soft_coulomb_t)
747 call spec%set_softening(softening)
748 class default
749 call messages_input_error(namespace, 'Species', 'softening can only be used with species_soft_coulomb', &
750 row=row, column=icol+1)
751 end select
752
753 case (option__species__file)
754 call check_duplication(option__species__file)
755 call parse_block_string(blk, row, icol + 1, filename)
756 call spec%set_filename(filename)
757
758 case (option__species__db_file)
759 call messages_write("The 'db_file' option for 'Species' block is obsolete. Please use", new_line = .true.)
760 call messages_write("the option 'set' or the variable 'PseudopotentialSet' instead.")
761 call messages_fatal(namespace=namespace)
762
763 case (option__species__set)
764 call check_duplication(option__species__set)
765 call parse_block_integer(blk, row, icol + 1, pseudopotential_set_id)
766
767 select type(spec)
768 type is(pseudopotential_t)
769 spec%pseudopotential_set_initialized = .true.
770 spec%pseudopotential_set_id = pseudopotential_set_id
771 call pseudo_set_init(spec%pseudopotential_set, get_set_directory(spec%pseudopotential_set_id), ierr)
772 class default
773 call messages_input_error(namespace, 'Species', 'set can only be used with species_pseudo', &
774 row=row, column=icol+1)
775 end select
776
777 case (option__species__potential_formula)
778 call check_duplication(option__species__potential_formula)
779 select type(spec)
781 call parse_block_string(blk, row, icol + 1, spec%potential_formula)
782 call conv_to_c_string(spec%potential_formula)
783 class default
784 call messages_input_error(namespace, 'Species', 'potential_formula can only be used with species_user_defined', &
785 row=row, column=icol+1)
786 end select
787
788 case (option__species__volume)
789 call check_duplication(option__species__volume)
790
791 select type(spec)
792 type is(jellium_charge_t)
793 call parse_block_string(blk, row, icol + 1, spec%density_formula)
794 call conv_to_c_string(spec%density_formula)
795
796 class default
797 call messages_input_error(namespace, 'Species', 'volume can only be used with species_jellium_charge_density', &
798 row=row, column=icol+1)
799 end select
800
801 case (option__species__density_formula)
802 call check_duplication(option__species__density_formula)
803
804 select type(spec)
806 call parse_block_string(blk, row, icol + 1, spec%density_formula)
807 call conv_to_c_string(spec%density_formula)
808
809 class default
810 call messages_input_error(namespace, 'Species', 'density_formula can only be used with species_charge_density', &
811 row=row, column=icol+1)
812 end select
813
814 case (option__species__thickness)
815 call check_duplication(option__species__thickness)
816 call parse_block_float(blk, row, icol + 1, jthick) ! thickness of the jellium slab
817
818 if (jthick <= m_zero) then
819 call messages_input_error(namespace, 'Species', 'the value of the thickness parameter in species '&
820 //trim(spec%get_label())//' must be positive.', row=row, column=icol+1)
821 end if
822
823 select type(spec)
824 type is(jellium_slab_t)
825 call spec%set_thickness(jthick)
826 class default
827 call messages_input_error(namespace, 'Species', 'thickness can only be used with species_jellium_slab', &
828 row=row, column=icol+1)
829 end select
830
831 case (option__species__vdw_radius)
832 call check_duplication(option__species__vdw_radius)
833 call parse_block_float(blk, row, icol + 1, vdw_radius, unit = units_inp%length)
834 call spec%set_vdw_radius(vdw_radius)
835
836 case default
837 call messages_input_error(namespace, 'Species', "Unknown parameter in species '"//trim(spec%get_label())//"'", &
838 row=row, column=icol)
839
840 end select
841
842 icol = icol + 2
843 end do
844 ! CHECK THAT WHAT WE PARSED MAKES SENSE
845
846 select type(spec)
847 type is(soft_coulomb_t)
848 if (.not. parameter_defined(option__species__softening)) then
849 call messages_input_error(namespace, 'Species', &
850 "The 'softening' parameter is missing for species "//trim(spec%get_label()))
851 end if
852
854 if (.not. parameter_defined(option__species__potential_formula)) then
855 call messages_input_error(namespace, 'Species', &
856 "The 'potential_formula' parameter is missing for species '"//trim(spec%get_label())//"'")
857 end if
858
860 if (.not. parameter_defined(option__species__density_formula)) then
861 call messages_input_error(namespace, 'Species', &
862 "The 'density_formula' parameter is missing for species '"//trim(spec%get_label())//"'")
863 end if
864
865 type is(species_from_file_t)
866 if( .not. (parameter_defined(option__species__file) .or. parameter_defined(option__species__db_file))) then
867 call messages_input_error(namespace, 'Species', &
868 "The 'file' or 'db_file' parameter is missing for species '"//trim(spec%get_label())//"'")
869 end if
870
871 type is(jellium_slab_t)
872 if (.not. parameter_defined(option__species__thickness)) then
873 call messages_input_error(namespace, 'Species', &
874 "The 'thickness' parameter is missing for species '"//trim(spec%get_label())//"'")
875 end if
876
877 type is(jellium_charge_t)
878 if (.not. parameter_defined(option__species__volume)) then
879 call messages_input_error(namespace, 'Species', &
880 "The 'volume' parameter is missing for species '"//trim(spec%get_label())//"'")
881 end if
882
883 type is(pseudopotential_t)
884 if (parameter_defined(option__species__lmax) .and. parameter_defined(option__species__lloc)) then
885 if (spec%get_user_lloc() > spec%get_user_lmax()) then
886 call messages_input_error(namespace, 'Species', &
887 "the 'lloc' parameter cannot be larger than the 'lmax' parameter in species "//trim(spec%get_label()))
888 end if
889 end if
890
891 if(.not. (parameter_defined(option__species__file) .or. parameter_defined(option__species__db_file))) then
892 ! we need to read the species from the pseudopotential set
893
894 !if the set was not defined, use the default set
895 if (.not. parameter_defined(option__species__set)) then
896 spec%pseudopotential_set_id = factory%default_pseudopotential_set_id
897 spec%pseudopotential_set = factory%default_pseudopotential_set
898 end if
899
900 call read_from_set(spec, spec%pseudopotential_set_id, spec%pseudopotential_set, set_read_data)
901
902 if (set_read_data == 0) then
903 call messages_write('Species '//trim(spec%get_label())//' is not defined in the requested pseudopotential set.', &
904 new_line=.true.)
905 call messages_write('( '//trim(get_set_directory(spec%pseudopotential_set_id))//' )')
906 call messages_fatal(namespace=namespace)
907 end if
908 end if
909
910 end select
911
912 select type (spec)
913 class is(pseudopotential_t)
915
916 class is(full_anc_t)
918
919 ! If z_val was not specified, we set it to be z
920 if (spec%get_zval() < m_zero) then
921 call spec%set_zval(spec%get_z())
922 end if
923
924 class is(full_gaussian_t)
926
927 ! If z_val was not specified, we set it to be z
928 if (spec%get_zval() < m_zero) then
929 call spec%set_zval(spec%get_z())
930 end if
931
932 class is(full_delta_t)
934
935 ! If z_val was not specified, we set it to be z
936 if (spec%get_zval() < m_zero) then
937 call spec%set_zval(spec%get_z())
938 end if
939
940 class default
941 if (.not. parameter_defined(option__species__mass)) then
942 call spec%set_mass(m_one)
943 call messages_write('Info: default mass for species '//trim(spec%get_label())//':')
944 call messages_write(spec%get_mass())
945 call messages_write(' amu.')
946 call messages_info(namespace=namespace)
947 end if
948
949 if (.not. parameter_defined(option__species__vdw_radius)) then
950 call spec%set_vdw_radius(m_zero)
951 call messages_write('Info: default vdW radius for species '//trim(spec%get_label())//':')
952 call messages_write(spec%get_vdw_radius())
953 call messages_write(' [b]')
954 call messages_info(namespace=namespace)
955 end if
956
957 if (.not. parameter_defined(option__species__valence)) then
958 if (spec%is_user_defined()) then
959 call spec%set_zval(m_zero)
960 else
961 call messages_input_error(namespace, 'Species', &
962 "The 'valence' parameter is missing for species '"//trim(spec%get_label())//"'")
963 end if
964 end if
965
966 end select
967
968 call iihash_end(read_parameters)
969
970 pop_sub(read_from_block)
971
972 contains
973
974 logical function parameter_defined(param) result(defined)
975 integer(int64), intent(in) :: param
976
977 integer :: tmp
978
980
981 tmp = iihash_lookup(read_parameters, int(-param), defined)
982
984 end function parameter_defined
985
986 !------------------------------------------------------
987 subroutine check_duplication(param)
988 integer(int64), intent(in) :: param
989
991
992 if (parameter_defined(param)) then
993 call messages_input_error(namespace, 'Species', "Duplicated parameter in species '"//trim(spec%get_label())//"'")
994 end if
995
996 call iihash_insert(read_parameters, int(-param), 1)
997
999 end subroutine check_duplication
1000
1001
1002 !------------------------------------------------------
1003 subroutine check_real_atom_species()
1004
1005 call element_init(element, get_symbol(spec%get_label()))
1006
1007 if (.not. element_valid(element)) then
1008 call messages_write('Cannot determine the element for species '//trim(spec%get_label())//'.')
1009 call messages_fatal(namespace=namespace)
1010 end if
1011
1012 call spec%set_z(real(element_atomic_number(element), real64))
1013
1014 if (spec%get_mass() < m_zero) then
1015 call spec%set_mass(element_mass(element))
1016 call messages_write('Info: default mass for species '//trim(spec%get_label())//':')
1017 call messages_write(spec%get_mass())
1018 call messages_write(' amu.')
1019 call messages_info(namespace=namespace)
1020 end if
1021
1022 if (spec%get_vdw_radius() < m_zero) then
1023 call spec%set_vdw_radius(element_vdw_radius(element))
1024 if (spec%get_vdw_radius() < m_zero) then
1025 call spec%set_vdw_radius(m_zero)
1026 call messages_write("The default vdW radius for species '"//trim(spec%get_label())//"' is not defined.", &
1027 new_line = .true.)
1028 call messages_write("You can specify the vdW radius in %Species block.")
1029 call messages_warning(namespace=namespace)
1030 end if
1031 call messages_write('Info: default vdW radius for species '//trim(spec%get_label())//':')
1032 call messages_write(spec%get_vdw_radius())
1033 call messages_write(' [b]')
1034 call messages_info(namespace=namespace)
1035 end if
1036
1037 call element_end(element)
1038
1039
1040 end subroutine check_real_atom_species
1041 end function read_from_block
1042
1043end module species_factory_oct_m
1044
1045!! Local Variables:
1046!! mode: f90
1047!! coding: utf-8
1048!! End:
subroutine check_duplication(param)
logical function parameter_defined(param)
Delete the C++ object.
Definition: pseudo_set.F90:152
Nullify the C++ pointer.
Definition: pseudo_set.F90:142
logical function, public element_valid(self)
Definition: element.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_epsilon
Definition: global.F90:203
type(conf_t), public conf
Global instance of Octopus configuration.
Definition: global.F90:177
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
This module implements a simple hash table for non-negative integer keys and integer values.
Definition: iihash.F90:125
subroutine, public iihash_end(h)
Free a hash table.
Definition: iihash.F90:184
subroutine, public iihash_insert(h, key, val)
Insert a (key, val) pair into the hash table h.
Definition: iihash.F90:206
integer function, public iihash_lookup(h, key, found)
Look up a value in the hash table h. If found is present, it indicates if key could be found in the t...
Definition: iihash.F90:231
subroutine, public iihash_init(h)
Initialize a hash table h.
Definition: iihash.F90:161
Definition: io.F90:114
integer, parameter, public species_charge_density
user-defined function for charge density
Definition: jellium.F90:139
integer, parameter, public species_jellium_charge_density
jellium volume read from file
Definition: jellium.F90:139
integer, parameter, public species_jellium
jellium sphere.
Definition: jellium.F90:139
integer, parameter, public species_from_file
Definition: jellium.F90:139
integer, parameter, public species_usdef
user-defined function for local potential
Definition: jellium.F90:139
integer, parameter, public species_jellium_slab
jellium slab.
Definition: jellium.F90:139
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public pseudo_set_init(pseudo_set, dirname, ierr)
Definition: pseudo_set.F90:187
character(len=max_path_len) function, public get_set_directory(set_id)
subroutine, public read_from_set(spec, set_id, set, read_data)
Creates a pseudopotential type from a set.
integer, parameter, public species_pseudo
pseudopotential
integer, parameter, public species_pspio
pseudopotential parsed by pspio library
class(species_t) function, pointer species_factory_create_from_input(factory, namespace, label, index)
Reads the information (from the input file) about a species_t variable, initializing part of it (it h...
subroutine, public species_factory_init(factory, namespace)
class(species_t) function, pointer read_from_block(factory, namespace, blk, row, label, index, read_data)
Parses the species block for a given species.
subroutine, public species_factory_end(factory)
character(len=label_len) function, public get_symbol(label)
Definition: species.F90:518
subroutine, public conv_to_c_string(str)
converts to c string
Definition: string.F90:252
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
subroutine check_real_atom_species()
An abstract class for species. Derived classes include jellium, all electron, and pseudopotential spe...
Definition: species.F90:143
int true(void)