55#include "../base/base_uses.f90"
59 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'ps_implicit_methods'
67 INTERFACE ps_implicit_compute_ehartree
68 MODULE PROCEDURE compute_ehartree_periodic_bc, &
69 compute_ehartree_mixed_bc
70 END INTERFACE ps_implicit_compute_ehartree
72 REAL(dp),
PRIVATE,
PARAMETER :: large_error = 1.0e4_dp
95 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_create'
97 INTEGER :: boundary_condition, handle, j, &
98 n_contacts, neumann_directions
101 CALL timeset(routinen, handle)
105 IF (.NOT.
ASSOCIATED(ps_implicit_env))
THEN
106 ALLOCATE (ps_implicit_env)
108 ps_implicit_env%do_dbc_cube = poisson_params%dbc_params%do_dbc_cube
109 boundary_condition = poisson_params%ps_implicit_params%boundary_condition
110 neumann_directions = poisson_params%ps_implicit_params%neumann_directions
113 NULLIFY (ps_implicit_env%dielectric)
114 SELECT CASE (boundary_condition)
116 CALL dielectric_create(ps_implicit_env%dielectric, pw_pool, poisson_params%dielectric_params)
119 CALL dielectric_create(ps_implicit_env%dielectric, pw_pool_xpndd, poisson_params%dielectric_params)
124 NULLIFY (ps_implicit_env%initial_guess)
127 NULLIFY (ps_implicit_env%v_eps)
128 ALLOCATE (ps_implicit_env%v_eps)
129 CALL pw_pool%create_pw(ps_implicit_env%v_eps)
130 CALL pw_zero(ps_implicit_env%v_eps)
133 NULLIFY (ps_implicit_env%cstr_charge)
134 SELECT CASE (boundary_condition)
136 ALLOCATE (ps_implicit_env%cstr_charge)
137 CALL pw_pool%create_pw(ps_implicit_env%cstr_charge)
138 CALL pw_zero(ps_implicit_env%cstr_charge)
141 ALLOCATE (ps_implicit_env%cstr_charge)
142 CALL pw_pool_xpndd%create_pw(ps_implicit_env%cstr_charge)
143 CALL pw_zero(ps_implicit_env%cstr_charge)
148 ps_implicit_env%ehartree = 0.0_dp
149 ps_implicit_env%electric_enthalpy = 0.0_dp
151 ps_implicit_env%times_called = 0
155 CALL dct_type_init(pw_pool%pw_grid, neumann_directions, ps_implicit_env%dct_env)
160 CALL ps_implicit_prepare_blocks(pw_pool, dct_pw_grid, green, poisson_params, ps_implicit_env)
163 (.NOT. poisson_params%dbc_params%do_dbc_cube))
THEN
164 n_contacts =
SIZE(ps_implicit_env%contacts)
172 CALL timestop(handle)
191 REAL(dp),
INTENT(OUT),
OPTIONAL :: ehartree
193 CHARACTER(LEN=*),
PARAMETER :: routinen =
'implicit_poisson_solver_periodic'
195 INTEGER :: handle, iter, max_iter, outp_unit, &
197 LOGICAL :: reached_max_iter, reached_tol, &
198 use_zero_initial_guess
199 REAL(dp) :: nabs_error, omega, pres_error, tol
207 CALL timeset(routinen, handle)
209 pw_pool => poisson_env%pw_pools(poisson_env%pw_level)%pool
210 dielectric => poisson_env%implicit_env%dielectric
211 green => poisson_env%green_fft
212 ps_implicit_env => poisson_env%implicit_env
214 tol = poisson_env%parameters%ps_implicit_params%tol
215 omega = poisson_env%parameters%ps_implicit_params%omega
216 max_iter = poisson_env%parameters%ps_implicit_params%max_iter
217 use_zero_initial_guess = poisson_env%parameters%ps_implicit_params%zero_initial_guess
218 times_called = ps_implicit_env%times_called
221 IF (times_called .EQ. 0)
CALL ps_implicit_initial_guess_create(ps_implicit_env, pw_pool)
223 CALL pw_pool%create_pw(g)
224 CALL pw_pool%create_pw(v_old)
225 CALL pw_pool%create_pw(res_old)
226 CALL pw_pool%create_pw(res_new)
227 CALL pw_pool%create_pw(qainvxres)
228 CALL pw_pool%create_pw(pxqainvxres)
230 IF (use_zero_initial_guess)
THEN
233 CALL pw_copy(ps_implicit_env%initial_guess, v_old)
236 g%array =
fourpi*density%array/dielectric%eps%array
239 CALL apply_poisson_operator_fft(pw_pool, green, dielectric, v_old, res_old)
244 CALL apply_inv_laplace_operator_fft(pw_pool, green, res_old, qainvxres)
256 CALL apply_p_operator(pw_pool, dielectric, qainvxres, pxqainvxres)
257 CALL pw_copy(pxqainvxres, res_new)
259 CALL pw_axpy(res_old, res_new, 1.0_dp - omega)
262 CALL ps_implicit_compute_error_fft(pw_pool, green, res_new, v_old, v_new, qainvxres, &
263 pres_error, nabs_error)
265 CALL ps_implicit_output(iter, pres_error, nabs_error, outp_unit)
266 IF (
PRESENT(ehartree))
THEN
267 CALL ps_implicit_compute_ehartree(density, v_new, ehartree)
268 CALL ps_implicit_report_ehartree(ps_implicit_env, outp_unit, ehartree)
269 ps_implicit_env%ehartree = ehartree
271 IF (outp_unit .GT. 0)
WRITE (outp_unit,
'(A1,/)')
275 reached_max_iter = iter .GT. max_iter
276 reached_tol = pres_error .LE. tol
277 IF (pres_error .GT. large_error) &
278 cpabort(
"Poisson solver did not converge.")
279 ps_implicit_env%times_called = ps_implicit_env%times_called + 1
280 IF (reached_max_iter .OR. reached_tol)
EXIT
287 CALL ps_implicit_print_convergence_msg(iter, max_iter, outp_unit)
289 IF ((times_called .NE. 0) .AND. (.NOT. use_zero_initial_guess)) &
290 CALL pw_copy(v_new, ps_implicit_env%initial_guess)
292 IF (
PRESENT(ehartree)) ehartree = ps_implicit_env%ehartree
296 v_eps%pw_grid => ps_implicit_env%v_eps%pw_grid
297 v_eps%array => ps_implicit_env%v_eps%array
298 CALL ps_implicit_compute_veps(pw_pool, dielectric, v_new, v_eps)
301 CALL pw_pool%give_back_pw(g)
302 CALL pw_pool%give_back_pw(v_old)
303 CALL pw_pool%give_back_pw(res_old)
304 CALL pw_pool%give_back_pw(res_new)
305 CALL pw_pool%give_back_pw(qainvxres)
306 CALL pw_pool%give_back_pw(pxqainvxres)
308 CALL timestop(handle)
326 TYPE(pw_poisson_type),
INTENT(IN) :: poisson_env
327 TYPE(pw_r3d_rs_type),
INTENT(IN) :: density
328 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: v_new
329 REAL(dp),
INTENT(OUT),
OPTIONAL :: ehartree
331 CHARACTER(LEN=*),
PARAMETER :: routinen =
'implicit_poisson_solver_neumann'
333 INTEGER :: handle, iter, max_iter, &
334 neumann_directions, outp_unit, &
336 LOGICAL :: reached_max_iter, reached_tol, &
337 use_zero_initial_guess
338 REAL(dp) :: nabs_error, omega, pres_error, tol, &
340 TYPE(dct_type),
POINTER :: dct_env
341 TYPE(dielectric_type),
POINTER :: dielectric
342 TYPE(greens_fn_type),
POINTER :: green
343 TYPE(ps_implicit_type),
POINTER :: ps_implicit_env
344 TYPE(pw_pool_type),
POINTER :: pw_pool, pw_pool_xpndd
345 TYPE(pw_r3d_rs_type) :: density_xpndd, g, pxqainvxres, &
346 qainvxres, res_new, res_old, &
347 v_eps_xpndd, v_new_xpndd, v_old
349 CALL timeset(routinen, handle)
351 pw_pool => poisson_env%pw_pools(poisson_env%pw_level)%pool
352 dielectric => poisson_env%implicit_env%dielectric
353 green => poisson_env%green_fft
354 ps_implicit_env => poisson_env%implicit_env
355 dct_env => ps_implicit_env%dct_env
357 tol = poisson_env%parameters%ps_implicit_params%tol
358 omega = poisson_env%parameters%ps_implicit_params%omega
359 max_iter = poisson_env%parameters%ps_implicit_params%max_iter
360 use_zero_initial_guess = poisson_env%parameters%ps_implicit_params%zero_initial_guess
361 neumann_directions = poisson_env%parameters%ps_implicit_params%neumann_directions
362 times_called = ps_implicit_env%times_called
364 SELECT CASE (neumann_directions)
367 CASE (neumannxy, neumannxz, neumannyz)
369 CASE (neumannx, neumanny, neumannz)
373 CALL pw_pool_create(pw_pool_xpndd, pw_grid=poisson_env%dct_pw_grid)
376 IF (times_called .EQ. 0)
CALL ps_implicit_initial_guess_create(ps_implicit_env, pw_pool_xpndd)
378 CALL pw_pool_xpndd%create_pw(g)
379 CALL pw_pool_xpndd%create_pw(v_old)
380 CALL pw_pool_xpndd%create_pw(res_old)
381 CALL pw_pool_xpndd%create_pw(res_new)
382 CALL pw_pool_xpndd%create_pw(qainvxres)
383 CALL pw_pool_xpndd%create_pw(pxqainvxres)
384 CALL pw_pool_xpndd%create_pw(density_xpndd)
385 CALL pw_pool_xpndd%create_pw(v_new_xpndd)
386 CALL pw_pool_xpndd%create_pw(v_eps_xpndd)
388 IF (use_zero_initial_guess)
THEN
391 CALL pw_copy(ps_implicit_env%initial_guess, v_old)
394 CALL pw_expand(neumann_directions, &
395 dct_env%recv_msgs_bnds, dct_env%dests_expand, dct_env%srcs_expand, &
396 dct_env%flipg_stat, dct_env%bounds_shftd, density, density_xpndd)
397 CALL pw_expand(neumann_directions, &
398 dct_env%recv_msgs_bnds, dct_env%dests_expand, dct_env%srcs_expand, &
399 dct_env%flipg_stat, dct_env%bounds_shftd, v_new, v_new_xpndd)
401 g%array = fourpi*density_xpndd%array/dielectric%eps%array
404 CALL apply_poisson_operator_dct(pw_pool_xpndd, green, dielectric, v_old, res_old)
405 CALL pw_scale(res_old, -1.0_dp)
406 CALL pw_axpy(g, res_old)
409 CALL apply_inv_laplace_operator_dct(pw_pool_xpndd, green, res_old, qainvxres)
415 CALL pw_scale(qainvxres, omega)
416 CALL pw_copy(qainvxres, v_new_xpndd)
417 CALL pw_axpy(v_old, v_new_xpndd)
421 CALL apply_p_operator(pw_pool_xpndd, dielectric, qainvxres, pxqainvxres)
422 CALL pw_copy(pxqainvxres, res_new)
423 CALL pw_scale(res_new, -1.0_dp)
424 CALL pw_axpy(res_old, res_new, 1.0_dp - omega)
427 CALL ps_implicit_compute_error_dct(pw_pool_xpndd, green, res_new, v_old, v_new_xpndd, qainvxres, &
428 pres_error, nabs_error)
430 CALL ps_implicit_output(iter, pres_error, nabs_error, outp_unit)
431 IF (
PRESENT(ehartree))
THEN
432 CALL ps_implicit_compute_ehartree(density_xpndd, v_new_xpndd, ehartree)
433 CALL ps_implicit_report_ehartree(ps_implicit_env, outp_unit, ehartree/vol_scfac)
434 ps_implicit_env%ehartree = ehartree/vol_scfac
436 IF (outp_unit .GT. 0)
WRITE (outp_unit,
'(A1,/)')
440 reached_max_iter = iter .GT. max_iter
441 reached_tol = pres_error .LE. tol
442 IF (pres_error .GT. large_error) &
443 cpabort(
"Poisson solver did not converge.")
444 ps_implicit_env%times_called = ps_implicit_env%times_called + 1
445 IF (reached_max_iter .OR. reached_tol)
EXIT
448 CALL pw_copy(v_new_xpndd, v_old)
449 CALL pw_copy(res_new, res_old)
452 CALL ps_implicit_print_convergence_msg(iter, max_iter, outp_unit)
454 CALL pw_shrink(neumann_directions, dct_env%dests_shrink, dct_env%srcs_shrink, &
455 dct_env%bounds_local_shftd, v_new_xpndd, v_new)
457 IF ((times_called .NE. 0) .AND. (.NOT. use_zero_initial_guess)) &
458 CALL pw_copy(v_new_xpndd, ps_implicit_env%initial_guess)
460 IF (
PRESENT(ehartree)) ehartree = ps_implicit_env%ehartree
463 CALL ps_implicit_compute_veps(pw_pool_xpndd, dielectric, v_new_xpndd, v_eps_xpndd)
465 TYPE(pw_r3d_rs_type) :: v_eps
466 v_eps%pw_grid => ps_implicit_env%v_eps%pw_grid
467 v_eps%array => ps_implicit_env%v_eps%array
468 CALL pw_shrink(neumann_directions, dct_env%dests_shrink, dct_env%srcs_shrink, &
469 dct_env%bounds_local_shftd, v_eps_xpndd, v_eps)
472 CALL pw_pool_xpndd%give_back_pw(g)
473 CALL pw_pool_xpndd%give_back_pw(v_old)
474 CALL pw_pool_xpndd%give_back_pw(res_old)
475 CALL pw_pool_xpndd%give_back_pw(res_new)
476 CALL pw_pool_xpndd%give_back_pw(qainvxres)
477 CALL pw_pool_xpndd%give_back_pw(pxqainvxres)
478 CALL pw_pool_xpndd%give_back_pw(density_xpndd)
479 CALL pw_pool_xpndd%give_back_pw(v_new_xpndd)
480 CALL pw_pool_xpndd%give_back_pw(v_eps_xpndd)
481 CALL pw_pool_release(pw_pool_xpndd)
483 CALL timestop(handle)
499 TYPE(pw_poisson_type),
INTENT(IN) :: poisson_env
500 TYPE(pw_r3d_rs_type),
INTENT(IN) :: density
501 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: v_new
502 REAL(dp),
INTENT(OUT),
OPTIONAL :: electric_enthalpy
504 CHARACTER(LEN=*),
PARAMETER :: routinen =
'implicit_poisson_solver_mixed_periodic'
506 INTEGER :: data_size, handle, iter, j, lb1, lb2, lb3, max_iter, n_contacts, n_tiles_tot, ng, &
507 ngpts_local, nt, nt_tot, outp_unit, times_called, ub1, ub2, ub3
508 INTEGER(KIND=int_8) :: ngpts
509 INTEGER,
DIMENSION(2, 3) :: bounds_local
510 INTEGER,
DIMENSION(3) :: npts_local
511 LOGICAL :: reached_max_iter, reached_tol, &
512 use_zero_initial_guess
513 REAL(dp) :: axvbar_avg, ehartree, eta, g_avg, &
514 gminusaxvbar_avg, nabs_error, omega, &
516 REAL(dp),
ALLOCATABLE,
DIMENSION(:) :: btxlambda_new, btxlambda_old, bxv_bar, bxv_new, &
517 lambda0, lambda_new, lambda_newneta, lambda_old, qsxlambda, v_bar1d, v_d, v_new1d, w
518 REAL(dp),
ALLOCATABLE,
DIMENSION(:, :) :: b, bt, qs, rinv
519 REAL(dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: btxlambda_new3d, btxlambda_old3d
520 TYPE(dielectric_type),
POINTER :: dielectric
521 TYPE(greens_fn_type),
POINTER :: green
522 TYPE(ps_implicit_type),
POINTER :: ps_implicit_env
523 TYPE(pw_grid_type),
POINTER :: pw_grid
524 TYPE(pw_pool_type),
POINTER :: pw_pool
525 TYPE(pw_r3d_rs_type) :: axvbar, g, pxqainvxres, qainvxres, &
526 res_new, res_old, v_old
528 CALL timeset(routinen, handle)
530 pw_pool => poisson_env%pw_pools(poisson_env%pw_level)%pool
531 pw_grid => pw_pool%pw_grid
532 dielectric => poisson_env%implicit_env%dielectric
533 green => poisson_env%green_fft
534 ps_implicit_env => poisson_env%implicit_env
536 ngpts_local = pw_grid%ngpts_local
537 ngpts = pw_grid%ngpts
538 npts_local = pw_grid%npts_local
539 bounds_local = pw_grid%bounds_local
540 tol = poisson_env%parameters%ps_implicit_params%tol
541 omega = poisson_env%parameters%ps_implicit_params%omega
542 max_iter = poisson_env%parameters%ps_implicit_params%max_iter
543 use_zero_initial_guess = poisson_env%parameters%ps_implicit_params%zero_initial_guess
544 times_called = ps_implicit_env%times_called
546 n_contacts =
SIZE(ps_implicit_env%contacts)
549 n_tiles_tot = n_tiles_tot + ps_implicit_env%contacts(j)%dirichlet_bc%n_tiles
552 IF (pw_grid%para%blocked)
THEN
553 data_size = product(npts_local)
554 ELSE IF (pw_grid%para%ray_distribution)
THEN
555 data_size = ngpts_local
557 data_size = product(npts_local)
561 IF (times_called .EQ. 0)
CALL ps_implicit_initial_guess_create(ps_implicit_env, pw_pool)
563 ALLOCATE (b(n_tiles_tot, data_size))
564 ALLOCATE (bt(data_size, n_tiles_tot))
565 ALLOCATE (qs(n_tiles_tot, n_tiles_tot))
566 ALLOCATE (rinv(n_tiles_tot + 1, n_tiles_tot + 1))
568 b(:, :) = ps_implicit_env%B
569 bt(:, :) = ps_implicit_env%Bt
570 qs(:, :) = ps_implicit_env%QS
571 rinv(:, :) = ps_implicit_env%Rinv
572 CALL get_voltage(poisson_env%parameters%dbc_params%time, ps_implicit_env%v_D, ps_implicit_env%osc_frac, &
573 ps_implicit_env%frequency, ps_implicit_env%phase, v_d)
575 lb1 = bounds_local(1, 1); ub1 = bounds_local(2, 1)
576 lb2 = bounds_local(1, 2); ub2 = bounds_local(2, 2)
577 lb3 = bounds_local(1, 3); ub3 = bounds_local(2, 3)
579 ALLOCATE (lambda0(n_tiles_tot), lambda_old(n_tiles_tot), lambda_new(n_tiles_tot))
580 ALLOCATE (btxlambda_old(data_size), btxlambda_new(data_size))
581 ALLOCATE (btxlambda_old3d(lb1:ub1, lb2:ub2, lb3:ub3), btxlambda_new3d(lb1:ub1, lb2:ub2, lb3:ub3))
582 ALLOCATE (qsxlambda(n_tiles_tot))
583 ALLOCATE (w(n_tiles_tot + 1))
584 ALLOCATE (lambda_newneta(n_tiles_tot + 1))
585 ALLOCATE (v_bar1d(data_size))
586 ALLOCATE (bxv_bar(n_tiles_tot))
588 ALLOCATE (v_new1d(data_size))
589 ALLOCATE (bxv_new(n_tiles_tot))
591 CALL pw_pool%create_pw(g)
592 CALL pw_pool%create_pw(v_old)
593 CALL pw_pool%create_pw(res_old)
594 CALL pw_pool%create_pw(res_new)
595 CALL pw_pool%create_pw(qainvxres)
596 CALL pw_pool%create_pw(pxqainvxres)
597 CALL pw_pool%create_pw(axvbar)
599 IF (use_zero_initial_guess)
THEN
603 CALL pw_copy(ps_implicit_env%initial_guess, v_old)
604 lambda0(:) = ps_implicit_env%initial_lambda
607 g%array = fourpi*density%array/dielectric%eps%array
608 g_avg = accurate_sum(g%array)/ngpts
610 lambda_old(:) = lambda0
613 CALL apply_poisson_operator_fft(pw_pool, green, dielectric, v_old, res_old)
614 CALL pw_scale(res_old, -1.0_dp)
615 CALL pw_axpy(g, res_old)
616 IF (data_size .NE. 0)
THEN
617 CALL dgemv(
'N', data_size, n_tiles_tot, 1.0_dp, bt, data_size, lambda_old, 1, 0.0_dp, btxlambda_old, 1)
619 CALL convert_1dto3d(ps_implicit_env%idx_1dto3d, btxlambda_old, btxlambda_old3d)
620 res_old%array = res_old%array - btxlambda_old3d
623 CALL apply_inv_laplace_operator_fft(pw_pool, green, res_old, qainvxres)
629 CALL pw_scale(qainvxres, omega)
630 CALL pw_copy(qainvxres, v_new)
631 CALL pw_axpy(v_old, v_new)
635 CALL apply_p_operator(pw_pool, dielectric, v_new, axvbar)
636 axvbar_avg = accurate_sum(axvbar%array)/ngpts
637 gminusaxvbar_avg = g_avg - axvbar_avg
638 CALL pw_grid%para%group%sum(gminusaxvbar_avg)
641 CALL dgemv(
'N', n_tiles_tot, n_tiles_tot, 1.0_dp, qs, n_tiles_tot, lambda_old, 1, 0.0_dp, qsxlambda, 1)
642 v_bar1d(ps_implicit_env%idx_1dto3d) = reshape(v_new%array, (/data_size/))
643 CALL dgemv(
'N', n_tiles_tot, data_size, 1.0_dp, b, n_tiles_tot, v_bar1d, 1, 0.0_dp, bxv_bar, 1)
644 CALL pw_grid%para%group%sum(bxv_bar)
647 w(:) = (/qsxlambda + v_d - bxv_bar, gminusaxvbar_avg/)
648 CALL dgemv(
'N', n_tiles_tot + 1, n_tiles_tot + 1, 1.0_dp, rinv, n_tiles_tot + 1, w, 1, 0.0_dp, lambda_newneta, 1)
649 lambda_new(:) = lambda_newneta(1:n_tiles_tot)
650 eta = lambda_newneta(n_tiles_tot + 1)
653 v_new%array = v_new%array + eta/ngpts
656 IF (data_size .NE. 0)
THEN
657 CALL dgemv(
'N', data_size, n_tiles_tot, 1.0_dp, bt, data_size, lambda_new, 1, 0.0_dp, btxlambda_new, 1)
659 CALL convert_1dto3d(ps_implicit_env%idx_1dto3d, btxlambda_new, btxlambda_new3d)
663 CALL pw_zero(res_new)
664 CALL apply_p_operator(pw_pool, dielectric, qainvxres, pxqainvxres)
665 CALL pw_axpy(pxqainvxres, res_new, -1.0_dp)
666 CALL pw_axpy(res_old, res_new, 1.0_dp - omega)
667 res_new%array = res_new%array + btxlambda_old3d - btxlambda_new3d
670 CALL ps_implicit_compute_error_fft(pw_pool, green, res_new, v_old, v_new, qainvxres, &
671 pres_error, nabs_error)
673 CALL ps_implicit_output(iter, pres_error, nabs_error, outp_unit)
674 IF (
PRESENT(electric_enthalpy))
THEN
675 CALL ps_implicit_compute_ehartree(dielectric, density, btxlambda_new3d, v_new, ehartree, electric_enthalpy)
676 CALL ps_implicit_report_ehartree(ps_implicit_env, outp_unit, ehartree)
677 ps_implicit_env%ehartree = ehartree
678 ps_implicit_env%electric_enthalpy = electric_enthalpy
680 IF (outp_unit .GT. 0)
WRITE (outp_unit,
'(A1,/)')
684 IF (poisson_env%parameters%dbc_params%verbose_output)
THEN
685 v_new1d(ps_implicit_env%idx_1dto3d) = reshape(v_new%array, (/data_size/))
686 CALL dgemv(
'N', n_tiles_tot, data_size, 1.0_dp, b, n_tiles_tot, v_new1d, 1, 0.0_dp, bxv_new, 1)
687 CALL pw_grid%para%group%sum(bxv_new)
688 IF (outp_unit .GT. 0)
THEN
689 WRITE (outp_unit,
'(T3,A,A)')
"======== verbose ", repeat(
'=', 61)
690 WRITE (outp_unit,
'(T20,A)')
"Drgn tile vhartree lambda "
691 WRITE (outp_unit,
'(T19,A)') repeat(
'-', 46)
693 DO ng = 1, n_contacts
694 DO nt = 1, ps_implicit_env%contacts(ng)%dirichlet_bc%n_tiles
695 WRITE (outp_unit,
'(T17,I6,5X,I6,3X,E13.4,E13.4)') ng, nt, bxv_new(nt_tot), lambda_new(nt_tot)
699 WRITE (outp_unit,
'(T3,A)') repeat(
'=', 78)
705 reached_max_iter = iter .GT. max_iter
706 reached_tol = pres_error .LE. tol
707 ps_implicit_env%times_called = ps_implicit_env%times_called + 1
708 IF (pres_error .GT. large_error) &
709 cpabort(
"Poisson solver did not converge.")
710 IF (reached_max_iter .OR. reached_tol)
EXIT
713 CALL pw_copy(v_new, v_old)
714 lambda_old(:) = lambda_new
715 CALL pw_copy(res_new, res_old)
716 btxlambda_old3d(:, :, :) = btxlambda_new3d
719 CALL ps_implicit_print_convergence_msg(iter, max_iter, outp_unit)
721 IF ((times_called .NE. 0) .AND. (.NOT. use_zero_initial_guess))
THEN
722 CALL pw_copy(v_new, ps_implicit_env%initial_guess)
723 ps_implicit_env%initial_lambda(:) = lambda_new
726 ps_implicit_env%cstr_charge%array = btxlambda_new3d
727 IF (
PRESENT(electric_enthalpy)) electric_enthalpy = ps_implicit_env%electric_enthalpy
730 TYPE(pw_r3d_rs_type) :: tmp
731 tmp%pw_grid => ps_implicit_env%v_eps%pw_grid
732 tmp%array => ps_implicit_env%v_eps%array
733 CALL ps_implicit_compute_veps(pw_pool, dielectric, v_new, tmp)
736 CALL pw_pool%give_back_pw(g)
737 CALL pw_pool%give_back_pw(v_old)
738 CALL pw_pool%give_back_pw(res_old)
739 CALL pw_pool%give_back_pw(res_new)
740 CALL pw_pool%give_back_pw(qainvxres)
741 CALL pw_pool%give_back_pw(pxqainvxres)
742 CALL pw_pool%give_back_pw(axvbar)
744 CALL timestop(handle)
761 TYPE(pw_poisson_type),
INTENT(IN) :: poisson_env
762 TYPE(pw_r3d_rs_type),
INTENT(IN) :: density
763 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: v_new
764 REAL(dp),
INTENT(OUT),
OPTIONAL :: electric_enthalpy
766 CHARACTER(LEN=*),
PARAMETER :: routinen =
'implicit_poisson_solver_mixed'
768 INTEGER :: data_size, handle, iter, j, lb1, lb2, lb3, max_iter, n_contacts, n_tiles_tot, &
769 neumann_directions, ng, ngpts_local, nt, nt_tot, outp_unit, times_called, ub1, ub2, ub3
770 INTEGER(KIND=int_8) :: ngpts
771 INTEGER,
DIMENSION(2, 3) :: bounds_local
772 INTEGER,
DIMENSION(3) :: npts_local
773 LOGICAL :: reached_max_iter, reached_tol, &
774 use_zero_initial_guess
775 REAL(dp) :: axvbar_avg, ehartree, eta, g_avg, &
776 gminusaxvbar_avg, nabs_error, omega, &
777 pres_error, tol, vol_scfac
778 REAL(dp),
ALLOCATABLE,
DIMENSION(:) :: btxlambda_new, btxlambda_old, bxv_bar, bxv_new, &
779 lambda0, lambda_new, lambda_newneta, lambda_old, qsxlambda, v_bar1d, v_d, v_new1d, w
780 REAL(dp),
ALLOCATABLE,
DIMENSION(:, :) :: b, bt, qs, rinv
781 REAL(dp),
ALLOCATABLE,
DIMENSION(:, :, :) :: btxlambda_new3d, btxlambda_old3d
782 TYPE(dct_type),
POINTER :: dct_env
783 TYPE(dielectric_type),
POINTER :: dielectric
784 TYPE(greens_fn_type),
POINTER :: green
785 TYPE(ps_implicit_type),
POINTER :: ps_implicit_env
786 TYPE(pw_grid_type),
POINTER :: dct_pw_grid, pw_grid
787 TYPE(pw_pool_type),
POINTER :: pw_pool, pw_pool_xpndd
788 TYPE(pw_r3d_rs_type) :: axvbar, density_xpndd, g, pxqainvxres, &
789 qainvxres, res_new, res_old, &
790 v_eps_xpndd, v_new_xpndd, v_old
792 CALL timeset(routinen, handle)
794 pw_pool => poisson_env%pw_pools(poisson_env%pw_level)%pool
795 pw_grid => pw_pool%pw_grid
796 dielectric => poisson_env%implicit_env%dielectric
797 green => poisson_env%green_fft
798 ps_implicit_env => poisson_env%implicit_env
799 dct_env => ps_implicit_env%dct_env
801 dct_pw_grid => poisson_env%dct_pw_grid
802 ngpts_local = dct_pw_grid%ngpts_local
803 ngpts = dct_pw_grid%ngpts
804 npts_local = dct_pw_grid%npts_local
805 bounds_local = dct_pw_grid%bounds_local
806 tol = poisson_env%parameters%ps_implicit_params%tol
807 omega = poisson_env%parameters%ps_implicit_params%omega
808 max_iter = poisson_env%parameters%ps_implicit_params%max_iter
809 use_zero_initial_guess = poisson_env%parameters%ps_implicit_params%zero_initial_guess
810 neumann_directions = poisson_env%parameters%ps_implicit_params%neumann_directions
811 times_called = ps_implicit_env%times_called
813 SELECT CASE (neumann_directions)
816 CASE (neumannxy, neumannxz, neumannyz)
818 CASE (neumannx, neumanny, neumannz)
822 n_contacts =
SIZE(ps_implicit_env%contacts)
825 n_tiles_tot = n_tiles_tot + ps_implicit_env%contacts(j)%dirichlet_bc%n_tiles
828 IF (dct_pw_grid%para%blocked)
THEN
829 data_size = product(npts_local)
830 ELSE IF (dct_pw_grid%para%ray_distribution)
THEN
831 data_size = ngpts_local
833 data_size = product(npts_local)
836 CALL pw_pool_create(pw_pool_xpndd, pw_grid=dct_pw_grid)
839 IF (times_called .EQ. 0)
CALL ps_implicit_initial_guess_create(ps_implicit_env, pw_pool_xpndd)
841 ALLOCATE (b(n_tiles_tot, data_size))
842 ALLOCATE (bt(data_size, n_tiles_tot))
843 ALLOCATE (qs(n_tiles_tot, n_tiles_tot))
844 ALLOCATE (rinv(n_tiles_tot + 1, n_tiles_tot + 1))
846 b(:, :) = ps_implicit_env%B
847 bt(:, :) = ps_implicit_env%Bt
848 qs(:, :) = ps_implicit_env%QS
849 rinv(:, :) = ps_implicit_env%Rinv
850 CALL get_voltage(poisson_env%parameters%dbc_params%time, ps_implicit_env%v_D, ps_implicit_env%osc_frac, &
851 ps_implicit_env%frequency, ps_implicit_env%phase, v_d)
853 lb1 = bounds_local(1, 1); ub1 = bounds_local(2, 1)
854 lb2 = bounds_local(1, 2); ub2 = bounds_local(2, 2)
855 lb3 = bounds_local(1, 3); ub3 = bounds_local(2, 3)
857 ALLOCATE (lambda0(n_tiles_tot), lambda_old(n_tiles_tot), lambda_new(n_tiles_tot))
858 ALLOCATE (btxlambda_old(data_size), btxlambda_new(data_size))
859 ALLOCATE (btxlambda_old3d(lb1:ub1, lb2:ub2, lb3:ub3), btxlambda_new3d(lb1:ub1, lb2:ub2, lb3:ub3))
860 ALLOCATE (qsxlambda(n_tiles_tot))
861 ALLOCATE (w(n_tiles_tot + 1))
862 ALLOCATE (lambda_newneta(n_tiles_tot + 1))
863 ALLOCATE (v_bar1d(data_size))
864 ALLOCATE (bxv_bar(n_tiles_tot))
866 ALLOCATE (v_new1d(data_size))
867 ALLOCATE (bxv_new(n_tiles_tot))
869 CALL pw_pool_xpndd%create_pw(g)
870 CALL pw_pool_xpndd%create_pw(v_old)
871 CALL pw_pool_xpndd%create_pw(res_old)
872 CALL pw_pool_xpndd%create_pw(res_new)
873 CALL pw_pool_xpndd%create_pw(qainvxres)
874 CALL pw_pool_xpndd%create_pw(pxqainvxres)
875 CALL pw_pool_xpndd%create_pw(axvbar)
876 CALL pw_pool_xpndd%create_pw(density_xpndd)
877 CALL pw_pool_xpndd%create_pw(v_new_xpndd)
878 CALL pw_pool_xpndd%create_pw(v_eps_xpndd)
880 IF (use_zero_initial_guess)
THEN
884 CALL pw_copy(ps_implicit_env%initial_guess, v_old)
885 lambda0(:) = ps_implicit_env%initial_lambda
888 CALL pw_expand(neumann_directions, &
889 dct_env%recv_msgs_bnds, dct_env%dests_expand, dct_env%srcs_expand, &
890 dct_env%flipg_stat, dct_env%bounds_shftd, density, density_xpndd)
891 CALL pw_expand(neumann_directions, &
892 dct_env%recv_msgs_bnds, dct_env%dests_expand, dct_env%srcs_expand, &
893 dct_env%flipg_stat, dct_env%bounds_shftd, v_new, v_new_xpndd)
895 g%array = fourpi*density_xpndd%array/dielectric%eps%array
896 g_avg = accurate_sum(g%array)/ngpts
898 lambda_old(:) = lambda0
901 CALL apply_poisson_operator_dct(pw_pool_xpndd, green, dielectric, v_old, res_old)
902 CALL pw_scale(res_old, -1.0_dp)
903 CALL pw_axpy(g, res_old)
904 IF (data_size .NE. 0)
THEN
905 CALL dgemv(
'N', data_size, n_tiles_tot, 1.0_dp, bt, data_size, lambda_old, 1, 0.0_dp, btxlambda_old, 1)
907 CALL convert_1dto3d(ps_implicit_env%idx_1dto3d, btxlambda_old, btxlambda_old3d)
908 res_old%array = res_old%array - btxlambda_old3d
911 CALL apply_inv_laplace_operator_dct(pw_pool_xpndd, green, res_old, qainvxres)
917 CALL pw_scale(qainvxres, omega)
918 CALL pw_copy(qainvxres, v_new_xpndd)
919 CALL pw_axpy(v_old, v_new_xpndd)
923 CALL apply_p_operator(pw_pool_xpndd, dielectric, v_new_xpndd, axvbar)
924 axvbar_avg = accurate_sum(axvbar%array)/ngpts
925 gminusaxvbar_avg = g_avg - axvbar_avg
926 CALL dct_pw_grid%para%group%sum(gminusaxvbar_avg)
929 CALL dgemv(
'N', n_tiles_tot, n_tiles_tot, 1.0_dp, qs, n_tiles_tot, lambda_old, 1, 0.0_dp, qsxlambda, 1)
930 v_bar1d(ps_implicit_env%idx_1dto3d) = reshape(v_new_xpndd%array, (/data_size/))
931 CALL dgemv(
'N', n_tiles_tot, data_size, 1.0_dp, b, n_tiles_tot, v_bar1d, 1, 0.0_dp, bxv_bar, 1)
932 CALL dct_pw_grid%para%group%sum(bxv_bar)
935 w(:) = (/qsxlambda + v_d - bxv_bar, gminusaxvbar_avg/)
936 CALL dgemv(
'N', n_tiles_tot + 1, n_tiles_tot + 1, 1.0_dp, rinv, n_tiles_tot + 1, w, 1, 0.0_dp, lambda_newneta, 1)
937 lambda_new(:) = lambda_newneta(1:n_tiles_tot)
938 eta = lambda_newneta(n_tiles_tot + 1)
941 v_new_xpndd%array = v_new_xpndd%array + eta/ngpts
944 IF (data_size .NE. 0)
THEN
945 CALL dgemv(
'N', data_size, n_tiles_tot, 1.0_dp, bt, data_size, lambda_new, 1, 0.0_dp, btxlambda_new, 1)
947 CALL convert_1dto3d(ps_implicit_env%idx_1dto3d, btxlambda_new, btxlambda_new3d)
951 CALL pw_zero(res_new)
952 CALL apply_p_operator(pw_pool_xpndd, dielectric, qainvxres, pxqainvxres)
953 CALL pw_axpy(pxqainvxres, res_new, -1.0_dp)
954 CALL pw_axpy(res_old, res_new, 1.0_dp - omega)
955 res_new%array = res_new%array - btxlambda_new3d + btxlambda_old3d
958 CALL ps_implicit_compute_error_dct(pw_pool_xpndd, green, res_new, v_old, v_new_xpndd, qainvxres, &
959 pres_error, nabs_error)
961 CALL ps_implicit_output(iter, pres_error, nabs_error, outp_unit)
962 IF (
PRESENT(electric_enthalpy))
THEN
963 CALL ps_implicit_compute_ehartree(dielectric, density_xpndd, btxlambda_new3d, v_new_xpndd, &
964 ehartree, electric_enthalpy)
965 CALL ps_implicit_report_ehartree(ps_implicit_env, outp_unit, ehartree/vol_scfac)
966 ps_implicit_env%ehartree = ehartree/vol_scfac
967 ps_implicit_env%electric_enthalpy = electric_enthalpy/vol_scfac
969 IF (outp_unit .GT. 0)
WRITE (outp_unit,
'(A1,/)')
973 IF (poisson_env%parameters%dbc_params%verbose_output)
THEN
974 v_new1d(ps_implicit_env%idx_1dto3d) = reshape(v_new_xpndd%array, (/data_size/))
975 CALL dgemv(
'N', n_tiles_tot, data_size, 1.0_dp, b, n_tiles_tot, v_new1d, 1, 0.0_dp, bxv_new, 1)
976 CALL pw_grid%para%group%sum(bxv_new)
977 IF (outp_unit .GT. 0)
THEN
978 WRITE (outp_unit,
'(T3,A)')
"======== verbose "//repeat(
'=', 61)
979 WRITE (outp_unit,
'(T20,A)')
"Drgn tile vhartree lambda "
980 WRITE (outp_unit,
'(T19,A)') repeat(
'-', 46)
982 DO ng = 1, n_contacts
983 DO nt = 1, ps_implicit_env%contacts(ng)%dirichlet_bc%n_tiles
984 WRITE (outp_unit,
'(T17,I6,5X,I6,3X,E13.4,E13.4)') ng, nt, bxv_new(nt_tot), lambda_new(nt_tot)
988 WRITE (outp_unit,
'(T3,A)') repeat(
'=', 78)
994 reached_max_iter = iter .GT. max_iter
995 reached_tol = pres_error .LE. tol
996 ps_implicit_env%times_called = ps_implicit_env%times_called + 1
997 IF (pres_error .GT. large_error) &
998 cpabort(
"Poisson solver did not converge.")
999 IF (reached_max_iter .OR. reached_tol)
EXIT
1002 CALL pw_copy(v_new_xpndd, v_old)
1003 lambda_old(:) = lambda_new
1004 CALL pw_copy(res_new, res_old)
1005 btxlambda_old3d(:, :, :) = btxlambda_new3d
1008 CALL ps_implicit_print_convergence_msg(iter, max_iter, outp_unit)
1010 CALL pw_shrink(neumann_directions, dct_env%dests_shrink, dct_env%srcs_shrink, &
1011 dct_env%bounds_local_shftd, v_new_xpndd, v_new)
1013 IF ((times_called .NE. 0) .AND. (.NOT. use_zero_initial_guess))
THEN
1014 CALL pw_copy(v_new_xpndd, ps_implicit_env%initial_guess)
1015 ps_implicit_env%initial_lambda(:) = lambda_new
1018 ps_implicit_env%cstr_charge%array = btxlambda_new3d
1019 IF (
PRESENT(electric_enthalpy)) electric_enthalpy = ps_implicit_env%electric_enthalpy
1021 CALL ps_implicit_compute_veps(pw_pool_xpndd, dielectric, v_new_xpndd, v_eps_xpndd)
1022 CALL pw_shrink(neumann_directions, dct_env%dests_shrink, dct_env%srcs_shrink, &
1023 dct_env%bounds_local_shftd, v_eps_xpndd, ps_implicit_env%v_eps)
1025 CALL pw_pool_xpndd%give_back_pw(g)
1026 CALL pw_pool_xpndd%give_back_pw(v_old)
1027 CALL pw_pool_xpndd%give_back_pw(res_old)
1028 CALL pw_pool_xpndd%give_back_pw(res_new)
1029 CALL pw_pool_xpndd%give_back_pw(qainvxres)
1030 CALL pw_pool_xpndd%give_back_pw(pxqainvxres)
1031 CALL pw_pool_xpndd%give_back_pw(axvbar)
1032 CALL pw_pool_xpndd%give_back_pw(density_xpndd)
1033 CALL pw_pool_xpndd%give_back_pw(v_new_xpndd)
1034 CALL pw_pool_xpndd%give_back_pw(v_eps_xpndd)
1035 CALL pw_pool_release(pw_pool_xpndd)
1037 CALL timestop(handle)
1049 SUBROUTINE ps_implicit_initial_guess_create(ps_implicit_env, pw_pool)
1051 TYPE(ps_implicit_type),
INTENT(INOUT),
POINTER :: ps_implicit_env
1052 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1054 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_initial_guess_create'
1056 INTEGER :: handle, n_tiles_tot
1058 CALL timeset(routinen, handle)
1060 n_tiles_tot =
SIZE(ps_implicit_env%v_D)
1061 NULLIFY (ps_implicit_env%initial_guess)
1062 ALLOCATE (ps_implicit_env%initial_guess)
1063 CALL pw_pool%create_pw(ps_implicit_env%initial_guess)
1064 CALL pw_zero(ps_implicit_env%initial_guess)
1065 ALLOCATE (ps_implicit_env%initial_lambda(n_tiles_tot))
1066 ps_implicit_env%initial_lambda = 0.0_dp
1068 CALL timestop(handle)
1070 END SUBROUTINE ps_implicit_initial_guess_create
1083 SUBROUTINE ps_implicit_prepare_blocks(pw_pool_orig, dct_pw_grid, green, &
1084 poisson_params, ps_implicit_env)
1086 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool_orig
1087 TYPE(pw_grid_type),
INTENT(IN),
POINTER :: dct_pw_grid
1088 TYPE(greens_fn_type),
INTENT(IN) :: green
1089 TYPE(pw_poisson_parameter_type),
INTENT(IN) :: poisson_params
1090 TYPE(ps_implicit_type),
INTENT(INOUT),
POINTER :: ps_implicit_env
1092 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_prepare_blocks'
1094 INTEGER :: data_size, handle, i, indx1, indx2, info, j, k, l, lb1, lb2, lb3, n_contacts, &
1095 n_tiles, n_tiles_tot, neumann_directions, ngpts_local, ub1, ub2, ub3, unit_nr
1096 INTEGER(KIND=int_8) :: ngpts
1097 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: ipiv
1098 INTEGER,
DIMENSION(2, 3) :: bounds, bounds_local
1099 INTEGER,
DIMENSION(3) :: npts, npts_local
1100 LOGICAL :: done_preparing
1101 REAL(dp) :: tile_volume, vol_scfac
1102 REAL(dp),
ALLOCATABLE,
DIMENSION(:) :: bxunit_vec, test_vec, work_arr
1103 REAL(dp),
ALLOCATABLE,
DIMENSION(:, :) :: qainvxbt, r
1104 TYPE(cp_logger_type),
POINTER :: logger
1105 TYPE(dct_type),
POINTER :: dct_env
1106 TYPE(pw_grid_type),
POINTER :: pw_grid_orig
1107 TYPE(pw_pool_type),
POINTER :: pw_pool_xpndd
1108 TYPE(pw_r3d_rs_type) :: pw_in, pw_out
1110 CALL timeset(routinen, handle)
1112 pw_grid_orig => pw_pool_orig%pw_grid
1114 logger => cp_get_default_logger()
1115 IF (logger%para_env%is_source())
THEN
1116 unit_nr = cp_logger_get_default_unit_nr(logger, local=.true.)
1121 SELECT CASE (poisson_params%ps_implicit_params%boundary_condition)
1124 ngpts_local = dct_pw_grid%ngpts_local
1125 ngpts = dct_pw_grid%ngpts
1126 npts_local = dct_pw_grid%npts_local
1127 npts = dct_pw_grid%npts
1128 bounds_local = dct_pw_grid%bounds_local
1129 bounds = dct_pw_grid%bounds
1130 dct_env => ps_implicit_env%dct_env
1132 neumann_directions = poisson_params%ps_implicit_params%neumann_directions
1134 SELECT CASE (neumann_directions)
1137 CASE (neumannxy, neumannxz, neumannyz)
1139 CASE (neumannx, neumanny, neumannz)
1144 lb1 = bounds_local(1, 1); ub1 = bounds_local(2, 1)
1145 lb2 = bounds_local(1, 2); ub2 = bounds_local(2, 2)
1146 lb3 = bounds_local(1, 3); ub3 = bounds_local(2, 3)
1148 IF (dct_pw_grid%para%blocked)
THEN
1149 data_size = product(npts_local)
1150 ELSE IF (dct_pw_grid%para%ray_distribution)
THEN
1151 data_size = ngpts_local
1153 data_size = product(npts_local)
1156 ALLOCATE (ps_implicit_env%idx_1dto3d(data_size))
1166 ps_implicit_env%idx_1dto3d(l) = (i - lb1 + 1) + &
1167 (j - lb2)*npts_local(1) + &
1168 (k - lb3)*npts_local(1)*npts_local(2)
1176 n_contacts =
SIZE(ps_implicit_env%contacts)
1178 DO j = 1, n_contacts
1179 n_tiles_tot = n_tiles_tot + ps_implicit_env%contacts(j)%dirichlet_bc%n_tiles
1182 ALLOCATE (ps_implicit_env%B(n_tiles_tot, data_size))
1183 ALLOCATE (ps_implicit_env%Bt(data_size, n_tiles_tot))
1184 ALLOCATE (ps_implicit_env%QS(n_tiles_tot, n_tiles_tot))
1185 ALLOCATE (ps_implicit_env%Rinv(n_tiles_tot + 1, n_tiles_tot + 1))
1186 ALLOCATE (ps_implicit_env%v_D(n_tiles_tot))
1187 ALLOCATE (ps_implicit_env%osc_frac(n_tiles_tot))
1188 ALLOCATE (ps_implicit_env%frequency(n_tiles_tot))
1189 ALLOCATE (ps_implicit_env%phase(n_tiles_tot))
1191 ALLOCATE (qainvxbt(data_size, n_tiles_tot))
1192 ALLOCATE (bxunit_vec(n_tiles_tot))
1193 ALLOCATE (test_vec(n_tiles_tot))
1194 ALLOCATE (r(n_tiles_tot + 1, n_tiles_tot + 1))
1195 ALLOCATE (work_arr(n_tiles_tot + 1), ipiv(n_tiles_tot + 1))
1198 CALL pw_pool_create(pw_pool_xpndd, pw_grid=dct_pw_grid)
1202 DO j = 1, n_contacts
1203 n_tiles = ps_implicit_env%contacts(j)%dirichlet_bc%n_tiles
1204 indx2 = indx1 + n_tiles - 1
1207 CALL pw_pool_xpndd%create_pw(pw_in)
1208 CALL pw_expand(neumann_directions, &
1209 dct_env%recv_msgs_bnds, dct_env%dests_expand, dct_env%srcs_expand, &
1210 dct_env%flipg_stat, dct_env%bounds_shftd, &
1211 ps_implicit_env%contacts(j)%dirichlet_bc%tiles(i)%tile%tile_pw, pw_in)
1213 tile_volume = ps_implicit_env%contacts(j)%dirichlet_bc%tiles(i)%tile%volume
1214 CALL pw_scale(pw_in, 1.0_dp/(vol_scfac*tile_volume))
1215 ps_implicit_env%Bt(ps_implicit_env%idx_1dto3d, indx1 + i - 1) = reshape(pw_in%array, (/data_size/))
1217 CALL pw_pool_xpndd%create_pw(pw_out)
1218 CALL apply_inv_laplace_operator_dct(pw_pool_xpndd, green, pw_in, pw_out)
1219 qainvxbt(ps_implicit_env%idx_1dto3d, indx1 + i - 1) = reshape(pw_out%array, (/data_size/))
1221 ps_implicit_env%v_D(indx1 + i - 1) = -1.0_dp*ps_implicit_env%contacts(j)%dirichlet_bc%v_D
1222 ps_implicit_env%osc_frac(indx1 + i - 1) = ps_implicit_env%contacts(j)%dirichlet_bc%osc_frac
1223 ps_implicit_env%frequency(indx1 + i - 1) = ps_implicit_env%contacts(j)%dirichlet_bc%frequency
1224 ps_implicit_env%phase(indx1 + i - 1) = ps_implicit_env%contacts(j)%dirichlet_bc%phase
1226 CALL pw_pool_xpndd%give_back_pw(pw_in)
1227 CALL pw_pool_xpndd%give_back_pw(pw_out)
1231 ps_implicit_env%B(:, :) = transpose(ps_implicit_env%Bt)
1234 IF (data_size .NE. 0)
THEN
1235 CALL dgemm(
'N',
'N', n_tiles_tot, n_tiles_tot, data_size, &
1236 -1.0_dp, ps_implicit_env%B, n_tiles_tot, qainvxbt, &
1237 data_size, 0.0_dp, ps_implicit_env%QS, n_tiles_tot)
1239 CALL pw_grid_orig%para%group%sum(ps_implicit_env%QS)
1242 bxunit_vec(:) = sum(ps_implicit_env%B, 2)/ngpts
1243 CALL pw_grid_orig%para%group%sum(bxunit_vec)
1246 r(1:n_tiles_tot, 1:n_tiles_tot) = ps_implicit_env%QS
1247 r(1:n_tiles_tot, n_tiles_tot + 1) = bxunit_vec
1248 r(n_tiles_tot + 1, 1:n_tiles_tot) = bxunit_vec
1250 ps_implicit_env%Rinv(:, :) = r
1251 CALL dgetrf(n_tiles_tot + 1, n_tiles_tot + 1, ps_implicit_env%Rinv, n_tiles_tot + 1, ipiv, info)
1253 CALL cp_abort(__location__, &
1254 "R is (nearly) singular! Either two Dirichlet constraints are identical or "// &
1255 "you need to reduce the number of tiles.")
1256 CALL dgetri(n_tiles_tot + 1, ps_implicit_env%Rinv, n_tiles_tot + 1, ipiv, work_arr, n_tiles_tot + 1, info)
1258 cpabort(
"Inversion of R failed!")
1260 DEALLOCATE (qainvxbt, bxunit_vec, r, work_arr, ipiv)
1261 CALL pw_pool_release(pw_pool_xpndd)
1263 done_preparing = .true.
1264 CALL pw_grid_orig%para%group%sum(done_preparing)
1265 IF ((unit_nr .GT. 0) .AND. done_preparing)
THEN
1266 WRITE (unit_nr,
"(T3,A,/,T3,A,/,A)")
"POISSON| ... Done. ", repeat(
'-', 78)
1269 CASE (mixed_periodic_bc)
1271 ngpts_local = pw_grid_orig%ngpts_local
1272 ngpts = pw_grid_orig%ngpts
1273 npts_local = pw_grid_orig%npts_local
1274 npts = pw_grid_orig%npts
1275 bounds_local = pw_grid_orig%bounds_local
1276 bounds = pw_grid_orig%bounds
1277 dct_env => ps_implicit_env%dct_env
1280 lb1 = bounds_local(1, 1); ub1 = bounds_local(2, 1)
1281 lb2 = bounds_local(1, 2); ub2 = bounds_local(2, 2)
1282 lb3 = bounds_local(1, 3); ub3 = bounds_local(2, 3)
1284 IF (pw_grid_orig%para%blocked)
THEN
1285 data_size = product(npts_local)
1286 ELSE IF (pw_grid_orig%para%ray_distribution)
THEN
1287 data_size = ngpts_local
1289 data_size = product(npts_local)
1292 ALLOCATE (ps_implicit_env%idx_1dto3d(data_size))
1302 ps_implicit_env%idx_1dto3d(l) = (i - lb1 + 1) + &
1303 (j - lb2)*npts_local(1) + &
1304 (k - lb3)*npts_local(1)*npts_local(2)
1312 n_contacts =
SIZE(ps_implicit_env%contacts)
1314 DO j = 1, n_contacts
1315 n_tiles_tot = n_tiles_tot + ps_implicit_env%contacts(j)%dirichlet_bc%n_tiles
1318 ALLOCATE (ps_implicit_env%B(n_tiles_tot, data_size))
1319 ALLOCATE (ps_implicit_env%Bt(data_size, n_tiles_tot))
1320 ALLOCATE (ps_implicit_env%QS(n_tiles_tot, n_tiles_tot))
1321 ALLOCATE (ps_implicit_env%Rinv(n_tiles_tot + 1, n_tiles_tot + 1))
1322 ALLOCATE (ps_implicit_env%v_D(n_tiles_tot))
1323 ALLOCATE (ps_implicit_env%osc_frac(n_tiles_tot))
1324 ALLOCATE (ps_implicit_env%frequency(n_tiles_tot))
1325 ALLOCATE (ps_implicit_env%phase(n_tiles_tot))
1327 ALLOCATE (qainvxbt(data_size, n_tiles_tot))
1328 ALLOCATE (bxunit_vec(n_tiles_tot))
1329 ALLOCATE (test_vec(n_tiles_tot))
1330 ALLOCATE (r(n_tiles_tot + 1, n_tiles_tot + 1))
1331 ALLOCATE (work_arr(n_tiles_tot + 1), ipiv(n_tiles_tot + 1))
1335 DO j = 1, n_contacts
1336 n_tiles = ps_implicit_env%contacts(j)%dirichlet_bc%n_tiles
1337 indx2 = indx1 + n_tiles - 1
1339 CALL pw_pool_orig%create_pw(pw_in)
1340 CALL pw_copy(ps_implicit_env%contacts(j)%dirichlet_bc%tiles(i)%tile%tile_pw, pw_in)
1342 tile_volume = ps_implicit_env%contacts(j)%dirichlet_bc%tiles(i)%tile%volume
1343 CALL pw_scale(pw_in, 1.0_dp/tile_volume)
1344 ps_implicit_env%Bt(ps_implicit_env%idx_1dto3d, indx1 + i - 1) = reshape(pw_in%array, (/data_size/))
1346 CALL pw_pool_orig%create_pw(pw_out)
1347 CALL apply_inv_laplace_operator_fft(pw_pool_orig, green, pw_in, pw_out)
1348 qainvxbt(ps_implicit_env%idx_1dto3d, indx1 + i - 1) = reshape(pw_out%array, (/data_size/))
1350 ps_implicit_env%v_D(indx1 + i - 1) = -1.0_dp*ps_implicit_env%contacts(j)%dirichlet_bc%v_D
1351 ps_implicit_env%osc_frac(indx1 + i - 1) = ps_implicit_env%contacts(j)%dirichlet_bc%osc_frac
1352 ps_implicit_env%frequency(indx1 + i - 1) = ps_implicit_env%contacts(j)%dirichlet_bc%frequency
1353 ps_implicit_env%phase(indx1 + i - 1) = ps_implicit_env%contacts(j)%dirichlet_bc%phase
1355 CALL pw_pool_orig%give_back_pw(pw_in)
1356 CALL pw_pool_orig%give_back_pw(pw_out)
1360 ps_implicit_env%B(:, :) = transpose(ps_implicit_env%Bt)
1363 IF (data_size .NE. 0)
THEN
1364 CALL dgemm(
'N',
'N', n_tiles_tot, n_tiles_tot, data_size, &
1365 -1.0_dp, ps_implicit_env%B, n_tiles_tot, qainvxbt, &
1366 data_size, 0.0_dp, ps_implicit_env%QS, n_tiles_tot)
1368 CALL pw_grid_orig%para%group%sum(ps_implicit_env%QS)
1371 bxunit_vec(:) = sum(ps_implicit_env%B, 2)/ngpts
1372 CALL pw_grid_orig%para%group%sum(bxunit_vec)
1375 r(1:n_tiles_tot, 1:n_tiles_tot) = ps_implicit_env%QS
1376 r(1:n_tiles_tot, n_tiles_tot + 1) = bxunit_vec
1377 r(n_tiles_tot + 1, 1:n_tiles_tot) = bxunit_vec
1379 ps_implicit_env%Rinv(:, :) = r
1380 CALL dgetrf(n_tiles_tot + 1, n_tiles_tot + 1, ps_implicit_env%Rinv, n_tiles_tot + 1, ipiv, info)
1382 CALL cp_abort(__location__, &
1383 "R is (nearly) singular! Either two Dirichlet constraints are identical or "// &
1384 "you need to reduce the number of tiles.")
1385 CALL dgetri(n_tiles_tot + 1, ps_implicit_env%Rinv, n_tiles_tot + 1, ipiv, work_arr, n_tiles_tot + 1, info)
1387 cpabort(
"Inversion of R failed!")
1389 DEALLOCATE (qainvxbt, bxunit_vec, r, work_arr, ipiv)
1391 done_preparing = .true.
1392 CALL pw_grid_orig%para%group%sum(done_preparing)
1393 IF ((unit_nr .GT. 0) .AND. done_preparing)
THEN
1394 WRITE (unit_nr,
"(T3,A,/,T3,A,/,A)")
"POISSON| ... Done. ", repeat(
'-', 78)
1397 CASE (periodic_bc, neumann_bc)
1399 ALLOCATE (ps_implicit_env%idx_1dto3d(1))
1400 ALLOCATE (ps_implicit_env%B(1, 1))
1401 ALLOCATE (ps_implicit_env%Bt(1, 1))
1402 ALLOCATE (ps_implicit_env%QS(1, 1))
1403 ALLOCATE (ps_implicit_env%Rinv(1, 1))
1404 ALLOCATE (ps_implicit_env%v_D(1))
1405 ALLOCATE (ps_implicit_env%osc_frac(1))
1406 ALLOCATE (ps_implicit_env%frequency(1))
1407 ALLOCATE (ps_implicit_env%phase(1))
1409 ps_implicit_env%idx_1dto3d = 1
1410 ps_implicit_env%B = 0.0_dp
1411 ps_implicit_env%Bt = 0.0_dp
1412 ps_implicit_env%QS = 0.0_dp
1413 ps_implicit_env%Rinv = 0.0_dp
1414 ps_implicit_env%v_D = 0.0_dp
1417 CALL cp_abort(__location__, &
1418 "Please specify the type of boundary conditions using the "// &
1419 "input file keyword BOUNDARY_CONDITIONS.")
1422 CALL timestop(handle)
1424 END SUBROUTINE ps_implicit_prepare_blocks
1437 SUBROUTINE apply_p_operator(pw_pool, dielectric, v, Pxv)
1439 TYPE(pw_pool_type),
POINTER :: pw_pool
1440 TYPE(dielectric_type),
INTENT(IN),
POINTER :: dielectric
1441 TYPE(pw_r3d_rs_type),
INTENT(IN) :: v
1442 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: pxv
1444 CHARACTER(LEN=*),
PARAMETER :: routinen =
'apply_P_operator'
1446 INTEGER :: handle, i
1447 TYPE(pw_r3d_rs_type),
DIMENSION(3) :: dv
1449 CALL timeset(routinen, handle)
1452 CALL pw_pool%create_pw(dv(i))
1455 CALL derive_fft(v, dv, pw_pool)
1456 associate(dln_eps => dielectric%dln_eps)
1457 pxv%array = -(dv(1)%array*dln_eps(1)%array + &
1458 dv(2)%array*dln_eps(2)%array + &
1459 dv(3)%array*dln_eps(3)%array)
1463 CALL pw_pool%give_back_pw(dv(i))
1466 CALL timestop(handle)
1468 END SUBROUTINE apply_p_operator
1480 SUBROUTINE apply_inv_laplace_operator_fft(pw_pool, green, pw_in, pw_out)
1482 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1483 TYPE(greens_fn_type),
INTENT(IN) :: green
1484 TYPE(pw_r3d_rs_type),
INTENT(IN) :: pw_in
1485 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: pw_out
1487 CHARACTER(LEN=*),
PARAMETER :: routinen =
'apply_inv_laplace_operator_fft'
1489 INTEGER :: handle, ig, ng
1490 REAL(dp) :: prefactor
1491 TYPE(pw_c1d_gs_type) :: pw_in_gs
1492 TYPE(pw_grid_type),
POINTER :: pw_grid
1494 CALL timeset(routinen, handle)
1497 prefactor = 1.0_dp/fourpi
1499 pw_grid => pw_pool%pw_grid
1500 ng =
SIZE(pw_grid%gsq)
1502 CALL pw_pool%create_pw(pw_in_gs)
1504 CALL pw_transfer(pw_in, pw_in_gs)
1506 pw_in_gs%array(ig) = prefactor*pw_in_gs%array(ig)*green%influence_fn%array(ig)
1508 CALL pw_transfer(pw_in_gs, pw_out)
1510 CALL pw_pool%give_back_pw(pw_in_gs)
1512 CALL timestop(handle)
1514 END SUBROUTINE apply_inv_laplace_operator_fft
1528 SUBROUTINE apply_inv_laplace_operator_dct(pw_pool, green, pw_in, pw_out)
1530 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1531 TYPE(greens_fn_type),
INTENT(IN) :: green
1532 TYPE(pw_r3d_rs_type),
INTENT(IN) :: pw_in
1533 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: pw_out
1535 CHARACTER(LEN=*),
PARAMETER :: routinen =
'apply_inv_laplace_operator_dct'
1537 INTEGER :: handle, ig, ng
1538 REAL(dp) :: prefactor
1539 TYPE(pw_c1d_gs_type) :: pw_in_gs
1540 TYPE(pw_grid_type),
POINTER :: pw_grid
1542 CALL timeset(routinen, handle)
1545 prefactor = 1.0_dp/fourpi
1547 pw_grid => pw_pool%pw_grid
1548 ng =
SIZE(pw_grid%gsq)
1550 CALL pw_pool%create_pw(pw_in_gs)
1552 CALL pw_transfer(pw_in, pw_in_gs)
1554 pw_in_gs%array(ig) = prefactor*pw_in_gs%array(ig)*green%dct_influence_fn%array(ig)
1556 CALL pw_transfer(pw_in_gs, pw_out)
1558 CALL pw_pool%give_back_pw(pw_in_gs)
1560 CALL timestop(handle)
1562 END SUBROUTINE apply_inv_laplace_operator_dct
1574 SUBROUTINE apply_laplace_operator_fft(pw_pool, green, pw_in, pw_out)
1576 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1577 TYPE(greens_fn_type),
INTENT(IN) :: green
1578 TYPE(pw_r3d_rs_type),
INTENT(IN) :: pw_in
1579 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: pw_out
1581 CHARACTER(LEN=*),
PARAMETER :: routinen =
'apply_laplace_operator_fft'
1583 INTEGER :: g0_index, handle, ig, ng
1585 REAL(dp) :: prefactor
1586 TYPE(pw_c1d_gs_type) :: pw_in_gs
1587 TYPE(pw_grid_type),
POINTER :: pw_grid
1589 CALL timeset(routinen, handle)
1594 pw_grid => pw_pool%pw_grid
1595 ng =
SIZE(pw_in%pw_grid%gsq)
1596 have_g0 = green%influence_fn%pw_grid%have_g0
1598 CALL pw_pool%create_pw(pw_in_gs)
1600 CALL pw_transfer(pw_in, pw_in_gs)
1603 g0_index = green%influence_fn%pw_grid%first_gne0 - 1
1604 pw_in_gs%array(g0_index) = 0.0_dp
1606 DO ig = green%influence_fn%pw_grid%first_gne0, ng
1607 pw_in_gs%array(ig) = prefactor*(pw_in_gs%array(ig)/green%influence_fn%array(ig))
1610 CALL pw_transfer(pw_in_gs, pw_out)
1612 CALL pw_pool%give_back_pw(pw_in_gs)
1614 CALL timestop(handle)
1616 END SUBROUTINE apply_laplace_operator_fft
1629 SUBROUTINE apply_laplace_operator_dct(pw_pool, green, pw_in, pw_out)
1631 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1632 TYPE(greens_fn_type),
INTENT(IN) :: green
1633 TYPE(pw_r3d_rs_type),
INTENT(IN) :: pw_in
1634 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: pw_out
1636 CHARACTER(LEN=*),
PARAMETER :: routinen =
'apply_laplace_operator_dct'
1638 INTEGER :: g0_index, handle, ig, ng
1640 REAL(dp) :: prefactor
1641 TYPE(pw_c1d_gs_type) :: pw_in_gs
1642 TYPE(pw_grid_type),
POINTER :: pw_grid
1644 CALL timeset(routinen, handle)
1649 pw_grid => pw_pool%pw_grid
1650 ng =
SIZE(pw_in%pw_grid%gsq)
1651 have_g0 = green%dct_influence_fn%pw_grid%have_g0
1653 CALL pw_pool%create_pw(pw_in_gs)
1655 CALL pw_transfer(pw_in, pw_in_gs)
1658 g0_index = green%dct_influence_fn%pw_grid%first_gne0 - 1
1659 pw_in_gs%array(g0_index) = 0.0_dp
1661 DO ig = green%dct_influence_fn%pw_grid%first_gne0, ng
1662 pw_in_gs%array(ig) = prefactor*(pw_in_gs%array(ig)/green%dct_influence_fn%array(ig))
1665 CALL pw_transfer(pw_in_gs, pw_out)
1667 CALL pw_pool%give_back_pw(pw_in_gs)
1669 CALL timestop(handle)
1671 END SUBROUTINE apply_laplace_operator_dct
1684 SUBROUTINE apply_poisson_operator_fft(pw_pool, green, dielectric, v, density)
1686 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1687 TYPE(greens_fn_type),
INTENT(IN) :: green
1688 TYPE(dielectric_type),
INTENT(IN),
POINTER :: dielectric
1689 TYPE(pw_r3d_rs_type),
INTENT(IN) :: v
1690 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: density
1692 CHARACTER(LEN=*),
PARAMETER :: routinen =
'apply_poisson_operator_fft'
1695 TYPE(pw_r3d_rs_type) :: pxv
1697 CALL timeset(routinen, handle)
1699 CALL pw_pool%create_pw(pxv)
1701 CALL apply_p_operator(pw_pool, dielectric, v, pxv)
1702 CALL apply_laplace_operator_fft(pw_pool, green, v, density)
1703 CALL pw_axpy(pxv, density)
1705 CALL pw_pool%give_back_pw(pxv)
1707 CALL timestop(handle)
1709 END SUBROUTINE apply_poisson_operator_fft
1724 SUBROUTINE apply_poisson_operator_dct(pw_pool, green, dielectric, v, density)
1726 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1727 TYPE(greens_fn_type),
INTENT(IN) :: green
1728 TYPE(dielectric_type),
INTENT(IN),
POINTER :: dielectric
1729 TYPE(pw_r3d_rs_type),
INTENT(IN) :: v
1730 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: density
1732 CHARACTER(LEN=*),
PARAMETER :: routinen =
'apply_poisson_operator_dct'
1735 TYPE(pw_r3d_rs_type) :: pxv
1737 CALL timeset(routinen, handle)
1739 CALL pw_pool%create_pw(pxv)
1741 CALL apply_p_operator(pw_pool, dielectric, v, pxv)
1742 CALL apply_laplace_operator_dct(pw_pool, green, v, density)
1743 CALL pw_axpy(pxv, density)
1745 CALL pw_pool%give_back_pw(pxv)
1747 CALL timestop(handle)
1749 END SUBROUTINE apply_poisson_operator_dct
1765 SUBROUTINE ps_implicit_compute_veps(pw_pool, dielectric, v, v_eps)
1767 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1768 TYPE(dielectric_type),
INTENT(IN),
POINTER :: dielectric
1769 TYPE(pw_r3d_rs_type),
INTENT(IN) :: v
1770 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: v_eps
1772 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_compute_veps'
1774 INTEGER :: handle, i
1776 TYPE(pw_r3d_rs_type) :: dv2
1777 TYPE(pw_r3d_rs_type),
DIMENSION(3) :: dv
1779 CALL timeset(routinen, handle)
1783 CALL pw_pool%create_pw(dv2)
1785 CALL pw_pool%create_pw(dv(i))
1788 CALL derive_fft(v, dv, pw_pool)
1791 dv2%array = dv(1)%array**2 + dv(2)%array**2 + dv(3)%array**2
1793 v_eps%array = -(1.0_dp/eightpi)*(dv2%array*dielectric%deps_drho%array)
1795 CALL pw_pool%give_back_pw(dv2)
1797 CALL pw_pool%give_back_pw(dv(i))
1800 CALL timestop(handle)
1802 END SUBROUTINE ps_implicit_compute_veps
1813 SUBROUTINE compute_ehartree_periodic_bc(density, v, ehartree)
1815 TYPE(pw_r3d_rs_type),
INTENT(IN) :: density, v
1816 REAL(dp),
INTENT(OUT) :: ehartree
1818 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_ehartree_periodic_bc'
1822 CALL timeset(routinen, handle)
1825 ehartree = 0.5_dp*pw_integral_ab(density, v)
1827 CALL timestop(handle)
1829 END SUBROUTINE compute_ehartree_periodic_bc
1844 SUBROUTINE compute_ehartree_mixed_bc(dielectric, density, Btxlambda, v, ehartree, electric_enthalpy)
1846 TYPE(dielectric_type),
INTENT(IN),
POINTER :: dielectric
1847 TYPE(pw_r3d_rs_type),
INTENT(IN) :: density
1848 REAL(dp),
ALLOCATABLE,
DIMENSION(:, :, :), &
1849 INTENT(IN) :: btxlambda
1850 TYPE(pw_r3d_rs_type),
INTENT(IN) :: v
1851 REAL(dp),
INTENT(OUT) :: ehartree, electric_enthalpy
1853 CHARACTER(LEN=*),
PARAMETER :: routinen =
'compute_ehartree_mixed_bc'
1856 REAL(dp) :: dvol, ehartree_rho, ehartree_rho_cstr
1857 TYPE(pw_grid_type),
POINTER :: pw_grid
1859 CALL timeset(routinen, handle)
1861 pw_grid => v%pw_grid
1868 ehartree_rho = accurate_sum(density%array*v%array)
1869 ehartree_rho_cstr = accurate_sum(dielectric%eps%array*btxlambda*v%array/fourpi)
1870 ehartree_rho = 0.5_dp*ehartree_rho*dvol
1871 ehartree_rho_cstr = 0.5_dp*ehartree_rho_cstr*dvol
1872 CALL pw_grid%para%group%sum(ehartree_rho)
1873 CALL pw_grid%para%group%sum(ehartree_rho_cstr)
1874 electric_enthalpy = ehartree_rho + ehartree_rho_cstr
1875 ehartree = ehartree_rho - ehartree_rho_cstr
1877 CALL timestop(handle)
1879 END SUBROUTINE compute_ehartree_mixed_bc
1896 SUBROUTINE ps_implicit_compute_error_fft(pw_pool, green, res_new, v_old, v_new, &
1897 QAinvxres_new, pres_error, nabs_error)
1899 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1900 TYPE(greens_fn_type),
INTENT(IN) :: green
1901 TYPE(pw_r3d_rs_type),
INTENT(IN) :: res_new, v_old, v_new
1902 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: qainvxres_new
1903 REAL(dp),
INTENT(OUT) :: pres_error, nabs_error
1905 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_compute_error_fft'
1910 CALL timeset(routinen, handle)
1912 vol = pw_pool%pw_grid%vol
1915 CALL apply_inv_laplace_operator_fft(pw_pool, green, res_new, qainvxres_new)
1917 pres_error = accurate_sum(qainvxres_new%array(:, :, :)**2)
1918 CALL pw_pool%pw_grid%para%group%sum(pres_error)
1919 pres_error = sqrt(pres_error)/vol
1923 nabs_error = accurate_sum(abs(v_old%array - v_new%array)**2)
1924 CALL pw_pool%pw_grid%para%group%sum(nabs_error)
1925 nabs_error = sqrt(nabs_error)/vol
1927 CALL timestop(handle)
1929 END SUBROUTINE ps_implicit_compute_error_fft
1947 SUBROUTINE ps_implicit_compute_error_dct(pw_pool, green, res_new, v_old, v_new, &
1948 QAinvxres_new, pres_error, nabs_error)
1950 TYPE(pw_pool_type),
INTENT(IN),
POINTER :: pw_pool
1951 TYPE(greens_fn_type),
INTENT(IN) :: green
1952 TYPE(pw_r3d_rs_type),
INTENT(IN) :: res_new, v_old, v_new
1953 TYPE(pw_r3d_rs_type),
INTENT(INOUT) :: qainvxres_new
1954 REAL(dp),
INTENT(OUT) :: pres_error, nabs_error
1956 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_compute_error_dct'
1961 CALL timeset(routinen, handle)
1963 vol = pw_pool%pw_grid%vol
1966 CALL apply_inv_laplace_operator_dct(pw_pool, green, res_new, qainvxres_new)
1968 pres_error = accurate_sum(qainvxres_new%array(:, :, :)**2)
1969 CALL pw_pool%pw_grid%para%group%sum(pres_error)
1970 pres_error = sqrt(pres_error)/vol
1974 nabs_error = accurate_sum(abs(v_old%array - v_new%array)**2)
1975 CALL pw_pool%pw_grid%para%group%sum(nabs_error)
1976 nabs_error = sqrt(nabs_error)/vol
1978 CALL timestop(handle)
1980 END SUBROUTINE ps_implicit_compute_error_dct
1992 SUBROUTINE ps_implicit_output(iter, pres_error, nabs_error, outp_unit)
1994 INTEGER,
INTENT(IN) :: iter
1995 REAL(dp),
INTENT(IN) :: pres_error, nabs_error
1996 INTEGER,
INTENT(OUT) :: outp_unit
1998 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_output'
1999 INTEGER,
PARAMETER :: low_print_level = 1
2002 TYPE(cp_logger_type),
POINTER :: logger
2004 CALL timeset(routinen, handle)
2006 logger => cp_get_default_logger()
2007 IF (logger%para_env%is_source())
THEN
2008 outp_unit = cp_logger_get_default_unit_nr(logger, local=.true.)
2013 IF (logger%iter_info%print_level .GT. low_print_level)
THEN
2014 IF ((outp_unit .GT. 0) .AND. (iter .EQ. 1))
THEN
2015 WRITE (outp_unit,
'(T3,A)') &
2016 "POISSON| iter pres error nabs error E_hartree delta E"
2019 IF (outp_unit .GT. 0)
THEN
2020 WRITE (outp_unit,
'(T3,A,I6,5X,E13.4,3X,E13.4)', advance=
'NO') &
2021 "POISSON| ", iter, pres_error, nabs_error
2025 CALL timestop(handle)
2027 END SUBROUTINE ps_implicit_output
2038 SUBROUTINE ps_implicit_report_ehartree(ps_implicit_env, outp_unit, ehartree)
2040 TYPE(ps_implicit_type) :: ps_implicit_env
2041 INTEGER,
INTENT(IN) :: outp_unit
2042 REAL(dp),
INTENT(IN) :: ehartree
2044 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_report_ehartree'
2045 INTEGER,
PARAMETER :: low_print_level = 1
2048 TYPE(cp_logger_type),
POINTER :: logger
2050 logger => cp_get_default_logger()
2051 CALL timeset(routinen, handle)
2052 IF (logger%iter_info%print_level .GT. low_print_level)
THEN
2053 IF (outp_unit .GT. 0)
WRITE (outp_unit,
'(F19.10,E10.2)') &
2054 ehartree, ehartree - ps_implicit_env%ehartree
2056 CALL timestop(handle)
2058 END SUBROUTINE ps_implicit_report_ehartree
2069 SUBROUTINE ps_implicit_print_convergence_msg(iter, max_iter, outp_unit)
2071 INTEGER,
INTENT(IN) :: iter, max_iter, outp_unit
2073 CHARACTER(LEN=*),
PARAMETER :: routinen =
'ps_implicit_print_convergence_msg'
2075 CHARACTER(LEN=12) :: msg
2076 INTEGER :: handle, last_iter
2078 CALL timeset(routinen, handle)
2080 last_iter = iter - 1
2082 IF (outp_unit .GT. 0)
THEN
2083 IF (last_iter .EQ. max_iter)
THEN
2084 WRITE (outp_unit,
'(T3,A)') &
2085 "POISSON| No convergence achieved within the maximum number of iterations."
2087 IF (last_iter .LT. max_iter)
THEN
2088 IF (last_iter .EQ. 1)
THEN
2091 msg =
" iterations."
2093 WRITE (outp_unit,
'(T3,A,I0,A)') &
2094 "POISSON| Poisson solver converged in ", last_iter, msg
2097 CALL timestop(handle)
2099 END SUBROUTINE ps_implicit_print_convergence_msg
2107 SUBROUTINE convert_1dto3d(idx_1dto3d, arr1d, arr3d)
2109 INTEGER,
DIMENSION(:),
INTENT(IN) :: idx_1dto3d
2110 REAL(dp),
DIMENSION(:),
INTENT(IN) :: arr1d
2111 REAL(dp),
ALLOCATABLE,
DIMENSION(:, :, :), &
2112 INTENT(INOUT) :: arr3d
2114 CHARACTER(LEN=*),
PARAMETER :: routinen =
'convert_1dto3d'
2116 INTEGER :: handle, i, j, k, l, lb1, lb2, lb3, &
2117 npts1, npts2, npts3, ub1, ub2, ub3
2119 CALL timeset(routinen, handle)
2121 lb1 = lbound(arr3d, 1); ub1 = ubound(arr3d, 1)
2122 lb2 = lbound(arr3d, 2); ub2 = ubound(arr3d, 2)
2123 lb3 = lbound(arr3d, 3); ub3 = ubound(arr3d, 3)
2125 npts1 = ub1 - lb1 + 1
2126 npts2 = ub2 - lb2 + 1
2127 npts3 = ub3 - lb3 + 1
2129 DO l = 1,
SIZE(idx_1dto3d)
2130 k = ((idx_1dto3d(l) - 1)/(npts1*npts2)) + lb3
2131 j = ((idx_1dto3d(l) - 1) - (k - lb3)*npts1*npts2)/npts1 + lb2
2132 i = idx_1dto3d(l) - ((j - lb2)*npts1 + (k - lb3)*npts1*npts2) + lb1 - 1
2133 arr3d(i, j, k) = arr1d(l)
2136 CALL timestop(handle)
2138 END SUBROUTINE convert_1dto3d
2149 SUBROUTINE get_voltage(time, v_D, osc_frac, frequency, phase, v_D_new)
2151 REAL(dp),
INTENT(IN) :: time
2152 REAL(dp),
DIMENSION(:),
INTENT(IN) :: v_d, osc_frac, frequency, phase
2153 REAL(dp),
ALLOCATABLE,
DIMENSION(:),
INTENT(OUT) :: v_d_new
2155 CHARACTER(LEN=*),
PARAMETER :: routinen =
'get_voltage'
2157 INTEGER :: handle, i
2159 CALL timeset(routinen, handle)
2161 ALLOCATE (v_d_new(
SIZE(v_d)))
2164 v_d_new(i) = v_d(i)*(1 - osc_frac(i)) + &
2165 v_d(i)*osc_frac(i)*cos(2*pi*time*frequency(i) + phase(i))
2168 CALL timestop(handle)
2170 END SUBROUTINE get_voltage
static void dgemm(const char transa, const char transb, const int m, const int n, const int k, const double alpha, const double *a, const int lda, const double *b, const int ldb, const double beta, double *c, const int ldc)
Convenient wrapper to hide Fortran nature of dgemm_, swapping a and b.
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public banihashemian2016
various routines to log and control the output. The idea is that decisions about where to log should ...
recursive integer function, public cp_logger_get_default_unit_nr(logger, local, skip_not_ionode)
asks the default unit number of the given logger. try to use cp_logger_get_unit_nr
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
the type I Discrete Cosine Transform (DCT-I)
integer, parameter, public neumannx
integer, parameter, public neumannxy
subroutine, public dct_type_init(pw_grid, neumann_directions, dct_env)
Initializes a dct_type.
integer, parameter, public neumannxz
subroutine, public pw_shrink(neumann_directions, dests_shrink, srcs_shrink, bounds_local_shftd, pw_in, pw_shrinked)
shrinks an evenly symmetric pw_r3d_rs_type data to a pw_r3d_rs_type data that is 8 times smaller (the...
integer, parameter, public neumannxyz
integer, parameter, public neumannz
integer, parameter, public neumannyz
integer, parameter, public neumanny
subroutine, public pw_expand(neumann_directions, recv_msgs_bnds, dests_expand, srcs_expand, flipg_stat, bounds_shftd, pw_in, pw_expanded)
expands a pw_r3d_rs_type data to an evenly symmetric pw_r3d_rs_type data that is 8 times larger than ...
methods for evaluating the dielectric constant
subroutine, public dielectric_create(dielectric, pw_pool, dielectric_params)
allocates memory for a dielectric data type
subroutine, public derive_fft(f, df, pw_pool)
computes the derivative of a function using FFT
dielectric constant data type
subroutines for defining and creating Dirichlet type subdomains
subroutine, public dirichlet_boundary_region_setup(pw_pool, poisson_params, dbcs)
Sets up the Dirichlet boundary condition.
Dirichlet boundary condition data types.
subroutine, public dbc_tile_release(dbc, pw_pool)
releases tiles
sums arrays of real/complex numbers with much reduced round-off as compared to a naive implementation...
Defines the basic variable types.
integer, parameter, public int_8
integer, parameter, public dp
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
real(kind=dp), parameter, public fourpi
The implicit (generalized) Poisson solver.
subroutine, public implicit_poisson_solver_periodic(poisson_env, density, v_new, ehartree)
implicit Poisson solver for periodic boundary conditions
subroutine, public implicit_poisson_solver_mixed(poisson_env, density, v_new, electric_enthalpy)
implicit Poisson solver for mixed boundary conditions (Neumann + Dirichlet)
subroutine, public implicit_poisson_solver_neumann(poisson_env, density, v_new, ehartree)
implicit Poisson solver: zero-average solution of the Poisson equation subject to homogeneous Neumann...
subroutine, public implicit_poisson_solver_mixed_periodic(poisson_env, density, v_new, electric_enthalpy)
implicit Poisson solver for mixed-periodic boundary conditions (periodic + Dirichlet)
subroutine, public ps_implicit_create(pw_pool, poisson_params, dct_pw_grid, green, ps_implicit_env)
Creates implicit Poisson solver environment.
Types containing essential information for running implicit (iterative) Poisson solver.
integer, parameter, public neumann_bc
integer, parameter, public mixed_bc
integer, parameter, public mixed_periodic_bc
integer, parameter, public periodic_bc
functions related to the poisson solver on regular grids
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
subroutine, public pw_pool_release(pool)
releases the given pool (see cp2k/doc/ReferenceCounting.html)
subroutine, public pw_pool_create(pool, pw_grid, max_cache)
creates a pool for pw
type of a logger, at the moment it contains just a print level starting at which level it should be l...
contains all the informations needed by the fft based poisson solvers
parameters for the poisson solver independet of input_section
environment for the poisson solver
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...