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