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
|
/**********************************************************************
FT_ProExpn_VNA.c:
FT_ProExpn_VNA.c is a subroutine to Fourier transform
VNA separable projectors.
Log of FT_ProExpn_VNA.c:
7/Apr/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
#ifdef noomp
#include "mimic_omp.h"
#else
#include <omp.h>
#endif
void FT_ProExpn_VNA()
{
int numprocs,myid,ID,tag=999;
int count,NumSpe;
int L,i,kj;
int Lspe,spe,GL,Mul;
double Sr,Dr;
double norm_k,h,dum0;
double rmin,rmax,r,sum;
double kmin,kmax,Sk,Dk;
double RGL[GL_Mesh + 2];
double *SumTmp;
double tmp0,tmp1;
double **SphB;
double *tmp_SphB,*tmp_SphBp;
double TStime, TEtime;
/* for MPI */
MPI_Status stat;
MPI_Request request;
/* for OpenMP */
int OMPID,Nthrds,Nprocs;
dtime(&TStime);
/* MPI */
MPI_Comm_size(mpi_comm_level1,&numprocs);
MPI_Comm_rank(mpi_comm_level1,&myid);
if (myid==Host_ID) printf("<FT_ProExpn_VNA> Fourier transform of VNA separable projectors\n");
for (Lspe=0; Lspe<MSpeciesNum; Lspe++){
spe = Species_Top[myid] + Lspe;
/* initalize */
/* tabulation on Gauss-Legendre radial grid */
rmin = Spe_VPS_RV[spe][0];
rmax = Spe_Atom_Cut1[spe] + 0.5;
Sr = rmax + rmin;
Dr = rmax - rmin;
for (i=0; i<GL_Mesh; i++){
RGL[i] = 0.50*(Dr*GL_Abscissae[i] + Sr);
}
kmin = Radial_kmin;
kmax = PAO_Nkmax;
Sk = kmax + kmin;
Dk = kmax - kmin;
/* loop for kj */
#pragma omp parallel shared(List_YOUSO,GL_Weight,GL_Abscissae,Dr,Dk,Sk,RGL,Projector_VNA,Spe_VPS_RV,Spe_Num_Mesh_VPS,Spe_VNA_Bessel) private(SumTmp,SphB,tmp_SphB,tmp_SphBp,OMPID,Nthrds,Nprocs,kj,norm_k,i,r,L,Mul,tmp0,dum0)
{
/* allocate arrays */
SumTmp = (double*)malloc(sizeof(double)*List_YOUSO[34]);
SphB = (double**)malloc(sizeof(double*)*(List_YOUSO[35]+3));
for(L=0; L<(List_YOUSO[35]+3); L++){
SphB[L] = (double*)malloc(sizeof(double)*GL_Mesh);
}
tmp_SphB = (double*)malloc(sizeof(double)*(List_YOUSO[35]+3));
tmp_SphBp = (double*)malloc(sizeof(double)*(List_YOUSO[35]+3));
/* get info. on OpenMP */
OMPID = omp_get_thread_num();
Nthrds = omp_get_num_threads();
Nprocs = omp_get_num_procs();
for ( kj=OMPID; kj<GL_Mesh; kj+=Nthrds ){
norm_k = 0.50*(Dk*GL_Abscissae[kj] + Sk);
/* calculate SphB */
for (i=0; i<GL_Mesh; i++){
r = RGL[i];
Spherical_Bessel(norm_k*r,List_YOUSO[35],tmp_SphB,tmp_SphBp);
for(L=0; L<=List_YOUSO[35]; L++){
SphB[L][i] = tmp_SphB[L];
}
}
/* loop for L */
for (L=0; L<=List_YOUSO[35]; L++){
/****************************************************
\int jL(k*r)RL r^2 dr
****************************************************/
for (Mul=0; Mul<List_YOUSO[34]; Mul++) SumTmp[Mul] = 0.0;
/* Gauss-Legendre quadrature */
for (i=0; i<GL_Mesh; i++){
r = RGL[i];
tmp0 = r*r*GL_Weight[i]*SphB[L][i];
for (Mul=0; Mul<List_YOUSO[34]; Mul++){
dum0 = PhiF(r, Projector_VNA[spe][L][Mul], Spe_VPS_RV[spe], Spe_Num_Mesh_VPS[spe]);
SumTmp[Mul] += dum0*tmp0;
}
}
for (Mul=0; Mul<List_YOUSO[34]; Mul++){
Spe_VNA_Bessel[spe][L][Mul][kj] = 0.5*Dr*SumTmp[Mul];
}
} /* L */
} /* kj */
/* free arrays */
free(SumTmp);
for(L=0; L<(List_YOUSO[35]+3); L++){
free(SphB[L]);
}
free(SphB);
free(tmp_SphB);
free(tmp_SphBp);
#pragma omp flush(Spe_VNA_Bessel)
} /* #pragma omp parallel */
} /* Lspe */
/****************************************************
regenerate radial grids in the k-space
for the MPI calculation
****************************************************/
for (kj=0; kj<GL_Mesh; kj++){
kmin = Radial_kmin;
kmax = PAO_Nkmax;
Sk = kmax + kmin;
Dk = kmax - kmin;
norm_k = 0.50*(Dk*GL_Abscissae[kj] + Sk);
GL_NormK[kj] = norm_k;
}
/***********************************************************
sending and receiving of Spe_VNA_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;
for (L=0; L<=List_YOUSO[35]; L++){
for (Mul=0; Mul<List_YOUSO[34]; Mul++){
MPI_Bcast(&Spe_VNA_Bessel[spe][L][Mul][0],
GL_Mesh,MPI_DOUBLE,ID,mpi_comm_level1);
}
}
}
}
/***********************************************************
elapsed time
***********************************************************/
dtime(&TEtime);
/*
printf("myid=%2d Elapsed Time (s) = %15.12f\n",myid,TEtime-TStime);
MPI_Finalize();
exit(0);
*/
}
|