File: Mixing_DM.c

package info (click to toggle)
openmx 3.2.4.dfsg-3
  • links: PTS, VCS
  • area: main
  • in suites: lenny, squeeze
  • size: 62,572 kB
  • ctags: 2,684
  • sloc: ansic: 130,666; python: 876; makefile: 560; xml: 63; perl: 18; sh: 4
file content (183 lines) | stat: -rw-r--r-- 5,851 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <time.h>
#include "openmx_common.h"

#ifdef nompi
#include "mimic_mpi.h"
#else
#include "mpi.h"
#endif

double Mixing_DM(int MD_iter,
                 int SCF_iter,
                 int SCF_iter0,
                 int SucceedReadingDMfile,
                 double *****ReRhok,
                 double *****ImRhok,
                 double ****ReBestRhok,
                 double ****ImBestRhok,
                 double *****Residual_ReRhok,
                 double *****Residual_ImRhok,
                 double ***ReV1,
                 double ***ImV1,
                 double ***ReV2,
                 double ***ImV2,
                 double ***ReRhoAtomk,
                 double ***ImRhoAtomk)
{
  int pSCF_iter,NumMix,NumSlide;
  int spin,ct_AN,wan1,TNO1,h_AN,Gh_AN,wan2;
  int TNO2,i,j,ian,jan,m,n;
  int spinmax,k1,k2,k3;
  double time0;
  double TStime,TEtime;
  int numprocs,myid,ID;

  /* MPI */
  MPI_Comm_size(mpi_comm_level1,&numprocs);
  MPI_Comm_rank(mpi_comm_level1,&myid);

  MPI_Barrier(mpi_comm_level1);
  dtime(&TStime);

  /*******************************************************
                     Simple Mixing
  *******************************************************/

  if (Mixing_switch==0){
    if (MD_iter==1 && SCF_iter==1){
      Simple_Mixing_DM(0,1.00,DM[0],DM[1],DM[2],iDM[0],iDM[1],iDM[2],ResidualDM[2],iResidualDM[2]);
    }
    else{
      Simple_Mixing_DM(1,Mixing_weight,DM[0],DM[1],DM[2],iDM[0],iDM[1],iDM[2],ResidualDM[2],iResidualDM[2]);
    }
  }

  /*******************************************************
    Pulay's method

    Residual Minimazation Method (RMM) using
    Direct Inversion in the Iterative Subspace (DIIS)
  *******************************************************/

  else if (Mixing_switch==1){
    if (SCF_iter==1){
      DIIS_Mixing_DM(1,ResidualDM,iResidualDM);
    }
    else if (SCF_iter<=(Pulay_SCF-1)){
      Simple_Mixing_DM(1,Mixing_weight,DM[0],DM[1],DM[2],iDM[0],iDM[1],iDM[2],ResidualDM[2],iResidualDM[2]);
    }
    else if (SCF_iter==Pulay_SCF)
      DIIS_Mixing_DM(2,ResidualDM,iResidualDM);
    else 
      DIIS_Mixing_DM(SCF_iter-(Pulay_SCF-2),ResidualDM,iResidualDM);
  }

  /*********************************************************
     Guaranteed-Reduction Pulay's method (GR-Pulay method)
  *********************************************************/

  else if (Mixing_switch==2){
    if (SCF_iter==1){
      GR_Pulay_DM(1,ResidualDM);
    }
    else if (SCF_iter<=(Pulay_SCF-1))
      Simple_Mixing_DM(1,Mixing_weight,DM[0],DM[1],DM[2],iDM[0],iDM[1],iDM[2],ResidualDM[2],iResidualDM[2]);
    else if (SCF_iter==Pulay_SCF)
      GR_Pulay_DM(2,ResidualDM);
    else 
      GR_Pulay_DM(SCF_iter-(Pulay_SCF-2),ResidualDM);
  }

  /*********************************************************
               Kerker simple mixing in k-space
  *********************************************************/

  else if (Mixing_switch==3){

    if (MD_iter==1 && SCF_iter0==1 && SucceedReadingDMfile==0){
      Kerker_Mixing_Rhok(0, Mixing_weight,
                         ReRhok,ImRhok,ReBestRhok,ImBestRhok,
                         Residual_ReRhok,Residual_ImRhok,
                         ReV1,ImV1,ReV2,ImV2,ReRhoAtomk,ImRhoAtomk);
    }

    else if (MD_iter==1 && SCF_iter==1){
      Kerker_Mixing_Rhok(0, 1.0,
                         ReRhok,ImRhok,ReBestRhok,ImBestRhok,
                         Residual_ReRhok,Residual_ImRhok,
                         ReV1,ImV1,ReV2,ImV2,ReRhoAtomk,ImRhoAtomk);
    }

    else if ( (SCF_iter%30)==1) {
      Kerker_Mixing_Rhok(1,Max_Mixing_weight,
                         ReRhok,ImRhok,ReBestRhok,ImBestRhok,
                         Residual_ReRhok,Residual_ImRhok,
                         ReV1,ImV1,ReV2,ImV2,ReRhoAtomk,ImRhoAtomk);
    }

    else{
      Kerker_Mixing_Rhok(1,Mixing_weight,
                         ReRhok,ImRhok,ReBestRhok,ImBestRhok,
                         Residual_ReRhok,Residual_ImRhok,
                         ReV1,ImV1,ReV2,ImV2,ReRhoAtomk,ImRhoAtomk);
    }
  }

  /*********************************************************
   RMM-DIIS mixing with Kerker's weighting factor in k-space
  *********************************************************/

  else if (Mixing_switch==4){

    if (MD_iter==1 && SCF_iter0==1 && SucceedReadingDMfile==0){
      Kerker_Mixing_Rhok(0, Mixing_weight,
                         ReRhok,ImRhok,ReBestRhok,ImBestRhok,
                         Residual_ReRhok,Residual_ImRhok,
                         ReV1,ImV1,ReV2,ImV2,ReRhoAtomk,ImRhoAtomk);
    }

    else if (MD_iter==1 && SCF_iter==1){
      Kerker_Mixing_Rhok(0,1.0,
                         ReRhok,ImRhok,ReBestRhok,ImBestRhok,
                         Residual_ReRhok,Residual_ImRhok,
                         ReV1,ImV1,ReV2,ImV2,ReRhoAtomk,ImRhoAtomk);
    }

    else if (SCF_iter<Pulay_SCF){
      Kerker_Mixing_Rhok(1,Mixing_weight,
                         ReRhok,ImRhok,ReBestRhok,ImBestRhok,
                         Residual_ReRhok,Residual_ImRhok,
                         ReV1,ImV1,ReV2,ImV2,ReRhoAtomk,ImRhoAtomk);
    }

    else{
      DIIS_Mixing_Rhok(SCF_iter-(Pulay_SCF-2),Mixing_weight,
                       ReRhok,ImRhok,ReBestRhok,ImBestRhok,
                       Residual_ReRhok,Residual_ImRhok,
                       ReV1,ImV1,ReV2,ImV2,ReRhoAtomk,ImRhoAtomk);
    }

  }

  /*******************************************************
                       not supported 
  *******************************************************/

  else {
    printf("Mixing_switch=%i is not supported.\n",Mixing_switch);
    MPI_Finalize();
    exit(0);
  }

  /* if SCF_iter0==1, then NormRD[0]=1 */

  if (SCF_iter0==1) NormRD[0] = 1.0;

  MPI_Barrier(mpi_comm_level1);
  dtime(&TEtime);
  time0 = TEtime - TStime;
  return time0;
}