Octopus
target_userdefined.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
24 use global_oct_m
25 use grid_oct_m
27 use io_oct_m
28 use ions_oct_m
30 use, intrinsic :: iso_fortran_env
34 use output_oct_m
36 use parser_oct_m
39 use td_oct_m
40 use space_oct_m
42 use string_oct_m
43
44 implicit none
45
46 private
47 public :: &
52
53
54contains
55
56 ! ----------------------------------------------------------------------
58 subroutine target_init_userdefined(gr, namespace, tg, td)
59 type(grid_t), intent(in) :: gr
60 type(namespace_t), intent(in) :: namespace
61 type(target_t), intent(inout) :: tg
62 type(td_t), intent(in) :: td
63
64 integer :: no_states, ib, ip, idim, inst, inik, id, ist, ik
65 type(block_t) :: blk
66 real(real64) :: xx(1:gr%box%dim), rr, psi_re, psi_im
67 complex(real64), allocatable :: zpsi(:, :)
68
70
71 message(1) = 'Info: Target is a user-defined state.'
72 call messages_info(1, namespace=namespace)
73
74 tg%move_ions = ion_dynamics_ions_move(td%ions_dyn)
75 tg%dt = td%dt
76
77 safe_allocate(zpsi(gr%np, 1:tg%st%d%dim))
78
79 !%Variable OCTTargetUserdefined
80 !%Type block
81 !%Section Calculation Modes::Optimal Control
82 !%Description
83 !% Define a target state. Syntax follows the one of the <tt>UserDefinedStates</tt> block.
84 !% Example:
85 !%
86 !% <tt>%OCTTargetUserdefined
87 !% <br>&nbsp;&nbsp; 1 | 1 | 1 | "exp(-r^2)*exp(-i*0.2*x)"
88 !% <br>%</tt>
89 !%
90 !%End
91 if (parse_block(namespace, 'OCTTargetUserdefined', blk) == 0) then
92
93 no_states = parse_block_n(blk)
94 do ib = 1, no_states
95 call parse_block_integer(blk, ib - 1, 0, idim)
96 call parse_block_integer(blk, ib - 1, 1, inst)
97 call parse_block_integer(blk, ib - 1, 2, inik)
98
99 ! read formula strings and convert to C strings
100 do id = 1, tg%st%d%dim
101 do ist = 1, tg%st%nst
102 do ik = 1, tg%st%nik
103
104 ! does the block entry match and is this node responsible?
105 if (.not. (id == idim .and. ist == inst .and. ik == inik &
106 .and. tg%st%st_start <= ist .and. tg%st%st_end >= ist)) cycle
107
108 ! parse formula string
109 call parse_block_string( &
110 blk, ib - 1, 3, tg%st%user_def_states(id, ist, ik))
111 ! convert to C string
112 call conv_to_c_string(tg%st%user_def_states(id, ist, ik))
113
114 do ip = 1, gr%np
115 xx = gr%x(ip, :)
116 rr = norm2(xx)
117
118 ! parse user-defined expressions
119 call parse_expression(psi_re, psi_im, &
120 gr%box%dim, xx, rr, m_zero, tg%st%user_def_states(id, ist, ik))
121 ! fill state
122 zpsi(ip, id) = cmplx(psi_re, psi_im, real64)
123 end do
124
125 ! normalize orbital
126 call zmf_normalize(gr, tg%st%d%dim, zpsi)
127
128 call states_elec_set_state(tg%st, gr, ist, ik, zpsi)
129
130 end do
131 end do
132 end do
133 end do
134 call parse_block_end(blk)
135 call density_calc(tg%st, gr, tg%st%rho)
136 else
137 call messages_variable_is_block(namespace, 'OCTTargetUserdefined')
138 end if
139
140 safe_deallocate_a(zpsi)
141
143 end subroutine target_init_userdefined
144
145
146 ! ----------------------------------------------------------------------
147 subroutine target_output_userdefined(tg, namespace, space, gr, dir, ions, hm, outp)
148 type(target_t), intent(in) :: tg
149 type(namespace_t), intent(in) :: namespace
150 class(space_t), intent(in) :: space
151 type(grid_t), intent(in) :: gr
152 character(len=*), intent(in) :: dir
153 type(ions_t), intent(in) :: ions
154 type(hamiltonian_elec_t), intent(in) :: hm
155 type(output_t), intent(in) :: outp
157
158 call io_mkdir(trim(dir), namespace)
159 call output_states(outp, namespace, space, trim(dir), tg%st, gr, ions, hm, -1)
160
162 end subroutine target_output_userdefined
163 ! ----------------------------------------------------------------------
164
165
166 ! ----------------------------------------------------------------------
168 real(real64) function target_j1_userdefined(tg, gr, psi) result(j1)
169 type(target_t), intent(in) :: tg
170 type(grid_t), intent(in) :: gr
171 type(states_elec_t), intent(in) :: psi
172
173 integer :: ik, ist
174 complex(real64), allocatable :: zpsi(:, :), zst(:, :)
175
176 push_sub(target_j1_userdefined)
177
178 safe_allocate(zpsi(1:gr%np, 1:tg%st%d%dim))
179 safe_allocate(zst(1:gr%np, 1:tg%st%d%dim))
180
181 j1 = m_zero
182 do ik = 1, psi%nik
183 do ist = psi%st_start, psi%st_end
184
185 call states_elec_get_state(psi, gr, ist, ik, zpsi)
186 call states_elec_get_state(tg%st, gr, ist, ik, zst)
187
188 j1 = j1 + psi%occ(ist, ik)*abs(zmf_dotp(gr, psi%d%dim, zpsi, zst))**2
189 end do
190 end do
191
192 safe_deallocate_a(zpsi)
193 safe_deallocate_a(zst)
194
195 pop_sub(target_j1_userdefined)
196 end function target_j1_userdefined
197
198
199 ! ----------------------------------------------------------------------
201 subroutine target_chi_userdefined(tg, gr, psi_in, chi_out)
202 type(target_t), intent(in) :: tg
203 type(grid_t), intent(in) :: gr
204 type(states_elec_t), intent(in) :: psi_in
205 type(states_elec_t), intent(inout) :: chi_out
206
207 integer :: ik, ist
208 complex(real64) :: olap
209 complex(real64), allocatable :: zpsi(:, :), zst(:, :), zchi(:, :)
210
211 push_sub(target_chi_userdefined)
212
213 safe_allocate(zpsi(1:gr%np, 1:tg%st%d%dim))
214 safe_allocate(zst(1:gr%np, 1:tg%st%d%dim))
215 safe_allocate(zchi(1:gr%np, 1:tg%st%d%dim))
216
217 do ik = 1, psi_in%nik
218 do ist = psi_in%st_start, psi_in%st_end
219
220 call states_elec_get_state(psi_in, gr, ist, ik, zpsi)
221 call states_elec_get_state(tg%st, gr, ist, ik, zst)
222
223 olap = zmf_dotp(gr, zst(:, 1), zpsi(:, 1))
224 zchi(1:gr%np, 1:tg%st%d%dim) = olap*zst(1:gr%np, 1:tg%st%d%dim)
225
226 call states_elec_set_state(chi_out, gr, ist, ik, zchi)
227
228 end do
229 end do
230
231 safe_deallocate_a(zpsi)
232 safe_deallocate_a(zst)
233 safe_deallocate_a(zchi)
234
236 end subroutine target_chi_userdefined
237
239
240!! Local Variables:
241!! mode: f90
242!! coding: utf-8
243!! End:
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:609
real(real64), parameter, public m_zero
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
Definition: io.F90:114
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
logical pure function, public ion_dynamics_ions_move(this)
This module defines various routines, operating on mesh functions.
subroutine, public zmf_normalize(mesh, dim, psi, norm)
Normalize a mesh function psi.
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_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
this module contains the low-level part of the output system
Definition: output_low.F90:115
this module contains the output system
Definition: output.F90:115
subroutine, public output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
Definition: output.F90:1888
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public conv_to_c_string(str)
converts to c string
Definition: string.F90:252
subroutine, public target_chi_userdefined(tg, gr, psi_in, chi_out)
real(real64) function, public target_j1_userdefined(tg, gr, psi)
subroutine, public target_init_userdefined(gr, namespace, tg, td)
subroutine, public target_output_userdefined(tg, namespace, space, gr, dir, ions, hm, outp)
Definition: td.F90:114
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
output handler class
Definition: output_low.F90:164
The states_elec_t class contains all electronic wave functions.