Octopus
orbitalbasis.F90
Go to the documentation of this file.
1!! Copyright (C) 2016 N. Tancogne-Dejean
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, 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#include "global.h"
19
21 use accel_oct_m
24 use debug_oct_m
26 use global_oct_m
27 use io_oct_m
28 use ions_oct_m
29 use, intrinsic :: iso_fortran_env
31 use math_oct_m
33 use mesh_oct_m
38 use parser_oct_m
43 use types_oct_m
44
45 implicit none
46
47 private
48
49 public :: &
57
59 private
60 type(orbitalset_t), allocatable, public :: orbsets(:)
61
62 integer, public :: norbsets = 0
63 integer, public :: maxnorbs = 0
64 integer, public :: max_np = 0
65 integer, public :: size = 0
66 integer, allocatable, public :: global2os(:,:)
67 integer, allocatable :: os2global(:,:)
68
69 integer(int64) :: truncation
70 real(real64) :: threshold = 0.01_real64
71
72 logical, public :: normalize = .true.
73 logical, public :: use_submesh = .false.
74 logical, public :: orthogonalization = .false.
75
76 character(len=256), public :: debugdir
77 end type orbitalbasis_t
78
79contains
80
81 subroutine orbitalbasis_init(this, namespace, periodic_dim)
82 type(orbitalbasis_t), intent(inout) :: this
83 type(namespace_t), intent(in) :: namespace
84 integer, intent(in) :: periodic_dim
85
86 push_sub(orbitalbasis_init)
87
88 !%Variable AOTruncation
89 !%Type flag
90 !%Default ao_full
91 !%Section Atomic Orbitals
92 !%Description
93 !% This option determines how Octopus will truncate the orbitals used for DFT+U.
94 !% Except for the full method, the other options are only there to get a quick idea.
95 !%Option ao_full bit(0)
96 !% The full size of the orbitals used. The radius is controled by variable AOThreshold.
97 !%Option ao_box bit(1)
98 !% The radius of the orbitals are restricted to the size of the simulation box.
99 !% This reduces the number of points used to discretize the orbitals.
100 !% This is mostly a debug option, and one should be aware that changing the size of the simulation box
101 !% will affect the result of the calculation. It is recommended to use ao_nlradius instead.
102 !%Option ao_nlradius bit(2)
103 !% The radius of the orbitals are restricted to the radius of the non-local part of the pseudopotential
104 !% of the corresponding atom.
105 !%End
106 call parse_variable(namespace, 'AOTruncation', option__aotruncation__ao_full, this%truncation)
107 call messages_print_var_option('AOTruncation', this%truncation, namespace=namespace)
108
109 !%Variable AOThreshold
110 !%Type float
111 !%Default 0.01
112 !%Section Atomic Orbitals
113 !%Description
114 !% Determines the threshold used to compute the radius of the atomic orbitals for DFT+U and for Wannier90.
115 !% This radius is computed by making sure that the
116 !% absolute value of the radial part of the atomic orbital is below the specified threshold.
117 !% This value should be converged to be sure that results do not depend on this value.
118 !% However increasing this value increases the number of grid points covered by the orbitals and directly affect performances.
119 !%End
120 call parse_variable(namespace, 'AOThreshold', 0.01_real64, this%threshold)
121 if (this%threshold <= m_zero) call messages_input_error(namespace, 'AOThreshold')
122 call messages_print_var_value('AOThreshold', this%threshold, namespace=namespace)
123
124 !%Variable AONormalize
125 !%Type logical
126 !%Default yes
127 !%Section Atomic Orbitals
128 !%Description
129 !% If set to yes, Octopus will normalize the atomic orbitals individually.
130 !% This variable is ignored is <tt>AOLoewdin</tt> is set to yes.
131 !%End
132 call parse_variable(namespace, 'AONormalize', .true., this%normalize)
133 call messages_print_var_value('AONormalize', this%normalize, namespace=namespace)
134
135 !%Variable AOSubmesh
136 !%Type logical
137 !%Section Atomic Orbitals
138 !%Description
139 !% If set to yes, Octopus will use submeshes to internally store the orbitals with
140 !% their phase instead of storing them on the mesh. This is usually slower for small
141 !% periodic systems, but becomes advantageous for large supercells.
142 !% Submeshes are not compatible with Loewdin orthogonalization.
143 !% For periodic systems, the default is set to no, whereas it is set to yes for isolated systems.
144 !%End
145 if (periodic_dim > 0) then
146 call parse_variable(namespace, 'AOSubmesh', .false., this%use_submesh)
147 else
148 call parse_variable(namespace, 'AOSubmesh', .true., this%use_submesh)
149 end if
150 call messages_print_var_value('AOSubmesh', this%use_submesh, namespace=namespace)
152 !%Variable AOLoewdin
153 !%Type logical
154 !%Default no
155 !%Section Atomic Orbitals
156 !%Description
157 !% This option determines if the atomic orbital basis is orthogonalized or not.
158 !% This is done for using the Loewdin orthogonalization scheme.
159 !% The default is set to no for the moment as this option is
160 !% not yet implemented for isolated systems, and seems to lead to important egg-box effect
161 !%End
162 call parse_variable(namespace, 'AOLoewdin', .false., this%orthogonalization)
163 call messages_print_var_value('AOLoewdin', this%orthogonalization, namespace=namespace)
164 if (this%orthogonalization) call messages_experimental("AOLoewdin")
165 if (this%orthogonalization) this%normalize = .false.
167 if (this%orthogonalization .and. this%use_submesh) then
168 call messages_not_implemented("AOLoewdin=yes with AOSubmesh=yes.")
169 end if
170
171 if (debug%info) then
172 write(this%debugdir, '(a)') 'debug/ao_basis'
173 call io_mkdir(this%debugdir, namespace)
174 end if
175
176 pop_sub(orbitalbasis_init)
177 end subroutine orbitalbasis_init
178
179
180 subroutine orbitalbasis_end(this)
181 type(orbitalbasis_t), intent(inout) :: this
182
183 integer :: ios
184
185 push_sub(orbitalbasis_end)
186
187 do ios = 1, this%norbsets
188 call orbitalset_end(this%orbsets(ios))
189 end do
190 safe_deallocate_a(this%orbsets)
191
192 safe_deallocate_a(this%global2os)
193 safe_deallocate_a(this%os2global)
194
195 pop_sub(orbitalbasis_end)
196 end subroutine orbitalbasis_end
197
198#include "undef.F90"
199#include "real.F90"
200#include "orbitalbasis_inc.F90"
201
202#include "undef.F90"
203#include "complex.F90"
204#include "orbitalbasis_inc.F90"
205
206end module orbitalbasis_oct_m
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
type(debug_t), save, public debug
Definition: debug.F90:151
real(real64), parameter, public m_zero
Definition: global.F90:187
Definition: io.F90:114
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:354
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
subroutine, public orbitalbasis_end(this)
subroutine, public zorbitalbasis_build(this, namespace, ions, mesh, kpt, ndim, skip_s_orb, use_all_orb, verbose)
This routine is an interface for constructing the orbital basis.
subroutine, public dorbitalbasis_build(this, namespace, ions, mesh, kpt, ndim, skip_s_orb, use_all_orb, verbose)
This routine is an interface for constructing the orbital basis.
subroutine, public zorbitalbasis_build_empty(this, mesh, kpt, ndim, norbsets, map_os, verbose)
This routine constructd an empty orbital basis.
subroutine, public dorbitalbasis_build_empty(this, mesh, kpt, ndim, norbsets, map_os, verbose)
This routine constructd an empty orbital basis.
subroutine, public orbitalbasis_init(this, namespace, periodic_dim)
subroutine, public orbitalset_end(this)
Definition: orbitalset.F90:225
int true(void)