romtools.vector_space.utils.orthogonalizer
The OrthogonalizerClass is used to orthogonalize a basis at the end of the construction of a vector space. Specifically, given a basis $$\boldsymbol \Phi \in \mathbb{R}^{N \times K},$$ the orthogonalizer will compute a new, orthogonalized basis $\boldsymbol \Phi_{*}$ where $$\boldsymbol \Phi_{*}^T \mathbf{W} \boldsymbol \Phi_{*} = \mathbf{I}.$$ In the above, $\mathbf{W}$ is a weighting matrix (typically the cell volumes).
1# 2# ************************************************************************ 3# 4# ROM Tools and Workflows 5# Copyright 2019 National Technology & Engineering Solutions of Sandia,LLC 6# (NTESS) 7# 8# Under the terms of Contract DE-NA0003525 with NTESS, the 9# U.S. Government retains certain rights in this software. 10# 11# ROM Tools and Workflows is licensed under BSD-3-Clause terms of use: 12# 13# Redistribution and use in source and binary forms, with or without 14# modification, are permitted provided that the following conditions 15# are met: 16# 17# 1. Redistributions of source code must retain the above copyright 18# notice, this list of conditions and the following disclaimer. 19# 20# 2. Redistributions in binary form must reproduce the above copyright 21# notice, this list of conditions and the following disclaimer in the 22# documentation and/or other materials provided with the distribution. 23# 24# 3. Neither the name of the copyright holder nor the names of its 25# contributors may be used to endorse or promote products derived 26# from this software without specific prior written permission. 27# 28# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS 29# "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT 30# LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS 31# FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE 32# COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, 33# INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES 34# (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR 35# SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) 36# HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, 37# STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING 38# IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE 39# POSSIBILITY OF SUCH DAMAGE. 40# 41# Questions? Contact Eric Parish (ejparis@sandia.gov) 42# 43# ************************************************************************ 44# 45 46''' 47The OrthogonalizerClass is used to orthogonalize a basis at the end of the 48construction of a vector space. Specifically, given a basis 49$$\\boldsymbol \\Phi \\in \\mathbb{R}^{N \\times K},$$ 50the orthogonalizer will compute a new, orthogonalized basis $\\boldsymbol \\Phi_{\\*}$ 51where 52$$\\boldsymbol \\Phi_{\\*}^T \\mathbf{W} \\boldsymbol \\Phi_{\\*} = \\mathbf{I}.$$ 53In the above, $\\mathbf{W}$ is a weighting matrix (typically the cell volumes). 54''' 55 56from typing import Protocol 57import numpy as np 58import scipy.sparse 59 60 61class Orthogonalizer(Protocol): 62 ''' 63 Interface for the Orthogonalizer class. 64 ''' 65 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 66 ... 67 68 69class NoOpOrthogonalizer: 70 ''' 71 No op class (doesn't do anything) 72 73 This class conforms to `Orthogonalizer` protocol. 74 ''' 75 def __init__(self) -> None: 76 pass 77 78 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 79 '''Returns the input array.''' 80 return my_array 81 82 83class EuclideanL2Orthogonalizer: 84 ''' 85 Orthogonalizes the basis in the standard Euclidean L2 inner product, i.e., 86 the output basis will satisfy 87 $$\\boldsymbol \\Phi_{\\*}^T \\boldsymbol \\Phi_{\\*} = \\mathbf{I}.$$ 88 89 This class conforms to `Orthogonalizer` protocol. 90 ''' 91 def __init__(self, qrFnc=None) -> None: 92 ''' 93 Constructor 94 Args: 95 96 qrFnc: a callable to use for computing the QR decomposition. 97 IMPORTANT: must conform to the API of 98 [np.linalg.qr](https://numpy.org/doc/stable/reference/generated/numpy.linalg.qr.html). 99 If `None`, internally we use `np.linalg.qr`. 100 Note: this is useful when you want to use a custom qr, for example when your snapshots are 101 distributed with MPI, or maybe you have a fancy qr function that you can use. 102 103 ''' 104 self.__qr_picked = np.linalg.qr if qrFnc is None else qrFnc 105 106 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 107 ''' 108 Orthogonalizes the input matrix and returns the orthogonalized matrix. 109 110 Args: 111 my_array (np.ndarray): The input matrix to be orthogonalized. 112 113 Returns: 114 np.ndarray: The orthogonalized matrix. 115 ''' 116 my_array, _ = self.__qr_picked(my_array, mode='reduced') 117 return my_array 118 119 120class EuclideanVectorWeightedL2Orthogonalizer: 121 ''' 122 Orthogonalizes the basis in vector-weighted Euclidean L2 inner product, 123 i.e., the output basis will satisfy 124 $$\\boldsymbol \\Phi_{\\*}^T \\mathrm{diag}(\\mathbf{w})\\boldsymbol \\Phi_{\\*} = \\mathbf{I},$$ 125 where $\\mathbf{w}$ is the weighting vector. Typically, this inner product 126 is used for orthogonalizing with respect to cell volumes 127 128 This class conforms to `Orthogonalizer` protocol. 129 ''' 130 def __init__(self, weighting_vector: np.ndarray, qrFnc=None) -> None: 131 ''' 132 Constructor for the EuclideanVectorWeightedL2Orthogonalizer that 133 initializes the orthogonalizer with the provided weighting vector and 134 an optional custom QR decomposition function. 135 Args: 136 weighting_vector (np.ndarray): a 1-D NumPy array that the matrix will be orthogonalized against. The 137 length of the array must match the number of rows in the matrix that will be orthogonalized. 138 qrFnc: a callable to use for computing the QR decomposition. 139 IMPORTANT: must conform to the API of 140 [np.linalg.qr](https://numpy.org/doc/stable/reference/generated/numpy.linalg.qr.html). 141 If `None`, internally we use `np.linalg.qr`. 142 Note: this is useful when you want to use a custom qr, for example when your snapshots are 143 distributed with MPI, or maybe you have a fancy qr function that you can use. 144 ''' 145 146 self.__weighting_vector = weighting_vector 147 self.__qr_picked = np.linalg.qr if qrFnc is None else qrFnc 148 149 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 150 ''' 151 Orthogonalizes the input matrix using the specified weighting vector and returns the orthogonalized matrix. 152 153 Args: 154 my_array (np.ndarray): The input matrix to be orthogonalized. 155 156 Returns: 157 np.ndarray: The orthogonalized matrix. 158 ''' 159 assert my_array.shape[0] == self.__weighting_vector.size, "Weighting vector does not match basis size" 160 tmp = scipy.sparse.diags(np.sqrt(self.__weighting_vector)) @ my_array 161 my_array, _ = self.__qr_picked(tmp, mode='reduced') 162 my_array = scipy.sparse.diags(np.sqrt(1./self.__weighting_vector)) @ my_array 163 return my_array
62class Orthogonalizer(Protocol): 63 ''' 64 Interface for the Orthogonalizer class. 65 ''' 66 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 67 ...
Interface for the Orthogonalizer class.
1771def _no_init_or_replace_init(self, *args, **kwargs): 1772 cls = type(self) 1773 1774 if cls._is_protocol: 1775 raise TypeError('Protocols cannot be instantiated') 1776 1777 # Already using a custom `__init__`. No need to calculate correct 1778 # `__init__` to call. This can lead to RecursionError. See bpo-45121. 1779 if cls.__init__ is not _no_init_or_replace_init: 1780 return 1781 1782 # Initially, `__init__` of a protocol subclass is set to `_no_init_or_replace_init`. 1783 # The first instantiation of the subclass will call `_no_init_or_replace_init` which 1784 # searches for a proper new `__init__` in the MRO. The new `__init__` 1785 # replaces the subclass' old `__init__` (ie `_no_init_or_replace_init`). Subsequent 1786 # instantiation of the protocol subclass will thus use the new 1787 # `__init__` and no longer call `_no_init_or_replace_init`. 1788 for base in cls.__mro__: 1789 init = base.__dict__.get('__init__', _no_init_or_replace_init) 1790 if init is not _no_init_or_replace_init: 1791 cls.__init__ = init 1792 break 1793 else: 1794 # should not happen 1795 cls.__init__ = object.__init__ 1796 1797 cls.__init__(self, *args, **kwargs)
70class NoOpOrthogonalizer: 71 ''' 72 No op class (doesn't do anything) 73 74 This class conforms to `Orthogonalizer` protocol. 75 ''' 76 def __init__(self) -> None: 77 pass 78 79 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 80 '''Returns the input array.''' 81 return my_array
No op class (doesn't do anything)
This class conforms to Orthogonalizer
protocol.
84class EuclideanL2Orthogonalizer: 85 ''' 86 Orthogonalizes the basis in the standard Euclidean L2 inner product, i.e., 87 the output basis will satisfy 88 $$\\boldsymbol \\Phi_{\\*}^T \\boldsymbol \\Phi_{\\*} = \\mathbf{I}.$$ 89 90 This class conforms to `Orthogonalizer` protocol. 91 ''' 92 def __init__(self, qrFnc=None) -> None: 93 ''' 94 Constructor 95 Args: 96 97 qrFnc: a callable to use for computing the QR decomposition. 98 IMPORTANT: must conform to the API of 99 [np.linalg.qr](https://numpy.org/doc/stable/reference/generated/numpy.linalg.qr.html). 100 If `None`, internally we use `np.linalg.qr`. 101 Note: this is useful when you want to use a custom qr, for example when your snapshots are 102 distributed with MPI, or maybe you have a fancy qr function that you can use. 103 104 ''' 105 self.__qr_picked = np.linalg.qr if qrFnc is None else qrFnc 106 107 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 108 ''' 109 Orthogonalizes the input matrix and returns the orthogonalized matrix. 110 111 Args: 112 my_array (np.ndarray): The input matrix to be orthogonalized. 113 114 Returns: 115 np.ndarray: The orthogonalized matrix. 116 ''' 117 my_array, _ = self.__qr_picked(my_array, mode='reduced') 118 return my_array
Orthogonalizes the basis in the standard Euclidean L2 inner product, i.e., the output basis will satisfy $$\boldsymbol \Phi_{*}^T \boldsymbol \Phi_{*} = \mathbf{I}.$$
This class conforms to Orthogonalizer
protocol.
92 def __init__(self, qrFnc=None) -> None: 93 ''' 94 Constructor 95 Args: 96 97 qrFnc: a callable to use for computing the QR decomposition. 98 IMPORTANT: must conform to the API of 99 [np.linalg.qr](https://numpy.org/doc/stable/reference/generated/numpy.linalg.qr.html). 100 If `None`, internally we use `np.linalg.qr`. 101 Note: this is useful when you want to use a custom qr, for example when your snapshots are 102 distributed with MPI, or maybe you have a fancy qr function that you can use. 103 104 ''' 105 self.__qr_picked = np.linalg.qr if qrFnc is None else qrFnc
Constructor
Arguments:
- qrFnc: a callable to use for computing the QR decomposition.
IMPORTANT: must conform to the API of
np.linalg.qr.
If
None
, internally we usenp.linalg.qr
. Note: this is useful when you want to use a custom qr, for example when your snapshots are distributed with MPI, or maybe you have a fancy qr function that you can use.
107 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 108 ''' 109 Orthogonalizes the input matrix and returns the orthogonalized matrix. 110 111 Args: 112 my_array (np.ndarray): The input matrix to be orthogonalized. 113 114 Returns: 115 np.ndarray: The orthogonalized matrix. 116 ''' 117 my_array, _ = self.__qr_picked(my_array, mode='reduced') 118 return my_array
Orthogonalizes the input matrix and returns the orthogonalized matrix.
Arguments:
- my_array (np.ndarray): The input matrix to be orthogonalized.
Returns:
np.ndarray: The orthogonalized matrix.
121class EuclideanVectorWeightedL2Orthogonalizer: 122 ''' 123 Orthogonalizes the basis in vector-weighted Euclidean L2 inner product, 124 i.e., the output basis will satisfy 125 $$\\boldsymbol \\Phi_{\\*}^T \\mathrm{diag}(\\mathbf{w})\\boldsymbol \\Phi_{\\*} = \\mathbf{I},$$ 126 where $\\mathbf{w}$ is the weighting vector. Typically, this inner product 127 is used for orthogonalizing with respect to cell volumes 128 129 This class conforms to `Orthogonalizer` protocol. 130 ''' 131 def __init__(self, weighting_vector: np.ndarray, qrFnc=None) -> None: 132 ''' 133 Constructor for the EuclideanVectorWeightedL2Orthogonalizer that 134 initializes the orthogonalizer with the provided weighting vector and 135 an optional custom QR decomposition function. 136 Args: 137 weighting_vector (np.ndarray): a 1-D NumPy array that the matrix will be orthogonalized against. The 138 length of the array must match the number of rows in the matrix that will be orthogonalized. 139 qrFnc: a callable to use for computing the QR decomposition. 140 IMPORTANT: must conform to the API of 141 [np.linalg.qr](https://numpy.org/doc/stable/reference/generated/numpy.linalg.qr.html). 142 If `None`, internally we use `np.linalg.qr`. 143 Note: this is useful when you want to use a custom qr, for example when your snapshots are 144 distributed with MPI, or maybe you have a fancy qr function that you can use. 145 ''' 146 147 self.__weighting_vector = weighting_vector 148 self.__qr_picked = np.linalg.qr if qrFnc is None else qrFnc 149 150 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 151 ''' 152 Orthogonalizes the input matrix using the specified weighting vector and returns the orthogonalized matrix. 153 154 Args: 155 my_array (np.ndarray): The input matrix to be orthogonalized. 156 157 Returns: 158 np.ndarray: The orthogonalized matrix. 159 ''' 160 assert my_array.shape[0] == self.__weighting_vector.size, "Weighting vector does not match basis size" 161 tmp = scipy.sparse.diags(np.sqrt(self.__weighting_vector)) @ my_array 162 my_array, _ = self.__qr_picked(tmp, mode='reduced') 163 my_array = scipy.sparse.diags(np.sqrt(1./self.__weighting_vector)) @ my_array 164 return my_array
Orthogonalizes the basis in vector-weighted Euclidean L2 inner product, i.e., the output basis will satisfy $$\boldsymbol \Phi_{*}^T \mathrm{diag}(\mathbf{w})\boldsymbol \Phi_{*} = \mathbf{I},$$ where $\mathbf{w}$ is the weighting vector. Typically, this inner product is used for orthogonalizing with respect to cell volumes
This class conforms to Orthogonalizer
protocol.
131 def __init__(self, weighting_vector: np.ndarray, qrFnc=None) -> None: 132 ''' 133 Constructor for the EuclideanVectorWeightedL2Orthogonalizer that 134 initializes the orthogonalizer with the provided weighting vector and 135 an optional custom QR decomposition function. 136 Args: 137 weighting_vector (np.ndarray): a 1-D NumPy array that the matrix will be orthogonalized against. The 138 length of the array must match the number of rows in the matrix that will be orthogonalized. 139 qrFnc: a callable to use for computing the QR decomposition. 140 IMPORTANT: must conform to the API of 141 [np.linalg.qr](https://numpy.org/doc/stable/reference/generated/numpy.linalg.qr.html). 142 If `None`, internally we use `np.linalg.qr`. 143 Note: this is useful when you want to use a custom qr, for example when your snapshots are 144 distributed with MPI, or maybe you have a fancy qr function that you can use. 145 ''' 146 147 self.__weighting_vector = weighting_vector 148 self.__qr_picked = np.linalg.qr if qrFnc is None else qrFnc
Constructor for the EuclideanVectorWeightedL2Orthogonalizer that initializes the orthogonalizer with the provided weighting vector and an optional custom QR decomposition function.
Arguments:
- weighting_vector (np.ndarray): a 1-D NumPy array that the matrix will be orthogonalized against. The length of the array must match the number of rows in the matrix that will be orthogonalized.
- qrFnc: a callable to use for computing the QR decomposition.
IMPORTANT: must conform to the API of
np.linalg.qr.
If
None
, internally we usenp.linalg.qr
. Note: this is useful when you want to use a custom qr, for example when your snapshots are distributed with MPI, or maybe you have a fancy qr function that you can use.
150 def orthogonalize(self, my_array: np.ndarray) -> np.ndarray: 151 ''' 152 Orthogonalizes the input matrix using the specified weighting vector and returns the orthogonalized matrix. 153 154 Args: 155 my_array (np.ndarray): The input matrix to be orthogonalized. 156 157 Returns: 158 np.ndarray: The orthogonalized matrix. 159 ''' 160 assert my_array.shape[0] == self.__weighting_vector.size, "Weighting vector does not match basis size" 161 tmp = scipy.sparse.diags(np.sqrt(self.__weighting_vector)) @ my_array 162 my_array, _ = self.__qr_picked(tmp, mode='reduced') 163 my_array = scipy.sparse.diags(np.sqrt(1./self.__weighting_vector)) @ my_array 164 return my_array
Orthogonalizes the input matrix using the specified weighting vector and returns the orthogonalized matrix.
Arguments:
- my_array (np.ndarray): The input matrix to be orthogonalized.
Returns:
np.ndarray: The orthogonalized matrix.