{ "info": { "author": "Radovan Bast", "author_email": "radovan.bast@uit.no", "bugtrack_url": null, "classifiers": [ "Development Status :: 3 - Alpha", "Intended Audience :: Science/Research", "Programming Language :: Python :: 3.6" ], "description": ".. image:: https://travis-ci.org/dftlibs/numgrid.svg?branch=master\n :target: https://travis-ci.org/dftlibs/numgrid/builds\n.. image:: https://coveralls.io/repos/dftlibs/numgrid/badge.png?branch=master\n :target: https://coveralls.io/r/dftlibs/numgrid?branch=master\n.. image:: https://img.shields.io/badge/license-%20MPL--v2.0-blue.svg\n :target: LICENSE\n.. image:: https://badge.fury.io/py/numgrid.svg\n :target: https://badge.fury.io/py/numgrid\n.. image:: https://zenodo.org/badge/DOI/10.5281/zenodo.1470277.svg\n :target: https://doi.org/10.5281/zenodo.1470277\n\n- `Changelog `__\n- `Build and test\n history `__\n- `Code coverage `__\n- Licensed under `MPL v2.0 `__ (except John\n Burkardt\u2019s Lebedev code which is redistributed under LGPL v3.0)\n- Built with `Autocmake `__\n\n\nAbout\n=====\n\nnumgrid is a library that produces numerical integration grid for\nmolecules based on atom coordinates, atom types, and basis set\ninformation. This library can be built with C, Fortran, and Python bindings.\n\n\nWho are the people behind this code?\n====================================\n\nAuthors\n-------\n\n- Radovan Bast\n\n\nContributors\n------------\n\n- Roberto Di Remigio (OS X testing, streamlined Travis testing, better\n C++, error handling)\n\nFor a list of all the contributions see\nhttps://github.com/dftlibs/numgrid/contributors.\n\n\nAcknowledgements\n----------------\n\n- Simon Neville (reporting issues)\n- Jaime Axel Rosal Sandberg (reporting issues)\n\nThis tool uses SPHERE_LEBEDEV_RULE, a C library written by John Burkardt which\ncomputes a Lebedev quadrature rule over the surface of the unit sphere in 3D,\nsee also:\nhttp://people.sc.fsu.edu/~jburkardt/c_src/sphere_lebedev_rule/sphere_lebedev_rule.html\n\nThis library uses and acknowledges the\nMolSSI BSE (https://molssi-bse.github.io/basis_set_exchange/),\nwhich is a rewrite of the Basis Set Exchange\n(https://bse.pnl.gov/bse/portal) and is a collaboration between the Molecular\nSciences Software Institute (http://www.molssi.org) and the Environmental\nMolecular Sciences Laboratory (https://www.emsl.pnl.gov).\n\n\nCitation\n--------\n\nIf you use this tool in a program or publication, please acknowledge its\nauthor(s). The metadata necessary for citing this tool can be found in the\n`CITATION.cff `__ file. For more information CITATION.cff files, see\nhttps://citation-file-format.github.io.\n\n\nWould you like to contribute?\n-----------------------------\n\nYes please! Please follow `this excellent\nguide `__. We do not require any\nformal copyright assignment or contributor license agreement. Any\ncontributions intentionally sent upstream are presumed to be offered\nunder terms of the Mozilla Public License Version 2.0.\n\n\nRequirements\n============\n\n- CMake\n- C and C++ compiler\n- Fortran compiler (to build the optional Fortran interface)\n- `CFFI `__ (to access the optional\n Python interface)\n- `pytest `__ (to test the optional Python\n interface)\n\n\nInstallation\n============\n\nInstalling via pip\n------------------\n\n::\n\n pip install numgrid\n\n\nBuilding and testing from sources\n---------------------------------\n\nFetch the code::\n\n git clone https://github.com/dftlibs/numgrid.git\n\nInstall Python dependencies (optional)::\n\n pipenv install\n\nBuild the code::\n\n cd numgrid\n ./setup --fc=gfortran --cc=gcc --cxx=g++\n cd build\n make\n make test\n\nThe Python interface is automatically tested by Travis CI:\nhttps://github.com/dftlibs/numgrid/blob/master/.travis.yml\n\n\nAPI\n===\n\nThe library provides a context-aware C interface. In addition it also\nprovides a Fortran and Python interfaces as thin layers on top of the C\ninterface::\n\n Python: numgrid/__init__.py\n \\\n \\ Fortran: numgrid/numgrid.f90\n \\ /\n C interface: numgrid/numgrid.h\n |\n implementation\n\n\nUnits\n-----\n\nCoordinates are in bohr.\n\n\nOverview\n--------\n\nGrid computation is done per atom/basis type and proceeds in five steps:\n\n- Create atom\n- Get number of points (depends on basis set range)\n- Allocate memory to hold the grid\n- Compute grid on this atom in a molecular environment\n- Free atom and its memory\n\nThe Python interface takes care of the allocation and deallocation part\nbut the essential point is that memory management is happening on the\nclient side.\n\nIf you have many atom centers that have the same atom type and same\nbasis set, it will make sense to create only one atom object and then\nreuse this object to compute the grid on all atoms with the same basis\ntype.\n\nIt is no problem to create several atom objects at the same time.\n\n\nPython example\n--------------\n\nThe Python interface is generated using\n`CFFI `__.\n\nAs an example let us generate a grid for the water molecule:\n\n.. code:: python\n\n import numgrid\n\n radial_precision = 1.0e-12\n min_num_angular_points = 86\n max_num_angular_points = 302\n\n num_centers = 3\n proton_charges = [8, 1, 1]\n\n x_coordinates_bohr = [0.0, 1.43, -1.43]\n y_coordinates_bohr = [0.0, 0.0, 0.0]\n z_coordinates_bohr = [0.0, 1.1, 1.1]\n\n # cc-pVDZ basis\n alpha_max = [11720.0, 13.01, 13.01] # O, H, H\n max_l_quantum_numbers = [2, 1, 1] # O, H, H\n alpha_min = [[0.3023, 0.2753, 1.185], # O\n [0.122, 0.727], # H\n [0.122, 0.727]] # H\n\n for center_index in range(num_centers):\n context = numgrid.new_atom_grid(radial_precision,\n min_num_angular_points,\n max_num_angular_points,\n proton_charges[center_index],\n alpha_max[center_index],\n max_l_quantum_numbers[center_index],\n alpha_min[center_index])\n\n num_points = numgrid.get_num_grid_points(context)\n\n # generate an atomic grid in the molecular environment\n x, y, z, w = numgrid.get_grid(context,\n num_centers,\n center_index,\n x_coordinates_bohr,\n y_coordinates_bohr,\n z_coordinates_bohr,\n proton_charges)\n\n num_radial_points = numgrid.get_num_radial_grid_points(context)\n\n # generate an isolated radial grid\n r, w = numgrid.get_radial_grid(context)\n\n numgrid.free_atom_grid(context)\n\n\n # generate an isolated angular grid\n x, y, z, w = numgrid.get_angular_grid(num_angular_grid_points=14)\n\n\nAvoiding explicit exponent ranges\n---------------------------------\n\nUsing the Python interface you can choose to not provide\nexplicit exponent ranges and instead specify the basis\nset which is then fetched directly from\nhttps://github.com/MolSSI-BSE/basis_set_exchange\nusing the wonderful\n`MolSSI BSE `__:\n\n.. code:: python\n\n context = numgrid.new_atom_grid_bse(radial_precision=1.0e-12,\n min_num_angular_points=86,\n max_num_angular_points=302,\n proton_charge=8,\n basis_set='cc-pVDZ')\n\n\nSaving grid in Numpy format\n---------------------------\n\nThe current API makes is relatively easy to export the computed grid in Numpy format.\n\nIn this example we save the radial grid positions and weights to two separate files\nin Numpy format:\n\n.. code:: python\n\n import numgrid\n import numpy as np\n\n # we assume the context is created\n # ...\n r, w = numgrid.get_radial_grid(context)\n\n np.save('radial_grid_r.npy', r)\n np.save('radial_grid_w.npy', w)\n\n\nC API\n-----\n\nTo see a real example, have a look at the `C++ test\ncase `__.\n\n\nCreating a new atom grid\n~~~~~~~~~~~~~~~~~~~~~~~~\n\n.. code:: c\n\n context_t *numgrid_new_atom_grid(const double radial_precision,\n const int min_num_angular_points,\n const int max_num_angular_points,\n const int proton_charge,\n const double alpha_max,\n const int max_l_quantum_number,\n const double alpha_min[]);\n\nThe smaller the ``radial_precision``, the better grid.\n\nFor ``min_num_angular_points`` and ``max_num_angular_points``, see\n\u201cAngular grid\u201d below.\n\n``alpha_max`` is the steepest basis set exponent.\n\n``alpha_min`` is an array of the size ``max_l_quantum_number`` + 1 and\nholds the smallest exponents for each angular momentum. If an angular\nmomentum set is missing \u201cin the middle\u201d, provide 0.0. In other words,\nimagine that you have a basis set which only contains *s* and *d*\nfunctions and no *p* functions and let us assume that the most diffuse\n*s* function has the exponent 0.1 and the most diffuse *d* function has\nthe exponent 0.2, then ``alpha_min`` would be an array of three numbers\nholding {0.1, 0.0, 0.2}.\n\n\nGet number of grid points on current atom\n~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n\nThe following two functions are probably self-explaining. We need to\nprovide the context which refers to a specific atom object.\n\n.. code:: c\n\n int numgrid_get_num_grid_points(const context_t *context);\n\n int numgrid_get_num_radial_grid_points(const context_t *context);\n\n\nGet grid on current atom, scaled by Becke partitioning\n~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n\nWe assume that ``grid_x_bohr``, ``grid_y_bohr``, ``grid_z_bohr``, and\n``grid_w`` are allocated by the caller and have the length that equals\nthe number of grid points.\n\n``x_coordinates_bohr``, ``y_coordinates_bohr``, ``z_coordinates_bohr``,\nand ``proton_charges`` refer to the molecular environment and have the\nsize ``num_centers``.\n\nUsing ``center_index`` we tell the code which of the atom centers is the\none we have computed the grid for.\n\n.. code:: c\n\n void numgrid_get_grid(const context_t *context,\n const int num_centers,\n const int center_index,\n const double x_coordinates_bohr[],\n const double y_coordinates_bohr[],\n const double z_coordinates_bohr[],\n const int proton_charges[],\n double grid_x_bohr[],\n double grid_y_bohr[],\n double grid_z_bohr[],\n double grid_w[]);\n\n\nGet radial grid on current atom\n~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n\nWe assume that ``radial_grid_r_bohr`` and ``radial_grid_w`` are\nallocated by the caller and have both the length that equals the number\nof radial grid points.\n\n.. code:: c\n\n void numgrid_get_radial_grid(const context_t *context,\n double radial_grid_r_bohr[],\n double radial_grid_w[]);\n\n\nGet angular grid\n~~~~~~~~~~~~~~~~\n\nThis does not refer to any specific atom and does not require any\ncontext.\n\n``num_angular_grid_points`` has to be one of the many supported Lebedev\ngrids (see table on the bottom of this page) and the code will assume\nthat the grid arrays are allocated by the caller and have at least the\nsize ``num_angular_grid_points``.\n\n.. code:: c\n\n void numgrid_get_angular_grid(const int num_angular_grid_points,\n double angular_grid_x_bohr[],\n double angular_grid_y_bohr[],\n double angular_grid_z_bohr[],\n double angular_grid_w[]);\n\n\nDestroy the atom and deallocate all data\n~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~\n\n.. code:: c\n\n void numgrid_free_atom_grid(context_t *context);\n\n\nFortran API\n-----------\n\nClosely follows the C API. To see a real example, have a look at the\n`Fortran test case `__.\n\n\nParallelization\n===============\n\nThe design decision was to not parallelize the library but rather\nparallelize over the atom/basis types by the caller. This simplifies\nmodularity and code reuse.\n\n\nSpace partitioning\n==================\n\nThe molecular integration grid is generated from atom-centered grids by\nscaling the grid weights according to the Becke partitioning scheme,\n`JCP 88, 2547 (1988) `__. The\ndefault Becke hardness is 3.\n\n\nRadial grid\n===========\n\nThe radial grid is generated according to Lindh, Malmqvist, and\nGagliardi, `TCA 106, 178\n(2001) `__.\n\nThe motivation for this choice is the nice feature of the above scheme\nthat the range of the radial grid is basis set dependent. The precision\ncan be tuned with one single radial precision parameter. The smaller the\nradial precision, the better quality grid you obtain.\n\nThe basis set (more precisely the Gaussian primitives/exponents) are\nused to generate the atomic radial grid range. This means that a more\ndiffuse basis set generates a more diffuse radial grid.\n\nIf you need a grid but you do not have a basis set or choose not to use\na specific one, then you can feed the library with a fantasy basis set\nconsisting of just two primitives. You can then adjust the range by\nmaking the exponents more steep or more diffuse.\n\n\nAngular grid\n============\n\nThe angular grid is generated according to Lebedev and Laikov [A\nquadrature formula for the sphere of the 131st algebraic order of\naccuracy, Russian Academy of Sciences Doklady Mathematics, Volume 59,\nNumber 3, 1999, pages 477-481].\n\nThe angular grid is pruned. The pruning is a primitive linear\ninterpolation between the minimum number and the maximum number of\nangular points per radial shell. The maximum number is reached at 0.2\ntimes the Bragg radius of the center.\n\nThe higher the values for minimum and maximum number of angular points,\nthe better.\n\nFor the minimum and maximum number of angular points the code will use\nthe following table and select the closest number with at least the\ndesired precision::\n\n {6, 14, 26, 38, 50, 74, 86, 110, 146,\n 170, 194, 230, 266, 302, 350, 434, 590, 770,\n 974, 1202, 1454, 1730, 2030, 2354, 2702, 3074, 3470,\n 3890, 4334, 4802, 5294, 5810}\n\nTaking the same number for the minimum and maximum number of angular\npoints switches off pruning.\n\n\nHow to include numgrid in a CMake project\n=========================================\n\nThere are multiple ways to achieve this. Here is how to include\nthe library using ``FetchContent``:\n\n.. code:: cmake\n\n cmake_minimum_required(VERSION 3.11 FATAL_ERROR)\n\n project(example LANGUAGES CXX)\n\n include(FetchContent)\n\n FetchContent_Declare(\n numgrid\n GIT_REPOSITORY https://github.com/dftlibs/numgrid.git\n GIT_TAG e14bf969d68e7847f5e40f36816f61f245211a9b\n )\n\n FetchContent_GetProperties(numgrid)\n\n if(NOT numgrid_POPULATED)\n FetchContent_Populate(numgrid)\n add_subdirectory(\n ${numgrid_SOURCE_DIR}\n ${numgrid_BINARY_DIR}\n )\n endif()\n\n add_executable(example \"\")\n\n target_sources(\n example\n PRIVATE\n main.cpp\n )\n\n target_link_libraries(\n example\n PRIVATE\n numgrid-objects\n )", "description_content_type": "", "docs_url": null, "download_url": "", "downloads": { "last_day": -1, "last_month": -1, "last_week": -1 }, "home_page": "https://github.com/dftlibs/numgrid", "keywords": "", "license": "MPL-v2.0", "maintainer": "", "maintainer_email": "", "name": "numgrid", "package_url": "https://pypi.org/project/numgrid/", "platform": "", "project_url": "https://pypi.org/project/numgrid/", "project_urls": { "Homepage": "https://github.com/dftlibs/numgrid" }, "release_url": "https://pypi.org/project/numgrid/1.1.0/", "requires_dist": null, "requires_python": "", "summary": "Numerical integration grid for molecules.", "version": "1.1.0" }, "last_serial": 5220940, "releases": { "1.0.2": [ { "comment_text": "", "digests": { "md5": "0daf81c5e9273bc5edf1640e56c81c44", "sha256": "b51b8a71317c1cb9c7b64192051a8aa198ba38a8987575f42585b818a4219993" }, "downloads": -1, "filename": "numgrid-1.0.2.tar.gz", "has_sig": false, "md5_digest": "0daf81c5e9273bc5edf1640e56c81c44", "packagetype": "sdist", "python_version": "source", "requires_python": null, "size": 87002, "upload_time": "2018-10-24T09:03:50", "url": "https://files.pythonhosted.org/packages/3d/33/e1e6e49166ef5e52b4045819a5708c05c3e15022b0b134be2d86449b13b1/numgrid-1.0.2.tar.gz" } ], "1.1.0": [ { "comment_text": "", "digests": { "md5": "29ad96312e49791c9523e5da45219284", "sha256": "ce6c359aee96ca3e544599e4662752185be8451a4c8708074e726edc9faad914" }, "downloads": -1, "filename": "numgrid-1.1.0.tar.gz", "has_sig": false, "md5_digest": "29ad96312e49791c9523e5da45219284", "packagetype": "sdist", "python_version": "source", "requires_python": null, "size": 89130, "upload_time": "2019-05-03T08:51:59", "url": "https://files.pythonhosted.org/packages/46/94/2099ab1ed4593c677f2fc0095704390001fc77d9943877074464d45a0d3f/numgrid-1.1.0.tar.gz" } ] }, "urls": [ { "comment_text": "", "digests": { "md5": "29ad96312e49791c9523e5da45219284", "sha256": "ce6c359aee96ca3e544599e4662752185be8451a4c8708074e726edc9faad914" }, "downloads": -1, "filename": "numgrid-1.1.0.tar.gz", "has_sig": false, "md5_digest": "29ad96312e49791c9523e5da45219284", "packagetype": "sdist", "python_version": "source", "requires_python": null, "size": 89130, "upload_time": "2019-05-03T08:51:59", "url": "https://files.pythonhosted.org/packages/46/94/2099ab1ed4593c677f2fc0095704390001fc77d9943877074464d45a0d3f/numgrid-1.1.0.tar.gz" } ] }