Octopus
magnetic_constrain.F90
Go to the documentation of this file.
1!! Copyright (C) 2021 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
19#include "global.h"
20
22 use debug_oct_m
24 use global_oct_m
25 use, intrinsic :: iso_fortran_env
29 use mesh_oct_m
31 use parser_oct_m
33 use space_oct_m
36
37 implicit none
38
39 private
40 public :: &
45
46 integer, public, parameter :: &
47 CONSTRAIN_NONE = 0, &
48 constrain_dir = 1, &
50
51
53 private
54 integer, public :: level
55 real(real64) :: lambda
56 real(real64), allocatable :: constrain(:,:)
57 real(real64), allocatable, public :: pot(:,:)
58 real(real64), public :: energy
59
60 real(real64) :: lmm_r
61
62 real(real64), pointer :: rho(:,:)
64
65contains
66
67 subroutine magnetic_constrain_init(this, namespace, mesh, std, rho, natoms, min_dist)
68 type(magnetic_constrain_t), intent(inout) :: this
69 type(namespace_t), intent(in) :: namespace
70 class(mesh_t), intent(in) :: mesh
71 type(states_elec_dim_t), intent(in) :: std
72 real(real64), target, intent(in) :: rho(:,:)
73 integer, intent(in) :: natoms
74 real(real64), intent(in) :: min_dist
75
76 integer :: ia, idir
77 real(real64) :: rmin, norm
78 type(block_t) :: blk
79
81
82 !%Variable MagneticConstrain
83 !%Type integer
84 !%Default no
85 !%Section Hamiltonian
86 !%Description
87 !% This variable selects which magnetic constrain expression is added to the Hamiltonian.
88 !%Option constrain_none 0
89 !% No constrain is added to the Hamiltonian.
90 !%Option constrain_dir 1
91 !% We are adding a constrain for the direction of the magnetic moments only,
92 !% see PRB 91, 054420 (2015).
93 !%Option constrain_full 2
94 !% We are adding a constrain for the direction and norm of the magnetic moments only,
95 !% see PRB 91, 054420 (2015).
96 !%End
97 call parse_variable(namespace, 'MagneticConstrain', constrain_none, this%level)
98 call messages_print_var_value('MagneticConstrain', this%level, namespace=namespace)
99 if (this%level == constrain_none) then
101 return
102 end if
103
104 call messages_experimental('MagneticConstrain', namespace=namespace)
105
106 !%Variable MagneticConstrainStrength
107 !%Type float
108 !%Default 0.01
109 !%Section Hamiltonian
110 !%Description
111 !% This variable determines the value of the Lagrange multiplier used for the constrain term.
112 !%End
113 call parse_variable(namespace, 'MagneticConstrainStrength', 0.01_real64, this%lambda)
114 call messages_print_var_value('MagneticConstrainStrength', this%lambda, namespace=namespace)
115
116
117 !We are using the AtomsMagnetDirection block for the contrain
118 if (parse_block(namespace, 'AtomsMagnetDirection', blk) < 0) then
119 message(1) = "AtomsMagnetDirection block is not defined."
120 message(2) = "Magnetic constrained DFT cannot be used without it."
121 call messages_fatal(2, namespace=namespace)
122 end if
123
124 if (parse_block_n(blk) /= natoms) then
125 message(1) = "AtomsMagnetDirection block has the wrong number of rows."
126 call messages_fatal(1, namespace=namespace)
127 end if
128
129 if(std%ispin == unpolarized) then
130 message(1) = "Magnetic constrains can only be used for spin-polized and spinor calculations."
131 call messages_fatal(1, namespace=namespace)
132 end if
133
134 safe_allocate(this%constrain(1:3, natoms))
135 do ia = 1, natoms
136 !Read from AtomsMagnetDirection block
137 if (std%ispin == spin_polarized) then
138 call parse_block_float(blk, ia-1, 0, this%constrain(3, ia))
139 this%constrain(1:2, ia) = m_zero
140 elseif (std%ispin == spinors) then
141 do idir = 1, 3
142 call parse_block_float(blk, ia-1, idir-1, this%constrain(idir, ia))
143 if (abs(this%constrain(idir, ia)) < 1.0e-20_real64) this%constrain(idir, ia) = m_zero
144 end do
145 end if
146
147 if(this%level == constrain_dir) then
148 !Renormalization for constrained directions
149 norm = norm2(this%constrain(1:3, ia))
150 this%constrain(1:3, ia) = this%constrain(1:3, ia) / norm
151 end if
152 end do
154
155 safe_allocate(this%pot(1:mesh%np, std%nspin))
156 this%pot = m_zero
157
158 rmin = min_dist
159 call parse_variable(namespace, 'LocalMagneticMomentsSphereRadius', min(m_half*rmin, 100.0_real64), &
160 this%lmm_r)
161
162 this%rho => rho
163
165 end subroutine magnetic_constrain_init
166
167 subroutine magnetic_constrain_end(this)
168 type(magnetic_constrain_t), intent(inout) :: this
169
170 push_sub(magnetic_constrain_end)
171
172 safe_deallocate_a(this%constrain)
173 safe_deallocate_a(this%pot)
174 nullify(this%rho)
175
177 end subroutine magnetic_constrain_end
178
179! ---------------------------------------------------------
180 subroutine magnetic_constrain_update(this, mesh, std, space, latt, pos)
181 type(magnetic_constrain_t), intent(inout) :: this
182 class(mesh_t), intent(in) :: mesh
183 type(states_elec_dim_t), intent(in) :: std
184 class(space_t), intent(in) :: space
185 type(lattice_vectors_t), intent(in) :: latt
186 real(real64), intent(in) :: pos(:,:)
187
188 integer :: ia, idir, ip
189 real(real64) :: bb(3), b2, lmm(3), dotp, norm, xx
190 real(real64), allocatable :: md(:,:), mdf(:), mask(:)
191 type(submesh_t) :: sphere
192 integer :: natoms
193
194 if (this%level == constrain_none) return
195
197
198 call profiling_in(tostring(magnetic_constrain))
199
200 this%pot = m_zero
201 this%energy = m_zero
202
203 safe_allocate(md(1:mesh%np, 1:max(space%dim, 3)))
204 safe_allocate(mdf(1:mesh%np))
205 call magnetic_density(mesh, std, this%rho, md)
206
207 natoms = size(pos, dim=2)
208
209 do ia = 1, natoms
210 call submesh_init(sphere, space, mesh, latt, pos(:, ia), this%lmm_r)
211
212 safe_allocate(mask(1:sphere%np))
213 mask = m_zero
214 do ip = 1, sphere%np
215 xx = sphere%r(ip) * m_pi / this%lmm_r
216 !if(xx < 0.1) then
217 ! mask(ip) = M_ONE-xx*xx/6.0_real64
218 !else
219 ! mask(ip) = sin(xx)/xx
220 !end if
221 mask(ip) = m_one
222 end do
223
224 ! We compute the local moment here
225 ! We multiply here by a function that decreases monotonically to zero
226 ! at the boundary of the atomic sphere.
227 mdf = m_zero
228 lmm = m_zero
229
230 do idir = 1, max(space%dim, 3)
231 do ip = 1, sphere%np
232 mdf(sphere%map(ip)) = md(sphere%map(ip), idir) * mask(ip)
233 end do
234 lmm(idir) = dsm_integrate_frommesh(mesh, sphere, mdf)
235 end do
236
237 ! See PRB 91, 054420 (2015)
238 ! See for instance https://www.vasp.at/wiki/index.php/I_CONSTRAINED_M
239 select case(this%level)
240 case(constrain_dir)
241 norm = max(norm2(lmm),1e-10_real64)
242 dotp = dot_product(lmm, this%constrain(1:3, ia))
243 this%energy = this%energy + this%lambda * (norm - dotp)
244 bb = lmm/norm - this%constrain(:,ia)
245 case(constrain_full)
246 this%energy = this%energy + this%lambda * sum((lmm - this%constrain(:,ia))**2)
247 bb = m_two * (lmm - this%constrain(:,ia))
248 end select
249
250 ! Here we need to be carefull because we might end up in the case where
251 ! the constrained is fullfilled in one direction, say along z
252 ! and small noise in the x-y plane would then break the symmetry and
253 ! if we have independent particles, this determines the magnetization direction
254 !
255 ! To prevent problems with the indenpendent particles,
256 ! we just set a small fraction of the constrain instead of zero, to break
257 ! the symmetries
258 do idir = 1, max(space%dim, 3)
259 if(abs(bb(idir)) < 1.e-3_real64) then
260 bb(idir) = -this%constrain(idir, ia) * 1.e-2_real64
261 else
262 bb(idir) = bb(idir) * this%lambda
263 end if
264 end do
265
266 ! We are adding an effective Zeeman term within the atomic sphere
267 select case(std%ispin)
268 case (spin_polarized)
269 b2 = norm2(bb(1:max(space%dim, 3)))
270 do ip = 1, sphere%np
271 this%pot(sphere%map(ip), 1) = this%pot(sphere%map(ip), 1) + b2 * mask(ip)
272 this%pot(sphere%map(ip), 2) = this%pot(sphere%map(ip), 2) - b2 * mask(ip)
273 end do
274
275 case (spinors)
276 do ip = 1, sphere%np
277 this%pot(sphere%map(ip), 1) = this%pot(sphere%map(ip), 1) + bb(3) * mask(ip)
278 this%pot(sphere%map(ip), 2) = this%pot(sphere%map(ip), 2) - bb(3) * mask(ip)
279 this%pot(sphere%map(ip), 3) = this%pot(sphere%map(ip), 3) + bb(1) * mask(ip)
280 this%pot(sphere%map(ip), 4) = this%pot(sphere%map(ip), 4) - bb(2) * mask(ip)
281 end do
282
283 end select
284
285 safe_deallocate_a(mask)
286 call submesh_end(sphere)
287 end do
288
289 safe_deallocate_a(md)
290 safe_deallocate_a(mdf)
291
292 call profiling_out(tostring(magnetic_constrain))
293
295 end subroutine magnetic_constrain_update
296
298
299!! Local Variables:
300!! mode: f90
301!! coding: utf-8
302!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:180
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
subroutine, public magnetic_constrain_end(this)
integer, parameter, public constrain_full
integer, parameter, public constrain_dir
subroutine, public magnetic_constrain_update(this, mesh, std, space, latt, pos)
subroutine, public magnetic_constrain_init(this, namespace, mesh, std, rho, natoms, min_dist)
subroutine, public magnetic_density(mesh, std, rho, md)
Definition: magnetic.F90:156
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
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_experimental(name, namespace)
Definition: messages.F90:1097
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
This module handles spin dimensions of the states and the k-point distribution.
real(real64) function, public dsm_integrate_frommesh(mesh, sm, ff, reduce)
Definition: submesh.F90:1203
subroutine, public submesh_end(this)
Definition: submesh.F90:737
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:282
Describes mesh distribution to nodes.
Definition: mesh.F90:186
class for organizing spins and k-points
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
Definition: submesh.F90:177