Octopus
target_spin.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
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
23 use global_oct_m
24 use grid_oct_m
25 use, intrinsic :: iso_fortran_env
29 use parser_oct_m
33
34 implicit none
35
36 private
37 public :: &
41
42
43contains
44
45 ! ----------------------------------------------------------------------
47 subroutine target_init_spin(tg, namespace)
48 type(target_t), intent(inout) :: tg
49 type(namespace_t), intent(in) :: namespace
50
51
52 type(block_t) :: blk
53 integer :: jst
54 complex(real64) :: alpha(3)
55
56 push_sub(target_init_spin)
57
58
59 message(1) = 'Info: Using a spin target'
60 call messages_info(1, namespace=namespace)
61
62 !%Variable OCTTargetSpin
63 !%Type block
64 !%Section Calculation Modes::Optimal Control
65 !%Description
66 !% (Experimental) Specify the targeted spin as a 3-component vector. It will be normalized.
67 !%End
68 if (parse_is_defined(namespace, 'OCTTargetSpin')) then
69 if (parse_block(namespace, 'OCTTargetSpin', blk) == 0) then
70 alpha = m_z0
71 do jst = 1, parse_block_cols(blk, 0)
72 call parse_block_cmplx(blk, 0, jst - 1, alpha(jst))
73 end do
74 call parse_block_end(blk)
75
76 alpha = alpha/norm2(abs(alpha))
77
78 tg%spin_matrix(1, 1) = alpha(3)
79 tg%spin_matrix(2, 2) = -alpha(3)
80 tg%spin_matrix(1, 2) = alpha(1) - m_zi * alpha(2)
81 tg%spin_matrix(2, 1) = alpha(1) + m_zi * alpha(2)
82
83 else
84 call messages_variable_is_block(namespace, 'OCTTargetSpin')
85 end if
86
87 else
88 message(1) = 'If "OCTTargetOperator = oct_tg_spin", then you must'
89 message(2) = 'supply a "OCTTargetSpin" block.'
90 call messages_fatal(2, namespace=namespace)
91 end if
92
93
94 pop_sub(target_init_spin)
95 end subroutine target_init_spin
96
97
98
99 ! ----------------------------------------------------------------------
101 real(real64) function target_j1_spin(tg, gr, psi) result(j1)
102 type(target_t), intent(in) :: tg
103 type(grid_t), intent(in) :: gr
104 type(states_elec_t), intent(in) :: psi
105
106 integer :: i, j
107 complex(real64), allocatable :: zpsi(:, :)
108
109 push_sub(target_j1_spin)
110
111 safe_allocate(zpsi(1:gr%np, 1:tg%st%d%dim))
112
113 call states_elec_get_state(psi, gr, 1, 1, zpsi)
115 j1 = m_zero
116 do i = 1, 2
117 do j = 1, 2
118 j1 = j1 + real(tg%spin_matrix(i,j)*zmf_dotp(gr, zpsi(:, i), zpsi(:, j)))
119 end do
120 end do
121
122 safe_deallocate_a(zpsi)
123
124 pop_sub(target_j1_spin)
125 end function target_j1_spin
126
127
128
129 ! ----------------------------------------------------------------------
131 subroutine target_chi_spin(tg, gr, psi_in, chi_out)
132 type(target_t), intent(in) :: tg
133 type(grid_t), intent(in) :: gr
134 type(states_elec_t), intent(in) :: psi_in
135 type(states_elec_t), intent(inout) :: chi_out
136
137 integer :: i, j
138 complex(real64), allocatable :: zpsi(:, :), zchi(:, :)
139
141
142 safe_allocate(zpsi(1:gr%np, 1:tg%st%d%dim))
143 safe_allocate(zchi(1:gr%np, 1:tg%st%d%dim))
144
145 call states_elec_get_state(psi_in, gr, 1, 1, zpsi)
146
147 zchi(1:gr%np, 1:tg%st%d%dim) = 0.0_real64
148
149 do i = 1, 2
150 do j = 1, 2
151 zchi(1:gr%np, i) = zchi(1:gr%np, i) + tg%spin_matrix(i, j)*zpsi(1:gr%np, j)
152 end do
153 end do
154
155 call states_elec_set_state(chi_out, gr, 1, 1, zchi)
156
157 safe_deallocate_a(zpsi)
158 safe_deallocate_a(zchi)
159
160 pop_sub(target_chi_spin)
161 end subroutine target_chi_spin
162
163end module target_spin_oct_m
164
165!! Local Variables:
166!! mode: f90
167!! coding: utf-8
168!! End:
real(real64), parameter, public m_zero
Definition: global.F90:188
complex(real64), parameter, public m_z0
Definition: global.F90:198
complex(real64), parameter, public m_zi
Definition: global.F90:202
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines various routines, operating on mesh functions.
subroutine, public messages_variable_is_block(namespace, name)
Definition: messages.F90:1069
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:414
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public target_init_spin(tg, namespace)
subroutine, public target_chi_spin(tg, gr, psi_in, chi_out)
real(real64) function, public target_j1_spin(tg, gr, psi)