TABLE OF CONTENTS
ABINIT/int_ang [ Functions ]
NAME
int_ang
FUNCTION
Evaluate angular part for <phi_i|nabla|phi_j> and <tphi_i|nabla|tphi_j>
COPYRIGHT
Copyright (C) 2005-2018 ABINIT group (VR,FJ,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 .
INPUTS
mpsang=1+ max. angular momentum
OUTPUT
ang_phipphj :: angular part for <phi_i|nabla|phi_j> and <tphi_i|nabla|tphi_j> ang_phipphj(i,j,1)=\int sin\theta cos\phi Si Sj d\omega ang_phipphj(i,j,2)=\int cos\theta cos\phi Si \frac{d}{d\theta}Sj d\Omega ang_phipphj(i,j,3)=\int -sin\phi Si \frac{d}{d\phi}Sj d\Omega ang_phipphj(i,j,4)=\int sin\theta sin\phi Si Sj d\Omega ang_phipphj(i,j,5)=\int cos\theta sin\phi Si \frac{d}{d\theta}Sj d\Omega ang_phipphj(i,j,6)=\int cos\phi Si \frac{d}{d\phi}Sj d\Omega ang_phipphj(i,j,7)=\int cos\theta Si Sj d\Omega ang_phipphj(i,j,8)=\int -sin\theta Si \frac{d}{d\theta}Sj d\Omega
PARENTS
optics_paw,optics_paw_core,pawnabla_init
CHILDREN
SOURCE
36 #if defined HAVE_CONFIG_H 37 #include "config.h" 38 #endif 39 40 #include "abi_common.h" 41 42 subroutine int_ang(ang_phipphj,mpsang) 43 44 use defs_basis 45 use m_profiling_abi 46 use m_errors 47 48 !This section has been created automatically by the script Abilint (TD). 49 !Do not modify the following lines by hand. 50 #undef ABI_FUNC 51 #define ABI_FUNC 'int_ang' 52 !End of the abilint section 53 54 implicit none 55 56 !Arguments ------------------------------------ 57 !scalars 58 integer,intent(in) :: mpsang 59 !arrays 60 real(dp),intent(out) :: ang_phipphj(mpsang**2,mpsang**2,8) 61 62 !Local variables------------------------------- 63 character(len=500) :: message 64 real(dp) :: ang_phipphj_tmp(16,16,8) 65 66 ! ************************************************************************ 67 68 if (mpsang>4) then 69 message = ' Not designed for angular momentum greater than 3 !' 70 MSG_ERROR(message) 71 end if 72 73 !8 angular integrals for l=0..3, m=-l..+l 74 !ang_phipphj(1,4,1)=\frac{1}{\sqrt{3}} 75 !ang_phipphj(2,5,1)=\frac{1}{\sqrt{5}} 76 !ang_phipphj(3,8,1)=\frac{1}{\sqrt{5}} 77 !ang_phipphj(4,1,1)=\frac{1}{\sqrt{3}} 78 !ang_phipphj(4,7,1)=-\frac{1}{\sqrt{15}} 79 !ang_phipphj(4,9,1)=\frac{1}{\sqrt{5}} 80 !ang_phipphj(5,2,1)=\frac{1}{\sqrt{5}} 81 !ang_phipphj(5,10,1)=\sqrt{\frac{3}{14}} 82 !ang_phipphj(5,12,1)=-\frac{1}{\sqrt{70}} 83 !ang_phipphj(6,11,1)=\frac{1}{\sqrt{7}} 84 !ang_phipphj(7,4,1)=-\frac{1}{\sqrt{15}} 85 !ang_phipphj(7,14,1)=\sqrt{\frac{6}{35}} 86 !ang_phipphj(8,3,1)=\frac{1}{\sqrt{5}} 87 !ang_phipphj(8,13,1)=-\sqrt{\frac{3}{35}} 88 !ang_phipphj(8,15,1)=\frac{1}{\sqrt{7}} 89 !ang_phipphj(9,4,1)=\frac{1}{\sqrt{5}} 90 !ang_phipphj(9,14,1)=-\frac{1}{\sqrt{70}} 91 !ang_phipphj(9,16,1)=\sqrt{\frac{3}{14}} 92 !ang_phipphj(10,5,1)=\sqrt{\frac{3}{14}} 93 !ang_phipphj(11,6,1)=\frac{1}{\sqrt{7}} 94 !ang_phipphj(12,5,1)=-\frac{1}{\sqrt{70}} 95 !ang_phipphj(13,8,1)=-\sqrt{\frac{3}{35}} 96 !ang_phipphj(14,7,1)=\sqrt{\frac{6}{35}} 97 !ang_phipphj(14,9,1)=-\frac{1}{\sqrt{70}} 98 !ang_phipphj(15,8,1)=\frac{1}{\sqrt{7}} 99 !ang_phipphj(16,9,1)=\sqrt{\frac{3}{14}} 100 !ang_phipphj(1,4,2)=\frac{1}{2 \sqrt{3}} 101 !ang_phipphj(1,14,2)=-\frac{\sqrt{\frac{7}{6}}}{2} 102 !ang_phipphj(2,5,2)=\frac{1}{2 \sqrt{5}} 103 !ang_phipphj(3,8,2)=\frac{1}{2 \sqrt{5}} 104 !ang_phipphj(4,7,2)=-\sqrt{\frac{3}{5}} 105 !ang_phipphj(4,9,2)=\frac{1}{2 \sqrt{5}} 106 !ang_phipphj(5,2,2)=\frac{1}{4 \sqrt{5}} 107 !ang_phipphj(5,10,2)=\frac{\sqrt{\frac{3}{14}}}{2} 108 !ang_phipphj(5,12,2)=-2 \sqrt{\frac{2}{35}} 109 !ang_phipphj(6,11,2)=\frac{1}{2 \sqrt{7}} 110 !ang_phipphj(7,4,2)=\frac{1}{\sqrt{15}} 111 !ang_phipphj(7,14,2)=\frac{13}{2 \sqrt{210}} 112 !ang_phipphj(8,3,2)=-\frac{1}{\sqrt{5}} 113 !ang_phipphj(8,13,2)=-4 \sqrt{\frac{3}{35}} 114 !ang_phipphj(8,15,2)=\frac{1}{2 \sqrt{7}} 115 !ang_phipphj(9,4,2)=\frac{1}{4 \sqrt{5}} 116 !ang_phipphj(9,14,2)=-2 \sqrt{\frac{2}{35}} 117 !ang_phipphj(9,16,2)=\frac{\sqrt{\frac{3}{14}}}{2} 118 !ang_phipphj(10,5,2)=\frac{1}{\sqrt{42}} 119 !ang_phipphj(11,6,2)=-\frac{1}{4 \sqrt{7}} 120 !ang_phipphj(12,5,2)=\sqrt{\frac{2}{35}} 121 !ang_phipphj(13,8,2)=2 \sqrt{\frac{3}{35}} 122 !ang_phipphj(14,7,2)=-2 \sqrt{\frac{6}{35}} 123 !ang_phipphj(14,9,2)=\sqrt{\frac{2}{35}} 124 !ang_phipphj(15,8,2)=-\frac{1}{4 \sqrt{7}} 125 !ang_phipphj(16,9,2)=\frac{1}{\sqrt{42}} 126 !ang_phipphj(1,4,3)=\frac{\sqrt{3}}{2} 127 !ang_phipphj(1,14,3)=\frac{\sqrt{\frac{7}{6}}}{2} 128 !ang_phipphj(2,5,3)=\frac{\sqrt{5}}{2} 129 !ang_phipphj(3,8,3)=\frac{\sqrt{5}}{2} 130 !ang_phipphj(4,9,3)=\frac{\sqrt{5}}{2} 131 !ang_phipphj(5,2,3)=-\frac{\sqrt{5}}{4} 132 !ang_phipphj(5,10,3)=\frac{\sqrt{\frac{21}{2}}}{2} 133 !ang_phipphj(6,11,3)=\frac{\sqrt{7}}{2} 134 !ang_phipphj(7,14,3)=\frac{\sqrt{\frac{35}{6}}}{2} 135 !ang_phipphj(8,15,3)=\frac{\sqrt{7}}{2} 136 !ang_phipphj(9,4,3)=-\frac{\sqrt{5}}{4} 137 !ang_phipphj(9,16,3)=\frac{\sqrt{\frac{21}{2}}}{2} 138 !ang_phipphj(10,5,3)=-\sqrt{\frac{7}{6}} 139 !ang_phipphj(11,6,3)=-\frac{\sqrt{7}}{4} 140 !ang_phipphj(15,8,3)=-\frac{\sqrt{7}}{4} 141 !ang_phipphj(16,9,3)=-\sqrt{\frac{7}{6}} 142 !ang_phipphj(1,2,4)=\frac{1}{\sqrt{3}} 143 !ang_phipphj(2,1,4)=\frac{1}{\sqrt{3}} 144 !ang_phipphj(2,7,4)=-\frac{1}{\sqrt{15}} 145 !ang_phipphj(2,9,4)=-\frac{1}{\sqrt{5}} 146 !ang_phipphj(3,6,4)=\frac{1}{\sqrt{5}} 147 !ang_phipphj(4,5,4)=\frac{1}{\sqrt{5}} 148 !ang_phipphj(5,4,4)=\frac{1}{\sqrt{5}} 149 !ang_phipphj(5,14,4)=-\frac{1}{\sqrt{70}} 150 !ang_phipphj(5,16,4)=-\sqrt{\frac{3}{14}} 151 !ang_phipphj(6,3,4)=\frac{1}{\sqrt{5}} 152 !ang_phipphj(6,13,4)=-\sqrt{\frac{3}{35}} 153 !ang_phipphj(6,15,4)=-\frac{1}{\sqrt{7}} 154 !ang_phipphj(7,2,4)=-\frac{1}{\sqrt{15}} 155 !ang_phipphj(7,12,4)=\sqrt{\frac{6}{35}} 156 !ang_phipphj(8,11,4)=\frac{1}{\sqrt{7}} 157 !ang_phipphj(9,2,4)=-\frac{1}{\sqrt{5}} 158 !ang_phipphj(9,10,4)=\sqrt{\frac{3}{14}} 159 !ang_phipphj(9,12,4)=\frac{1}{\sqrt{70}} 160 !ang_phipphj(10,9,4)=\sqrt{\frac{3}{14}} 161 !ang_phipphj(11,8,4)=\frac{1}{\sqrt{7}} 162 !ang_phipphj(12,7,4)=\sqrt{\frac{6}{35}} 163 !ang_phipphj(12,9,4)=\frac{1}{\sqrt{70}} 164 !ang_phipphj(13,6,4)=-\sqrt{\frac{3}{35}} 165 !ang_phipphj(14,5,4)=-\frac{1}{\sqrt{70}} 166 !ang_phipphj(15,6,4)=-\frac{1}{\sqrt{7}} 167 !ang_phipphj(16,5,4)=-\sqrt{\frac{3}{14}} 168 !ang_phipphj(1,2,5)=\frac{1}{2 \sqrt{3}} 169 !ang_phipphj(1,12,5)=-\frac{\sqrt{\frac{7}{6}}}{2} 170 !ang_phipphj(2,7,5)=-\sqrt{\frac{3}{5}} 171 !ang_phipphj(2,9,5)=-\frac{1}{2 \sqrt{5}} 172 !ang_phipphj(3,6,5)=\frac{1}{2 \sqrt{5}} 173 !ang_phipphj(4,5,5)=\frac{1}{2 \sqrt{5}} 174 !ang_phipphj(5,4,5)=\frac{1}{4 \sqrt{5}} 175 !ang_phipphj(5,14,5)=-2 \sqrt{\frac{2}{35}} 176 !ang_phipphj(5,16,5)=-\frac{\sqrt{\frac{3}{14}}}{2} 177 !ang_phipphj(6,3,5)=-\frac{1}{\sqrt{5}} 178 !ang_phipphj(6,13,5)=-4 \sqrt{\frac{3}{35}} 179 !ang_phipphj(6,15,5)=-\frac{1}{2 \sqrt{7}} 180 !ang_phipphj(7,2,5)=\frac{1}{\sqrt{15}} 181 !ang_phipphj(7,12,5)=\frac{13}{2 \sqrt{210}} 182 !ang_phipphj(8,11,5)=\frac{1}{2 \sqrt{7}} 183 !ang_phipphj(9,2,5)=-\frac{1}{4 \sqrt{5}} 184 !ang_phipphj(9,10,5)=\frac{\sqrt{\frac{3}{14}}}{2} 185 !ang_phipphj(9,12,5)=2 \sqrt{\frac{2}{35}} 186 !ang_phipphj(10,9,5)=\frac{1}{\sqrt{42}} 187 !ang_phipphj(11,8,5)=-\frac{1}{4 \sqrt{7}} 188 !ang_phipphj(12,7,5)=-2 \sqrt{\frac{6}{35}} 189 !ang_phipphj(12,9,5)=-\sqrt{\frac{2}{35}} 190 !ang_phipphj(13,6,5)=2 \sqrt{\frac{3}{35}} 191 !ang_phipphj(14,5,5)=\sqrt{\frac{2}{35}} 192 !ang_phipphj(15,6,5)=\frac{1}{4 \sqrt{7}} 193 !ang_phipphj(16,5,5)=-\frac{1}{\sqrt{42}} 194 !ang_phipphj(1,2,6)=\frac{\sqrt{3}}{2} 195 !ang_phipphj(1,12,6)=\frac{\sqrt{\frac{7}{6}}}{2} 196 !ang_phipphj(2,9,6)=-\frac{\sqrt{5}}{2} 197 !ang_phipphj(3,6,6)=\frac{\sqrt{5}}{2} 198 !ang_phipphj(4,5,6)=\frac{\sqrt{5}}{2} 199 !ang_phipphj(5,4,6)=-\frac{\sqrt{5}}{4} 200 !ang_phipphj(5,16,6)=-\frac{\sqrt{\frac{21}{2}}}{2} 201 !ang_phipphj(6,15,6)=-\frac{\sqrt{7}}{2} 202 !ang_phipphj(7,12,6)=\frac{\sqrt{\frac{35}{6}}}{2} 203 !ang_phipphj(8,11,6)=\frac{\sqrt{7}}{2} 204 !ang_phipphj(9,2,6)=\frac{\sqrt{5}}{4} 205 !ang_phipphj(9,10,6)=\frac{\sqrt{\frac{21}{2}}}{2} 206 !ang_phipphj(10,9,6)=-\sqrt{\frac{7}{6}} 207 !ang_phipphj(11,8,6)=-\frac{\sqrt{7}}{4} 208 !ang_phipphj(15,6,6)=\frac{\sqrt{7}}{4} 209 !ang_phipphj(16,5,6)=\sqrt{\frac{7}{6}} 210 !ang_phipphj(1,3,7)=\frac{1}{\sqrt{3}} 211 !ang_phipphj(2,6,7)=\frac{1}{\sqrt{5}} 212 !ang_phipphj(3,1,7)=\frac{1}{\sqrt{3}} 213 !ang_phipphj(3,7,7)=\frac{2}{\sqrt{15}} 214 !ang_phipphj(4,8,7)=\frac{1}{\sqrt{5}} 215 !ang_phipphj(5,11,7)=\frac{1}{\sqrt{7}} 216 !ang_phipphj(6,2,7)=\frac{1}{\sqrt{5}} 217 !ang_phipphj(6,12,7)=2 \sqrt{\frac{2}{35}} 218 !ang_phipphj(7,3,7)=\frac{2}{\sqrt{15}} 219 !ang_phipphj(7,13,7)=\frac{3}{\sqrt{35}} 220 !ang_phipphj(8,4,7)=\frac{1}{\sqrt{5}} 221 !ang_phipphj(8,14,7)=2 \sqrt{\frac{2}{35}} 222 !ang_phipphj(9,15,7)=\frac{1}{\sqrt{7}} 223 !ang_phipphj(11,5,7)=\frac{1}{\sqrt{7}} 224 !ang_phipphj(12,6,7)=2 \sqrt{\frac{2}{35}} 225 !ang_phipphj(13,7,7)=\frac{3}{\sqrt{35}} 226 !ang_phipphj(14,8,7)=2 \sqrt{\frac{2}{35}} 227 !ang_phipphj(15,9,7)=\frac{1}{\sqrt{7}} 228 !ang_phipphj(1,3,8)=\frac{2}{\sqrt{3}} 229 !ang_phipphj(2,6,8)=\frac{3}{\sqrt{5}} 230 !ang_phipphj(3,7,8)=2 \sqrt{\frac{3}{5}} 231 !ang_phipphj(4,8,8)=\frac{3}{\sqrt{5}} 232 !ang_phipphj(5,11,8)=\frac{4}{\sqrt{7}} 233 !ang_phipphj(6,2,8)=-\frac{1}{\sqrt{5}} 234 !ang_phipphj(6,12,8)=8 \sqrt{\frac{2}{35}} 235 !ang_phipphj(7,3,8)=-\frac{2}{\sqrt{15}} 236 !ang_phipphj(7,13,8)=\frac{12}{\sqrt{35}} 237 !ang_phipphj(8,4,8)=-\frac{1}{\sqrt{5}} 238 !ang_phipphj(8,14,8)=8 \sqrt{\frac{2}{35}} 239 !ang_phipphj(9,15,8)=\frac{4}{\sqrt{7}} 240 !ang_phipphj(11,5,8)=-\frac{2}{\sqrt{7}} 241 !ang_phipphj(12,6,8)=-4 \sqrt{\frac{2}{35}} 242 !ang_phipphj(13,7,8)=-\frac{6}{\sqrt{35}} 243 !ang_phipphj(14,8,8)=-4 \sqrt{\frac{2}{35}} 244 !ang_phipphj(15,9,8)=-\frac{2}{\sqrt{7}} 245 246 247 ang_phipphj_tmp=zero 248 ! 249 ang_phipphj_tmp(1,4,1)=0.57735026918962576451_dp 250 ang_phipphj_tmp(2,5,1)=0.44721359549995793928_dp 251 ang_phipphj_tmp(3,8,1)=0.44721359549995793928_dp 252 ang_phipphj_tmp(4,1,1)=0.57735026918962576451_dp 253 ang_phipphj_tmp(4,7,1)=-0.25819888974716112568_dp 254 ang_phipphj_tmp(4,9,1)=0.44721359549995793928_dp 255 ang_phipphj_tmp(5,2,1)=0.44721359549995793928_dp 256 ang_phipphj_tmp(5,10,1)=0.46291004988627573078_dp 257 ang_phipphj_tmp(5,12,1)=-0.11952286093343936400_dp 258 ang_phipphj_tmp(6,11,1)=0.37796447300922722721_dp 259 ang_phipphj_tmp(7,4,1)=-0.25819888974716112568_dp 260 ang_phipphj_tmp(7,14,1)=0.41403933560541253068_dp 261 ang_phipphj_tmp(8,3,1)=0.44721359549995793928_dp 262 ang_phipphj_tmp(8,13,1)=-0.29277002188455995381_dp 263 ang_phipphj_tmp(8,15,1)=0.37796447300922722721_dp 264 ang_phipphj_tmp(9,4,1)=0.44721359549995793928_dp 265 ang_phipphj_tmp(9,14,1)=-0.11952286093343936400_dp 266 ang_phipphj_tmp(9,16,1)=0.46291004988627573078_dp 267 ang_phipphj_tmp(10,5,1)=0.46291004988627573078_dp 268 ang_phipphj_tmp(11,6,1)=0.37796447300922722721_dp 269 ang_phipphj_tmp(12,5,1)=-0.11952286093343936400_dp 270 ang_phipphj_tmp(13,8,1)=-0.29277002188455995381_dp 271 ang_phipphj_tmp(14,7,1)=0.41403933560541253068_dp 272 ang_phipphj_tmp(14,9,1)=-0.11952286093343936400_dp 273 ang_phipphj_tmp(15,8,1)=0.37796447300922722721_dp 274 ang_phipphj_tmp(16,9,1)=0.46291004988627573078_dp 275 ! 276 ang_phipphj_tmp(1,4,2)=0.28867513459481288225_dp 277 ang_phipphj_tmp(1,14,2)=-0.54006172486732168591_dp 278 ang_phipphj_tmp(2,5,2)=0.22360679774997896964_dp 279 ang_phipphj_tmp(3,8,2)=0.22360679774997896964_dp 280 ang_phipphj_tmp(4,7,2)=-0.77459666924148337704_dp 281 ang_phipphj_tmp(4,9,2)=0.22360679774997896964_dp 282 ang_phipphj_tmp(5,2,2)=0.11180339887498948482_dp 283 ang_phipphj_tmp(5,10,2)=0.23145502494313786539_dp 284 ang_phipphj_tmp(5,12,2)=-0.47809144373375745599_dp 285 ang_phipphj_tmp(6,11,2)=0.18898223650461361361_dp 286 ang_phipphj_tmp(7,4,2)=0.25819888974716112568_dp 287 ang_phipphj_tmp(7,14,2)=0.44854261357253024157_dp 288 ang_phipphj_tmp(8,3,2)=-0.44721359549995793928_dp 289 ang_phipphj_tmp(8,13,2)=-1.1710800875382398152_dp 290 ang_phipphj_tmp(8,15,2)=0.18898223650461361361_dp 291 ang_phipphj_tmp(9,4,2)=0.11180339887498948482_dp 292 ang_phipphj_tmp(9,14,2)=-0.47809144373375745599_dp 293 ang_phipphj_tmp(9,16,2)=0.23145502494313786539_dp 294 ang_phipphj_tmp(10,5,2)=0.15430334996209191026_dp 295 ang_phipphj_tmp(11,6,2)=-0.094491118252306806804_dp 296 ang_phipphj_tmp(12,5,2)=0.23904572186687872799_dp 297 ang_phipphj_tmp(13,8,2)=0.58554004376911990761_dp 298 ang_phipphj_tmp(14,7,2)=-0.82807867121082506136_dp 299 ang_phipphj_tmp(14,9,2)=0.23904572186687872799_dp 300 ang_phipphj_tmp(15,8,2)=-0.094491118252306806804_dp 301 ang_phipphj_tmp(16,9,2)=0.15430334996209191026_dp 302 ! 303 ang_phipphj_tmp(1,4,3)=0.86602540378443864676_dp 304 ang_phipphj_tmp(1,14,3)=0.54006172486732168591_dp 305 ang_phipphj_tmp(2,5,3)=1.1180339887498948482_dp 306 ang_phipphj_tmp(3,8,3)=1.1180339887498948482_dp 307 ang_phipphj_tmp(4,9,3)=1.1180339887498948482_dp 308 ang_phipphj_tmp(5,2,3)=-0.55901699437494742410_dp 309 ang_phipphj_tmp(5,10,3)=1.6201851746019650577_dp 310 ang_phipphj_tmp(6,11,3)=1.3228756555322952953_dp 311 ang_phipphj_tmp(7,14,3)=1.2076147288491198811_dp 312 ang_phipphj_tmp(8,15,3)=1.3228756555322952953_dp 313 ang_phipphj_tmp(9,4,3)=-0.55901699437494742410_dp 314 ang_phipphj_tmp(9,16,3)=1.6201851746019650577_dp 315 ang_phipphj_tmp(10,5,3)=-1.0801234497346433718_dp 316 ang_phipphj_tmp(11,6,3)=-0.66143782776614764763_dp 317 ang_phipphj_tmp(15,8,3)=-0.66143782776614764763_dp 318 ang_phipphj_tmp(16,9,3)=-1.0801234497346433718_dp 319 ! 320 ang_phipphj_tmp(1,2,4)=0.57735026918962576451_dp 321 ang_phipphj_tmp(2,1,4)=0.57735026918962576451_dp 322 ang_phipphj_tmp(2,7,4)=-0.25819888974716112568_dp 323 ang_phipphj_tmp(2,9,4)=-0.44721359549995793928_dp 324 ang_phipphj_tmp(3,6,4)=0.44721359549995793928_dp 325 ang_phipphj_tmp(4,5,4)=0.44721359549995793928_dp 326 ang_phipphj_tmp(5,4,4)=0.44721359549995793928_dp 327 ang_phipphj_tmp(5,14,4)=-0.11952286093343936400_dp 328 ang_phipphj_tmp(5,16,4)=-0.46291004988627573078_dp 329 ang_phipphj_tmp(6,3,4)=0.44721359549995793928_dp 330 ang_phipphj_tmp(6,13,4)=-0.29277002188455995381_dp 331 ang_phipphj_tmp(6,15,4)=-0.37796447300922722721_dp 332 ang_phipphj_tmp(7,2,4)=-0.25819888974716112568_dp 333 ang_phipphj_tmp(7,12,4)=0.41403933560541253068_dp 334 ang_phipphj_tmp(8,11,4)=0.37796447300922722721_dp 335 ang_phipphj_tmp(9,2,4)=-0.44721359549995793928_dp 336 ang_phipphj_tmp(9,10,4)=0.46291004988627573078_dp 337 ang_phipphj_tmp(9,12,4)=0.11952286093343936400_dp 338 ang_phipphj_tmp(10,9,4)=0.46291004988627573078_dp 339 ang_phipphj_tmp(11,8,4)=0.37796447300922722721_dp 340 ang_phipphj_tmp(12,7,4)=0.41403933560541253068_dp 341 ang_phipphj_tmp(12,9,4)=0.11952286093343936400_dp 342 ang_phipphj_tmp(13,6,4)=-0.29277002188455995381_dp 343 ang_phipphj_tmp(14,5,4)=-0.11952286093343936400_dp 344 ang_phipphj_tmp(15,6,4)=-0.37796447300922722721_dp 345 ang_phipphj_tmp(16,5,4)=-0.46291004988627573078_dp 346 ! 347 ang_phipphj_tmp(1,2,5)=0.28867513459481288225_dp 348 ang_phipphj_tmp(1,12,5)=-0.54006172486732168591_dp 349 ang_phipphj_tmp(2,7,5)=-0.77459666924148337704_dp 350 ang_phipphj_tmp(2,9,5)=-0.22360679774997896964_dp 351 ang_phipphj_tmp(3,6,5)=0.22360679774997896964_dp 352 ang_phipphj_tmp(4,5,5)=0.22360679774997896964_dp 353 ang_phipphj_tmp(5,4,5)=0.11180339887498948482_dp 354 ang_phipphj_tmp(5,14,5)=-0.47809144373375745599_dp 355 ang_phipphj_tmp(5,16,5)=-0.23145502494313786539_dp 356 ang_phipphj_tmp(6,3,5)=-0.44721359549995793928_dp 357 ang_phipphj_tmp(6,13,5)=-1.1710800875382398152_dp 358 ang_phipphj_tmp(6,15,5)=-0.18898223650461361361_dp 359 ang_phipphj_tmp(7,2,5)=0.25819888974716112568_dp 360 ang_phipphj_tmp(7,12,5)=0.44854261357253024157_dp 361 ang_phipphj_tmp(8,11,5)=0.18898223650461361361_dp 362 ang_phipphj_tmp(9,2,5)=-0.11180339887498948482_dp 363 ang_phipphj_tmp(9,10,5)=0.23145502494313786539_dp 364 ang_phipphj_tmp(9,12,5)=0.47809144373375745599_dp 365 ang_phipphj_tmp(10,9,5)=0.15430334996209191026_dp 366 ang_phipphj_tmp(11,8,5)=-0.094491118252306806804_dp 367 ang_phipphj_tmp(12,7,5)=-0.82807867121082506136_dp 368 ang_phipphj_tmp(12,9,5)=-0.23904572186687872799_dp 369 ang_phipphj_tmp(13,6,5)=0.58554004376911990761_dp 370 ang_phipphj_tmp(14,5,5)=0.23904572186687872799_dp 371 ang_phipphj_tmp(15,6,5)=0.094491118252306806804_dp 372 ang_phipphj_tmp(16,5,5)=-0.15430334996209191026_dp 373 ! 374 ang_phipphj_tmp(1,2,6)=0.86602540378443864676_dp 375 ang_phipphj_tmp(1,12,6)=0.54006172486732168591_dp 376 ang_phipphj_tmp(2,9,6)=-1.1180339887498948482_dp 377 ang_phipphj_tmp(3,6,6)=1.1180339887498948482_dp 378 ang_phipphj_tmp(4,5,6)=1.1180339887498948482_dp 379 ang_phipphj_tmp(5,4,6)=-0.55901699437494742410_dp 380 ang_phipphj_tmp(5,16,6)=-1.6201851746019650577_dp 381 ang_phipphj_tmp(6,15,6)=-1.3228756555322952953_dp 382 ang_phipphj_tmp(7,12,6)=1.2076147288491198811_dp 383 ang_phipphj_tmp(8,11,6)=1.3228756555322952953_dp 384 ang_phipphj_tmp(9,2,6)=0.55901699437494742410_dp 385 ang_phipphj_tmp(9,10,6)=1.6201851746019650577_dp 386 ang_phipphj_tmp(10,9,6)=-1.0801234497346433718_dp 387 ang_phipphj_tmp(11,8,6)=-0.66143782776614764763_dp 388 ang_phipphj_tmp(15,6,6)=0.66143782776614764763_dp 389 ang_phipphj_tmp(16,5,6)=1.0801234497346433718_dp 390 ! 391 ang_phipphj_tmp(1,3,7)=0.57735026918962576451_dp 392 ang_phipphj_tmp(2,6,7)=0.44721359549995793928_dp 393 ang_phipphj_tmp(3,1,7)=0.57735026918962576451_dp 394 ang_phipphj_tmp(3,7,7)=0.51639777949432225136_dp 395 ang_phipphj_tmp(4,8,7)=0.44721359549995793928_dp 396 ang_phipphj_tmp(5,11,7)=0.37796447300922722721_dp 397 ang_phipphj_tmp(6,2,7)=0.44721359549995793928_dp 398 ang_phipphj_tmp(6,12,7)=0.47809144373375745599_dp 399 ang_phipphj_tmp(7,3,7)=0.51639777949432225136_dp 400 ang_phipphj_tmp(7,13,7)=0.50709255283710994651_dp 401 ang_phipphj_tmp(8,4,7)=0.44721359549995793928_dp 402 ang_phipphj_tmp(8,14,7)=0.47809144373375745599_dp 403 ang_phipphj_tmp(9,15,7)=0.37796447300922722721_dp 404 ang_phipphj_tmp(11,5,7)=0.37796447300922722721_dp 405 ang_phipphj_tmp(12,6,7)=0.47809144373375745599_dp 406 ang_phipphj_tmp(13,7,7)=0.50709255283710994651_dp 407 ang_phipphj_tmp(14,8,7)=0.47809144373375745599_dp 408 ang_phipphj_tmp(15,9,7)=0.37796447300922722721_dp 409 ! 410 ang_phipphj_tmp(1,3,8)=1.1547005383792515290_dp 411 ang_phipphj_tmp(2,6,8)=1.3416407864998738178_dp 412 ang_phipphj_tmp(3,7,8)=1.5491933384829667541_dp 413 ang_phipphj_tmp(4,8,8)=1.3416407864998738178_dp 414 ang_phipphj_tmp(5,11,8)=1.5118578920369089089_dp 415 ang_phipphj_tmp(6,2,8)=-0.44721359549995793928_dp 416 ang_phipphj_tmp(6,12,8)=1.9123657749350298240_dp 417 ang_phipphj_tmp(7,3,8)=-0.51639777949432225136_dp 418 ang_phipphj_tmp(7,13,8)=2.0283702113484397860_dp 419 ang_phipphj_tmp(8,4,8)=-0.44721359549995793928_dp 420 ang_phipphj_tmp(8,14,8)=1.9123657749350298240_dp 421 ang_phipphj_tmp(9,15,8)=1.5118578920369089089_dp 422 ang_phipphj_tmp(11,5,8)=-0.75592894601845445443_dp 423 ang_phipphj_tmp(12,6,8)=-0.95618288746751491198_dp 424 ang_phipphj_tmp(13,7,8)=-1.0141851056742198930_dp 425 ang_phipphj_tmp(14,8,8)=-0.95618288746751491198_dp 426 ang_phipphj_tmp(15,9,8)=-0.75592894601845445443_dp 427 428 ang_phipphj(:,:,:)=ang_phipphj_tmp(1:mpsang**2,1:mpsang**2,:) 429 430 431 end subroutine int_ang