• R/O
  • HTTP
  • SSH
  • HTTPS

Commit

Tags
No Tags

Frequently used words (click to add to your profile)

javac++androidlinuxc#windowsobjective-ccocoa誰得qtpythonphprubygameguibathyscaphec計画中(planning stage)翻訳omegatframeworktwitterdomtestvb.netdirectxゲームエンジンbtronarduinopreviewer

Commit MetaInfo

Revision665c8e35162062dab0dec240ae4990d1493ec623 (tree)
Time2013-08-05 18:49:09
AuthorMikiya Fujii <mikiya.fujii@gmai...>
CommiterMikiya Fujii

Log Message

refactoring Mndo::CalcDiatomicTwoElecTwoCore2ndDerivatives also. #31814

git-svn-id: https://svn.sourceforge.jp/svnroot/molds/trunk@1455 1136aad2-a195-0410-b898-f5ea1d11b9d8

Change Summary

Incremental Difference

--- a/src/mndo/Mndo.cpp
+++ b/src/mndo/Mndo.cpp
@@ -975,7 +975,9 @@ void Mndo::MallocTempMatricesEachThreadCalcHessianSCF(double***** diatomicOve
975975 double******** diatomicTwoElecTwoCore2ndDerivs,
976976 double*** tmpRotMat,
977977 double**** tmpRotMat1stDerivs,
978- double***** tmpDiatomicTwoElecTwoCore) const{
978+ double***** tmpRotMat2ndDerivs,
979+ double***** tmpDiatomicTwoElecTwoCore,
980+ double****** tmpDiatomicTwoElecTwoCore1stDerivs) const{
979981 MallocerFreer::GetInstance()->Malloc<double>(diatomicOverlapAOs1stDerivs,
980982 this->molecule->GetNumberAtoms(),
981983 OrbitalType_end,
@@ -1009,11 +1011,22 @@ void Mndo::MallocTempMatricesEachThreadCalcHessianSCF(double***** diatomicOve
10091011 OrbitalType_end,
10101012 OrbitalType_end,
10111013 CartesianType_end);
1014+ MallocerFreer::GetInstance()->Malloc<double>(tmpRotMat2ndDerivs,
1015+ OrbitalType_end,
1016+ OrbitalType_end,
1017+ CartesianType_end,
1018+ CartesianType_end);
10121019 MallocerFreer::GetInstance()->Malloc<double>(tmpDiatomicTwoElecTwoCore,
10131020 dxy,
10141021 dxy,
10151022 dxy,
10161023 dxy);
1024+ MallocerFreer::GetInstance()->Malloc<double>(tmpDiatomicTwoElecTwoCore1stDerivs,
1025+ dxy,
1026+ dxy,
1027+ dxy,
1028+ dxy,
1029+ CartesianType_end);
10171030 }
10181031
10191032 void Mndo::FreeTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverlapAOs1stDerivs,
@@ -1022,7 +1035,9 @@ void Mndo::FreeTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverl
10221035 double******** diatomicTwoElecTwoCore2ndDerivs,
10231036 double*** tmpRotMat,
10241037 double**** tmpRotMat1stDerivs,
1025- double***** tmpDiatomicTwoElecTwoCore) const{
1038+ double***** tmpRotMat2ndDerivs,
1039+ double***** tmpDiatomicTwoElecTwoCore,
1040+ double****** tmpDiatomicTwoElecTwoCore1stDerivs) const{
10261041 MallocerFreer::GetInstance()->Free<double>(diatomicOverlapAOs1stDerivs,
10271042 this->molecule->GetNumberAtoms(),
10281043 OrbitalType_end,
@@ -1056,11 +1071,22 @@ void Mndo::FreeTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverl
10561071 OrbitalType_end,
10571072 OrbitalType_end,
10581073 CartesianType_end);
1074+ MallocerFreer::GetInstance()->Free<double>(tmpRotMat2ndDerivs,
1075+ OrbitalType_end,
1076+ OrbitalType_end,
1077+ CartesianType_end,
1078+ CartesianType_end);
10591079 MallocerFreer::GetInstance()->Free<double>(tmpDiatomicTwoElecTwoCore,
10601080 dxy,
10611081 dxy,
10621082 dxy,
10631083 dxy);
1084+ MallocerFreer::GetInstance()->Free<double>(tmpDiatomicTwoElecTwoCore1stDerivs,
1085+ dxy,
1086+ dxy,
1087+ dxy,
1088+ dxy,
1089+ CartesianType_end);
10641090 }
10651091
10661092 // mu and nu is included in atomA' AO.
@@ -1602,13 +1628,15 @@ void Mndo::CalcHessianSCF(double** hessianSCF, bool isMassWeighted) const{
16021628 stringstream ompErrors;
16031629 #pragma omp parallel
16041630 {
1605- double**** diatomicOverlapAOs1stDerivs = NULL;
1606- double***** diatomicOverlapAOs2ndDerivs = NULL;
1607- double****** diatomicTwoElecTwoCore1stDerivs = NULL;
1608- double******* diatomicTwoElecTwoCore2ndDerivs = NULL;
1609- double** tmpRotMat = NULL;
1610- double*** tmpRotMat1stDerivs = NULL;
1611- double**** tmpDiatomicTwoElecTwoCore = NULL;
1631+ double**** diatomicOverlapAOs1stDerivs = NULL;
1632+ double***** diatomicOverlapAOs2ndDerivs = NULL;
1633+ double****** diatomicTwoElecTwoCore1stDerivs = NULL;
1634+ double******* diatomicTwoElecTwoCore2ndDerivs = NULL;
1635+ double** tmpRotMat = NULL;
1636+ double*** tmpRotMat1stDerivs = NULL;
1637+ double**** tmpRotMat2ndDerivs = NULL;
1638+ double**** tmpDiatomicTwoElecTwoCore = NULL;
1639+ double***** tmpDiatomicTwoElecTwoCore1stDerivs = NULL;
16121640 try{
16131641 this->MallocTempMatricesEachThreadCalcHessianSCF(&diatomicOverlapAOs1stDerivs,
16141642 &diatomicOverlapAOs2ndDerivs,
@@ -1616,7 +1644,9 @@ void Mndo::CalcHessianSCF(double** hessianSCF, bool isMassWeighted) const{
16161644 &diatomicTwoElecTwoCore2ndDerivs,
16171645 &tmpRotMat,
16181646 &tmpRotMat1stDerivs,
1619- &tmpDiatomicTwoElecTwoCore);
1647+ &tmpRotMat2ndDerivs,
1648+ &tmpDiatomicTwoElecTwoCore,
1649+ &tmpDiatomicTwoElecTwoCore1stDerivs);
16201650 #pragma omp for schedule(auto)
16211651 for(int indexAtomA=0; indexAtomA<this->molecule->GetNumberAtoms(); indexAtomA++){
16221652 const Atom& atomA = *this->molecule->GetAtom(indexAtomA);
@@ -1640,6 +1670,11 @@ void Mndo::CalcHessianSCF(double** hessianSCF, bool isMassWeighted) const{
16401670 indexAtomA,
16411671 indexAtomB);
16421672 this->CalcDiatomicTwoElecTwoCore2ndDerivatives(diatomicTwoElecTwoCore2ndDerivs[indexAtomB],
1673+ tmpRotMat,
1674+ tmpRotMat1stDerivs,
1675+ tmpRotMat2ndDerivs,
1676+ tmpDiatomicTwoElecTwoCore,
1677+ tmpDiatomicTwoElecTwoCore1stDerivs,
16431678 indexAtomA,
16441679 indexAtomB);
16451680 }
@@ -1701,7 +1736,9 @@ void Mndo::CalcHessianSCF(double** hessianSCF, bool isMassWeighted) const{
17011736 &diatomicTwoElecTwoCore2ndDerivs,
17021737 &tmpRotMat,
17031738 &tmpRotMat1stDerivs,
1704- &tmpDiatomicTwoElecTwoCore);
1739+ &tmpRotMat2ndDerivs,
1740+ &tmpDiatomicTwoElecTwoCore,
1741+ &tmpDiatomicTwoElecTwoCore1stDerivs);
17051742 }// end of omp-region
17061743 // Exception throwing for omp-region
17071744 if(!ompErrors.str().empty()){
@@ -3393,8 +3430,13 @@ void Mndo::CalcDiatomicTwoElecTwoCore1stDerivatives(double***** matrix,
33933430 // Note taht d-orbital cannot be treated,
33943431 // that is, matrix[dxy][dxy][dxy][dxy][CartesianType_end][CartesianType_end] cannot be treatable.
33953432 void Mndo::CalcDiatomicTwoElecTwoCore2ndDerivatives(double****** matrix,
3396- int indexAtomA,
3397- int indexAtomB) const{
3433+ double** tmpRotMat,
3434+ double*** tmpRotMat1stDerivs,
3435+ double**** tmpRotMat2ndDerivs,
3436+ double**** tmpDiatomicTwoElecTwoCore,
3437+ double***** tmpDiatomicTwoElecTwoCore1stDerivs,
3438+ int indexAtomA,
3439+ int indexAtomB) const{
33983440 const Atom& atomA = *this->molecule->GetAtom(indexAtomA);
33993441 const Atom& atomB = *this->molecule->GetAtom(indexAtomB);
34003442 if(indexAtomA == indexAtomB){
@@ -3420,148 +3462,57 @@ void Mndo::CalcDiatomicTwoElecTwoCore2ndDerivatives(double****** matrix,
34203462 CartesianType_end,
34213463 CartesianType_end);
34223464
3423- double** rotatingMatrix = NULL;
3424- double*** rotMat1stDerivatives = NULL;
3425- double**** rotMat2ndDerivatives = NULL;
3426- double**** diatomicTwoElecTwoCore = NULL;
3427- double***** diatomicTwoElecTwoCore1stDerivatives = NULL;
3428- try{
3429- this->MallocDiatomicTwoElecTwoCore2ndDeriTemps(&rotatingMatrix,
3430- &rotMat1stDerivatives,
3431- &rotMat2ndDerivatives,
3432- &diatomicTwoElecTwoCore,
3433- &diatomicTwoElecTwoCore1stDerivatives);
3434- // calclation in diatomic frame
3435- for(int mu=0; mu<atomA.GetValenceSize(); mu++){
3436- for(int nu=0; nu<atomA.GetValenceSize(); nu++){
3437- for(int lambda=0; lambda<atomB.GetValenceSize(); lambda++){
3438- for(int sigma=0; sigma<atomB.GetValenceSize(); sigma++){
3439- for(int dimA1=0; dimA1<CartesianType_end; dimA1++){
3440- for(int dimA2=0; dimA2<CartesianType_end; dimA2++){
3441- matrix[mu][nu][lambda][sigma][dimA1][dimA2]
3442- = this->GetNddoRepulsionIntegral2ndDerivative(
3443- atomA,
3444- atomA.GetValence(mu),
3445- atomA.GetValence(nu),
3446- atomB,
3447- atomB.GetValence(lambda),
3448- atomB.GetValence(sigma),
3449- static_cast<CartesianType>(dimA1),
3450- static_cast<CartesianType>(dimA2));
3451- }
3452- diatomicTwoElecTwoCore1stDerivatives[mu][nu][lambda][sigma][dimA1]
3453- = this->GetNddoRepulsionIntegral1stDerivative(
3465+ // calclation in diatomic frame
3466+ for(int mu=0; mu<atomA.GetValenceSize(); mu++){
3467+ for(int nu=0; nu<atomA.GetValenceSize(); nu++){
3468+ for(int lambda=0; lambda<atomB.GetValenceSize(); lambda++){
3469+ for(int sigma=0; sigma<atomB.GetValenceSize(); sigma++){
3470+ for(int dimA1=0; dimA1<CartesianType_end; dimA1++){
3471+ for(int dimA2=0; dimA2<CartesianType_end; dimA2++){
3472+ matrix[mu][nu][lambda][sigma][dimA1][dimA2]
3473+ = this->GetNddoRepulsionIntegral2ndDerivative(
34543474 atomA,
34553475 atomA.GetValence(mu),
34563476 atomA.GetValence(nu),
34573477 atomB,
34583478 atomB.GetValence(lambda),
34593479 atomB.GetValence(sigma),
3460- static_cast<CartesianType>(dimA1));
3461- }
3462- diatomicTwoElecTwoCore[mu][nu][lambda][sigma]
3463- = this->GetNddoRepulsionIntegral(
3480+ static_cast<CartesianType>(dimA1),
3481+ static_cast<CartesianType>(dimA2));
3482+ }
3483+ tmpDiatomicTwoElecTwoCore1stDerivs[mu][nu][lambda][sigma][dimA1]
3484+ = this->GetNddoRepulsionIntegral1stDerivative(
34643485 atomA,
34653486 atomA.GetValence(mu),
34663487 atomA.GetValence(nu),
34673488 atomB,
34683489 atomB.GetValence(lambda),
3469- atomB.GetValence(sigma));
3470- }
3490+ atomB.GetValence(sigma),
3491+ static_cast<CartesianType>(dimA1));
3492+ }
3493+ tmpDiatomicTwoElecTwoCore[mu][nu][lambda][sigma]
3494+ = this->GetNddoRepulsionIntegral(
3495+ atomA,
3496+ atomA.GetValence(mu),
3497+ atomA.GetValence(nu),
3498+ atomB,
3499+ atomB.GetValence(lambda),
3500+ atomB.GetValence(sigma));
34713501 }
34723502 }
34733503 }
3474-
3475- // rotate matirix into the space frame
3476- this->CalcRotatingMatrix(rotatingMatrix, atomA, atomB);
3477- this->CalcRotatingMatrix1stDerivatives(rotMat1stDerivatives, atomA, atomB);
3478- this->CalcRotatingMatrix2ndDerivatives(rotMat2ndDerivatives, atomA, atomB);
3479- this->RotateDiatomicTwoElecTwoCore2ndDerivativesToSpaceFrame(matrix,
3480- diatomicTwoElecTwoCore,
3481- diatomicTwoElecTwoCore1stDerivatives,
3482- rotatingMatrix,
3483- rotMat1stDerivatives,
3484- rotMat2ndDerivatives);
3485- }
3486- catch(MolDSException ex){
3487- this->FreeDiatomicTwoElecTwoCore2ndDeriTemps(&rotatingMatrix,
3488- &rotMat1stDerivatives,
3489- &rotMat2ndDerivatives,
3490- &diatomicTwoElecTwoCore,
3491- &diatomicTwoElecTwoCore1stDerivatives);
3492- throw ex;
34933504 }
3494- this->FreeDiatomicTwoElecTwoCore2ndDeriTemps(&rotatingMatrix,
3495- &rotMat1stDerivatives,
3496- &rotMat2ndDerivatives,
3497- &diatomicTwoElecTwoCore,
3498- &diatomicTwoElecTwoCore1stDerivatives);
3499-}
3500-
3501-void Mndo::MallocDiatomicTwoElecTwoCore1stDeriTemps(double*** rotatingMatrix,
3502- double**** rotMat1stDerivatives,
3503- double***** diatomicTwoElecTwoCore) const{
3504- MallocerFreer::GetInstance()->Malloc<double>(rotatingMatrix,
3505- OrbitalType_end, OrbitalType_end);
3506- MallocerFreer::GetInstance()->Malloc<double>(rotMat1stDerivatives,
3507- OrbitalType_end,
3508- OrbitalType_end,
3509- CartesianType_end);
3510- MallocerFreer::GetInstance()->Malloc<double>(diatomicTwoElecTwoCore, dxy, dxy, dxy, dxy);
3511-}
3512-
3513-void Mndo::MallocDiatomicTwoElecTwoCore2ndDeriTemps(double*** rotatingMatrix,
3514- double**** rotMat1stDerivatives,
3515- double***** rotMat2ndDerivatives,
3516- double***** diatomicTwoElecTwoCore,
3517- double****** diatomicTwoElecTwoCore1stDerivatives) const{
3518- this->MallocDiatomicTwoElecTwoCore1stDeriTemps(rotatingMatrix,
3519- rotMat1stDerivatives,
3520- diatomicTwoElecTwoCore);
3521- MallocerFreer::GetInstance()->Malloc<double>(rotMat2ndDerivatives,
3522- OrbitalType_end,
3523- OrbitalType_end,
3524- CartesianType_end,
3525- CartesianType_end);
3526- MallocerFreer::GetInstance()->Malloc<double>(diatomicTwoElecTwoCore1stDerivatives,
3527- dxy,
3528- dxy,
3529- dxy,
3530- dxy,
3531- CartesianType_end);
3532-}
35333505
3534-void Mndo::FreeDiatomicTwoElecTwoCore1stDeriTemps(double*** rotatingMatrix,
3535- double**** rotMat1stDerivatives,
3536- double***** diatomicTwoElecTwoCore) const{
3537- MallocerFreer::GetInstance()->Free<double>(rotatingMatrix,
3538- OrbitalType_end, OrbitalType_end);
3539- MallocerFreer::GetInstance()->Free<double>(rotMat1stDerivatives,
3540- OrbitalType_end,
3541- OrbitalType_end,
3542- CartesianType_end);
3543- MallocerFreer::GetInstance()->Free<double>(diatomicTwoElecTwoCore, dxy, dxy, dxy, dxy);
3544-}
3545-
3546-void Mndo::FreeDiatomicTwoElecTwoCore2ndDeriTemps(double*** rotatingMatrix,
3547- double**** rotMat1stDerivatives,
3548- double***** rotMat2ndDerivatives,
3549- double***** diatomicTwoElecTwoCore,
3550- double****** diatomicTwoElecTwoCore1stDerivatives) const{
3551- this->FreeDiatomicTwoElecTwoCore1stDeriTemps(rotatingMatrix,
3552- rotMat1stDerivatives,
3553- diatomicTwoElecTwoCore);
3554- MallocerFreer::GetInstance()->Free<double>(rotMat2ndDerivatives,
3555- OrbitalType_end,
3556- OrbitalType_end,
3557- CartesianType_end,
3558- CartesianType_end);
3559- MallocerFreer::GetInstance()->Free<double>(diatomicTwoElecTwoCore1stDerivatives,
3560- dxy,
3561- dxy,
3562- dxy,
3563- dxy,
3564- CartesianType_end);
3506+ // rotate matirix into the space frame
3507+ this->CalcRotatingMatrix(tmpRotMat, atomA, atomB);
3508+ this->CalcRotatingMatrix1stDerivatives(tmpRotMat1stDerivs, atomA, atomB);
3509+ this->CalcRotatingMatrix2ndDerivatives(tmpRotMat2ndDerivs, atomA, atomB);
3510+ this->RotateDiatomicTwoElecTwoCore2ndDerivativesToSpaceFrame(matrix,
3511+ tmpDiatomicTwoElecTwoCore,
3512+ tmpDiatomicTwoElecTwoCore1stDerivs,
3513+ tmpRotMat,
3514+ tmpRotMat1stDerivs,
3515+ tmpRotMat2ndDerivs);
35653516 }
35663517
35673518 // Rotate 4-dimensional matrix from diatomic frame to space frame
--- a/src/mndo/Mndo.h
+++ b/src/mndo/Mndo.h
@@ -153,14 +153,18 @@ private:
153153 double******** diatomicTwoElecTwoCore2ndDerivs,
154154 double*** tmpRotMat,
155155 double**** tmpRotMat1stDerivs,
156- double***** tmpDiatomicTwoElecTwo) const;
156+ double***** tmpRotMat2ndDerivs,
157+ double***** tmpDiatomicTwoElecTwo,
158+ double****** tmpDiatomicTwoElecTwo1stDerivs) const;
157159 void FreeTempMatricesEachThreadCalcHessianSCF(double***** diatomicOverlapAOs1stDerivs,
158160 double****** diatomicOverlapAOs2ndDerivs,
159161 double******* diatomicTwoElecTwoCore1stDerivs,
160162 double******** diatomicTwoElecTwoCore2ndDerivs,
161163 double*** tmpRotMat,
162164 double**** tmpRotMat1stDerivs,
163- double***** tmpDiatomicTwoElecTwo) const;
165+ double***** tmpRotMat2ndDerivs,
166+ double***** tmpDiatomicTwoElecTwo,
167+ double****** tmpDiatomicTwoElecTwo1stDerivs) const;
164168 double GetAuxiliaryHessianElement1(int mu,
165169 int nu,
166170 int indexAtomA,
@@ -284,24 +288,13 @@ private:
284288 int indexAtomA,
285289 int indexAtomB) const;
286290 void CalcDiatomicTwoElecTwoCore2ndDerivatives(double****** matrix,
291+ double** tmpRotMat,
292+ double*** tmpRotMat1stDerivs,
293+ double**** tmpRotMat2ndDerivs,
294+ double**** tmpDiatomicTwoElecTwoCore,
295+ double***** tmpDiatomicTwoElecTwoCore1stDerivs,
287296 int indexAtomA,
288297 int indexAtomB) const;
289- void MallocDiatomicTwoElecTwoCore1stDeriTemps(double*** rotatingMatrix,
290- double**** rotMat1stDerivatives,
291- double***** diatomicTwoElecTwoCore) const;
292- void MallocDiatomicTwoElecTwoCore2ndDeriTemps(double*** rotatingMatrix,
293- double**** rotMat1stDerivatives,
294- double***** rotMat2ndDerivatives,
295- double***** diatomicTwoElecTwoCore,
296- double****** diatomicTwoElecTwoCore1stDerivatives) const;
297- void FreeDiatomicTwoElecTwoCore1stDeriTemps(double*** rotatingMatrix,
298- double**** rotMat1stDerivatives,
299- double***** diatomicTwoElecTwoCore) const;
300- void FreeDiatomicTwoElecTwoCore2ndDeriTemps(double*** rotatingMatrix,
301- double**** rotMat1stDerivatives,
302- double***** rotMat2ndDerivatives,
303- double***** diatomicTwoElecTwoCore,
304- double****** diatomicTwoElecTwoCore1stDerivatives) const;
305298 void RotateDiatomicTwoElecTwoCoreToSpaceFrame(double**** matrix,
306299 double const* const* rotatingMatrix) const;
307300 void RotateDiatomicTwoElecTwoCore1stDerivativesToSpaceFrame(double***** matrix,