pygbe package

Submodules

pygbe.class_initialization module

It contains the necessary functions to set up the surface to be solved.

pygbe.class_initialization.initialize_field(filename, param, field=None)[source]

Initialize all the regions in the surface to be solved.

Parameters:
  • filename (name of the file that contains the surface information.) –
  • param (class, parameters related to the surface.) –
  • field (dictionary with preloaded field values for programmatic) – interaction with PyGBe
Returns:

field_array

Return type:

array, contains the Field classes of each region on the surface.

pygbe.class_initialization.initialize_surface(field_array, filename, param)[source]

Initialize the surface of the molecule.

Parameters:
  • field_array (array, contains the Field classes of each region on the surface.) –
  • filename (name of the file that contains the surface information.) –
  • param (class, parameters related to the surface.) –
Returns:

surf_array – surface.

Return type:

array, contains the surface classes of each region on the

pygbe.classes module

class pygbe.classes.Event[source]

Bases: object

Class for logging like in pycuda’s cuda.Event()

record()[source]
synchronize()[source]
time_till(toc)[source]
class pygbe.classes.Field(LorY, kappa, E, coulomb, pot)[source]

Bases: object

Field class. Information about each region in the molecule.

parent

list, Pointer to “parent” surface.

child

list, Pointer to “children” surfaces.

LorY

int, 1: Laplace, 2: Yukawa.

kappa

float, inverse of Debye length.

E

float, dielectric constant.

xq

list, position of charges.

q

list, value of charges.

pot

int, 1: calculate energy on this field 0: ignore

coulomb

int, 1: perform Coulomb interaction calculation, 0: don’t do Coulomb.

# Device data
xq_gpu

list, x position of charges on GPU.

yq_gpu

list, y position of charges on GPU.

zq_gpu

list, z position of charges on GPU.

q_gpu

list, value of charges on GPU.

load_charges(qfile, REAL)[source]
class pygbe.classes.IndexConstant[source]

Bases: object

Precompute indices required for the treecode computation.

II

list, multipole order in the x-direction for the treecode.

JJ

list, multipole order in the y-direction for the treecode.

KK

list, multipole order in the z-direction for the treecode.

index_large

list, pointers to the position of multipole order i, j, k – in the multipole array, organized in a 1D array of size P*P*P. Index is given by index[i*P*P+j*P+k]

index_small

list, pointers to the position of multipole order i, j, k – in the multipole array, organized in a 1D array which is compressed with respect to index_large (does not consider combinations of i,j,k which do not have a multipole).

index

list, copy of index_small

index_ptr

list, pointer to index_small. Data in index_small is organized – in a i-major fashion (i,j,k), and index_ptr points at the position in index_small where the order i changes.

combII

array, combinatory of (I, i) where I is the maximum i multipole. – Used in coefficients of M2M.

combJJ

array, combinatory of (J, j) where J is the maximum j multipole. – Used in coefficients of M2M.

combKK

array, combinatory of (K, k) where K is the maximum k multipole. – Used in coefficients of M2M.

IImii

array, I-i where I is the maximum i multipole. – Used in exponents of M2M.

JJmjj

array, J-j where J is the maximum j multipole. – Used in exponents of M2M.

KKmkk

array, K-k where K is the maximum k multipole. – Used in exponents of M2M.

# Device data
indexDev

list, index_large on GPU.

class pygbe.classes.Parameters[source]

Bases: object

Parameters class. It contains the information of the parameters needed to run the code.

kappa

float, inverse of Debye length.

restart

int, Restart of GMRES.

tol

float, Tolerance of GMRES.

max_iter

int, Max number of GMRES iterations.

P

int, Order of Taylor expansion.

eps

int, Epsilon machine.

Nm

int, Number of terms in Taylor expansion.

NCRIT

int, Max number of targets per twig box.

theta

float, MAC criterion for treecode.

K

int, Number of Gauss points per element.

K_fine

int, Number of Gauss points per element for near singular integrals.

threshold

float, L/d criterion for semi-analytic intergrals.

Nk

int, Gauss points per side for semi-analytical integrals.

BSZ

int, CUDA block size.

Nround

int, Max size of sorted target array.

BlocksPerTwig

