TABLE OF CONTENTS


ABINIT/frskerker2 [ Modules ]

[ Top ] [ Modules ]

NAME

 frskerker2

FUNCTION

 provide the ability to compute the
 penalty function and its first derivative associated
 with some residuals and a real space dielectric function

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DCA, XG, MT)
 This file is distributed under the terms of the
 GNU General Public License, see ~ABINIT/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~ABINIT/Infos/contributors .

NOTES

 this is neither a function nor a subroutine. This is a module
 It is made of two functions and one init subroutine

PARENTS

CHILDREN

SOURCE

28 #if defined HAVE_CONFIG_H
29 #include "config.h"
30 #endif
31 
32 #include "abi_common.h"
33 
34 module frskerker2 
35 
36   use defs_basis
37   use defs_abitypes
38   use m_profiling_abi
39   use interfaces_32_util        ! THIS IS MANDATORY TO CALL dotproduct
40   use interfaces_53_spacepar
41   use interfaces_56_recipspace  ! THIS IS MANDATORY TO CALL LAPLACIAN
42 
43   implicit none
44 
45   !! common variables copied from input
46   integer,save,private                  :: nfft,nspden,ngfft(18)
47   real(dp),save,allocatable,private    :: deltaW(:,:),mat(:,:),rdielng(:)
48   real(dp),save,private                :: gprimd(3,3)
49   type(dataset_type),pointer,save,private  :: dtset_ptr
50   type(MPI_type),pointer,save,private  :: mpi_enreg_ptr
51   !! common variables computed
52   logical,save,private :: ok=.false.
53 ! *************************************************************************
54 
55 contains

frskerker2/frskerker2__dpf [ Functions ]

[ Top ] [ frskerker2 ] [ Functions ]

NAME

 frskerker2__dpf

FUNCTION

 derivative of the penalty function
 actually not the derivative but something allowing minimization
 at constant density
 formula from the work of rackowski,canning and wang
 H*phi - int(phi**2H d3r)phi
 that is the simple projection of the change on a direction
 normal to the density changes

INPUTS

OUTPUT

PARENTS

  Will be filled automatically by the parent script

CHILDREN

  Will be filled automatically by the parent script

SOURCE

288   function frskerker2__dpf(nv1,nv2,vrespc)
289 
290 
291 !This section has been created automatically by the script Abilint (TD).
292 !Do not modify the following lines by hand.
293 #undef ABI_FUNC
294 #define ABI_FUNC 'frskerker2__dpf'
295  use interfaces_56_recipspace
296 !End of the abilint section
297 
298     implicit none
299 !Arguments ------------------------------------
300     integer,intent(in) :: nv1,nv2
301     real(dp),intent(in)::vrespc(nv1,nv2)
302     real(dp)           :: frskerker2__dpf(nv1,nv2)
303 !Local variables-------------------------------
304 
305     real(dp):: buffer1(nv1,nv2),buffer2(nv1,nv2)
306     integer :: ispden
307 
308 ! *************************************************************************
309 
310   if(ok) then
311    buffer1=vrespc
312    call laplacian(gprimd,mpi_enreg_ptr,nfft,nspden,ngfft,dtset_ptr%paral_kgb,rdfuncr=buffer1,laplacerdfuncr=buffer2)
313    do ispden=1,nspden
314     frskerker2__dpf(:,ispden)= vrespc(:,ispden)-deltaW(:,ispden)-((rdielng(:))**2)*buffer2(:,ispden)
315    end do
316   else
317    frskerker2__dpf = zero
318   end if
319  end function frskerker2__dpf
320 
321 end module frskerker2

frskerker2/frskerker2__end [ Functions ]

[ Top ] [ frskerker2 ] [ Functions ]

NAME

 frskerker2__end

FUNCTION

 ending subroutine
 deallocate memory areas

INPUTS

OUTPUT

PARENTS

      prcrskerker2

CHILDREN

      laplacian

SOURCE

137   subroutine frskerker2__end()
138     use defs_basis
139     use defs_abitypes
140 
141 !This section has been created automatically by the script Abilint (TD).
142 !Do not modify the following lines by hand.
143 #undef ABI_FUNC
144 #define ABI_FUNC 'frskerker2__end'
145 !End of the abilint section
146 
147     implicit none
148 !Arguments ------------------------------------
149 
150 !Local variables-------------------------------
151 
152 ! *************************************************************************
153   if(ok) then
154 !  ! set ok to false which prevent using the pf and dpf
155    ok = .false.
156 !  ! free memory
157    ABI_DEALLOCATE(deltaW)
158    ABI_DEALLOCATE(mat)
159    ABI_DEALLOCATE(rdielng)
160   end if
161  end subroutine frskerker2__end

frskerker2/frskerker2__init [ Functions ]

[ Top ] [ frskerker2 ] [ Functions ]

NAME

 frskerker2__init

FUNCTION

 initialisation subroutine
 Copy every variables required for the energy calculation
 Allocate the required memory

INPUTS

OUTPUT

PARENTS

      prcrskerker2

CHILDREN

      laplacian

