File: FT_VNA.c

package info (click to toggle)
openmx 3.7.6-1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd, stretch
  • size: 325,856 kB
  • ctags: 3,575
  • sloc: ansic: 152,655; f90: 2,080; python: 876; makefile: 675; sh: 25; perl: 18
file content (259 lines) | stat: -rw-r--r-- 6,299 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
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
/**********************************************************************
  FT_VNA.c:

     FT_VNA.c is a subroutine to Fourier transform VNA potentials

  Log of FT_VNA.c:

     18/May/2004  Released by T.Ozaki

***********************************************************************/

#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <time.h>
/*  stat section */
#include <sys/types.h>
#include <sys/stat.h>
#include <unistd.h>
/*  end stat section */
#include "openmx_common.h"
#include "mpi.h"
#include <omp.h>





void FT_VNA()
{
  int numprocs,myid,ID,tag=999;
  int count,NumSpe;
  int L,i,j;
  int Lspe,spe,GL,Mul;
  int RestartRead_Succeed;
  double Sr,Dr,Sk,Dk,kmin,kmax;
  double norm_k,h,dum0;
  double xmin,xmax,x,r,sum;
  double sj;
  double tmp0,tmp1;
  double *tmp_SphB,*tmp_SphBp;
  double TStime, TEtime;
  /* for MPI */
  MPI_Status stat;
  MPI_Request request;
  /* for OpenMP */
  int OMPID,Nthrds,Nprocs;

  char fileFT[YOUSO10];
  char operate[300];
  FILE *fp;
  size_t size; 

  dtime(&TStime);

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

  if (myid==Host_ID && 0<level_stdout) printf("<FT_VNA>          Fourier transform of VNA potentials\n");

  RestartRead_Succeed = 0;

  /***********************************************************
   In case of Scf_RestartFromFile==1, read Spe_CrudeVNA_Bessel
  ***********************************************************/

  if (Scf_RestartFromFile){

    /****************************************************
         regenerate radial grids in the k-space
         for the MPI calculation
    ****************************************************/

    for (j=0; j<GL_Mesh; j++){
      kmin = Radial_kmin;
      kmax = PAO_Nkmax;
      Sk = kmax + kmin;
      Dk = kmax - kmin;
      norm_k = 0.50*(Dk*GL_Abscissae[j] + Sk);
      GL_NormK[j] = norm_k;
    }

    /***********************************************************
                        read Spe_VNA_Bessel
    ***********************************************************/

    sprintf(fileFT,"%s%s_rst/%s.ftCvna",filepath,filename,filename);

    if ((fp = fopen(fileFT,"rb")) != NULL){

      RestartRead_Succeed = 1;

      for (spe=0; spe<SpeciesNum; spe++){
        for (L=0; L<=List_YOUSO[35]; L++){
	  for (Mul=0; Mul<List_YOUSO[34]; Mul++){

	    size = fread(&Spe_CrudeVNA_Bessel[spe][0],sizeof(double),GL_Mesh,fp);
	    if (size!=GL_Mesh) RestartRead_Succeed = 0;
	  }
        }
      }

      fclose(fp);
    }
    else{
      printf("Could not open a file %s in FT_VNA\n",fileFT);
    }
  }

  /***********************************************************
   if (RestartRead_Succeed==0), calculate Spe_CrudeVNA_Bessel
  ***********************************************************/

  if (RestartRead_Succeed==0){

    /* loop for Lspe */

    for (Lspe=0; Lspe<MSpeciesNum; Lspe++){

      spe = Species_Top[myid] + Lspe;

      /****************************************************
                   \int jL(k*r)RL r^2 dr 
      ****************************************************/

      /* tabulation on Gauss-Legendre radial grid */

      kmin = Radial_kmin;
      kmax = PAO_Nkmax;
      Sk = kmax + kmin;
      Dk = kmax - kmin;
      xmin = sqrt(Spe_PAO_RV[spe][0]);
      xmax = sqrt(Spe_Atom_Cut1[spe] + 0.5);
      h = (xmax - xmin)/(double)OneD_Grid;

      /* loop for j */

#pragma omp parallel shared(spe,Dk,Sk,GL_Abscissae,xmin,xmax,h,OneD_Grid,Spe_CrudeVNA_Bessel)  private(OMPID,Nthrds,Nprocs,j,norm_k,sum,x,r,i,tmp_SphB,tmp_SphBp,sj)
      {

	/* allocate arrays */

	tmp_SphB  = (double*)malloc(sizeof(double)*3);
	tmp_SphBp = (double*)malloc(sizeof(double)*3);

	/* get info. on OpenMP */ 

	OMPID = omp_get_thread_num();
	Nthrds = omp_get_num_threads();
	Nprocs = omp_get_num_procs();

	for ( j=OMPID; j<GL_Mesh; j+=Nthrds ){

	  norm_k = 0.50*(Dk*GL_Abscissae[j] + Sk);

	  /**************************
           trapezoidal rule

            grid: r = x^2
                  dr = 2*x*dx 
	  ***************************/

	  sum = 0.0;

	  for (i=0; i<=OneD_Grid; i++){
	    x = xmin + (double)i*h;
	    r = x*x; 

	    Spherical_Bessel(norm_k*r,0,tmp_SphB,tmp_SphBp);
	    sj = tmp_SphB[0];

	    if (i==0 || i==OneD_Grid)
	      sum += r*r*x*sj*VNAF(spe,r);
	    else 
	      sum += 2.0*r*r*x*sj*VNAF(spe,r);
	  }
	  sum = sum*h;

	  Spe_CrudeVNA_Bessel[spe][j] = sum;

	} /* j */

	/* free arrays */

	free(tmp_SphB);
	free(tmp_SphBp);

#pragma omp flush(Spe_CrudeVNA_Bessel)

      } /* #pragma omp parallel */
    } /* Lspe */

    /****************************************************
         regenerate radial grids in the k-space
         for the MPI calculation
    ****************************************************/

    for (j=0; j<GL_Mesh; j++){
      kmin = Radial_kmin;
      kmax = PAO_Nkmax;
      Sk = kmax + kmin;
      Dk = kmax - kmin;
      norm_k = 0.50*(Dk*GL_Abscissae[j] + Sk);
      GL_NormK[j] = norm_k;
    }

    /***********************************************************
      sending and receiving of Spe_CrudeVNA_Bessel by MPI
    ***********************************************************/

    for (ID=0; ID<Num_Procs2; ID++){
      NumSpe = Species_End[ID] - Species_Top[ID] + 1;
      for (Lspe=0; Lspe<NumSpe; Lspe++){
	spe = Species_Top[ID] + Lspe;
	MPI_Bcast(&Spe_CrudeVNA_Bessel[spe][0],
		  GL_Mesh,MPI_DOUBLE,ID,mpi_comm_level1);
      }
    }

    /***********************************************************
                      save Spe_CrudeVNA_Bessel
    ***********************************************************/

    if (myid==Host_ID){

      sprintf(fileFT,"%s%s_rst/%s.ftCvna",filepath,filename,filename);

      if ((fp = fopen(fileFT,"wb")) != NULL){

	for (spe=0; spe<SpeciesNum; spe++){
	  fwrite(&Spe_CrudeVNA_Bessel[spe][0],sizeof(double),GL_Mesh,fp);
	}

	fclose(fp);
      }
      else{
	printf("Could not open a file %s in FT_VNA\n",fileFT);
      }
    }

  } /* if (RestartRead_Succeed==0) */

  /***********************************************************
                         elapsed time
  ***********************************************************/

  dtime(&TEtime);

  /*
  printf("myid=%2d Elapsed Time (s) = %15.12f\n",myid,TEtime-TStime);
  MPI_Finalize();
  exit(0);
  */
}