int, Number of CUDA blocks that fit per tree twig.

N

int, Total number of elements.

Neq

int, Total number of equations.

qe

float, Charge of an electron (1.60217646e-19).

Na

float, Avogadro’s number (6.0221415e23).

E_0

float, Vacuum dielectric constant (8.854187818e-12).

REAL

Data type.

E_field

list, Regions where energy will be calculated.

GPU

int, =1: with GPU, =0: no GPU.

class pygbe.classes.Surface(Nsurf, surf_type, phi0_file)[source]

Bases: object

Surface class. Information about the solvent excluded surface.

triangle

list, indices to triangle vertices.

vertex

list, position of vertices.

XinV

list, weights input for single layer potential.

XinK

list, weights input for double layer potential.

Xout_int

list, output vector of interior operators.

Xout_ext

list, output vector of exterior operators.

xi

list, x component of center.

yi

list, y component of center.

zi

list, z component of center.

xj

list, x component of gauss nodes.

yj

list, y component of gauss nodes.

zj

list, z component of gauss nodes.

area

list, areas of triangles.

normal

list, normal of triangles.

sglInt_int

list, singular integrals for V for internal equation.

sglInt_ext

list, singular integrals for V for external equation.

xk

list, position of gauss points on edges.

wk

list, weight of gauss points on edges.

Xsk

list, position of gauss points for near singular integrals.

Wsk

list, weight of gauss points for near singular integrals.

tree

list, tree structure.

twig

list, tree twigs.

xiSort

list, sorted x component of center.

yiSort

list, sorted y component of center.

ziSort

list, sorted z component of center.

xjSort

list, sorted x component of gauss nodes.

yjSort

list, sorted y component of gauss nodes.

zjSort

list, sorted z component of gauss nodes.

xcSort

list, sorted x component of the box centers according to – M2P_list array.

ycSort

list, sorted y component of the box centers according to – M2P_list array.

zcSort

list, sorted z component of the box centers according to – M2P_list array.

areaSort

list, sorted array of areas.

sglInt_intSort

list, sorted array of singular integrals for V for internal – equation.

sglInt_extSort

list, sorted array of singular integrals for V for external – equation.

unsort

list, array of indices to unsort targets.

triangleSort

list, sorted array of triangles.

sortTarget

list, indices to sort targets.

sortSource

list, indices to sort sources.

offsetSource

list, offsets to sorted source array.

offsetTarget

list, offsets to sorted target array.

sizeTarget

list, number of targets per twig.

offsetTwigs

list, offset to twig in P2P list array.

P2P_list

list, pointers to twigs for P2P interaction list.

offsetMlt

list, offset to multipoles in M2P list array.

M2P_list

list, pointers to boxes for M2P interaction list.

Precond

list, sparse representation of preconditioner for self – interaction block.

Ein

float, permitivitty inside surface.

Eout

float, permitivitty outside surface.

E_hat

float, ratio of Ein/Eout.

kappa_in

float, kappa inside surface.

kappa_out

float, kappa inside surface.

LorY_in

int, Laplace (1) or Yukawa (2) in inner region.

LorY_out

int, Laplace (1) or Yukawa (2) in outer region.

surf_type

int, Surface type: internal_cavity (=0), stern or – dielecric_interface (=1).

phi0

list, known surface potential (dirichlet) or derivative of – potential (neumann).

phi

list, potential on surface.

dphi

list, derivative of potential on surface.

# Device data
xiDev

list, sorted x component of center (on the GPU).

yiDev

list, sorted y component of center (on the GPU).

ziDev

list, sorted z component of center (on the GPU).

xjDev

list, sorted x component of gauss nodes (on the GPU).

yjDev

list, sorted y component of gauss nodes (on the GPU).

zjDev

list, sorted z component of gauss nodes (on the GPU).

xcDev

list, sorted x component of the box centers according to – M2P_list array (on the GPU).

ycDev

list, sorted y component of the box centers according to – M2P_list array (on the GPU).

zcDev

list, sorted z component of the box centers according to – M2P_list array (on the GPU).

areaDev

list, areas of triangles (on the GPU).

sglInt_intDev

list, singular integrals for V for internal equation (on the – GPU).

sglInt_extDev

