Octopus
target_local.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
20#include "global.h"
21
23 use debug_oct_m
24 use global_oct_m
25 use grid_oct_m
26 use ions_oct_m
28 use io_oct_m
30 use, intrinsic :: iso_fortran_env
31 use mesh_oct_m
36 use parser_oct_m
38 use space_oct_m
40 use string_oct_m
42 use td_oct_m
43 use unit_oct_m
45
46
47 implicit none
48
49 private
50 public :: &
56
57
58contains
59
60
61 ! ----------------------------------------------------------------------
63 subroutine target_init_local(gr, namespace, tg, td)
64 type(grid_t), intent(in) :: gr
65 type(namespace_t), intent(in) :: namespace
66 type(target_t), intent(inout) :: tg
67 type(td_t), intent(in) :: td
68
69 integer :: ip
70 real(real64) :: xx(1:gr%box%dim), rr, psi_re, psi_im
71 character(len=1024) :: expression
72 push_sub(target_init_local)
73
74 tg%move_ions = td%ions_dyn%ions_move()
75 tg%dt = td%dt
76
77 !%Variable OCTLocalTarget
78 !%Type string
79 !%Section Calculation Modes::Optimal Control
80 !%Description
81 !% If <tt>OCTTargetOperator = oct_tg_local</tt>, then one must supply a function
82 !% that defines the target. This should be done by defining it through a string, using
83 !% the variable <tt>OCTLocalTarget</tt>.
84 !%End
85 if (parse_is_defined(namespace, 'OCTLocalTarget')) then
86 safe_allocate(tg%rho(1:gr%np))
87 tg%rho = m_zero
88 call parse_variable(namespace, 'OCTLocalTarget', "0", expression)
89 call conv_to_c_string(expression)
90 do ip = 1, gr%np
91 call mesh_r(gr, ip, rr, coords = xx)
92 ! parse user-defined expression
93 call parse_expression(psi_re, psi_im, gr%box%dim, xx, rr, m_zero, expression)
94 tg%rho(ip) = psi_re
95 end do
96 else
97 message(1) = 'If OCTTargetOperator = oct_tg_local, then you must give the shape'
98 message(2) = 'of this target in variable "OCTLocalTarget".'
99 call messages_fatal(2, namespace=namespace)
100 end if
101
102 pop_sub(target_init_local)
103 end subroutine target_init_local
104
105
106 ! ----------------------------------------------------------------------
108 subroutine target_end_local(tg)
109 type(target_t), intent(inout) :: tg
110
111 push_sub(target_end_local)
112
113 safe_deallocate_a(tg%rho)
114
116 end subroutine target_end_local
117
118
119 ! ----------------------------------------------------------------------
120 subroutine target_output_local(tg, namespace, space, mesh, dir, ions, outp)
121 type(target_t), intent(in) :: tg
122 type(namespace_t), intent(in) :: namespace
123 class(space_t), intent(in) :: space
124 class(mesh_t), intent(in) :: mesh
125 character(len=*), intent(in) :: dir
126 type(ions_t), intent(in) :: ions
127 type(output_t), intent(in) :: outp
128
129 integer :: ierr
130 push_sub(target_output_local)
131
132 call io_mkdir(trim(dir), namespace)
133 call dio_function_output(outp%how(0), trim(dir), 'local_target', namespace, space, mesh, &
134 tg%rho, units_out%length**(-space%dim), ierr, pos=ions%pos, atoms=ions%atom)
135
136
137 pop_sub(target_output_local)
138 end subroutine target_output_local
139 ! ----------------------------------------------------------------------
140
141
142 ! ----------------------------------------------------------------------
144 real(real64) function target_j1_local(mesh, tg, psi) result(j1)
145 class(mesh_t), intent(in) :: mesh
146 type(target_t), intent(in) :: tg
147 type(states_elec_t), intent(in) :: psi
148
149 integer :: is
150 push_sub(target_j1_local)
151
152 j1 = m_zero
153 do is = 1, psi%d%spin_channels
154 j1 = j1 + dmf_dotp(mesh, tg%rho, psi%rho(:, is))
155 end do
157 pop_sub(target_j1_local)
158 end function target_j1_local
159
160
161 ! ----------------------------------------------------------------------
163 subroutine target_chi_local(tg, mesh, psi_in, chi_out)
164 type(target_t), intent(in) :: tg
165 class(mesh_t), intent(in) :: mesh
166 type(states_elec_t), intent(in) :: psi_in
167 type(states_elec_t), intent(inout) :: chi_out
168
169 integer :: ik, idim, ist, ip
170 complex(real64), allocatable :: zpsi(:, :)
171
172 push_sub(target_chi_local)
173
174 safe_allocate(zpsi(1:mesh%np, 1:psi_in%d%dim))
175
176 do ik = 1, psi_in%nik
177 do idim = 1, psi_in%d%dim
178 do ist = psi_in%st_start, psi_in%st_end
179 call states_elec_get_state(psi_in, mesh, ist, ik, zpsi)
180 do ip = 1, mesh%np
181 zpsi(ip, idim) = psi_in%occ(ist, ik)*tg%rho(ip)*zpsi(ip, idim)
182 end do
183 call states_elec_set_state(chi_out, mesh, ist, ik, zpsi)
184 end do
185 end do
186 end do
187
188 safe_deallocate_a(zpsi)
189
190 pop_sub(target_chi_local)
191 end subroutine target_chi_local
192
193end module target_local_oct_m
194
195!! Local Variables:
196!! mode: f90
197!! coding: utf-8
198!! End:
real(real64), parameter, public m_zero
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
Definition: io.F90:114
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
pure subroutine, public mesh_r(mesh, ip, rr, origin, coords)
return the distance to the origin for a given grid point
Definition: mesh.F90:336
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
this module contains the low-level part of the output system
Definition: output_low.F90:115
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:502
subroutine, public conv_to_c_string(str)
converts to c string
Definition: string.F90:252
real(real64) function, public target_j1_local(mesh, tg, psi)
subroutine, public target_init_local(gr, namespace, tg, td)
subroutine, public target_output_local(tg, namespace, space, mesh, dir, ions, outp)
subroutine, public target_chi_local(tg, mesh, psi_in, chi_out)
subroutine, public target_end_local(tg)
Definition: td.F90:114
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out