using FEMBase
From the beginning of a project we had a clear concept in our mind: "everything is a field". That is, everything can vary temporally and spatially. We think that constant is just a special case of field which does not vary in temporal nor spatial direction. Fields can vary in spatial direction, i.e. can be either constant or variable, and in temporal direction, i.e. can be time variant or time invariant. From this pondering we can think that there exists four kind of (discrete) fields:
- discrete, constant, time invariant (DCTI)
- discrete, variable, time invariant (DVTI)
- discrete, constant, time variant (DCTV)
- discrete, variable, time variant (DVTV)
Discrete, in this context, means that field is defined in point-wise in
Then we have discrete, variable, time invariant fields DVTI
. Notice
the use of Tuple
when defining field.
b = DVTI( (1.0, 2.0) )
Discrete, constant, time variant field DCTV
is constant in spatial
direction =>
syntax is used. New values can be
added to field using function update!
. If there already exists a
value for that particular time, it will be overridden. It is assumed that
content of field in time direction is monotonically increasing, i.e.
For the sake of clarity let's also mention that update!
works for time
invariant fields as well if content needs to be updated.
c = DCTV(0.0 => 1.0, 1.0 => 2.0)
update!(c, 2.0 => 3.0)
Discrete, variable, time variant DVTV
field is the most general
one, allowing values of field to vary in both spatial and time direction.
d = DVTV(0.0 => (1.0, 2.0), 1.0 => (2.0, 3.0))
update!(d, 2.0 => (3.0, 4.0))
In examples above, all fields defined was scalar fields. Defining vector or tensor fields goes in the same way. The only difference is that now we define vectors and tensors instead of a single scalar value. They can vary in spatial and time direction in the same way than scalar fields. Here is example of defining the abovementioned vector fields:
a = DCTI([1.0, 2.0])
b = DVTI(([1.0, 2.0], [2.0, 3.0]))
c = DCTV(0.0 => [1.0, 2.0], 1.0 => [2.0, 3.0])
d = DVTV(0.0 => ([1.0, 2.0], [2.0, 3.0]), 1.0 => ([2.0, 3.0], [3.0, 4.0]))
Accessing fields in time direction is done using a function interpolate
.
For example, if we have (constant)
c = DCTV(0.0 => [1.0,2.0], 1.0 => [3.0,4.0])
interpolate(c, 0.5)
If field is spatially varying, a Tuple
will be returned, having one value
for each "node". This can then be interpolated in spatial direction, typically
using basis functions defined on element, i.e.
d = DVTV(0.0 => (1.0,2.0), 1.0 => (3.0,4.0))
interpolate(d, 0.5)
Although the two fields above looks quite same, the key difference is that in
DCTV field we have a constant vectorial value (defined using square brackets [])
and in latter DVTV field we have a scalar value (defined using round brackets)
changing in spatial direction from 1.0 to 2.0 at time
If a field is time invariant, interpolation in time direction returns a trivial solution:
interpolate(DCTI(1.0), 0.5)
interpolate(DVTI((1.0,2.0)), 0.5)
For spatially varying fields, one can access to ith element using getindex:
a = DVTI((1.0,2.0))
getindex(a, 1)
For field varying both temporally and spatially, one has first to interpolate in time direction to get iterable tuple:
d = DVTV(0.0 => (1.0,2.0), 1.0 => (3.0,4.0))
result = interpolate(d, 0.5)
getindex(result, 1)
Internally, each field is a subtype of AbstractField
having a field
data
, which be accessed directly for hacking purposes.
d.data
Continuous fields may be useful when defining analytical boundary conditions.
For that we have CVTV
, where "C" stands for continuous.
f(xi,t) = xi[1]*xi[2]*t
g = CVTV(f)
g((1.0,2.0), 3.0)
Usually it is assumed that size of length of discrete field matches to the number of basis functions of a single element, typically something small like 1-10.
However, there might be cases where it is more practical to define field in a
sense that indexing is not continuous or starting from 1. For example, we might
want to define field common for a set of elements. In that case it's natural to
think that each index in field corresponds to the certain id-number of node. For
example, if we have a triangle element connecting nodes 1, 1000 and 100000, we
still want to access that field naturally using getindex
, e.g. f[1]
, f[1000]
and f[100000]
. In that case, more appropriate internal structure for field is
based on a dictionary, not tuple.
It only makes sense to define dictionary fields for spatially varying fields.
Two new fields are introduced: DVTId
and DVTVd
, where the
last character "d" stands for "dictionary".
Keep on mind, that this type of field has one restriction. If and when this field is typically defined on nodes of several elements, field must be continuous between elements. That is, if field value in node 1000 is for example 1.0, then it's 1.0 in all elements connecting to that node. To define jumps on field, one must define field element-wise.
Define eg. "geometry" for nodes 1,1000,100000:
X = Dict(1=>[0.0,0.0], 1000=>[1.0,0.0], 100000=>[1.0,1.0])
G = DVTId(X)
G[1], G[1000], G[100000]
Interpolation in time directions works in a same way than with other fields depends from time.
Y = Dict(1=>[1.0,1.0], 1000=>[2.0,1.0], 100000=>[2.0,2.0])
F = DVTVd(0.0 => X, 1.0 => Y)
interpolate(F,0.5)[100000]
Now we have introduced total of 7 fields: DCTI, DCTV, DVTI, DVTV, CVTV, DVTId,
DVTVd. A good question arises that how to remember all this stuff and is it
even necessary? Luckily not, because one can use a single constructor called
field
to create all kind of fields. Type of field is inspected from
data type. It's not necessary to remember all this technical stuff, just declare
new field using more of less intuitive syntax and field
-function.
f1 = field(1)
f2 = field(1, 2)
f3 = field(0.0 => 1)
f4 = field(0.0 => (1, 2), 1.0 => (2, 3))
f5 = field((xi,t) -> xi[1]*t)
f6 = field(1 => 1, 2 => 2)
f7 = field(0.0 => (1=>1, 10=>2), 1.0 => (1=>2,10=>3))
If the FEMBase ones are not enough, it's always possible to define new ones.
Minimum requirements is that field is a subtype of AbstractField
and
interpolate
, getindex
, has been defined to it. Field can, for example
fetch data from random.org or market stocks, read data from hard drive or add
some stochastics behavior to it.
CurrentModule = FEMBase
DocTestSetup = quote
using FEMBase
end
AbstractField
DCTI
DVTI
DCTV
DVTV
CVTV
DVTId
DVTVd
These functions needs to be defined when developing new fields:
new_field
update_field!
interpolate_field
field(x)
update!(field::F, data) where {F<:AbstractField}
interpolate(field::F, time) where {F<:AbstractField}