TABLE OF CONTENTS
ABINIT/dfptnl_doutput [ Functions ]
NAME
dfptnl_doutput
FUNCTION
Write the matrix of third-order derivatives to the output file and the DDB
COPYRIGHT
Copyright (C) 1999-2018 ABINIT group (MVeithen) 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/doc/developers/contributors.txt .
INPUTS
blkflg(3,mpert,3,mpert,3,mpert)= ( 1 if the element of the 3dte has been calculated ; 0 otherwise ) d3(2,3,mpert,3,mpert,3,mpert)= matrix of the 3DTE mpert =maximum number of ipert natom=Number of atoms ntypat=Number of type of atoms unddb = unit number for DDB output
NOTES
d3 holds the third-order derivatives before computing the permutations of the perturbations.
PARENTS
nonlinear
CHILDREN
ddb_free,ddb_malloc,ddb_write_blok,wrtout
SOURCE
37 #if defined HAVE_CONFIG_H 38 #include "config.h" 39 #endif 40 41 #include "abi_common.h" 42 43 44 subroutine dfptnl_doutput(blkflg,d3,mband,mpert,nkpt,natom,ntypat,unddb) 45 46 use defs_basis 47 use m_profiling_abi 48 use m_ddb 49 50 !This section has been created automatically by the script Abilint (TD). 51 !Do not modify the following lines by hand. 52 #undef ABI_FUNC 53 #define ABI_FUNC 'dfptnl_doutput' 54 use interfaces_14_hidewrite 55 !End of the abilint section 56 57 implicit none 58 59 !Arguments ------------------------------- 60 !scalars 61 integer,intent(in) :: mband,mpert,nkpt,unddb,natom,ntypat 62 !arrays 63 integer,intent(in) :: blkflg(3,mpert,3,mpert,3,mpert) 64 real(dp),intent(in) :: d3(2,3,mpert,3,mpert,3,mpert) 65 66 !Local variables ------------------------- 67 !scalars 68 integer :: choice,i1dir,i1pert,i2dir,i2pert,i3dir,i3pert,index,msize 69 character(len=500) :: message 70 type(ddb_type) :: ddb 71 72 !************************************************************************* 73 74 msize = 27*mpert*mpert*mpert 75 call ddb_malloc(ddb,msize,1,natom,ntypat) 76 77 choice = 2 78 79 ddb%typ = 3 80 ddb%nrm = one 81 ddb%qpt = zero ! this has to be changed in case anharmonic 82 !force constants have been computed 83 84 85 !Write blok of third-order derivatives to ouput file 86 87 write(message,'(a,a,a,a,a)')ch10,& 88 & ' Matrix of third-order derivatives (reduced coordinates)',ch10,& 89 & ' before computing the permutations of the perturbations',ch10 90 call wrtout(ab_out,message,'COLL') 91 92 write(ab_out,*)' j1 j2 j3 matrix element' 93 write(ab_out,*)' dir pert dir pert dir pert real part imaginary part' 94 95 do i1pert=1,mpert 96 do i1dir=1,3 97 98 do i2pert=1,mpert 99 do i2dir=1,3 100 101 do i3pert=1,mpert 102 do i3dir=1,3 103 104 index = i1dir + & 105 & 3*((i1pert-1)+mpert*((i2dir-1) + & 106 & 3*((i2pert-1)+mpert*((i3dir-1) + 3*(i3pert-1))))) 107 ddb%flg(index,1) = blkflg(i1dir,i1pert,i2dir,i2pert,i3dir,i3pert) 108 ddb%val(:,index,1)= d3(:,i1dir,i1pert,i2dir,i2pert,i3dir,i3pert) 109 110 if (blkflg(i1dir,i1pert,i2dir,i2pert,i3dir,i3pert)/=0) then 111 112 write(ab_out,'(3(i4,i5),2f22.10)')& 113 & i1dir,i1pert,i2dir,i2pert,i3dir,i3pert,& 114 & d3(:,i1dir,i1pert,i2dir,i2pert,i3dir,i3pert) 115 116 end if 117 118 end do 119 end do 120 121 end do 122 end do 123 124 end do 125 end do 126 127 !Write blok of third-order derivatives to DDB 128 call ddb_write_blok(ddb,1,choice,mband,mpert,msize,nkpt,unddb) 129 130 call ddb_free(ddb) 131 132 end subroutine dfptnl_doutput