49#include "./base/base_uses.f90"
54 LOGICAL,
PRIVATE,
PARAMETER :: debug_this_module = .true.
55 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'qmmm_per_elpot'
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
87 TYPE(
cell_type),
POINTER :: qm_cell_small, mm_cell
88 LOGICAL,
INTENT(IN) :: compatibility
90 REAL(kind=
dp),
INTENT(IN) :: eps_mm_rspace, maxchrg
91 INTEGER,
INTENT(IN) :: ncp(3), ncpl(3)
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, &
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), &
102 REAL(kind=
dp),
DIMENSION(:),
POINTER :: gx, gy, gz, lg
106 NULLIFY (lg, gx, gy, gz)
107 ncoarset = product(ncp)
108 ncoarsel = product(ncpl)
111 rmax = sqrt(mm_cell%hmat(1, 1)**2 + &
112 mm_cell%hmat(2, 2)**2 + &
113 mm_cell%hmat(3, 3)**2)
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)* &
121 ndim = (kmax(1) + 1)*(2*kmax(2) + 1)*(2*kmax(3) + 1)
123 n_rep_real = n_rep_real_val
124 IF (compatibility .AND. (qmmm_coupl_type ==
do_qmmm_gauss)) ig_start = 2
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)
131 rc = pgfs(k)%pgf%Elp_Radius
132 ALLOCATE (per_potentials(k)%Pot)
133 SELECT CASE (qmmm_coupl_type)
149 SELECT CASE (qmmm_coupl_type)
157 DO iy = -kmax(2), kmax(2)
158 DO iz = -kmax(3), kmax(3)
160 IF (ix == 0 .AND. iy == 0 .AND. iz == 0)
THEN
161 DO ig = ig_start, pgf%number_of_gaussians
163 ak = pgf%Ak(ig)*
pi**(3.0_dp/2.0_dp)*gk**3.0_dp
164 lg(idim) = lg(idim) - ak
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)
173 lg(idim) = 4.0_dp*
pi/g2*exp(-(g2*rc2)/4.0_dp)
176 lg(idim) = 4.0_dp*
pi*tmp**2/(g2*(g2 + tmp)**2)
178 DO ig = ig_start, pgf%number_of_gaussians
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)
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)
192 IF (all(n_rep_real == -1))
THEN
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)
198 prefactor = pgf%Ak(i)*maxchrg
199 mymaxradius = max(mymaxradius,
exp_radius(0, alpha, eps_mm_rspace, prefactor, rlow=mymaxradius))
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
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))
214 DEALLOCATE (per_potentials(k)%Pot)
215 IF (output_unit > 0)
WRITE (output_unit,
'(A)')
" QMMM Periodic Potential - not Initialized!"
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
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, &
229 qmmm_per_section=qmmm_periodic, print_section=print_section)
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)
251 END DO potential_type
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)
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
286 INTEGER,
DIMENSION(:),
POINTER :: ngrids
287 REAL(kind=
dp) :: hmat(3, 3)
294 pot%mm_atom_index => mm_atom_index
297 pot%n_rep_real = n_rep_real
302 NULLIFY (pot%pw_grid)
303 NULLIFY (pot%pw_pool)
304 NULLIFY (pot%TabLR, ngrids)
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)
314 END SUBROUTINE qmmm_per_pot_type_create
331 qmmm_periodic, print_section)
334 INTEGER,
INTENT(IN) :: qmmm_coupl_type
339 INTEGER :: ewald_type, gmax(3), iw, o_spline, ounit
340 LOGICAL :: do_multipoles
341 REAL(kind=
dp) :: alpha, rcut
348 cpassert(.NOT.
ASSOCIATED(ewald_env))
349 cpassert(.NOT.
ASSOCIATED(ewald_pw))
355 CALL ewald_env_set(ewald_env, poisson_section=poisson_section)
360 CALL ewald_pw_create(ewald_pw, ewald_env, mm_cell, mm_cell, print_section=ewald_print_section)
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)
365 cpabort(
"No multipole force fields allowed in QM-QM Ewald long range correction")
367 SELECT CASE (qmmm_coupl_type)
369 cpabort(
"QM-QM long range correction not possible with COULOMB coupling")
373 cpabort(
"QM-QM long range correction not possible with GAUSS/SWAVE coupling")
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)
387 cpabort(
"QM-QM long range correction not compatible with Ewald=NONE")
389 cpabort(
"QM-QM long range correction not possible with Ewald=PME")
391 cpabort(
"QM-QM long range correction not possible with Ewald method")
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,
"-"
399 WRITE (unit=iw, fmt=
"(T2,A)") repeat(
"-", 79)
All kind of helpful little routines.
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**...
Handles all functions related to the CELL.
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.
Defines the basic variable types.
integer, parameter, public dp
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.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
to build arrays of pointers
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