list, singular integrals for V for external equation (on the – GPU).

vertexDev

list, sorted vertex of the triangles.

sizeTarDev

list, number of targets per twig (on the GPU).

offSrcDev

list, offsets to sorted source array (on the GPU).

offMltDev

list, offset to multipoles in M2P list array (on the GPU).

offTwgDev

list, offset to twig in P2P list array (on the GPU).

M2P_lstDev

list, pointers to boxes for M2P interaction list (on the GPU).

P2P_lstDev

list, pointers to twigs for P2P interaction list (on the GPU).

xkDev

list, position of gauss points on edges (on the GPU).

wkDev

list, weight of gauss points on edges (on the GPU).

XskDev

list, position of gauss points for near singular integrals – (on the GPU).

WskDev

list, weight of gauss points for near singular integrals (on – the GPU).

kDev

list, quadrature number of each quadrature point, in order. – (on the GPU)

calc_centers()[source]

Calculate the centers of each element of the surface

calc_distance(param)[source]

Calculate the radius spanned by the points on the surface

calc_norms()[source]

Calculate the surface normal vector

define_regions(field_array, i)[source]

Look for regions inside/outside

define_surface(files, param)[source]

Load the vertices and triangles that define the molecule surface

fill_phi(phi, s_start=0)[source]

Place the result vector on surface structure.

Parameters:
  • phi (array, result vector.) –
  • s_start (int, offset to grab the corresponding section of phi) –
fill_surface(param)[source]

Fill the surface with all the necessary information to solve it.

-Set the Gauss points. -Generate tree, compute the indices and precompute terms for M2M. -Generate preconditioner. -Compute the diagonal integral for internal and external equations.

Parameters:param (class, parameters related to the surface we are studying.) –
generate_preconditioner()[source]

Generate preconditioner

Notes

Uses block-diagonal preconditioner [3]

[3]Altman, M. D., Bardhan, J. P., White, J. K., & Tidor, B. (2009). Accurate solution of multi‐region continuum biomolecule electrostatic problems using the linearized Poisson–Boltzmann equation with curved boundary elements. Journal of computational chemistry, 30(1), 132-153.
get_gauss_points(n)[source]

Get the Gauss points for far away integrals.

Parameters:n (int (1,3,4,7), desired Gauss points per element.) –
zero_areas(triangle_raw, area_null)[source]

Looks for “zero-areas”, areas that are really small, almost zero. It appends them to area_null list.

Parameters:
  • s (class, surface where we whan to look for zero areas.) –
  • triangle_raw (list, triangles of the surface.) –
  • area_null (list, contains the zero areas.) –
Returns:

area_null

Return type:

list, indices of the triangles with zero-areas.

class pygbe.classes.Timing[source]

Bases: object

Timing class. Timing information for different parts of the code.

time_an

float, time spent in compute the near singular integrals.

time_P2P

float, time spent in compute the P2P part of the treecode.

time_P2M

float, time spent in compute the P2M part of the treecode.

time_M2M

float, time spent in compute the M2M part of the treecode.

time_M2P

float, time spent in compute the M2P part of the treecode.

time_trans

float, time spent in transfer data to and from the GPU.

time_sort

float, time spent in sorting data to send to the GPU.

time_mass

float, time spent in compute the mass of the sources in treecode.

AI_int

int, counter of the amount of near singular integrals solved.

pygbe.gmres module

Generalized Minimum Residual Method (GMRES).

GMRES iteratively refines the initial solution guess to the system Ax=b.

This implementation was based mainly on the gmres_mgs from PyAMG, where modified Gram-Schmidt is used to orthogonalize the Krylov Space and Givens Rotations are used to provide the residual norm each iteration.

Guidance code:

pygbe.gmres.apply_givens(Q, v, k)[source]

Apply the first k Givens rotations in Q to the vector v.

Parameters:
  • Q (list, list of consecutive 2x2 Givens rotations) –
  • v (array, vector to apply the rotations to) –
  • k (int, number of rotations to apply) –
Returns:

v

Return type:

array, that is changed in place.

pygbe.gmres.gmres_mgs(surf_array, field_array, X, b, param, ind0, timing, kernel)[source]

GMRES solver.

