ergo
|
#include <stdlib.h>
#include <math.h>
#include <stdio.h>
#include <errno.h>
#include <memory.h>
#include <time.h>
#include <stdarg.h>
#include <assert.h>
#include <vector>
#include "integrals_1el.h"
#include "integrals_1el_potential.h"
#include "memorymanag.h"
#include "pi.h"
#include "output.h"
#include "utilities.h"
#include "boysfunction.h"
#include "integral_info.h"
#include "integrals_general.h"
#include "box_system.h"
#include "multipole.h"
#include "integrals_2el_single.h"
#include "integrals_1el_single.h"
#include "integrals_hermite.h"
#include "matrix_norm.h"
#include "mm_limit_table.h"
Classes | |
struct | atom_box_struct |
struct | DistributionSpecStructWithPairIdx |
struct | group_struct |
Functions | |
static ergo_real | get_distance_3d (const ergo_real *x, const ergo_real *y) |
static int | create_nuclei_mm_tree (int nAtoms, const Atom *atomList, ergo_real boxSize, atom_box_struct **return_boxList, int *return_numberOfLevels, Atom **return_atomListReordered) |
static int | do_interaction_recursive (const IntegralInfo &integralInfo, ergo_real *V_list, int noOfBasisFuncIndexPairs, const DistributionSpecStructWithPairIdx *list, int nDistrs, const multipole_struct_small *multipoleList, const ergo_real *maxMomentVectorNormForDistrsList, int maxNoOfMomentsForDistrs, int maxDegreeForDistrs, ergo_real distrExtent, const Atom *atomListReordered, ergo_real threshold, const atom_box_struct *boxList, MMInteractor &interactor, int boxIndex, int currLevel, int numberOfLevels) |
Take care of interaction between list of distrs and box. | |
static int | get_list_of_distrs_for_V (const BasisInfoStruct &basisInfo, const basis_func_index_pair_struct_1el *basisFuncIndexPairList, int noOfBasisFuncIndexPairs, ergo_real threshold, ergo_real maxCharge, DistributionSpecStructWithPairIdx *resultList, int maxCountResult) |
static int | compare_distrs (const void *p1, const void *p2) |
static int | sort_distr_list (DistributionSpecStructWithPairIdx *list, int n) |
static ergo_real | getSafeMaxDistance (const BasisInfoStruct &basisInfo, const Molecule &molecule) |
int | compute_V_linear (const BasisInfoStruct &basisInfo, const IntegralInfo &integralInfo, const Molecule &molecule, ergo_real threshold, ergo_real boxSize, const basis_func_index_pair_struct_1el *basisFuncIndexPairList, ergo_real *V_list, int noOfBasisFuncIndexPairs) |
static ergo_real | simplePrimVintegralSingle (DistributionSpecStruct *prim, const Atom *atom, const IntegralInfo *integralInfo) |
static ergo_real | simplePrimVintegral_list (DistributionSpecStruct *list, int nPrims, const Atom *atom, ergo_real threshold, const IntegralInfo *integralInfo) |
int | compute_V_matrix_full (const BasisInfoStruct &basisInfo, const IntegralInfo &integralInfo, int nAtoms, const Atom *atomList, ergo_real threshold, ergo_real *result) |
Variables | |
const ergo_real | tolernance_dist = 1e-10 |
const ergo_real | tolernance_exponent = 1e-11 |
static int compare_distrs | ( | const void * | p1, |
const void * | p2 | ||
) | [static] |
References tolernance_dist, tolernance_exponent, DistributionSpecStructWithPairIdx::distr, DistributionSpecStruct_::centerCoords, and DistributionSpecStruct_::exponent.
Referenced by sort_distr_list(), and compute_V_linear().
int compute_V_linear | ( | const BasisInfoStruct & | basisInfo, |
const IntegralInfo & | integralInfo, | ||
const Molecule & | molecule, | ||
ergo_real | threshold, | ||
ergo_real | boxSize, | ||
const basis_func_index_pair_struct_1el * | basisFuncIndexPairList, | ||
ergo_real * | V_list, | ||
int | noOfBasisFuncIndexPairs | ||
) |
References do_output(), LOG_CAT_INFO, LOG_AREA_INTEGRALS, Molecule::noOfAtoms, Molecule::atoms, Atom::charge, get_list_of_distrs_for_V(), LOG_CAT_ERROR, sort_distr_list(), compare_distrs(), init_multipole_code(), getSafeMaxDistance(), mm_limits_init(), create_nuclei_mm_tree(), compute_multipole_moments(), MAX_MULTIPOLE_DEGREE_BASIC, multipole_struct_small::noOfMoments, multipole_struct_small::degree, A, multipole_struct_small::momentList, DistributionSpecStructWithPairIdx::distr, DistributionSpecStruct_::exponent, and do_interaction_recursive().
Referenced by compute_V_sparse().
int compute_V_matrix_full | ( | const BasisInfoStruct & | basisInfo, |
const IntegralInfo & | integralInfo, | ||
int | nAtoms, | ||
const Atom * | atomList, | ||
ergo_real | threshold, | ||
ergo_real * | result | ||
) |
References A, BasisInfoStruct::noOfBasisFuncs, BasisInfoStruct::basisFuncList, BasisFuncStruct_::noOfSimplePrimitives, BasisFuncStruct_::simplePrimitiveIndex, BasisInfoStruct::simplePrimitiveList, get_product_simple_prims(), do_output(), LOG_CAT_ERROR, LOG_AREA_INTEGRALS, and simplePrimVintegral_list().
Referenced by compute_h_core_matrix_full(), savePotential(), and test_V_by_explicit_comparison().
static int create_nuclei_mm_tree | ( | int | nAtoms, |
const Atom * | atomList, | ||
ergo_real | boxSize, | ||
atom_box_struct ** | return_boxList, | ||
int * | return_numberOfLevels, | ||
Atom ** | return_atomListReordered | ||
) | [static] |
References box_item_struct::originalIndex, BoxSystem::create_box_system(), do_output(), LOG_CAT_ERROR, LOG_AREA_INTEGRALS, BoxSystem::noOfLevels, BoxSystem::totNoOfBoxes, BoxSystem::boxList, BoxSystem::levelList, box_level_struct::startIndexInBoxList, box_level_struct::noOfBoxes, atom_box_struct::basicBox, box_struct_basic::firstItemIndex, box_struct_basic::noOfItems, Atom::charge, Atom::coords, atom_box_struct::centerOfChargeCoords, box_struct_basic::centerCoords, multipole_struct_large::centerCoords, multipole_struct_large::degree, MAX_MULTIPOLE_DEGREE, multipole_struct_large::noOfMoments, MAX_NO_OF_MOMENTS_PER_MULTIPOLE, multipole_struct_small::degree, multipole_struct_small::noOfMoments, multipole_struct_small::momentList, MMTranslator::getTranslationMatrix(), A, B, multipole_struct_large::momentList, atom_box_struct::multipole, box_struct_basic::noOfChildBoxes, box_struct_basic::firstChildBoxIndex, and setup_multipole_maxAbsMomentList().
Referenced by compute_V_linear().
static int do_interaction_recursive | ( | const IntegralInfo & | integralInfo, |
ergo_real * | V_list, | ||
int | noOfBasisFuncIndexPairs, | ||
const DistributionSpecStructWithPairIdx * | list, | ||
int | nDistrs, | ||
const multipole_struct_small * | multipoleList, | ||
const ergo_real * | maxMomentVectorNormForDistrsList, | ||
int | maxNoOfMomentsForDistrs, | ||
int | maxDegreeForDistrs, | ||
ergo_real | distrExtent, | ||
const Atom * | atomListReordered, | ||
ergo_real | threshold, | ||
const atom_box_struct * | boxList, | ||
MMInteractor & | interactor, | ||
int | boxIndex, | ||
int | currLevel, | ||
int | numberOfLevels | ||
) | [static] |
Take care of interaction between list of distrs and box.
References atom_box_struct::multipole, distance(), get_distance_3d(), atom_box_struct::basicBox, box_struct_basic::centerCoords, box_struct_basic::width, multipole_struct_large::centerCoords, mm_limits_get_minimum_multipole_degree_needed(), MAX_MULTIPOLE_DEGREE, DistributionSpecStructWithPairIdx::distr, DistributionSpecStruct_::centerCoords, MMInteractor::getInteractionMatrix(), MAX_NO_OF_MOMENTS_PER_MULTIPOLE, A, B, multipole_struct_large::momentList, multipole_struct_small::noOfMoments, DistributionSpecStructWithPairIdx::pairIdx, DistributionSpecStruct_::monomialInts, box_struct_basic::firstItemIndex, box_struct_basic::noOfItems, IntegralInfo::monomial_info, monomial_info_struct::no_of_monomials_list, DistributionSpecStruct_::exponent, Atom::coords, pi, get_related_integrals_hermite(), IntegralInfo::hermite_conversion_info, hermite_conversion_info_struct::multiply_by_hermite_conversion_matrix_from_right(), monomial_info_struct::monomial_index_list, Atom::charge, DistributionSpecStruct_::coeff, box_struct_basic::noOfChildBoxes, and box_struct_basic::firstChildBoxIndex.
Referenced by compute_V_linear().
Referenced by do_interaction_recursive().
static int get_list_of_distrs_for_V | ( | const BasisInfoStruct & | basisInfo, |
const basis_func_index_pair_struct_1el * | basisFuncIndexPairList, | ||
int | noOfBasisFuncIndexPairs, | ||
ergo_real | threshold, | ||
ergo_real | maxCharge, | ||
DistributionSpecStructWithPairIdx * | resultList, | ||
int | maxCountResult | ||
) | [static] |
References basis_func_index_pair_struct_1el::index_1, basis_func_index_pair_struct_1el::index_2, POLY_PRODUCT_MAX_DISTRS, get_product_simple_primitives(), do_output(), LOG_CAT_ERROR, LOG_AREA_INTEGRALS, pi, DistributionSpecStruct_::coeff, DistributionSpecStruct_::exponent, DistributionSpecStructWithPairIdx::distr, and DistributionSpecStructWithPairIdx::pairIdx.
Referenced by compute_V_linear().
static ergo_real getSafeMaxDistance | ( | const BasisInfoStruct & | basisInfo, |
const Molecule & | molecule | ||
) | [static] |
static ergo_real simplePrimVintegral_list | ( | DistributionSpecStruct * | list, |
int | nPrims, | ||
const Atom * | atom, | ||
ergo_real | threshold, | ||
const IntegralInfo * | integralInfo | ||
) | [static] |
References simplePrimVintegralSingle().
Referenced by compute_V_matrix_full().
static ergo_real simplePrimVintegralSingle | ( | DistributionSpecStruct * | prim, |
const Atom * | atom, | ||
const IntegralInfo * | integralInfo | ||
) | [static] |
References do_1e_repulsion_integral_using_symb_info(), Atom::charge, and Atom::coords.
Referenced by simplePrimVintegral_list().
static int sort_distr_list | ( | DistributionSpecStructWithPairIdx * | list, |
int | n | ||
) | [static] |
References compare_distrs().
Referenced by compute_V_linear().
const ergo_real tolernance_dist = 1e-10 |
Referenced by compare_distrs().
const ergo_real tolernance_exponent = 1e-11 |
Referenced by compare_distrs().