(git:ed6f26b)
Loading...
Searching...
No Matches
qmmm_per_elpot.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Setting up the potential for QM/MM periodic boundary conditions calculations
10!> \par History
11!> 07.2005 created [tlaino]
12!> \author Teodoro Laino
13! **************************************************************************************************
15 USE ao_util, ONLY: exp_radius
16 USE cell_types, ONLY: cell_type
37 USE kinds, ONLY: dp
38 USE mathconstants, ONLY: pi
49#include "./base/base_uses.f90"
50
51 IMPLICIT NONE
52 PRIVATE
53
54 LOGICAL, PRIVATE, PARAMETER :: debug_this_module = .true.
55 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'qmmm_per_elpot'
57
58CONTAINS
59
60! **************************************************************************************************
61!> \brief Initialize the QMMM potential stored on vector,
62!> according the qmmm_coupl_type
63!> \param qmmm_coupl_type ...
64!> \param per_potentials ...
65!> \param potentials ...
66!> \param pgfs ...
67!> \param qm_cell_small ...
68!> \param mm_cell ...
69!> \param compatibility ...
70!> \param qmmm_periodic ...
71!> \param print_section ...
72!> \param eps_mm_rspace ...
73!> \param maxchrg ...
74!> \param ncp ...
75!> \param ncpl ...
76!> \par History
77!> 09.2004 created [tlaino]
78!> \author Teodoro Laino
79! **************************************************************************************************
80 SUBROUTINE qmmm_per_potential_init(qmmm_coupl_type, per_potentials, potentials, &
81 pgfs, qm_cell_small, mm_cell, compatibility, qmmm_periodic, print_section, &
82 eps_mm_rspace, maxchrg, ncp, ncpl)
83 INTEGER, INTENT(IN) :: qmmm_coupl_type
84 TYPE(qmmm_per_pot_p_type), DIMENSION(:), POINTER :: per_potentials
85 TYPE(qmmm_pot_p_type), DIMENSION(:), POINTER :: potentials
86 TYPE(qmmm_gaussian_p_type), DIMENSION(:), POINTER :: pgfs
87 TYPE(cell_type), POINTER :: qm_cell_small, mm_cell
88 LOGICAL, INTENT(IN) :: compatibility
89 TYPE(section_vals_type), POINTER :: qmmm_periodic, print_section
90 REAL(kind=dp), INTENT(IN) :: eps_mm_rspace, maxchrg
91 INTEGER, INTENT(IN) :: ncp(3), ncpl(3)
92
93 INTEGER :: i, idim, ig, ig_start, iw, ix, iy, iz, &
94 k, kmax(3), n_rep_real(3), &
95 n_rep_real_val, ncoarsel, ncoarset, &
96 ndim, output_unit
97 INTEGER, DIMENSION(:), POINTER :: mm_atom_index
98 REAL(kind=dp) :: ak, alpha, box(3), fac(3), fs, g, g2, &
99 gk, gmax, mymaxradius, npl, npt, &
100 prefactor, rc, rc2, rmax, tmp, vec(3), &
101 vol
102 REAL(kind=dp), DIMENSION(:), POINTER :: gx, gy, gz, lg
103 TYPE(cp_logger_type), POINTER :: logger
104 TYPE(qmmm_gaussian_type), POINTER :: pgf
105
106 NULLIFY (lg, gx, gy, gz)
107 ncoarset = product(ncp)
108 ncoarsel = product(ncpl)
109 logger => cp_get_default_logger()
110 output_unit = cp_logger_get_default_io_unit(logger)
111 rmax = sqrt(mm_cell%hmat(1, 1)**2 + &
112 mm_cell%hmat(2, 2)**2 + &
113 mm_cell%hmat(3, 3)**2)
114 CALL section_vals_val_get(qmmm_periodic, "GMAX", r_val=gmax)
115 CALL section_vals_val_get(qmmm_periodic, "REPLICA", i_val=n_rep_real_val)
116 fac = 2.0e0_dp*pi/(/mm_cell%hmat(1, 1), mm_cell%hmat(2, 2), mm_cell%hmat(3, 3)/)
117 kmax = ceiling(gmax/fac)
118 vol = mm_cell%hmat(1, 1)* &
119 mm_cell%hmat(2, 2)* &
120 mm_cell%hmat(3, 3)
121 ndim = (kmax(1) + 1)*(2*kmax(2) + 1)*(2*kmax(3) + 1)
122 ig_start = 1
123 n_rep_real = n_rep_real_val
124 IF (compatibility .AND. (qmmm_coupl_type == do_qmmm_gauss)) ig_start = 2
125
126 cpassert(.NOT. ASSOCIATED(per_potentials))
127 ALLOCATE (per_potentials(SIZE(pgfs)))
128 cpassert(SIZE(pgfs) == SIZE(potentials))
129 potential_type: DO k = 1, SIZE(pgfs)
130
131 rc = pgfs(k)%pgf%Elp_Radius
132 ALLOCATE (per_potentials(k)%Pot)
133 SELECT CASE (qmmm_coupl_type)
135 ! Not yet implemented for this case
136 cpabort("")
138 ALLOCATE (lg(ndim))
139 ALLOCATE (gx(ndim))
140 ALLOCATE (gy(ndim))
141 ALLOCATE (gz(ndim))
142 END SELECT
143
144 lg = 0.0_dp
145 gx = 0.0_dp
146 gy = 0.0_dp
147 gz = 0.0_dp
148
149 SELECT CASE (qmmm_coupl_type)
151 ! Not yet implemented for this case
152 cpabort("")
154 pgf => pgfs(k)%pgf
155 idim = 0
156 DO ix = 0, kmax(1)
157 DO iy = -kmax(2), kmax(2)
158 DO iz = -kmax(3), kmax(3)
159 idim = idim + 1
160 IF (ix == 0 .AND. iy == 0 .AND. iz == 0) THEN
161 DO ig = ig_start, pgf%number_of_gaussians
162 gk = pgf%Gk(ig)
163 ak = pgf%Ak(ig)*pi**(3.0_dp/2.0_dp)*gk**3.0_dp
164 lg(idim) = lg(idim) - ak
165 END DO
166 ELSE
167 fs = 2.0_dp; IF (ix == 0) fs = 1.0_dp
168 vec = fac*(/real(ix, kind=dp), real(iy, kind=dp), real(iz, kind=dp)/)
169 g2 = dot_product(vec, vec)
170 rc2 = rc*rc
171 g = sqrt(g2)
172 IF (qmmm_coupl_type == do_qmmm_gauss) THEN
173 lg(idim) = 4.0_dp*pi/g2*exp(-(g2*rc2)/4.0_dp)
174 ELSEIF (qmmm_coupl_type == do_qmmm_swave) THEN
175 tmp = 4.0_dp/rc2
176 lg(idim) = 4.0_dp*pi*tmp**2/(g2*(g2 + tmp)**2)
177 END IF
178 DO ig = ig_start, pgf%number_of_gaussians
179 gk = pgf%Gk(ig)
180 ak = pgf%Ak(ig)*pi**(3.0_dp/2.0_dp)*gk**3.0_dp
181 lg(idim) = lg(idim) - ak*exp(-(g*gk)**2.0_dp/4.0_dp)
182 END DO
183 END IF
184 lg(idim) = fs*lg(idim)*1.0_dp/vol
185 gx(idim) = fac(1)*real(ix, kind=dp)
186 gy(idim) = fac(2)*real(iy, kind=dp)
187 gz(idim) = fac(3)*real(iz, kind=dp)
188 END DO
189 END DO
190 END DO
191
192 IF (all(n_rep_real == -1)) THEN
193 mymaxradius = 0.0_dp
194 DO i = 1, pgf%number_of_gaussians
195 IF (pgf%Gk(i) /= 0.0_dp) THEN
196 alpha = 1.0_dp/pgf%Gk(i)
197 alpha = alpha*alpha
198 prefactor = pgf%Ak(i)*maxchrg
199 mymaxradius = max(mymaxradius, exp_radius(0, alpha, eps_mm_rspace, prefactor, rlow=mymaxradius))
200 END IF
201 END DO
202 box(1) = (qm_cell_small%hmat(1, 1) - mm_cell%hmat(1, 1))/2.0_dp
203 box(2) = (qm_cell_small%hmat(2, 2) - mm_cell%hmat(2, 2))/2.0_dp
204 box(3) = (qm_cell_small%hmat(3, 3) - mm_cell%hmat(3, 3))/2.0_dp
205 IF (any(box > 0.0_dp)) THEN
206 cpabort("")
207 END IF
208 n_rep_real(1) = ceiling((box(1) + mymaxradius)/mm_cell%hmat(1, 1))
209 n_rep_real(2) = ceiling((box(2) + mymaxradius)/mm_cell%hmat(2, 2))
210 n_rep_real(3) = ceiling((box(3) + mymaxradius)/mm_cell%hmat(3, 3))
211 END IF
212
213 CASE DEFAULT
214 DEALLOCATE (per_potentials(k)%Pot)
215 IF (output_unit > 0) WRITE (output_unit, '(A)') " QMMM Periodic Potential - not Initialized!"
216 cycle potential_type
217 END SELECT
218
219 NULLIFY (mm_atom_index)
220 ALLOCATE (mm_atom_index(SIZE(potentials(k)%pot%mm_atom_index)))
221 mm_atom_index = potentials(k)%pot%mm_atom_index
222
223 NULLIFY (per_potentials(k)%Pot%LG, per_potentials(k)%Pot%mm_atom_index, &
224 per_potentials(k)%Pot%gx, per_potentials(k)%Pot%gy, per_potentials(k)%Pot%gz)
225 CALL qmmm_per_pot_type_create(per_potentials(k)%Pot, lg=lg, gx=gx, gy=gy, gz=gz, &
226 gmax=gmax, kmax=kmax, n_rep_real=n_rep_real, &
227 fac=fac, mm_atom_index=mm_atom_index, &
228 mm_cell=mm_cell, &
229 qmmm_per_section=qmmm_periodic, print_section=print_section)
230
231 iw = cp_print_key_unit_nr(logger, print_section, "PERIODIC_INFO", &
232 extension=".log")
233 IF (iw > 0) THEN
234 npt = real(ncoarset, kind=dp)*real(ndim, kind=dp)*real(SIZE(mm_atom_index), kind=dp)
235 npl = real(ncoarsel, kind=dp)*real(ndim, kind=dp)*real(SIZE(mm_atom_index), kind=dp)
236 WRITE (unit=iw, fmt="(/,T2,A)") repeat("-", 79)
237 WRITE (unit=iw, fmt="(T2,A,T20,A,T80,A)") "-", "QMMM PERIODIC BOUNDARY CONDITION INFO", "-"
238 WRITE (unit=iw, fmt="(T2,A)") repeat("-", 79)
239 WRITE (unit=iw, fmt="(T2,A,T10,A,F15.6,T50,A,3I5,T80,A)") "-", "RADIUS =", rc, "REPLICA =", n_rep_real, "-"
240 WRITE (unit=iw, fmt="(T2,A,T10,A,F15.6,T50,A,I15,T80,A)") "-", "MINGVAL =", minval(abs(lg)), &
241 "GPOINTS =", ndim, "-"
242 WRITE (unit=iw, fmt="(T2,A,T10,A,3I5,T50,A,3I5,T80,A)") "-", "NCOARSL =", ncpl, &
243 "NCOARST =", ncp, "-"
244 WRITE (unit=iw, fmt="(T2,A,T10,A,F15.0,T50,A,F15.0,T80,A)") "-", "NFLOP-L ~", npl, &
245 "NFLOP-T ~", npt, "-"
246 WRITE (unit=iw, fmt="(T2,A)") repeat("-", 79)
247 END IF
248 CALL cp_print_key_finished_output(iw, logger, print_section, &
249 "PERIODIC_INFO")
250
251 END DO potential_type
252
253 END SUBROUTINE qmmm_per_potential_init
254
255! **************************************************************************************************
256!> \brief Creates the qmmm_pot_type structure
257!> \param Pot ...
258!> \param LG ...
259!> \param gx ...
260!> \param gy ...
261!> \param gz ...
262!> \param GMax ...
263!> \param Kmax ...
264!> \param n_rep_real ...
265!> \param Fac ...
266!> \param mm_atom_index ...
267!> \param mm_cell ...
268!> \param qmmm_per_section ...
269!> \param print_section ...
270!> \par History
271!> 09.2004 created [tlaino]
272!> \author Teodoro Laino
273! **************************************************************************************************
274 SUBROUTINE qmmm_per_pot_type_create(Pot, LG, gx, gy, gz, GMax, Kmax, n_rep_real, &
275 Fac, mm_atom_index, mm_cell, qmmm_per_section, print_section)
276 TYPE(qmmm_per_pot_type), POINTER :: pot
277 REAL(kind=dp), DIMENSION(:), POINTER :: lg, gx, gy, gz
278 REAL(kind=dp), INTENT(IN) :: gmax
279 INTEGER, INTENT(IN) :: kmax(3), n_rep_real(3)
280 REAL(kind=dp), INTENT(IN) :: fac(3)
281 INTEGER, DIMENSION(:), POINTER :: mm_atom_index
282 TYPE(cell_type), POINTER :: mm_cell
283 TYPE(section_vals_type), POINTER :: qmmm_per_section, print_section
284
285 INTEGER :: npts(3)
286 INTEGER, DIMENSION(:), POINTER :: ngrids
287 REAL(kind=dp) :: hmat(3, 3)
288 TYPE(section_vals_type), POINTER :: grid_print_section
289
290 pot%LG => lg
291 pot%gx => gx
292 pot%gy => gy
293 pot%gz => gz
294 pot%mm_atom_index => mm_atom_index
295 pot%Gmax = gmax
296 pot%Kmax = kmax
297 pot%n_rep_real = n_rep_real
298 pot%Fac = fac
299 !
300 ! Setting Up Fit Procedure
301 !
302 NULLIFY (pot%pw_grid)
303 NULLIFY (pot%pw_pool)
304 NULLIFY (pot%TabLR, ngrids)
305 CALL section_vals_val_get(qmmm_per_section, "ngrids", i_vals=ngrids)
306 npts = ngrids
307 hmat = mm_cell%hmat
308
309 grid_print_section => section_vals_get_subs_vals(print_section, "GRID_INFORMATION")
310 CALL setup_ewald_spline(pw_grid=pot%pw_grid, pw_pool=pot%pw_pool, coeff=pot%TabLR, &
311 lg=lg, gx=gx, gy=gy, gz=gz, hmat=hmat, npts=npts, param_section=qmmm_per_section, &
312 tag="qmmm", print_section=grid_print_section)
313
314 END SUBROUTINE qmmm_per_pot_type_create
315
316! **************************************************************************************************
317!> \brief Initialize the QMMM Ewald potential needed for QM-QM Coupling using
318!> point charges
319!> \param ewald_env ...
320!> \param ewald_pw ...
321!> \param qmmm_coupl_type ...
322!> \param mm_cell ...
323!> \param para_env ...
324!> \param qmmm_periodic ...
325!> \param print_section ...
326!> \par History
327!> 10.2014 created [JGH]
328!> \author JGH
329! **************************************************************************************************
330 SUBROUTINE qmmm_ewald_potential_init(ewald_env, ewald_pw, qmmm_coupl_type, mm_cell, para_env, &
331 qmmm_periodic, print_section)
332 TYPE(ewald_environment_type), POINTER :: ewald_env
333 TYPE(ewald_pw_type), POINTER :: ewald_pw
334 INTEGER, INTENT(IN) :: qmmm_coupl_type
335 TYPE(cell_type), POINTER :: mm_cell
336 TYPE(mp_para_env_type), POINTER :: para_env
337 TYPE(section_vals_type), POINTER :: qmmm_periodic, print_section
338
339 INTEGER :: ewald_type, gmax(3), iw, o_spline, ounit
340 LOGICAL :: do_multipoles
341 REAL(kind=dp) :: alpha, rcut
342 TYPE(cp_logger_type), POINTER :: logger
343 TYPE(section_vals_type), POINTER :: ewald_print_section, ewald_section, &
344 poisson_section
345
346 logger => cp_get_default_logger()
347 ounit = cp_logger_get_default_io_unit(logger)
348 cpassert(.NOT. ASSOCIATED(ewald_env))
349 cpassert(.NOT. ASSOCIATED(ewald_pw))
350
351 ! Create Ewald environments
352 poisson_section => section_vals_get_subs_vals(qmmm_periodic, "POISSON")
353 ALLOCATE (ewald_env)
354 CALL ewald_env_create(ewald_env, para_env)
355 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
356 ewald_section => section_vals_get_subs_vals(poisson_section, "EWALD")
357 CALL read_ewald_section(ewald_env, ewald_section)
358 ewald_print_section => section_vals_get_subs_vals(print_section, "GRID_INFORMATION")
359 ALLOCATE (ewald_pw)
360 CALL ewald_pw_create(ewald_pw, ewald_env, mm_cell, mm_cell, print_section=ewald_print_section)
361
362 CALL ewald_env_get(ewald_env, ewald_type=ewald_type, do_multipoles=do_multipoles, &
363 gmax=gmax, o_spline=o_spline, alpha=alpha, rcut=rcut)
364 IF (do_multipoles) &
365 cpabort("No multipole force fields allowed in QM-QM Ewald long range correction")
366
367 SELECT CASE (qmmm_coupl_type)
368 CASE (do_qmmm_coulomb)
369 cpabort("QM-QM long range correction not possible with COULOMB coupling")
370 CASE (do_qmmm_pcharge)
371 ! OK
373 cpabort("QM-QM long range correction not possible with GAUSS/SWAVE coupling")
374 CASE DEFAULT
375 ! We should never get to this point
376 cpabort("")
377 END SELECT
378
379 iw = cp_print_key_unit_nr(logger, print_section, "PERIODIC_INFO", extension=".log")
380 IF (iw > 0) THEN
381 WRITE (unit=iw, fmt="(/,T2,A)") repeat("-", 79)
382 WRITE (unit=iw, fmt="(T2,A,T20,A,T80,A)") "-", "QMMM PERIODIC BOUNDARY CONDITION INFO", "-"
383 WRITE (unit=iw, fmt="(T2,A,T25,A,T80,A)") "-", "QM-QM Long Range Correction", "-"
384 WRITE (unit=iw, fmt="(T2,A)") repeat("-", 79)
385 SELECT CASE (ewald_type)
386 CASE (do_ewald_none)
387 cpabort("QM-QM long range correction not compatible with Ewald=NONE")
388 CASE (do_ewald_pme)
389 cpabort("QM-QM long range correction not possible with Ewald=PME")
390 CASE (do_ewald_ewald)
391 cpabort("QM-QM long range correction not possible with Ewald method")
392 CASE (do_ewald_spme)
393 WRITE (unit=iw, fmt="(T2,A,T35,A,T75,A,T80,A)") "-", "Ewald type", "SPME", "-"
394 WRITE (unit=iw, fmt="(T2,A,T35,A,T61,3I6,T80,A)") "-", "GMAX values", gmax, "-"
395 WRITE (unit=iw, fmt="(T2,A,T35,A,T73,I6,T80,A)") "-", "Spline order", o_spline, "-"
396 WRITE (unit=iw, fmt="(T2,A,T35,A,T67,F12.4,T80,A)") "-", "Alpha value", alpha, "-"
397 WRITE (unit=iw, fmt="(T2,A,T35,A,T67,F12.4,T80,A)") "-", "Real space cutoff value", rcut, "-"
398 END SELECT
399 WRITE (unit=iw, fmt="(T2,A)") repeat("-", 79)
400 END IF
401 CALL cp_print_key_finished_output(iw, logger, print_section, "PERIODIC_INFO")
402
403 END SUBROUTINE qmmm_ewald_potential_init
404
405END MODULE qmmm_per_elpot
All kind of helpful little routines.
Definition ao_util.F:14
real(kind=dp) function, public exp_radius(l, alpha, threshold, prefactor, epsabs, epsrel, rlow)
The radius of a primitive Gaussian function for a given threshold is calculated. g(r) = prefactor*r**...
Definition ao_util.F:96
Handles all functions related to the CELL.
Definition cell_types.F:15
various routines to log and control the output. The idea is that decisions about where to log should ...
integer function, public cp_logger_get_default_io_unit(logger)
returns the unit nr for the ionode (-1 on all other processors) skips as well checks if the procs cal...
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,...
subroutine, public ewald_env_set(ewald_env, ewald_type, alpha, epsilon, eps_pol, gmax, ns_max, precs, o_spline, para_env, poisson_section, interaction_cutoffs, cell_hmat)
Purpose: Set the EWALD environment.
subroutine, public ewald_env_create(ewald_env, para_env)
allocates and intitializes a ewald_env
subroutine, public read_ewald_section(ewald_env, ewald_section)
Purpose: read the EWALD section.
subroutine, public ewald_env_get(ewald_env, ewald_type, alpha, eps_pol, epsilon, gmax, ns_max, o_spline, group, para_env, poisson_section, precs, rcut, do_multipoles, max_multipole, do_ipol, max_ipol_iter, interaction_cutoffs, cell_hmat)
Purpose: Get the EWALD environment.
subroutine, public ewald_pw_create(ewald_pw, ewald_env, cell, cell_ref, print_section)
creates the structure ewald_pw_type
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)
Setup of the G-space Ewald Term Spline Coefficients.
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public do_qmmm_pcharge
integer, parameter, public do_qmmm_coulomb
integer, parameter, public do_qmmm_swave
integer, parameter, public do_qmmm_gauss
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
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
real(kind=dp), dimension(0:maxfac), parameter, public fac
Interface to the message passing library MPI.
functions related to the poisson solver on regular grids
integer, parameter, public do_ewald_pme
integer, parameter, public do_ewald_ewald
integer, parameter, public do_ewald_none
integer, parameter, public do_ewald_spme
Sets the typo for the gaussian treatment of the qm/mm interaction.
Setting up the potential for QM/MM periodic boundary conditions calculations.
subroutine, public qmmm_per_potential_init(qmmm_coupl_type, per_potentials, potentials, pgfs, qm_cell_small, mm_cell, compatibility, qmmm_periodic, print_section, eps_mm_rspace, maxchrg, ncp, ncpl)
Initialize the QMMM potential stored on vector, according the qmmm_coupl_type.
subroutine, public qmmm_ewald_potential_init(ewald_env, ewald_pw, qmmm_coupl_type, mm_cell, para_env, qmmm_periodic, print_section)
Initialize the QMMM Ewald potential needed for QM-QM Coupling using point charges.
Type defining parameters related to the simulation cell.
Definition cell_types.F:55
type of a logger, at the moment it contains just a print level starting at which level it should be l...
stores all the informations relevant to an mpi environment
represent a pointer to a qmmm_gaussian_type, to be able to create arrays of pointers