Parameters:
  • surf_array (array, contains the surface classes of each region on the) – surface.
  • field_array (array, contains the Field classes of each region on the surface.) –
  • X (array, initial guess.) –
  • b (array, right hand side.) –
  • param (class, parameters related to the surface.) –
  • ind0 (class, it contains the indices related to the treecode) – computation.
  • timing (class, it contains timing information for different parts of) – the code.
  • kernel (pycuda source module.) –
Returns:

  • X (array, an updated guess to the solution.)
  • iteration (int, number of outer iterations for convergence)

References

[1]Yousef Saad, “Iterative Methods for Sparse Linear Systems, Second Edition”, SIAM, pp. 151-172, pp. 272-275, 2003 http://www-users.cs.umn.edu/~saad/books.html
[2]
    1. Kelley, http://www4.ncsu.edu/~ctk/matlab_roots.html

pygbe.gpuio module

It contains the function in charge of the data transfer to the GPU.

pygbe.gpuio.dataTransfer(surf_array, field_array, ind, param, kernel)[source]

It manages the data transfer to the GPU.

Parameters:
  • surf_array (array, contains the surface classes of each region on the) – surface.
  • field_array (array, contains the Field classes of each region on the surface.) –
  • ind (class, it contains the indices related to the treecode) – computation.
  • param (class, parameters related to the surface.) –
  • kernel (pycuda source module.) –

pygbe.main module

This is the main function of the program. We use a boundary element method (BEM) to perform molecular electrostatics calculations with a continuum approach. It calculates solvation energies for proteins modeled with any number of dielectric regions.

class pygbe.main.Logger(filename='Default.log')[source]

Bases: object

Allow writing both to STDOUT on screen and sending text to file in conjunction with the command sys.stdout = Logger(“desired_log_file.txt”)

flush()[source]

Required for Python 3

write(message)[source]
pygbe.main.check_file_exists(filename)[source]

Try to open the file filename and return True if it’s valid

pygbe.main.check_for_nvcc()[source]

Check system PATH for nvcc, exit if not found

pygbe.main.find_config_files(cliargs)[source]

Check that .config and .param files exist and can be opened. If either file isn’t found, PyGBe exits (and should print which file was not found). Otherwise return the path to the config and param files

Parameters:cliargs (parser) – parser containing cli arguments passed to PyGBe
Returns:
  • cliargs.config (string) – path to config file
  • cliargs.param (string) – path to param file
pygbe.main.main(argv=['/home/travis/miniconda/envs/pygbe/bin/sphinx-build', '-b', 'html', '-d', '_build/doctrees', '.', '_build/html'], log_output=True, return_output_fname=False, return_results_dict=False, field=None)[source]

Run a PyGBe problem, write outputs to STDOUT and to log file in problem directory

Parameters:
  • log_output (Bool, default True.) – If False, output is written only to STDOUT and not to a log file.
  • return_output_fname (Bool, default False.) – If True, function main() returns the name of the output log file. This is used for the regression tests.
  • return_results_dict (Bool, default False.) – If True, function main() returns the results of the run packed in a dictionary. Used in testing and programmatic use of PyGBe
  • field (Dictionary, defaults to None.) – If passed, this dictionary will supercede any config file found, useful in programmatically stepping through slight changes in a problem
Returns:

output_fname – The name of the log file containing problem output

Return type:

str, if kwarg is True.

pygbe.main.read_inputs(args)[source]

Parse command-line arguments to determine which config and param files to run Assumes that in the absence of specific command line arguments that pygbe problem folder resembles the following structure

lys/ - lys.param - lys.config - built_parse.pqr - geometry/Lys1.face - geometry/Lys1.vert - output/

pygbe.main.resolve_relative_config_file(config_file, full_path)[source]

Does its level-headed best to find the config files specified by the user

Parameters:
  • config_file (str) – the given path to a .param or .config file from the command line
  • full_path (str) – the full path to the problem folder

pygbe.matrixfree module

Matrix free formulation of the matrix vector product in the GMRES.

pygbe.matrixfree.calc_aux(field, surface)[source]

Helper function to calculate aux

Parameters:
  • field (Field object, current field being evaluated) –
  • surface (Surface object, current surface being evaluated) –
