Révision | 545330eb4cc38ddca93a31569d467d2d79e8b701 (tree) |
---|---|
l'heure | 2013-08-05 04:12:21 |
Auteur | Mikiya Fujii <mikiya.fujii@gmai...> |
Commiter | Mikiya Fujii |
Mndo::CalcDiatomicTwoElecTwoCore1stDerivatives is refactored. #31814
git-svn-id: https://svn.sourceforge.jp/svnroot/molds/trunk@1454 1136aad2-a195-0410-b898-f5ea1d11b9d8
@@ -893,7 +893,7 @@ void Mndo::CalcCISMatrix(double** matrixCIS) const{ | ||
893 | 893 | } |
894 | 894 | } // end of k-loop |
895 | 895 | |
896 | - // communication to collect all matrix data on rank 0 | |
896 | + // communication to collect all matrix data on head-rank | |
897 | 897 | int mpiHeadRank = MolDS_mpi::MpiProcess::GetInstance()->GetHeadRank(); |
898 | 898 | if(mpiRank == mpiHeadRank){ |
899 | 899 | // receive the matrix data from other ranks |
@@ -905,7 +905,7 @@ void Mndo::CalcCISMatrix(double** matrixCIS) const{ | ||
905 | 905 | } |
906 | 906 | } |
907 | 907 | else{ |
908 | - // send the matrix data to rank-0 | |
908 | + // send the matrix data to head-rank | |
909 | 909 | for(int k=0; k<this->matrixCISdimension; k++){ |
910 | 910 | if(k%mpiSize != mpiRank){continue;} |
911 | 911 | int dest = mpiHeadRank; |
@@ -969,10 +969,13 @@ double Mndo::GetCISCoefficientTwoElecIntegral(int k, | ||
969 | 969 | return value; |
970 | 970 | } |
971 | 971 | |
972 | -void Mndo::MallocTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverlapAOs1stDerivs, | |
973 | - double****** diatomicOverlapAOs2ndDerivs, | |
974 | - double******* diatomicTwoElecTwoCore1stDerivs, | |
975 | - double******** diatomicTwoElecTwoCore2ndDerivs) const{ | |
972 | +void Mndo::MallocTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverlapAOs1stDerivs, | |
973 | + double****** diatomicOverlapAOs2ndDerivs, | |
974 | + double******* diatomicTwoElecTwoCore1stDerivs, | |
975 | + double******** diatomicTwoElecTwoCore2ndDerivs, | |
976 | + double*** tmpRotMat, | |
977 | + double**** tmpRotMat1stDerivs, | |
978 | + double***** tmpDiatomicTwoElecTwoCore) const{ | |
976 | 979 | MallocerFreer::GetInstance()->Malloc<double>(diatomicOverlapAOs1stDerivs, |
977 | 980 | this->molecule->GetNumberAtoms(), |
978 | 981 | OrbitalType_end, |
@@ -999,12 +1002,27 @@ void Mndo::MallocTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverla | ||
999 | 1002 | dxy, |
1000 | 1003 | CartesianType_end, |
1001 | 1004 | CartesianType_end); |
1005 | + MallocerFreer::GetInstance()->Malloc<double>(tmpRotMat, | |
1006 | + OrbitalType_end, | |
1007 | + OrbitalType_end); | |
1008 | + MallocerFreer::GetInstance()->Malloc<double>(tmpRotMat1stDerivs, | |
1009 | + OrbitalType_end, | |
1010 | + OrbitalType_end, | |
1011 | + CartesianType_end); | |
1012 | + MallocerFreer::GetInstance()->Malloc<double>(tmpDiatomicTwoElecTwoCore, | |
1013 | + dxy, | |
1014 | + dxy, | |
1015 | + dxy, | |
1016 | + dxy); | |
1002 | 1017 | } |
1003 | 1018 | |
1004 | 1019 | void Mndo::FreeTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverlapAOs1stDerivs, |
1005 | 1020 | double****** diatomicOverlapAOs2ndDerivs, |
1006 | 1021 | double******* diatomicTwoElecTwoCore1stDerivs, |
1007 | - double******** diatomicTwoElecTwoCore2ndDerivs) const{ | |
1022 | + double******** diatomicTwoElecTwoCore2ndDerivs, | |
1023 | + double*** tmpRotMat, | |
1024 | + double**** tmpRotMat1stDerivs, | |
1025 | + double***** tmpDiatomicTwoElecTwoCore) const{ | |
1008 | 1026 | MallocerFreer::GetInstance()->Free<double>(diatomicOverlapAOs1stDerivs, |
1009 | 1027 | this->molecule->GetNumberAtoms(), |
1010 | 1028 | OrbitalType_end, |
@@ -1031,6 +1049,18 @@ void Mndo::FreeTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverl | ||
1031 | 1049 | dxy, |
1032 | 1050 | CartesianType_end, |
1033 | 1051 | CartesianType_end); |
1052 | + MallocerFreer::GetInstance()->Free<double>(tmpRotMat, | |
1053 | + OrbitalType_end, | |
1054 | + OrbitalType_end); | |
1055 | + MallocerFreer::GetInstance()->Free<double>(tmpRotMat1stDerivs, | |
1056 | + OrbitalType_end, | |
1057 | + OrbitalType_end, | |
1058 | + CartesianType_end); | |
1059 | + MallocerFreer::GetInstance()->Free<double>(tmpDiatomicTwoElecTwoCore, | |
1060 | + dxy, | |
1061 | + dxy, | |
1062 | + dxy, | |
1063 | + dxy); | |
1034 | 1064 | } |
1035 | 1065 | |
1036 | 1066 | // mu and nu is included in atomA' AO. |
@@ -1572,15 +1602,21 @@ void Mndo::CalcHessianSCF(double** hessianSCF, bool isMassWeighted) const{ | ||
1572 | 1602 | stringstream ompErrors; |
1573 | 1603 | #pragma omp parallel |
1574 | 1604 | { |
1575 | - double**** diatomicOverlapAOs1stDerivs = NULL; | |
1576 | - double***** diatomicOverlapAOs2ndDerivs = NULL; | |
1605 | + double**** diatomicOverlapAOs1stDerivs = NULL; | |
1606 | + double***** diatomicOverlapAOs2ndDerivs = NULL; | |
1577 | 1607 | double****** diatomicTwoElecTwoCore1stDerivs = NULL; |
1578 | 1608 | double******* diatomicTwoElecTwoCore2ndDerivs = NULL; |
1579 | - this->MallocTempMatricesEachThreadCalcHessianSCF(&diatomicOverlapAOs1stDerivs, | |
1580 | - &diatomicOverlapAOs2ndDerivs, | |
1581 | - &diatomicTwoElecTwoCore1stDerivs, | |
1582 | - &diatomicTwoElecTwoCore2ndDerivs); | |
1609 | + double** tmpRotMat = NULL; | |
1610 | + double*** tmpRotMat1stDerivs = NULL; | |
1611 | + double**** tmpDiatomicTwoElecTwoCore = NULL; | |
1583 | 1612 | try{ |
1613 | + this->MallocTempMatricesEachThreadCalcHessianSCF(&diatomicOverlapAOs1stDerivs, | |
1614 | + &diatomicOverlapAOs2ndDerivs, | |
1615 | + &diatomicTwoElecTwoCore1stDerivs, | |
1616 | + &diatomicTwoElecTwoCore2ndDerivs, | |
1617 | + &tmpRotMat, | |
1618 | + &tmpRotMat1stDerivs, | |
1619 | + &tmpDiatomicTwoElecTwoCore); | |
1584 | 1620 | #pragma omp for schedule(auto) |
1585 | 1621 | for(int indexAtomA=0; indexAtomA<this->molecule->GetNumberAtoms(); indexAtomA++){ |
1586 | 1622 | const Atom& atomA = *this->molecule->GetAtom(indexAtomA); |
@@ -1598,6 +1634,9 @@ void Mndo::CalcHessianSCF(double** hessianSCF, bool isMassWeighted) const{ | ||
1598 | 1634 | indexAtomA, |
1599 | 1635 | indexAtomB); |
1600 | 1636 | this->CalcDiatomicTwoElecTwoCore1stDerivatives(diatomicTwoElecTwoCore1stDerivs[indexAtomB], |
1637 | + tmpRotMat, | |
1638 | + tmpRotMat1stDerivs, | |
1639 | + tmpDiatomicTwoElecTwoCore, | |
1601 | 1640 | indexAtomA, |
1602 | 1641 | indexAtomB); |
1603 | 1642 | this->CalcDiatomicTwoElecTwoCore2ndDerivatives(diatomicTwoElecTwoCore2ndDerivs[indexAtomB], |
@@ -1659,7 +1698,10 @@ void Mndo::CalcHessianSCF(double** hessianSCF, bool isMassWeighted) const{ | ||
1659 | 1698 | this->FreeTempMatricesEachThreadCalcHessianSCF(&diatomicOverlapAOs1stDerivs, |
1660 | 1699 | &diatomicOverlapAOs2ndDerivs, |
1661 | 1700 | &diatomicTwoElecTwoCore1stDerivs, |
1662 | - &diatomicTwoElecTwoCore2ndDerivs); | |
1701 | + &diatomicTwoElecTwoCore2ndDerivs, | |
1702 | + &tmpRotMat, | |
1703 | + &tmpRotMat1stDerivs, | |
1704 | + &tmpDiatomicTwoElecTwoCore); | |
1663 | 1705 | }// end of omp-region |
1664 | 1706 | // Exception throwing for omp-region |
1665 | 1707 | if(!ompErrors.str().empty()){ |
@@ -1841,10 +1883,17 @@ void Mndo::CalcStaticFirstOrderFock(double* staticFirstOrderFock, | ||
1841 | 1883 | MallocerFreer::GetInstance()->Initialize<double>(staticFirstOrderFock, |
1842 | 1884 | nonRedundantQIndeces.size()+redundantQIndeces.size()); |
1843 | 1885 | double***** diatomicTwoElecTwoCore1stDerivs = NULL; |
1844 | - double*** diatomicOverlapAOs1stDerivs = NULL; | |
1886 | + double*** diatomicOverlapAOs1stDerivs = NULL; | |
1887 | + double** tmpRotMat = NULL; | |
1888 | + double*** tmpRotMat1stDerivs = NULL; | |
1889 | + double**** tmpDiatomicTwoElecTwoCore = NULL; | |
1845 | 1890 | |
1846 | 1891 | try{ |
1847 | - this->MallocTempMatricesStaticFirstOrderFock(&diatomicTwoElecTwoCore1stDerivs, &diatomicOverlapAOs1stDerivs); | |
1892 | + this->MallocTempMatricesStaticFirstOrderFock(&diatomicTwoElecTwoCore1stDerivs, | |
1893 | + &diatomicOverlapAOs1stDerivs, | |
1894 | + &tmpRotMat, | |
1895 | + &tmpRotMat1stDerivs, | |
1896 | + &tmpDiatomicTwoElecTwoCore); | |
1848 | 1897 | const Atom& atomA = *molecule->GetAtom(indexAtomA); |
1849 | 1898 | int firstAOIndexA = atomA.GetFirstAOIndex(); |
1850 | 1899 | int lastAOIndexA = atomA.GetLastAOIndex(); |
@@ -1857,7 +1906,11 @@ void Mndo::CalcStaticFirstOrderFock(double* staticFirstOrderFock, | ||
1857 | 1906 | int coreChargeB = atomB.GetCoreCharge(); |
1858 | 1907 | |
1859 | 1908 | // calc. first derivative of two elec two core interaction |
1860 | - this->CalcDiatomicTwoElecTwoCore1stDerivatives(diatomicTwoElecTwoCore1stDerivs, indexAtomA, indexAtomB); | |
1909 | + this->CalcDiatomicTwoElecTwoCore1stDerivatives(diatomicTwoElecTwoCore1stDerivs, | |
1910 | + tmpRotMat, | |
1911 | + tmpRotMat1stDerivs, | |
1912 | + tmpDiatomicTwoElecTwoCore, | |
1913 | + indexAtomA, indexAtomB); | |
1861 | 1914 | // calc. first derivative of overlapAOs. |
1862 | 1915 | this->CalcDiatomicOverlapAOs1stDerivatives(diatomicOverlapAOs1stDerivs, atomA, atomB); |
1863 | 1916 |
@@ -1950,10 +2003,18 @@ void Mndo::CalcStaticFirstOrderFock(double* staticFirstOrderFock, | ||
1950 | 2003 | } |
1951 | 2004 | } |
1952 | 2005 | catch(MolDSException ex){ |
1953 | - this->FreeTempMatricesStaticFirstOrderFock(&diatomicTwoElecTwoCore1stDerivs, &diatomicOverlapAOs1stDerivs); | |
2006 | + this->FreeTempMatricesStaticFirstOrderFock(&diatomicTwoElecTwoCore1stDerivs, | |
2007 | + &diatomicOverlapAOs1stDerivs, | |
2008 | + &tmpRotMat, | |
2009 | + &tmpRotMat1stDerivs, | |
2010 | + &tmpDiatomicTwoElecTwoCore); | |
1954 | 2011 | throw ex; |
1955 | 2012 | } |
1956 | - this->FreeTempMatricesStaticFirstOrderFock(&diatomicTwoElecTwoCore1stDerivs, &diatomicOverlapAOs1stDerivs); | |
2013 | + this->FreeTempMatricesStaticFirstOrderFock(&diatomicTwoElecTwoCore1stDerivs, | |
2014 | + &diatomicOverlapAOs1stDerivs, | |
2015 | + &tmpRotMat, | |
2016 | + &tmpRotMat1stDerivs, | |
2017 | + &tmpDiatomicTwoElecTwoCore); | |
1957 | 2018 | |
1958 | 2019 | /* |
1959 | 2020 | printf("staticFirstOrderFock(atomA:%d axis:%s)\n",indexAtomA,CartesianTypeStr(axisA)); |
@@ -1964,7 +2025,10 @@ void Mndo::CalcStaticFirstOrderFock(double* staticFirstOrderFock, | ||
1964 | 2025 | } |
1965 | 2026 | |
1966 | 2027 | void Mndo::MallocTempMatricesStaticFirstOrderFock(double****** diatomicTwoElecTwoCore1stDeriv, |
1967 | - double**** diatomicOverlapAOs1stDeriv)const{ | |
2028 | + double**** diatomicOverlapAOs1stDeriv, | |
2029 | + double*** tmpRotMat, | |
2030 | + double**** tmpRotMat1stDerivs, | |
2031 | + double***** tmpDiatomicTwoElecTwoCore)const{ | |
1968 | 2032 | MallocerFreer::GetInstance()->Malloc<double>(diatomicTwoElecTwoCore1stDeriv, |
1969 | 2033 | dxy, |
1970 | 2034 | dxy, |
@@ -1975,10 +2039,25 @@ void Mndo::MallocTempMatricesStaticFirstOrderFock(double****** diatomicTwoElecTw | ||
1975 | 2039 | OrbitalType_end, |
1976 | 2040 | OrbitalType_end, |
1977 | 2041 | CartesianType_end); |
2042 | + MallocerFreer::GetInstance()->Malloc<double>(tmpRotMat, | |
2043 | + OrbitalType_end, | |
2044 | + OrbitalType_end); | |
2045 | + MallocerFreer::GetInstance()->Malloc<double>(tmpRotMat1stDerivs, | |
2046 | + OrbitalType_end, | |
2047 | + OrbitalType_end, | |
2048 | + CartesianType_end); | |
2049 | + MallocerFreer::GetInstance()->Malloc<double>(tmpDiatomicTwoElecTwoCore, | |
2050 | + dxy, | |
2051 | + dxy, | |
2052 | + dxy, | |
2053 | + dxy); | |
1978 | 2054 | } |
1979 | 2055 | |
1980 | 2056 | void Mndo::FreeTempMatricesStaticFirstOrderFock(double****** diatomicTwoElecTwoCore1stDeriv, |
1981 | - double**** diatomicOverlapAOs1stDeriv)const{ | |
2057 | + double**** diatomicOverlapAOs1stDeriv, | |
2058 | + double*** tmpRotMat, | |
2059 | + double**** tmpRotMat1stDerivs, | |
2060 | + double***** tmpDiatomicTwoElecTwoCore)const{ | |
1982 | 2061 | MallocerFreer::GetInstance()->Free<double>(diatomicTwoElecTwoCore1stDeriv, |
1983 | 2062 | dxy, |
1984 | 2063 | dxy, |
@@ -1989,6 +2068,18 @@ void Mndo::FreeTempMatricesStaticFirstOrderFock(double****** diatomicTwoElecTwoC | ||
1989 | 2068 | OrbitalType_end, |
1990 | 2069 | OrbitalType_end, |
1991 | 2070 | CartesianType_end); |
2071 | + MallocerFreer::GetInstance()->Free<double>(tmpRotMat, | |
2072 | + OrbitalType_end, | |
2073 | + OrbitalType_end); | |
2074 | + MallocerFreer::GetInstance()->Free<double>(tmpRotMat1stDerivs, | |
2075 | + OrbitalType_end, | |
2076 | + OrbitalType_end, | |
2077 | + CartesianType_end); | |
2078 | + MallocerFreer::GetInstance()->Free<double>(tmpDiatomicTwoElecTwoCore, | |
2079 | + dxy, | |
2080 | + dxy, | |
2081 | + dxy, | |
2082 | + dxy); | |
1992 | 2083 | } |
1993 | 2084 | |
1994 | 2085 | // see (40) - (46) in [PT_1996]. |
@@ -2287,8 +2378,15 @@ void Mndo::CalcForce(const vector<int>& elecStates){ | ||
2287 | 2378 | { |
2288 | 2379 | double***** diatomicTwoElecTwoCore1stDerivs = NULL; |
2289 | 2380 | double*** diatomicOverlapAOs1stDerivs = NULL; |
2381 | + double** tmpRotMat = NULL; | |
2382 | + double*** tmpRotMat1stDerivs = NULL; | |
2383 | + double**** tmpDiatomicTwoElecTwoCore = NULL; | |
2290 | 2384 | try{ |
2291 | - this->MallocTempMatricesCalcForce(&diatomicOverlapAOs1stDerivs, &diatomicTwoElecTwoCore1stDerivs); | |
2385 | + this->MallocTempMatricesCalcForce(&diatomicOverlapAOs1stDerivs, | |
2386 | + &diatomicTwoElecTwoCore1stDerivs, | |
2387 | + &tmpRotMat, | |
2388 | + &tmpRotMat1stDerivs, | |
2389 | + &tmpDiatomicTwoElecTwoCore); | |
2292 | 2390 | |
2293 | 2391 | #pragma omp for schedule(auto) |
2294 | 2392 | for(int b=0; b<this->molecule->GetNumberAtoms(); b++){ |
@@ -2300,7 +2398,11 @@ void Mndo::CalcForce(const vector<int>& elecStates){ | ||
2300 | 2398 | // calc. first derivative of overlapAOs. |
2301 | 2399 | this->CalcDiatomicOverlapAOs1stDerivatives(diatomicOverlapAOs1stDerivs, atomA, atomB); |
2302 | 2400 | // calc. first derivative of two elec two core interaction |
2303 | - this->CalcDiatomicTwoElecTwoCore1stDerivatives(diatomicTwoElecTwoCore1stDerivs, a, b); | |
2401 | + this->CalcDiatomicTwoElecTwoCore1stDerivatives(diatomicTwoElecTwoCore1stDerivs, | |
2402 | + tmpRotMat, | |
2403 | + tmpRotMat1stDerivs, | |
2404 | + tmpDiatomicTwoElecTwoCore, | |
2405 | + a, b); | |
2304 | 2406 | |
2305 | 2407 | // core repulsion part |
2306 | 2408 | double coreRepulsion[CartesianType_end] = {0.0,0.0,0.0}; |
@@ -2407,7 +2509,11 @@ void Mndo::CalcForce(const vector<int>& elecStates){ | ||
2407 | 2509 | #pragma omp critical |
2408 | 2510 | ex.Serialize(ompErrors); |
2409 | 2511 | } |
2410 | - this->FreeTempMatricesCalcForce(&diatomicOverlapAOs1stDerivs, &diatomicTwoElecTwoCore1stDerivs); | |
2512 | + this->FreeTempMatricesCalcForce(&diatomicOverlapAOs1stDerivs, | |
2513 | + &diatomicTwoElecTwoCore1stDerivs, | |
2514 | + &tmpRotMat, | |
2515 | + &tmpRotMat1stDerivs, | |
2516 | + &tmpDiatomicTwoElecTwoCore); | |
2411 | 2517 | } // end of omp-parallelized region |
2412 | 2518 | // Exception throwing for omp-region |
2413 | 2519 | if(!ompErrors.str().empty()){ |
@@ -2418,10 +2524,13 @@ void Mndo::CalcForce(const vector<int>& elecStates){ | ||
2418 | 2524 | // communication to reduce thsi->matrixForce on all node (namely, all_reduce) |
2419 | 2525 | int numTransported = elecStates.size()*this->molecule->GetNumberAtoms()*CartesianType_end; |
2420 | 2526 | MolDS_mpi::MpiProcess::GetInstance()->AllReduce(&this->matrixForce[0][0][0], numTransported, std::plus<double>()); |
2421 | - | |
2422 | 2527 | } |
2423 | 2528 | |
2424 | -void Mndo::MallocTempMatricesCalcForce(double**** diatomicOverlapAOs1stDerivs, double****** diatomicTwoElecTwoCore1stDerivs) const{ | |
2529 | +void Mndo::MallocTempMatricesCalcForce(double**** diatomicOverlapAOs1stDerivs, | |
2530 | + double****** diatomicTwoElecTwoCore1stDerivs, | |
2531 | + double*** tmpRotMat, | |
2532 | + double**** tmpRotMat1stDerivs, | |
2533 | + double***** tmpDiatomicTwoElecTwoCore) const{ | |
2425 | 2534 | MallocerFreer::GetInstance()->Malloc<double>(diatomicOverlapAOs1stDerivs, |
2426 | 2535 | OrbitalType_end, |
2427 | 2536 | OrbitalType_end, |
@@ -2432,9 +2541,25 @@ void Mndo::MallocTempMatricesCalcForce(double**** diatomicOverlapAOs1stDerivs, d | ||
2432 | 2541 | dxy, |
2433 | 2542 | dxy, |
2434 | 2543 | CartesianType_end); |
2544 | + MallocerFreer::GetInstance()->Malloc<double>(tmpRotMat, | |
2545 | + OrbitalType_end, | |
2546 | + OrbitalType_end); | |
2547 | + MallocerFreer::GetInstance()->Malloc<double>(tmpRotMat1stDerivs, | |
2548 | + OrbitalType_end, | |
2549 | + OrbitalType_end, | |
2550 | + CartesianType_end); | |
2551 | + MallocerFreer::GetInstance()->Malloc<double>(tmpDiatomicTwoElecTwoCore, | |
2552 | + dxy, | |
2553 | + dxy, | |
2554 | + dxy, | |
2555 | + dxy); | |
2435 | 2556 | } |
2436 | 2557 | |
2437 | -void Mndo::FreeTempMatricesCalcForce(double**** diatomicOverlapAOs1stDerivs, double****** diatomicTwoElecTwoCore1stDerivs) const{ | |
2558 | +void Mndo::FreeTempMatricesCalcForce(double**** diatomicOverlapAOs1stDerivs, | |
2559 | + double****** diatomicTwoElecTwoCore1stDerivs, | |
2560 | + double*** tmpRotMat, | |
2561 | + double**** tmpRotMat1stDerivs, | |
2562 | + double***** tmpDiatomicTwoElecTwoCore) const{ | |
2438 | 2563 | MallocerFreer::GetInstance()->Free<double>(diatomicOverlapAOs1stDerivs, |
2439 | 2564 | OrbitalType_end, |
2440 | 2565 | OrbitalType_end, |
@@ -2445,6 +2570,18 @@ void Mndo::FreeTempMatricesCalcForce(double**** diatomicOverlapAOs1stDerivs, dou | ||
2445 | 2570 | dxy, |
2446 | 2571 | dxy, |
2447 | 2572 | CartesianType_end); |
2573 | + MallocerFreer::GetInstance()->Free<double>(tmpRotMat, | |
2574 | + OrbitalType_end, | |
2575 | + OrbitalType_end); | |
2576 | + MallocerFreer::GetInstance()->Free<double>(tmpRotMat1stDerivs, | |
2577 | + OrbitalType_end, | |
2578 | + OrbitalType_end, | |
2579 | + CartesianType_end); | |
2580 | + MallocerFreer::GetInstance()->Free<double>(tmpDiatomicTwoElecTwoCore, | |
2581 | + dxy, | |
2582 | + dxy, | |
2583 | + dxy, | |
2584 | + dxy); | |
2448 | 2585 | } |
2449 | 2586 | |
2450 | 2587 | // see (18) in [PT_1997] |
@@ -3174,6 +3311,9 @@ void Mndo::CalcDiatomicTwoElecTwoCore(double**** matrix, int indexAtomA, int ind | ||
3174 | 3311 | // Note taht d-orbital cannot be treated, |
3175 | 3312 | // that is, matrix[dxy][dxy][dxy][dxy][CartesianType_end] cannot be treatable. |
3176 | 3313 | void Mndo::CalcDiatomicTwoElecTwoCore1stDerivatives(double***** matrix, |
3314 | + double** tmpRotMat, | |
3315 | + double*** tmpRotMat1stDerivs, | |
3316 | + double**** tmpDiatomicTwoElecTwoCore, | |
3177 | 3317 | int indexAtomA, |
3178 | 3318 | int indexAtomB) const{ |
3179 | 3319 | const Atom& atomA = *this->molecule->GetAtom(indexAtomA); |
@@ -3200,65 +3340,48 @@ void Mndo::CalcDiatomicTwoElecTwoCore1stDerivatives(double***** matrix, | ||
3200 | 3340 | dxy, |
3201 | 3341 | CartesianType_end); |
3202 | 3342 | |
3203 | - double** rotatingMatrix = NULL; | |
3204 | - double*** rotMat1stDerivatives = NULL; | |
3205 | - double**** diatomicTwoElecTwoCore = NULL; | |
3206 | - try{ | |
3207 | - this->MallocDiatomicTwoElecTwoCore1stDeriTemps(&rotatingMatrix, | |
3208 | - &rotMat1stDerivatives, | |
3209 | - &diatomicTwoElecTwoCore); | |
3210 | - // calclation in diatomic frame | |
3211 | - for(int mu=0; mu<atomA.GetValenceSize(); mu++){ | |
3212 | - for(int nu=mu; nu<atomA.GetValenceSize(); nu++){ | |
3213 | - for(int lambda=0; lambda<atomB.GetValenceSize(); lambda++){ | |
3214 | - for(int sigma=lambda; sigma<atomB.GetValenceSize(); sigma++){ | |
3215 | - for(int dimA=0; dimA<CartesianType_end; dimA++){ | |
3216 | - matrix[mu][nu][lambda][sigma][dimA] | |
3217 | - = this->GetNddoRepulsionIntegral1stDerivative( | |
3218 | - atomA, | |
3219 | - atomA.GetValence(mu), | |
3220 | - atomA.GetValence(nu), | |
3221 | - atomB, | |
3222 | - atomB.GetValence(lambda), | |
3223 | - atomB.GetValence(sigma), | |
3224 | - static_cast<CartesianType>(dimA)); | |
3225 | - matrix[nu][mu][lambda][sigma][dimA] = matrix[mu][nu][lambda][sigma][dimA]; | |
3226 | - matrix[nu][mu][sigma][lambda][dimA] = matrix[mu][nu][lambda][sigma][dimA]; | |
3227 | - matrix[mu][nu][sigma][lambda][dimA] = matrix[mu][nu][lambda][sigma][dimA]; | |
3228 | - } | |
3229 | - diatomicTwoElecTwoCore[mu][nu][lambda][sigma] | |
3230 | - = this->GetNddoRepulsionIntegral( | |
3343 | + // calclation in diatomic frame | |
3344 | + for(int mu=0; mu<atomA.GetValenceSize(); mu++){ | |
3345 | + for(int nu=mu; nu<atomA.GetValenceSize(); nu++){ | |
3346 | + for(int lambda=0; lambda<atomB.GetValenceSize(); lambda++){ | |
3347 | + for(int sigma=lambda; sigma<atomB.GetValenceSize(); sigma++){ | |
3348 | + for(int dimA=0; dimA<CartesianType_end; dimA++){ | |
3349 | + matrix[mu][nu][lambda][sigma][dimA] | |
3350 | + = this->GetNddoRepulsionIntegral1stDerivative( | |
3231 | 3351 | atomA, |
3232 | 3352 | atomA.GetValence(mu), |
3233 | 3353 | atomA.GetValence(nu), |
3234 | 3354 | atomB, |
3235 | 3355 | atomB.GetValence(lambda), |
3236 | - atomB.GetValence(sigma)); | |
3237 | - diatomicTwoElecTwoCore[nu][mu][lambda][sigma] = diatomicTwoElecTwoCore[mu][nu][lambda][sigma]; | |
3238 | - diatomicTwoElecTwoCore[nu][mu][sigma][lambda] = diatomicTwoElecTwoCore[mu][nu][lambda][sigma]; | |
3239 | - diatomicTwoElecTwoCore[mu][nu][sigma][lambda] = diatomicTwoElecTwoCore[mu][nu][lambda][sigma]; | |
3240 | - } | |
3356 | + atomB.GetValence(sigma), | |
3357 | + static_cast<CartesianType>(dimA)); | |
3358 | + matrix[nu][mu][lambda][sigma][dimA] = matrix[mu][nu][lambda][sigma][dimA]; | |
3359 | + matrix[nu][mu][sigma][lambda][dimA] = matrix[mu][nu][lambda][sigma][dimA]; | |
3360 | + matrix[mu][nu][sigma][lambda][dimA] = matrix[mu][nu][lambda][sigma][dimA]; | |
3361 | + } | |
3362 | + tmpDiatomicTwoElecTwoCore[mu][nu][lambda][sigma] | |
3363 | + = this->GetNddoRepulsionIntegral( | |
3364 | + atomA, | |
3365 | + atomA.GetValence(mu), | |
3366 | + atomA.GetValence(nu), | |
3367 | + atomB, | |
3368 | + atomB.GetValence(lambda), | |
3369 | + atomB.GetValence(sigma)); | |
3370 | + tmpDiatomicTwoElecTwoCore[nu][mu][lambda][sigma] = tmpDiatomicTwoElecTwoCore[mu][nu][lambda][sigma]; | |
3371 | + tmpDiatomicTwoElecTwoCore[nu][mu][sigma][lambda] = tmpDiatomicTwoElecTwoCore[mu][nu][lambda][sigma]; | |
3372 | + tmpDiatomicTwoElecTwoCore[mu][nu][sigma][lambda] = tmpDiatomicTwoElecTwoCore[mu][nu][lambda][sigma]; | |
3241 | 3373 | } |
3242 | 3374 | } |
3243 | 3375 | } |
3244 | - | |
3245 | - // rotate matirix into the space frame | |
3246 | - this->CalcRotatingMatrix(rotatingMatrix, atomA, atomB); | |
3247 | - this->CalcRotatingMatrix1stDerivatives(rotMat1stDerivatives, atomA, atomB); | |
3248 | - this->RotateDiatomicTwoElecTwoCore1stDerivativesToSpaceFrame(matrix, | |
3249 | - diatomicTwoElecTwoCore, | |
3250 | - rotatingMatrix, | |
3251 | - rotMat1stDerivatives); | |
3252 | 3376 | } |
3253 | - catch(MolDSException ex){ | |
3254 | - this->FreeDiatomicTwoElecTwoCore1stDeriTemps(&rotatingMatrix, | |
3255 | - &rotMat1stDerivatives, | |
3256 | - &diatomicTwoElecTwoCore); | |
3257 | - throw ex; | |
3258 | - } | |
3259 | - this->FreeDiatomicTwoElecTwoCore1stDeriTemps(&rotatingMatrix, | |
3260 | - &rotMat1stDerivatives, | |
3261 | - &diatomicTwoElecTwoCore); | |
3377 | + | |
3378 | + // rotate matirix into the space frame | |
3379 | + this->CalcRotatingMatrix(tmpRotMat, atomA, atomB); | |
3380 | + this->CalcRotatingMatrix1stDerivatives(tmpRotMat1stDerivs, atomA, atomB); | |
3381 | + this->RotateDiatomicTwoElecTwoCore1stDerivativesToSpaceFrame(matrix, | |
3382 | + tmpDiatomicTwoElecTwoCore, | |
3383 | + tmpRotMat, | |
3384 | + tmpRotMat1stDerivs); | |
3262 | 3385 | } |
3263 | 3386 | |
3264 | 3387 | // Calculation of second derivatives of the two electrons two cores integral in space fixed frame, |
@@ -150,11 +150,17 @@ private: | ||
150 | 150 | void MallocTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverlapAOs1stDerivs, |
151 | 151 | double****** diatomicOverlapAOs2ndDerivs, |
152 | 152 | double******* diatomicTwoElecTwoCore1stDerivs, |
153 | - double******** diatomicTwoElecTwoCore2ndDerivs) const; | |
153 | + double******** diatomicTwoElecTwoCore2ndDerivs, | |
154 | + double*** tmpRotMat, | |
155 | + double**** tmpRotMat1stDerivs, | |
156 | + double***** tmpDiatomicTwoElecTwo) const; | |
154 | 157 | void FreeTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverlapAOs1stDerivs, |
155 | 158 | double****** diatomicOverlapAOs2ndDerivs, |
156 | - double******* diatomicTwoElecTwoCore1stDeriv, | |
157 | - double******** diatomicTwoElecTwoCore2ndDerivs) const; | |
159 | + double******* diatomicTwoElecTwoCore1stDerivs, | |
160 | + double******** diatomicTwoElecTwoCore2ndDerivs, | |
161 | + double*** tmpRotMat, | |
162 | + double**** tmpRotMat1stDerivs, | |
163 | + double***** tmpDiatomicTwoElecTwo) const; | |
158 | 164 | double GetAuxiliaryHessianElement1(int mu, |
159 | 165 | int nu, |
160 | 166 | int indexAtomA, |
@@ -241,9 +247,15 @@ private: | ||
241 | 247 | int indexAtomA, |
242 | 248 | MolDS_base::CartesianType axisA) const; |
243 | 249 | void MallocTempMatricesStaticFirstOrderFock(double****** diatomicTwoElecTwoCore1stDeriv, |
244 | - double**** diatomicOverlapAOs1stDeriv) const; | |
250 | + double**** diatomicOverlapAOs1stDeriv, | |
251 | + double*** tmpRotMat, | |
252 | + double**** tmpRotMat1stDerivs, | |
253 | + double***** tmpDiatomicTwoElecTwo) const; | |
245 | 254 | void FreeTempMatricesStaticFirstOrderFock(double****** diatomicTwoElecTwoCore1stDeriv, |
246 | - double**** diatomicOverlapAOs1stDeriv) const; | |
255 | + double**** diatomicOverlapAOs1stDeriv, | |
256 | + double*** tmpRotMat, | |
257 | + double**** tmpRotMat1stDerivs, | |
258 | + double***** tmpDiatomicTwoElecTwo) const; | |
247 | 259 | void CalcMatrixCPHF(double** matrixCPHF, |
248 | 260 | const std::vector<MoIndexPair>& nonRedundantQIndeces, |
249 | 261 | const std::vector<MoIndexPair>& redundantQIndeces) const; |
@@ -266,6 +278,9 @@ private: | ||
266 | 278 | MolDS_base::CartesianType axisA) const; |
267 | 279 | void CalcDiatomicTwoElecTwoCore(double**** matrix, int indexAtomA, int indexAtomB) const; |
268 | 280 | void CalcDiatomicTwoElecTwoCore1stDerivatives(double***** matrix, |
281 | + double** tmpRotMat, | |
282 | + double*** tmpRotMat1stDerivs, | |
283 | + double**** tmpDiatomicTwoElecTwoCore, | |
269 | 284 | int indexAtomA, |
270 | 285 | int indexAtomB) const; |
271 | 286 | void CalcDiatomicTwoElecTwoCore2ndDerivatives(double****** matrix, |
@@ -356,9 +371,15 @@ private: | ||
356 | 371 | MolDS_base::MultipoleType multipoleB, |
357 | 372 | double rAB) const; |
358 | 373 | void MallocTempMatricesCalcForce(double**** diatomicOverlapAOs1stDerivs, |
359 | - double****** diatomiTwoElecTwoCore1stDerivs) const; | |
374 | + double****** diatomiTwoElecTwoCore1stDerivs, | |
375 | + double*** tmpRotMat, | |
376 | + double**** tmpRotMat1stDerivs, | |
377 | + double***** tmpDiatomicTwoElecTwoCore) const; | |
360 | 378 | void FreeTempMatricesCalcForce(double**** diatomicOverlapAOs1stDerivs, |
361 | - double****** diatomiTwoElecTwoCore1stDerivs) const; | |
379 | + double****** diatomiTwoElecTwoCore1stDerivs, | |
380 | + double*** tmpRotMat, | |
381 | + double**** tmpRotMat1stDerivs, | |
382 | + double***** tmpDiatomicTwoElecTwoCore) const; | |
362 | 383 | void CalcForceSCFElecCoreAttractionPart(double* force, |
363 | 384 | int indexAtomA, |
364 | 385 | int indexAtomB, |