LCOV - code coverage report
Current view: top level - src - ewald_spline_util.F (source / functions) Hit Total Coverage
Test: CP2K Regtests (git:4dc10b3) Lines: 178 182 97.8 %
Date: 2024-11-21 06:45:46 Functions: 3 3 100.0 %

          Line data    Source code
       1             : !--------------------------------------------------------------------------------------------------!
       2             : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3             : !   Copyright 2000-2024 CP2K developers group <https://cp2k.org>                                   !
       4             : !                                                                                                  !
       5             : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6             : !--------------------------------------------------------------------------------------------------!
       7             : 
       8             : ! **************************************************************************************************
       9             : !> \brief Setting up the Spline coefficients used to Interpolate the G-Term
      10             : !>      in Ewald sums
      11             : !> \par History
      12             : !>      12.2005 created [tlaino]
      13             : !> \author Teodoro Laino
      14             : ! **************************************************************************************************
      15             : MODULE ewald_spline_util
      16             :    USE cell_methods,                    ONLY: cell_create
      17             :    USE cell_types,                      ONLY: cell_release,&
      18             :                                               cell_type
      19             :    USE cp_log_handling,                 ONLY: cp_get_default_logger,&
      20             :                                               cp_logger_type
      21             :    USE cp_output_handling,              ONLY: cp_print_key_finished_output,&
      22             :                                               cp_print_key_unit_nr
      23             :    USE input_section_types,             ONLY: section_vals_get_subs_vals,&
      24             :                                               section_vals_type,&
      25             :                                               section_vals_val_get
      26             :    USE kinds,                           ONLY: dp
      27             :    USE message_passing,                 ONLY: mp_comm_self
      28             :    USE pw_grid_types,                   ONLY: HALFSPACE,&
      29             :                                               pw_grid_type
      30             :    USE pw_grids,                        ONLY: pw_grid_create
      31             :    USE pw_methods,                      ONLY: pw_zero
      32             :    USE pw_pool_types,                   ONLY: pw_pool_create,&
      33             :                                               pw_pool_type
      34             :    USE pw_spline_utils,                 ONLY: &
      35             :         Eval_Interp_Spl3_pbc, Eval_d_Interp_Spl3_pbc, find_coeffs, pw_spline_do_precond, &
      36             :         pw_spline_precond_create, pw_spline_precond_release, pw_spline_precond_set_kind, &
      37             :         pw_spline_precond_type, spl3_pbc
      38             :    USE pw_types,                        ONLY: pw_r3d_rs_type
      39             : !NB parallelization
      40             : #include "./base/base_uses.f90"
      41             : 
      42             :    IMPLICIT NONE
      43             :    PRIVATE
      44             : 
      45             :    LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .TRUE.
      46             :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'ewald_spline_util'
      47             :    PUBLIC :: Setup_Ewald_Spline
      48             : 
      49             : CONTAINS
      50             : 
      51             : ! **************************************************************************************************
      52             : !> \brief Setup of the G-space Ewald Term Spline Coefficients
      53             : !> \param pw_grid ...
      54             : !> \param pw_pool ...
      55             : !> \param coeff ...
      56             : !> \param LG ...
      57             : !> \param gx ...
      58             : !> \param gy ...
      59             : !> \param gz ...
      60             : !> \param hmat ...
      61             : !> \param npts ...
      62             : !> \param param_section ...
      63             : !> \param tag ...
      64             : !> \param print_section ...
      65             : !> \par History
      66             : !>      12.2005 created [tlaino]
      67             : !> \author Teodoro Laino
      68             : ! **************************************************************************************************
      69         102 :    SUBROUTINE Setup_Ewald_Spline(pw_grid, pw_pool, coeff, LG, gx, gy, gz, hmat, npts, &
      70             :                                  param_section, tag, print_section)
      71             :       TYPE(pw_grid_type), POINTER                        :: pw_grid
      72             :       TYPE(pw_pool_type), POINTER                        :: pw_pool
      73             :       TYPE(pw_r3d_rs_type), POINTER                      :: coeff
      74             :       REAL(KIND=dp), DIMENSION(:), POINTER               :: LG, gx, gy, gz
      75             :       REAL(KIND=dp), INTENT(IN)                          :: hmat(3, 3)
      76             :       INTEGER, INTENT(IN)                                :: npts(3)
      77             :       TYPE(section_vals_type), POINTER                   :: param_section
      78             :       CHARACTER(LEN=*), INTENT(IN)                       :: tag
      79             :       TYPE(section_vals_type), POINTER                   :: print_section
      80             : 
      81             :       INTEGER                                            :: bo(2, 3), iounit
      82             :       TYPE(cell_type), POINTER                           :: cell
      83             :       TYPE(cp_logger_type), POINTER                      :: logger
      84             :       TYPE(pw_r3d_rs_type)                               :: pw
      85             : 
      86             : !
      87             : ! Setting Up Fit Procedure
      88             : !
      89             : 
      90           0 :       CPASSERT(.NOT. ASSOCIATED(pw_grid))
      91         102 :       CPASSERT(.NOT. ASSOCIATED(pw_pool))
      92         102 :       CPASSERT(.NOT. ASSOCIATED(coeff))
      93         102 :       NULLIFY (cell)
      94             : 
      95         102 :       CALL cell_create(cell, hmat=hmat, periodic=(/1, 1, 1/))
      96         102 :       logger => cp_get_default_logger()
      97             :       iounit = cp_print_key_unit_nr(logger, print_section, "", &
      98         102 :                                     extension=".Log")
      99         408 :       bo(1, 1:3) = 0
     100         408 :       bo(2, 1:3) = npts(1:3) - 1
     101         102 :       CALL pw_grid_create(pw_grid, mp_comm_self, cell%hmat, grid_span=HALFSPACE, bounds=bo, iounit=iounit)
     102             : 
     103             :       CALL cp_print_key_finished_output(iounit, logger, print_section, &
     104         102 :                                         "")
     105             :       ! pw_pool initialized
     106         102 :       CALL pw_pool_create(pw_pool, pw_grid=pw_grid)
     107         102 :       ALLOCATE (coeff)
     108         102 :       CALL pw_pool%create_pw(pw)
     109         102 :       CALL pw_pool%create_pw(coeff)
     110             :       ! Evaluate function on grid
     111             :       CALL eval_pw_TabLR(pw, pw_pool, coeff, Lg, gx, gy, gz, hmat_mm=hmat, &
     112         102 :                          param_section=param_section, tag=tag)
     113         102 :       CALL pw_pool%give_back_pw(pw)
     114         102 :       CALL cell_release(cell)
     115             : 
     116         102 :    END SUBROUTINE Setup_Ewald_Spline
     117             : 
     118             : ! **************************************************************************************************
     119             : !> \brief Evaluates the function G-Term in reciprocal space on the grid
     120             : !>      and find the coefficients of the Splines
     121             : !> \param grid ...
     122             : !> \param pw_pool ...
     123             : !> \param TabLR ...
     124             : !> \param Lg ...
     125             : !> \param gx ...
     126             : !> \param gy ...
     127             : !> \param gz ...
     128             : !> \param hmat_mm ...
     129             : !> \param param_section ...
     130             : !> \param tag ...
     131             : !> \par History
     132             : !>      12.2005 created [tlaino]
     133             : !> \author Teodoro Laino
     134             : ! **************************************************************************************************
     135         102 :    SUBROUTINE eval_pw_TabLR(grid, pw_pool, TabLR, Lg, gx, gy, gz, hmat_mm, &
     136             :                             param_section, tag)
     137             :       TYPE(pw_r3d_rs_type), INTENT(INOUT)                :: grid
     138             :       TYPE(pw_pool_type), POINTER                        :: pw_pool
     139             :       TYPE(pw_r3d_rs_type), POINTER                      :: TabLR
     140             :       REAL(KIND=dp), DIMENSION(:), POINTER               :: Lg, gx, gy, gz
     141             :       REAL(KIND=dp), DIMENSION(3, 3)                     :: hmat_mm
     142             :       TYPE(section_vals_type), POINTER                   :: param_section
     143             :       CHARACTER(LEN=*), INTENT(IN)                       :: tag
     144             : 
     145             :       CHARACTER(len=*), PARAMETER                        :: routineN = 'eval_pw_TabLR'
     146             : 
     147             :       INTEGER :: act_nx, act_ny, aint_precond, handle, i, iii, is, j, js, k, ks, Lg_loc_max, &
     148             :          Lg_loc_min, max_iter, my_i, my_j, my_k, n1, n2, n3, n_extra, NLg_loc, nxlim, nylim, &
     149             :          nzlim, precond_kind
     150             :       INTEGER, DIMENSION(2, 3)                           :: gbo
     151             :       LOGICAL                                            :: success
     152             :       REAL(KIND=dp)                                      :: dr1, dr2, dr3, eps_r, eps_x, xs1, xs2, &
     153             :                                                             xs3
     154         102 :       REAL(KIND=dp), ALLOCATABLE                         :: cos_gx(:, :), cos_gy(:, :), &
     155         102 :                                                             cos_gz(:, :), lhs(:, :), rhs(:, :), &
     156         102 :                                                             sin_gx(:, :), sin_gy(:, :), &
     157         102 :                                                             sin_gz(:, :)
     158             :       TYPE(pw_spline_precond_type)                       :: precond
     159             :       TYPE(section_vals_type), POINTER                   :: interp_section
     160             : 
     161             : !NB pull expensive Cos() out of inner looop
     162             : !NB temporaries for holding stuff so that dgemm can be used
     163             : 
     164             :       EXTERNAL :: DGEMM
     165             : 
     166         102 :       CALL timeset(routineN, handle)
     167         102 :       n1 = grid%pw_grid%npts(1)
     168         102 :       n2 = grid%pw_grid%npts(2)
     169         102 :       n3 = grid%pw_grid%npts(3)
     170         102 :       dr1 = grid%pw_grid%dr(1)
     171         102 :       dr2 = grid%pw_grid%dr(2)
     172         102 :       dr3 = grid%pw_grid%dr(3)
     173        1020 :       gbo = grid%pw_grid%bounds
     174         102 :       nxlim = FLOOR(REAL(n1, KIND=dp)/2.0_dp)
     175         102 :       nylim = FLOOR(REAL(n2, KIND=dp)/2.0_dp)
     176         102 :       nzlim = FLOOR(REAL(n3, KIND=dp)/2.0_dp)
     177         102 :       is = 0
     178         102 :       js = 0
     179         102 :       ks = 0
     180         102 :       IF (2*nxlim /= n1) is = 1
     181         102 :       IF (2*nylim /= n2) js = 1
     182         102 :       IF (2*nzlim /= n3) ks = 1
     183         102 :       CALL pw_zero(grid)
     184             : 
     185             :       ! Used the full symmetry to reduce the evaluation to 1/64th
     186             :       !NB parallelization
     187         102 :       iii = 0
     188             :       !NB allocate temporaries for Cos refactoring
     189         408 :       ALLOCATE (cos_gx(SIZE(Lg), gbo(1, 1):gbo(2, 1)))
     190         306 :       ALLOCATE (sin_gx(SIZE(Lg), gbo(1, 1):gbo(2, 1)))
     191         408 :       ALLOCATE (cos_gy(SIZE(Lg), gbo(1, 2):gbo(2, 2)))
     192         306 :       ALLOCATE (sin_gy(SIZE(Lg), gbo(1, 2):gbo(2, 2)))
     193         408 :       ALLOCATE (cos_gz(SIZE(Lg), gbo(1, 3):gbo(2, 3)))
     194         306 :       ALLOCATE (sin_gz(SIZE(Lg), gbo(1, 3):gbo(2, 3)))
     195             :       !NB precalculate Cos(gx*xs1) etc for Cos refactoring
     196        5094 :       DO k = gbo(1, 3), gbo(2, 3)
     197        4992 :          my_k = k - gbo(1, 3)
     198        4992 :          xs3 = REAL(my_k, dp)*dr3
     199        4992 :          IF (k > nzlim) CYCLE
     200     1824444 :          cos_gz(1:SIZE(Lg), k) = COS(gz(1:SIZE(Lg))*xs3)
     201     1826942 :          sin_gz(1:SIZE(Lg), k) = SIN(gz(1:SIZE(Lg))*xs3)
     202             :       END DO ! k
     203             :       xs2 = 0.0_dp
     204        5118 :       DO j = gbo(1, 2), gbo(2, 2)
     205        5016 :          IF (j > nylim) CYCLE
     206     1920156 :          cos_gy(1:SIZE(Lg), j) = COS(gy(1:SIZE(Lg))*xs2)
     207     1920156 :          sin_gy(1:SIZE(Lg), j) = SIN(gy(1:SIZE(Lg))*xs2)
     208        5118 :          xs2 = xs2 + dr2
     209             :       END DO ! j
     210             :       xs1 = 0.0_dp
     211        5070 :       DO i = gbo(1, 1), gbo(2, 1)
     212        4968 :          IF (i > nxlim) CYCLE
     213     1728732 :          cos_gx(1:SIZE(Lg), i) = COS(gx(1:SIZE(Lg))*xs1)
     214     1728732 :          sin_gx(1:SIZE(Lg), i) = SIN(gx(1:SIZE(Lg))*xs1)
     215        5070 :          xs1 = xs1 + dr1
     216             :       END DO ! i
     217             : 
     218             :       !NB use DGEMM to compute sum over kg for each i, j, k
     219             :       ! number of elements per node, round down
     220         102 :       NLg_loc = SIZE(Lg)/grid%pw_grid%para%group%num_pe
     221             :       ! number of extra elements not yet accounted for
     222         102 :       n_extra = MOD(SIZE(Lg), grid%pw_grid%para%group%num_pe)
     223             :       ! first n_extra nodes get NLg_loc+1, remaining get NLg_loc
     224         102 :       IF (grid%pw_grid%para%group%mepos < n_extra) THEN
     225           0 :          Lg_loc_min = (NLg_loc + 1)*grid%pw_grid%para%group%mepos + 1
     226           0 :          Lg_loc_max = Lg_loc_min + (NLg_loc + 1) - 1
     227             :       ELSE
     228         102 :          Lg_loc_min = (NLg_loc + 1)*n_extra + NLg_loc*(grid%pw_grid%para%group%mepos - n_extra) + 1
     229         102 :          Lg_loc_max = Lg_loc_min + NLg_loc - 1
     230             :       END IF
     231             :       ! shouldn't be necessary
     232         102 :       Lg_loc_max = MIN(SIZE(Lg), Lg_loc_max)
     233         102 :       NLg_loc = Lg_loc_max - Lg_loc_min + 1
     234             : 
     235         102 :       IF (NLg_loc > 0) THEN ! some work for this node
     236             : 
     237         102 :          act_nx = MIN(gbo(2, 1), nxlim) - gbo(1, 1) + 1
     238         102 :          act_ny = MIN(gbo(2, 2), nylim) - gbo(1, 2) + 1
     239             :          !NB temporaries for DGEMM use
     240         408 :          ALLOCATE (lhs(act_nx, NLg_loc))
     241         408 :          ALLOCATE (rhs(act_ny, NLg_loc))
     242             : 
     243             :          ! do cos(gx) cos(gy+gz) term
     244        5070 :          DO i = gbo(1, 1), gbo(2, 1)
     245        4968 :             IF (i > nxlim) CYCLE
     246     1728834 :             lhs(i - gbo(1, 1) + 1, 1:NLg_loc) = lg(Lg_loc_min:Lg_loc_max)*cos_gx(Lg_loc_min:Lg_loc_max, i)
     247             :          END DO
     248        5094 :          DO k = gbo(1, 3), gbo(2, 3)
     249        4992 :             IF (k > nzlim) CYCLE
     250      131388 :             DO j = gbo(1, 2), gbo(2, 2)
     251      128792 :                IF (j > nylim) CYCLE
     252             :                rhs(j - gbo(1, 2) + 1, 1:NLg_loc) = cos_gy(Lg_loc_min:Lg_loc_max, j)*cos_gz(Lg_loc_min:Lg_loc_max, k) - &
     253    43226812 :                                                    sin_gy(Lg_loc_min:Lg_loc_max, j)*sin_gz(Lg_loc_min:Lg_loc_max, k)
     254             :             END DO
     255             :             CALL DGEMM('N', 'T', act_nx, act_ny, NLg_loc, 1.0D0, lhs(1, 1), act_nx, rhs(1, 1), act_ny, 0.0D0, &
     256        5094 :                        grid%array(gbo(1, 1), gbo(1, 2), k), SIZE(grid%array, 1))
     257             :          END DO
     258             : 
     259             :          ! do sin(gx) sin(gy+gz) term
     260        5070 :          DO i = gbo(1, 1), gbo(2, 1)
     261        4968 :             IF (i > nxlim) CYCLE
     262     1728834 :             lhs(i - gbo(1, 1) + 1, 1:NLg_loc) = -lg(Lg_loc_min:Lg_loc_max)*sin_gx(Lg_loc_min:Lg_loc_max, i)
     263             :          END DO
     264        5094 :          DO k = gbo(1, 3), gbo(2, 3)
     265        4992 :             IF (k > nzlim) CYCLE
     266      131388 :             DO j = gbo(1, 2), gbo(2, 2)
     267      128792 :                IF (j > nylim) CYCLE
     268             :                rhs(j - gbo(1, 2) + 1, 1:NLg_loc) = cos_gy(Lg_loc_min:Lg_loc_max, j)*sin_gz(Lg_loc_min:Lg_loc_max, k) + &
     269    43226812 :                                                    sin_gy(Lg_loc_min:Lg_loc_max, j)*cos_gz(Lg_loc_min:Lg_loc_max, k)
     270             :             END DO
     271             :             CALL DGEMM('N', 'T', act_nx, act_ny, NLg_loc, 1.0D0, lhs(1, 1), act_nx, rhs(1, 1), act_ny, 1.0D0, &
     272        5094 :                        grid%array(gbo(1, 1), gbo(1, 2), k), SIZE(grid%array, 1))
     273             :          END DO
     274             : 
     275             :          !NB deallocate temporaries for DGEMM use
     276         102 :          DEALLOCATE (lhs)
     277         102 :          DEALLOCATE (rhs)
     278             :          !NB deallocate temporaries for Cos refactoring
     279         102 :          DEALLOCATE (cos_gx)
     280         102 :          DEALLOCATE (sin_gx)
     281         102 :          DEALLOCATE (cos_gy)
     282         102 :          DEALLOCATE (sin_gy)
     283         102 :          DEALLOCATE (cos_gz)
     284         102 :          DEALLOCATE (sin_gz)
     285             :          !NB parallelization
     286             :       ELSE ! no work for this node, just zero contribution
     287           0 :          grid%array(gbo(1, 1):nxlim, gbo(1, 2):nylim, gbo(1, 3):nzlim) = 0.0_dp
     288             :       END IF ! NLg_loc > 0
     289             : 
     290     3597086 :       CALL grid%pw_grid%para%group%sum(grid%array(gbo(1, 1):nxlim, gbo(1, 2):nylim, gbo(1, 3):nzlim))
     291             : 
     292        5094 :       Fake_LoopOnGrid: DO k = gbo(1, 3), gbo(2, 3)
     293        4992 :          my_k = k
     294        4992 :          IF (k > nzlim) my_k = nzlim - ABS(nzlim - k) + ks
     295      252762 :          DO j = gbo(1, 2), gbo(2, 2)
     296      247668 :             my_j = j
     297      247668 :             IF (j > nylim) my_j = nylim - ABS(nylim - j) + js
     298    12548016 :             DO i = gbo(1, 1), gbo(2, 1)
     299    12295356 :                my_i = i
     300    12295356 :                IF (i > nxlim) my_i = nxlim - ABS(nxlim - i) + is
     301    12543024 :                grid%array(i, j, k) = grid%array(my_i, my_j, my_k)
     302             :             END DO
     303             :          END DO
     304             :       END DO Fake_LoopOnGrid
     305             :       !
     306             :       ! Solve for spline coefficients
     307             :       !
     308         102 :       interp_section => section_vals_get_subs_vals(param_section, "INTERPOLATOR")
     309         102 :       CALL section_vals_val_get(interp_section, "aint_precond", i_val=aint_precond)
     310         102 :       CALL section_vals_val_get(interp_section, "precond", i_val=precond_kind)
     311         102 :       CALL section_vals_val_get(interp_section, "max_iter", i_val=max_iter)
     312         102 :       CALL section_vals_val_get(interp_section, "eps_r", r_val=eps_r)
     313         102 :       CALL section_vals_val_get(interp_section, "eps_x", r_val=eps_x)
     314             :       !
     315             :       ! Solve for spline coefficients
     316             :       !
     317             :       CALL pw_spline_precond_create(precond, precond_kind=aint_precond, &
     318         102 :                                     pool=pw_pool, pbc=.TRUE., transpose=.FALSE.)
     319         102 :       CALL pw_spline_do_precond(precond, grid, TabLR)
     320         102 :       CALL pw_spline_precond_set_kind(precond, precond_kind)
     321             :       success = find_coeffs(values=grid, coeffs=TabLR, &
     322             :                             linOp=spl3_pbc, preconditioner=precond, pool=pw_pool, &
     323             :                             eps_r=eps_r, eps_x=eps_x, &
     324         102 :                             max_iter=max_iter)
     325         102 :       CPASSERT(success)
     326         102 :       CALL pw_spline_precond_release(precond)
     327             :       !
     328             :       ! Check for the interpolation Spline
     329             :       !
     330             :       CALL check_spline_interp_TabLR(hmat_mm, Lg, gx, gy, gz, TabLR, param_section, &
     331         102 :                                      tag)
     332         102 :       CALL timestop(handle)
     333        1020 :    END SUBROUTINE eval_pw_TabLR
     334             : 
     335             : ! **************************************************************************************************
     336             : !> \brief Routine to check the accuracy for the Spline Interpolation
     337             : !> \param hmat_mm ...
     338             : !> \param Lg ...
     339             : !> \param gx ...
     340             : !> \param gy ...
     341             : !> \param gz ...
     342             : !> \param TabLR ...
     343             : !> \param param_section ...
     344             : !> \param tag ...
     345             : !> \par History
     346             : !>      12.2005 created [tlaino]
     347             : !> \author Teodoro Laino
     348             : ! **************************************************************************************************
     349         204 :    SUBROUTINE check_spline_interp_TabLR(hmat_mm, Lg, gx, gy, gz, TabLR, &
     350             :                                         param_section, tag)
     351             :       REAL(KIND=dp), DIMENSION(3, 3)                     :: hmat_mm
     352             :       REAL(KIND=dp), DIMENSION(:), POINTER               :: Lg, gx, gy, gz
     353             :       TYPE(pw_r3d_rs_type), POINTER                      :: TabLR
     354             :       TYPE(section_vals_type), POINTER                   :: param_section
     355             :       CHARACTER(LEN=*), INTENT(IN)                       :: tag
     356             : 
     357             :       CHARACTER(len=*), PARAMETER :: routineN = 'check_spline_interp_TabLR'
     358             : 
     359             :       INTEGER                                            :: handle, i, iw, kg, npoints
     360             :       REAL(KIND=dp)                                      :: dn(3), dr1, dr2, dr3, dxTerm, dyTerm, &
     361             :                                                             dzTerm, errd, errf, Fterm, maxerrord, &
     362             :                                                             maxerrorf, Na, Nn, Term, tmp1, tmp2, &
     363             :                                                             vec(3), xs1, xs2, xs3
     364             :       TYPE(cp_logger_type), POINTER                      :: logger
     365             : 
     366         102 :       NULLIFY (logger)
     367         102 :       logger => cp_get_default_logger()
     368             :       iw = cp_print_key_unit_nr(logger, param_section, "check_spline", &
     369         102 :                                 extension="."//TRIM(tag)//"Log")
     370         102 :       CALL timeset(routineN, handle)
     371         102 :       IF (iw > 0) THEN
     372          46 :          npoints = 100
     373          46 :          errf = 0.0_dp
     374          46 :          maxerrorf = 0.0_dp
     375          46 :          errd = 0.0_dp
     376          46 :          maxerrord = 0.0_dp
     377          46 :          dr1 = hmat_mm(1, 1)/REAL(npoints, KIND=dp)
     378          46 :          dr2 = hmat_mm(2, 2)/REAL(npoints, KIND=dp)
     379          46 :          dr3 = hmat_mm(3, 3)/REAL(npoints, KIND=dp)
     380          46 :          xs1 = 0.0_dp
     381          46 :          xs2 = 0.0_dp
     382          46 :          xs3 = 0.0_dp
     383             :          WRITE (iw, '(A,T5,A15,4X,A17,T50,4X,A,5X,A,T80,A,T85,A15,4X,A17,T130,4X,A,5X,A)') &
     384          46 :             "#", "Analytical Term", "Interpolated Term", "Error", "MaxError", &
     385          92 :             "*", " Analyt Deriv  ", "Interp Deriv Mod ", "Error", "MaxError"
     386        4692 :          DO i = 1, npoints + 1
     387        4646 :             Term = 0.0_dp
     388        4646 :             dxTerm = 0.0_dp
     389        4646 :             dyTerm = 0.0_dp
     390        4646 :             dzTerm = 0.0_dp
     391             :             ! Sum over k vectors
     392     4670947 :             DO kg = 1, SIZE(Lg)
     393    18665204 :                vec = (/REAL(gx(kg), KIND=dp), REAL(gy(kg), KIND=dp), REAL(gz(kg), KIND=dp)/)
     394     4666301 :                Term = Term + lg(kg)*COS(vec(1)*xs1 + vec(2)*xs2 + vec(3)*xs3)
     395     4666301 :                dxTerm = dxTerm - lg(kg)*SIN(vec(1)*xs1 + vec(2)*xs2 + vec(3)*xs3)*vec(1)
     396     4666301 :                dyTerm = dyTerm - lg(kg)*SIN(vec(1)*xs1 + vec(2)*xs2 + vec(3)*xs3)*vec(2)
     397     4670947 :                dzTerm = dzTerm - lg(kg)*SIN(vec(1)*xs1 + vec(2)*xs2 + vec(3)*xs3)*vec(3)
     398             :             END DO
     399        4646 :             Na = SQRT(dxTerm*dxTerm + dyTerm*dyTerm + dzTerm*dzTerm)
     400       18584 :             dn = Eval_d_Interp_Spl3_pbc((/xs1, xs2, xs3/), TabLR)
     401       18584 :             Nn = SQRT(DOT_PRODUCT(dn, dn))
     402       18584 :             Fterm = Eval_Interp_Spl3_pbc((/xs1, xs2, xs3/), TabLR)
     403        4646 :             tmp1 = ABS(Term - Fterm)
     404       18584 :             tmp2 = SQRT(DOT_PRODUCT(dn - (/dxTerm, dyTerm, dzTerm/), dn - (/dxTerm, dyTerm, dzTerm/)))
     405        4646 :             errf = errf + tmp1
     406        4646 :             maxerrorf = MAX(maxerrorf, tmp1)
     407        4646 :             errd = errd + tmp2
     408        4646 :             maxerrord = MAX(maxerrord, tmp2)
     409             :             WRITE (iw, '(T5,F15.10,5X,F15.10,T50,2F12.9,T80,A,T85,F15.10,5X,F15.10,T130,2F12.9)') &
     410        4646 :                Term, Fterm, tmp1, maxerrorf, "*", Na, Nn, tmp2, maxerrord
     411        4646 :             xs1 = xs1 + dr1
     412        4646 :             xs2 = xs2 + dr2
     413        4692 :             xs3 = xs3 + dr3
     414             :          END DO
     415          46 :          WRITE (iw, '(A,T5,A,T50,F12.9,T130,F12.9)') "#", "Averages", errf/REAL(npoints, kind=dp), &
     416          92 :             errd/REAL(npoints, kind=dp)
     417             :       END IF
     418         102 :       CALL timestop(handle)
     419         102 :       CALL cp_print_key_finished_output(iw, logger, param_section, "check_spline")
     420             : 
     421         102 :    END SUBROUTINE check_spline_interp_TabLR
     422             : 
     423             : END MODULE ewald_spline_util
     424             : 

Generated by: LCOV version 1.15