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
25 use debug_oct_m
27 use global_oct_m
28 use, intrinsic :: iso_fortran_env
32 use mesh_oct_m
34 use parser_oct_m
36 use space_oct_m
39
40 implicit none
41
42 private
43 public :: &
49
50 integer, public, parameter :: &
51 CONSTRAIN_NONE = 0, &
52 constrain_dir = 1, &
54
55
58 private
59 integer, public :: level = constrain_none
60 real(real64) :: lambda = m_zero
61 real(real64), allocatable :: constrain(:,:)
62 real(real64), allocatable, public :: pot(:,:)
63 real(real64), public :: energy
64
65 real(real64) :: lmm_r
67
68 ! Threshold for considering zero magnetic moments
69 real(real64), parameter :: tol_mag_norm = 1.0e-6_real64
70
71contains
72
73 ! ---------------------------------------------------------
75 subroutine magnetic_constrain_init(this, namespace, mesh, std, natoms, min_dist)
76 type(magnetic_constrain_t), intent(inout) :: this
77 type(namespace_t), intent(in) :: namespace
78 class(mesh_t), intent(in) :: mesh
79 type(states_elec_dim_t), intent(in) :: std
80 integer, intent(in) :: natoms
81 real(real64), intent(in) :: min_dist
82
83 integer :: ia, idir
84 real(real64) :: rmin, norm
85 type(block_t) :: blk
86
88
89 !%Variable MagneticConstrain
90 !%Type integer
91 !%Default no
92 !%Section Hamiltonian
93 !%Description
94 !% This variable selects which magnetic constrain expression is added to the Hamiltonian.
95 !% The <tt>AtomsMagnetDirection</tt> block is used for determining the constrained direction.
96 !%
97 !% This is ignored for independent particles, as in this case, there is nothing that
98 !% determines the magnetic direction.
99 !%
100 !%Option constrain_none 0
101 !% No constrain is added to the Hamiltonian.
102 !%Option constrain_dir 1
103 !% We are adding a constrain for the direction and sign of the magnetic moments only,
104 !% see Ma and Dudarev, PRB 91, 054420 (2015).
105 !%Option constrain_full 2
106 !% We are adding a constrain for the direction and norm of the magnetic moments only,
107 !% see Ma and Dudarev, PRB 91, 054420 (2015).
108 !%End
109 call parse_variable(namespace, 'MagneticConstrain', constrain_none, this%level)
110 call messages_print_var_value('MagneticConstrain', this%level, namespace=namespace)
111 if (this%level == constrain_none) then
113 return
114 end if
115
116 call messages_experimental('MagneticConstrain', namespace=namespace)
117
118 !%Variable MagneticConstrainStrength
119 !%Type float
120 !%Default 0.01
121 !%Section Hamiltonian
122 !%Description
123 !% This variable determines the value of the Lagrange multiplier used for the constrain term.
124 !%End
125 call parse_variable(namespace, 'MagneticConstrainStrength', 0.01_real64, this%lambda)
126 call messages_print_var_value('MagneticConstrainStrength', this%lambda, namespace=namespace)
127
128
129 !We are using the AtomsMagnetDirection block for the contrain
130 if (parse_block(namespace, 'AtomsMagnetDirection', blk) < 0) then
131 message(1) = "AtomsMagnetDirection block is not defined."
132 message(2) = "Magnetic constrained DFT cannot be used without it."
133 call messages_fatal(2, namespace=namespace)
134 end if
135
136 if (parse_block_n(blk) /= natoms) then
137 message(1) = "AtomsMagnetDirection block has the wrong number of rows."
138 call messages_fatal(1, namespace=namespace)
139 end if
140
141 if(std%ispin == unpolarized) then
142 message(1) = "Magnetic constrains can only be used for spin-polized and spinor calculations."
143 call messages_fatal(1, namespace=namespace)
144 end if
146 safe_allocate(this%constrain(1:3, natoms))
147 do ia = 1, natoms
148 !Read from AtomsMagnetDirection block
149 if (std%ispin == spin_polarized) then
150 call parse_block_float(blk, ia-1, 0, this%constrain(3, ia))
151 this%constrain(1:2, ia) = m_zero
152 elseif (std%ispin == spinors) then
153 do idir = 1, 3
154 call parse_block_float(blk, ia-1, idir-1, this%constrain(idir, ia))
155 if (abs(this%constrain(idir, ia)) < tol_mag_norm) this%constrain(idir, ia) = m_zero
156 end do
157 end if
159 if(this%level == constrain_dir) then
160 !Renormalization for constrained directions
161 norm = norm2(this%constrain(1:3, ia))
162 this%constrain(1:3, ia) = this%constrain(1:3, ia) / norm
163 end if
164 end do
165 call parse_block_end(blk)
166
167
168 safe_allocate(this%pot(1:mesh%np, 1:std%nspin))
169 this%pot = m_zero
171 rmin = min_dist
172 call parse_variable(namespace, 'LocalMagneticMomentsSphereRadius', min(m_half*rmin, 100.0_real64), &
173 this%lmm_r)
174
176 end subroutine magnetic_constrain_init
177
178 ! ---------------------------------------------------------
180 subroutine magnetic_constrain_end(this)
181 type(magnetic_constrain_t), intent(inout) :: this
182
183 push_sub(magnetic_constrain_end)
184
185 safe_deallocate_a(this%constrain)
186 safe_deallocate_a(this%pot)
187
189 end subroutine magnetic_constrain_end
190
191 ! ---------------------------------------------------------
193 subroutine magnetic_constrain_copy(this_out, this_in)
194 type(magnetic_constrain_t), intent(out) :: this_out
195 type(magnetic_constrain_t), intent(in) :: this_in
196
198
199 this_out%level = this_in%level
200
201 ! Keep a deterministic state when the source object was never initialized
202 ! through magnetic_constrain_init.
203 this_out%lambda = m_zero
204 this_out%energy = m_zero
205 this_out%lmm_r = m_zero
206 if (this_in%level /= constrain_none) then
207 this_out%lambda = this_in%lambda
208 this_out%energy = this_in%energy
209 this_out%lmm_r = this_in%lmm_r
210 else if (allocated(this_in%pot) .or. allocated(this_in%constrain)) then
211 this_out%energy = this_in%energy
212 end if
213
214 if (allocated(this_in%constrain)) then
215 safe_allocate(this_out%constrain(1:size(this_in%constrain, dim=1), 1:size(this_in%constrain, dim=2)))
216 this_out%constrain = this_in%constrain
217 end if
218
219 if (allocated(this_in%pot)) then
220 safe_allocate(this_out%pot(1:size(this_in%pot, dim=1), 1:size(this_in%pot, dim=2)))
221 this_out%pot = this_in%pot
222 end if
223
225 end subroutine magnetic_constrain_copy
226
227 ! ---------------------------------------------------------
229 subroutine magnetic_constrain_update(this, mesh, std, space, latt, pos, rho)
230 type(magnetic_constrain_t), intent(inout) :: this
231 class(mesh_t), intent(in) :: mesh
232 type(states_elec_dim_t), intent(in) :: std
233 class(space_t), intent(in) :: space
234 type(lattice_vectors_t), intent(in) :: latt
235 real(real64), intent(in) :: pos(:,:)
236 real(real64), intent(in) :: rho(:,:)
237
238 integer :: ia, idir, ip
239 real(real64) :: bb(3), b2, lmm(3), dotp, norm, xx, max_r
240 real(real64), allocatable :: md(:,:), mdf(:), mask(:)
241 type(submesh_t) :: sphere
242 integer :: natoms
243 ! Threshold for computes sinc(x) with a Taylor expansion
244 ! The error in the Taylor expansion is below 8e-11 for this value
245 real(real64), parameter :: tol_sinc = 0.01_real64
246
247 if (this%level == constrain_none) return
248
250
251 call profiling_in(tostring(magnetic_constrain))
252
253 this%pot = m_zero
254 this%energy = m_zero
255
256 safe_allocate(md(1:mesh%np, 1:max(space%dim, 3)))
257 call magnetic_density(mesh, std, rho, md)
258
259 natoms = size(pos, dim=2)
260
261 do ia = 1, natoms
262 call submesh_init(sphere, space, mesh, latt, pos(:, ia), this%lmm_r)
263
264 ! Create a mask to be applied to the magnetic moment, see
265 ! Ma and Dudarev, PRB 91, 054420 (2015)
266 ! This helps the convergence in avoiding jitter at the edge of the sphere,
267 ! and produces a continuous potential
268 max_r = maxval(sphere%r)
269 safe_allocate(mask(1:sphere%np))
270 mask = m_zero
271 !$omp parallel do private(xx)
272 do ip = 1, sphere%np
273 xx = sphere%r(ip) * m_pi / max_r
274 if(xx < tol_sinc) then
275 mask(ip) = m_one-xx*xx/6.0_real64
276 else
277 mask(ip) = sin(xx)/xx
278 end if
279 end do
280
281 ! We compute the local moment here
282 ! We multiply here by a function that decreases monotonically to zero
283 ! at the boundary of the atomic sphere.
284 lmm = m_zero
285 safe_allocate(mdf(1:sphere%np))
286 do idir = 1, max(space%dim, 3)
287 !$omp parallel do
288 do ip = 1, sphere%np
289 mdf(ip) = md(sphere%map(ip), idir) * mask(ip)
290 end do
291 lmm(idir) = dsm_integrate(mesh, sphere, mdf)
292 end do
293 safe_deallocate_a(mdf)
294
295 ! See Ma and Dudarev, PRB 91, 054420 (2015)
296 ! See for instance https://www.vasp.at/wiki/index.php/I_CONSTRAINED_M
297 select case(this%level)
298 case(constrain_dir)
299 norm = norm2(lmm)
300 if (norm > tol_mag_norm) then
301 dotp = dot_product(lmm, this%constrain(1:3, ia))
302 this%energy = this%energy + this%lambda * (norm - dotp)
303 bb = lmm/norm - this%constrain(:,ia)
304 else
305 cycle
306 end if
307 case(constrain_full)
308 this%energy = this%energy + this%lambda * sum((lmm - this%constrain(:,ia))**2)
309 bb = m_two * (lmm - this%constrain(:,ia))
310 end select
311 bb = bb * this%lambda
312
313 ! We are adding an effective Zeeman term within the atomic sphere
314 select case(std%ispin)
315 case (spin_polarized)
316 b2 = norm2(bb(1:max(space%dim, 3)))
317 do ip = 1, sphere%np
318 this%pot(sphere%map(ip), 1) = this%pot(sphere%map(ip), 1) + b2 * mask(ip)
319 this%pot(sphere%map(ip), 2) = this%pot(sphere%map(ip), 2) - b2 * mask(ip)
320 end do
321
322 case (spinors)
323 do ip = 1, sphere%np
324 this%pot(sphere%map(ip), 1) = this%pot(sphere%map(ip), 1) + bb(3) * mask(ip)
325 this%pot(sphere%map(ip), 2) = this%pot(sphere%map(ip), 2) - bb(3) * mask(ip)
326 this%pot(sphere%map(ip), 3) = this%pot(sphere%map(ip), 3) + bb(1) * mask(ip)
327 this%pot(sphere%map(ip), 4) = this%pot(sphere%map(ip), 4) - bb(2) * mask(ip)
328 end do
329
330 end select
331
332 safe_deallocate_a(mask)
333 call submesh_end(sphere)
334 end do
335
336 safe_deallocate_a(md)
337 safe_deallocate_a(mdf)
338
339 call profiling_out(tostring(magnetic_constrain))
340
342 end subroutine magnetic_constrain_update
343
345
346!! Local Variables:
347!! mode: f90
348!! coding: utf-8
349!! End:
Prints out to iunit a message in the form: ["InputVariable" = value] where "InputVariable" is given b...
Definition: messages.F90:182
double sin(double __x) __attribute__((__nothrow__
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
This modules implements the routines for doing constrain DFT for noncollinear magnetism.
subroutine, public magnetic_constrain_copy(this_out, this_in)
Deep-copy magnetic constrain data.
subroutine, public magnetic_constrain_end(this)
Releases memory of the magnetic constrain.
subroutine, public magnetic_constrain_update(this, mesh, std, space, latt, pos, rho)
Recomputes the magnetic contraining potential.
integer, parameter, public constrain_full
integer, parameter, public constrain_dir
subroutine, public magnetic_constrain_init(this, namespace, mesh, std, natoms, min_dist)
Initilializes the magnetic_constrain_t object.
subroutine, public magnetic_density(mesh, std, rho, md)
Definition: magnetic.F90:164
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:410
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1063
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
This module handles spin dimensions of the states and the k-point distribution.
real(real64) function, public dsm_integrate(mesh, sm, ff, reduce)
Definition: submesh.F90:1163
subroutine, public submesh_end(this)
Definition: submesh.F90:733
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:282
Datatype containing the magnetic constrain information.
Describes mesh distribution to nodes.
Definition: mesh.F90:187
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