agd.Metrics.isotropic
1# Copyright 2020 Jean-Marie Mirebeau, University Paris-Sud, CNRS, University Paris-Saclay 2# Distributed WITHOUT ANY WARRANTY. Licensed under the Apache License, Version 2.0, see http://www.apache.org/licenses/LICENSE-2.0 3 4from .base import Base 5from .. import AutomaticDifferentiation as ad 6import numpy as np 7 8class Isotropic(Base): 9 r""" 10 An Isotropic norm takes the form 11 $$ 12 F(x) = cost * \sqrt{< x,x>}, 13 $$ 14 where cost is a given positive scalar. 15 16 Member fields and __init__ arguments : 17 - cost : an array of arbirary shape (n1,..,nk). 18 - vdim (optional) : the ambient space dimension 19 """ 20 21 def __init__(self,cost,vdim=None): 22 self.cost = ad.asarray(cost) 23 self._vdim = vdim 24 25 @classmethod 26 def from_speed(cls,speed,vdim=None): 27 """ 28 Produces a metric whose cost equals 1/speed. 29 """ 30 return cls(1./speed,vdim) 31 32 def dual(self): 33 other = self.from_speed(self.cost) 34 if hasattr(self,'_vdim'): other._vdim = self._vdim 35 return other 36 37 def norm(self,v): 38 return self.cost*ad.Optimization.norm(v,ord=2,axis=0) 39 40 def gradient(self,v): 41 return self.cost*v/ad.Optimization.norm(v,ord=2,axis=0) 42 43 def is_definite(self): 44 return self.cost>0. 45 46 def anisotropy(self): 47 return 1. 48 49 def cost_bound(self): 50 return self.cost 51 52 @property 53 def vdim(self): 54 if self._vdim is not None: return self._vdim 55 if self.cost.ndim>0: return self.cost.ndim 56 raise ValueError("Could not determine dimension of isotropic metric") 57 58 @vdim.setter 59 def vdim(self,vdim): 60 if self.cost.ndim>0: 61 assert(self.cost.ndim==vdim) 62 else: 63 self._vdim = vdim 64 65 @property 66 def shape(self): return self.cost.shape 67 68 def rotate(self,a): return self 69 def with_cost(self,cost): return Isotropic(cost*self.cost) 70 def with_costs(self,costs): 71 if len(costs)>1 and not np.allclose(costs[1:],costs[0]): 72 raise ValueError("Costs must be identical along all axes for Metrics.Isotropic. Consider using Metrics.Diagonal class") 73 return self.with_cost(costs[0]) 74 75 def flatten(self): return self.cost 76 @classmethod 77 def expand(cls,arr): return cls(arr) 78 79 def to_HFM(self): return self.cost 80 @classmethod 81 def from_HFM(cls,arr): return cls(arr) 82 83 def model_HFM(self): 84 return "Isotropic"+str(self.vdim) 85 86 @classmethod 87 def from_cast(cls,metric): 88 if isinstance(metric,cls): return metric 89 raise ValueError("Isotropic.from_cast error : cannot cast an isotropic metric from ",type(metric)) 90 91 def __iter__(self): 92 yield self.cost 93 94 def make_proj_dual(self): 95 """Returns the projection operator onto the unit ball of the dual metric""" 96 a=1./self.cost 97 norm = np.linalg.norm 98 return lambda x : x / np.maximum(1.,a*norm(x,axis=0))
9class Isotropic(Base): 10 r""" 11 An Isotropic norm takes the form 12 $$ 13 F(x) = cost * \sqrt{< x,x>}, 14 $$ 15 where cost is a given positive scalar. 16 17 Member fields and __init__ arguments : 18 - cost : an array of arbirary shape (n1,..,nk). 19 - vdim (optional) : the ambient space dimension 20 """ 21 22 def __init__(self,cost,vdim=None): 23 self.cost = ad.asarray(cost) 24 self._vdim = vdim 25 26 @classmethod 27 def from_speed(cls,speed,vdim=None): 28 """ 29 Produces a metric whose cost equals 1/speed. 30 """ 31 return cls(1./speed,vdim) 32 33 def dual(self): 34 other = self.from_speed(self.cost) 35 if hasattr(self,'_vdim'): other._vdim = self._vdim 36 return other 37 38 def norm(self,v): 39 return self.cost*ad.Optimization.norm(v,ord=2,axis=0) 40 41 def gradient(self,v): 42 return self.cost*v/ad.Optimization.norm(v,ord=2,axis=0) 43 44 def is_definite(self): 45 return self.cost>0. 46 47 def anisotropy(self): 48 return 1. 49 50 def cost_bound(self): 51 return self.cost 52 53 @property 54 def vdim(self): 55 if self._vdim is not None: return self._vdim 56 if self.cost.ndim>0: return self.cost.ndim 57 raise ValueError("Could not determine dimension of isotropic metric") 58 59 @vdim.setter 60 def vdim(self,vdim): 61 if self.cost.ndim>0: 62 assert(self.cost.ndim==vdim) 63 else: 64 self._vdim = vdim 65 66 @property 67 def shape(self): return self.cost.shape 68 69 def rotate(self,a): return self 70 def with_cost(self,cost): return Isotropic(cost*self.cost) 71 def with_costs(self,costs): 72 if len(costs)>1 and not np.allclose(costs[1:],costs[0]): 73 raise ValueError("Costs must be identical along all axes for Metrics.Isotropic. Consider using Metrics.Diagonal class") 74 return self.with_cost(costs[0]) 75 76 def flatten(self): return self.cost 77 @classmethod 78 def expand(cls,arr): return cls(arr) 79 80 def to_HFM(self): return self.cost 81 @classmethod 82 def from_HFM(cls,arr): return cls(arr) 83 84 def model_HFM(self): 85 return "Isotropic"+str(self.vdim) 86 87 @classmethod 88 def from_cast(cls,metric): 89 if isinstance(metric,cls): return metric 90 raise ValueError("Isotropic.from_cast error : cannot cast an isotropic metric from ",type(metric)) 91 92 def __iter__(self): 93 yield self.cost 94 95 def make_proj_dual(self): 96 """Returns the projection operator onto the unit ball of the dual metric""" 97 a=1./self.cost 98 norm = np.linalg.norm 99 return lambda x : x / np.maximum(1.,a*norm(x,axis=0))
An Isotropic norm takes the form $$ F(x) = cost * \sqrt{< x,x>}, $$ where cost is a given positive scalar.
Member fields and __init__ arguments :
- cost : an array of arbirary shape (n1,..,nk).
- vdim (optional) : the ambient space dimension
26 @classmethod 27 def from_speed(cls,speed,vdim=None): 28 """ 29 Produces a metric whose cost equals 1/speed. 30 """ 31 return cls(1./speed,vdim)
Produces a metric whose cost equals 1/speed.
33 def dual(self): 34 other = self.from_speed(self.cost) 35 if hasattr(self,'_vdim'): other._vdim = self._vdim 36 return other
Dual norm
, mathematically defined by
$N^*(x) = max\{ < x, y> ; N(y)\leq 1 \}$
Norm or quasi-norm defined by the class, often denoted $N$ in mathematical formulas. Unless incorrect data is provided, this member function obeys, for all vectors $u,v\in R^d$ and all $\alpha \geq 0$
- $N(u+v) \leq N(u)+N(v)$
- $N(\alpha u) = \alpha N(u)$
- $N(u)\geq 0$ with equality iff $u=0$.
Broadcasting will occur depending on the shape of $v$ and of the class data.
Anisotropy ratio of the norm
denoted $N$.
Defined as
$$
\max_{|u| = |v| = 1} \frac {N(u)}{N(v)}.
$$
Upper bound on $N(u)$, for any unit vector $u$, where $N$ is the norm
defined by the class.
53 @property 54 def vdim(self): 55 if self._vdim is not None: return self._vdim 56 if self.cost.ndim>0: return self.cost.ndim 57 raise ValueError("Could not determine dimension of isotropic metric")
The ambient vector space dimension, often denoted $d$ in mathematical formulas.
Dimensions of the underlying domain.
Expected to be the empty tuple, or a tuple of length vdim
.
69 def rotate(self,a): return self
Rotation of the norm, by a given rotation matrix. The new unit ball is the direct image of the previous one.
70 def with_cost(self,cost): return Isotropic(cost*self.cost)
Produces a norm $N'$ obeying $N'(x) = N(cost*x)$.
71 def with_costs(self,costs): 72 if len(costs)>1 and not np.allclose(costs[1:],costs[0]): 73 raise ValueError("Costs must be identical along all axes for Metrics.Isotropic. Consider using Metrics.Diagonal class") 74 return self.with_cost(costs[0])
Produces a norm $N'$ defined by $$ N'(x) = N(costs * x) $$ where the multiplication is elementwise.
76 def flatten(self): return self.cost
Flattens and concatenate the member fields into a single array.
Inverse of the flatten member function. Turns a suitable array into a metric.
80 def to_HFM(self): return self.cost
Formats a metric for the HFM library. This may include flattening some symmetric matrices, concatenating with vector fields, and moving the first axis last.
Inverse of the to_HFM member function. Turns a suitable array into a metric.
87 @classmethod 88 def from_cast(cls,metric): 89 if isinstance(metric,cls): return metric 90 raise ValueError("Isotropic.from_cast error : cannot cast an isotropic metric from ",type(metric))
Produces a metric by casting another metric of a compatible type.
95 def make_proj_dual(self): 96 """Returns the projection operator onto the unit ball of the dual metric""" 97 a=1./self.cost 98 norm = np.linalg.norm 99 return lambda x : x / np.maximum(1.,a*norm(x,axis=0))
Returns the projection operator onto the unit ball of the dual metric