SOURCE

 78   subroutine frskerker2__init(dtset_in,mpi_enreg_in,nfft_in,ngfft_in,nspden_in,rdielng_in,deltaW_in,gprimd_in,mat_in )
 79 
 80 
 81 !This section has been created automatically by the script Abilint (TD).
 82 !Do not modify the following lines by hand.
 83 #undef ABI_FUNC
 84 #define ABI_FUNC 'frskerker2__init'
 85 !End of the abilint section
 86 
 87     implicit none
 88 !Arguments ------------------------------------
 89     type(dataset_type),target,intent(in) :: dtset_in
 90     integer,intent(in)  :: nfft_in,ngfft_in(18),nspden_in
 91     real(dp),intent(in) :: deltaW_in(nfft_in,nspden_in),mat_in(nfft_in,nspden_in)
 92     real(dp),intent(in) :: rdielng_in(nfft_in)
 93     real(dp),intent(in)  :: gprimd_in(3,3)
 94     type(MPI_type),target,intent(in)  :: mpi_enreg_in
 95 
 96 ! *************************************************************************
 97 ! !allocation and data transfer
 98 ! !Thought it would have been more logical to use the privates intrinsic of the module as
 99 ! !input variables it seems that it is not possible...
100   if(.not.ok) then
101    dtset_ptr => dtset_in
102    mpi_enreg_ptr => mpi_enreg_in
103    nspden=nspden_in
104    ngfft=ngfft_in
105    nfft=nfft_in
106    ABI_ALLOCATE(deltaW,(size(deltaW_in,1),size(deltaW_in,2)))
107    ABI_ALLOCATE(mat,(size(mat_in,1),size(mat_in,2)))
108    ABI_ALLOCATE(rdielng,(size(rdielng_in)))
109    deltaW=deltaW_in
110    rdielng=rdielng_in
111    mat=mat_in
112    gprimd=gprimd_in
113    ok = .true.
114   end if
115  end subroutine frskerker2__init

frskerker2/frskerker2__newvres2 [ Functions ]

[ Top ] [ frskerker2 ] [ Functions ]

NAME

 frskerker2__newvres2

FUNCTION

 affectation subroutine
 do the required renormalisation when providing a new value for
 the density after application of the gradient

INPUTS

OUTPUT

PARENTS

CHILDREN

      laplacian

SOURCE

183   subroutine frskerker2__newvres2(nv1,nv2,x, grad, vrespc)
184 
185 
186 !This section has been created automatically by the script Abilint (TD).
187 !Do not modify the following lines by hand.
188 #undef ABI_FUNC
189 #define ABI_FUNC 'frskerker2__newvres2'
190 !End of the abilint section
191 
192     implicit none
193 !Arguments ------------------------------------
194 
195     integer,intent(in)    :: nv1,nv2
196     real(dp),intent(in)   :: x
197     real(dp),intent(inout)::grad(nv1,nv2)
198     real(dp),intent(inout)::vrespc(nv1,nv2)
199 !Local variables-------------------------------
200 
201 ! *************************************************************************
202   grad(:,:)=x*grad(:,:)
203   vrespc(:,:)=vrespc(:,:)+grad(:,:)
204  end subroutine frskerker2__newvres2

frskerker2/frskerker2__pf [ Functions ]

[ Top ] [ frskerker2 ] [ Functions ]

NAME

 frskerker2__pf

FUNCTION

 penalty function associated with the preconditionned residuals    !!

INPUTS

OUTPUT

PARENTS

  Will be filled automatically by the parent script

CHILDREN

  Will be filled automatically by the parent script

SOURCE

225   function frskerker2__pf(nv1,nv2,vrespc)
226 
227 
228 !This section has been created automatically by the script Abilint (TD).
229 !Do not modify the following lines by hand.
230 #undef ABI_FUNC
231 #define ABI_FUNC 'frskerker2__pf'
232  use interfaces_56_recipspace
233  use interfaces_62_cg_noabirule
234 !End of the abilint section
235 
236     implicit none
237 !Arguments ------------------------------------
238     integer,intent(in)    :: nv1,nv2
239     real(dp),intent(in) ::vrespc(nv1,nv2)
240     real(dp)            ::frskerker2__pf
241 !Local variables-------------------------------
242 
243     real(dp)            :: buffer1(nv1,nv2),buffer2(nv1,nv2)
244     integer             :: ispden
245 ! *************************************************************************
246 
247 ! Added by YP [HAVE_FORTRAN_INTERFACES]
248 
249   if(ok) then
250    buffer1=vrespc
251    call laplacian(gprimd,mpi_enreg_ptr,nfft,nspden,ngfft,dtset_ptr%paral_kgb,rdfuncr=buffer1,laplacerdfuncr=buffer2)
252    do ispden=1,nspden
253     buffer2(:,ispden)=(vrespc(:,ispden)-((rdielng(:))**2)*buffer2(:,ispden))  &
254 &    *half  -  deltaW(:,ispden)
255    end do
256 !  pf_rscgres=dotproduct(vrespc,buffer2)*half-dotproduct(vrespc,deltaW)
257    frskerker2__pf=dotproduct(nv1,nv2,vrespc,buffer2)
258   else
259    frskerker2__pf=zero
260   end if
261  end function frskerker2__pf