• 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

Revisión89b871e47633cbd2df42ccb670abef94c37c813f (tree)
Tiempo2012-12-05 16:55:13
AutorKatsuhiko Nishimra <ktns.87@gmai...>
CommiterKatsuhiko Nishimra

Log Message

Reimplement BFGS::ShiftHessian using BLAS functions. #28764

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

Cambiar Resumen

Diferencia incremental

--- a/src/optimization/BFGS.cpp
+++ b/src/optimization/BFGS.cpp
@@ -407,16 +407,20 @@ void BFGS::UpdateHessian(double **matrixHessian,
407407
408408 void BFGS::ShiftHessianRedundantMode(double** matrixHessian,
409409 const Molecule& molecule) const{
410- const double largeEigenvalue = 1.0e3;
411- const int numAtoms = molecule.GetNumberAtoms();
412- const int dimension = numAtoms *CartesianType_end;
413- const int numTranslationalModes = 3;
414- const int numRotationalModes = 3;
415- int numRedundantModes = numTranslationalModes + numRotationalModes;
416- double** vectorsHessianModes = NULL;
417- double* vectorHessianEigenValues = NULL;
418- double** matrixesRedundantModes[] = {NULL, NULL, NULL, NULL, NULL, NULL};
419- double* vectorsRedundantModes[] = {NULL, NULL, NULL, NULL, NULL, NULL};
410+ const double one = 1;
411+ const double largeEigenvalue = 1.0e3;
412+ const int numAtoms = molecule.GetNumberAtoms();
413+ const int dimension = numAtoms *CartesianType_end;
414+ const int numTranslationalModes = 3;
415+ const int numRotationalModes = 3;
416+ int numRedundantModes = numTranslationalModes + numRotationalModes;
417+ double** vectorsHessianModes = NULL;
418+ double* vectorHessianEigenValues = NULL;
419+ double** matrixesRedundantModes[] = {NULL, NULL, NULL, NULL, NULL, NULL};
420+ double* vectorsRedundantModes[] = {NULL, NULL, NULL, NULL, NULL, NULL};
421+ double** matrixProjection = NULL;
422+ double* vectorProjectedRedundantMode = NULL;
423+ double** matrixShiftedHessianBuffer = NULL;
420424 const double matrixesRotationalModeGenerators[numRotationalModes]
421425 [CartesianType_end]
422426 [CartesianType_end]
@@ -455,42 +459,17 @@ void BFGS::ShiftHessianRedundantMode(double** matrixHessian,
455459 }
456460 }
457461 }
458- // Orthonormalize redundant modes
459- // using Gram-Schmidt process
460- for(int c=0; c<numRedundantModes;c++){
461- for(int d=0;d<c;d++){
462- double dotproduct = 0.0;
463-#pragma omp parallel for schedule(auto) reduction(+:dotproduct)
464- for(int i=0;i<dimension;i++){
465- dotproduct += vectorsRedundantModes[d][i] * vectorsRedundantModes[c][i];
466- }
467-#pragma omp parallel for schedule(auto)
468- for(int i=0;i<dimension;i++){
469- vectorsRedundantModes[c][i] -= dotproduct * vectorsRedundantModes[d][i];
470- }
471- }
472- double norm = 0.0;
473-#pragma omp parallel for schedule(auto) reduction(+:norm)
474- for(int i=0;i<dimension;i++){
475- norm += vectorsRedundantModes[c][i] * vectorsRedundantModes[c][i];
476- }
477- norm = sqrt(norm);
478- // Eliminate a linear dependent mode
479- if(norm < 1e-5){
480- numRedundantModes--;
481- for(int d=c;d<numRedundantModes;d++){
482-#pragma omp parallel for schedule(auto)
483- for(int i=0;i<dimension;i++){
484- vectorsRedundantModes[d][i] = vectorsRedundantModes[d+1][i];
485- }
486- }
487- }
488- else{
489-#pragma omp parallel for schedule(auto)
490- for(int i=0;i<dimension;i++){
491- vectorsRedundantModes[c][i] /= norm;
492- }
493- }
462+ //Prepare the identity matrix I
463+ MallocerFreer::GetInstance()->Malloc(&matrixProjection, dimension, dimension);
464+ MolDS_wrappers::Blas::GetInstance()->Dcopy(dimension, &one, 0, &matrixProjection[0][0], dimension+1);
465+
466+ //Prepare the projection matrix P = I - sum_i u_i u_i^T, that projects redundant modes to 0 vector
467+ MallocerFreer::GetInstance()->Malloc(&vectorProjectedRedundantMode, dimension);
468+ for(int c=0; c<numRedundantModes; c++){
469+ double normSquare = 0;
470+ MolDS_wrappers::Blas::GetInstance()->Dsymv(dimension, matrixProjection, vectorsRedundantModes[c], vectorProjectedRedundantMode);
471+ normSquare = MolDS_wrappers::Blas::GetInstance()->Ddot(dimension, vectorProjectedRedundantMode, vectorProjectedRedundantMode);
472+ MolDS_wrappers::Blas::GetInstance()->Dsyr(dimension, -1.0/normSquare, vectorProjectedRedundantMode, matrixProjection);
494473 }
495474
496475 // Diagonalize hessian
@@ -503,7 +482,7 @@ void BFGS::ShiftHessianRedundantMode(double** matrixHessian,
503482 dimension,
504483 calcEigenVectors);
505484
506- // Output eigenvalues of the raw Hessianto the log
485+ // Output eigenvalues of the raw Hessian to the log
507486 this->OutputLog(this->messageRawHessianEigenvalues);
508487 for(int i=0;i<dimension;i++){
509488 if((i%6) == 0){
@@ -515,46 +494,27 @@ void BFGS::ShiftHessianRedundantMode(double** matrixHessian,
515494 }
516495 this->OutputLog("\n");
517496
518- // Orthogonalize vectorsHessianModes against vectorsRedundantModes
519-#pragma omp parallel for schedule(auto)
520- for(int i=0;i<dimension-numRedundantModes;i++){
521- for(int c=0; c<numRedundantModes;c++){
522- double dotproduct = 0.0;
523- for(int j=0;j<dimension;j++){
524- dotproduct += vectorsHessianModes[i][j] * vectorsRedundantModes[c][j];
525- }
526- for(int j=0;j<dimension;j++){
527- vectorsHessianModes[i][j] -= dotproduct * vectorsRedundantModes[c][j];
528- }
529- double norm = 0.0;
530- for(int j=0;j<dimension;j++){
531- norm += vectorsHessianModes[i][j] * vectorsHessianModes[i][j];
532- }
533- norm = sqrt(norm);
534- for(int j=0;j<dimension;j++){
535- vectorsHessianModes[i][j] /= norm;
536- }
537- }
538- }
497+ // Project Hessian H' = P H P
498+ MallocerFreer::GetInstance()->Malloc(&matrixShiftedHessianBuffer, dimension, dimension);
499+ // TODO: Use dsymm instead of dgemm
500+ MolDS_wrappers::Blas::GetInstance()->Dgemm(dimension, dimension, dimension,
501+ matrixProjection , matrixHessian , matrixShiftedHessianBuffer);
502+ MolDS_wrappers::Blas::GetInstance()->Dgemm(dimension, dimension, dimension,
503+ matrixShiftedHessianBuffer, matrixProjection, matrixHessian);
539504
540- // Calculate projected eigenvalues of new modes
541- // h_i' = x_i'^T * H * x_i'
542-#pragma omp parallel for schedule(auto)
543- for(int i=0;i<dimension-numRedundantModes;i++){
544- double tmp = 0.0;
545- for(int j=0;j<dimension;j++){
546- for(int k=0;k<dimension;k++){
547- tmp += vectorsHessianModes[i][j] * matrixHessian[j][k] * vectorsHessianModes[i][k];
548- }
549- }
550- vectorHessianEigenValues[i] = tmp;
551- }
552-#pragma omp parallel for schedule(auto)
553- for(int i=dimension-numRedundantModes;i<dimension;i++){
554- vectorHessianEigenValues[i] = largeEigenvalue;
555- }
505+ // Shift eigenvalues for redundant mode by adding L sum_i u_i*u_i^T = L ( I - P )= -L (P - I)
506+ MolDS_wrappers::Blas::GetInstance()->Daxpy(dimension, -1, &one, 0, &matrixProjection[0][0], dimension + 1); // (P - I)
507+ MolDS_wrappers::Blas::GetInstance()->Daxpy(dimension*dimension, -largeEigenvalue, &matrixProjection[0][0], &matrixHessian[0][0]);
508+
509+ // Diagonalize shifted hessian
510+ MolDS_wrappers::Blas::GetInstance()->Dcopy(dimension*dimension, &matrixHessian[0][0], &vectorsHessianModes[0][0]);
511+ calcEigenVectors = true;
512+ MolDS_wrappers::Lapack::GetInstance()->Dsyevd(&vectorsHessianModes[0],
513+ &vectorHessianEigenValues[0],
514+ dimension,
515+ calcEigenVectors);
556516
557- // Output eigenvalues of the raw Hessianto the log
517+ // Output eigenvalues of the shifted Hessian to the log
558518 this->OutputLog(this->messageShiftedHessianEigenvalues);
559519 for(int i=0;i<dimension;i++){
560520 if((i%6) == 0){
@@ -565,34 +525,24 @@ void BFGS::ShiftHessianRedundantMode(double** matrixHessian,
565525 }
566526 }
567527 this->OutputLog("\n");
568-
569- // Calculate shifted Hessian from eigenvalues and modes
570- // H' = sum x_i' h_i' x_i^T
571-#pragma omp parallel for schedule(auto)
572- for(int i=0;i<dimension;i++){
573- for(int j=0;j<dimension;j++){
574- double tmp = 0.0;
575- for(int k=0;k<dimension-numRedundantModes;k++){
576- tmp += vectorsHessianModes[k][i] * vectorsHessianModes[k][j] * vectorHessianEigenValues[k];
577- }
578- for(int k=0;k<numRedundantModes;k++){
579- tmp += vectorsRedundantModes[k][i] * vectorsRedundantModes[k][j] * largeEigenvalue;
580- }
581- matrixHessian[i][j] = tmp;
582- }
583- }
584528 }
585529 catch(MolDSException ex)
586530 {
587531 for(int i=0;i<numRedundantModes;i++){
588532 MallocerFreer::GetInstance()->Free(&matrixesRedundantModes[i], numAtoms, CartesianType_end);
589533 }
534+ MallocerFreer::GetInstance()->Free(&matrixProjection, dimension, dimension);
535+ MallocerFreer::GetInstance()->Free(&vectorProjectedRedundantMode, dimension);
536+ MallocerFreer::GetInstance()->Free(&matrixShiftedHessianBuffer, dimension, dimension);
590537 MallocerFreer::GetInstance()->Free(&vectorHessianEigenValues, dimension);
591538 MallocerFreer::GetInstance()->Free(&vectorsHessianModes, dimension, dimension);
592539 }
593540 for(int i=0;i<numRedundantModes;i++){
594541 MallocerFreer::GetInstance()->Free(&matrixesRedundantModes[i], numAtoms, CartesianType_end);
595542 }
543+ MallocerFreer::GetInstance()->Free(&matrixProjection, dimension, dimension);
544+ MallocerFreer::GetInstance()->Free(&vectorProjectedRedundantMode, dimension);
545+ MallocerFreer::GetInstance()->Free(&matrixShiftedHessianBuffer, dimension, dimension);
596546 MallocerFreer::GetInstance()->Free(&vectorHessianEigenValues, dimension);
597547 MallocerFreer::GetInstance()->Free(&vectorsHessianModes, dimension, dimension);
598548 }