(git:34ef472)
ewald_spline_util.F
Go to the documentation of this file.
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 ! **************************************************************************************************
16  USE cell_methods, ONLY: cell_create
17  USE cell_types, ONLY: cell_release,&
18  cell_type
20  cp_logger_type
24  section_vals_type,&
26  USE kinds, ONLY: dp
27  USE message_passing, ONLY: mp_para_env_type
28  USE pw_grid_types, ONLY: halfspace,&
29  pw_grid_type
30  USE pw_grids, ONLY: pw_grid_create,&
32  USE pw_methods, ONLY: pw_zero
33  USE pw_pool_types, ONLY: pw_pool_create,&
34  pw_pool_type
35  USE pw_spline_utils, ONLY: &
38  pw_spline_precond_type, spl3_pbc
39  USE pw_types, ONLY: pw_r3d_rs_type
40 !NB parallelization
41 #include "./base/base_uses.f90"
42 
43  IMPLICIT NONE
44  PRIVATE
45 
46  LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .true.
47  CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'ewald_spline_util'
48  PUBLIC :: setup_ewald_spline
49 
50 CONTAINS
51 
52 ! **************************************************************************************************
53 !> \brief Setup of the G-space Ewald Term Spline Coefficients
54 !> \param pw_grid ...
55 !> \param pw_pool ...
56 !> \param coeff ...
57 !> \param LG ...
58 !> \param gx ...
59 !> \param gy ...
60 !> \param gz ...
61 !> \param hmat ...
62 !> \param npts ...
63 !> \param param_section ...
64 !> \param tag ...
65 !> \param print_section ...
66 !> \param para_env ...
67 !> \par History
68 !> 12.2005 created [tlaino]
69 !> \author Teodoro Laino
70 ! **************************************************************************************************
71  SUBROUTINE setup_ewald_spline(pw_grid, pw_pool, coeff, LG, gx, gy, gz, hmat, npts, &
72  param_section, tag, print_section, para_env)
73  TYPE(pw_grid_type), POINTER :: pw_grid
74  TYPE(pw_pool_type), POINTER :: pw_pool
75  TYPE(pw_r3d_rs_type), POINTER :: coeff
76  REAL(kind=dp), DIMENSION(:), POINTER :: lg, gx, gy, gz
77  REAL(kind=dp), INTENT(IN) :: hmat(3, 3)
78  INTEGER, INTENT(IN) :: npts(3)
79  TYPE(section_vals_type), POINTER :: param_section
80  CHARACTER(LEN=*), INTENT(IN) :: tag
81  TYPE(section_vals_type), POINTER :: print_section
82  TYPE(mp_para_env_type), POINTER :: para_env
83 
84  INTEGER :: bo(2, 3), iounit
85  TYPE(cell_type), POINTER :: cell
86  TYPE(cp_logger_type), POINTER :: logger
87  TYPE(pw_r3d_rs_type) :: pw
88 
89 !
90 ! Setting Up Fit Procedure
91 !
92 
93  cpassert(.NOT. ASSOCIATED(pw_grid))
94  cpassert(.NOT. ASSOCIATED(pw_pool))
95  cpassert(.NOT. ASSOCIATED(coeff))
96  NULLIFY (cell)
97 
98  CALL cell_create(cell, hmat=hmat, periodic=(/1, 1, 1/))
99  CALL pw_grid_create(pw_grid, para_env, local=.true.)
100  logger => cp_get_default_logger()
101  iounit = cp_print_key_unit_nr(logger, print_section, "", &
102  extension=".Log")
103  bo(1, 1:3) = 0
104  bo(2, 1:3) = npts(1:3) - 1
105  CALL pw_grid_setup(cell%hmat, pw_grid, grid_span=halfspace, bounds=bo, iounit=iounit)
106 
107  CALL cp_print_key_finished_output(iounit, logger, print_section, &
108  "")
109  ! pw_pool initialized
110  CALL pw_pool_create(pw_pool, pw_grid=pw_grid)
111  ALLOCATE (coeff)
112  CALL pw_pool%create_pw(pw)
113  CALL pw_pool%create_pw(coeff)
114  ! Evaluate function on grid
115  CALL eval_pw_tablr(pw, pw_pool, coeff, lg, gx, gy, gz, hmat_mm=hmat, &
116  param_section=param_section, tag=tag)
117  CALL pw_pool%give_back_pw(pw)
118  CALL cell_release(cell)
119 
120  END SUBROUTINE setup_ewald_spline
121 
122 ! **************************************************************************************************
123 !> \brief Evaluates the function G-Term in reciprocal space on the grid
124 !> and find the coefficients of the Splines
125 !> \param grid ...
126 !> \param pw_pool ...
127 !> \param TabLR ...
128 !> \param Lg ...
129 !> \param gx ...
130 !> \param gy ...
131 !> \param gz ...
132 !> \param hmat_mm ...
133 !> \param param_section ...
134 !> \param tag ...
135 !> \par History
136 !> 12.2005 created [tlaino]
137 !> \author Teodoro Laino
138 ! **************************************************************************************************
139  SUBROUTINE eval_pw_tablr(grid, pw_pool, TabLR, Lg, gx, gy, gz, hmat_mm, &
140  param_section, tag)
141  TYPE(pw_r3d_rs_type), INTENT(INOUT) :: grid
142  TYPE(pw_pool_type), POINTER :: pw_pool
143  TYPE(pw_r3d_rs_type), POINTER :: tablr
144  REAL(kind=dp), DIMENSION(:), POINTER :: lg, gx, gy, gz
145  REAL(kind=dp), DIMENSION(3, 3) :: hmat_mm
146  TYPE(section_vals_type), POINTER :: param_section
147  CHARACTER(LEN=*), INTENT(IN) :: tag
148 
149  CHARACTER(len=*), PARAMETER :: routinen = 'eval_pw_TabLR'
150 
151  INTEGER :: act_nx, act_ny, aint_precond, handle, i, iii, is, j, js, k, ks, lg_loc_max, &
152  lg_loc_min, max_iter, my_i, my_j, my_k, n1, n2, n3, n_extra, nlg_loc, nxlim, nylim, &
153  nzlim, precond_kind
154  INTEGER, DIMENSION(2, 3) :: gbo
155  LOGICAL :: success
156  REAL(kind=dp) :: dr1, dr2, dr3, eps_r, eps_x, xs1, xs2, &
157  xs3
158  REAL(kind=dp), ALLOCATABLE :: cos_gx(:, :), cos_gy(:, :), &
159  cos_gz(:, :), lhs(:, :), rhs(:, :), &
160  sin_gx(:, :), sin_gy(:, :), &
161  sin_gz(:, :)
162  TYPE(pw_spline_precond_type) :: precond
163  TYPE(section_vals_type), POINTER :: interp_section
164 
165 !NB pull expensive Cos() out of inner looop
166 !NB temporaries for holding stuff so that dgemm can be used
167 
168  EXTERNAL :: dgemm
169 
170  CALL timeset(routinen, handle)
171  n1 = grid%pw_grid%npts(1)
172  n2 = grid%pw_grid%npts(2)
173  n3 = grid%pw_grid%npts(3)
174  dr1 = grid%pw_grid%dr(1)
175  dr2 = grid%pw_grid%dr(2)
176  dr3 = grid%pw_grid%dr(3)
177  gbo = grid%pw_grid%bounds
178  nxlim = floor(real(n1, kind=dp)/2.0_dp)
179  nylim = floor(real(n2, kind=dp)/2.0_dp)
180  nzlim = floor(real(n3, kind=dp)/2.0_dp)
181  is = 0
182  js = 0
183  ks = 0
184  IF (2*nxlim /= n1) is = 1
185  IF (2*nylim /= n2) js = 1
186  IF (2*nzlim /= n3) ks = 1
187  CALL pw_zero(grid)
188 
189  ! Used the full symmetry to reduce the evaluation to 1/64th
190  !NB parallelization
191  iii = 0
192  !NB allocate temporaries for Cos refactoring
193  ALLOCATE (cos_gx(SIZE(lg), gbo(1, 1):gbo(2, 1)))
194  ALLOCATE (sin_gx(SIZE(lg), gbo(1, 1):gbo(2, 1)))
195  ALLOCATE (cos_gy(SIZE(lg), gbo(1, 2):gbo(2, 2)))
196  ALLOCATE (sin_gy(SIZE(lg), gbo(1, 2):gbo(2, 2)))
197  ALLOCATE (cos_gz(SIZE(lg), gbo(1, 3):gbo(2, 3)))
198  ALLOCATE (sin_gz(SIZE(lg), gbo(1, 3):gbo(2, 3)))
199  !NB precalculate Cos(gx*xs1) etc for Cos refactoring
200  DO k = gbo(1, 3), gbo(2, 3)
201  my_k = k - gbo(1, 3)
202  xs3 = real(my_k, dp)*dr3
203  IF (k > nzlim) cycle
204  cos_gz(1:SIZE(lg), k) = cos(gz(1:SIZE(lg))*xs3)
205  sin_gz(1:SIZE(lg), k) = sin(gz(1:SIZE(lg))*xs3)
206  END DO ! k
207  xs2 = 0.0_dp
208  DO j = gbo(1, 2), gbo(2, 2)
209  IF (j > nylim) cycle
210  cos_gy(1:SIZE(lg), j) = cos(gy(1:SIZE(lg))*xs2)
211  sin_gy(1:SIZE(lg), j) = sin(gy(1:SIZE(lg))*xs2)
212  xs2 = xs2 + dr2
213  END DO ! j
214  xs1 = 0.0_dp
215  DO i = gbo(1, 1), gbo(2, 1)
216  IF (i > nxlim) cycle
217  cos_gx(1:SIZE(lg), i) = cos(gx(1:SIZE(lg))*xs1)
218  sin_gx(1:SIZE(lg), i) = sin(gx(1:SIZE(lg))*xs1)
219  xs1 = xs1 + dr1
220  END DO ! i
221 
222  !NB use DGEMM to compute sum over kg for each i, j, k
223  ! number of elements per node, round down
224  nlg_loc = SIZE(lg)/grid%pw_grid%para%group_size
225  ! number of extra elements not yet accounted for
226  n_extra = mod(SIZE(lg), grid%pw_grid%para%group_size)
227  ! first n_extra nodes get NLg_loc+1, remaining get NLg_loc
228  IF (grid%pw_grid%para%my_pos < n_extra) THEN
229  lg_loc_min = (nlg_loc + 1)*grid%pw_grid%para%my_pos + 1
230  lg_loc_max = lg_loc_min + (nlg_loc + 1) - 1
231  ELSE
232  lg_loc_min = (nlg_loc + 1)*n_extra + nlg_loc*(grid%pw_grid%para%my_pos - n_extra) + 1
233  lg_loc_max = lg_loc_min + nlg_loc - 1
234  END IF
235  ! shouldn't be necessary
236  lg_loc_max = min(SIZE(lg), lg_loc_max)
237  nlg_loc = lg_loc_max - lg_loc_min + 1
238 
239  IF (nlg_loc > 0) THEN ! some work for this node
240 
241  act_nx = min(gbo(2, 1), nxlim) - gbo(1, 1) + 1
242  act_ny = min(gbo(2, 2), nylim) - gbo(1, 2) + 1
243  !NB temporaries for DGEMM use
244  ALLOCATE (lhs(act_nx, nlg_loc))
245  ALLOCATE (rhs(act_ny, nlg_loc))
246 
247  ! do cos(gx) cos(gy+gz) term
248  DO i = gbo(1, 1), gbo(2, 1)
249  IF (i > nxlim) cycle
250  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)
251  END DO
252  DO k = gbo(1, 3), gbo(2, 3)
253  IF (k > nzlim) cycle
254  DO j = gbo(1, 2), gbo(2, 2)
255  IF (j > nylim) cycle
256  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) - &
257  sin_gy(lg_loc_min:lg_loc_max, j)*sin_gz(lg_loc_min:lg_loc_max, k)
258  END DO
259  CALL dgemm('N', 'T', act_nx, act_ny, nlg_loc, 1.0d0, lhs(1, 1), act_nx, rhs(1, 1), act_ny, 0.0d0, &
260  grid%array(gbo(1, 1), gbo(1, 2), k), SIZE(grid%array, 1))
261  END DO
262 
263  ! do sin(gx) sin(gy+gz) term
264  DO i = gbo(1, 1), gbo(2, 1)
265  IF (i > nxlim) cycle
266  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)
267  END DO
268  DO k = gbo(1, 3), gbo(2, 3)
269  IF (k > nzlim) cycle
270  DO j = gbo(1, 2), gbo(2, 2)
271  IF (j > nylim) cycle
272  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) + &
273  sin_gy(lg_loc_min:lg_loc_max, j)*cos_gz(lg_loc_min:lg_loc_max, k)
274  END DO
275  CALL dgemm('N', 'T', act_nx, act_ny, nlg_loc, 1.0d0, lhs(1, 1), act_nx, rhs(1, 1), act_ny, 1.0d0, &
276  grid%array(gbo(1, 1), gbo(1, 2), k), SIZE(grid%array, 1))
277  END DO
278 
279  !NB deallocate temporaries for DGEMM use
280  DEALLOCATE (lhs)
281  DEALLOCATE (rhs)
282  !NB deallocate temporaries for Cos refactoring
283  DEALLOCATE (cos_gx)
284  DEALLOCATE (sin_gx)
285  DEALLOCATE (cos_gy)
286  DEALLOCATE (sin_gy)
287  DEALLOCATE (cos_gz)
288  DEALLOCATE (sin_gz)
289  !NB parallelization
290  ELSE ! no work for this node, just zero contribution
291  grid%array(gbo(1, 1):nxlim, gbo(1, 2):nylim, gbo(1, 3):nzlim) = 0.0_dp
292  END IF ! NLg_loc > 0
293 
294  CALL grid%pw_grid%para%group%sum(grid%array(gbo(1, 1):nxlim, gbo(1, 2):nylim, gbo(1, 3):nzlim))
295 
296  fake_loopongrid: DO k = gbo(1, 3), gbo(2, 3)
297  my_k = k
298  IF (k > nzlim) my_k = nzlim - abs(nzlim - k) + ks
299  DO j = gbo(1, 2), gbo(2, 2)
300  my_j = j
301  IF (j > nylim) my_j = nylim - abs(nylim - j) + js
302  DO i = gbo(1, 1), gbo(2, 1)
303  my_i = i
304  IF (i > nxlim) my_i = nxlim - abs(nxlim - i) + is
305  grid%array(i, j, k) = grid%array(my_i, my_j, my_k)
306  END DO
307  END DO
308  END DO fake_loopongrid
309  !
310  ! Solve for spline coefficients
311  !
312  interp_section => section_vals_get_subs_vals(param_section, "INTERPOLATOR")
313  CALL section_vals_val_get(interp_section, "aint_precond", i_val=aint_precond)
314  CALL section_vals_val_get(interp_section, "precond", i_val=precond_kind)
315  CALL section_vals_val_get(interp_section, "max_iter", i_val=max_iter)
316  CALL section_vals_val_get(interp_section, "eps_r", r_val=eps_r)
317  CALL section_vals_val_get(interp_section, "eps_x", r_val=eps_x)
318  !
319  ! Solve for spline coefficients
320  !
321  CALL pw_spline_precond_create(precond, precond_kind=aint_precond, &
322  pool=pw_pool, pbc=.true., transpose=.false.)
323  CALL pw_spline_do_precond(precond, grid, tablr)
324  CALL pw_spline_precond_set_kind(precond, precond_kind)
325  success = find_coeffs(values=grid, coeffs=tablr, &
326  linop=spl3_pbc, preconditioner=precond, pool=pw_pool, &
327  eps_r=eps_r, eps_x=eps_x, &
328  max_iter=max_iter)
329  cpassert(success)
330  CALL pw_spline_precond_release(precond)
331  !
332  ! Check for the interpolation Spline
333  !
334  CALL check_spline_interp_tablr(hmat_mm, lg, gx, gy, gz, tablr, param_section, &
335  tag)
336  CALL timestop(handle)
337  END SUBROUTINE eval_pw_tablr
338 
339 ! **************************************************************************************************
340 !> \brief Routine to check the accuracy for the Spline Interpolation
341 !> \param hmat_mm ...
342 !> \param Lg ...
343 !> \param gx ...
344 !> \param gy ...
345 !> \param gz ...
346 !> \param TabLR ...
347 !> \param param_section ...
348 !> \param tag ...
349 !> \par History
350 !> 12.2005 created [tlaino]
351 !> \author Teodoro Laino
352 ! **************************************************************************************************
353  SUBROUTINE check_spline_interp_tablr(hmat_mm, Lg, gx, gy, gz, TabLR, &
354  param_section, tag)
355  REAL(kind=dp), DIMENSION(3, 3) :: hmat_mm
356  REAL(kind=dp), DIMENSION(:), POINTER :: lg, gx, gy, gz
357  TYPE(pw_r3d_rs_type), POINTER :: tablr
358  TYPE(section_vals_type), POINTER :: param_section
359  CHARACTER(LEN=*), INTENT(IN) :: tag
360 
361  CHARACTER(len=*), PARAMETER :: routinen = 'check_spline_interp_TabLR'
362 
363  INTEGER :: handle, i, iw, kg, npoints
364  REAL(kind=dp) :: dn(3), dr1, dr2, dr3, dxterm, dyterm, &
365  dzterm, errd, errf, fterm, maxerrord, &
366  maxerrorf, na, nn, term, tmp1, tmp2, &
367  vec(3), xs1, xs2, xs3
368  TYPE(cp_logger_type), POINTER :: logger
369 
370  NULLIFY (logger)
371  logger => cp_get_default_logger()
372  iw = cp_print_key_unit_nr(logger, param_section, "check_spline", &
373  extension="."//trim(tag)//"Log")
374  CALL timeset(routinen, handle)
375  IF (iw > 0) THEN
376  npoints = 100
377  errf = 0.0_dp
378  maxerrorf = 0.0_dp
379  errd = 0.0_dp
380  maxerrord = 0.0_dp
381  dr1 = hmat_mm(1, 1)/real(npoints, kind=dp)
382  dr2 = hmat_mm(2, 2)/real(npoints, kind=dp)
383  dr3 = hmat_mm(3, 3)/real(npoints, kind=dp)
384  xs1 = 0.0_dp
385  xs2 = 0.0_dp
386  xs3 = 0.0_dp
387  WRITE (iw, '(A,T5,A15,4X,A17,T50,4X,A,5X,A,T80,A,T85,A15,4X,A17,T130,4X,A,5X,A)') &
388  "#", "Analytical Term", "Interpolated Term", "Error", "MaxError", &
389  "*", " Analyt Deriv ", "Interp Deriv Mod ", "Error", "MaxError"
390  DO i = 1, npoints + 1
391  term = 0.0_dp
392  dxterm = 0.0_dp
393  dyterm = 0.0_dp
394  dzterm = 0.0_dp
395  ! Sum over k vectors
396  DO kg = 1, SIZE(lg)
397  vec = (/real(gx(kg), kind=dp), real(gy(kg), kind=dp), real(gz(kg), kind=dp)/)
398  term = term + lg(kg)*cos(vec(1)*xs1 + vec(2)*xs2 + vec(3)*xs3)
399  dxterm = dxterm - lg(kg)*sin(vec(1)*xs1 + vec(2)*xs2 + vec(3)*xs3)*vec(1)
400  dyterm = dyterm - lg(kg)*sin(vec(1)*xs1 + vec(2)*xs2 + vec(3)*xs3)*vec(2)
401  dzterm = dzterm - lg(kg)*sin(vec(1)*xs1 + vec(2)*xs2 + vec(3)*xs3)*vec(3)
402  END DO
403  na = sqrt(dxterm*dxterm + dyterm*dyterm + dzterm*dzterm)
404  dn = eval_d_interp_spl3_pbc((/xs1, xs2, xs3/), tablr)
405  nn = sqrt(dot_product(dn, dn))
406  fterm = eval_interp_spl3_pbc((/xs1, xs2, xs3/), tablr)
407  tmp1 = abs(term - fterm)
408  tmp2 = sqrt(dot_product(dn - (/dxterm, dyterm, dzterm/), dn - (/dxterm, dyterm, dzterm/)))
409  errf = errf + tmp1
410  maxerrorf = max(maxerrorf, tmp1)
411  errd = errd + tmp2
412  maxerrord = max(maxerrord, tmp2)
413  WRITE (iw, '(T5,F15.10,5X,F15.10,T50,2F12.9,T80,A,T85,F15.10,5X,F15.10,T130,2F12.9)') &
414  term, fterm, tmp1, maxerrorf, "*", na, nn, tmp2, maxerrord
415  xs1 = xs1 + dr1
416  xs2 = xs2 + dr2
417  xs3 = xs3 + dr3
418  END DO
419  WRITE (iw, '(A,T5,A,T50,F12.9,T130,F12.9)') "#", "Averages", errf/real(npoints, kind=dp), &
420  errd/real(npoints, kind=dp)
421  END IF
422  CALL timestop(handle)
423  CALL cp_print_key_finished_output(iw, logger, param_section, "check_spline")
424 
425  END SUBROUTINE check_spline_interp_tablr
426 
427 END MODULE ewald_spline_util
428 
subroutine pbc(r, r_pbc, s, s_pbc, a, b, c, alpha, beta, gamma, debug, info, pbc0, h, hinv)
...
Definition: dumpdcd.F:1203
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.
Handles all functions related to the CELL.
Definition: cell_methods.F:15
subroutine, public cell_create(cell, hmat, periodic, tag)
allocates and initializes a cell
Definition: cell_methods.F:85
Handles all functions related to the CELL.
Definition: cell_types.F:15
subroutine, public cell_release(cell)
releases the given cell (see doc/ReferenceCounting.html)
Definition: cell_types.F:559
various routines to log and control the output. The idea is that decisions about where to log should ...
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
Setting up the Spline coefficients used to Interpolate the G-Term in Ewald sums.
subroutine, public setup_ewald_spline(pw_grid, pw_pool, coeff, LG, gx, gy, gz, hmat, npts, param_section, tag, print_section, para_env)
Setup of the G-space Ewald Term Spline Coefficients.
objects that represent the structure of input sections and the data contained in an input section
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
Defines the basic variable types.
Definition: kinds.F:23
integer, parameter, public dp
Definition: kinds.F:34
Interface to the message passing library MPI.
computes preconditioners, and implements methods to apply them currently used in qs_ot
integer, parameter, public halfspace
Definition: pw_grid_types.F:28
This module defines the grid data type and some basic operations on it.
Definition: pw_grids.F:36
subroutine, public pw_grid_setup(cell_hmat, pw_grid, grid_span, cutoff, bounds, bounds_local, npts, spherical, odd, fft_usage, ncommensurate, icommensurate, blocked, ref_grid, rs_dims, iounit)
sets up a pw_grid
Definition: pw_grids.F:286
subroutine, public pw_grid_create(pw_grid, pe_group, local)
Initialize a PW grid with all defaults.
Definition: pw_grids.F:93
Manages a pool of grids (to be used for example as tmp objects), but can also be used to instantiate ...
Definition: pw_pool_types.F:24
subroutine, public pw_pool_create(pool, pw_grid, max_cache)
creates a pool for pw
different utils that are useful to manipulate splines on the regular grid of a pw
logical function, public find_coeffs(values, coeffs, linOp, preconditioner, pool, eps_r, eps_x, max_iter, sumtype)
solves iteratively (CG) a systmes of linear equations linOp(coeffs)=values (for example those needed ...
subroutine, public pw_spline_precond_release(preconditioner)
releases the preconditioner
subroutine, public pw_spline_precond_create(preconditioner, precond_kind, pool, pbc, transpose)
...
subroutine, public pw_spline_do_precond(preconditioner, in_v, out_v)
applies the preconditioner to the system of equations to find the coefficients of the spline
subroutine, public pw_spline_precond_set_kind(preconditioner, precond_kind, pbc, transpose)
switches the types of precoditioner to use
real(kind=dp) function, dimension(3), public eval_d_interp_spl3_pbc(vec, pw)
Evaluates the derivatives of the PBC interpolated Spline (pw) function on the generic input vector (v...
real(kind=dp) function, public eval_interp_spl3_pbc(vec, pw)
Evaluates the PBC interpolated Spline (pw) function on the generic input vector (vec)
subroutine, public spl3_pbc(pw_in, pw_out)
...