Returns:

aux

Return type:

numpy array

pygbe.matrixfree.calculate_solvation_energy(surf_array, field_array, param, kernel)[source]

It calculates the solvation energy.

Parameters:
  • surf_array (array, contains the surface classes of each region on the) – surface.
  • field_array (array, contains the Field classes of each region on the surface.) –
  • param (class, parameters related to the surface.) –
  • kernel (pycuda source module.) –
Returns:

E_solv

Return type:

float, solvation energy.

pygbe.matrixfree.calculate_surface_energy(surf_array, field_array, param, kernel)[source]

It calculates the surface energy

Parameters:
  • surf_array (array, contains the surface classes of each region on the) – surface.
  • field_array (array, contains the Field classes of each region on the surface.) –
  • param (class, parameters related to the surface.) –
  • kernel (pycuda source module.) –
Returns:

E_surf

Return type:

float, surface energy.

pygbe.matrixfree.coulomb_energy(f, param)[source]

It calculates the Coulomb energy .

Parameters:
  • f (class, region in the field array where we desire to calculate the) – coloumb energy.
  • param (class, parameters related to the surface.) –
Returns:

Ecoul

Return type:

float, coloumb energy.

pygbe.matrixfree.generateRHS(field_array, surf_array, param, kernel, timing, ind0)[source]

It generate the right hand side (RHS) for the GMRES.

Parameters:
  • field_array (array, contains the Field classes of each region on the surface.) –
  • surf_array (array, contains the surface classes of each region on the) – surface.
  • param (class, parameters related to the surface.) –
  • kernel (pycuda source module.) –
  • timing (class, it contains timing information for different parts of) – the code.
  • ind0 (class, it contains the indices related to the treecode computation.) –
Returns:

F

Return type:

array, RHS.

pygbe.matrixfree.generateRHS_gpu(field_array, surf_array, param, kernel, timing, ind0)[source]

It generate the right hand side (RHS) for the GMRES suitable for the GPU.

Parameters:
  • field_array (array, contains the Field classes of each region on the surface.) –
  • surf_array (array, contains the surface classes of each region on the) – surface.
  • param (class, parameters related to the surface.) –
  • kernel (pycuda source module.) –
  • timing (class, it contains timing information for different parts of) – the code.
  • ind0 (class, it contains the indices related to the treecode) – computation.
Returns:

F

Return type:

array, RHS suitable for the GPU.

pygbe.matrixfree.gmres_dot(X, surf_array, field_array, ind0, param, timing, kernel)[source]

It computes the matrix-vector product in the GMRES.

Parameters:
  • X (array, initial vector guess.) –
  • surf_array (array, contains the surface classes of each region on the) – surface.
  • field_array (array, contains the Field classes of each region on the) – surface.
  • ind0 (class, it contains the indices related to the treecode) – computation.
  • param (class, parameters related to the surface.) –
  • timing (class, it contains timing information for different parts of) – the code.
  • kernel (pycuda source module.) –
Returns:

MV

Return type:

array, resulting matrix-vector multiplication.

pygbe.matrixfree.locate_s_in_RHS(surf_index, surf_array)[source]

Find starting index for current surface on RHS

This is assembling the RHS of the block matrix. Needs to go through all the previous surfaces to find out where on the RHS vector it belongs. If any previous surfaces were dirichlet, neumann or asc, then they take up half the number of spots in the RHS, so we act accordingly.

Parameters:
  • surf_index (int, index of surface in question) –
  • surf_array (list, list of all surfaces in problem) –
Returns:

s_start

Return type:

int, index to insert values on RHS

pygbe.matrixfree.nonselfExterior(surf, src, tar, LorY, param, ind0, timing, kernel)[source]

Non-self exterior operator:

The evaluation point resides in a surface that is outside the region enclosed by the surface with the strength of the single and double layer potentials. If both surfaces share a common external region, we add the result to the external equation. If they don’t share a common external region we add the result to internal equation.

Parameters:
  • surf (array, contains all the classes of the surface.) –
  • src (int, position in the array of the source surface (surface that) – contains the gauss points).
  • tar (int, position in the array of the target surface (surface that) – contains the collocation points).
  • LorY (int, Laplace (1) or Yukawa (2)) –
  • param (class, parameters related to the surface.) –
  • ind0 (class, it contains the indices related to the treecode computation.) –
  • timing (class, it contains timing information for different parts of the) – code.
  • kernel (pycuda source module.) –
