Package solvcon :: Package kerpak :: Module elaslin
[hide private]
[frames] | no frames]

Source Code for Module solvcon.kerpak.elaslin

  1  # -*- coding: UTF-8 -*- 
  2  # 
  3  # Copyright (C) 2010-2011 Yung-Yu Chen <yyc@solvcon.net>. 
  4  # 
  5  # This program is free software; you can redistribute it and/or modify 
  6  # it under the terms of the GNU General Public License as published by 
  7  # the Free Software Foundation; either version 2 of the License, or 
  8  # (at your option) any later version. 
  9  # 
 10  # This program is distributed in the hope that it will be useful, 
 11  # but WITHOUT ANY WARRANTY; without even the implied warranty of 
 12  # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the 
 13  # GNU General Public License for more details. 
 14  # 
 15  # You should have received a copy of the GNU General Public License along 
 16  # with this program; if not, write to the Free Software Foundation, Inc., 
 17  # 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. 
 18   
 19  """ 
 20  Anisotropic velocity-stress equations solver using the linear CESE method. 
 21  """ 
 22   
 23  from solvcon.gendata import TypeNameRegistry 
 24  from solvcon.anchor import Anchor 
 25  from solvcon.hook import BlockHook 
 26  from solvcon.kerpak.cese import CeseBC 
 27  from solvcon.kerpak.lincese import (LinceseSolver, LinceseCase, 
 28      PlaneWaveSolution) 
