| Title: | The Exterior Calculus |
|---|---|
| Description: | Provides functionality for working with tensors, alternating forms, wedge products, Stokes's theorem, and related concepts from the exterior calculus. Uses 'disordR' discipline (Hankin, 2022, <doi:10.48550/arXiv.2210.03856>). The canonical reference would be M. Spivak (1965, ISBN:0-8053-9021-9) "Calculus on Manifolds". To cite the package in publications please use Hankin (2022) <doi:10.48550/arXiv.2210.17008>. |
| Authors: | Robin K. S. Hankin [aut, cre] (ORCID: <https://orcid.org/0000-0001-5982-0415>) |
| Maintainer: | Robin K. S. Hankin <[email protected]> |
| License: | GPL-2 |
| Version: | 1.2-4 |
| Built: | 2026-05-10 20:12:38 UTC |
| Source: | https://github.com/robinhankin/stokes |
Provides functionality for working with tensors, alternating forms, wedge products, Stokes's theorem, and related concepts from the exterior calculus. Uses 'disordR' discipline (Hankin, 2022, <doi:10.48550/arXiv.2210.03856>). The canonical reference would be M. Spivak (1965, ISBN:0-8053-9021-9) "Calculus on Manifolds". To cite the package in publications please use Hankin (2022) <doi:10.48550/arXiv.2210.17008>.
The DESCRIPTION file:
| Package: | stokes |
| Type: | Package |
| Title: | The Exterior Calculus |
| Version: | 1.2-4 |
| Depends: | R (>= 4.2.0) |
| Suggests: | covr, Deriv, knitr, magrittr, markdown, quadform, rmarkdown, testthat, bookdown |
| VignetteBuilder: | knitr |
| Imports: | disordR (>= 0.9-7), methods, partitions, permutations (>= 1.1-2), spray (>= 1.0-26) |
| Authors@R: | person( given=c("Robin", "K. S."), family="Hankin", role = c("aut","cre"), email="[email protected]", comment = c(ORCID = "0000-0001-5982-0415")) |
| Maintainer: | Robin K. S. Hankin <[email protected]> |
| Description: | Provides functionality for working with tensors, alternating forms, wedge products, Stokes's theorem, and related concepts from the exterior calculus. Uses 'disordR' discipline (Hankin, 2022, <doi:10.48550/arXiv.2210.03856>). The canonical reference would be M. Spivak (1965, ISBN:0-8053-9021-9) "Calculus on Manifolds". To cite the package in publications please use Hankin (2022) <doi:10.48550/arXiv.2210.17008>. |
| License: | GPL-2 |
| LazyData: | yes |
| URL: | https://github.com/RobinHankin/stokes, https://robinhankin.github.io/stokes/ |
| BugReports: | https://github.com/RobinHankin/stokes/issues |
| RoxygenNote: | 7.3.3 |
| Encoding: | UTF-8 |
| Config/pak/sysreqs: | libgmp3-dev libicu-dev |
| Repository: | https://robinhankin.r-universe.dev |
| Date/Publication: | 2026-05-10 12:11:22 UTC |
| RemoteUrl: | https://github.com/robinhankin/stokes |
| RemoteRef: | HEAD |
| RemoteSha: | 332ce2a540f75f2d146242323c5283305a6b2eff |
| Author: | Robin K. S. Hankin [aut, cre] (ORCID: <https://orcid.org/0000-0001-5982-0415>) |
Index of help topics:
Alt Alternating multilinear forms
as.1form Coerce vectors to 1-forms
coeffs Extract and manipulate coefficients
consolidate Various low-level helper functions
contract Contractions of k-forms
dovs Dimension of the underlying vector space
dx Elementary forms in three-dimensional space
ex Basis vectors in three-dimensional space
hodge Hodge star operator
inner Inner product operator
issmall Is a form zero to within numerical precision?
keep Keep or drop variables
kform k-forms
kinner Inner product of two kforms
ktensor k-tensors
Ops.kform Arithmetic Ops Group Methods for 'kform' and
'ktensor' objects
phi Elementary tensors
print.stokes Print methods for k-tensors and k-forms
rform Random kforms and ktensors
scalar Scalars and losing attributes
stokes-package The Exterior Calculus
summary.stokes Summaries of tensors and alternating forms
symbolic Symbolic form
tensorprod Tensor products of k-tensors
transform Linear transforms of k-forms
vector_cross_product The Vector cross product
volume The volume element
wedge Wedge products
zap Zap small values in k-forms and k-tensors
zero Zero tensors and zero forms
Generally in the package, arguments that are -forms are
denoted K, -tensors by U, and spray objects by
S. Multilinear maps (which may be either -forms or
-tensors) are denoted by M.
Robin K. S. Hankin [aut, cre] (ORCID: <https://orcid.org/0000-0001-5982-0415>)
Maintainer: Robin K. S. Hankin <[email protected]>
M. Spivak 1971. Calculus on manifolds, Addison-Wesley.
R. K. S. Hankin 2022. “Disordered vectors in R: introducing the disordR package.” https://arxiv.org/abs/2210.03856.
R. K. S. Hankin 2022. “Sparse arrays in R: the spray package. https://arxiv.org/abs/2210.03856.”
## Some k-tensors: U1 <- as.ktensor(matrix(1:15,5,3)) U2 <- as.ktensor(cbind(1:3,2:4),1:3) ## Coerce a tensor to functional form, here mapping V^3 -> R (here V=R^15): as.function(U1)(matrix(rnorm(45),15,3)) ## Tensor product is tensorprod() or %X%: U1 %X% U2 ## A k-form is an alternating k-tensor: K1 <- as.kform(cbind(1:5,2:6),rnorm(5)) K2 <- kform_general(3:6,2,1:6) K3 <- rform(9,3,9) ## The distributive law is true (K1 + K2) ^ K3 == K1 ^ K3 + K2 ^ K3 # TRUE to numerical precision ## Wedge product is associative (non-trivial): (K1 ^ K2) ^ K3 K1 ^ (K2 ^ K3) ## k-forms can be coerced to a function and wedge product: f <- as.function(K1 ^ K2 ^ K3) ## E is a a random point in V^k: E <- matrix(rnorm(63),9,7) ## f() is alternating: f(E) f(E[,7:1]) ## The package blurs the distinction between symbolic and numeric computing: dx <- as.kform(1) dy <- as.kform(2) dz <- as.kform(3) dx ^ dy ^ dz K3 ^ dx ^ dy ^ dz## Some k-tensors: U1 <- as.ktensor(matrix(1:15,5,3)) U2 <- as.ktensor(cbind(1:3,2:4),1:3) ## Coerce a tensor to functional form, here mapping V^3 -> R (here V=R^15): as.function(U1)(matrix(rnorm(45),15,3)) ## Tensor product is tensorprod() or %X%: U1 %X% U2 ## A k-form is an alternating k-tensor: K1 <- as.kform(cbind(1:5,2:6),rnorm(5)) K2 <- kform_general(3:6,2,1:6) K3 <- rform(9,3,9) ## The distributive law is true (K1 + K2) ^ K3 == K1 ^ K3 + K2 ^ K3 # TRUE to numerical precision ## Wedge product is associative (non-trivial): (K1 ^ K2) ^ K3 K1 ^ (K2 ^ K3) ## k-forms can be coerced to a function and wedge product: f <- as.function(K1 ^ K2 ^ K3) ## E is a a random point in V^k: E <- matrix(rnorm(63),9,7) ## f() is alternating: f(E) f(E[,7:1]) ## The package blurs the distinction between symbolic and numeric computing: dx <- as.kform(1) dy <- as.kform(2) dz <- as.kform(3) dx ^ dy ^ dz K3 ^ dx ^ dy ^ dz
Converts a -tensor to alternating form
Alt(S,give_kform)Alt(S,give_kform)
S |
A multilinear form, an object of class |
give_kform |
Boolean, with default |
Given a -tensor , we have
Thus for example if :
and it is reasonably easy to see that
is alternating, in the sense that
Function Alt() is intended to take and return an object of
class ktensor; but if given a kform object, it just
returns its argument unchanged.
A short vignette is provided with the package: type
vignette("Alt") at the commandline.
Returns an alternating -tensor. To work with -forms,
which are a much more efficient representation of alternating tensors,
use as.kform().
Robin K. S. Hankin
(X <- ktensor(spray(rbind(1:3),6))) Alt(X) Alt(X,give_kform=TRUE) S <- as.ktensor(expand.grid(1:3,1:3),rnorm(9)) S Alt(S) issmall(Alt(S) - Alt(Alt(S))) # should be TRUE; Alt() is idempotent a <- rtensor() V <- matrix(rnorm(21),ncol=3) LHS <- as.function(Alt(a))(V) RHS <- as.function(Alt(a,give_kform=TRUE))(V) c(LHS=LHS,RHS=RHS,diff=LHS-RHS)(X <- ktensor(spray(rbind(1:3),6))) Alt(X) Alt(X,give_kform=TRUE) S <- as.ktensor(expand.grid(1:3,1:3),rnorm(9)) S Alt(S) issmall(Alt(S) - Alt(Alt(S))) # should be TRUE; Alt() is idempotent a <- rtensor() V <- matrix(rnorm(21),ncol=3) LHS <- as.function(Alt(a))(V) RHS <- as.function(Alt(a,give_kform=TRUE))(V) c(LHS=LHS,RHS=RHS,diff=LHS-RHS)
Given a vector, return the corresponding 1-form; the exterior
derivative of a 0-form (that is, a scalar function). Function
grad() is a synonym.
as.1form(v) grad(v)as.1form(v) grad(v)
v |
A vector with element |
The exterior derivative of a -form is a
-form given by
We can use the facts that
and
to calculate differentials of general -forms. Specifically, if
then
The entry in square brackets is given by grad(). See the
examples for appropriate R idiom.
A one-form
Robin K. S. Hankin
as.1form(1:9) # note ordering of terms as.1form(rnorm(20)) grad(c(4,7)) ^ grad(1:4)as.1form(1:9) # note ordering of terms as.1form(rnorm(20)) grad(c(4,7)) ^ grad(1:4)
Extract and manipulate coefficients of ktensor and kform objects; this using the methods of the spray package.
Functions as.spray() and nterms() are imported from
spray.
To see the coefficients of a kform or ktensor object,
use coeffs(), which returns a disord object (this is
actually spray::coeffs()). Replacement methods also use the
methods of the spray package. Note that disordR
discipline is enforced.
Experimental functionality for “pure” extraction and
replacement is provided, following spray version 1.0-25 or
above. Thus idiom such as a[abs(coeffs(a)) > 0.1] or indeed
a[coeffs(a) < 1] <- 0 should work as expected. Compare the old
idiom which would be something like coeffs(a)[coeffs(a) < 0.1]
<- 0.
Robin K. S. Hankin
(a <- kform_general(5,2,1:10)) coeffs(a) # a disord object coeffs(a)[coeffs(a)%%2==1] <- 100 # replace every odd coeff with 100 a coeffs(a*0) a <- rform() a[coeffs(a) < 5] # experimental a[coeffs(a) > 3] <- 99 # experimental(a <- kform_general(5,2,1:10)) coeffs(a) # a disord object coeffs(a)[coeffs(a)%%2==1] <- 100 # replace every odd coeff with 100 a coeffs(a*0) a <- rform() a[coeffs(a) < 5] # experimental a[coeffs(a) > 3] <- 99 # experimental
Various low-level helper functions used in Alt() and
kform()
consolidate(S) kill_trivial_rows(S) include_perms(S) kform_to_ktensor(S)consolidate(S) kill_trivial_rows(S) include_perms(S) kform_to_ktensor(S)
S |
Object of class |
Low-level helper functions, used in Alt() and kform().
Function consolidate() takes a spray object, and
combines any rows that are identical up to a permutation, respecting
the sign of the permutation.
Function kill_trivial_rows() takes a spray object and
deletes any rows with a repeated entry (which have -forms
identically zero)
Function include_perms() replaces each row of a
spray object with all its permutations, respecting the sign
of the permutation
Function ktensor_to_kform() coerces a -form
to a -tensor
The functions documented here all return a spray object.
Robin K. S. Hankin
(S <- spray(matrix(c(1,1,2,2,1,3,3,1,3,5),ncol=2,byrow=TRUE),5:1)) kill_trivial_rows(S) # (rows "1 1; 5" and and "2 2; 4" killed - identically zero) consolidate(S) # (merges rows "1 3; 3" and "3 1; 2") include_perms(S) # returns a spray object, not an alternating tensor. include_perms(spray(t(matrix(1:4))))(S <- spray(matrix(c(1,1,2,2,1,3,3,1,3,5),ncol=2,byrow=TRUE),5:1)) kill_trivial_rows(S) # (rows "1 1; 5" and and "2 2; 4" killed - identically zero) consolidate(S) # (merges rows "1 3; 3" and "3 1; 2") include_perms(S) # returns a spray object, not an alternating tensor. include_perms(spray(t(matrix(1:4))))
-formsA contraction is a natural linear map from -forms to -forms.
contract(K, v, lose=TRUE) contract_elementary(o, v)contract(K, v, lose=TRUE) contract_elementary(o, v)
K |
A |
o |
Integer-valued vector corresponding to one row of an index matrix |
lose |
Boolean, with default |
v |
A vector; in function |
Given a -form and a vector ,
the contraction of
and is a -form with
provided ; if we specify
.
Function contract_elementary() is a low-level helper function
that translates elementary -forms with coefficient 1 (in the form
of an integer vector corresponding to one row of an index matrix) into
its contraction with .
There is an extensive vignette in the package,
vignette("contract").
Returns an object of class kform.
Robin K. S. Hankin
Steven H. Weintraub 2014. “Differential forms: theory and practice”, Elsevier (Definition 2.2.23, chapter 2, page 77).
contract(as.kform(1:5), 1:8) contract(as.kform(1), 3) # 0-form contract_elementary(c(1,2,5), c(1,2,10,11,71)) ## Now some verification [takes ~10s to run]: #o <- kform(spray(t(replicate(2, sample(9,4))), runif(2))) #V <- matrix(rnorm(36),ncol=4) #jj <- c( # as.function(o)(V), # as.function(contract(o,V[,1,drop=TRUE]))(V[,-1]), # scalar # as.function(contract(o,V[,1:2]))(V[,-(1:2),drop=FALSE]), # as.function(contract(o,V[,1:3]))(V[,-(1:3),drop=FALSE]), # as.function(contract(o,V[,1:4],lose=FALSE))(V[,-(1:4),drop=FALSE]) #) #print(jj) #max(jj) - min(jj) # zero to numerical precisioncontract(as.kform(1:5), 1:8) contract(as.kform(1), 3) # 0-form contract_elementary(c(1,2,5), c(1,2,10,11,71)) ## Now some verification [takes ~10s to run]: #o <- kform(spray(t(replicate(2, sample(9,4))), runif(2))) #V <- matrix(rnorm(36),ncol=4) #jj <- c( # as.function(o)(V), # as.function(contract(o,V[,1,drop=TRUE]))(V[,-1]), # scalar # as.function(contract(o,V[,1:2]))(V[,-(1:2),drop=FALSE]), # as.function(contract(o,V[,1:3]))(V[,-(1:3),drop=FALSE]), # as.function(contract(o,V[,1:4],lose=FALSE))(V[,-(1:4),drop=FALSE]) #) #print(jj) #max(jj) - min(jj) # zero to numerical precision
A -form maps
to the reals, where .
Function dovs() returns , the dimensionality of the
underlying vector space.
dovs(K)dovs(K)
K |
A |
The function itself is almost trivial, returning the maximum of the index matrix. Special dispensation is given for zero-forms and zero tensors, which return zero.
Vignette dovs provides more discussion.
Returns a non-negative integer
Robin K. S. Hankin
dovs(rform()) table(replicate(20, dovs(rform(3))))dovs(rform()) table(replicate(20, dovs(rform(3))))
Objects dx, dy and dz are the three elementary
one-forms on three-dimensional space. These objects can be generated by
running script ‘vignettes/dx.Rmd’, which includes some further
discussion and technical documentation and creates file ‘dx.rda’
which resides in the data/ directory.
data(dx)data(dx)
See vignettes dx and ex for an extended discussion; a
use-case is given in vector_cross_product.
The default print method is a little opaque for these objects. To print them more intuitively, use
options(kform_symbolic_print = "dx")
which is documented at print.Rd.
Robin K. S. Hankin
M. Spivak 1971. Calculus on manifolds, Addison-Wesley
dx hodge(dx) hodge(dx, 3) dx # default print method, not particularly intelligible options(kform_symbolic_print = 'dx') # shows dx dy dz dx dx^dz hodge(dx, 3) as.function(dx)(ex) options(kform_symbolic_print = NULL) # revert to defaultdx hodge(dx) hodge(dx, 3) dx # default print method, not particularly intelligible options(kform_symbolic_print = 'dx') # shows dx dy dz dx dx^dz hodge(dx, 3) as.function(dx)(ex) options(kform_symbolic_print = NULL) # revert to default
Objects ex, ey and ez are the three elementary
one-forms on three-dimensional space, sometimes denoted
.
See vignettes dx and ex for an extended discussion; a
use-case is given in vector_cross_product.
These objects can be generated by running script
‘vignettes/ex.Rmd’, which includes some further discussion and
technical documentation and creates file ‘exeyez.rda’ which resides
in the data/ directory.
Robin K. S. Hankin
M. Spivak 1971. Calculus on manifolds, Addison-Wesley
as.function(dx)(ex) (X <- as.kform(matrix(1:12, nrow=4), c(1,2,7,11))) as.function(X)(cbind(e(2,12), e(6,12), e(10,12)))as.function(dx)(ex) (X <- as.kform(matrix(1:12, nrow=4), c(1,2,7,11))) as.function(X)(cbind(e(2,12), e(6,12), e(10,12)))
Given a -form, return its Hodge dual
hodge(K, n=dovs(K), g, lose=TRUE)hodge(K, n=dovs(K), g, lose=TRUE)
K |
Object of class |
n |
Dimensionality of space, defaulting to the largest element of the index |
g |
Diagonal of the metric tensor, with missing default being the
standard metric of the identity matrix. Currently, only entries of
|
lose |
Boolean, with default |
Given a -form, in an -dimensional space,
return a -form.
Most authors write the Hodge dual of as
or , but Weintraub
uses .
Robin K. S. Hankin
(o <- kform_general(5,2,1:10)) hodge(o) o == hodge(hodge(o)) Faraday <- kform_general(4,2,runif(6)) # Faraday electromagnetic tensor mink <- c(-1,1,1,1) # Minkowski metric hodge(Faraday,g=mink) Faraday == Faraday |> hodge(g=mink) |> hodge(g=mink) |> hodge(g=mink) |> hodge(g=mink) hodge(dx,3) == dy^dz ## Some edge-cases: hodge(scalar(1),2) hodge(zeroform(5),9) hodge(volume(5)) hodge(volume(5),lose=TRUE) hodge(scalar(7),n=9)(o <- kform_general(5,2,1:10)) hodge(o) o == hodge(hodge(o)) Faraday <- kform_general(4,2,runif(6)) # Faraday electromagnetic tensor mink <- c(-1,1,1,1) # Minkowski metric hodge(Faraday,g=mink) Faraday == Faraday |> hodge(g=mink) |> hodge(g=mink) |> hodge(g=mink) |> hodge(g=mink) hodge(dx,3) == dy^dz ## Some edge-cases: hodge(scalar(1),2) hodge(zeroform(5),9) hodge(volume(5)) hodge(volume(5),lose=TRUE) hodge(scalar(7),n=9)
The inner product
inner(M)inner(M)
M |
square matrix |
The inner product of two vectors and
is usually written
or
, but the most general form would
be where is a matrix.
Noting that inner products are multilinear, that is
and ,
we see that the inner product is indeed a multilinear map, that is, a
tensor.
Given a square matrix , function inner(M) returns the
-form that maps to
. Non-square matrices are
effectively padded with zeros.
A short vignette is provided with the package: type
vignette("inner") at the commandline.
Returns a -tensor, an inner product
Robin K. S. Hankin
inner(diag(7)) inner(matrix(1:9,3,3)) ## Compare the following two: Alt(inner(matrix(1:9,3,3))) # An alternating k tensor as.kform(inner(matrix(1:9,3,3))) # Same thing coerced to a kform f <- as.function(inner(diag(7))) X <- matrix(rnorm(14),ncol=2) # random element of (R^7)^2 f(X) - sum(X[,1]*X[,2]) # zero to numerical precision ## verify positive-definiteness: g <- as.function(inner(crossprod(matrix(rnorm(56),8,7)))) stopifnot(g(kronecker(rnorm(7),t(c(1,1))))>0)inner(diag(7)) inner(matrix(1:9,3,3)) ## Compare the following two: Alt(inner(matrix(1:9,3,3))) # An alternating k tensor as.kform(inner(matrix(1:9,3,3))) # Same thing coerced to a kform f <- as.function(inner(diag(7))) X <- matrix(rnorm(14),ncol=2) # random element of (R^7)^2 f(X) - sum(X[,1]*X[,2]) # zero to numerical precision ## verify positive-definiteness: g <- as.function(inner(crossprod(matrix(rnorm(56),8,7)))) stopifnot(g(kronecker(rnorm(7),t(c(1,1))))>0)
Given a -form, return TRUE if it is “small”
issmall(M, tol=1e-8)issmall(M, tol=1e-8)
M |
Object of class |
tol |
Small tolerance, defaulting to |
Returns a logical
Robin K. S. Hankin
o <- kform_general(3,2,runif(3)) M <- matrix(rnorm(9),3,3) discrepancy <- o - pullback(pullback(o,M),solve(M)) discrepancy # print method might imply coefficients are zeros issmall(discrepancy) # should be TRUE is.zero(discrepancy) # might be FALSEo <- kform_general(3,2,runif(3)) M <- matrix(rnorm(9),3,3) discrepancy <- o - pullback(pullback(o,M),solve(M)) discrepancy # print method might imply coefficients are zeros issmall(discrepancy) # should be TRUE is.zero(discrepancy) # might be FALSE
Keep or drop variables
keep(K, yes) discard(K, no)keep(K, yes) discard(K, no)
K |
Object of class |
yes, no
|
Specification of dimensions to either keep (yes) or discard (no) |
Given a -form , we have
, where
. Now, discarding
dimension is equivalent to asserting (or guaranteeing) that
for .
Alternatively, we may say that is
independent of for
. If this is the case, we may ignore any
row in which an appears.
For -forms, discarding (and its dual, keeping) is carried out in
functions discard() and keep(). Function
keep(omega, yes) keeps the terms specified and
discard(omega, no) discards the terms specified.
The functions documented here all return a kform object.
Robin K. S. Hankin
(o <- kform_general(7, 3, seq_len(choose(7, 3)))) keep(o, 1:4) # keeps only terms with dimensions 1-4 discard(o, 1:2) # loses any term with "1" or "2" in the index(o <- kform_general(7, 3, seq_len(choose(7, 3)))) keep(o, 1:4) # keeps only terms with dimensions 1-4 discard(o, 1:2) # loses any term with "1" or "2" in the index
Functionality for dealing with -forms
kform(S) as.kform(M, coeffs, lose=TRUE) kform_basis(n, k) kform_general(W, k, coeffs, lose=TRUE) is.kform(x) d(i) e(i, n) ## S3 method for class 'kform' as.function(x, ...)kform(S) as.kform(M, coeffs, lose=TRUE) kform_basis(n, k) kform_general(W, k, coeffs, lose=TRUE) is.kform(x) d(i) e(i, n) ## S3 method for class 'kform' as.function(x, ...)
n |
Dimension of the vector space |
i |
Integer |
k |
A |
W |
Integer vector of dimensions |
M, coeffs
|
Index matrix and coefficients for a |
S |
Object of class |
lose |
Boolean, with default |
x |
Object of class |
... |
Further arguments, currently ignored |
A -form is an alternating -tensor. In the
package, -forms are represented as sparse arrays
(spray objects), but with a class of c("kform",
"spray"). The constructor function kform() takes a
spray object and returns a kform object: it ensures that
rows of the index matrix are strictly nonnegative integers, have no
repeated entries, and are strictly increasing. Function
as.kform() is more user-friendly.
kform() is the constructor function. It takes a
spray object and returns a kform.
as.kform() also returns a kform but is a bit more
user-friendly than kform().
kform_basis() is a low-level helper function that
returns a matrix whose rows constitute a basis for the vector space
of -forms.
kform_general() returns a kform object with terms
that span the space of alternating tensors.
is.kform() returns TRUE if its argument is a
kform object.
d() is an easily-typed synonym for
as.kform(). The idea is that d(1) = dx,
d(2)=dy, d(5)=dx^5, etc. Also note that, for
example, d(1:3)=dx^dy^dz, the volume form.
Recall that a -tensor is a multilinear map from
to the reals, where is a vector space.
A multilinear -tensor is alternating if it
satisfies
In the package, an object of class kform is an efficient
representation of an alternating tensor.
Function kform_basis() is a low-level helper function that
returns a matrix whose rows constitute a basis for the vector space
of -forms:
and indeed we have:
where is a basis for
.
All functions documented here return a kform object except
as.function.kform(), which returns a function, and
is.kform(), which returns a Boolean, and e(), which
returns a conjugate basis to that of d().
Hubbard and Hubbard use the term “-form”, but Spivak
does not.
Robin K. S. Hankin
Hubbard and Hubbard; Spivak
as.kform(cbind(1:5,2:6), rnorm(5)) kform_general(1:4, 2, coeffs=1:6) # used in electromagnetism K1 <- as.kform(cbind(1:5, 2:6), rnorm(5)) K2 <- kform_general(5:8, 2, 1:6) K1^K2 # or wedge(K1, K2) d(1:3) dx^dy^dz # same thing d(sample(9)) # coeff is +/-1 depending on even/odd permutation of 1:9 f <- as.function(wedge(K1, K2)) E <- matrix(rnorm(32), 8, 4) f(E) + f(E[,c(1,3,2,4)]) # should be zero by alternating property options(kform_symbolic_print = 'd') (d(5)+d(7)) ^ (d(2)^d(5) + 6*d(4)^d(7)) options(kform_symbolic_print = NULL) # revert to defaultas.kform(cbind(1:5,2:6), rnorm(5)) kform_general(1:4, 2, coeffs=1:6) # used in electromagnetism K1 <- as.kform(cbind(1:5, 2:6), rnorm(5)) K2 <- kform_general(5:8, 2, 1:6) K1^K2 # or wedge(K1, K2) d(1:3) dx^dy^dz # same thing d(sample(9)) # coeff is +/-1 depending on even/odd permutation of 1:9 f <- as.function(wedge(K1, K2)) E <- matrix(rnorm(32), 8, 4) f(E) + f(E[,c(1,3,2,4)]) # should be zero by alternating property options(kform_symbolic_print = 'd') (d(5)+d(7)) ^ (d(2)^d(5) + 6*d(4)^d(7)) options(kform_symbolic_print = NULL) # revert to default
Given two -forms and ,
return the inner product
. Here our
underlying vector space is .
The inner product is a symmetric bilinear form defined in two stages.
First, we specify its behaviour on decomposable -forms
and
as
and secondly, we extend to the whole of
through linearity.
kinner(o1, o2, M)kinner(o1, o2, M)
o1, o2
|
Objects of class |
M |
Matrix |
Returns a real number
There is a vignette available: type vignette("kinner") at
the command line.
Robin K. S. Hankin
a <- (2*dx)^(3*dy) b <- (5*dx)^(7*dy) kinner(a, b) det(matrix(c(2*5,0,0,3*7), 2, 2)) # mathematically identical, slight numerical mismatcha <- (2*dx)^(3*dy) b <- (5*dx)^(7*dy) kinner(a, b) det(matrix(c(2*5,0,0,3*7), 2, 2)) # mathematically identical, slight numerical mismatch
Functionality for -tensors
ktensor(S) as.ktensor(M, coeffs) is.ktensor(x) ## S3 method for class 'ktensor' as.function(x, ...)ktensor(S) as.ktensor(M, coeffs) is.ktensor(x) ## S3 method for class 'ktensor' as.function(x, ...)
M, coeffs
|
Matrix of indices and coefficients, as in
|
S |
Object of class |
x |
Object of class |
... |
Further arguments, currently ignored |
A -tensor object is a map from
to the reals , where is a vector space (here
) that satisfies multilinearity:
and
Note that this is not equivalent to linearity over
(see examples).
In the stokes package, -tensors are represented as
sparse arrays (spray objects), but with a class of
c("ktensor", "spray"). This is a natural and efficient
representation for tensors that takes advantage of sparsity using
spray package features.
Function as.ktensor() will coerce a -form to a
-tensor via kform_to_ktensor().
All functions documented here return a ktensor object
except as.function.ktensor(), which returns a function.
Robin K. S. Hankin
Spivak 1961
as.ktensor(cbind(1:4,2:5,3:6),1:4) ## Test multilinearity: k <- 4 n <- 5 u <- 3 ## Define a randomish k-tensor: S <- ktensor(spray(matrix(1+sample(u*k)%%n,u,k),seq_len(u))) ## And a random point in V^k: E <- matrix(rnorm(n*k),n,k) E1 <- E2 <- E3 <- E x1 <- rnorm(n) x2 <- rnorm(n) r1 <- rnorm(1) r2 <- rnorm(1) # change one column: E1[,2] <- x1 E2[,2] <- x2 E3[,2] <- r1*x1 + r2*x2 f <- as.function(S) r1*f(E1) + r2*f(E2) -f(E3) # should be small ## Note that multilinearity is different from linearity: r1*f(E1) + r2*f(E2) - f(r1*E1 + r2*E2) # not small!as.ktensor(cbind(1:4,2:5,3:6),1:4) ## Test multilinearity: k <- 4 n <- 5 u <- 3 ## Define a randomish k-tensor: S <- ktensor(spray(matrix(1+sample(u*k)%%n,u,k),seq_len(u))) ## And a random point in V^k: E <- matrix(rnorm(n*k),n,k) E1 <- E2 <- E3 <- E x1 <- rnorm(n) x2 <- rnorm(n) r1 <- rnorm(1) r2 <- rnorm(1) # change one column: E1[,2] <- x1 E2[,2] <- x2 E3[,2] <- r1*x1 + r2*x2 f <- as.function(S) r1*f(E1) + r2*f(E2) -f(E3) # should be small ## Note that multilinearity is different from linearity: r1*f(E1) + r2*f(E2) - f(r1*E1 + r2*E2) # not small!
kform and ktensor
objectsAllows arithmetic operators to be used for -forms and
-tensors such as addition, multiplication, etc, where
defined.
## S3 method for class 'kform' Ops(e1, e2 = NULL) ## S3 method for class 'ktensor' Ops(e1, e2 = NULL)## S3 method for class 'kform' Ops(e1, e2 = NULL) ## S3 method for class 'ktensor' Ops(e1, e2 = NULL)
e1, e2
|
Objects of class |
The functions Ops.kform() and Ops.ktensor() pass unary
and binary arithmetic operators (“+”, “-”,
“*”, “/” and “^”) to the
appropriate specialist function by coercing to spray objects.
Binary addition and subtraction, as in x+y or x-y, are
performed by spray::spray_plus_spray().
For wedge products of -forms, use wedge() or
%^% or ^; and for tensor products of
-tensors, use tensorprod() or %X%.
All functions documented here return an object of class
kform or ktensor.
A plain asterisk, “*” behaves differently for ktensors
and kforms.
Given two ktensors T1 and T2, then
“T1*T2” will return their tensor product,
tensorprod(T1,T2). This on the grounds that the idiom has only
one natural interpretation. But its use is discouraged (use
%X% or tensorprod() instead). An asterisk can also be
used to multiply a tensor by a scalar, as in T1*5, in which the
scalar is interpreted as a -tensor.
Given two kforms K1 and K2, asterisks behave
differently. An asterisk cannot be used to multiply K1 and
K2: K1*K2 will return an error. Multiplication by
scalars is OK, as in K1*6.
Powers simply do not make sense for alternating forms: if we wish to
define “” (that is, the “square” of ),
then the only natural interpretation of this is S^S [that is,
wedge(S,S)], but this is zero identically. In the package,
therefore, the caret (“^”) exclusively evaluates the
wedge product; note that %^% is also acceptable.
Here the caret is interpreted consistently as a wedge product, and if
one of the factors is numeric it is interpreted as a zero-form (that
is, a scalar). Thus S^2 = wedge(S,2) = 2^S = S*2 = S+S, and
indeed S^n == S*n. Caveat emptor! If S is a kform
object, it is very tempting [but incorrect] to interpret
“S^3” as something like “S to the power
3”. See also the note at Ops.clifford in the
clifford package.
Powers are not currently implemented for ktensors, but this might be implemented in the future: a ktensor to the power zero is the 0-tensor with unit coefficient.
Note that one has to take care with order of operations if we mix
^ with *. For example, dx ^ (6*dy) is perfectly
acceptable; but (dx ^ 6)*dy) will return an error, as will the
unbracketed form dx ^ 6 * dy. In the second case we attempt to
use an asterisk to multiply two k-forms, which triggers the error.
Robin K. S. Hankin
## dx_1 ^ dx_2 + 6dx_5 ^ dx_6: as.kform(1) ^ as.kform(2) + 6*as.kform(5) ^ as.kform(6) k1 <- kform_general(4, 2, rnorm(6)) k2 <- kform_general(4, 2, rnorm(6)) E <- matrix(rnorm(8), 4, 2) as.function(k1 + k2)(E) ## verify linearity, here 2*k1 + 3*k2: as.function(2*k1 + 3*k2)(E)-(2*as.function(k1)(E) + 3*as.function(k2)(E)) ## should be small## dx_1 ^ dx_2 + 6dx_5 ^ dx_6: as.kform(1) ^ as.kform(2) + 6*as.kform(5) ^ as.kform(6) k1 <- kform_general(4, 2, rnorm(6)) k2 <- kform_general(4, 2, rnorm(6)) E <- matrix(rnorm(8), 4, 2) as.function(k1 + k2)(E) ## verify linearity, here 2*k1 + 3*k2: as.function(2*k1 + 3*k2)(E)-(2*as.function(k1)(E) + 3*as.function(k2)(E)) ## should be small
Creates the elementary tensors or tensor products of elementary tensors
phi(n)phi(n)
n |
Vector of strictly non-negative integers |
If is the standard basis for
then is defined so
that . phi(n) returns
.
If n is a vector of strictly positive integers, then
phi(n) returns the tensor cross product of
applied to the individual elements of n [which is a lot
easier and more obvious than it sounds].
There is a vignette, phi
Robin K. S. Hankin
phi(6) phi(6:8) v <- sample(9) phi(v) == Reduce("%X%",sapply(v,phi))phi(6) phi(6:8) v <- sample(9) phi(v) == Reduce("%X%",sapply(v,phi))
-tensors and -formsPrint methods for objects with options for printing in matrix form or multivariate polynomial form
## S3 method for class 'kform' print(x, ...) ## S3 method for class 'ktensor' print(x, ...)## S3 method for class 'kform' print(x, ...) ## S3 method for class 'ktensor' print(x, ...)
x |
|
... |
Further arguments (currently ignored) |
Printing is dispatched to print.ktensor() and
print.kform() depending on its argument. Special dispensation is
given for the zero object.
Although -forms are alternating tensors and thus mathematically
are tensors, they are handled differently.
The default print method uses the spray print methods,
and as such respects the polyform option. However, setting
polyform to TRUE can give misleading output, because
spray objects are interpreted as multivariate polynomials not
differential forms (and in particular uses the caret to signify
powers).
It is much better to use options ktensor_symbolic_print or
kform_symbolic_print instead: the bespoke print methods
print.kform() and print.ktensor() are sensitive to these
options.
For kform objects, if option kform_symbolic_print is
non-null, the print method uses as.symbolic() to give an
alternate way of displaying -tensors and -forms. The
generic non-null value for this option would be “x”
which gives output like “dx1 ^ dx2”. However, it has
two special values: set kform_symbolic_print to
“dx” for output like “dx ^ dz” and
“txyz” for output like “dt ^ dx”, useful
in relativistic physics with a Minkowski metric. See the examples.
For ktensor objects, if option ktensor_symbolic_print is
TRUE, a different system is used. Given a tensor
, for example (where
), the method will give output that looks
like “+3 d4*d1 -5 d2*d2”. I am not entirely happy with
this and it might change in future.
More detail is given at symbolic.Rd and the dx vignette.
Returns its argument invisibly.
For both kform and ktensor objects, the print method
asserts that its argument is a map from to
with . Here,
is the largest element in the index matrix. However, such a map
naturally furnishes a map from to
, provided that via the
natural projection from to
. Formally this would be
. In the case of the zero -form or -tensor,
“n” is to be interpreted as “any ”. See also dovs().
Robin K. S. Hankin
a <- rform() a options(kform_symbolic_print = "x") a options(kform_symbolic_print = "dx") kform(spray(kform_basis(3,2),1:3)) kform(spray(kform_basis(4,2),1:6)) # runs out of symbols options(kform_symbolic_print = "txyz") kform(spray(kform_basis(4,2),1:6)) # standard notation options(kform_symbolic_print = NULL) # revert to default aa <- rform() a options(kform_symbolic_print = "x") a options(kform_symbolic_print = "dx") kform(spray(kform_basis(3,2),1:3)) kform(spray(kform_basis(4,2),1:6)) # runs out of symbols options(kform_symbolic_print = "txyz") kform(spray(kform_basis(4,2),1:6)) # standard notation options(kform_symbolic_print = NULL) # revert to default a
Random -form objects and -tensors,
intended as quick “get you going” examples
rform(terms=9, k=3, n=7, ensure=TRUE, integer=TRUE) rformm(terms=30, k=7, n=20, ensure=TRUE, integer=TRUE) rformmm(terms=90, k=15, n=30, ensure=TRUE, integer=TRUE) rtensor(terms=9, k=3, n=7, integer=TRUE)rform(terms=9, k=3, n=7, ensure=TRUE, integer=TRUE) rformm(terms=30, k=7, n=20, ensure=TRUE, integer=TRUE) rformmm(terms=90, k=15, n=30, ensure=TRUE, integer=TRUE) rtensor(terms=9, k=3, n=7, integer=TRUE)
terms |
Number of distinct terms |
k, n
|
A |
ensure |
Boolean with default |
integer |
Boolean specifying whether the coefficients are integers or not |
Random -form objects and -tensors.
By default, function rform() returns a simple -form;
rformm() and rformmm() return successively more
complicated objects. Note that argument terms is an upper
bound, as the index matrix might contain repeats which are combined.
Function rtensor() returns a random tensor.
All functions documented here return an object of class kform or
ktensor.
Robin K. S. Hankin
(a <- rform()) (b <- rform()) a ^ b a a ^ dx a ^ dx ^ dy (x <- rtensor()) x %X% x(a <- rform()) (b <- rform()) a ^ b a a ^ dx a ^ dx ^ dy (x <- rtensor()) x %X% x
Scalars: -forms and -tensors
scalar(s, kform=TRUE, lose=FALSE) is.scalar(M) `0form`(s=1, lose=FALSE) `0tensor`(s=1, lose=FALSE) ## S3 method for class 'kform' lose(M) ## S3 method for class 'ktensor' lose(M)scalar(s, kform=TRUE, lose=FALSE) is.scalar(M) `0form`(s=1, lose=FALSE) `0tensor`(s=1, lose=FALSE) ## S3 method for class 'kform' lose(M) ## S3 method for class 'ktensor' lose(M)
s |
A scalar value; a number |
kform |
Boolean with default |
M |
Object of class |
lose |
In function |
A -tensor (including -forms) maps vectors
to a scalar. If , then a -tensor maps no vectors
to a scalar, that is, mapping nothing at all to a scalar, or what normal
people would call a plain old scalar. Such forms are created by a
couple of constructions in the package, specifically scalar(),
kform_general(1,0) and contract(). These functions take a
lose argument that behaves much like the drop argument in
base extraction. Functions 0form() and 0tensor() are
wrappers for scalar().
Function lose() takes an object of class ktensor or
kform and, if of arity zero, returns the coefficient.
Note that function kform() always returns a kform
object, it never loses attributes.
There is a slight terminological problem. A -form maps
vectors to the reals: so a -form maps vectors to the
reals. This is what anyone on the planet would call a scalar.
Similarly, a -tensor maps vectors to the reals, and so it
too is a scalar. Mathematically, there is no difference between
-forms and -tensors, but the package print methods make a
distinction:
> scalar(5,kform=TRUE)
An alternating linear map from V^0 to R with V=R^0:
val
= 5
> scalar(5,kform=FALSE)
A linear map from V^0 to R with V=R^0:
val
= 5
>
The functions documented here return an object of class
kform or ktensor, except for is.scalar(), which
returns a Boolean.
Compare zero tensors and zero forms. These are different from
-tensors and -scalars. A zero tensor maps to
the real number zero, and a zero form is an alternating tensor mapping
to zero (so a zero tensor is necessarily alternating). See
zero.Rd.
Robin K. S. Hankin
o <- scalar(5) o lose(o) kform_general(1,0) kform_general(1,0,lose=FALSE)o <- scalar(5) o lose(o) kform_general(1,0) kform_general(1,0,lose=FALSE)
A summary method for tensors and alternating forms, and a print method for summaries.
## S3 method for class 'kform' summary(object, ...) ## S3 method for class 'ktensor' summary(object, ...) ## S3 method for class 'summary.kform' print(x, ...) ## S3 method for class 'summary.ktensor' print(x, ...)## S3 method for class 'kform' summary(object, ...) ## S3 method for class 'ktensor' summary(object, ...) ## S3 method for class 'summary.kform' print(x, ...) ## S3 method for class 'summary.ktensor' print(x, ...)
object, x
|
Object of class |
... |
Further arguments, passed to |
Summary methods for tensors and alternating forms. Uses spray::summary().
Robin K. S. Hankin
a <- rform(100) summary(a) options(kform_symbolic_print = TRUE) summary(a) options(kform_symbolic_print = NULL) # restore defaulta <- rform(100) summary(a) options(kform_symbolic_print = TRUE) summary(a) options(kform_symbolic_print = NULL) # restore default
Returns a character string representing -tensor and
-form objects in symbolic form. Used by the print method if
either option kform_symbolic_print or
ktensor_symbolic_print is non-null.
as.symbolic(M, symbols=letters, d="")as.symbolic(M, symbols=letters, d="")
M |
Object of class |
symbols |
A character vector giving the names of the symbols |
d |
String specifying the appearance of the differential operator |
Spivak (p89), in archetypically terse writing, states:
A function is considered to be a 0-form and
is also written
. If
is
differentiable, then
.
By a minor modification we therefore obtain a 1-form
, defined by
Let us consider in particular the 1-forms
. It is customary to let
denote the function (On we often denote ,
, and by , , and
). This standard notation has obvious disadvantages but it
allows many classical results to be expressed by formulas of equally
classical appearance. Since
, we see that
is just the dual basis to
. Thus every
k-form can be written
Function as.symbolic() uses this format. For completeness, we
add (p77) that -tensors may be expressed in the form
and this form is used for -tensors. The print method for
tensors, print.ktensor(), writes d1 for ,
d2 for [where ].
Returns a “noquote” character string.
Robin K. S. Hankin
(o <- kform_general(3, 2, 1:3)) as.symbolic(o, d="d", symbols=letters[23:26]) (a <- rform(n=50)) as.symbolic(a, symbols=state.abb)(o <- kform_general(3, 2, 1:3)) as.symbolic(o, d="d", symbols=letters[23:26]) (a <- rform(n=50)) as.symbolic(a, symbols=state.abb)
-tensorsTensor products of -tensors
tensorprod(U, ...) tensorprod2(U1, U2)tensorprod(U, ...) tensorprod2(U1, U2)
U, U1, U2
|
Object of class |
... |
Further arguments, currently ignored |
Given a -tensor and an -tensor
, we can form the tensor product , defined as
Package idiom for this includes tensorprod(S,T) and S %X%
T; note that the tensor product is not commutative. Function
tensorprod() can take any number of arguments (the result is
well-defined because the tensor product is associative); it uses
tensorprod2() as a low-level helper function.
The functions documented here all return a ktensor object.
The binary form %X% uses uppercase X to avoid clashing with
%x% which is the Kronecker product in base R.
Taking a tensor product of a -form will return an error.
Robin K. S. Hankin
Spivak 1961
(A <- ktensor(spray(matrix(c(1,1,2,2,3,3),2,3,byrow=TRUE),1:2))) (B <- ktensor(spray(10+matrix(4:9,3,2),5:7))) tensorprod(A,B) A %X% B - B %X% A Va <- matrix(rnorm(9),3,3) Vb <- matrix(rnorm(38),19,2) LHS <- as.function(A %X% B)(cbind(rbind(Va,matrix(0,19-3,3)),Vb)) RHS <- as.function(A)(Va) * as.function(B)(Vb) c(LHS=LHS,RHS=RHS,diff=LHS-RHS)(A <- ktensor(spray(matrix(c(1,1,2,2,3,3),2,3,byrow=TRUE),1:2))) (B <- ktensor(spray(10+matrix(4:9,3,2),5:7))) tensorprod(A,B) A %X% B - B %X% A Va <- matrix(rnorm(9),3,3) Vb <- matrix(rnorm(38),19,2) LHS <- as.function(A %X% B)(cbind(rbind(Va,matrix(0,19-3,3)),Vb)) RHS <- as.function(A)(Va) * as.function(B)(Vb) c(LHS=LHS,RHS=RHS,diff=LHS-RHS)
-formsExpress a -form in terms of linear combinations of the
pullback(K, M) stretch(K, d)pullback(K, M) stretch(K, d)
K |
Object of class |
M |
Matrix of transformation |
d |
Numeric vector representing the diagonal elements of a diagonal matrix |
Function pullback() calculates the pullback of a function. A
vignette is provided at ‘pullback.Rmd’.
Suppose we are given a two-form
and relationships
then we would have
The general situation would be a -form where we would have
giving
The transform() function does all this but it is slow. I am not
100% sure that there isn't a much more efficient way to do such a
transformation. There are a few tests in tests/testthat and a
discussion in the stokes vignette.
Function stretch() carries out the same operation but for
a diagonal matrix. It is much faster than transform(). It can
take a ktensor or a kform object.
The functions documented here return an object of class
kform.
Robin K. S. Hankin
S. H. Weintraub 2019. Differential forms: theory and practice. Elsevier. (Chapter 3)
# Example in the text: K <- as.kform(matrix(c(1,1,2,3),2,2),c(1,5)) M <- matrix(1:9,3,3) pullback(K,M) # Demonstrate that the result can be complicated: M <- matrix(rnorm(25),5,5) pullback(as.kform(1:2),M) # Numerical verification: o <- volume(3) o2 <- pullback(pullback(o,M),solve(M)) max(abs(coeffs(o-o2))) # zero to numerical precision # Following should be zero: pullback(as.kform(1),M)-as.kform(matrix(1:5),c(crossprod(M,c(1,rep(0,4))))) # Following should be TRUE: issmall(pullback(o,crossprod(matrix(rnorm(10),2,5)))) # Some stretch() use-cases: p <- rform() p stretch(p,seq_len(7)) stretch(p,c(1,0,0,1,1,1,1)) # kills dimensions 2 and 3# Example in the text: K <- as.kform(matrix(c(1,1,2,3),2,2),c(1,5)) M <- matrix(1:9,3,3) pullback(K,M) # Demonstrate that the result can be complicated: M <- matrix(rnorm(25),5,5) pullback(as.kform(1:2),M) # Numerical verification: o <- volume(3) o2 <- pullback(pullback(o,M),solve(M)) max(abs(coeffs(o-o2))) # zero to numerical precision # Following should be zero: pullback(as.kform(1),M)-as.kform(matrix(1:5),c(crossprod(M,c(1,rep(0,4))))) # Following should be TRUE: issmall(pullback(o,crossprod(matrix(rnorm(10),2,5)))) # Some stretch() use-cases: p <- rform() p stretch(p,seq_len(7)) stretch(p,c(1,0,0,1,1,1,1)) # kills dimensions 2 and 3
The vector cross product
for is defined
in elementary school as
Function vcp3() is a convenience wrapper for this. However, the
vector cross product may easily be generalized to a product of
-tuples of vectors in , given by
package function vector_cross_product().
Vignette vector_cross_product, supplied with the package, gives
an extensive discussion of vector cross products, including formal
definitions and verification of identities.
vector_cross_product(M) vcp3(u,v)vector_cross_product(M) vcp3(u,v)
M |
Matrix with one more row than column; columns are interpreted as vectors |
u, v
|
Vectors of length 3, representing vectors in |
A joint function profile for vector_cross_product() and
vcp3() is given with the package at
vignette("vector_cross_product").
Returns a vector
Robin K. S. Hankin
vector_cross_product(matrix(1:6,3,2)) M <- matrix(rnorm(30),6,5) LHS <- hodge(as.1form(M[,1])^as.1form(M[,2])^as.1form(M[,3])^as.1form(M[,4])^as.1form(M[,5])) RHS <- as.1form(vector_cross_product(M)) LHS-RHS # zero to numerical precision # Alternatively: hodge(Reduce(`^`,sapply(seq_len(5),function(i){as.1form(M[,i])},simplify=FALSE)))vector_cross_product(matrix(1:6,3,2)) M <- matrix(rnorm(30),6,5) LHS <- hodge(as.1form(M[,1])^as.1form(M[,2])^as.1form(M[,3])^as.1form(M[,4])^as.1form(M[,5])) RHS <- as.1form(vector_cross_product(M)) LHS-RHS # zero to numerical precision # Alternatively: hodge(Reduce(`^`,sapply(seq_len(5),function(i){as.1form(M[,i])},simplify=FALSE)))
The volume element in dimensions
volume(n) is.volume(K, n=dovs(K))volume(n) is.volume(K, n=dovs(K))
n |
Dimension of the space |
K |
Object of class |
Spivak phrases it well (theorem 4.6, page 82):
If has dimension , it follows that
has dimension 1. Thus all alternating
-tensors on are multiples of any non-zero one.
Since the determinant is an example of such a member of
it is not surprising to find it in
the following theorem:
Let be a basis for and
let . If then
(see the examples for numerical verification of this).
Neither the zero -form, nor scalars, are considered to be a
volume element.
Function volume() returns an object of class kform;
function is.volume() returns a Boolean.
Robin K. S. Hankin
M. Spivak 1971. Calculus on manifolds, Addison-Wesley
dx^dy^dz == volume(3) p <- 1 for(i in 1:7){p <- p ^ as.kform(i)} p p == volume(7) # should be TRUE o <- volume(5) M <- matrix(runif(25),5,5) det(M) - as.function(o)(M) # should be zero is.volume(d(1) ^ d(2) ^ d(3) ^ d(4)) is.volume(d(1:9))dx^dy^dz == volume(3) p <- 1 for(i in 1:7){p <- p ^ as.kform(i)} p p == volume(7) # should be TRUE o <- volume(5) M <- matrix(runif(25),5,5) det(M) - as.function(o)(M) # should be zero is.volume(d(1) ^ d(2) ^ d(3) ^ d(4)) is.volume(d(1:9))
Wedge products of -forms
wedge2(K1, K2) wedge(x, ...)wedge2(K1, K2) wedge(x, ...)
K1, K2, x, ...
|
|
Wedge product of -forms.
The functions documented here return an object of class
kform.
In general use, use wedge() or ^ or %^%, as
documented under Ops. Function wedge() uses low-level
helper function wedge2(), which takes only two arguments.
A short vignette is provided with the package: type
vignette("wedge") at the commandline.
Robin K. S. Hankin
k1 <- as.kform(cbind(1:5,2:6),1:5) k2 <- as.kform(cbind(5:7,6:8,7:9),1:3) k3 <- kform_general(1:6,2) a1 <- wedge2(k1,wedge2(k2,k3)) a2 <- wedge2(wedge2(k1,k2),k3) is.zero(a1-a2) # NB terms of a1, a2 in a different order! # This is why wedge(k1,k2,k3) is well-defined. Can also use ^: k1 ^ k2 ^ k3k1 <- as.kform(cbind(1:5,2:6),1:5) k2 <- as.kform(cbind(5:7,6:8,7:9),1:3) k3 <- kform_general(1:6,2) a1 <- wedge2(k1,wedge2(k2,k3)) a2 <- wedge2(wedge2(k1,k2),k3) is.zero(a1-a2) # NB terms of a1, a2 in a different order! # This is why wedge(k1,k2,k3) is well-defined. Can also use ^: k1 ^ k2 ^ k3
-forms and -tensorsEquivalent to zapsmall()
zap(X) ## S3 method for class 'kform' zap(X) ## S3 method for class 'ktensor' zap(X)zap(X) ## S3 method for class 'kform' zap(X) ## S3 method for class 'ktensor' zap(X)
X |
Tensor or |
Given an object of class ktensor or kform, coefficients
close to zero are ‘zapped’, i.e., replaced by ‘0’, using
base::zapsmall().
Note, zap() actually changes the numeric value, it is not just
a print method.
Returns an object of the same class
Robin K. S. Hankin
S <- rform(7) S == zap(S) # should be TRUE because the coeffs are integers (a <- rform()) (b <- rform()*1e-11) a+b zap(a+b)S <- rform(7) S == zap(S) # should be TRUE because the coeffs are integers (a <- rform()) (b <- rform()*1e-11) a+b zap(a+b)
Correct idiom for generating zero -tensors and -forms
zeroform(n) zerotensor(n) is.zero(x) is.empty(x)zeroform(n) zerotensor(n) is.zero(x) is.empty(x)
n |
Arity of the |
x |
Object to be tested for zero |
Returns an object of class kform or ktensor.
Idiom such as as.ktensor(rep(1,5),0) and
as.kform(rep(1,5),0) and indeed as.kform(1:5,0) will
return the zero tensor or -form (in earlier versions of the
package, these were held to be incorrect as the arity of the tensor
was lost).
A -form is not the same thing as a zero tensor. A
-form maps to the reals; a scalar. A zero
tensor maps to zero. Some discussion is given at
scalar.Rd.
Robin K. S. Hankin
zerotensor(5) zeroform(3) x <- rform(k=3) x*0 == zeroform(3) # should be true x == x + zeroform(3) # should be true y <- rtensor(k=3) y*0 == zerotensor(3) # should be true y == y+zerotensor(3) # should be true ## Following idiom is plausible but fails because as.ktensor(coeffs=0) ## and as.kform(coeffs=0) do not retain arity: ## as.ktensor(1+diag(5)) + as.ktensor(rep(1,5),0) # fails ## as.kform(matrix(1:6,2,3)) + as.kform(1:3,0) # also failszerotensor(5) zeroform(3) x <- rform(k=3) x*0 == zeroform(3) # should be true x == x + zeroform(3) # should be true y <- rtensor(k=3) y*0 == zerotensor(3) # should be true y == y+zerotensor(3) # should be true ## Following idiom is plausible but fails because as.ktensor(coeffs=0) ## and as.kform(coeffs=0) do not retain arity: ## as.ktensor(1+diag(5)) + as.ktensor(rep(1,5),0) # fails ## as.kform(matrix(1:6,2,3)) + as.kform(1:3,0) # also fails