Returns:

v – non-self exterior interaction.

Return type:

array, result of the matrix-vector product corresponding to the

pygbe.matrixfree.nonselfInterior(surf, src, tar, LorY, param, ind0, timing, kernel)[source]

Non-self interior operator:

The evaluation point resides in a surface that is inside the region enclosed by the surface with the strength of the single and double layer potentials, and the result has to be added to the exterior equation.

Parameters:
  • surf (array, contains all the classes of the surface.) –
  • src (int, position in the array of the source surface (surface that) – contains the gauss points).
  • tar (int, position in the array of the target surface (surface that) – contains the collocation points).
  • LorY (int, Laplace (1) or Yukawa (2)) –
  • param (class, parameters related to the surface.) –
  • ind0 (class, it contains the indices related to the treecode computation.) –
  • timing (class, it contains timing information for different parts of the) – code.
  • kernel (pycuda source module.) –
Returns:

v – non-self interior interaction.

Return type:

array, result of the matrix-vector product corresponding to the

pygbe.matrixfree.selfASC(surf, src, tar, LorY, param, ind0, timing, kernel)[source]

Self interaction for the aparent surface charge (ASC) formulation.

Parameters:
  • surf (array, contains all the classes of the surface.) –
  • src (int, position in the array of the source surface (surface that) – contains the gauss points).
  • tar (int, position in the array of the target surface (surface that) – contains the collocation points).
  • LorY (int, Laplace (1) or Yukawa (2)) –
  • param (class, parameters related to the surface.) –
  • ind0 (class, it contains the indices related to the treecode computation.) –
  • timing (class, it contains timing information for different parts of the) – code.
  • kernel (pycuda source module.) –
Returns:

v – self interaction in the ASC formulation.

Return type:

array, result of the matrix-vector product corresponding to the

pygbe.matrixfree.selfExterior(surf, s, LorY, param, ind0, timing, kernel)[source]

Self surface exterior operator:

The evaluation point and the strengths of both the single and double layer potential are on the same surface. When we take the limit, the evaluation point approaches the surface from the outside, then the result is added to the exterior equation.

Parameters:
  • surf (array, contains all the classes of the surface.) –
  • s (int, position of the surface in the surface array.) –
  • LorY (int, Laplace (1) or Yukawa (2)) –
  • param (class, parameters related to the surface.) –
  • ind0 (class, it contains the indices related to the treecode computation.) –
  • timing (class, it contains timing information for different parts of the) – code.
  • kernel (pycuda source module.) –
Returns:

  • v (array, result of the matrix-vector product corresponding to the self) – exterior interaction.
  • K_lyr (array, self exterior double layer potential.)
  • V_lyr (array, self exterior single layer potential.)

pygbe.matrixfree.selfInterior(surf, s, LorY, param, ind0, timing, kernel)[source]

Self surface interior operator:

The evaluation point and strengths of the single and double layer potential are on the same surface. When we take the limit, the evaluation point approaches the surface from the inside, then the result is added to the interior equation.

Parameters:
  • surf (array, contains all the classes of the surface.) –
  • s (int, position of the surface in the surface array.) –
  • LorY (int, Laplace (1) or Yukawa (2)) –
  • param (class, parameters related to the surface.) –
  • ind0 (class, it contains the indices related to the treecode computation.) –
  • timing (class, it contains timing information for different parts of the) – code.
  • kernel (pycuda source module.) –
Returns:

v – interior interaction.

Return type:

array, result of the matrix-vector product corresponding to the self

pygbe.output module

Prints output with the main information.

pygbe.output.print_summary(surf_array, field_array, param, results_dict)[source]

Prints a summary with the main information of the run.

Parameters:
  • surf_array (array, contains the surface classes of each region on the) – surface.
  • field_array (array, contains the Field classes of each region on the surface.) –
  • param (class, parameters related to the surface.) –

pygbe.projection module

It contains the functions to calculate the different potentials: -The single and double layer potential. -The adjoint double layer potential. -The reaction potential.

