TABLE OF CONTENTS


ABINIT/m_frskerker1 [ Modules ]

[ Top ] [ Modules ]

NAME

 m_frskerker1

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-2024 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

SOURCE

23 #if defined HAVE_CONFIG_H
24 #include "config.h"
25 #endif
26 
27 #include "abi_common.h"
28 
29 module m_frskerker1
30 
31   use defs_basis
32   use m_abicore
33   use m_dtset
34 
35   use defs_abitypes, only : MPI_type
36   use m_spacepar,  only : laplacian
37   use m_numeric_tools, only : dotproduct
38 
39   implicit none
40   !! common variables copied from input
41   integer,save,private                  :: nfft,nspden,ngfft(18)
42   real(dp),save,allocatable,private    :: deltaW(:,:),mat(:,:),g2cart(:)
43   real(dp),save,private                :: gprimd(3,3),dielng
44   type(dataset_type),pointer,save,private  :: dtset_ptr
45   type(MPI_type),save,private,pointer :: mpi_enreg_ptr
46   !! common variables computed
47   logical,save,private :: ok=.false.
48 
49 contains

m_frskerker1/frskerker1__dpf [ Functions ]

[ Top ] [ m_frskerker1 ] [ Functions ]

NAME

 frskerker1__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

INPUTS

OUTPUT

SOURCE

215 function frskerker1__dpf(nv1,nv2,vrespc)
216 
217 !Arguments ------------------------------------
218  integer,intent(in) :: nv1,nv2
219  real(dp),intent(in):: vrespc(nv1,nv2)
220  real(dp) ::frskerker1__dpf(nv1,nv2)
221 
222 !Local variables-------------------------------
223  real(dp) :: buffer1(nv1,nv2),buffer2(nv1,nv2)
224 ! *************************************************************************
225 
226   if(ok) then
227    buffer1=vrespc
228    call laplacian(gprimd,mpi_enreg_ptr,nfft,nspden,ngfft,&
229 &    rdfuncr=buffer1,laplacerdfuncr=buffer2,g2cart_in=g2cart)
230    frskerker1__dpf(:,:)= vrespc(:,:)-deltaW-((dielng)**2)*buffer2(:,:)
231   else
232    frskerker1__dpf = zero
233   end if
234 
235 end function frskerker1__dpf

m_frskerker1/frskerker1__end [ Functions ]

[ Top ] [ m_frskerker1 ] [ Functions ]

NAME

 frskerker1__end

FUNCTION

 ending subroutine
 deallocate memory areas

INPUTS

OUTPUT

SOURCE

115 subroutine frskerker1__end()
116 
117 ! *************************************************************************
118   if(ok) then
119 !  ! set ok to false which prevent using the pf and dpf
120    ok = .false.
121 !  ! free memory
122    ABI_FREE(deltaW)
123    ABI_FREE(mat)
124    ABI_FREE(g2cart)
125   end if
126 
127  end subroutine frskerker1__end

m_frskerker1/frskerker1__init [ Functions ]

[ Top ] [ m_frskerker1 ] [ Functions ]

NAME

 frskerker1__init

FUNCTION

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

INPUTS

OUTPUT

SOURCE

67 subroutine frskerker1__init(dtset_in,mpi_enreg_in,nfft_in,ngfft_in,nspden_in,dielng_in,deltaW_in,gprimd_in,mat_in,g2cart_in )
68 
69 !Arguments ------------------------------------
70  integer,intent(in)                  :: nfft_in,ngfft_in(18),nspden_in
71  real(dp),intent(in)                 :: deltaW_in(nfft_in,nspden_in),mat_in(nfft_in,nspden_in),g2cart_in(nfft_in)
72  real(dp),dimension(3,3),intent(in)  :: gprimd_in
73  real(dp),intent(in)                 :: dielng_in
74  type(dataset_type),target,intent(in) :: dtset_in
75  type(MPI_type),target,intent(in)  :: mpi_enreg_in
76 
77 ! *************************************************************************
78 
79 ! !allocation and data transfer
80 ! !Thought it would have been more logical to use the privates intrinsic of the module as
81 ! !input variables it seems that it is not possible...
82   if(.not.ok) then
83    dtset_ptr  => dtset_in
84    mpi_enreg_ptr => mpi_enreg_in
85    nspden=nspden_in
86    ngfft=ngfft_in
87    nfft=nfft_in
88    ABI_MALLOC(deltaW,(size(deltaW_in,1),size(deltaW_in,2)))
89    ABI_MALLOC(mat,(size(mat_in,1),size(mat_in,2)))
90    ABI_MALLOC(g2cart,(size(g2cart_in,1)))
91    deltaW=deltaW_in
92    dielng=dielng_in
93    mat=mat_in
94    gprimd=gprimd_in
95    g2cart=g2cart_in
96    ok = .true.
97   end if
98  end subroutine frskerker1__init

m_frskerker1/frskerker1__newvres [ Functions ]

[ Top ] [ m_frskerker1 ] [ Functions ]

NAME

 frskerker1__newvres

FUNCTION

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

INPUTS

OUTPUT

SOURCE

145 subroutine frskerker1__newvres(nv1,nv2,x, grad, vrespc)
146 
147 !Arguments ------------------------------------
148  integer,intent(in) :: nv1,nv2
149  real(dp),intent(in):: x
150  real(dp),intent(inout)::grad(nv1,nv2)
151  real(dp),intent(inout)::vrespc(nv1,nv2)
152 
153 ! *************************************************************************
154 
155   grad(:,:)=x*grad(:,:)
156   vrespc(:,:)=vrespc(:,:)+grad(:,:)
157 
158 end subroutine frskerker1__newvres

m_frskerker1/frskerker1__pf [ Functions ]

[ Top ] [ m_frskerker1 ] [ Functions ]

NAME

 frskerker1__pf

FUNCTION

 penalty function associated with the preconditionned residuals

INPUTS

OUTPUT

SOURCE

174 function frskerker1__pf(nv1,nv2,vrespc)
175 
176 !Arguments ------------------------------------
177  integer,intent(in) :: nv1,nv2
178  real(dp),intent(in)::vrespc(nv1,nv2)
179  real(dp)           ::frskerker1__pf
180 
181 !Local variables-------------------------------
182  real(dp):: buffer1(nv1,nv2),buffer2(nv1,nv2)
183 ! *************************************************************************
184 
185   if(ok) then
186    buffer1=vrespc
187    call laplacian(gprimd,mpi_enreg_ptr,nfft,nspden,ngfft,&
188 &    rdfuncr=buffer1,laplacerdfuncr=buffer2,g2cart_in=g2cart)
189    buffer2(:,:)=(vrespc(:,:)-((dielng)**2)*buffer2(:,:)) * half -deltaW
190    frskerker1__pf=dotproduct(nv1,nv2,vrespc,buffer2) !*half-dotproduct(nv1,nv2,vrespc,deltaW)
191   else
192    frskerker1__pf=zero
193   end if
194 
195  end function frskerker1__pf