Function reference¶
Throughout this document the sum \(\sum_\tau\) implies the periodic sum \(\sum_{p \in \mathbb Z^3}\) with \(\tau = (p_1L_1, p_2L_2, p_3L_3)\).
SE¶
Ewald summation for the 3-periodic electrostatic potential
-
SE.se3p_fourier_space(x, f, opt)¶ initialize time array
-
SE.se3p_precomp(x, opt)¶ SPECTRAL EWALD 3P, pre-computation of FGG vectors Fast Ewald method for electrostatic potential calculation, k-space part.
-
SE.spectral_ewald(eval_idx, x, q, xi, opt)¶ SPECTRAL EWALD Fast Ewald method for electrostatic potential calculation, k-space part.
phi = spectral_ewald(idx, x, q, xi, box, opt)
Returns: phi – k-space part of periodic potential
Parameters: - idx – evaluation indices (particle indices)
- x – particle positions (N-by-3)
- q – particle charge (N-by-1)
- xi – Ewald parameter
- opt – Options:
- opt.M – grid size, eg. [31 31 31] (required)
- opt.P – Gaussian support (required)
- opt.box – Periodic cell, eg. [1 1 1]
- opt.m – Gaussian shape function (default: m=m(P))
- opt.w – Gaussian width (default: w=h*P/2)
SE2P¶
Ewald summation for the 2-periodic electrostatic potential.
SE1P¶
Ewald summation for the 1-periodic electrostatic potential.
SE0P¶
Ewald summation for the free-space electrostatic (stokeslet, stresslet, rotlet) potential.
SE_kaiser¶
Ewald summation for the electrostatic potential with arbitrary periodicity using Barnett-Magland kernel.
SE2P_Stokes¶
Ewald summation for the 2-periodic stokeslet potential.
-
SE2P_Stokes.SE2P_Stokes(eval_idx, x, f, xi, opt)¶ SPECTRAL EWALD 2DP Fast Ewald method for electrostatic potential calculation, k-space part.
phi = spectral_ewald_2dp(idx, x, q, xi, opt)
Dag Lindbo, dag@kth.se, Jul. 2011
-
SE2P_Stokes.SE2P_Stokes_pre(x, xi, opt)¶ SPECTRAL EWALD 2DP, pre-computation of FGG vectors
Dag Lindbo, dag@kth.se, Jul. 2011
-
SE2P_Stokes.parse_params(opt)¶ check that we have all mandatory options
SE_fast_gridding¶
Routines for fast gaussian gridding (FGG) and fast Kaiser gridding (FKG), which is central to Spectral Ewald.
-
SE_fast_gridding.SE_FGG_precomp(x, xi, opt)¶ SPECTRAL EWALD 3P, pre-computation of FGG vectors Fast Ewald method for electrostatic potential calculation, k-space part.
-
SE_fast_gridding.SE_fg_extend_fcn(F, opt)¶ parameters and constants
-
SE_fast_gridding.SE_fg_int(x, U, opt)¶ - u = SE_fg_int(x, U, opt)
- Return integration at x of on-grid values U
- u = SE_fg_int(x, {U1, …, UN}, opt)
- Return integration at x of on-grid values U1,…,UN u has N columns
-
SE_fast_gridding.SE_fg_int_extended(x, U, opt)¶ parameters and constants
-
SE_fast_gridding.parse_params(opt)¶ Parse parameter struct for fast Gaussian gridding
Mandatory parameters:
Parameters: - opt.M – grid size (1 x 3)
- opt.P – Gaussian width
- opt.box – box size (1 x 3)
SE_Rotlet¶
Ewald summation for the 3-periodic rotlet potential
where
-
SE_Rotlet.SE_Rotlet(xe, x, f, xi, opt)¶ - u = SE_Rotlet(x_targets, x_sources, f_sources, xi, opt)
- Returns the rotlet potential
- [U1, U2, U3] = SE_Rotlet(x_targets, x_sources, f_sources, xi, opt)
- Returns the grid with Fourier coefficients, \(U_i\in\mathbb{C}^{M_1\times M_2 \times M_3}\)
To compute u from U1,U2,U3:
F{1} = real( ifftn( ifftshift( U1 ))); F{2} = real( ifftn( ifftshift( U2 ))); F{3} = real( ifftn( ifftshift( U3 ))); u = SE_fg_int(xe, F, opt);
opt must be identical to that passed to SE_Rotlet
TODO: Add possibility of passing precomputed structures
-
SE_Rotlet.rotlet_direct_fd(xe, x, f, xi, L, kmax)¶ Direct summation of rotlet k-space part.
Parameters: - xe – target positions, M-by-3
- x – source positions, N-by-3
- f – source strengths, N-by-3
- xi – Ewald parameter,
- L – box size, 1-by-3
- kmax – k-space truncation
-
SE_Rotlet.rotlet_direct_real(idx, x, f, xi, L, varargin)¶ Direct summation of rotlet real space part
- phi = rotlet_direct_real(idx, x, t, xi, box, ‘layers’, M, ‘tol’, TOL);
- Compute interactions from all particles within M layers, stop at tolerance TOL.
- phi = rotlet_direct_real(idx, x, t, xi, box, ‘mode’,’cutoff’, ‘rc’, rc);
- Compute interactions from all particles within cutoff radius rc.
-
SE_Rotlet.rotlet_direct_rsrc(xe, x, t, opt)¶ Rotlet real space direct summation with cutoff radius, MEX implementation.
Parameters: - xe – target locations
- x – source locaitons
- t – source strengths
- opt – Ewald options
- opt.box – box size
- opt.xi – Ewald parameter \(\xi\)
- opt.rc – Cutoff radius
-
SE_Rotlet.rotlet_ewald_sum(xe, x, t, opt)¶ Ewald summation for 3P rotlet.
Parameters: - xe – target locations
- x – source locaitons
- t – source strengths
- opt – Ewald options
- opt.box – box size
- opt.xi – Ewald parameter \(\xi\)
- opt.M – Fourier space grid size
- opt.P – Gaussian width
- opt.rc – Real space cutoff radius
SE_Stokes¶
Ewald summation for the 3-periodic stokeslet potential
where
-
SE_Stokes.SE_Stokes(eval_idx, x, f, xi, opt)¶ Compute Fourier space part of Ewald sum for periodic stokeslet potential.
- u = SE_Stokes(eval_idx,x,f,xi,opt)
- Return potential
- [U1, U2, U3] = SE_Stokes(eval_idx,x,f,xi,opt)
- Return Fourier coefficients
Parameters: - eval_idx – index of source locations where potential should be evaluated
- x – source locations (Nx3)
- f – source strengths (Nx3)
- xi – Ewald paramter
- opt – Ewald options
- opt.M – grid size (M1, M2, M3)
- opt.P – Gaussian width
- opt.box – Box size (L1, L2, L3)
Returns: phi – Fourier space potential
Returns: U1,U2,U3 – Fourier space potential
-
SE_Stokes.SE_Stokes_ext(xe, xs, f, xi, opt)¶ Evaluate Stokeslet at external points
-
SE_Stokes.SE_Stokes_par(eval_idx, x, f, xi, opt)¶ parfor:ed version
SE_Stresslet¶
Ewald summation for the 3-periodic stresslet potential
where
-
SE_Stresslet.SE_Stresslet(eval_idx, x, f, n, xi, opt, varargin)¶ Compute Fourier space part of Ewald sum for periodic stresslet potential.
- phi = SE_Stresslet(eval_idx,x,f,n,xi,opt)
- Return potential
- [U1, U2, U3] = SE_Stresslet(eval_idx,x,f,n,xi,opt)
- Return Fourier coefficients
Parameters: - eval_idx – index of source locations where potential should be evaluated
- x – source locations (Nx3)
- f – source strengths (Nx3)
- n – source normals (Nx3)
- xi – Ewald paramter
- opt – Ewald options:
- opt.M – grid size (M1, M2, M3)
- opt.P – Gaussian width
- opt.box – Box size (L1, L2, L3)
- varargin=SE_static – Precomputations from SE_Stresslet_pre() to use SIMD FGG code.
Returns: phi – Fourier space potential
Returns: U1,U2,U3 – Fourier space coefficients
-
SE_Stresslet.SE_Stresslet_pre(x, xi, opt)¶ Precompute data for fast Gaussian gridding to enable optimized kernels
Returns: SE_static
-
SE_Stresslet.generate_state(N, L, varargin)¶ Returning x: Nx3 matrix with coords inside box L, f: Nx3 matrix with point forces. nvec: Nx3 matrix with associated normal vectors.
-
SE_Stresslet.stresslet_direct_fd(idx, x, f, nvec, xi, L, kmax)¶ Fourier space Ewald sum for stresslet, using (very slow) direct summation.
Parameters: kmax – Fourier space truncation
-
SE_Stresslet.stresslet_direct_fd_zero(idx, xvec, fvec, nvec, L)¶ Component for \(k=0\) of stresslet Fourier sum.
Parameters: - idx – index of target positions (Nx1)
- xvec – source positions (Nx3)
- fvec – source strengths (Nx3)
- nvec – normal vectors (Nx3)
- L – box size (1x3)
Returns: phi_k0 – (Nx3)
-
SE_Stresslet.stresslet_direct_real(idx, x, f, nvec, xi, L, nbox, TOL, varargin)¶ Real space Ewald sum for stresslet, using (very slow) direct summation.
Parameters: - nbox – Number of periodic replications in each direction.
- TOL – Break when truncation error < TOL.
- varargin=eval_x – Evaluate at eval_x instead of x(idx,:)
-
SE_Stresslet.stresslet_direct_real_fast(idx, x, f, nvec, xi, L, nbox, rc, varargin)¶ Ewald summation for the stresslet – Real space part. Fast in the sense that it saves interactions as matrix for subsequent iterations Sums over layers
- [phi A] = stresslet_direct_real_fast( idx, x, f, nvec, xi, L, nbox, rc, [A, sing_sub, ns_quad_mtrx])
Parameters: - x – positions :param (N-by-3)
- f – source strengths (N-by-3)
- nvec – normal vectors (N-by-3)
- xi – Ewald parameter
- L – Box size (3)
- nbox – periodic replications
- rc – cutoff radius
- A – Precomputed matrix
- sing_sub – Use singularity subtraction
- ns_quad_mtrx – Correction matrices for nearly singular quadrature
-
SE_Stresslet.stresslet_op_fd(k, xi)¶ Returns imaginary part of Bi(k, xi) ie B = 1i*Bi;
-
SE_Stresslet.stresslet_real_rc(x, q, nvec, xi, box, rc, varargin)¶ Stresslet real space summation with truncation radius rc
Usage:
- [res] = stresslet_real_rc( x, f, nvec, xi, box, rc);
- [res AMAT] = stresslet_real_rc( x, f, nvec, xi, box, rc, [], sing_sub, ns_quad_mtrx);
- [res AMAT] = stresslet_real_rc( x, f, nvec, xi, box, rc, AMAT, sing_sub, ns_quad_mtrx);
- [res AMAT R C V PER] = stresslet_real_rc( x, f, nvec, xi, box, rc);
Parameters: - sing_sub – singularity subtraction, default 0. Only for matrix comp.
- rc – cutoff radius, must be <= min(box)/2