pygbe.projection.get_phir(XK, XV, surface, xq, Cells, par_reac, ind_reac)[source]

It computes the reaction potential. To compute this potential we need more terms in the Taylor expansion, that is the reason why we need fine parameters (par_reac class) and a different array of indices (ind_reac) than ind0.

Parameters:
  • XK (array, input for the double layer potential.) –
  • XV (array, input for the single layer potential.) –
  • surface (class, surface where we are computing the reaction potential.) –
  • xq (array, it contains the position of the charges.) –
  • Cells (array, it contains the tree cells.) –
  • par_reac (class, fine parameters related to the surface.) –
  • ind_reac (array, it contains the indices related to the treecode) – computation.
Returns:

  • phi_reac (array, reaction potential.)
  • AI_int (int, counter of the amount of near singular integrals solved.)

pygbe.projection.get_phir_gpu(XK, XV, surface, field, par_reac, kernel)[source]

It computes the reaction potential on the GPU and it brings the data to the cpu.

Parameters:
  • XK (array, input for the double layer potential.) –
  • XV (array, input for the single layer potential.) –
  • surface (class, surface where we are computing the reaction potential.) –
  • field (class, information about the different regions in the molecule.) –
  • par_reac (class, fine parameters related to the surface.) –
Returns:

  • phir_cpu (array, reaction potential brought from the GPU to the cpu.)
  • AI_int (int, counter of the amount of near singular integrals solved.)

pygbe.projection.project(XK, XV, LorY, surfSrc, surfTar, K_diag, V_diag, IorE, self, param, ind0, timing, kernel)[source]

It computes the single and double layer potentials.

Parameters:
  • XK (array, input for the double layer potential.) –
  • XV (array, input for the single layer potential.) –
  • LorY (int, Laplace (1) or Yukawa (2)) –
  • surfSrc (class, source surface, the one that contains the gauss points.) –
  • surfTar (class, target surface, the one that contains the collocation) – points.
  • K_diag (array, diagonal elements of the double layer integral operator.) –
  • V_diag (array, diagonal elements of the single layer integral operator.) –
  • IorE (int, internal (1) or external (2)) –
  • self (int, position in the surface array of the source surface.) –
  • param (class, parameters related to the surface.) –
  • ind0 (array, it contains the indices related to the treecode computation.) –
  • timing (class, it contains timing information for different parts of the) – code.
  • kernel (pycuda source module.) –
Returns:

  • K_lyr (array, double layer potential.)
  • V_lyr (array, single layer potential.)

pygbe.projection.project_Kt(XKt, LorY, surfSrc, surfTar, Kt_diag, self, param, ind0, timing, kernel)[source]

It computes the adjoint double layer potential.

Parameters:
  • XKt (array, input for the adjoint double layer potential.) –
  • LorY (int, Laplace (1) or Yukawa (2)) –
  • surfSrc (class, source surface, the one that contains the gauss points.) –
  • surfTar (class, target surface, the one that contains the collocation points.) –
  • Kt_diag (array, diagonal elements of the adjoint double layer integral) – operator.
  • self (int, position in the surface array of the source surface.) –
  • param (class, parameters related to the surface.) –
  • ind0 (array, it contains the indices related to the treecode computation.) –
  • timing (class, it contains timing information for different parts of the) – code.
  • kernel (pycuda source module.) –
Returns:

Kt_lyr

Return type:

array, adjoint double layer potential.

pygbe.quadrature module

It contains the functions to compute the fine Gaussian quadrature and the wights and gauss points for the regular Gauss quadrature.

pygbe.quadrature.getWeights(K)[source]

It gets the weights of the Gauss points.

Parameters:K (int, number of Gauss points per element. (1, 3, 4, and 7 are supported)) –
Returns:w
Return type:K-size array, weights of the Gauss points.
pygbe.quadrature.quadratureRule_fine(K)[source]

Fine quadrature rule, to solve the near singular integrals.

Parameters:K (int (1, 7, 13, 17, 19, 25, 37, 48, 52, 61, 79), number of Gauss points) – per element.
Returns:
  • X (array, position of the gauss quadrature points.)
  • W (array, gauss quadrature weights.)

Module contents