TABLE OF CONTENTS
ABINIT/mati3inv [ Functions ]
NAME
mati3inv
FUNCTION
Invert and transpose orthogonal 3x3 matrix of INTEGER elements.
COPYRIGHT
Copyright (C) 1998-2018 ABINIT group (DCA, XG, GMR) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
INPUTS
mm = integer matrix to be inverted
OUTPUT
mit = inverse of mm input matrix
NOTES
Used for symmetry operations. This routine applies to ORTHOGONAL matrices only. Since these form a group, inverses are also integer arrays. Returned array is TRANSPOSE of inverse, as needed. Note use of integer arithmetic.
PARENTS
cg_rotate,chkgrp,classify_bands,debug_tools,dfpt_nstdy,get_full_kgrid get_npert_rbz,getkgrid,ingeo,m_ab7_symmetry,m_crystal,m_ddb,m_dvdb m_dynmat,m_fft_mesh,m_fock,m_ptgroups,m_tdep_sym,matpointsym memory_eval,optic,read_gkk,setsym,strainsym,thmeig,wfconv
CHILDREN
SOURCE
38 #if defined HAVE_CONFIG_H 39 #include "config.h" 40 #endif 41 42 #include "abi_common.h" 43 44 subroutine mati3inv(mm,mit) 45 46 use defs_basis 47 use m_profiling_abi 48 use m_errors 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 'mati3inv' 54 !End of the abilint section 55 56 implicit none 57 58 !Arguments ------------------------------------ 59 !arrays 60 integer,intent(in) :: mm(3,3) 61 integer,intent(out) :: mit(3,3) 62 63 !Local variables------------------------------- 64 !scalars 65 integer :: dd 66 character(len=500) :: message 67 !arrays 68 integer :: tt(3,3) 69 70 ! ************************************************************************* 71 72 tt(1,1) = mm(2,2) * mm(3,3) - mm(3,2) * mm(2,3) 73 tt(2,1) = mm(3,2) * mm(1,3) - mm(1,2) * mm(3,3) 74 tt(3,1) = mm(1,2) * mm(2,3) - mm(2,2) * mm(1,3) 75 tt(1,2) = mm(3,1) * mm(2,3) - mm(2,1) * mm(3,3) 76 tt(2,2) = mm(1,1) * mm(3,3) - mm(3,1) * mm(1,3) 77 tt(3,2) = mm(2,1) * mm(1,3) - mm(1,1) * mm(2,3) 78 tt(1,3) = mm(2,1) * mm(3,2) - mm(3,1) * mm(2,2) 79 tt(2,3) = mm(3,1) * mm(1,2) - mm(1,1) * mm(3,2) 80 tt(3,3) = mm(1,1) * mm(2,2) - mm(2,1) * mm(1,2) 81 dd = mm(1,1) * tt(1,1) + mm(2,1) * tt(2,1) + mm(3,1) * tt(3,1) 82 83 !Make sure matrix is not singular 84 if (dd/=0) then 85 mit(:,:)=tt(:,:)/dd 86 else 87 write(message, '(2a,2x,9i5,a)' )& 88 & ' Attempting to invert integer array',ch10,& 89 & mm(:,:),' ==> determinant is zero.' 90 MSG_BUG(message) 91 end if 92 93 !If matrix is orthogonal, determinant must be 1 or -1 94 if (abs(dd)/=1) then 95 write(message, '(2a,2x,9i5,a)' )& 96 & ' Absolute value of determinant should be one',ch10,& 97 & ' but determinant=',dd 98 MSG_BUG(message) 99 end if 100 101 102 end subroutine mati3inv