29 30 ############################################################################### 31 # Solver. 32 ############################################################################### 33 34 -class ElaslinSolver(LinceseSolver):
35 """ 36 Basic elastic solver. 37 38 @ivar mtrldict: map from names to material objects. 39 @itype mtrldict: dict 40 @ivar mtrllist: list of all material objects. 41 @itype mtrllist: list 42 """ 43 from solvcon.dependency import getcdll 44 __clib_elaslin = { 45 2: getcdll('elaslin2d'), 46 3: getcdll('elaslin3d'), 47 } 48 del getcdll 49 @property
50 - def _clib_elaslin(self):
51 return self.__clib_elaslin[self.ndim]
52 @property
53 - def _gdlen_(self):
54 return 9 * 9 * self.ndim
55 - def __init__(self, *args, **kw):
56 super(ElaslinSolver, self).__init__(*args, **kw) 57 self.mtrldict = kw.pop('mtrldict', {}) 58 self.mtrllist = None
59 - def make_grpda(self):
60 self.mtrllist = self._build_mtrllist(self.grpnames, self.mtrldict) 61 for igrp in range(len(self.grpnames)): 62 mtrl = self.mtrllist[igrp] 63 jaco = self.grpda[igrp].reshape(self.neq, self.neq, self.ndim) 64 mjacos = mtrl.get_jacos() 65 for idm in range(self.ndim): 66 jaco[:,:,idm] = mjacos[idm,:,:]
67 @staticmethod
68 - def _build_mtrllist(grpnames, mtrldict):
69 """ 70 Build the material list out of the mapping dict. 71 72 @type grpnames: list 73 @param mtrldict: the map from names to material objects. 74 @type mtrldict: dict 75 @return: the list of material object. 76 @rtype: Material 77 """ 78 mtrllist = list() 79 default_mtuple = mtrldict.get(None, None) 80 for grpname in grpnames: 81 try: 82 mtrl = mtrldict.get(grpname, default_mtuple) 83 except KeyError, e: 84 args = e.args[:] 85 args.append('no material named %s in mtrldict'%grpname) 86 e.args = args 87 raise 88 mtrllist.append(mtrl) 89 return mtrllist
90
91 ############################################################################### 92 # Case. 93 ############################################################################### 94 95 -class ElaslinCase(LinceseCase):
96 """ 97 Case for anisotropic elastic solids. 98 """ 99 defdict = { 100 'execution.neq': 9, 101 'solver.solvertype': ElaslinSolver, 102 'solver.mtrldict': dict, 103 }
104 - def make_solver_keywords(self):
105 kw = super(ElaslinCase, self).make_solver_keywords() 106 # setup material mapper. 107 kw['mtrldict'] = self.solver.mtrldict 108 return kw
109
110 ############################################################################### 111 # Boundary conditions. 112 ############################################################################### 113 114 -class ElaslinBC(CeseBC):
115 """ 116 Basic BC class for elastic problems. 117 """ 118 from solvcon.dependency import getcdll 119 __clib_elaslinb = { 120 2: getcdll('elaslinb2d'), 121 3: getcdll('elaslinb3d'), 122 } 123 del getcdll 124 @property
125 - def _clib_elaslinb(self):
126 return self.__clib_elaslinb[self.svr.ndim]
127
128 -class ElaslinTraction(ElaslinBC):
129 vnames = [ 130 'bfcsys', 'tau1', 'tau2', 'tau3', 'freq', 'phase', 131 ] 132 vdefaults = { 133 'bfcsys': 0.0, 134 'tau1': 0.0, 'tau2': 0.0, 'tau3': 0.0, 'freq': 0.0, 'phase': 0.0, 135 } 136 _ghostgeom_ = 'compress'
137 - def soln(self):
138 from ctypes import byref, c_int 139 self._clib_boundcond.bound_traction_soln( 140 byref(self.svr.exd), 141 c_int(self.facn.shape[0]), 142 self.facn.ctypes._as_parameter_, 143 c_int(self.value.shape[1]), 144 self.value.ctypes._as_parameter_, 145 )
146 - def dsoln(self):
147 from ctypes import byref, c_int 148 self._clib_boundcond.bound_traction_dsoln( 149 byref(self.svr.exd), 150 c_int(self.facn.shape[0]), 151 self.facn.ctypes._as_parameter_, 152 )
153
154 -class ElaslinTractionFree(ElaslinBC):
155 _ghostgeom_ = 'mirror'
156 - def soln(self):
157 from ctypes import byref, c_int 158 self._clib_boundcond.bound_traction_free_soln( 159 byref(self.svr.exd), 160 c_int(self.facn.shape[0]), 161 self.facn.ctypes._as_parameter_, 162 )
163 - def dsoln(self):
164 from ctypes import byref, c_int 165 self._clib_boundcond.bound_traction_free_dsoln( 166 byref(self.svr.exd), 167 c_int(self.facn.shape[0]), 168 self.facn.ctypes._as_parameter_, 169 )
170
171 -class ElaslinTractionFree2(ElaslinBC):
172 _ghostgeom_ = 'mirror'
173 - def soln(self):
174 from ctypes import byref, c_int 175 self._clib_boundcond.bound_traction_free2_soln( 176 byref(self.svr.exd), 177 c_int(self.facn.shape[0]), 178 self.facn.ctypes._as_parameter_, 179 )
180 - def dsoln(self):
181 from ctypes import byref, c_int 182 self._clib_boundcond.bound_traction_free2_dsoln( 183 byref(self.svr.exd), 184 c_int(self.facn.shape[0]), 185 self.facn.ctypes._as_parameter_, 186 )
187
188 ################################################################################ 189 # Plane wave solution. 190 ################################################################################ 191 192 -class ElaslinPWSolution(PlaneWaveSolution):
193 - def _calc_eigen(self, **kw):
194 from numpy import sqrt 195 from numpy.linalg import eig 196 wvec = kw['wvec'] 197 mtrl = kw['mtrl'] 198 idx = kw['idx'] 199 nml = wvec/sqrt((wvec**2).sum()) 200 jacos = mtrl.get_jacos() 201 jaco = jacos[0] * nml[0] 202 for idm in range(1, len(nml)): 203 jaco += jacos[idm] * nml[idm] 204 evl, evc = eig(jaco) 205 srt = evl.argsort() 206 evl = evl[srt[idx]].real 207 evc = evc[:,srt[idx]].real 208 evc *= evc[0]/abs(evc[0]+1.e-200) 209 return evl, evc
210
211 ################################################################################ 212 # Anchor. 213 ################################################################################ 214 215 -class ElaslinOAnchor(Anchor):
216 """ 217 Calculate total energy, i.e., the summation of kinetic energy and strain 218 energy. 219 """
220 - def _calculate_physics(self):
221 from ctypes import byref 222 from numpy import empty 223 from numpy.linalg import inv 224 svr = self.svr 225 # input arrays. 226 rhos = empty(svr.ngroup, dtype=svr.fpdtype) 227 comps = empty((svr.ngroup, 6, 6), dtype=svr.fpdtype) 228 for igp in range(svr.ngroup): 229 mtrl = svr.mtrllist[igp] 230 rhos[igp] = mtrl.rho 231 comps[igp,:,:] = inv(mtrl.stiff).T 232 # output arrays. 233 svr._clib_elaslin.calc_energy( 234 byref(svr.exd), 235 rhos.ctypes._as_parameter_, 236 comps.ctypes._as_parameter_, 237 svr.der['energy'].ctypes._as_parameter_, 238 )
239 - def provide(self):
240 from numpy import empty 241 svr = self.svr 242 svr.der['energy'] = empty(svr.ngstcell+svr.ncell, dtype=svr.fpdtype) 243 self._calculate_physics()
244 - def postfull(self):
245 self._calculate_physics()
246 247 ################################################################################ 248 # Material definition. 249 ################################################################################ 250 251 mltregy = TypeNameRegistry() # registry singleton.
252 -class MaterialMeta(type):
253 """ 254 Meta class for material class. 255 """
256 - def __new__(cls, name, bases, namespace):
257 newcls = super(MaterialMeta, cls).__new__(cls, name, bases, namespace) 258 # register. 259 mltregy.register(newcls) 260 return newcls
261
262 -class Material(object):
263 """ 264 Material properties. The constitutive relation needs not be symmetric. 265 266 @cvar _zeropoints_: list of tuples for indices where the content should be 267 zero. 268 @ctype _zeropoints_: list 269 @ivar rho: density 270 @ivar al: alpha angle. 271 @ivar be: beta angle. 272 @ivar ga: gamma angle. 273 @ivar origstiff: stiffness matrix in the crystal coordinate. 274 @ivar stiff: stiffness matrix in the transformed global coordinate. 275 """ 276 277 __metaclass__ = MaterialMeta 278 279 _zeropoints_ = [] 280 281 from numpy import array 282 K = array([ [ 283 [1, 0, 0, 0, 0, 0], 284 [0, 0, 0, 0, 0, 1], 285 [0, 0, 0, 0, 1, 0], 286 ], [ 287 [0, 0, 0, 0, 0, 1], 288 [0, 1, 0, 0, 0, 0], 289 [0, 0, 0, 1, 0, 0], 290 ], [ 291 [0, 0, 0, 0, 1, 0], 292 [0, 0, 0, 1, 0, 0], 293 [0, 0, 1, 0, 0, 0], 294 ], ], dtype='float64') 295 del array 296
297 - def __init__(self, *args, **kw):
298 from numpy import empty, dot 299 self.rho = kw.pop('rho') 300 self.al = kw.pop('al') 301 self.be = kw.pop('be') 302 self.ga = kw.pop('ga') 303 # set stiffness matrix. 304 origstiff = empty((6,6), dtype='float64') 305 origstiff.fill(0.0) 306 for key in kw.keys(): # becaues I pop out the key. 307 if len(key) == 4 and key[:2] == 'co': 308 try: 309 i = int(key[2])-1 310 j = int(key[3])-1 311 except: 312 continue 313 assert i < origstiff.shape[0] 314 assert j < origstiff.shape[1] 315 val = kw.pop(key) 316 origstiff[i,j] = val 317 self.origstiff = origstiff 318 # check for zeros. 319 self._check_origstiffzero(self.origstiff) 320 # compute the stiffness matrix in the transformed global coordinate 321 # system. 322 bondmat = self.get_bondmat() 323 self.stiff = dot(bondmat, dot(self.origstiff, bondmat.T)) 324 super(Material, self).__init__(*args, **kw)
325
326 - def __getattr__(self, key):
327 if len(key) == 4 and key[:2] == 'co': 328 i = int(key[2]) 329 j = int(key[3]) 330 if 1 <= i <= 6 and 1 <= j <= 6: 331 return self.origstiff[i-1,j-1] 332 elif len(key) == 3 and key[0] == 'c': 333 i = int(key[1]) 334 j = int(key[2]) 335 if 1 <= i <= 6 and 1 <= j <= 6: 336 return self.stiff[i-1,j-1] 337 else: 338 raise AttributeError
339
340 - def __str__(self):
341 from math import pi 342 return '[%s: al=%.2f be=%.2f ga=%.2f (deg)]' % (self.__class__.__name__, 343 self.al/(pi/180), self.be/(pi/180), self.ga/(pi/180))
344 345 @classmethod
346 - def _check_origstiffzero(cls, origstiff):
347 """ 348 Check for zero in original stiffness matrix. 349 350 @note: no assumed symmetry. 351 """ 352 for i, j in cls._zeropoints_: 353 assert origstiff[i,j] == 0.0
354
355 - def get_rotmat(self):
356 """ 357 Coordinate transformation matrix for three successive rotations through 358 the Euler angles. 359 360 @return: the transformation matrix. 361 @rtype: numpy.ndarray 362 """ 363 from numpy import array, cos, sin, dot 364 al = self.al; be = self.be; ga = self.ga 365 almat = array([ 366 [cos(al), sin(al), 0], 367 [-sin(al), cos(al), 0], 368 [0, 0, 1], 369 ], dtype='float64') 370 bemat = array([ 371 [1, 0, 0], 372 [0, cos(be), sin(be)], 373 [0, -sin(be), cos(be)], 374 ], dtype='float64') 375 gamat = array([ 376 [cos(ga), sin(ga), 0], 377 [-sin(ga), cos(ga), 0], 378 [0, 0, 1], 379 ], dtype='float64') 380 return dot(gamat, dot(bemat, almat))
381
382 - def get_bondmat(self):
383 """ 384 The Bond's matrix M as a shorthand of coordinate transformation for the 385 6-component stress vector. 386 387 @return: the Bond's matrix. 388 @rtype: numpy.ndarray 389 """ 390 from numpy import empty 391 rotmat = self.get_rotmat() 392 bond = empty((6,6), dtype='float64') 393 # upper left. 394 bond[:3,:3] = rotmat[:,:]**2 395 # upper right. 396 bond[0,3] = 2*rotmat[0,1]*rotmat[0,2] 397 bond[0,4] = 2*rotmat[0,2]*rotmat[0,0] 398 bond[0,5] = 2*rotmat[0,0]*rotmat[0,1] 399 bond[1,3] = 2*rotmat[1,1]*rotmat[1,2] 400 bond[1,4] = 2*rotmat[1,2]*rotmat[1,0] 401 bond[1,5] = 2*rotmat[1,0]*rotmat[1,1] 402 bond[2,3] = 2*rotmat[2,1]*rotmat[2,2] 403 bond[2,4] = 2*rotmat[2,2]*rotmat[2,0] 404 bond[2,5] = 2*rotmat[2,0]*rotmat[2,1] 405 # lower left. 406 bond[3,0] = rotmat[1,0]*rotmat[2,0] 407 bond[3,1] = rotmat[1,1]*rotmat[2,1] 408 bond[3,2] = rotmat[1,2]*rotmat[2,2] 409 bond[4,0] = rotmat[2,0]*rotmat[0,0] 410 bond[4,1] = rotmat[2,1]*rotmat[0,1] 411 bond[4,2] = rotmat[2,2]*rotmat[0,2] 412 bond[5,0] = rotmat[0,0]*rotmat[1,0] 413 bond[5,1] = rotmat[0,1]*rotmat[1,1] 414 bond[5,2] = rotmat[0,2]*rotmat[1,2] 415 # lower right. 416 bond[3,3] = rotmat[1,1]*rotmat[2,2] + rotmat[1,2]*rotmat[2,1] 417 bond[3,4] = rotmat[1,0]*rotmat[2,2] + rotmat[1,2]*rotmat[2,0] 418 bond[3,5] = rotmat[1,1]*rotmat[2,0] + rotmat[1,0]*rotmat[2,1] 419 bond[4,3] = rotmat[0,1]*rotmat[2,2] + rotmat[0,2]*rotmat[2,1] 420 bond[4,4] = rotmat[0,0]*rotmat[2,2] + rotmat[0,2]*rotmat[2,0] 421 bond[4,5] = rotmat[0,1]*rotmat[2,0] + rotmat[0,0]*rotmat[2,1] 422 bond[5,3] = rotmat[0,1]*rotmat[1,2] + rotmat[0,2]*rotmat[1,1] 423 bond[5,4] = rotmat[0,0]*rotmat[1,2] + rotmat[0,2]*rotmat[1,0] 424 bond[5,5] = rotmat[0,1]*rotmat[1,0] + rotmat[0,0]*rotmat[1,1] 425 return bond
426
427 - def get_jacos(self):
428 """ 429 Obtain the Jacobian matrices for the solid. 430 431 @param K: the K matrix. 432 @type K: numpy.ndarray 433 @return: the Jacobian matrices 434 @rtype: numpy.ndarray 435 """ 436 from numpy import zeros, dot 437 rho = self.rho 438 sf = self.stiff 439 jacos = zeros((3,9,9), dtype='float64') 440 for idm in range(3): 441 K = self.K[idm] 442 jaco = jacos[idm] 443 jaco[:3,3:] = K/(-rho) # the upper right submatrix. 444 jaco[3:,:3] = -dot(sf, K.T) # the lower left submatrix. 445 return jacos
446
447 ################################################################################ 448 # Symmetry. 449 ################################################################################ 450 451 -class Triclinic(Material):
452 """ 453 The stiffness matrix has to be symmetric. 454 """ 455 _zeropoints_ = []
456 - def __init__(self, *args, **kw):
457 for key in kw.keys(): # becaues I modify the key. 458 if len(key) == 4 and key[:2] == 'co': 459 try: 460 i = int(key[2]) 461 j = int(key[3]) 462 except: 463 continue 464 symkey = 'co%d%d' % (j, i) 465 if i != j: 466 assert symkey not in kw 467 kw[symkey] = kw[key] 468 super(Triclinic, self).__init__(*args, **kw)
469 @classmethod
470 - def _check_origstiffzero(cls, origstiff):
471 for i, j in cls._zeropoints_: 472 assert origstiff[i,j] == origstiff[j,i] == 0.0
473
474 -class Monoclinic(Triclinic):
475 _zeropoints_ = [ 476 (0,3), (0,5), 477 (1,3), (1,5), 478 (2,3), (2,5), 479 (3,4), (4,5), 480 ]
481
482 -class Orthorhombic(Triclinic):
483 _zeropoints_ = [ 484 (0,3), (0,4), (0,5), 485 (1,3), (1,4), (1,5), 486 (2,3), (2,4), (2,5), 487 (3,4), (3,5), (4,5), 488 ]
489
490 -class Tetragonal(Triclinic):
491 _zeropoints_ = [ 492 (0,3), (0,4), 493 (1,3), (1,4), 494 (2,3), (2,4), (2,5), 495 (3,4), (3,5), (4,5), 496 ]
497 - def __init__(self, *args, **kw):
498 kw['co22'] = kw['co11'] 499 kw['co23'] = kw['co13'] 500 kw['co26'] = -kw.get('co16', 0.0) 501 kw['co55'] = kw['co44'] 502 super(Tetragonal, self).__init__(*args, **kw)
503
504 -class Trigonal(Triclinic):
505 _zeropoints_ = [ 506 (0,5), (1,5), 507 (2,3), (2,4), (2,5), 508 (3,4), 509 ]
510 - def __init__(self, *args, **kw):
511 kw['co15'] = -kw.get('co25', 0.0) 512 kw['co22'] = kw['co11'] 513 kw['co23'] = kw['co13'] 514 kw['co24'] = -kw.get('co14', 0.0) 515 kw['co46'] = kw.get('co25', 0.0) 516 kw['co55'] = kw['co44'] 517 kw['co56'] = kw.get('co14', 0.0) 518 kw['co66'] = (kw['co11'] - kw['co12'])/2 519 super(Trigonal, self).__init__(*args, **kw)
520
521 -class Hexagonal(Trigonal):
522 _zeropoints_ = [ 523 (0,3), (0,4), (0,5), 524 (1,3), (1,4), (1,5), 525 (2,3), (2,4), (2,5), 526 (3,4), (3,5), (4,5), 527 ]
528
529 -class Cubic(Triclinic):
530 _zeropoints_ = [ 531 (0,3), (0,4), (0,5), 532 (1,3), (1,4), (1,5), 533 (2,3), (2,4), (2,5), 534 (3,4), (3,5), (4,5), 535 ]
536 - def __init__(self, *args, **kw):
537 kw['co13'] = kw['co12'] 538 kw['co22'] = kw['co11'] 539 kw['co23'] = kw['co12'] 540 kw['co33'] = kw['co11'] 541 kw['co55'] = kw['co44'] 542 kw['co66'] = kw['co44'] 543 super(Cubic, self).__init__(*args, **kw)
544
545 -class Isotropic(Triclinic):
546 _zeropoints_ = [ 547 (0,3), (0,4), (0,5), 548 (1,3), (1,4), (1,5), 549 (2,3), (2,4), (2,5), 550 (3,4), (3,5), (4,5), 551 ]
552 - def __init__(self, *args, **kw):
553 kw['co12'] = kw['co11']-2*kw['co44'] 554 kw['co13'] = kw['co11']-2*kw['co44'] 555 kw['co22'] = kw['co11'] 556 kw['co23'] = kw['co11']-2*kw['co44'] 557 kw['co33'] = kw['co11'] 558 kw['co55'] = kw['co44'] 559 kw['co66'] = kw['co44'] 560 super(Isotropic, self).__init__(*args, **kw)
561
562 ################################################################################ 563 # Material properties. 564 ################################################################################ 565 566 -class GaAs(Cubic):
567 - def __init__(self, *args, **kw):
568 kw.setdefault('rho', 5307.0) 569 kw.setdefault('co11', 11.88e10) 570 kw.setdefault('co12', 5.38e10) 571 kw.setdefault('co44', 5.94e10) 572 super(GaAs, self).__init__(*args, **kw)
573
574 -class ZnO(Hexagonal):
575 - def __init__(self, *args, **kw):
576 kw.setdefault('rho', 5680.0) 577 kw.setdefault('co11', 20.97e10) 578 kw.setdefault('co12', 12.11e10) 579 kw.setdefault('co13', 10.51e10) 580 kw.setdefault('co33', 21.09e10) 581 kw.setdefault('co44', 4.247e10) 582 super(ZnO, self).__init__(*args, **kw)
583
584 -class CdS(Hexagonal):
585 - def __init__(self, *args, **kw):
586 kw.setdefault('rho', 4820.0) 587 kw.setdefault('co11', 9.07e10) 588 kw.setdefault('co12', 5.81e10) 589 kw.setdefault('co13', 5.1e10) 590 kw.setdefault('co33', 9.38e10) 591 kw.setdefault('co44', 1.504e10) 592 super(CdS, self).__init__(*args, **kw)
593
594 -class Zinc(Hexagonal):
595 - def __init__(self, *args, **kw):
596 kw.setdefault('rho', 7.1*1.e-3/(1.e-2**3)) 597 kw.setdefault('co11', 14.3e11*1.e-5/(1.e-2**2)) 598 kw.setdefault('co12', 1.7e11*1.e-5/(1.e-2**2)) 599 kw.setdefault('co13', 3.3e11*1.e-5/(1.e-2**2)) 600 kw.setdefault('co33', 5.0e11*1.e-5/(1.e-2**2)) 601 kw.setdefault('co44', 4.0e11*1.e-5/(1.e-2**2)) 602 super(Zinc, self).__init__(*args, **kw)
603
604 -class Beryl(Hexagonal):
605 - def __init__(self, *args, **kw):
606 kw.setdefault('rho', 2.7*1.e-3/(1.e-2**3)) 607 kw.setdefault('co11', 26.94e11*1.e-5/(1.e-2**2)) 608 kw.setdefault('co12', 9.61e11*1.e-5/(1.e-2**2)) 609 kw.setdefault('co13', 6.61e11*1.e-5/(1.e-2**2)) 610 kw.setdefault('co33', 23.63e11*1.e-5/(1.e-2**2)) 611 kw.setdefault('co44', 6.53e11*1.e-5/(1.e-2**2)) 612 super(Beryl, self).__init__(*args, **kw)
613
614 -class Albite(Triclinic):
615 - def __init__(self, *args, **kw):
616 #kw.setdefault('rho', ) 617 kw.setdefault('co11', 69.9e9) 618 kw.setdefault('co22', 183.5e9) 619 kw.setdefault('co33', 179.5e9) 620 kw.setdefault('co44', 24.9e9) 621 kw.setdefault('co55', 26.8e9) 622 kw.setdefault('co66', 33.5e9) 623 kw.setdefault('co12', 34.0e9) 624 kw.setdefault('co13', 30.8e9) 625 kw.setdefault('co14', 5.1e9) 626 kw.setdefault('co15', -2.4e9) 627 kw.setdefault('co16', -0.9e9) 628 kw.setdefault('co23', 5.5e9) 629 kw.setdefault('co24', -3.9e9) 630 kw.setdefault('co25', -7.7e9) 631 kw.setdefault('co26', -5.8e9) 632 kw.setdefault('co34', -8.7e9) 633 kw.setdefault('co35', 7.1e9) 634 kw.setdefault('co36', -9.8e9) 635 kw.setdefault('co45', -2.4e9) 636 kw.setdefault('co46', -7.2e9) 637 kw.setdefault('co56', 0.5e9) 638 super(Albite, self).__init__(*args, **kw)
639
640 -class Acmite(Monoclinic):
641 - def __init__(self, *args, **kw):
642 kw.setdefault('rho', 3.5e3) 643 kw.setdefault('co11', 185.8e9) 644 kw.setdefault('co22', 181.3e9) 645 kw.setdefault('co33', 234.4e9) 646 kw.setdefault('co44', 62.9e9) 647 kw.setdefault('co55', 51.0e9) 648 kw.setdefault('co66', 47.4e9) 649 kw.setdefault('co12', 68.5e9) 650 kw.setdefault('co13', 70.7e9) 651 kw.setdefault('co15', 9.8e9) 652 kw.setdefault('co23', 62.9e9) 653 kw.setdefault('co25', 9.4e9) 654 kw.setdefault('co35', 21.4e9) 655 kw.setdefault('co46', 7.7e9) 656 super(Acmite, self).__init__(*args, **kw)
657
658 -class AlphaUranium(Orthorhombic):
659 - def __init__(self, *args, **kw):
660 #kw.setdefault('rho', ) 661 kw.setdefault('rho', 8.2e3) # a false value. 662 kw.setdefault('co11', 215.e9) 663 kw.setdefault('co22', 199.e9) 664 kw.setdefault('co33', 267.e9) 665 kw.setdefault('co44', 124.e9) 666 kw.setdefault('co55', 73.e9) 667 kw.setdefault('co66', 74.e9) 668 kw.setdefault('co12', 46.e9) 669 kw.setdefault('co13', 22.e9) 670 kw.setdefault('co23', 107.e9) 671 super(AlphaUranium, self).__init__(*args, **kw)
672
673 -class BariumTitanate(Tetragonal):
674 - def __init__(self, *args, **kw):
675 kw.setdefault('rho', 6.2e3) 676 kw.setdefault('co11', 275.0e9) 677 kw.setdefault('co33', 165.0e9) 678 kw.setdefault('co44', 54.3e9) 679 kw.setdefault('co66', 113.0e9) 680 kw.setdefault('co12', 179.0e9) 681 kw.setdefault('co13', 151.0e9) 682 super(BariumTitanate, self).__init__(*args, **kw)
683
684 -class AlphaQuartz(Trigonal):
685 - def __init__(self, *args, **kw):
686 kw.setdefault('rho', 2.651e3) 687 kw.setdefault('co11', 87.6e9) 688 kw.setdefault('co33', 106.8e9) 689 kw.setdefault('co44', 57.2e9) 690 kw.setdefault('co12', 6.1e9) 691 kw.setdefault('co13', 13.3e9) 692 kw.setdefault('co14', 17.3e9) 693 super(AlphaQuartz, self).__init__(*args, **kw)
694
695 -class RickerSample(Isotropic):
696 - def __init__(self, *args, **kw):
697 kw.setdefault('rho', 2200.e0) 698 kw.setdefault('co11', 3200.e0**2*2200.e0) 699 kw.setdefault('co44', 1847.5e0**2*2200.e0) 700 super(RickerSample, self).__init__(*args, **kw)
701 -class RickerSampleLight(Isotropic):
702 - def __init__(self, *args, **kw):
703 scale = 1.e-3 704 kw.setdefault('rho', 2200.e0*scale) 705 kw.setdefault('co11', 3200.e0**2*2200.e0*scale) 706 kw.setdefault('co44', 1847.5e0**2*2200.e0*scale) 707 super(RickerSampleLight, self).__init__(*args, **kw)
708