TABLE OF CONTENTS
- ABINIT/frskerker1
- frskerker1/frskerker1__dpf
- frskerker1/frskerker1__end
- frskerker1/frskerker1__init
- frskerker1/frskerker1__newvres
- frskerker1/frskerker1__pf
ABINIT/frskerker1 [ Modules ]
NAME
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-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
29 #if defined HAVE_CONFIG_H 30 #include "config.h" 31 #endif 32 33 #include "abi_common.h" 34 35 module frskerker1 36 37 use defs_basis 38 use defs_abitypes 39 use m_profiling_abi 40 use interfaces_32_util ! THIS IS MANDATORY TO CALL dotproduct 41 use interfaces_53_spacepar 42 use interfaces_56_recipspace ! THIS IS MANDATORY TO CALL LAPLACIAN 43 44 implicit none 45 !! common variables copied from input 46 integer,save,private :: nfft,nspden,ngfft(18) 47 real(dp),save,allocatable,private :: deltaW(:,:),mat(:,:),g2cart(:) 48 real(dp),save,private :: gprimd(3,3),dielng 49 type(dataset_type),pointer,save,private :: dtset_ptr 50 type(MPI_type),save,private,pointer :: mpi_enreg_ptr 51 !! common variables computed 52 logical,save,private :: ok=.false. 53 54 contains
frskerker1/frskerker1__dpf [ Functions ]
[ Top ] [ 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
PARENTS
Will be filled automatically by the parent script
CHILDREN
Will be filled automatically by the parent script
SOURCE
283 function frskerker1__dpf(nv1,nv2,vrespc) 284 285 286 !This section has been created automatically by the script Abilint (TD). 287 !Do not modify the following lines by hand. 288 #undef ABI_FUNC 289 #define ABI_FUNC 'frskerker1__dpf' 290 use interfaces_56_recipspace 291 !End of the abilint section 292 293 implicit none 294 !Arguments ------------------------------------ 295 integer,intent(in) :: nv1,nv2 296 real(dp),intent(in):: vrespc(nv1,nv2) 297 real(dp) ::frskerker1__dpf(nv1,nv2) 298 !Local variables------------------------------- 299 300 real(dp) :: buffer1(nv1,nv2),buffer2(nv1,nv2) 301 ! ************************************************************************* 302 303 if(ok) then 304 buffer1=vrespc 305 call laplacian(gprimd,mpi_enreg_ptr,nfft,nspden,ngfft,dtset_ptr%paral_kgb,& 306 & rdfuncr=buffer1,laplacerdfuncr=buffer2,g2cart_in=g2cart) 307 frskerker1__dpf(:,:)= vrespc(:,:)-deltaW-((dielng)**2)*buffer2(:,:) 308 else 309 frskerker1__dpf = zero 310 end if 311 end function frskerker1__dpf 312 313 end module frskerker1
frskerker1/frskerker1__end [ Functions ]
[ Top ] [ frskerker1 ] [ Functions ]
NAME
frskerker1__end
FUNCTION
ending subroutine deallocate memory areas
INPUTS
OUTPUT
PARENTS
prcrskerker1
CHILDREN
laplacian
SOURCE
139 subroutine frskerker1__end() 140 use defs_basis 141 use defs_abitypes 142 143 !This section has been created automatically by the script Abilint (TD). 144 !Do not modify the following lines by hand. 145 #undef ABI_FUNC 146 #define ABI_FUNC 'frskerker1__end' 147 !End of the abilint section 148 149 implicit none 150 !Arguments ------------------------------------ 151 152 !Local variables------------------------------- 153 154 ! ************************************************************************* 155 if(ok) then 156 ! ! set ok to false which prevent using the pf and dpf 157 ok = .false. 158 ! ! free memory 159 ABI_DEALLOCATE(deltaW) 160 ABI_DEALLOCATE(mat) 161 ABI_DEALLOCATE(g2cart) 162 end if 163 end subroutine frskerker1__end
frskerker1/frskerker1__init [ Functions ]
[ Top ] [ frskerker1 ] [ Functions ]
NAME
frskerker1__init
FUNCTION
initialisation subroutine Copy every variables required for the energy calculation Allocate the required memory
INPUTS
OUTPUT
PARENTS
prcrskerker1
CHILDREN
laplacian
SOURCE
77 subroutine frskerker1__init(dtset_in,mpi_enreg_in,nfft_in,ngfft_in,nspden_in,dielng_in,deltaW_in,gprimd_in,mat_in,g2cart_in ) 78 79 80 !This section has been created automatically by the script Abilint (TD). 81 !Do not modify the following lines by hand. 82 #undef ABI_FUNC 83 #define ABI_FUNC 'frskerker1__init' 84 !End of the abilint section 85 86 implicit none 87 !Arguments ------------------------------------ 88 integer,intent(in) :: nfft_in,ngfft_in(18),nspden_in 89 real(dp),intent(in) :: deltaW_in(nfft_in,nspden_in),mat_in(nfft_in,nspden_in),g2cart_in(nfft_in) 90 real(dp),dimension(3,3),intent(in) :: gprimd_in 91 real(dp),intent(in) :: dielng_in 92 type(dataset_type),target,intent(in) :: dtset_in 93 type(MPI_type),target,intent(in) :: mpi_enreg_in 94 !Local variables------------------------------- 95 96 ! ************************************************************************* 97 98 ! !allocation and data transfer 99 ! !Thought it would have been more logical to use the privates intrinsic of the module as 100 ! !input variables it seems that it is not possible... 101 if(.not.ok) then 102 dtset_ptr => dtset_in 103 mpi_enreg_ptr => mpi_enreg_in 104 nspden=nspden_in 105 ngfft=ngfft_in 106 nfft=nfft_in 107 ABI_ALLOCATE(deltaW,(size(deltaW_in,1),size(deltaW_in,2))) 108 ABI_ALLOCATE(mat,(size(mat_in,1),size(mat_in,2))) 109 ABI_ALLOCATE(g2cart,(size(g2cart_in,1))) 110 deltaW=deltaW_in 111 dielng=dielng_in 112 mat=mat_in 113 gprimd=gprimd_in 114 g2cart=g2cart_in 115 ok = .true. 116 end if 117 end subroutine frskerker1__init
frskerker1/frskerker1__newvres [ Functions ]
[ Top ] [ 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
PARENTS
CHILDREN
laplacian
SOURCE
185 subroutine frskerker1__newvres(nv1,nv2,x, grad, vrespc) 186 187 188 !This section has been created automatically by the script Abilint (TD). 189 !Do not modify the following lines by hand. 190 #undef ABI_FUNC 191 #define ABI_FUNC 'frskerker1__newvres' 192 !End of the abilint section 193 194 implicit none 195 !Arguments ------------------------------------ 196 integer,intent(in) :: nv1,nv2 197 real(dp),intent(in):: x 198 real(dp),intent(inout)::grad(nv1,nv2) 199 real(dp),intent(inout)::vrespc(nv1,nv2) 200 !Local variables------------------------------- 201 202 ! ************************************************************************* 203 204 grad(:,:)=x*grad(:,:) 205 vrespc(:,:)=vrespc(:,:)+grad(:,:) 206 end subroutine frskerker1__newvres
frskerker1/frskerker1__pf [ Functions ]
[ Top ] [ frskerker1 ] [ Functions ]
NAME
frskerker1__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
227 function frskerker1__pf(nv1,nv2,vrespc) 228 229 230 !This section has been created automatically by the script Abilint (TD). 231 !Do not modify the following lines by hand. 232 #undef ABI_FUNC 233 #define ABI_FUNC 'frskerker1__pf' 234 use interfaces_56_recipspace 235 use interfaces_62_cg_noabirule 236 !End of the abilint section 237 238 implicit none 239 !Arguments ------------------------------------ 240 integer,intent(in) :: nv1,nv2 241 real(dp),intent(in)::vrespc(nv1,nv2) 242 real(dp) ::frskerker1__pf 243 !Local variables------------------------------- 244 245 real(dp):: buffer1(nv1,nv2),buffer2(nv1,nv2) 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,& 252 & rdfuncr=buffer1,laplacerdfuncr=buffer2,g2cart_in=g2cart) 253 buffer2(:,:)=(vrespc(:,:)-((dielng)**2)*buffer2(:,:)) * half -deltaW 254 frskerker1__pf=dotproduct(nv1,nv2,vrespc,buffer2) !*half-dotproduct(nv1,nv2,vrespc,deltaW) 255 else 256 frskerker1__pf=zero 257 end if 258 end function frskerker1__pf