TABLE OF CONTENTS


ABINIT/getgh2c [ Functions ]

[ Top ] [ Functions ]

NAME

 getgh2c

FUNCTION

 Compute <G|H^(2)|C> (or <G|H^(2)-Eps.S^(2)|C>) for input vector |C> expressed in reciprocal space.
 (H^(2) is the 2nd-order pertubed Hamiltonian, S^(2) is the 2nd-order perturbed overlap operator).
 Result is put in array gh2c.
 If required, part of <G|K(2)+Vnonlocal^(2)|C> not depending on VHxc^(2) is also returned in gvnl2.
 If required, <G|S^(2)|C> is returned in gs2c (S=overlap - PAW only)
 Available for the following cases :
  ipert = natom+10 (dkdk)   :       2nd derivative w.r.t wavevector
          natom+11 (dkdE)   : mixed 2nd derivative w.r.t wavector     and eletric field
  also if natom+12<=ipert<=2*natom+11 :
                   (dtaudE) : mixed 2nd derivative w.r.t atom. displ. and eletric field (nonlocal only)

INPUTS

  cwavef(2,npw*nspinor)=input wavefunction, in reciprocal space
  cwaveprj(natom,nspinor*usecprj)=<p_lmn|C> coefficients for wavefunction |C>
  gs_hamkq <type(gs_hamiltonian_type)>=all data for the Hamiltonian
  idir=direction of the perturbation
  ipert=type of the perturbation
  lambda=real use to apply H^(2)-lambda.S^(2)
  mpi_enreg=information about MPI parallelization
  optlocal=0: local part of H^(2) is not computed in gh2c=<G|H^(2)|C>
           1: local part of H^(2) is computed in gh2c=<G|H^(2)|C>
  optnl=0: non-local part of H^(2) is not computed in gh2c=<G|H^(2)|C>
        1: non-local part of H^(2) depending on VHxc^(2) is not computed in gh2c=<G|H^(2)|C>
        2: non-local part of H^(2) is totally computed in gh2c=<G|H^(2)|C>
  opt_gvnl2=option controlling the use of gvnl2 array:
            0: not used
            1: used as input:    - used only for PAW and ipert=natom+11/+12
               At input: contains the derivative w.r.t wavevector of cwavef (times i)
  rf_hamkq <type(rf_hamiltonian_type)>=all data for the 2nd-order Hamiltonian at k,k+q
  sij_opt= -PAW ONLY-  if  0, only matrix elements <G|H^(2)|C> have to be computed
     (S=overlap)       if  1, matrix elements <G|S^(2)|C> have to be computed in gs2c in addition to gh2c
                       if -1, matrix elements <G|H^(2)-lambda.S^(2)|C> have to be computed in gh2c (gs2c not used)
  tim_getgh2c=timing code of the calling subroutine (can be set to 0 if not attributed)
  usevnl=1 if gvnl2=(part of <G|K^(2)+Vnl^(2)-lambda.S^(2)|C> not depending on VHxc^(2)) has to be input/output

OUTPUT

 gh2c(2,npw1*nspinor)= <G|H^(2)|C> or  <G|H^(2)-lambda.S^(2)|C>
                     (only kinetic+non-local parts if optlocal=0)
 if (usevnl==1)
  gvnl2(2,npw1*nspinor)=  part of <G|K^(2)+Vnl^(2)|C> not depending on VHxc^(2)              (sij_opt/=-1)
                       or part of <G|K^(2)+Vnl^(2)-lambda.S^(2)|C> not depending on VHxc^(2) (sij_opt==-1)
 if (sij_opt=1)
  gs2c(2,npw1*nspinor)=<G|S^(2)|C> (S=overlap).

SOURCE

 95 subroutine getgh2c(cwavef,cwaveprj,gh2c,gs2c,gs_hamkq,gvnl2,idir,ipert,lambda,&
 96 &                  mpi_enreg,optlocal,optnl,opt_gvnl2,rf_hamkq,sij_opt,tim_getgh2c,usevnl,conj,enl,optkin)
 97 
 98 !Arguments ------------------------------------
 99 !scalars
