CheMPS2
HeffDiagrams1.cpp
1 /*
2  CheMPS2: a spin-adapted implementation of DMRG for ab initio quantum chemistry
3  Copyright (C) 2013-2016 Sebastian Wouters
4 
5  This program is free software; you can redistribute it and/or modify
6  it under the terms of the GNU General Public License as published by
7  the Free Software Foundation; either version 2 of the License, or
8  (at your option) any later version.
9 
10  This program is distributed in the hope that it will be useful,
11  but WITHOUT ANY WARRANTY; without even the implied warranty of
12  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  GNU General Public License for more details.
14 
15  You should have received a copy of the GNU General Public License along
16  with this program; if not, write to the Free Software Foundation, Inc.,
17  51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
18 */
19 
20 #include <math.h>
21 
22 #include "Heff.h"
23 #include "Lapack.h"
24 #include "MPIchemps2.h"
25 
26 void CheMPS2::Heff::addDiagram1A(const int ikappa, double * memS, double * memHeff, const Sobject * denS, TensorX * Xleft) const{
27  int dimL = denBK->gCurrentDim(denS->gIndex(), denS->gNL(ikappa), denS->gTwoSL(ikappa), denS->gIL(ikappa));
28  int dimR = denBK->gCurrentDim(denS->gIndex()+2, denS->gNR(ikappa), denS->gTwoSR(ikappa), denS->gIR(ikappa));
29  double * BlockX = Xleft->gStorage( denS->gNL(ikappa), denS->gTwoSL(ikappa), denS->gIL(ikappa), denS->gNL(ikappa), denS->gTwoSL(ikappa), denS->gIL(ikappa) );
30 
31  double one = 1.0;
32  char notr = 'N';
33  dgemm_(&notr,&notr,&dimL,&dimR,&dimL,&one,BlockX,&dimL,memS+denS->gKappa2index(ikappa),&dimL,&one,memHeff+denS->gKappa2index(ikappa),&dimL);
34 }
35 
36 void CheMPS2::Heff::addDiagram1B(const int ikappa, double * memS, double * memHeff, const Sobject * denS, TensorX * Xright) const{
37  int dimL = denBK->gCurrentDim(denS->gIndex(), denS->gNL(ikappa), denS->gTwoSL(ikappa), denS->gIL(ikappa));
38  int dimR = denBK->gCurrentDim(denS->gIndex()+2, denS->gNR(ikappa), denS->gTwoSR(ikappa), denS->gIR(ikappa));
39  double * BlockX = Xright->gStorage( denS->gNR(ikappa), denS->gTwoSR(ikappa), denS->gIR(ikappa), denS->gNR(ikappa), denS->gTwoSR(ikappa), denS->gIR(ikappa) );
40 
41  double one = 1.0;
42  char notr = 'N';
43  char trans = 'T';
44  dgemm_(&notr,&trans,&dimL,&dimR,&dimR,&one,memS+denS->gKappa2index(ikappa),&dimL,BlockX,&dimR,&one,memHeff+denS->gKappa2index(ikappa),&dimL);
45 }
46 
47 void CheMPS2::Heff::addDiagram1C(const int ikappa, double * memS, double * memHeff, const Sobject * denS, double Helem_links) const{
48  if (denS->gN1(ikappa)==2){
49  int inc = 1;
50  int ptr = denS->gKappa2index(ikappa);
51  int dim = denS->gKappa2index(ikappa+1) - ptr;
52  daxpy_(&dim,&Helem_links,memS+ptr,&inc,memHeff+ptr,&inc);
53  }
54 }
55 
56 void CheMPS2::Heff::addDiagram1D(const int ikappa, double * memS, double * memHeff, const Sobject * denS, double Helem_rechts) const{
57  if (denS->gN2(ikappa)==2){
58  int inc = 1;
59  int ptr = denS->gKappa2index(ikappa);
60  int dim = denS->gKappa2index(ikappa+1) - ptr;
61  daxpy_(&dim,&Helem_rechts,memS+ptr,&inc,memHeff+ptr,&inc);
62  }
63 }
64 
65 void CheMPS2::Heff::addDiagramExcitations(const int ikappa, double * memS, double * memHeff, const Sobject * denS, int nLower, double ** VeffTilde) const{
66 
67  int dimTotal = denS->gKappa2index(denS->gNKappa());
68  int ptr = denS->gKappa2index(ikappa);
69  int dimBlock = denS->gKappa2index(ikappa+1) - ptr;
70  int inc = 1;
71  #ifdef CHEMPS2_MPI_COMPILATION
72  const int MPIRANK = MPIchemps2::mpi_rank();
73  #endif
74 
75  for (int state=0; state<nLower; state++){
76  #ifdef CHEMPS2_MPI_COMPILATION
77  if ( MPIchemps2::owner_specific_excitation( Prob->gL(), state ) == MPIRANK )
78  #endif
79  {
80  double alpha = ddot_(&dimTotal, memS, &inc, VeffTilde[state], &inc);
81  daxpy_(&dimBlock,&alpha,VeffTilde[state]+ptr,&inc,memHeff+ptr,&inc);
82  }
83  }
84 
85 }
86 
int gCurrentDim(const int boundary, const int N, const int TwoS, const int irrep) const
Get the current virtual dimensions.
static int mpi_rank()
Get the rank of this MPI process.
Definition: MPIchemps2.h:131
int gL() const
Get the number of orbitals.
Definition: Problem.h:51