File: FT_VNA.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 (149 lines) | stat: -rw-r--r-- 3,431 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
/**********************************************************************
  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>
#include "openmx_common.h"

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

void FT_VNA()
{
  int numprocs,myid,ID,tag=999;
  int count,NumSpe;
  int L,i,j;
  int Lspe,spe,GL,Mul;
  double Sr,Dr,Sk,Dk,kmin,kmax;
  double norm_k,h,dum0;
  double xmin,xmax,x,r,sum;
  double sj,sy,sjp,syp;
  double RGL[GL_Mesh + 2];
  double SumTmp;
  double tmp0,tmp1;
  double *tmp_SphB,*tmp_SphBp;
  double TStime, TEtime;
  MPI_Status stat;
  MPI_Request request;

  dtime(&TStime);

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

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

  /* allocate arrays */

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

  /* 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;

    for (j=0; j<GL_Mesh; j++){

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

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

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

      xmin = sqrt(Spe_PAO_RV[spe][0]);
      xmax = sqrt(Spe_Atom_Cut1[spe] + 0.5);

      sum = 0.0;
      h = (xmax - xmin)/(double)OneD_Grid;

      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;
      GL_NormK[j] = norm_k;
    }
  }

  /****************************************************
         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);
    }
  }

  /* free arrays */

  free(tmp_SphB);
  free(tmp_SphBp);

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

  dtime(&TEtime);

}