69 class(mesh_t),
intent(in) :: mesh
70 type(states_elec_dim_t),
intent(in) :: std
71 real(real64),
intent(in) :: rho(:,:)
72 real(real64),
intent(out) :: md(:,:)
76 select case (std%ispin)
82 md(1:mesh%np, 3) = rho(1:mesh%np, 1) - rho(1:mesh%np, 2)
85 md(1:mesh%np, 1) =
m_two*rho(1:mesh%np, 3)
86 md(1:mesh%np, 2) = -
m_two*rho(1:mesh%np, 4)
87 md(1:mesh%np, 3) = rho(1:mesh%np, 1) - rho(1:mesh%np, 2)
96 class(mesh_t),
intent(in) :: mesh
97 type(states_elec_t),
intent(in) :: st
98 real(real64),
intent(in) :: rho(:,:)
99 real(real64),
intent(out) :: mm(3)
101 real(real64),
allocatable :: md(:,:)
103 push_sub(states_elec_magnetic_moment)
105 safe_allocate(md(1:mesh%np, 1:3))
110 safe_deallocate_a(md)
112 pop_sub(states_elec_magnetic_moment)
119 type(grid_t),
intent(in) :: gr
120 type(states_elec_t),
intent(in) :: st
121 type(phase_t),
intent(in) :: phase
122 type(epot_t),
intent(in) :: ep
123 type(ions_t),
intent(in) :: ions
124 real(real64),
intent(in) :: lmm_r
125 logical,
optional,
intent(in) :: calc_orb_moments
126 integer,
optional,
intent(in) :: iunit
127 type(namespace_t),
optional,
intent(in) :: namespace
130 real(real64) :: mm(max(gr%box%dim, 3))
131 real(real64),
allocatable :: lmm(:,:)
136 safe_allocate(lmm(1:max(gr%box%dim, 3), 1:ions%natoms))
139 message(1) =
'Total Spin Magnetic Moment:'
142 write(
message(1),
'(a,f10.6)')
' mz = ', mm(3)
144 else if (st%d%ispin ==
spinors)
then
145 write(
message(1),
'(1x,3(a,f10.6,3x))')
'mx = ', mm(1),
'my = ', mm(2),
'mz = ', mm(3)
149 write(
message(1),
'(a,a,a,f7.3,a)')
'Local Spin Magnetic Moments (sphere radius [', &
153 write(
message(1),
'(a,6x,14x,a)')
' Ion',
'mz'
155 do ia = 1, ions%natoms
156 write(
message(1),
'(i4,a10,f15.6)') ia, trim(ions%atom(ia)%species%get_label()), lmm(3, ia)
159 else if (st%d%ispin ==
spinors)
then
160 write(
message(1),
'(a,8x,13x,a,13x,a,13x,a)')
' Ion',
'mx',
'my',
'mz'
162 do ia = 1, ions%natoms
163 write(
message(1),
'(i4,a10,9f15.6)') ia, trim(ions%atom(ia)%species%get_label()), lmm(:, ia)
171 if (
mpi_world%is_root())
write(iunit,
'(1x)')
175 write(
message(1),
'(a,a,a,f7.3,a)')
'Local Orbital Moments (sphere radius [', &
178 write(
message(1),
'(a,8x,13x,a,13x,a,13x,a)')
' Ion',
'Lx',
'Ly',
'Lz'
180 do ia = 1, ions%natoms
181 write(
message(1),
'(i4,a10,9f15.6)') ia, trim(ions%atom(ia)%species%get_label()), lmm(:, ia)
187 safe_deallocate_a(lmm)
194 class(
mesh_t),
intent(in) :: mesh
196 type(
ions_t),
intent(in) :: ions
198 real(real64),
intent(in) :: rho(:,:)
199 real(real64),
intent(in) :: rr
200 real(real64),
intent(out) :: lmm(max(mesh%box%dim, 3), ions%natoms)
202 integer :: ia, idir, is
203 real(real64),
allocatable :: md(:, :)
205 real(real64) :: cosqr, sinqr
206 complex(real64),
allocatable :: phase_spiral(:)
210 safe_allocate(md(1:mesh%np, 1:max(mesh%box%dim, 3)))
214 do ia = 1, ions%natoms
215 call submesh_init(sphere, ions%space, mesh, ions%latt, ions%pos(:, ia), rr)
217 if (boundaries%spiral)
then
218 safe_allocate(phase_spiral(1:sphere%np))
220 phase_spiral(is) =
exp(+
m_zi*sum((sphere%rel_x(:,is) + sphere%center - mesh%x(:, sphere%map(is))) &
221 *boundaries%spiral_q(1:mesh%box%dim)))
224 if (mesh%box%dim>= 3)
then
230 cosqr = real(phase_spiral(is), real64)
231 sinqr = aimag(phase_spiral(is))
232 lmm(1,ia) = lmm(1,ia)+md(sphere%map(is),1)*cosqr - md(sphere%map(is),2)*sinqr
233 lmm(2,ia) = lmm(2,ia)+md(sphere%map(is),1)*sinqr + md(sphere%map(is),2)*cosqr
235 lmm(1,ia) = lmm(1,ia)*mesh%volume_element
236 lmm(2,ia) = lmm(2,ia)*mesh%volume_element
239 assert(.not. boundaries%spiral)
242 safe_deallocate_a(phase_spiral)
244 do idir = 1, max(mesh%box%dim, 3)
252 call mesh%allreduce(lmm)
254 safe_deallocate_a(md)
271 type(
grid_t),
intent(in) :: gr
273 type(
phase_t),
intent(in) :: phase
274 type(
epot_t),
intent(in) :: ep
275 type(
ions_t),
intent(in) :: ions
276 real(real64),
intent(in) :: rho(:,:)
277 real(real64),
intent(in) :: rr
278 real(real64),
intent(out) :: lom(3, ions%natoms)
280 integer :: ia, idir, is, map_ip, ibatch, idim, ist, ib, ik
285 complex(real64),
allocatable :: psi(:,:), gf(:,:,:), gf_nl(:,:,:)
286 complex(real64) :: acc(3), acc_nl(3)
287 real(real64) :: x1, x2, x3, factor
292 assert(ions%space%dim == 3)
293 assert(.not. gr%use_curvilinear)
297 safe_allocate(psi(1:gr%np_part, 1:st%d%dim))
298 safe_allocate(gf(1:gr%np, 1:st%d%dim, 1:gr%der%dim))
299 safe_allocate(gf_nl(1:gr%np, 1:st%d%dim, 1:gr%der%dim))
300 safe_allocate_type_array(
wfs_elec_t, gpsib, (1:gr%der%dim))
302 do ik = st%d%kpt%start, st%d%kpt%end
303 do ib = st%group%block_start, st%group%block_end
305 call phase%copy_and_set_phase(gr, st%d%kpt, st%group%psib(ib, ik), epsib)
310 do ibatch = 1, epsib%nst
311 ist = st%group%psib(ib, ik)%ist(ibatch)
312 factor = st%kweights(ik)*st%occ(ist, ik)
324 do ia = 1, ions%natoms
325 call submesh_init(sphere, ions%space, gr, ions%latt, ions%pos(:, ia), rr)
328 call zprojector_commute_r(ep%proj(ia), gr, gr%der%boundaries, st%d%dim, 1, ik, psi, gf_nl(:, :, 1))
329 call zprojector_commute_r(ep%proj(ia), gr, gr%der%boundaries, st%d%dim, 2, ik, psi, gf_nl(:, :, 2))
330 call zprojector_commute_r(ep%proj(ia), gr, gr%der%boundaries, st%d%dim, 3, ik, psi, gf_nl(:, :, 3))
334 do idim = 1, st%d%dim
336 map_ip = sphere%map(is)
337 x1 = sphere%rel_x(1, is)
338 x2 = sphere%rel_x(2, is)
339 x3 = sphere%rel_x(3, is)
341 acc(1) = acc(1) + conjg(psi(map_ip, idim)) * (x2 * gf(map_ip, idim, 3) - x3 * gf(map_ip, idim, 2))
342 acc(2) = acc(2) + conjg(psi(map_ip, idim)) * (x3 * gf(map_ip, idim, 1) - x1 * gf(map_ip, idim, 3))
343 acc(3) = acc(3) + conjg(psi(map_ip, idim)) * (x1 * gf(map_ip, idim, 2) - x2 * gf(map_ip, idim, 1))
345 acc_nl(:) = acc_nl(:) + conjg(psi(map_ip, idim)) * gf_nl(map_ip, idim, :)
351 lom(:, ia) = lom(:, ia) + factor * gr%volume_element * (aimag(acc) +
dcross_product(ions%pos(:, ia), aimag(acc_nl)))
355 do idir = 1, gr%der%dim
356 call gpsib(idir)%end()
363 call gr%allreduce(lom)
364 if (st%parallel_in_states .or. st%d%kpt%parallel)
then
368 safe_deallocate_a(psi)
369 safe_deallocate_a(gf)
370 safe_deallocate_a(gf_nl)
371 safe_deallocate_a(gpsib)
379 class(
mesh_t),
intent(in) :: mesh
381 real(real64),
intent(in) :: qq(:)
382 complex(real64),
intent(out) :: trans_mag(6)
385 complex(real64),
allocatable :: tmp(:,:)
386 real(real64),
allocatable :: md(:, :)
387 complex(real64) :: expqr
393 safe_allocate(tmp(1:mesh%np, 1:6))
394 safe_allocate(md(1:mesh%np, 1:max(mesh%box%dim, 3)))
398 expqr =
exp(-
m_zi*sum(mesh%x(1:mesh%box%dim, ip)*qq(1:mesh%box%dim)))
399 tmp(ip,1) = expqr*md(ip,1)
400 tmp(ip,2) = expqr*md(ip,2)
401 tmp(ip,3) = expqr*md(ip,3)
402 tmp(ip,4) = conjg(expqr)*md(ip,1)
403 tmp(ip,5) = conjg(expqr)*md(ip,2)
404 tmp(ip,6) = conjg(expqr)*md(ip,3)
409 safe_deallocate_a(md)
410 safe_deallocate_a(tmp)
421 subroutine magnetic_induced(namespace, gr, st, psolver, kpoints, a_ind, b_ind)
423 type(
grid_t),
intent(in) :: gr
427 real(real64),
contiguous,
intent(out) :: a_ind(:, :)
428 real(real64),
contiguous,
intent(out) :: b_ind(:, :)
433 real(real64),
allocatable :: jj(:, :, :)
446 safe_allocate(jj(1:gr%np_part, 1:gr%der%dim, 1:st%d%nspin))
450 if (st%d%nspin > 1)
then
451 do idir = 1, gr%der%dim
452 jj(:, idir, 1) = jj(:, idir, 1) + jj(:, idir, 2)
457 do idir = 1, gr%der%dim
458 call dpoisson_solve(psolver, namespace, a_ind(:, idir), jj(:, idir, 1))
463 a_ind(1:gr%np, 1:gr%der%dim) = - a_ind(1:gr%np, 1:gr%der%dim) /
p_c
467 safe_deallocate_a(jj)
472 integer,
intent(in) :: iunit
474 real(real64),
intent(in) :: vxc(:,:)
477 real(real64),
allocatable :: torque(:,:)
478 real(real64) :: tt(3)
482 safe_allocate(torque(1:mesh%np, 1:3))
489 write(iunit,
'(a)')
'Total xc torque:'
490 write(iunit,
'(1x,3(a,es10.3,3x))')
'Tx = ', tt(1),
'Ty = ', tt(2),
'Tz = ', tt(3)
493 safe_deallocate_a(torque)
500 class(
mesh_t),
intent(in) :: mesh
501 real(real64),
intent(in) :: vxc(:,:)
503 real(real64),
intent(inout) :: torque(:,:)
505 real(real64) :: mag(3), bxc(3)
513 mag(1) =
m_two * st%rho(ip, 3)
514 mag(2) = -
m_two * st%rho(ip, 4)
515 mag(3) = st%rho(ip, 1) - st%rho(ip, 2)
517 bxc(2) =
m_two * vxc(ip, 4)
518 bxc(3) = -(vxc(ip, 1) - vxc(ip, 2))
519 torque(ip, 1) = mag(2) * bxc(3) - mag(3) * bxc(2)
520 torque(ip, 2) = mag(3) * bxc(1) - mag(1) * bxc(3)
521 torque(ip, 3) = mag(1) * bxc(2) - mag(2) * bxc(1)
double exp(double __x) __attribute__((__nothrow__
This module implements common operations on batches of mesh functions.
Module implementing boundary conditions in Octopus.
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_curl(der, ff, op_ff, ghost_update, set_bc)
apply the curl operator to a vector of mesh functions
subroutine, public zderivatives_batch_grad(der, ffb, opffb, ghost_update, set_bc, to_cartesian, factor)
apply the gradient to a batch of mesh functions
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
real(real64), parameter, public m_two
real(real64), parameter, public m_zero
complex(real64), parameter, public m_zi
real(real64), parameter, public m_epsilon
real(real64), parameter, public p_c
Electron gyromagnetic ratio, see Phys. Rev. Lett. 130, 071801 (2023)
This module implements the underlying real-space grid.
subroutine, public magnetic_local_moments(mesh, st, ions, boundaries, rho, rr, lmm)
subroutine magnetic_orbital_moments(gr, st, phase, ep, ions, rho, rr, lom)
@bief Computes orbital moments around the atoms in the GIPAW gauge
subroutine, public calc_xc_torque(mesh, vxc, st, torque)
subroutine, public magnetic_moment(mesh, st, rho, mm)
subroutine, public magnetic_total_magnetization(mesh, st, qq, trans_mag)
subroutine, public compute_and_write_magnetic_moments(gr, st, phase, ep, ions, lmm_r, calc_orb_moments, iunit, namespace)
Computes and prints the global and local magnetic moments.
subroutine, public write_total_xc_torque(iunit, mesh, vxc, st)
subroutine, public magnetic_density(mesh, std, rho, md)
subroutine, public magnetic_induced(namespace, gr, st, psolver, kpoints, a_ind, b_ind)
This subroutine receives as input a current, and produces as an output the vector potential that it i...
This module is intended to contain "only mathematical" functions and procedures.
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
type(mpi_grp_t), public mpi_world
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
subroutine, public zprojector_commute_r(pj, mesh, bnd, dim, idir, ik, psi, cpsi)
This function calculates |cpsi> += [x, V_nl] |psi>
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine, public states_elec_calc_quantities(gr, st, kpoints, nlcc, kinetic_energy_density, paramagnetic_current, density_gradient, density_laplacian, gi_kinetic_energy_density, st_end)
calculated selected quantities
real(real64) function, public dsm_integrate_frommesh(mesh, sm, ff, reduce)
subroutine, public submesh_end(this)
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
character(len=20) pure function, public units_abbrev(this)
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
This class contains information about the boundary conditions.
Description of the grid, containing information on derivatives, stencil, and symmetries.
Describes mesh distribution to nodes.
A container for the phase.
The states_elec_t class contains all electronic wave functions.
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
batches of electronic states