_triinterpolate.py 61 KB

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