100  logical,intent(in),optional :: conj
101  integer,intent(in) :: idir,ipert,optlocal,optnl,opt_gvnl2,sij_opt,tim_getgh2c,usevnl
102  integer,intent(in),optional :: optkin
103  real(dp),intent(in) :: lambda
104  type(MPI_type),intent(in) :: mpi_enreg
105  type(gs_hamiltonian_type),intent(inout),target :: gs_hamkq
106  type(rf_hamiltonian_type),intent(inout),target :: rf_hamkq
107 !arrays
108  real(dp),intent(in),optional,target :: enl(gs_hamkq%dimekb1,gs_hamkq%dimekb2,gs_hamkq%nspinor**2,gs_hamkq%dimekbq)
109  real(dp),intent(inout) :: cwavef(:,:)
110  real(dp),intent(inout),target :: gvnl2(:,:)
111  real(dp),intent(out) :: gh2c(:,:),gs2c(:,:)
112  type(pawcprj_type),intent(inout),target :: cwaveprj(:,:)
113 
114 !Local variables-------------------------------
115 !scalars
116  integer,parameter :: tim_nonlop=0
117  integer :: choice,cpopt,iatm,idir1,idir2,idirc,ipw,ipws,ispinor,my_nspinor
118  integer :: natom,ncpgr,nnlout=1,npw,npw1,paw_opt,signs,usecprj
119  logical :: compute_conjugate,has_kin,has_vnl,pert_phon_elfd
120  real(dp) :: enlout_dum(1)
121  character(len=500) :: msg
122 !arrays
123 ! integer,parameter :: alpha(6)=(/1,2,3,3,3,2/),beta(6)=(/1,2,3,2,1,1/)
124  integer,parameter :: alpha(9)=(/1,2,3,2,1,1,3,3,2/),beta(9)=(/1,2,3,3,3,2,2,1,1/)
125 !real(dp) :: tsec(2)
126  real(dp) :: svectout_dum(1,1),vectout_dum(1,1)
127  real(dp),allocatable :: nonlop_out(:,:)
128  real(dp), pointer :: gvnl2_(:,:)
129  real(dp), pointer :: ddkinpw(:),kinpw1(:),enl_ptr(:,:,:,:)
130  real(dp),allocatable,target :: enl_temp(:,:,:,:)
131  type(pawcprj_type),allocatable,target :: cwaveprj_tmp(:,:)
132  type(pawcprj_type),pointer :: cwaveprj_ptr(:,:)
133 
134 ! *********************************************************************
135 
136  DBG_ENTER("COLL")
137  ABI_UNUSED(tim_getgh2c)
138 
139 !Keep track of total time spent in getgh2c
140 !call timab(196+tim_getgh2c,1,tsec)
141 
142 !======================================================================
143 !== Initialisations and compatibility tests
144 !======================================================================
145 
146  npw  =gs_hamkq%npw_k
147  npw1 =gs_hamkq%npw_kp
148  natom=gs_hamkq%natom
149 
150 !Compatibility tests
151  if(ipert/=natom+10.and.ipert/=natom+11.and.ipert>2*natom+11)then
152    msg='only ipert==natom+10/+11 and natom+11<=ipert<=2*natom+11 implemented!'
153    ABI_BUG(msg)
154  end if
155  pert_phon_elfd = .false.
156  if (ipert>natom+11.and.ipert<=2*natom+11) pert_phon_elfd = .true.
157  if (mpi_enreg%paral_spinor==1) then
158    msg='Not compatible with parallelization over spinorial components!'
159    ABI_BUG(msg)
160  end if
161  if (gs_hamkq%nvloc>1) then
162    msg='Not compatible with nvloc=4 (non-coll. magnetism)!'
163    ABI_BUG(msg)
164  end if
165  if((ipert==natom+11.or.pert_phon_elfd).and.gs_hamkq%usepaw==1.and.optnl>=1) then
166    if (gs_hamkq%nvloc>1) then
167      msg='Not compatible with nvloc=4 (non-coll. magnetism)!'
168      ABI_BUG(msg)
169    end if
170    if (present(enl)) then
171      enl_ptr => enl
172    else if (associated(rf_hamkq%e1kbfr).and.associated(rf_hamkq%e1kbsc).and.optnl==2) then
173      ABI_CHECK(size(rf_hamkq%e1kbfr,4)==1,'BUG in getgh2c: qphase>1!')
174      ABI_CHECK(size(rf_hamkq%e1kbsc,4)==1,'BUG in getgh2c: qphase>1!')
175      ABI_MALLOC(enl_temp,(gs_hamkq%dimekb1,gs_hamkq%dimekb2,gs_hamkq%nspinor**2,gs_hamkq%dimekbq))
176      enl_temp(:,:,:,:) = rf_hamkq%e1kbfr(:,:,:,:) + rf_hamkq%e1kbsc(:,:,:,:)
177      enl_ptr => enl_temp
178    else if (associated(rf_hamkq%e1kbfr)) then
179      ABI_CHECK(size(rf_hamkq%e1kbfr,4)==1,'BUG in getgh2c: qphase>1!')
180      enl_ptr => rf_hamkq%e1kbfr
181    else
182      msg='For ipert=natom+11/pert_phon_elfd : e1kbfr and/or e1kbsc must be associated or enl optional input must be present.'
183      ABI_BUG(msg)
184    end if
185    if (usevnl==0) then
186      msg='gvnl2 must be allocated for ipert=natom+11/pert_phon_elfd !'
187      ABI_BUG(msg)
188    end if
189    if(opt_gvnl2==0) then
190      msg='opt_gvnl2=0 not compatible with ipert=natom+11/pert_phon_elfd !'
191      ABI_BUG(msg)
192    end if
193  end if
194 
195 !Check sizes
196  my_nspinor=max(1,gs_hamkq%nspinor/mpi_enreg%nproc_spinor)
197  if (size(cwavef)<2*npw*my_nspinor) then
198    msg='wrong size for cwavef!'
199    ABI_BUG(msg)
200  end if
201  if (size(gh2c)<2*npw1*my_nspinor) then
202    msg='wrong size for gh2c!'
203    ABI_BUG(msg)
204  end if
205  if (usevnl/=0) then
206    if (size(gvnl2)<2*npw1*my_nspinor) then
207      msg='wrong size for gvnl2!'
208      ABI_BUG(msg)
209    end if
210  end if
211  if (sij_opt==1) then
212    if (size(gs2c)<2*npw1*my_nspinor) then
213      msg='wrong size for gs2c!'
214      ABI_BUG(msg)
215    end if
216  end if
217 
218 !PAW: specific treatment for usecprj input arg
219 !     force it to zero if cwaveprj is not allocated
220  usecprj=gs_hamkq%usecprj ; ncpgr=0
221  if(gs_hamkq%usepaw==1) then
222    if (size(cwaveprj)==0) usecprj=0
223    if (usecprj/=0) then
224      ncpgr=cwaveprj(1,1)%ncpgr
225      if (size(cwaveprj)<gs_hamkq%natom*my_nspinor) then
226        msg='wrong size for cwaveprj!'
227        ABI_BUG(msg)
228      end if
229    end if
230  else
231    if(usecprj==1)then
232      msg='usecprj==1 not allowed for NC psps !'
233      ABI_BUG(msg)
234    end if
235  end if
236 
237 ! tim_nonlop=8
238 ! if (tim_getgh2c==1.and.ipert<=natom) tim_nonlop=7
239 ! if (tim_getgh2c==2.and.ipert<=natom) tim_nonlop=5
240 ! if (tim_getgh2c==1.and.ipert> natom) tim_nonlop=8
241 ! if (tim_getgh2c==2.and.ipert> natom) tim_nonlop=5
242 ! if (tim_getgh2c==3                 ) tim_nonlop=0
243 
244  idir1=alpha(idir);idir2=beta(idir)
245 
246  compute_conjugate = .false.
247  if(present(conj)) compute_conjugate = conj
248 
249 !======================================================================
250 !== Apply the 2nd-order local potential to the wavefunction
251 !======================================================================
252 
253  if (ipert/=natom+10.and.ipert/=natom+11.and.optlocal>0) then
254    msg='local part not implemented'
255    ABI_BUG(msg)
256  else
257 !  In the case of ddk operator, no local contribution (also because no self-consistency)
258 !$OMP PARALLEL DO
259    do ipw=1,npw1*my_nspinor
260      gh2c(:,ipw)=zero
261    end do
262 
263  end if
264 
265 !======================================================================
266 !== Apply the 2st-order non-local potential to the wavefunction
267 !======================================================================
268 
269  has_vnl=(ipert==natom+10.or.ipert==natom+11.or.pert_phon_elfd)
270 
271 !Use of gvnl2 depends on usevnl
272  if (usevnl==1) then
273    gvnl2_ => gvnl2
274  else
275    ABI_MALLOC(gvnl2_,(2,npw1*my_nspinor))
276  end if
277 
278  if (has_vnl.and.(optnl>0.or.sij_opt/=0)) then
279 
280    idirc=3*(idir1-1)+idir2 !xx=1, xy=2, xz=3, yx=4, yy=5, yz=6, zx=7, zy=8, zz=9, (xyz,xyz)=(idir1,idir2)
281 
282 ! d^2[H_nl]/dk1dk2
283 !  -------------------------------------------
284    if (ipert==natom+10) then
285      if (gs_hamkq%usepaw==1) then
286        if (usecprj==1) then
287          cwaveprj_ptr => cwaveprj
288        else
289          ABI_MALLOC(cwaveprj_tmp,(natom,my_nspinor))
290          call pawcprj_alloc(cwaveprj_tmp,0,gs_hamkq%dimcprj)
291          cwaveprj_ptr => cwaveprj_tmp
292        end if
293        cpopt=-1+5*usecprj
294        choice=8; signs=2; paw_opt=1; if (sij_opt/=0) paw_opt=sij_opt+3
295        call nonlop(choice,cpopt,cwaveprj_ptr,enlout_dum,gs_hamkq,idirc,(/lambda/),mpi_enreg,1,nnlout,&
296 &       paw_opt,signs,gs2c,tim_nonlop,cwavef,gvnl2_)
297        if (usecprj==0) then
298          call pawcprj_free(cwaveprj_tmp)
299          ABI_FREE(cwaveprj_tmp)
300        end if
301        nullify(cwaveprj_ptr)
302      else
303        choice=8; signs=2; cpopt=-1 ; paw_opt=0
304        call nonlop(choice,cpopt,cwaveprj,enlout_dum,gs_hamkq,idirc,(/zero/),mpi_enreg,1,nnlout,&
305 &       paw_opt,signs,svectout_dum,tim_nonlop,cwavef,gvnl2_)
306      end if
307 
308 ! d^2[H_nl]/dk1dE2 : Non-zero only in PAW
309 !  -------------------------------------------
310    else if (ipert==natom+11.and.gs_hamkq%usepaw==1) then
311 
312      ABI_MALLOC(nonlop_out,(2,npw1*my_nspinor))
313 
314      if (usecprj==1) then
315        cwaveprj_ptr => cwaveprj
316      else
317        ABI_MALLOC(cwaveprj_tmp,(natom,my_nspinor))
318        call pawcprj_alloc(cwaveprj_tmp,2,gs_hamkq%dimcprj)
319        cwaveprj_ptr => cwaveprj_tmp
320      end if
321 
322      if (opt_gvnl2==1.and.optnl>=1) then
323 
324 !      Compute application of dS/dk1 to i*d[cwavef]/dk2
325 !      sum_{i,j} s_ij d(|p_i><p_j|)/dk(idir1) | i*psi^(k(idir2)) >
326        cpopt=-1 ; choice=5 ; paw_opt=3 ; signs=2
327        call nonlop(choice,cpopt,cwaveprj_ptr,enlout_dum,gs_hamkq,idir1,(/zero/),mpi_enreg,1,nnlout,&
328 &       paw_opt,signs,nonlop_out,tim_nonlop,gvnl2_,vectout_dum)
329 
330 !$OMP PARALLEL DO
331        do ipw=1,npw1*my_nspinor
332          gvnl2_(:,ipw)=nonlop_out(:,ipw)
333        end do
334 
335 !      Compute part of H^(2) due to derivative of projectors (idir1) and derivative of Dij (idir2)
336 !      sum_{i,j} chi_ij(idir2) d(|p_i><p_j|)/dk(idir1) | psi^(0) >
337        cpopt=4*usecprj ; choice=5 ; paw_opt=1 ; signs=2
338        call nonlop(choice,cpopt,cwaveprj_ptr,enlout_dum,gs_hamkq,idir1,(/zero/),mpi_enreg,1,nnlout,&
339 &       paw_opt,signs,svectout_dum,tim_nonlop,cwavef,nonlop_out,enl=enl_ptr)
340 
341 !$OMP PARALLEL DO
342        do ipw=1,npw1*my_nspinor
343          gvnl2_(:,ipw)=gvnl2_(:,ipw)+nonlop_out(:,ipw)
344        end do
345 
346      else
347 
348 !$OMP PARALLEL DO
349        do ipw=1,npw1*my_nspinor
350          gvnl2_(:,ipw)=zero
351        end do
352 
353      end if ! opt_gvnl2==1
354 
355 !    Compute derivatives due to projectors |d^2[p_i]/dk1dk2>,|d[p_i]/dk1>,|d[p_i]/dk2>
356 !    i * sum_{i,j} (d(|p_i><dp_j/dk(idir2)|)/dk(idir1) | psi^(0) >
357      cpopt=-1+5*usecprj ; choice=81 ; paw_opt=3 ; signs=2
358      call nonlop(choice,cpopt,cwaveprj_ptr,enlout_dum,gs_hamkq,idirc,(/zero/),mpi_enreg,1,nnlout,&
359 &     paw_opt,signs,nonlop_out,tim_nonlop,cwavef,vectout_dum)
360 
361      if(compute_conjugate) then
362 !$OMP PARALLEL DO
363        do ipw=1,npw1*my_nspinor ! Note the multiplication by -i
364          gvnl2_(1,ipw)=gvnl2_(1,ipw)+nonlop_out(2,ipw)
365          gvnl2_(2,ipw)=gvnl2_(2,ipw)-nonlop_out(1,ipw)
366        end do
367      else
368 !$OMP PARALLEL DO
369        do ipw=1,npw1*my_nspinor ! Note the multiplication by i
370          gvnl2_(1,ipw)=gvnl2_(1,ipw)-nonlop_out(2,ipw)
371          gvnl2_(2,ipw)=gvnl2_(2,ipw)+nonlop_out(1,ipw)
372        end do
373      end if
374 
375      ABI_FREE(nonlop_out)
376      if (sij_opt==1) gs2c=zero
377      if (usecprj==0) then
378        call pawcprj_free(cwaveprj_tmp)
379        ABI_FREE(cwaveprj_tmp)
380      end if
381      nullify(cwaveprj_ptr)
382 
383 ! d^2[H_nl]/dtau1dE2 : Non-zero only in PAW
384 !  -------------------------------------------
385    else if (pert_phon_elfd.and.gs_hamkq%usepaw==1) then
386 
387      iatm = ipert-(natom+11)
388      if (iatm<1.or.iatm>natom) then
389        ABI_BUG(" iatm must be between 1 and natom")
390      end if
391 
392      ABI_MALLOC(nonlop_out,(2,npw1*my_nspinor))
393 
394      if (usecprj==1) then
395        cwaveprj_ptr => cwaveprj
396      else
397        ABI_MALLOC(cwaveprj_tmp,(natom,my_nspinor))
398        call pawcprj_alloc(cwaveprj_tmp,2,gs_hamkq%dimcprj)
399        cwaveprj_ptr => cwaveprj_tmp
400      end if
401 
402      if (opt_gvnl2==1) then
403 
404 !      Compute application of dS/dtau1 to i*d[cwavef]/dk2
405 !      sum_{i,j} s_ij d(|p_i><p_j|)/dtau(idir1) | i*psi^(k(idir2)) >
406        cpopt=-1 ; choice=2 ; paw_opt=3 ; signs=2
407        call nonlop(choice,cpopt,cwaveprj_ptr,enlout_dum,gs_hamkq,idir1,(/zero/),mpi_enreg,1,nnlout,&
408 &       paw_opt,signs,nonlop_out,tim_nonlop,gvnl2_,vectout_dum,iatom_only=iatm)
409 
410 !$OMP PARALLEL DO
411        do ipw=1,npw1*my_nspinor
412          gvnl2_(:,ipw)=nonlop_out(:,ipw)
413        end do
414 
415 !      Compute part of H^(2) due to derivative of projectors (idir1) and derivative of Dij (idir2)
416 !      sum_{i,j} chi_ij(idir2) d(|p_i><p_j|)/dtau(idir1) | psi^(0) >
417        cpopt=4*usecprj ; choice=2 ; paw_opt=1 ; signs=2
418        call nonlop(choice,cpopt,cwaveprj_ptr,enlout_dum,gs_hamkq,idir1,(/zero/),mpi_enreg,1,nnlout,&
419 &       paw_opt,signs,svectout_dum,tim_nonlop,cwavef,nonlop_out,enl=enl_ptr,iatom_only=iatm)
420 
421 !$OMP PARALLEL DO
422        do ipw=1,npw1*my_nspinor
423          gvnl2_(:,ipw)=gvnl2_(:,ipw)+nonlop_out(:,ipw)
424        end do
425 
426      else
427 
428 !$OMP PARALLEL DO
429        do ipw=1,npw1*my_nspinor
430          gvnl2_(:,ipw)=zero
431        end do
432 
433      end if ! opt_gvnl2==1
434 
435 !    Compute derivatives due to projectors |d^2[p_i]/dtau1dk2>,|d[p_i]/dtau1>,|d[p_i]/dk2>
436 !    i * sum_{i,j} (d(|p_i><dp_j/dk(idir2)|)/dtau(idir1) | psi^(0) >
437      cpopt=-1+5*usecprj ; choice=54 ; paw_opt=3 ; signs=2
438      call nonlop(choice,cpopt,cwaveprj_ptr,enlout_dum,gs_hamkq,idirc,(/zero/),mpi_enreg,1,nnlout,&
439 &     paw_opt,signs,nonlop_out,tim_nonlop,cwavef,vectout_dum,iatom_only=iatm)
440 
441      if(compute_conjugate) then
442 !$OMP PARALLEL DO
443        do ipw=1,npw1*my_nspinor ! Note the multiplication by -i
444          gvnl2_(1,ipw)=gvnl2_(1,ipw)+nonlop_out(2,ipw)
445          gvnl2_(2,ipw)=gvnl2_(2,ipw)-nonlop_out(1,ipw)
446        end do
447      else
448 !$OMP PARALLEL DO
449        do ipw=1,npw1*my_nspinor ! Note the multiplication by i
450          gvnl2_(1,ipw)=gvnl2_(1,ipw)-nonlop_out(2,ipw)
451          gvnl2_(2,ipw)=gvnl2_(2,ipw)+nonlop_out(1,ipw)
452        end do
453      end if
454 
455      ABI_FREE(nonlop_out)
456      if (sij_opt==1) gs2c=zero
457      if (usecprj==0) then
458        call pawcprj_free(cwaveprj_tmp)
459        ABI_FREE(cwaveprj_tmp)
460      end if
461      nullify(cwaveprj_ptr)
462 
463    end if
464 
465 !No non-local part
466 !-------------------------------------------
467  else
468 
469    if (optnl>=1) then
470  !$OMP PARALLEL DO
471      do ipw=1,npw1*my_nspinor
472        gvnl2_(:,ipw)=zero
473      end do
474    end if
475    if (sij_opt/=0) then
476  !$OMP PARALLEL DO
477      do ipw=1,npw1*my_nspinor
478        gs2c(:,ipw)=zero
479      end do
480    end if
481 
482  end if
483 
484  if (associated(enl_ptr)) then
485    nullify(enl_ptr)
486  end if
487  if (allocated(enl_temp)) then
488    ABI_FREE(enl_temp)
489  end if
490 
491 !======================================================================
492 !== Apply the 2nd-order kinetic operator to the wavefunction
493 !======================================================================
494 
495  if (present(optkin)) then
496    has_kin=(optkin/=0.and.ipert==natom+10)
497  else
498    has_kin=(ipert==natom+10)
499  end if
500 
501 !k-point perturbation
502 !-------------------------------------------
503  if (associated(gs_hamkq%kinpw_kp)) then
504    kinpw1 => gs_hamkq%kinpw_kp
505  else if (optnl>=1.or.has_kin) then
506    msg='need kinpw1 allocated!'
507    ABI_BUG(msg)
508  end if
509  if (associated(rf_hamkq%ddkinpw_k)) then
510    ddkinpw => rf_hamkq%ddkinpw_k
511  else if (has_kin) then
512    msg='need ddkinpw allocated!'
513    ABI_BUG(msg)
514  end if
515 
516  if (has_kin) then
517    do ispinor=1,my_nspinor
518  !$OMP PARALLEL DO PRIVATE(ipw,ipws) SHARED(cwavef,ispinor,gvnl2_,ddkinpw,kinpw1,npw,my_nspinor)
519      do ipw=1,npw
520        ipws=ipw+npw*(ispinor-1)
521        if(kinpw1(ipw)<huge(zero)*1.d-11)then
522          gvnl2_(1,ipws)=gvnl2_(1,ipws)+ddkinpw(ipw)*cwavef(1,ipws)
523          gvnl2_(2,ipws)=gvnl2_(2,ipws)+ddkinpw(ipw)*cwavef(2,ipws)
524        else
525          gvnl2_(1,ipws)=zero
526          gvnl2_(2,ipws)=zero
527        end if
528      end do
529    end do
530  end if
531 
532 !======================================================================
533 !== Sum contributions to get the application of H^(2) to the wf
534 !======================================================================
535 !Also filter the wavefunctions for large modified kinetic energy
536 
537 !Add non-local+kinetic to local part
538  if (optnl>=1.or.has_kin) then
539    do ispinor=1,my_nspinor
540      ipws=(ispinor-1)*npw1
541  !$OMP PARALLEL DO PRIVATE(ipw) SHARED(gh2c,gvnl2_,kinpw1,ipws,npw1)
542      do ipw=1+ipws,npw1+ipws
543        if(kinpw1(ipw-ipws)<huge(zero)*1.d-11)then
544          gh2c(1,ipw)=gh2c(1,ipw)+gvnl2_(1,ipw)
545          gh2c(2,ipw)=gh2c(2,ipw)+gvnl2_(2,ipw)
546        else
547          gh2c(1,ipw)=zero
548          gh2c(2,ipw)=zero
549        end if
550      end do
551    end do
552  end if
553 
554  if (usevnl==1) then
555    nullify(gvnl2_)
556  else
557    ABI_FREE(gvnl2_)
558  end if
559 
560 !call timab(196+tim_getgh2c,2,tsec)
561 
562  DBG_EXIT("COLL")
563 
564 end subroutine getgh2c

ABINIT/m_getgh2c [ Modules ]

[ Top ] [ Modules ]

NAME

  m_getgh2c

FUNCTION

COPYRIGHT

  Copyright (C) 2015-2024 ABINIT group (MT,JLJ)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

SOURCE

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 module m_getgh2c
23 
24  use defs_basis
25  use m_abicore
26  use m_errors
27 
28  use defs_abitypes, only : mpi_type
29  use m_pawcprj,     only : pawcprj_type,pawcprj_alloc,pawcprj_free
30  use m_hamiltonian, only : gs_hamiltonian_type,rf_hamiltonian_type
31  use m_nonlop,      only : nonlop
32 
33  implicit none
34 
35  private