Octopus
potential_interpolation.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
25 use math_oct_m
26 use mesh_oct_m
31 use space_oct_m
32
33 implicit none
34
35 private
36 public :: &
48
50 private
51 real(real64), allocatable, public :: v_old(:, :, :)
52 real(real64), allocatable, public :: vtau_old(:, :, :)
53 logical :: mgga_with_exc
55
56 integer :: interpolation_steps
57
58contains
59
60 ! ---------------------------------------------------------
61 subroutine potential_interpolation_copy(vkso, vksi)
62 type(potential_interpolation_t), intent(inout) :: vkso
63 type(potential_interpolation_t), intent(in) :: vksi
64
66
67 safe_allocate_source_a(vkso%v_old, vksi%v_old)
68 safe_allocate_source_a(vkso%vtau_old, vksi%vtau_old)
69 vkso%mgga_with_exc = vksi%mgga_with_exc
70
72 end subroutine potential_interpolation_copy
73 ! ---------------------------------------------------------
74
75 ! ---------------------------------------------------------
76 subroutine potential_interpolation_init(potential_interpolation, np, nspin, mgga_with_exc, order)
77 type(potential_interpolation_t), intent(inout) :: potential_interpolation
78 integer, intent(in) :: np, nspin
79 logical, intent(in) :: mgga_with_exc
80 integer, optional, intent(in) :: order
81
83
84 interpolation_steps = optional_default(order, 3)
85
86 safe_allocate(potential_interpolation%v_old(1:np, 1:nspin, 0:interpolation_steps))
87 potential_interpolation%v_old(:, :, :) = m_zero
88
89 potential_interpolation%mgga_with_exc = mgga_with_exc
90
91 if (potential_interpolation%mgga_with_exc) then
92
93 safe_allocate(potential_interpolation%vtau_old(1:np, 1:nspin, 0:interpolation_steps))
94 potential_interpolation%vtau_old(:, :, :) = m_zero
95
96 end if
97
99 end subroutine potential_interpolation_init
100 ! ---------------------------------------------------------
101
102 ! ---------------------------------------------------------
103 subroutine potential_interpolation_end(potential_interpolation)
104 type(potential_interpolation_t), intent(inout) :: potential_interpolation
105
107
108 assert(allocated(potential_interpolation%v_old))
109 safe_deallocate_a(potential_interpolation%v_old)
110 safe_deallocate_a(potential_interpolation%vtau_old)
111
113 end subroutine potential_interpolation_end
114 ! ---------------------------------------------------------
115
116 ! ---------------------------------------------------------
117 subroutine potential_interpolation_run_zero_iter(potential_interpolation, np, nspin, vhxc, vtau)
118 type(potential_interpolation_t), intent(inout) :: potential_interpolation
119 integer, intent(in) :: np, nspin
120 real(real64), intent(in) :: vhxc(:, :)
121 real(real64), optional, intent(in) :: vtau(:, :)
122
123 integer :: i, ispin, ip
124
126
127 do i = 1, interpolation_steps
128 do ispin = 1, nspin
129 do ip = 1, np
130 potential_interpolation%v_old(ip, ispin, i) = vhxc(ip, ispin)
131 end do
132 end do
133 end do
134
135 if (potential_interpolation%mgga_with_exc) then
136 assert(present(vtau))
137 do i = 1, interpolation_steps
138 do ispin = 1, nspin
139 do ip = 1, np
140 potential_interpolation%vtau_old(ip, ispin, i) = vtau(ip, ispin)
141 end do
142 end do
143 end do
144 end if
148 ! ---------------------------------------------------------
150 ! ---------------------------------------------------------
151 subroutine potential_interpolation_new(potential_interpolation, np, nspin, time, dt, vhxc, vtau)
152 type(potential_interpolation_t), intent(inout) :: potential_interpolation
153 integer, intent(in) :: np, nspin
154 real(real64), intent(in) :: time, dt
155 real(real64), contiguous, intent(in) :: vhxc(:, :)
156 real(real64), contiguous, optional, intent(in) :: vtau(:, :)
157
158 real(real64), allocatable :: times(:)
159 integer :: j
160
162
163 safe_allocate(times(1:interpolation_steps))
164 do j = 1, interpolation_steps
165 times(j) = time - j*dt
166 end do
167
168 do j = interpolation_steps, 2, -1
169 call lalg_copy(np, nspin, potential_interpolation%v_old(:, :, j-1), potential_interpolation%v_old(:, :, j))
170 end do
171 call lalg_copy(np, nspin, vhxc(:, :), potential_interpolation%v_old(:, :, 1))
172 call interpolate( times, potential_interpolation%v_old(:, :, 1:interpolation_steps), &
173 time, potential_interpolation%v_old(:, :, 0))
174
175 if (potential_interpolation%mgga_with_exc) then
176 assert(present(vtau))
177 do j = interpolation_steps, 2, -1
178 call lalg_copy(np, nspin, potential_interpolation%vtau_old(:, :, j-1), potential_interpolation%vtau_old(:, :, j))
179 end do
180 call lalg_copy(np, nspin, vtau(:, :), potential_interpolation%vtau_old(:, :, 1))
181 call interpolate( times, potential_interpolation%vtau_old(:, :, 1:interpolation_steps), &
182 time, potential_interpolation%vtau_old(:, :, 0))
183 end if
184
185 safe_deallocate_a(times)
187 end subroutine potential_interpolation_new
188 ! ---------------------------------------------------------
189
190 ! ---------------------------------------------------------
191 subroutine potential_interpolation_set(potential_interpolation, np, nspin, i, vhxc, vtau)
192 type(potential_interpolation_t), intent(inout) :: potential_interpolation
193 integer, intent(in) :: np, nspin
194 integer, intent(in) :: i
195 real(real64), contiguous, intent(inout) :: vhxc(:, :)
196 real(real64), contiguous, optional,intent(in) :: vtau(:, :)
197
199
200 call lalg_copy(np, nspin, vhxc, potential_interpolation%v_old(:, :, i))
201
202 if (potential_interpolation%mgga_with_exc) then
203 assert(present(vtau))
204 call lalg_copy(np, nspin, vtau, potential_interpolation%vtau_old(:, :, i))
205 end if
206
208 end subroutine potential_interpolation_set
209 ! ---------------------------------------------------------
211 ! ---------------------------------------------------------
212 subroutine potential_interpolation_get(potential_interpolation, np, nspin, i, vhxc, vtau)
213 type(potential_interpolation_t), intent(inout) :: potential_interpolation
214 integer, intent(in) :: np, nspin
215 integer, intent(in) :: i
216 real(real64), contiguous, intent(inout) :: vhxc(:, :)
217 real(real64), optional, contiguous, intent(inout) :: vtau(:, :)
218
220
221 call lalg_copy(np, nspin, potential_interpolation%v_old(:, :, i), vhxc)
222
223 if (potential_interpolation%mgga_with_exc) then
224 assert(present(vtau))
225 call lalg_copy(np, nspin, potential_interpolation%vtau_old(:, :, i), vtau)
226 end if
227
229 end subroutine potential_interpolation_get
230 ! ---------------------------------------------------------
231
232
233 ! ---------------------------------------------------------
234 subroutine potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
235 type(potential_interpolation_t), intent(inout) :: potential_interpolation
236 integer, intent(in) :: order
237 real(real64), intent(in) :: time, dt, t
238 real(real64), contiguous, intent(inout) :: vhxc(:, :)
239 real(real64), contiguous, optional, intent(inout) :: vtau(:, :)
240
241 integer :: j
242 real(real64), allocatable :: times(:)
243
245
246 safe_allocate(times(1:interpolation_steps))
247 do j = 1, interpolation_steps
248 times(j) = time - (j-1)*dt
249 end do
250
251 call interpolate( times(1:order), potential_interpolation%v_old(:, :, 0:order-1), t, vhxc(:, :))
252
253 if (potential_interpolation%mgga_with_exc) then
254 assert(present(vtau))
255 call interpolate( times(1:order), potential_interpolation%vtau_old(:, :, 0:order-1), t, vtau(:, :))
256 end if
257
258 safe_deallocate_a(times)
261 ! ---------------------------------------------------------
262
263 ! ---------------------------------------------------------
264 subroutine potential_interpolation_dump(potential_interpolation, space, restart, mesh, nspin, err2)
265 type(potential_interpolation_t), intent(in) :: potential_interpolation
266 type(restart_t), intent(in) :: restart
267 class(space_t), intent(in) :: space
268 class(mesh_t), intent(in) :: mesh
269 integer, intent(in) :: nspin
270 integer, intent(out) :: err2
271
272 integer :: ii, is, err
273 character(len=256) :: filename
274 !type(mpi_grp_t) :: mpi_grp
275
277
278 err2 = 0
279 do ii = 1, interpolation_steps - 1
280 do is = 1, nspin
281 write(filename,'(a6,i2.2,i3.3)') 'vprev_', ii, is
282 call restart_write_mesh_function(restart, space, filename, mesh, &
283 potential_interpolation%v_old(1:mesh%np, is, ii), err)
284 ! the unit is energy actually, but this only for restart, and can be kept in atomic units
285 ! for simplicity
286 if (err /= 0) err2 = err2 + 1
287 end do
288 end do
289
290
291 if (potential_interpolation%mgga_with_exc) then
292 err2 = 0
293 do ii = 1, interpolation_steps - 1
294 do is = 1, nspin
295 write(filename,'(a9,i2.2,i3.3)') 'vtauprev_', ii, is
296 call restart_write_mesh_function(restart, space, filename, mesh, &
297 potential_interpolation%vtau_old(1:mesh%np, is, ii), err)
298 ! the unit is energy actually, but this only for restart, and can be kept in atomic units
299 ! for simplicity
300 if (err /= 0) err2 = err2 + 1
301 end do
302 end do
303 end if
304
306 end subroutine potential_interpolation_dump
307 ! ---------------------------------------------------------
308
309 ! ---------------------------------------------------------
310 subroutine potential_interpolation_load(potential_interpolation, namespace, space, restart, mesh, nspin, err2)
311 type(potential_interpolation_t), intent(inout) :: potential_interpolation
312 type(namespace_t), intent(in) :: namespace
313 class(space_t), intent(in) :: space
314 type(restart_t), intent(in) :: restart
315 class(mesh_t), intent(in) :: mesh
316 integer, intent(in) :: nspin
317 integer, intent(out) :: err2
318
319 integer :: ii, is, err
320 character(len=256) :: filename
321
323
324 err2 = 0
325 do ii = 1, interpolation_steps - 1
326 do is = 1, nspin
327 write(filename,'(a,i2.2,i3.3)') 'vprev_', ii, is
328 call restart_read_mesh_function(restart, space, filename, mesh, &
329 potential_interpolation%v_old(1:mesh%np, is, ii), err)
330 if (err /= 0) then
331 err2 = err2 + 1
332 message(1) = "Unable to read VKS restart file '" // trim(filename) // "'"
333 call messages_warning(1, namespace=namespace)
334 end if
335 end do
336 end do
337
338 if (potential_interpolation%mgga_with_exc) then
339 err2 = 0
340 do ii = 1, interpolation_steps - 1
341 do is = 1, nspin
342 write(filename,'(a,i2.2,i3.3)') 'vtauprev_', ii, is
343 call restart_read_mesh_function(restart, space, filename, mesh, &
344 potential_interpolation%vtau_old(1:mesh%np, is, ii), err)
345 if (err /= 0) then
346 err2 = err2 + 1
347 message(1) = "Unable to read VKS restart file '" // trim(filename) // "'"
348 call messages_warning(1, namespace=namespace)
349 end if
350 end do
351 end do
352 end if
353
355 end subroutine potential_interpolation_load
356 ! ---------------------------------------------------------
359
360!! Local Variables:
361!! mode: f90
362!! coding: utf-8
363!! End:
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:185
This is the common interface to a simple-minded polynomical interpolation procedure (simple use of th...
Definition: math.F90:186
real(real64), parameter, public m_zero
Definition: global.F90:187
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public potential_interpolation_copy(vkso, vksi)
subroutine, public potential_interpolation_end(potential_interpolation)
subroutine, public potential_interpolation_interpolate(potential_interpolation, order, time, dt, t, vhxc, vtau)
subroutine, public potential_interpolation_set(potential_interpolation, np, nspin, i, vhxc, vtau)
subroutine, public potential_interpolation_new(potential_interpolation, np, nspin, time, dt, vhxc, vtau)
subroutine, public potential_interpolation_run_zero_iter(potential_interpolation, np, nspin, vhxc, vtau)
subroutine, public potential_interpolation_get(potential_interpolation, np, nspin, i, vhxc, vtau)
subroutine, public potential_interpolation_load(potential_interpolation, namespace, space, restart, mesh, nspin, err2)
subroutine, public potential_interpolation_dump(potential_interpolation, space, restart, mesh, nspin, err2)
subroutine, public potential_interpolation_init(potential_interpolation, np, nspin, mgga_with_exc, order)
Describes mesh distribution to nodes.
Definition: mesh.F90:186