Octopus
curv_modine.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
27
30 use debug_oct_m
31 use global_oct_m
32 use, intrinsic :: iso_fortran_env
36 use parser_oct_m
39 use unit_oct_m
41
42 implicit none
43
44 private
45 public :: &
48
49 type, extends(coordinate_system_t) :: curv_modine_t
50 private
51 real(real64), allocatable :: lsize(:)
52 real(real64) :: xbar
53 real(real64) :: Jbar
54
55 real(real64), allocatable :: Jlocal(:)
56 real(real64), allocatable :: Jrange(:)
57
58 real(real64), allocatable :: chi_atoms(:,:)
59 real(real64), allocatable :: csi(:,:)
60
61 integer :: npos
62 type(root_solver_t) :: rs
63 contains
64 procedure :: to_cartesian => curv_modine_to_cartesian
65 procedure :: from_cartesian => curv_modine_from_cartesian
66 procedure :: det_jac => curv_modine_det_jac
67 procedure :: write_info => curv_modine_write_info
68 procedure :: surface_element => curv_modine_surface_element
70 end type curv_modine_t
71
72 interface curv_modine_t
73 procedure curv_modine_constructor
74 end interface curv_modine_t
75
76 ! Auxiliary variables for the root solver.
77 class(curv_modine_t), pointer :: modine_p
78 real(real64), allocatable :: x_p(:)
79
80contains
81
82 ! ---------------------------------------------------------
83 function curv_modine_constructor(namespace, dim, npos, pos, lsize, spacing) result(modine)
84 type(namespace_t), intent(in) :: namespace
85 integer, intent(in) :: dim
86 integer, intent(in) :: npos
87 real(real64), intent(in) :: pos(1:dim,1:npos)
88 real(real64), intent(in) :: lsize(1:dim)
89 real(real64), intent(in) :: spacing(1:dim)
90 class(curv_modine_t), pointer :: modine
91
93
94 safe_allocate(modine)
95
96 modine%dim = dim
97 modine%local_basis = .true.
98 modine%orthogonal = .true. ! This needs to be checked.
99
100 modine%npos = npos
101
102 !%Variable CurvModineXBar
103 !%Type float
104 !%Default 1/3
105 !%Section Mesh::Curvilinear::Modine
106 !%Description
107 !% Size of central flat region (in units of <tt>Lsize</tt>). Must be between 0 and 1.
108 !% See N. A. Modine, G. Zumbach, and E. Kaxiras, <i>Phys. Rev. B</i> <b>55</b>, 10289-10301 (1997).
109 !%End
110 call parse_variable(namespace, 'CurvModineXBar', m_one/m_three, modine%xbar)
111
112 !%Variable CurvModineJBar
113 !%Type float
114 !%Default 1/2
115 !%Section Mesh::Curvilinear::Modine
116 !%Description
117 !% Increase in density of points is inverse of this parameter.
118 !% See N. A. Modine, G. Zumbach, and E. Kaxiras, <i>Phys. Rev. B</i> <b>55</b>, 10289-10301 (1997).
119 !%End
120 call parse_variable(namespace, 'CurvModineJBar', m_half, modine%Jbar)
122 safe_allocate(modine%lsize(1:dim))
123 modine%lsize = lsize / modine%Jbar
124
125 if (modine%xbar < m_zero .or. modine%xbar > m_one) then
126 message(1) = 'The parameter "CurvModineXBar" must lie between 0 and 1.'
127 call messages_fatal(1, namespace=namespace)
128 end if
129
130 safe_allocate(modine%Jlocal(1:npos))
131 safe_allocate(modine%Jrange(1:npos))
132
133 ! \warning: the reading has to be done for each atom kind
134
135 !%Variable CurvModineJlocal
136 !%Type float
137 !%Default 0.25
138 !%Section Mesh::Curvilinear::Modine
139 !%Description
140 !% Local refinement around the atoms. Must be between 0 and 1.
141 !% See N. A. Modine, G. Zumbach, and E. Kaxiras, <i>Phys. Rev. B</i> <b>55</b>, 10289-10301 (1997).
142 !%End
143 call parse_variable(namespace, 'CurvModineJlocal', 0.25_real64, modine%Jlocal(1))
145 !%Variable CurvModineJrange
146 !%Type float
147 !%Default 2 b
148 !%Section Mesh::Curvilinear::Modine
149 !%Description
150 !% Local refinement range (a length).
151 !% See N. A. Modine, G. Zumbach, and E. Kaxiras, <i>Phys. Rev. B</i> <b>55</b>, 10289-10301 (1997).
152 !%End
153 call parse_variable(namespace, 'CurvModineJrange', m_two, modine%Jrange(1), units_inp%length)
155 if (modine%Jlocal(1)<m_zero.or.modine%Jlocal(1)>m_one) then
156 message(1) = 'The parameter "CurvModineJlocal" must lie between 0 and 1.'
157 call messages_fatal(1, namespace=namespace)
158 end if
160 modine%Jlocal(:) = modine%Jlocal(1)
161 modine%Jrange(:) = modine%Jrange(1)
163 ! initialize root solver for the optimization
164 call root_solver_init(modine%rs, namespace, dim, &
165 solver_type = root_newton, maxiter = 500, abs_tolerance = 1.0e-10_real64)
166
167 call find_atom_points()
168 call optimize()
169
170 modine%min_mesh_scaling_product = modine%Jbar**modine%dim
173 contains
174
175 subroutine find_atom_points()
176 integer :: iatom, jj
177
179
180 ! Initialize csi
181 safe_allocate(modine%csi(1:modine%dim, 1:modine%npos))
182 modine%csi = npos
183
184 ! get first estimate for chi_atoms
185 safe_allocate(modine%chi_atoms(1:modine%dim, 1:modine%npos))
186 do jj = 1, 10 ! \warning: make something better
187 do iatom = 1, modine%npos
188 modine%chi_atoms(:, iatom) = modine%from_cartesian(pos(:, iatom))
189 end do
190 modine%csi(:,:) = modine%chi_atoms(:,:)
191 end do
192
193 do iatom = 1, modine%npos
194 ! These are the chi positions where we want the atoms.
195 modine%chi_atoms(:, iatom) = nint(modine%chi_atoms(:, iatom) / spacing(:)) * spacing(:)
196 end do
197
199 end subroutine find_atom_points
200
201 subroutine optimize()
202 logical :: conv
203 integer :: iatom, idim, index
204 real(real64), allocatable :: my_csi(:), start_csi(:)
205
207
208 modine_p => modine
209
210 safe_allocate(x_p(1:modine%dim*modine%npos))
211 safe_allocate(my_csi(1:modine%dim*modine%npos))
212 safe_allocate(start_csi(1:modine%dim*modine%npos))
213
214 do iatom = 1, modine%npos
215 do idim = 1, modine%dim
216 index = (iatom-1)*dim + idim
217 x_p(index) = pos(idim, iatom)
218 start_csi(index) = modine%chi_atoms(idim, iatom)
219 end do
220 end do
221
222 call droot_solver_run(modine%rs, getf2, my_csi, conv, startval=start_csi)
223
224 if (.not. conv) then
225 message(1) = "During the construction of the adaptive grid, the Newton-Raphson"
226 message(2) = "method did not converge."
227 call messages_fatal(2, namespace=namespace)
228 end if
229
230 ! Now set csi to the new values
231 do iatom = 1, modine%npos
232 do idim = 1, modine%dim
233 index = (iatom-1)*modine%dim + idim
234 modine_p%csi(idim, iatom) = my_csi(index)
235 end do
236 end do
237
238 safe_deallocate_a(x_p)
239 safe_deallocate_a(my_csi)
240 safe_deallocate_a(start_csi)
241
242 nullify(modine_p)
243
245 end subroutine optimize
246
247 end function curv_modine_constructor
248
249 ! ---------------------------------------------------------
250 subroutine curv_modine_copy(this_out, this_in)
251 type(curv_modine_t), intent(inout) :: this_out
252 type(curv_modine_t), intent(in) :: this_in
253
254 push_sub(curv_modine_copy)
255
256 this_out%dim = this_in%dim
257 this_out%local_basis = this_in%local_basis
258 safe_allocate_source_a(this_out%lsize, this_in%lsize)
259 this_out%xbar = this_in%xbar
260 this_out%Jbar = this_in%Jbar
261 safe_allocate_source_a(this_out%Jlocal, this_in%Jlocal)
262 safe_allocate_source_a(this_out%Jrange, this_in%Jrange)
263 safe_allocate_source_a(this_out%chi_atoms, this_in%chi_atoms)
264 safe_allocate_source_a(this_out%csi, this_in%csi)
265 this_out%npos = this_in%npos
266
267 pop_sub(curv_modine_copy)
268 end subroutine curv_modine_copy
269
270 ! ---------------------------------------------------------
271 subroutine curv_modine_finalize(this)
272 type(curv_modine_t), intent(inout) :: this
273
274 push_sub(curv_modine_finalize)
275
276 safe_deallocate_a(this%lsize)
277 safe_deallocate_a(this%Jlocal)
278 safe_deallocate_a(this%Jrange)
279 safe_deallocate_a(this%chi_atoms)
280 safe_deallocate_a(this%csi)
281
282 pop_sub(curv_modine_finalize)
283 end subroutine curv_modine_finalize
284
285 ! ---------------------------------------------------------
286 function curv_modine_to_cartesian(this, chi) result(xx)
287 class(curv_modine_t), target, intent(in) :: this
288 real(real64), intent(in) :: chi(:)
289 real(real64) :: xx(1:this%dim)
290
291 real(real64) :: chi2(this%dim), rr2, dd
292 integer :: iatom
293
294 ! no PUSH_SUB, called too often
295
296 call curv_modine_chi2chi2(this, chi, chi2)
297
298 xx(:) = chi2(:)
299 do iatom = 1, this%npos
300 rr2 = sum((chi2(:) - this%csi(:, iatom))**2)
301 dd = exp(-rr2/(m_two*this%Jrange(iatom)**2))
302
303 xx(:) = xx(:) - this%Jlocal(iatom)*(chi2(:) - this%csi(:, iatom)) * dd
304 end do
305
306 end function curv_modine_to_cartesian
307
308 ! ---------------------------------------------------------
309 function curv_modine_from_cartesian(this, xx) result(chi)
310 class(curv_modine_t), target, intent(in) :: this
311 real(real64), intent(in) :: xx(:)
312 real(real64) :: chi(1:this%dim)
313
314 logical :: conv
315
316 ! no PUSH_SUB, called too often
317
318 modine_p => this
319 safe_allocate(x_p(1:this%dim))
320 x_p(:) = xx(:)
321
322 call droot_solver_run(this%rs, getf, chi, conv, startval = xx)
323
324 safe_deallocate_a(x_p)
325 nullify(modine_p)
326
327 if (.not. conv) then
328 message(1) = "During the construction of the adaptive grid, the Newton-Raphson"
329 message(2) = "method did not converge for point:"
330 write(message(3),'(3f14.6)') xx(1:this%dim)
331 call messages_fatal(3)
332 end if
333
334 end function curv_modine_from_cartesian
335
336 ! ---------------------------------------------------------
337 real(real64) function curv_modine_det_jac(this, xx, chi) result(jdet)
338 class(curv_modine_t), intent(in) :: this
339 real(real64), intent(in) :: xx(:)
340 real(real64), intent(in) :: chi(:)
341
342 real(real64) :: dummy(this%dim), jac(1:this%dim, 1:this%dim)
344 ! no PUSH_SUB, called too often
345
346 call curv_modine_jacobian_inv(this, chi, dummy, jac)
347 jdet = lalg_determinant(this%dim, jac, preserve_mat = .false.)
348
349 end function curv_modine_det_jac
350
351 ! ---------------------------------------------------------
352 subroutine curv_modine_write_info(this, iunit, namespace)
353 class(curv_modine_t), intent(in) :: this
354 integer, optional, intent(in) :: iunit
355 type(namespace_t), optional, intent(in) :: namespace
356
357 push_sub(curv_modine_write_info)
358
359 write(message(1), '(a)') ' Curvilinear Method = modine'
360 call messages_info(1, iunit=iunit, namespace=namespace)
361
363 end subroutine curv_modine_write_info
365 ! ---------------------------------------------------------
366 real(real64) function curv_modine_surface_element(this, idir) result(ds)
367 class(curv_modine_t), intent(in) :: this
368 integer, intent(in) :: idir
369
370 ds = m_zero
371 message(1) = 'Surface element with modine curvilinear coordinates not implemented'
372 call messages_fatal(1)
373
374 end function curv_modine_surface_element
375
376 ! ---------------------------------------------------------
377 pure subroutine curv_modine_chi2chi2(this, chi_, chi2, Jac)
378 class(curv_modine_t), intent(in) :: this
379 real(real64), intent(in) :: chi_(:)
380 real(real64), intent(out) :: chi2(:)
381 real(real64), optional, intent(out) :: jac(:)
382
383 integer, parameter :: qq = 3
384 real(real64) :: chibar(this%dim), rr, chi
385 logical :: neg
386 integer :: i
387
388 ! no PUSH_SUB, called too often
389
390 chibar = this%xbar*this%lsize
391
392 do i = 1, this%dim
393 neg = (chi_(i) < 0)
394 chi = abs(chi_(i))
395
396 chi2(i) = this%Jbar * chi
397 if (present(jac)) jac(i) = this%Jbar
398
399 if (chi > chibar(i)) then
400 rr = (chi - chibar(i))/(this%lsize(i) - chibar(i))
401
402 chi2(i) = chi2(i) + this%lsize(i)/m_two*(1 - this%Jbar) * rr**qq * &
403 (qq + m_one - (qq - m_one)*rr)
404
405 if (present(jac)) then
406 jac(i) = jac(i) + this%lsize(i)/m_two*(m_one - this%Jbar) * rr**(qq - 1)/(this%lsize(i) - chibar(i)) * &
407 (qq*(qq + m_one) - (qq**2 - m_one)*rr)
408 end if
409 end if
410
411 if (neg) chi2(i) = -chi2(i)
412 ! CHECK if Jacobian does not have to be negated!
413 end do
414
415 end subroutine curv_modine_chi2chi2
416
417 ! ---------------------------------------------------------
418 pure subroutine curv_modine_jacobian_inv(this, chi, xx, Jac)
419 type(curv_modine_t), intent(in) :: this
420 real(real64), intent(in) :: chi(:)
421 real(real64), intent(out) :: xx(:)
422 real(real64), intent(out) :: jac(:, :)
423
424 real(real64) :: chi2(this%dim), rr2, dd, j2(this%dim)
425 integer :: iatom, idim, idim2
426
427 ! no PUSH_SUB, called too often
428
429 call curv_modine_chi2chi2(this, chi, chi2, j2)
431 ! initialize both xx and the Jacobian
432 xx(:) = chi2(:)
433 jac(:,:) = m_zero
434 do idim = 1, this%dim
435 jac(idim, idim) = m_one
436 end do
437
438 do iatom = 1, this%npos
439 rr2 = sum((chi2(:) - this%csi(:, iatom))**2)
440 dd = exp(-rr2/(m_two*this%Jrange(iatom)**2))
441
442 xx(:) = xx(:) - this%Jlocal(iatom)*(chi2(:) - this%csi(:, iatom)) * dd
443
444 do idim = 1, this%dim
445 jac(idim, idim) = jac(idim, idim) - this%Jlocal(iatom) * dd
446 do idim2 = 1, this%dim
447 jac(idim, idim2) = jac(idim, idim2) + &
448 this%Jlocal(iatom)*(chi2(idim) - this%csi(idim, iatom))*(chi2(idim2) - this%csi(idim2, iatom)) * &
449 m_two/(m_two*this%Jrange(iatom)**2) * dd
450 end do
451 end do
452 end do
453
454 do idim = 1, this%dim
455 jac(idim, :) = jac(idim, :) * j2(:)
456 end do
457
458 end subroutine curv_modine_jacobian_inv
460 ! ---------------------------------------------------------
461 pure subroutine getf(yy, ff, jf)
462 real(real64), intent(in) :: yy(:)
463 real(real64), intent(out) :: ff(:), jf(:, :)
464
465 ! no PUSH_SUB, called too often
466
467 call curv_modine_jacobian_inv(modine_p, yy, ff, jf)
468 ff(:) = ff(:) - x_p(:)
469
470 end subroutine getf
471
472 ! ---------------------------------------------------------
473 subroutine getf2(csi, ff, jf)
474 real(real64), intent(in) :: csi(:)
475 real(real64), intent(out) :: ff(:), jf(:, :)
476
477 integer :: i1, j1, i2, j2, index1, index2
478 real(real64) :: rr2, dd, dd2
479 real(real64), allocatable :: xx(:), chi2(:)
480
481 ! no PUSH_SUB, called too often
482
483 safe_allocate(xx(1:modine_p%dim))
484 safe_allocate(chi2(1:modine_p%dim))
485
486 ! first we fill in coord_system%csi with the values we have
487 index1 = 1
488 do i1 = 1, modine_p%npos
489 do j1 = 1, modine_p%dim
490 modine_p%csi(j1, i1) = csi(index1)
491 index1 = index1 + 1
492 end do
493 end do
494
495 ! get ff and jf
496 jf(:,:) = m_zero
497 do i1 = 1, modine_p%npos
498 call curv_modine_chi2chi2(modine_p, modine_p%chi_atoms(:,i1), chi2)
499 xx = chi2
500
501 do i2 = 1, modine_p%npos
502 rr2 = sum((chi2 - modine_p%csi(:,i2))**2)
503 dd = exp(-rr2/(m_two*modine_p%Jrange(i2)**2))
504
505 xx = xx - modine_p%Jlocal(i2)*(chi2 - modine_p%csi(:,i2)) * dd
506 end do
507
508 do j1 = 1, modine_p%dim
509 index1 = (i1 - 1)*modine_p%dim + j1
510 ff(index1) = xx(j1) - x_p(index1)
512 do i2 = 1, modine_p%npos
513 rr2 = sum((chi2 - modine_p%csi(:,i2))**2)
514 dd = exp(-rr2/(m_two*modine_p%Jrange(i2)**2))
515 dd2 = -m_two/(m_two*modine_p%Jrange(i2)**2)*dd
516
517 index2 = (i2 - 1)*modine_p%dim + j1
518 jf(index1, index2) = modine_p%Jlocal(i2) * dd
519
520 do j2 = 1, modine_p%dim
521 index2 = (i2 - 1)*modine_p%dim + j2
522
523 jf(index1, index2) = jf(index1, index2) + modine_p%Jlocal(i2) * dd2 * &
524 (chi2(j1) - modine_p%csi(j1,i2))*(chi2(j2) - modine_p%csi(j2,i2))
525 end do
526 end do
527 end do
528 end do
529
530 safe_deallocate_a(xx)
531 safe_deallocate_a(chi2)
532
533 end subroutine getf2
534
535end module curv_modine_oct_m
536
537!! Local Variables:
538!! mode: f90
539!! coding: utf-8
540!! End:
subroutine optimize()
subroutine find_atom_points()
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine.
Definition: lalg_adv.F90:187
double exp(double __x) __attribute__((__nothrow__
This module implements the curvilinear coordinates given in N. A. Modine, G. Zumbach,...
real(real64) function curv_modine_det_jac(this, xx, chi)
class(curv_modine_t), pointer modine_p
real(real64) function, dimension(1:this%dim) curv_modine_from_cartesian(this, xx)
real(real64) function curv_modine_surface_element(this, idir)
pure subroutine curv_modine_jacobian_inv(this, chi, xx, Jac)
class(curv_modine_t) function, pointer curv_modine_constructor(namespace, dim, npos, pos, lsize, spacing)
subroutine getf2(csi, ff, jf)
pure subroutine getf(yy, ff, jf)
subroutine curv_modine_finalize(this)
real(real64), dimension(:), allocatable x_p
pure subroutine curv_modine_chi2chi2(this, chi_, chi2, Jac)
subroutine, public curv_modine_copy(this_out, this_in)
subroutine curv_modine_write_info(this, iunit, namespace)
real(real64) function, dimension(1:this%dim) curv_modine_to_cartesian(this, chi)
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
real(real64), parameter, public m_three
Definition: global.F90:190
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:420
integer, parameter, public root_newton
subroutine, public root_solver_init(rs, namespace, dimensionality, solver_type, maxiter, rel_tolerance, abs_tolerance)
subroutine, public droot_solver_run(rs, func, root, success, startval)
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_inp
the units systems for reading and writing
abstract class to describe coordinate systems
int true(void)