triinterpolate.py 63 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475147614771478147914801481148214831484148514861487148814891490149114921493149414951496149714981499150015011502150315041505150615071508150915101511151215131514151515161517151815191520152115221523152415251526152715281529153015311532153315341535153615371538153915401541154215431544154515461547154815491550155115521553155415551556155715581559156015611562156315641565156615671568156915701571157215731574157515761577157815791580158115821583158415851586158715881589159015911592159315941595159615971598159916001601160216031604160516061607160816091610161116121613
  1. """
  2. Interpolation inside triangular grids.
  3. """
  4. import numpy as np
  5. from matplotlib import cbook
  6. from matplotlib.tri import Triangulation
  7. from matplotlib.tri.trifinder import TriFinder
  8. from matplotlib.tri.tritools import TriAnalyzer
  9. __all__ = ('TriInterpolator', 'LinearTriInterpolator', 'CubicTriInterpolator')
  10. class TriInterpolator:
  11. """
  12. Abstract base class for classes used to perform interpolation on
  13. triangular grids.
  14. Derived classes implement the following methods:
  15. - ``__call__(x, y)`` ,
  16. where x, y are array-like point coordinates of the same shape, and
  17. that returns a masked array of the same shape containing the
  18. interpolated z-values.
  19. - ``gradient(x, y)`` ,
  20. where x, y are array-like point coordinates of the same
  21. shape, and that returns a list of 2 masked arrays of the same shape
  22. containing the 2 derivatives of the interpolator (derivatives of
  23. interpolated z values with respect to x and y).
  24. """
  25. def __init__(self, triangulation, z, trifinder=None):
  26. cbook._check_isinstance(Triangulation, triangulation=triangulation)
  27. self._triangulation = triangulation
  28. self._z = np.asarray(z)
  29. if self._z.shape != self._triangulation.x.shape:
  30. raise ValueError("z array must have same length as triangulation x"
  31. " and y arrays")
  32. cbook._check_isinstance((TriFinder, None), trifinder=trifinder)
  33. self._trifinder = trifinder or self._triangulation.get_trifinder()
  34. # Default scaling factors : 1.0 (= no scaling)
  35. # Scaling may be used for interpolations for which the order of
  36. # magnitude of x, y has an impact on the interpolant definition.
  37. # Please refer to :meth:`_interpolate_multikeys` for details.
  38. self._unit_x = 1.0
  39. self._unit_y = 1.0
  40. # Default triangle renumbering: None (= no renumbering)
  41. # Renumbering may be used to avoid unnecessary computations
  42. # if complex calculations are done inside the Interpolator.
  43. # Please refer to :meth:`_interpolate_multikeys` for details.
  44. self._tri_renum = None
  45. # __call__ and gradient docstrings are shared by all subclasses
  46. # (except, if needed, relevant additions).
  47. # However these methods are only implemented in subclasses to avoid
  48. # confusion in the documentation.
  49. _docstring__call__ = """
  50. Returns a masked array containing interpolated values at the specified
  51. (x, y) points.
  52. Parameters
  53. ----------
  54. x, y : array-like
  55. x and y coordinates of the same shape and any number of
  56. dimensions.
  57. Returns
  58. -------
  59. z : np.ma.array
  60. Masked array of the same shape as *x* and *y*; values corresponding
  61. to (*x*, *y*) points outside of the triangulation are masked out.
  62. """
  63. _docstringgradient = r"""
  64. Returns a list of 2 masked arrays containing interpolated derivatives
  65. at the specified (x, y) points.
  66. Parameters
  67. ----------
  68. x, y : array-like
  69. x and y coordinates of the same shape and any number of
  70. dimensions.
  71. Returns
  72. -------
  73. dzdx, dzdy : np.ma.array
  74. 2 masked arrays of the same shape as *x* and *y*; values
  75. corresponding to (x, y) points outside of the triangulation
  76. are masked out.
  77. The first returned array contains the values of
  78. :math:`\frac{\partial z}{\partial x}` and the second those of
  79. :math:`\frac{\partial z}{\partial y}`.
  80. """
  81. def _interpolate_multikeys(self, x, y, tri_index=None,
  82. return_keys=('z',)):
  83. """
  84. Versatile (private) method defined for all TriInterpolators.
  85. :meth:`_interpolate_multikeys` is a wrapper around method
  86. :meth:`_interpolate_single_key` (to be defined in the child
  87. subclasses).
  88. :meth:`_interpolate_single_key actually performs the interpolation,
  89. but only for 1-dimensional inputs and at valid locations (inside
  90. unmasked triangles of the triangulation).
  91. The purpose of :meth:`_interpolate_multikeys` is to implement the
  92. following common tasks needed in all subclasses implementations:
  93. - calculation of containing triangles
  94. - dealing with more than one interpolation request at the same
  95. location (e.g., if the 2 derivatives are requested, it is
  96. unnecessary to compute the containing triangles twice)
  97. - scaling according to self._unit_x, self._unit_y
  98. - dealing with points outside of the grid (with fill value np.nan)
  99. - dealing with multi-dimensional *x*, *y* arrays: flattening for
  100. :meth:`_interpolate_params` call and final reshaping.
  101. (Note that np.vectorize could do most of those things very well for
  102. you, but it does it by function evaluations over successive tuples of
  103. the input arrays. Therefore, this tends to be more time consuming than
  104. using optimized numpy functions - e.g., np.dot - which can be used
  105. easily on the flattened inputs, in the child-subclass methods
  106. :meth:`_interpolate_single_key`.)
  107. It is guaranteed that the calls to :meth:`_interpolate_single_key`
  108. will be done with flattened (1-d) array-like input parameters *x*, *y*
  109. and with flattened, valid `tri_index` arrays (no -1 index allowed).
  110. Parameters
  111. ----------
  112. x, y : array-like
  113. x and y coordinates indicating where interpolated values are
  114. requested.
  115. tri_index : array-like of int, optional
  116. Array of the containing triangle indices, same shape as
  117. *x* and *y*. Defaults to None. If None, these indices
  118. will be computed by a TriFinder instance.
  119. (Note: For point outside the grid, tri_index[ipt] shall be -1).
  120. return_keys : tuple of keys from {'z', 'dzdx', 'dzdy'}
  121. Defines the interpolation arrays to return, and in which order.
  122. Returns
  123. -------
  124. ret : list of arrays
  125. Each array-like contains the expected interpolated values in the
  126. order defined by *return_keys* parameter.
  127. """
  128. # Flattening and rescaling inputs arrays x, y
  129. # (initial shape is stored for output)
  130. x = np.asarray(x, dtype=np.float64)
  131. y = np.asarray(y, dtype=np.float64)
  132. sh_ret = x.shape
  133. if x.shape != y.shape:
  134. raise ValueError("x and y shall have same shapes."
  135. " Given: {0} and {1}".format(x.shape, y.shape))
  136. x = np.ravel(x)
  137. y = np.ravel(y)
  138. x_scaled = x/self._unit_x
  139. y_scaled = y/self._unit_y
  140. size_ret = np.size(x_scaled)
  141. # Computes & ravels the element indexes, extract the valid ones.
  142. if tri_index is None:
  143. tri_index = self._trifinder(x, y)
  144. else:
  145. if tri_index.shape != sh_ret:
  146. raise ValueError(
  147. "tri_index array is provided and shall"
  148. " have same shape as x and y. Given: "
  149. "{0} and {1}".format(tri_index.shape, sh_ret))
  150. tri_index = np.ravel(tri_index)
  151. mask_in = (tri_index != -1)
  152. if self._tri_renum is None:
  153. valid_tri_index = tri_index[mask_in]
  154. else:
  155. valid_tri_index = self._tri_renum[tri_index[mask_in]]
  156. valid_x = x_scaled[mask_in]
  157. valid_y = y_scaled[mask_in]
  158. ret = []
  159. for return_key in return_keys:
  160. # Find the return index associated with the key.
  161. try:
  162. return_index = {'z': 0, 'dzdx': 1, 'dzdy': 2}[return_key]
  163. except KeyError:
  164. raise ValueError("return_keys items shall take values in"
  165. " {'z', 'dzdx', 'dzdy'}")
  166. # Sets the scale factor for f & df components
  167. scale = [1., 1./self._unit_x, 1./self._unit_y][return_index]
  168. # Computes the interpolation
  169. ret_loc = np.empty(size_ret, dtype=np.float64)
  170. ret_loc[~mask_in] = np.nan
  171. ret_loc[mask_in] = self._interpolate_single_key(
  172. return_key, valid_tri_index, valid_x, valid_y) * scale
  173. ret += [np.ma.masked_invalid(ret_loc.reshape(sh_ret), copy=False)]
  174. return ret
  175. def _interpolate_single_key(self, return_key, tri_index, x, y):
  176. """
  177. Performs the interpolation at points belonging to the triangulation
  178. (inside an unmasked triangles).
  179. Parameters
  180. ----------
  181. return_index : {'z', 'dzdx', 'dzdy'}
  182. Identifies the requested values (z or its derivatives)
  183. tri_index : 1d integer array
  184. Valid triangle index (-1 prohibited)
  185. x, y : 1d arrays, same shape as `tri_index`
  186. Valid locations where interpolation is requested.
  187. Returns
  188. -------
  189. ret : 1-d array
  190. Returned array of the same size as *tri_index*
  191. """
  192. raise NotImplementedError("TriInterpolator subclasses" +
  193. "should implement _interpolate_single_key!")
  194. class LinearTriInterpolator(TriInterpolator):
  195. """
  196. A LinearTriInterpolator performs linear interpolation on a triangular grid.
  197. Each triangle is represented by a plane so that an interpolated value at
  198. point (x, y) lies on the plane of the triangle containing (x, y).
  199. Interpolated values are therefore continuous across the triangulation, but
  200. their first derivatives are discontinuous at edges between triangles.
  201. Parameters
  202. ----------
  203. triangulation : :class:`~matplotlib.tri.Triangulation` object
  204. The triangulation to interpolate over.
  205. z : array-like of shape (npoints,)
  206. Array of values, defined at grid points, to interpolate between.
  207. trifinder : :class:`~matplotlib.tri.TriFinder` object, optional
  208. If this is not specified, the Triangulation's default TriFinder will
  209. be used by calling
  210. :func:`matplotlib.tri.Triangulation.get_trifinder`.
  211. Methods
  212. -------
  213. `__call__` (x, y) : Returns interpolated values at (x, y) points.
  214. `gradient` (x, y) : Returns interpolated derivatives at (x, y) points.
  215. """
  216. def __init__(self, triangulation, z, trifinder=None):
  217. TriInterpolator.__init__(self, triangulation, z, trifinder)
  218. # Store plane coefficients for fast interpolation calculations.
  219. self._plane_coefficients = \
  220. self._triangulation.calculate_plane_coefficients(self._z)
  221. def __call__(self, x, y):
  222. return self._interpolate_multikeys(x, y, tri_index=None,
  223. return_keys=('z',))[0]
  224. __call__.__doc__ = TriInterpolator._docstring__call__
  225. def gradient(self, x, y):
  226. return self._interpolate_multikeys(x, y, tri_index=None,
  227. return_keys=('dzdx', 'dzdy'))
  228. gradient.__doc__ = TriInterpolator._docstringgradient
  229. def _interpolate_single_key(self, return_key, tri_index, x, y):
  230. if return_key == 'z':
  231. return (self._plane_coefficients[tri_index, 0]*x +
  232. self._plane_coefficients[tri_index, 1]*y +
  233. self._plane_coefficients[tri_index, 2])
  234. elif return_key == 'dzdx':
  235. return self._plane_coefficients[tri_index, 0]
  236. elif return_key == 'dzdy':
  237. return self._plane_coefficients[tri_index, 1]
  238. else:
  239. raise ValueError("Invalid return_key: " + return_key)
  240. class CubicTriInterpolator(TriInterpolator):
  241. r"""
  242. A CubicTriInterpolator performs cubic interpolation on triangular grids.
  243. In one-dimension - on a segment - a cubic interpolating function is
  244. defined by the values of the function and its derivative at both ends.
  245. This is almost the same in 2-d inside a triangle, except that the values
  246. of the function and its 2 derivatives have to be defined at each triangle
  247. node.
  248. The CubicTriInterpolator takes the value of the function at each node -
  249. provided by the user - and internally computes the value of the
  250. derivatives, resulting in a smooth interpolation.
  251. (As a special feature, the user can also impose the value of the
  252. derivatives at each node, but this is not supposed to be the common
  253. usage.)
  254. Parameters
  255. ----------
  256. triangulation : :class:`~matplotlib.tri.Triangulation` object
  257. The triangulation to interpolate over.
  258. z : array-like of shape (npoints,)
  259. Array of values, defined at grid points, to interpolate between.
  260. kind : {'min_E', 'geom', 'user'}, optional
  261. Choice of the smoothing algorithm, in order to compute
  262. the interpolant derivatives (defaults to 'min_E'):
  263. - if 'min_E': (default) The derivatives at each node is computed
  264. to minimize a bending energy.
  265. - if 'geom': The derivatives at each node is computed as a
  266. weighted average of relevant triangle normals. To be used for
  267. speed optimization (large grids).
  268. - if 'user': The user provides the argument *dz*, no computation
  269. is hence needed.
  270. trifinder : :class:`~matplotlib.tri.TriFinder` object, optional
  271. If not specified, the Triangulation's default TriFinder will
  272. be used by calling
  273. :func:`matplotlib.tri.Triangulation.get_trifinder`.
  274. dz : tuple of array-likes (dzdx, dzdy), optional
  275. Used only if *kind* ='user'. In this case *dz* must be provided as
  276. (dzdx, dzdy) where dzdx, dzdy are arrays of the same shape as *z* and
  277. are the interpolant first derivatives at the *triangulation* points.
  278. Methods
  279. -------
  280. `__call__` (x, y) : Returns interpolated values at (x, y) points.
  281. `gradient` (x, y) : Returns interpolated derivatives at (x, y) points.
  282. Notes
  283. -----
  284. This note is a bit technical and details the way a
  285. :class:`~matplotlib.tri.CubicTriInterpolator` computes a cubic
  286. interpolation.
  287. The interpolation is based on a Clough-Tocher subdivision scheme of
  288. the *triangulation* mesh (to make it clearer, each triangle of the
  289. grid will be divided in 3 child-triangles, and on each child triangle
  290. the interpolated function is a cubic polynomial of the 2 coordinates).
  291. This technique originates from FEM (Finite Element Method) analysis;
  292. the element used is a reduced Hsieh-Clough-Tocher (HCT)
  293. element. Its shape functions are described in [1]_.
  294. The assembled function is guaranteed to be C1-smooth, i.e. it is
  295. continuous and its first derivatives are also continuous (this
  296. is easy to show inside the triangles but is also true when crossing the
  297. edges).
  298. In the default case (*kind* ='min_E'), the interpolant minimizes a
  299. curvature energy on the functional space generated by the HCT element
  300. shape functions - with imposed values but arbitrary derivatives at each
  301. node. The minimized functional is the integral of the so-called total
  302. curvature (implementation based on an algorithm from [2]_ - PCG sparse
  303. solver):
  304. .. math::
  305. E(z) = \frac{1}{2} \int_{\Omega} \left(
  306. \left( \frac{\partial^2{z}}{\partial{x}^2} \right)^2 +
  307. \left( \frac{\partial^2{z}}{\partial{y}^2} \right)^2 +
  308. 2\left( \frac{\partial^2{z}}{\partial{y}\partial{x}} \right)^2
  309. \right) dx\,dy
  310. If the case *kind* ='geom' is chosen by the user, a simple geometric
  311. approximation is used (weighted average of the triangle normal
  312. vectors), which could improve speed on very large grids.
  313. References
  314. ----------
  315. .. [1] Michel Bernadou, Kamal Hassan, "Basis functions for general
  316. Hsieh-Clough-Tocher triangles, complete or reduced.",
  317. International Journal for Numerical Methods in Engineering,
  318. 17(5):784 - 789. 2.01.
  319. .. [2] C.T. Kelley, "Iterative Methods for Optimization".
  320. """
  321. def __init__(self, triangulation, z, kind='min_E', trifinder=None,
  322. dz=None):
  323. TriInterpolator.__init__(self, triangulation, z, trifinder)
  324. # Loads the underlying c++ _triangulation.
  325. # (During loading, reordering of triangulation._triangles may occur so
  326. # that all final triangles are now anti-clockwise)
  327. self._triangulation.get_cpp_triangulation()
  328. # To build the stiffness matrix and avoid zero-energy spurious modes
  329. # we will only store internally the valid (unmasked) triangles and
  330. # the necessary (used) points coordinates.
  331. # 2 renumbering tables need to be computed and stored:
  332. # - a triangle renum table in order to translate the result from a
  333. # TriFinder instance into the internal stored triangle number.
  334. # - a node renum table to overwrite the self._z values into the new
  335. # (used) node numbering.
  336. tri_analyzer = TriAnalyzer(self._triangulation)
  337. (compressed_triangles, compressed_x, compressed_y, tri_renum,
  338. node_renum) = tri_analyzer._get_compressed_triangulation(True, True)
  339. self._triangles = compressed_triangles
  340. self._tri_renum = tri_renum
  341. # Taking into account the node renumbering in self._z:
  342. node_mask = (node_renum == -1)
  343. self._z[node_renum[~node_mask]] = self._z
  344. self._z = self._z[~node_mask]
  345. # Computing scale factors
  346. self._unit_x = np.ptp(compressed_x)
  347. self._unit_y = np.ptp(compressed_y)
  348. self._pts = np.column_stack([compressed_x / self._unit_x,
  349. compressed_y / self._unit_y])
  350. # Computing triangle points
  351. self._tris_pts = self._pts[self._triangles]
  352. # Computing eccentricities
  353. self._eccs = self._compute_tri_eccentricities(self._tris_pts)
  354. # Computing dof estimations for HCT triangle shape function
  355. self._dof = self._compute_dof(kind, dz=dz)
  356. # Loading HCT element
  357. self._ReferenceElement = _ReducedHCT_Element()
  358. def __call__(self, x, y):
  359. return self._interpolate_multikeys(x, y, tri_index=None,
  360. return_keys=('z',))[0]
  361. __call__.__doc__ = TriInterpolator._docstring__call__
  362. def gradient(self, x, y):
  363. return self._interpolate_multikeys(x, y, tri_index=None,
  364. return_keys=('dzdx', 'dzdy'))
  365. gradient.__doc__ = TriInterpolator._docstringgradient
  366. def _interpolate_single_key(self, return_key, tri_index, x, y):
  367. tris_pts = self._tris_pts[tri_index]
  368. alpha = self._get_alpha_vec(x, y, tris_pts)
  369. ecc = self._eccs[tri_index]
  370. dof = np.expand_dims(self._dof[tri_index], axis=1)
  371. if return_key == 'z':
  372. return self._ReferenceElement.get_function_values(
  373. alpha, ecc, dof)
  374. elif return_key in ['dzdx', 'dzdy']:
  375. J = self._get_jacobian(tris_pts)
  376. dzdx = self._ReferenceElement.get_function_derivatives(
  377. alpha, J, ecc, dof)
  378. if return_key == 'dzdx':
  379. return dzdx[:, 0, 0]
  380. else:
  381. return dzdx[:, 1, 0]
  382. else:
  383. raise ValueError("Invalid return_key: " + return_key)
  384. def _compute_dof(self, kind, dz=None):
  385. """
  386. Computes and returns nodal dofs according to kind
  387. Parameters
  388. ----------
  389. kind : {'min_E', 'geom', 'user'}
  390. Choice of the _DOF_estimator subclass to perform the gradient
  391. estimation.
  392. dz : tuple of array-likes (dzdx, dzdy), optional
  393. Used only if *kind*=user; in this case passed to the
  394. :class:`_DOF_estimator_user`.
  395. Returns
  396. -------
  397. dof : array-like, shape (npts, 2)
  398. Estimation of the gradient at triangulation nodes (stored as
  399. degree of freedoms of reduced-HCT triangle elements).
  400. """
  401. if kind == 'user':
  402. if dz is None:
  403. raise ValueError("For a CubicTriInterpolator with "
  404. "*kind*='user', a valid *dz* "
  405. "argument is expected.")
  406. TE = _DOF_estimator_user(self, dz=dz)
  407. elif kind == 'geom':
  408. TE = _DOF_estimator_geom(self)
  409. elif kind == 'min_E':
  410. TE = _DOF_estimator_min_E(self)
  411. else:
  412. raise ValueError("CubicTriInterpolator *kind* proposed: {0}; "
  413. "should be one of: "
  414. "'user', 'geom', 'min_E'".format(kind))
  415. return TE.compute_dof_from_df()
  416. @staticmethod
  417. def _get_alpha_vec(x, y, tris_pts):
  418. """
  419. Fast (vectorized) function to compute barycentric coordinates alpha.
  420. Parameters
  421. ----------
  422. x, y : array-like of dim 1 (shape (nx,))
  423. Coordinates of the points whose points barycentric
  424. coordinates are requested
  425. tris_pts : array like of dim 3 (shape: (nx, 3, 2))
  426. Coordinates of the containing triangles apexes.
  427. Returns
  428. -------
  429. alpha : array of dim 2 (shape (nx, 3))
  430. Barycentric coordinates of the points inside the containing
  431. triangles.
  432. """
  433. ndim = tris_pts.ndim-2
  434. a = tris_pts[:, 1, :] - tris_pts[:, 0, :]
  435. b = tris_pts[:, 2, :] - tris_pts[:, 0, :]
  436. abT = np.stack([a, b], axis=-1)
  437. ab = _transpose_vectorized(abT)
  438. OM = np.stack([x, y], axis=1) - tris_pts[:, 0, :]
  439. metric = _prod_vectorized(ab, abT)
  440. # Here we try to deal with the colinear cases.
  441. # metric_inv is in this case set to the Moore-Penrose pseudo-inverse
  442. # meaning that we will still return a set of valid barycentric
  443. # coordinates.
  444. metric_inv = _pseudo_inv22sym_vectorized(metric)
  445. Covar = _prod_vectorized(ab, _transpose_vectorized(
  446. np.expand_dims(OM, ndim)))
  447. ksi = _prod_vectorized(metric_inv, Covar)
  448. alpha = _to_matrix_vectorized([
  449. [1-ksi[:, 0, 0]-ksi[:, 1, 0]], [ksi[:, 0, 0]], [ksi[:, 1, 0]]])
  450. return alpha
  451. @staticmethod
  452. def _get_jacobian(tris_pts):
  453. """
  454. Fast (vectorized) function to compute triangle jacobian matrix.
  455. Parameters
  456. ----------
  457. tris_pts : array like of dim 3 (shape: (nx, 3, 2))
  458. Coordinates of the containing triangles apexes.
  459. Returns
  460. -------
  461. J : array of dim 3 (shape (nx, 2, 2))
  462. Barycentric coordinates of the points inside the containing
  463. triangles.
  464. J[itri,:,:] is the jacobian matrix at apex 0 of the triangle
  465. itri, so that the following (matrix) relationship holds:
  466. [dz/dksi] = [J] x [dz/dx]
  467. with x: global coordinates
  468. ksi: element parametric coordinates in triangle first apex
  469. local basis.
  470. """
  471. a = np.array(tris_pts[:, 1, :] - tris_pts[:, 0, :])
  472. b = np.array(tris_pts[:, 2, :] - tris_pts[:, 0, :])
  473. J = _to_matrix_vectorized([[a[:, 0], a[:, 1]],
  474. [b[:, 0], b[:, 1]]])
  475. return J
  476. @staticmethod
  477. def _compute_tri_eccentricities(tris_pts):
  478. """
  479. Computes triangle eccentricities
  480. Parameters
  481. ----------
  482. tris_pts : array like of dim 3 (shape: (nx, 3, 2))
  483. Coordinates of the triangles apexes.
  484. Returns
  485. -------
  486. ecc : array like of dim 2 (shape: (nx, 3))
  487. The so-called eccentricity parameters [1] needed for
  488. HCT triangular element.
  489. """
  490. a = np.expand_dims(tris_pts[:, 2, :] - tris_pts[:, 1, :], axis=2)
  491. b = np.expand_dims(tris_pts[:, 0, :] - tris_pts[:, 2, :], axis=2)
  492. c = np.expand_dims(tris_pts[:, 1, :] - tris_pts[:, 0, :], axis=2)
  493. # Do not use np.squeeze, this is dangerous if only one triangle
  494. # in the triangulation...
  495. dot_a = _prod_vectorized(_transpose_vectorized(a), a)[:, 0, 0]
  496. dot_b = _prod_vectorized(_transpose_vectorized(b), b)[:, 0, 0]
  497. dot_c = _prod_vectorized(_transpose_vectorized(c), c)[:, 0, 0]
  498. # Note that this line will raise a warning for dot_a, dot_b or dot_c
  499. # zeros, but we choose not to support triangles with duplicate points.
  500. return _to_matrix_vectorized([[(dot_c-dot_b) / dot_a],
  501. [(dot_a-dot_c) / dot_b],
  502. [(dot_b-dot_a) / dot_c]])
  503. # FEM element used for interpolation and for solving minimisation
  504. # problem (Reduced HCT element)
  505. class _ReducedHCT_Element:
  506. """
  507. Implementation of reduced HCT triangular element with explicit shape
  508. functions.
  509. Computes z, dz, d2z and the element stiffness matrix for bending energy:
  510. E(f) = integral( (d2z/dx2 + d2z/dy2)**2 dA)
  511. *** Reference for the shape functions: ***
  512. [1] Basis functions for general Hsieh-Clough-Tocher _triangles, complete or
  513. reduced.
  514. Michel Bernadou, Kamal Hassan
  515. International Journal for Numerical Methods in Engineering.
  516. 17(5):784 - 789. 2.01
  517. *** Element description: ***
  518. 9 dofs: z and dz given at 3 apex
  519. C1 (conform)
  520. """
  521. # 1) Loads matrices to generate shape functions as a function of
  522. # triangle eccentricities - based on [1] p.11 '''
  523. M = np.array([
  524. [ 0.00, 0.00, 0.00, 4.50, 4.50, 0.00, 0.00, 0.00, 0.00, 0.00],
  525. [-0.25, 0.00, 0.00, 0.50, 1.25, 0.00, 0.00, 0.00, 0.00, 0.00],
  526. [-0.25, 0.00, 0.00, 1.25, 0.50, 0.00, 0.00, 0.00, 0.00, 0.00],
  527. [ 0.50, 1.00, 0.00, -1.50, 0.00, 3.00, 3.00, 0.00, 0.00, 3.00],
  528. [ 0.00, 0.00, 0.00, -0.25, 0.25, 0.00, 1.00, 0.00, 0.00, 0.50],
  529. [ 0.25, 0.00, 0.00, -0.50, -0.25, 1.00, 0.00, 0.00, 0.00, 1.00],
  530. [ 0.50, 0.00, 1.00, 0.00, -1.50, 0.00, 0.00, 3.00, 3.00, 3.00],
  531. [ 0.25, 0.00, 0.00, -0.25, -0.50, 0.00, 0.00, 0.00, 1.00, 1.00],
  532. [ 0.00, 0.00, 0.00, 0.25, -0.25, 0.00, 0.00, 1.00, 0.00, 0.50]])
  533. M0 = np.array([
  534. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  535. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  536. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  537. [-1.00, 0.00, 0.00, 1.50, 1.50, 0.00, 0.00, 0.00, 0.00, -3.00],
  538. [-0.50, 0.00, 0.00, 0.75, 0.75, 0.00, 0.00, 0.00, 0.00, -1.50],
  539. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  540. [ 1.00, 0.00, 0.00, -1.50, -1.50, 0.00, 0.00, 0.00, 0.00, 3.00],
  541. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  542. [ 0.50, 0.00, 0.00, -0.75, -0.75, 0.00, 0.00, 0.00, 0.00, 1.50]])
  543. M1 = np.array([
  544. [-0.50, 0.00, 0.00, 1.50, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  545. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  546. [-0.25, 0.00, 0.00, 0.75, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  547. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  548. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  549. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  550. [ 0.50, 0.00, 0.00, -1.50, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  551. [ 0.25, 0.00, 0.00, -0.75, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  552. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00]])
  553. M2 = np.array([
  554. [ 0.50, 0.00, 0.00, 0.00, -1.50, 0.00, 0.00, 0.00, 0.00, 0.00],
  555. [ 0.25, 0.00, 0.00, 0.00, -0.75, 0.00, 0.00, 0.00, 0.00, 0.00],
  556. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  557. [-0.50, 0.00, 0.00, 0.00, 1.50, 0.00, 0.00, 0.00, 0.00, 0.00],
  558. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  559. [-0.25, 0.00, 0.00, 0.00, 0.75, 0.00, 0.00, 0.00, 0.00, 0.00],
  560. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  561. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00],
  562. [ 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00, 0.00]])
  563. # 2) Loads matrices to rotate components of gradient & Hessian
  564. # vectors in the reference basis of triangle first apex (a0)
  565. rotate_dV = np.array([[ 1., 0.], [ 0., 1.],
  566. [ 0., 1.], [-1., -1.],
  567. [-1., -1.], [ 1., 0.]])
  568. rotate_d2V = np.array([[1., 0., 0.], [0., 1., 0.], [ 0., 0., 1.],
  569. [0., 1., 0.], [1., 1., 1.], [ 0., -2., -1.],
  570. [1., 1., 1.], [1., 0., 0.], [-2., 0., -1.]])
  571. # 3) Loads Gauss points & weights on the 3 sub-_triangles for P2
  572. # exact integral - 3 points on each subtriangles.
  573. # NOTE: as the 2nd derivative is discontinuous , we really need those 9
  574. # points!
  575. n_gauss = 9
  576. gauss_pts = np.array([[13./18., 4./18., 1./18.],
  577. [ 4./18., 13./18., 1./18.],
  578. [ 7./18., 7./18., 4./18.],
  579. [ 1./18., 13./18., 4./18.],
  580. [ 1./18., 4./18., 13./18.],
  581. [ 4./18., 7./18., 7./18.],
  582. [ 4./18., 1./18., 13./18.],
  583. [13./18., 1./18., 4./18.],
  584. [ 7./18., 4./18., 7./18.]], dtype=np.float64)
  585. gauss_w = np.ones([9], dtype=np.float64) / 9.
  586. # 4) Stiffness matrix for curvature energy
  587. E = np.array([[1., 0., 0.], [0., 1., 0.], [0., 0., 2.]])
  588. # 5) Loads the matrix to compute DOF_rot from tri_J at apex 0
  589. J0_to_J1 = np.array([[-1., 1.], [-1., 0.]])
  590. J0_to_J2 = np.array([[ 0., -1.], [ 1., -1.]])
  591. def get_function_values(self, alpha, ecc, dofs):
  592. """
  593. Parameters
  594. ----------
  595. alpha : is a (N x 3 x 1) array (array of column-matrices) of
  596. barycentric coordinates,
  597. ecc : is a (N x 3 x 1) array (array of column-matrices) of triangle
  598. eccentricities,
  599. dofs : is a (N x 1 x 9) arrays (arrays of row-matrices) of computed
  600. degrees of freedom.
  601. Returns
  602. -------
  603. Returns the N-array of interpolated function values.
  604. """
  605. subtri = np.argmin(alpha, axis=1)[:, 0]
  606. ksi = _roll_vectorized(alpha, -subtri, axis=0)
  607. E = _roll_vectorized(ecc, -subtri, axis=0)
  608. x = ksi[:, 0, 0]
  609. y = ksi[:, 1, 0]
  610. z = ksi[:, 2, 0]
  611. x_sq = x*x
  612. y_sq = y*y
  613. z_sq = z*z
  614. V = _to_matrix_vectorized([
  615. [x_sq*x], [y_sq*y], [z_sq*z], [x_sq*z], [x_sq*y], [y_sq*x],
  616. [y_sq*z], [z_sq*y], [z_sq*x], [x*y*z]])
  617. prod = _prod_vectorized(self.M, V)
  618. prod += _scalar_vectorized(E[:, 0, 0],
  619. _prod_vectorized(self.M0, V))
  620. prod += _scalar_vectorized(E[:, 1, 0],
  621. _prod_vectorized(self.M1, V))
  622. prod += _scalar_vectorized(E[:, 2, 0],
  623. _prod_vectorized(self.M2, V))
  624. s = _roll_vectorized(prod, 3*subtri, axis=0)
  625. return _prod_vectorized(dofs, s)[:, 0, 0]
  626. def get_function_derivatives(self, alpha, J, ecc, dofs):
  627. """
  628. Parameters
  629. ----------
  630. *alpha* is a (N x 3 x 1) array (array of column-matrices of
  631. barycentric coordinates)
  632. *J* is a (N x 2 x 2) array of jacobian matrices (jacobian matrix at
  633. triangle first apex)
  634. *ecc* is a (N x 3 x 1) array (array of column-matrices of triangle
  635. eccentricities)
  636. *dofs* is a (N x 1 x 9) arrays (arrays of row-matrices) of computed
  637. degrees of freedom.
  638. Returns
  639. -------
  640. Returns the values of interpolated function derivatives [dz/dx, dz/dy]
  641. in global coordinates at locations alpha, as a column-matrices of
  642. shape (N x 2 x 1).
  643. """
  644. subtri = np.argmin(alpha, axis=1)[:, 0]
  645. ksi = _roll_vectorized(alpha, -subtri, axis=0)
  646. E = _roll_vectorized(ecc, -subtri, axis=0)
  647. x = ksi[:, 0, 0]
  648. y = ksi[:, 1, 0]
  649. z = ksi[:, 2, 0]
  650. x_sq = x*x
  651. y_sq = y*y
  652. z_sq = z*z
  653. dV = _to_matrix_vectorized([
  654. [ -3.*x_sq, -3.*x_sq],
  655. [ 3.*y_sq, 0.],
  656. [ 0., 3.*z_sq],
  657. [ -2.*x*z, -2.*x*z+x_sq],
  658. [-2.*x*y+x_sq, -2.*x*y],
  659. [ 2.*x*y-y_sq, -y_sq],
  660. [ 2.*y*z, y_sq],
  661. [ z_sq, 2.*y*z],
  662. [ -z_sq, 2.*x*z-z_sq],
  663. [ x*z-y*z, x*y-y*z]])
  664. # Puts back dV in first apex basis
  665. dV = _prod_vectorized(dV, _extract_submatrices(
  666. self.rotate_dV, subtri, block_size=2, axis=0))
  667. prod = _prod_vectorized(self.M, dV)
  668. prod += _scalar_vectorized(E[:, 0, 0],
  669. _prod_vectorized(self.M0, dV))
  670. prod += _scalar_vectorized(E[:, 1, 0],
  671. _prod_vectorized(self.M1, dV))
  672. prod += _scalar_vectorized(E[:, 2, 0],
  673. _prod_vectorized(self.M2, dV))
  674. dsdksi = _roll_vectorized(prod, 3*subtri, axis=0)
  675. dfdksi = _prod_vectorized(dofs, dsdksi)
  676. # In global coordinates:
  677. # Here we try to deal with the simplest colinear cases, returning a
  678. # null matrix.
  679. J_inv = _safe_inv22_vectorized(J)
  680. dfdx = _prod_vectorized(J_inv, _transpose_vectorized(dfdksi))
  681. return dfdx
  682. def get_function_hessians(self, alpha, J, ecc, dofs):
  683. """
  684. Parameters
  685. ----------
  686. *alpha* is a (N x 3 x 1) array (array of column-matrices) of
  687. barycentric coordinates
  688. *J* is a (N x 2 x 2) array of jacobian matrices (jacobian matrix at
  689. triangle first apex)
  690. *ecc* is a (N x 3 x 1) array (array of column-matrices) of triangle
  691. eccentricities
  692. *dofs* is a (N x 1 x 9) arrays (arrays of row-matrices) of computed
  693. degrees of freedom.
  694. Returns
  695. -------
  696. Returns the values of interpolated function 2nd-derivatives
  697. [d2z/dx2, d2z/dy2, d2z/dxdy] in global coordinates at locations alpha,
  698. as a column-matrices of shape (N x 3 x 1).
  699. """
  700. d2sdksi2 = self.get_d2Sidksij2(alpha, ecc)
  701. d2fdksi2 = _prod_vectorized(dofs, d2sdksi2)
  702. H_rot = self.get_Hrot_from_J(J)
  703. d2fdx2 = _prod_vectorized(d2fdksi2, H_rot)
  704. return _transpose_vectorized(d2fdx2)
  705. def get_d2Sidksij2(self, alpha, ecc):
  706. """
  707. Parameters
  708. ----------
  709. *alpha* is a (N x 3 x 1) array (array of column-matrices) of
  710. barycentric coordinates
  711. *ecc* is a (N x 3 x 1) array (array of column-matrices) of triangle
  712. eccentricities
  713. Returns
  714. -------
  715. Returns the arrays d2sdksi2 (N x 3 x 1) Hessian of shape functions
  716. expressed in covariant coordinates in first apex basis.
  717. """
  718. subtri = np.argmin(alpha, axis=1)[:, 0]
  719. ksi = _roll_vectorized(alpha, -subtri, axis=0)
  720. E = _roll_vectorized(ecc, -subtri, axis=0)
  721. x = ksi[:, 0, 0]
  722. y = ksi[:, 1, 0]
  723. z = ksi[:, 2, 0]
  724. d2V = _to_matrix_vectorized([
  725. [ 6.*x, 6.*x, 6.*x],
  726. [ 6.*y, 0., 0.],
  727. [ 0., 6.*z, 0.],
  728. [ 2.*z, 2.*z-4.*x, 2.*z-2.*x],
  729. [2.*y-4.*x, 2.*y, 2.*y-2.*x],
  730. [2.*x-4.*y, 0., -2.*y],
  731. [ 2.*z, 0., 2.*y],
  732. [ 0., 2.*y, 2.*z],
  733. [ 0., 2.*x-4.*z, -2.*z],
  734. [ -2.*z, -2.*y, x-y-z]])
  735. # Puts back d2V in first apex basis
  736. d2V = _prod_vectorized(d2V, _extract_submatrices(
  737. self.rotate_d2V, subtri, block_size=3, axis=0))
  738. prod = _prod_vectorized(self.M, d2V)
  739. prod += _scalar_vectorized(E[:, 0, 0],
  740. _prod_vectorized(self.M0, d2V))
  741. prod += _scalar_vectorized(E[:, 1, 0],
  742. _prod_vectorized(self.M1, d2V))
  743. prod += _scalar_vectorized(E[:, 2, 0],
  744. _prod_vectorized(self.M2, d2V))
  745. d2sdksi2 = _roll_vectorized(prod, 3*subtri, axis=0)
  746. return d2sdksi2
  747. def get_bending_matrices(self, J, ecc):
  748. """
  749. Parameters
  750. ----------
  751. *J* is a (N x 2 x 2) array of jacobian matrices (jacobian matrix at
  752. triangle first apex)
  753. *ecc* is a (N x 3 x 1) array (array of column-matrices) of triangle
  754. eccentricities
  755. Returns
  756. -------
  757. Returns the element K matrices for bending energy expressed in
  758. GLOBAL nodal coordinates.
  759. K_ij = integral [ (d2zi/dx2 + d2zi/dy2) * (d2zj/dx2 + d2zj/dy2) dA]
  760. tri_J is needed to rotate dofs from local basis to global basis
  761. """
  762. n = np.size(ecc, 0)
  763. # 1) matrix to rotate dofs in global coordinates
  764. J1 = _prod_vectorized(self.J0_to_J1, J)
  765. J2 = _prod_vectorized(self.J0_to_J2, J)
  766. DOF_rot = np.zeros([n, 9, 9], dtype=np.float64)
  767. DOF_rot[:, 0, 0] = 1
  768. DOF_rot[:, 3, 3] = 1
  769. DOF_rot[:, 6, 6] = 1
  770. DOF_rot[:, 1:3, 1:3] = J
  771. DOF_rot[:, 4:6, 4:6] = J1
  772. DOF_rot[:, 7:9, 7:9] = J2
  773. # 2) matrix to rotate Hessian in global coordinates.
  774. H_rot, area = self.get_Hrot_from_J(J, return_area=True)
  775. # 3) Computes stiffness matrix
  776. # Gauss quadrature.
  777. K = np.zeros([n, 9, 9], dtype=np.float64)
  778. weights = self.gauss_w
  779. pts = self.gauss_pts
  780. for igauss in range(self.n_gauss):
  781. alpha = np.tile(pts[igauss, :], n).reshape(n, 3)
  782. alpha = np.expand_dims(alpha, 2)
  783. weight = weights[igauss]
  784. d2Skdksi2 = self.get_d2Sidksij2(alpha, ecc)
  785. d2Skdx2 = _prod_vectorized(d2Skdksi2, H_rot)
  786. K += weight * _prod_vectorized(_prod_vectorized(d2Skdx2, self.E),
  787. _transpose_vectorized(d2Skdx2))
  788. # 4) With nodal (not elem) dofs
  789. K = _prod_vectorized(_prod_vectorized(_transpose_vectorized(DOF_rot),
  790. K), DOF_rot)
  791. # 5) Need the area to compute total element energy
  792. return _scalar_vectorized(area, K)
  793. def get_Hrot_from_J(self, J, return_area=False):
  794. """
  795. Parameters
  796. ----------
  797. *J* is a (N x 2 x 2) array of jacobian matrices (jacobian matrix at
  798. triangle first apex)
  799. Returns
  800. -------
  801. Returns H_rot used to rotate Hessian from local basis of first apex,
  802. to global coordinates.
  803. if *return_area* is True, returns also the triangle area (0.5*det(J))
  804. """
  805. # Here we try to deal with the simplest colinear cases; a null
  806. # energy and area is imposed.
  807. J_inv = _safe_inv22_vectorized(J)
  808. Ji00 = J_inv[:, 0, 0]
  809. Ji11 = J_inv[:, 1, 1]
  810. Ji10 = J_inv[:, 1, 0]
  811. Ji01 = J_inv[:, 0, 1]
  812. H_rot = _to_matrix_vectorized([
  813. [Ji00*Ji00, Ji10*Ji10, Ji00*Ji10],
  814. [Ji01*Ji01, Ji11*Ji11, Ji01*Ji11],
  815. [2*Ji00*Ji01, 2*Ji11*Ji10, Ji00*Ji11+Ji10*Ji01]])
  816. if not return_area:
  817. return H_rot
  818. else:
  819. area = 0.5 * (J[:, 0, 0]*J[:, 1, 1] - J[:, 0, 1]*J[:, 1, 0])
  820. return H_rot, area
  821. def get_Kff_and_Ff(self, J, ecc, triangles, Uc):
  822. """
  823. Builds K and F for the following elliptic formulation:
  824. minimization of curvature energy with value of function at node
  825. imposed and derivatives 'free'.
  826. Builds the global Kff matrix in cco format.
  827. Builds the full Ff vec Ff = - Kfc x Uc
  828. Parameters
  829. ----------
  830. *J* is a (N x 2 x 2) array of jacobian matrices (jacobian matrix at
  831. triangle first apex)
  832. *ecc* is a (N x 3 x 1) array (array of column-matrices) of triangle
  833. eccentricities
  834. *triangles* is a (N x 3) array of nodes indexes.
  835. *Uc* is (N x 3) array of imposed displacements at nodes
  836. Returns
  837. -------
  838. (Kff_rows, Kff_cols, Kff_vals) Kff matrix in coo format - Duplicate
  839. (row, col) entries must be summed.
  840. Ff: force vector - dim npts * 3
  841. """
  842. ntri = np.size(ecc, 0)
  843. vec_range = np.arange(ntri, dtype=np.int32)
  844. c_indices = np.full(ntri, -1, dtype=np.int32) # for unused dofs, -1
  845. f_dof = [1, 2, 4, 5, 7, 8]
  846. c_dof = [0, 3, 6]
  847. # vals, rows and cols indices in global dof numbering
  848. f_dof_indices = _to_matrix_vectorized([[
  849. c_indices, triangles[:, 0]*2, triangles[:, 0]*2+1,
  850. c_indices, triangles[:, 1]*2, triangles[:, 1]*2+1,
  851. c_indices, triangles[:, 2]*2, triangles[:, 2]*2+1]])
  852. expand_indices = np.ones([ntri, 9, 1], dtype=np.int32)
  853. f_row_indices = _prod_vectorized(_transpose_vectorized(f_dof_indices),
  854. _transpose_vectorized(expand_indices))
  855. f_col_indices = _prod_vectorized(expand_indices, f_dof_indices)
  856. K_elem = self.get_bending_matrices(J, ecc)
  857. # Extracting sub-matrices
  858. # Explanation & notations:
  859. # * Subscript f denotes 'free' degrees of freedom (i.e. dz/dx, dz/dx)
  860. # * Subscript c denotes 'condensated' (imposed) degrees of freedom
  861. # (i.e. z at all nodes)
  862. # * F = [Ff, Fc] is the force vector
  863. # * U = [Uf, Uc] is the imposed dof vector
  864. # [ Kff Kfc ]
  865. # * K = [ ] is the laplacian stiffness matrix
  866. # [ Kcf Kff ]
  867. # * As F = K x U one gets straightforwardly: Ff = - Kfc x Uc
  868. # Computing Kff stiffness matrix in sparse coo format
  869. Kff_vals = np.ravel(K_elem[np.ix_(vec_range, f_dof, f_dof)])
  870. Kff_rows = np.ravel(f_row_indices[np.ix_(vec_range, f_dof, f_dof)])
  871. Kff_cols = np.ravel(f_col_indices[np.ix_(vec_range, f_dof, f_dof)])
  872. # Computing Ff force vector in sparse coo format
  873. Kfc_elem = K_elem[np.ix_(vec_range, f_dof, c_dof)]
  874. Uc_elem = np.expand_dims(Uc, axis=2)
  875. Ff_elem = - _prod_vectorized(Kfc_elem, Uc_elem)[:, :, 0]
  876. Ff_indices = f_dof_indices[np.ix_(vec_range, [0], f_dof)][:, 0, :]
  877. # Extracting Ff force vector in dense format
  878. # We have to sum duplicate indices - using bincount
  879. Ff = np.bincount(np.ravel(Ff_indices), weights=np.ravel(Ff_elem))
  880. return Kff_rows, Kff_cols, Kff_vals, Ff
  881. # :class:_DOF_estimator, _DOF_estimator_user, _DOF_estimator_geom,
  882. # _DOF_estimator_min_E
  883. # Private classes used to compute the degree of freedom of each triangular
  884. # element for the TriCubicInterpolator.
  885. class _DOF_estimator:
  886. """
  887. Abstract base class for classes used to perform estimation of a function
  888. first derivatives, and deduce the dofs for a CubicTriInterpolator using a
  889. reduced HCT element formulation.
  890. Derived classes implement compute_df(self,**kwargs), returning
  891. np.vstack([dfx,dfy]).T where : dfx, dfy are the estimation of the 2
  892. gradient coordinates.
  893. """
  894. def __init__(self, interpolator, **kwargs):
  895. cbook._check_isinstance(
  896. CubicTriInterpolator, interpolator=interpolator)
  897. self._pts = interpolator._pts
  898. self._tris_pts = interpolator._tris_pts
  899. self.z = interpolator._z
  900. self._triangles = interpolator._triangles
  901. (self._unit_x, self._unit_y) = (interpolator._unit_x,
  902. interpolator._unit_y)
  903. self.dz = self.compute_dz(**kwargs)
  904. self.compute_dof_from_df()
  905. def compute_dz(self, **kwargs):
  906. raise NotImplementedError
  907. def compute_dof_from_df(self):
  908. """
  909. Computes reduced-HCT elements degrees of freedom, knowing the
  910. gradient.
  911. """
  912. J = CubicTriInterpolator._get_jacobian(self._tris_pts)
  913. tri_z = self.z[self._triangles]
  914. tri_dz = self.dz[self._triangles]
  915. tri_dof = self.get_dof_vec(tri_z, tri_dz, J)
  916. return tri_dof
  917. @staticmethod
  918. def get_dof_vec(tri_z, tri_dz, J):
  919. """
  920. Computes the dof vector of a triangle, knowing the value of f, df and
  921. of the local Jacobian at each node.
  922. *tri_z*: array of shape (3,) of f nodal values
  923. *tri_dz*: array of shape (3, 2) of df/dx, df/dy nodal values
  924. *J*: Jacobian matrix in local basis of apex 0
  925. Returns dof array of shape (9,) so that for each apex iapex:
  926. dof[iapex*3+0] = f(Ai)
  927. dof[iapex*3+1] = df(Ai).(AiAi+)
  928. dof[iapex*3+2] = df(Ai).(AiAi-)]
  929. """
  930. npt = tri_z.shape[0]
  931. dof = np.zeros([npt, 9], dtype=np.float64)
  932. J1 = _prod_vectorized(_ReducedHCT_Element.J0_to_J1, J)
  933. J2 = _prod_vectorized(_ReducedHCT_Element.J0_to_J2, J)
  934. col0 = _prod_vectorized(J, np.expand_dims(tri_dz[:, 0, :], axis=2))
  935. col1 = _prod_vectorized(J1, np.expand_dims(tri_dz[:, 1, :], axis=2))
  936. col2 = _prod_vectorized(J2, np.expand_dims(tri_dz[:, 2, :], axis=2))
  937. dfdksi = _to_matrix_vectorized([
  938. [col0[:, 0, 0], col1[:, 0, 0], col2[:, 0, 0]],
  939. [col0[:, 1, 0], col1[:, 1, 0], col2[:, 1, 0]]])
  940. dof[:, 0:7:3] = tri_z
  941. dof[:, 1:8:3] = dfdksi[:, 0]
  942. dof[:, 2:9:3] = dfdksi[:, 1]
  943. return dof
  944. class _DOF_estimator_user(_DOF_estimator):
  945. """dz is imposed by user; accounts for scaling if any."""
  946. def compute_dz(self, dz):
  947. (dzdx, dzdy) = dz
  948. dzdx = dzdx * self._unit_x
  949. dzdy = dzdy * self._unit_y
  950. return np.vstack([dzdx, dzdy]).T
  951. class _DOF_estimator_geom(_DOF_estimator):
  952. """Fast 'geometric' approximation, recommended for large arrays."""
  953. def compute_dz(self):
  954. """
  955. self.df is computed as weighted average of _triangles sharing a common
  956. node. On each triangle itri f is first assumed linear (= ~f), which
  957. allows to compute d~f[itri]
  958. Then the following approximation of df nodal values is then proposed:
  959. f[ipt] = SUM ( w[itri] x d~f[itri] , for itri sharing apex ipt)
  960. The weighted coeff. w[itri] are proportional to the angle of the
  961. triangle itri at apex ipt
  962. """
  963. el_geom_w = self.compute_geom_weights()
  964. el_geom_grad = self.compute_geom_grads()
  965. # Sum of weights coeffs
  966. w_node_sum = np.bincount(np.ravel(self._triangles),
  967. weights=np.ravel(el_geom_w))
  968. # Sum of weighted df = (dfx, dfy)
  969. dfx_el_w = np.empty_like(el_geom_w)
  970. dfy_el_w = np.empty_like(el_geom_w)
  971. for iapex in range(3):
  972. dfx_el_w[:, iapex] = el_geom_w[:, iapex]*el_geom_grad[:, 0]
  973. dfy_el_w[:, iapex] = el_geom_w[:, iapex]*el_geom_grad[:, 1]
  974. dfx_node_sum = np.bincount(np.ravel(self._triangles),
  975. weights=np.ravel(dfx_el_w))
  976. dfy_node_sum = np.bincount(np.ravel(self._triangles),
  977. weights=np.ravel(dfy_el_w))
  978. # Estimation of df
  979. dfx_estim = dfx_node_sum/w_node_sum
  980. dfy_estim = dfy_node_sum/w_node_sum
  981. return np.vstack([dfx_estim, dfy_estim]).T
  982. def compute_geom_weights(self):
  983. """
  984. Builds the (nelems x 3) weights coeffs of _triangles angles,
  985. renormalized so that np.sum(weights, axis=1) == np.ones(nelems)
  986. """
  987. weights = np.zeros([np.size(self._triangles, 0), 3])
  988. tris_pts = self._tris_pts
  989. for ipt in range(3):
  990. p0 = tris_pts[:, (ipt) % 3, :]
  991. p1 = tris_pts[:, (ipt+1) % 3, :]
  992. p2 = tris_pts[:, (ipt-1) % 3, :]
  993. alpha1 = np.arctan2(p1[:, 1]-p0[:, 1], p1[:, 0]-p0[:, 0])
  994. alpha2 = np.arctan2(p2[:, 1]-p0[:, 1], p2[:, 0]-p0[:, 0])
  995. # In the below formula we could take modulo 2. but
  996. # modulo 1. is safer regarding round-off errors (flat triangles).
  997. angle = np.abs(((alpha2-alpha1) / np.pi) % 1)
  998. # Weight proportional to angle up np.pi/2; null weight for
  999. # degenerated cases 0 and np.pi (note that *angle* is normalized
  1000. # by np.pi).
  1001. weights[:, ipt] = 0.5 - np.abs(angle-0.5)
  1002. return weights
  1003. def compute_geom_grads(self):
  1004. """
  1005. Compute the (global) gradient component of f assumed linear (~f).
  1006. returns array df of shape (nelems, 2)
  1007. df[ielem].dM[ielem] = dz[ielem] i.e. df = dz x dM = dM.T^-1 x dz
  1008. """
  1009. tris_pts = self._tris_pts
  1010. tris_f = self.z[self._triangles]
  1011. dM1 = tris_pts[:, 1, :] - tris_pts[:, 0, :]
  1012. dM2 = tris_pts[:, 2, :] - tris_pts[:, 0, :]
  1013. dM = np.dstack([dM1, dM2])
  1014. # Here we try to deal with the simplest colinear cases: a null
  1015. # gradient is assumed in this case.
  1016. dM_inv = _safe_inv22_vectorized(dM)
  1017. dZ1 = tris_f[:, 1] - tris_f[:, 0]
  1018. dZ2 = tris_f[:, 2] - tris_f[:, 0]
  1019. dZ = np.vstack([dZ1, dZ2]).T
  1020. df = np.empty_like(dZ)
  1021. # With np.einsum : could be ej,eji -> ej
  1022. df[:, 0] = dZ[:, 0]*dM_inv[:, 0, 0] + dZ[:, 1]*dM_inv[:, 1, 0]
  1023. df[:, 1] = dZ[:, 0]*dM_inv[:, 0, 1] + dZ[:, 1]*dM_inv[:, 1, 1]
  1024. return df
  1025. class _DOF_estimator_min_E(_DOF_estimator_geom):
  1026. """
  1027. The 'smoothest' approximation, df is computed through global minimization
  1028. of the bending energy:
  1029. E(f) = integral[(d2z/dx2 + d2z/dy2 + 2 d2z/dxdy)**2 dA]
  1030. """
  1031. def __init__(self, Interpolator):
  1032. self._eccs = Interpolator._eccs
  1033. _DOF_estimator_geom.__init__(self, Interpolator)
  1034. def compute_dz(self):
  1035. """
  1036. Elliptic solver for bending energy minimization.
  1037. Uses a dedicated 'toy' sparse Jacobi PCG solver.
  1038. """
  1039. # Initial guess for iterative PCG solver.
  1040. dz_init = _DOF_estimator_geom.compute_dz(self)
  1041. Uf0 = np.ravel(dz_init)
  1042. reference_element = _ReducedHCT_Element()
  1043. J = CubicTriInterpolator._get_jacobian(self._tris_pts)
  1044. eccs = self._eccs
  1045. triangles = self._triangles
  1046. Uc = self.z[self._triangles]
  1047. # Building stiffness matrix and force vector in coo format
  1048. Kff_rows, Kff_cols, Kff_vals, Ff = reference_element.get_Kff_and_Ff(
  1049. J, eccs, triangles, Uc)
  1050. # Building sparse matrix and solving minimization problem
  1051. # We could use scipy.sparse direct solver; however to avoid this
  1052. # external dependency an implementation of a simple PCG solver with
  1053. # a simple diagonal Jacobi preconditioner is implemented.
  1054. tol = 1.e-10
  1055. n_dof = Ff.shape[0]
  1056. Kff_coo = _Sparse_Matrix_coo(Kff_vals, Kff_rows, Kff_cols,
  1057. shape=(n_dof, n_dof))
  1058. Kff_coo.compress_csc()
  1059. Uf, err = _cg(A=Kff_coo, b=Ff, x0=Uf0, tol=tol)
  1060. # If the PCG did not converge, we return the best guess between Uf0
  1061. # and Uf.
  1062. err0 = np.linalg.norm(Kff_coo.dot(Uf0) - Ff)
  1063. if err0 < err:
  1064. # Maybe a good occasion to raise a warning here ?
  1065. cbook._warn_external("In TriCubicInterpolator initialization, "
  1066. "PCG sparse solver did not converge after "
  1067. "1000 iterations. `geom` approximation is "
  1068. "used instead of `min_E`")
  1069. Uf = Uf0
  1070. # Building dz from Uf
  1071. dz = np.empty([self._pts.shape[0], 2], dtype=np.float64)
  1072. dz[:, 0] = Uf[::2]
  1073. dz[:, 1] = Uf[1::2]
  1074. return dz
  1075. # The following private :class:_Sparse_Matrix_coo and :func:_cg provide
  1076. # a PCG sparse solver for (symmetric) elliptic problems.
  1077. class _Sparse_Matrix_coo:
  1078. def __init__(self, vals, rows, cols, shape):
  1079. """
  1080. Creates a sparse matrix in coo format
  1081. *vals*: arrays of values of non-null entries of the matrix
  1082. *rows*: int arrays of rows of non-null entries of the matrix
  1083. *cols*: int arrays of cols of non-null entries of the matrix
  1084. *shape*: 2-tuple (n, m) of matrix shape
  1085. """
  1086. self.n, self.m = shape
  1087. self.vals = np.asarray(vals, dtype=np.float64)
  1088. self.rows = np.asarray(rows, dtype=np.int32)
  1089. self.cols = np.asarray(cols, dtype=np.int32)
  1090. def dot(self, V):
  1091. """
  1092. Dot product of self by a vector *V* in sparse-dense to dense format
  1093. *V* dense vector of shape (self.m,)
  1094. """
  1095. assert V.shape == (self.m,)
  1096. return np.bincount(self.rows,
  1097. weights=self.vals*V[self.cols],
  1098. minlength=self.m)
  1099. def compress_csc(self):
  1100. """
  1101. Compress rows, cols, vals / summing duplicates. Sort for csc format.
  1102. """
  1103. _, unique, indices = np.unique(
  1104. self.rows + self.n*self.cols,
  1105. return_index=True, return_inverse=True)
  1106. self.rows = self.rows[unique]
  1107. self.cols = self.cols[unique]
  1108. self.vals = np.bincount(indices, weights=self.vals)
  1109. def compress_csr(self):
  1110. """
  1111. Compress rows, cols, vals / summing duplicates. Sort for csr format.
  1112. """
  1113. _, unique, indices = np.unique(
  1114. self.m*self.rows + self.cols,
  1115. return_index=True, return_inverse=True)
  1116. self.rows = self.rows[unique]
  1117. self.cols = self.cols[unique]
  1118. self.vals = np.bincount(indices, weights=self.vals)
  1119. def to_dense(self):
  1120. """
  1121. Returns a dense matrix representing self.
  1122. Mainly for debugging purposes.
  1123. """
  1124. ret = np.zeros([self.n, self.m], dtype=np.float64)
  1125. nvals = self.vals.size
  1126. for i in range(nvals):
  1127. ret[self.rows[i], self.cols[i]] += self.vals[i]
  1128. return ret
  1129. def __str__(self):
  1130. return self.to_dense().__str__()
  1131. @property
  1132. def diag(self):
  1133. """
  1134. Returns the (dense) vector of the diagonal elements.
  1135. """
  1136. in_diag = (self.rows == self.cols)
  1137. diag = np.zeros(min(self.n, self.n), dtype=np.float64) # default 0.
  1138. diag[self.rows[in_diag]] = self.vals[in_diag]
  1139. return diag
  1140. def _cg(A, b, x0=None, tol=1.e-10, maxiter=1000):
  1141. """
  1142. Use Preconditioned Conjugate Gradient iteration to solve A x = b
  1143. A simple Jacobi (diagonal) preconditionner is used.
  1144. Parameters
  1145. ----------
  1146. A : _Sparse_Matrix_coo
  1147. *A* must have been compressed before by compress_csc or
  1148. compress_csr method.
  1149. b : array
  1150. Right hand side of the linear system.
  1151. Returns
  1152. -------
  1153. x : array
  1154. The converged solution.
  1155. err : float
  1156. The absolute error np.linalg.norm(A.dot(x) - b)
  1157. Other parameters
  1158. ----------------
  1159. x0 : array
  1160. Starting guess for the solution.
  1161. tol : float
  1162. Tolerance to achieve. The algorithm terminates when the relative
  1163. residual is below tol.
  1164. maxiter : integer
  1165. Maximum number of iterations. Iteration will stop
  1166. after maxiter steps even if the specified tolerance has not
  1167. been achieved.
  1168. """
  1169. n = b.size
  1170. assert A.n == n
  1171. assert A.m == n
  1172. b_norm = np.linalg.norm(b)
  1173. # Jacobi pre-conditioner
  1174. kvec = A.diag
  1175. # For diag elem < 1e-6 we keep 1e-6.
  1176. kvec = np.maximum(kvec, 1e-6)
  1177. # Initial guess
  1178. if x0 is None:
  1179. x = np.zeros(n)
  1180. else:
  1181. x = x0
  1182. r = b - A.dot(x)
  1183. w = r/kvec
  1184. p = np.zeros(n)
  1185. beta = 0.0
  1186. rho = np.dot(r, w)
  1187. k = 0
  1188. # Following C. T. Kelley
  1189. while (np.sqrt(abs(rho)) > tol*b_norm) and (k < maxiter):
  1190. p = w + beta*p
  1191. z = A.dot(p)
  1192. alpha = rho/np.dot(p, z)
  1193. r = r - alpha*z
  1194. w = r/kvec
  1195. rhoold = rho
  1196. rho = np.dot(r, w)
  1197. x = x + alpha*p
  1198. beta = rho/rhoold
  1199. #err = np.linalg.norm(A.dot(x) - b) # absolute accuracy - not used
  1200. k += 1
  1201. err = np.linalg.norm(A.dot(x) - b)
  1202. return x, err
  1203. # The following private functions:
  1204. # :func:`_safe_inv22_vectorized`
  1205. # :func:`_pseudo_inv22sym_vectorized`
  1206. # :func:`_prod_vectorized`
  1207. # :func:`_scalar_vectorized`
  1208. # :func:`_transpose_vectorized`
  1209. # :func:`_roll_vectorized`
  1210. # :func:`_to_matrix_vectorized`
  1211. # :func:`_extract_submatrices`
  1212. # provide fast numpy implementation of some standard operations on arrays of
  1213. # matrices - stored as (:, n_rows, n_cols)-shaped np.arrays.
  1214. # Development note: Dealing with pathologic 'flat' triangles in the
  1215. # CubicTriInterpolator code and impact on (2, 2)-matrix inversion functions
  1216. # :func:`_safe_inv22_vectorized` and :func:`_pseudo_inv22sym_vectorized`.
  1217. #
  1218. # Goals:
  1219. # 1) The CubicTriInterpolator should be able to handle flat or almost flat
  1220. # triangles without raising an error,
  1221. # 2) These degenerated triangles should have no impact on the automatic dof
  1222. # calculation (associated with null weight for the _DOF_estimator_geom and
  1223. # with null energy for the _DOF_estimator_min_E),
  1224. # 3) Linear patch test should be passed exactly on degenerated meshes,
  1225. # 4) Interpolation (with :meth:`_interpolate_single_key` or
  1226. # :meth:`_interpolate_multi_key`) shall be correctly handled even *inside*
  1227. # the pathologic triangles, to interact correctly with a TriRefiner class.
  1228. #
  1229. # Difficulties:
  1230. # Flat triangles have rank-deficient *J* (so-called jacobian matrix) and
  1231. # *metric* (the metric tensor = J x J.T). Computation of the local
  1232. # tangent plane is also problematic.
  1233. #
  1234. # Implementation:
  1235. # Most of the time, when computing the inverse of a rank-deficient matrix it
  1236. # is safe to simply return the null matrix (which is the implementation in
  1237. # :func:`_safe_inv22_vectorized`). This is because of point 2), itself
  1238. # enforced by:
  1239. # - null area hence null energy in :class:`_DOF_estimator_min_E`
  1240. # - angles close or equal to 0 or np.pi hence null weight in
  1241. # :class:`_DOF_estimator_geom`.
  1242. # Note that the function angle -> weight is continuous and maximum for an
  1243. # angle np.pi/2 (refer to :meth:`compute_geom_weights`)
  1244. # The exception is the computation of barycentric coordinates, which is done
  1245. # by inversion of the *metric* matrix. In this case, we need to compute a set
  1246. # of valid coordinates (1 among numerous possibilities), to ensure point 4).
  1247. # We benefit here from the symmetry of metric = J x J.T, which makes it easier
  1248. # to compute a pseudo-inverse in :func:`_pseudo_inv22sym_vectorized`
  1249. def _safe_inv22_vectorized(M):
  1250. """
  1251. Inversion of arrays of (2, 2) matrices, returns 0 for rank-deficient
  1252. matrices.
  1253. *M* : array of (2, 2) matrices to inverse, shape (n, 2, 2)
  1254. """
  1255. assert M.ndim == 3
  1256. assert M.shape[-2:] == (2, 2)
  1257. M_inv = np.empty_like(M)
  1258. prod1 = M[:, 0, 0]*M[:, 1, 1]
  1259. delta = prod1 - M[:, 0, 1]*M[:, 1, 0]
  1260. # We set delta_inv to 0. in case of a rank deficient matrix; a
  1261. # rank-deficient input matrix *M* will lead to a null matrix in output
  1262. rank2 = (np.abs(delta) > 1e-8*np.abs(prod1))
  1263. if np.all(rank2):
  1264. # Normal 'optimized' flow.
  1265. delta_inv = 1./delta
  1266. else:
  1267. # 'Pathologic' flow.
  1268. delta_inv = np.zeros(M.shape[0])
  1269. delta_inv[rank2] = 1./delta[rank2]
  1270. M_inv[:, 0, 0] = M[:, 1, 1]*delta_inv
  1271. M_inv[:, 0, 1] = -M[:, 0, 1]*delta_inv
  1272. M_inv[:, 1, 0] = -M[:, 1, 0]*delta_inv
  1273. M_inv[:, 1, 1] = M[:, 0, 0]*delta_inv
  1274. return M_inv
  1275. def _pseudo_inv22sym_vectorized(M):
  1276. """
  1277. Inversion of arrays of (2, 2) SYMMETRIC matrices; returns the
  1278. (Moore-Penrose) pseudo-inverse for rank-deficient matrices.
  1279. In case M is of rank 1, we have M = trace(M) x P where P is the orthogonal
  1280. projection on Im(M), and we return trace(M)^-1 x P == M / trace(M)**2
  1281. In case M is of rank 0, we return the null matrix.
  1282. *M* : array of (2, 2) matrices to inverse, shape (n, 2, 2)
  1283. """
  1284. assert M.ndim == 3
  1285. assert M.shape[-2:] == (2, 2)
  1286. M_inv = np.empty_like(M)
  1287. prod1 = M[:, 0, 0]*M[:, 1, 1]
  1288. delta = prod1 - M[:, 0, 1]*M[:, 1, 0]
  1289. rank2 = (np.abs(delta) > 1e-8*np.abs(prod1))
  1290. if np.all(rank2):
  1291. # Normal 'optimized' flow.
  1292. M_inv[:, 0, 0] = M[:, 1, 1] / delta
  1293. M_inv[:, 0, 1] = -M[:, 0, 1] / delta
  1294. M_inv[:, 1, 0] = -M[:, 1, 0] / delta
  1295. M_inv[:, 1, 1] = M[:, 0, 0] / delta
  1296. else:
  1297. # 'Pathologic' flow.
  1298. # Here we have to deal with 2 sub-cases
  1299. # 1) First sub-case: matrices of rank 2:
  1300. delta = delta[rank2]
  1301. M_inv[rank2, 0, 0] = M[rank2, 1, 1] / delta
  1302. M_inv[rank2, 0, 1] = -M[rank2, 0, 1] / delta
  1303. M_inv[rank2, 1, 0] = -M[rank2, 1, 0] / delta
  1304. M_inv[rank2, 1, 1] = M[rank2, 0, 0] / delta
  1305. # 2) Second sub-case: rank-deficient matrices of rank 0 and 1:
  1306. rank01 = ~rank2
  1307. tr = M[rank01, 0, 0] + M[rank01, 1, 1]
  1308. tr_zeros = (np.abs(tr) < 1.e-8)
  1309. sq_tr_inv = (1.-tr_zeros) / (tr**2+tr_zeros)
  1310. #sq_tr_inv = 1. / tr**2
  1311. M_inv[rank01, 0, 0] = M[rank01, 0, 0] * sq_tr_inv
  1312. M_inv[rank01, 0, 1] = M[rank01, 0, 1] * sq_tr_inv
  1313. M_inv[rank01, 1, 0] = M[rank01, 1, 0] * sq_tr_inv
  1314. M_inv[rank01, 1, 1] = M[rank01, 1, 1] * sq_tr_inv
  1315. return M_inv
  1316. def _prod_vectorized(M1, M2):
  1317. """
  1318. Matrix product between arrays of matrices, or a matrix and an array of
  1319. matrices (*M1* and *M2*)
  1320. """
  1321. sh1 = M1.shape
  1322. sh2 = M2.shape
  1323. assert len(sh1) >= 2
  1324. assert len(sh2) >= 2
  1325. assert sh1[-1] == sh2[-2]
  1326. ndim1 = len(sh1)
  1327. t1_index = [*range(ndim1-2), ndim1-1, ndim1-2]
  1328. return np.sum(np.transpose(M1, t1_index)[..., np.newaxis] *
  1329. M2[..., np.newaxis, :], -3)
  1330. def _scalar_vectorized(scalar, M):
  1331. """
  1332. Scalar product between scalars and matrices.
  1333. """
  1334. return scalar[:, np.newaxis, np.newaxis]*M
  1335. def _transpose_vectorized(M):
  1336. """
  1337. Transposition of an array of matrices *M*.
  1338. """
  1339. return np.transpose(M, [0, 2, 1])
  1340. def _roll_vectorized(M, roll_indices, axis):
  1341. """
  1342. Rolls an array of matrices along an axis according to an array of indices
  1343. *roll_indices*
  1344. *axis* can be either 0 (rolls rows) or 1 (rolls columns).
  1345. """
  1346. assert axis in [0, 1]
  1347. ndim = M.ndim
  1348. assert ndim == 3
  1349. ndim_roll = roll_indices.ndim
  1350. assert ndim_roll == 1
  1351. sh = M.shape
  1352. r, c = sh[-2:]
  1353. assert sh[0] == roll_indices.shape[0]
  1354. vec_indices = np.arange(sh[0], dtype=np.int32)
  1355. # Builds the rolled matrix
  1356. M_roll = np.empty_like(M)
  1357. if axis == 0:
  1358. for ir in range(r):
  1359. for ic in range(c):
  1360. M_roll[:, ir, ic] = M[vec_indices, (-roll_indices+ir) % r, ic]
  1361. elif axis == 1:
  1362. for ir in range(r):
  1363. for ic in range(c):
  1364. M_roll[:, ir, ic] = M[vec_indices, ir, (-roll_indices+ic) % c]
  1365. return M_roll
  1366. def _to_matrix_vectorized(M):
  1367. """
  1368. Builds an array of matrices from individuals np.arrays of identical
  1369. shapes.
  1370. *M*: ncols-list of nrows-lists of shape sh.
  1371. Returns M_res np.array of shape (sh, nrow, ncols) so that:
  1372. M_res[..., i, j] = M[i][j]
  1373. """
  1374. assert isinstance(M, (tuple, list))
  1375. assert all(isinstance(item, (tuple, list)) for item in M)
  1376. c_vec = np.asarray([len(item) for item in M])
  1377. assert np.all(c_vec-c_vec[0] == 0)
  1378. r = len(M)
  1379. c = c_vec[0]
  1380. M00 = np.asarray(M[0][0])
  1381. dt = M00.dtype
  1382. sh = [M00.shape[0], r, c]
  1383. M_ret = np.empty(sh, dtype=dt)
  1384. for irow in range(r):
  1385. for icol in range(c):
  1386. M_ret[:, irow, icol] = np.asarray(M[irow][icol])
  1387. return M_ret
  1388. def _extract_submatrices(M, block_indices, block_size, axis):
  1389. """
  1390. Extracts selected blocks of a matrices *M* depending on parameters
  1391. *block_indices* and *block_size*.
  1392. Returns the array of extracted matrices *Mres* so that:
  1393. M_res[...,ir,:] = M[(block_indices*block_size+ir), :]
  1394. """
  1395. assert block_indices.ndim == 1
  1396. assert axis in [0, 1]
  1397. r, c = M.shape
  1398. if axis == 0:
  1399. sh = [block_indices.shape[0], block_size, c]
  1400. elif axis == 1:
  1401. sh = [block_indices.shape[0], r, block_size]
  1402. dt = M.dtype
  1403. M_res = np.empty(sh, dtype=dt)
  1404. if axis == 0:
  1405. for ir in range(block_size):
  1406. M_res[:, ir, :] = M[(block_indices*block_size+ir), :]
  1407. elif axis == 1:
  1408. for ic in range(block_size):
  1409. M_res[:, :, ic] = M[:, (block_indices*block_size+ic)]
  1410. return M_res