Documentation ¶
Index ¶
- func ApplyModalFilter(filter ModalFilter, freq Frequency, data []complex128)
- func DefaultNonLinSolver() nonlin.NewtonKrylov
- func GetFieldName(term string, fieldNames []string) string
- func GetNonLinearFieldExpressions(pattern string, field string, fieldNames []string) string
- func GetPower(pattern string) float64
- func Indicator(x float64) float64
- func IndicatorDeriv(x float64) float64
- func LoadFloat64(fname string) []float64
- func RealPartAsUint8(data []complex128, min float64, max float64) []uint8
- func SaveCsv(fname string, data []CsvData, domainSize []int)
- func SaveFloat64(fname string, data []float64)
- func SortFactors(expr string) string
- func ValidName(name string, model *Model) bool
- func WriteXDMF(xdmfFile string, fields []string, prefix string, num int, domainSize []int)
- type Advection
- func (ad *Advection) AllFieldsExist(m *Model) bool
- func (ad *Advection) Construct(bricks map[string]Brick) Term
- func (ad *Advection) GetName() string
- func (ad *Advection) GradName(comp int) string
- func (ad *Advection) OnStepFinished(t float64)
- func (ad *Advection) PrepareModel(N int, m *Model, FT FourierTransform)
- type BoundTransform
- type Brick
- type BrickPlaceholder
- type ChargeTransport
- type ConservativeNoise
- func (cn *ConservativeNoise) Construct(bricks map[string]Brick) Term
- func (cn *ConservativeNoise) CurrentFieldsAreRegistered(bricks map[string]Brick) bool
- func (cn *ConservativeNoise) GetCurrentName(comp int) string
- func (cn *ConservativeNoise) OnStepFinished(t float64, bricks map[string]Brick)
- func (cn *ConservativeNoise) RequiredDerivedFields(numNodes int) []DerivedField
- type CsvData
- type CsvIO
- type DerivedField
- type DerivedFieldCalc
- type DiffOp
- type DisplacementGetter
- type DivGrad
- type Euler
- type ExplicitPairCorrelationTerm
- type Field
- type FieldDB
- func (fdb *FieldDB) Comment(comment string)
- func (fdb *FieldDB) Load(simID int, timestep int) []Field
- func (fdb *FieldDB) LoadLast(simID int) []Field
- func (fdb *FieldDB) SaveFields(s *Solver, epoch int)
- func (fdb *FieldDB) SetAttr(attr map[string]float64)
- func (fdb *FieldDB) SetTextAttr(attr map[string]string)
- func (fdb *FieldDB) TimeSeries(data map[string]float64, timestep int)
- type Float64IO
- type FourierTransform
- type Frequency
- type GenericFunction
- type GradientCalculator
- type HomogeneousModulusLinElast
- func (h *HomogeneousModulusLinElast) Construct(bricks map[string]Brick) Term
- func (h *HomogeneousModulusLinElast) Force(indicator []complex128) [][]complex128
- func (h *HomogeneousModulusLinElast) Freq(i int) []float64
- func (h *HomogeneousModulusLinElast) OnStepFinished(t float64, bricks map[string]Brick)
- type IdealMixtureTerm
- func (idt *IdealMixtureTerm) ConstructLinear(bricks map[string]Brick) Term
- func (idt *IdealMixtureTerm) ConstructNonLinear(bricks map[string]Brick) Term
- func (idt *IdealMixtureTerm) DerivedField(numNodes int, bricks map[string]Brick) DerivedField
- func (idt *IdealMixtureTerm) Eval(i int, bricks map[string]Brick) complex128
- func (idt *IdealMixtureTerm) GetEnergy(bricks map[string]Brick, nodes int) float64
- func (idt *IdealMixtureTerm) OnStepFinished(t float64, bricks map[string]Brick)
- type ImplicitEuler
- type LaplacianN
- type MixedTerm
- type ModalFilter
- type Model
- func (m *Model) AddEquation(eq string)
- func (m *Model) AddField(f Field)
- func (m *Model) AddScalar(s Scalar)
- func (m *Model) AddSource(eqNo int, s Source)
- func (m *Model) AllFieldNames() []string
- func (m *Model) AllVariableNames() []string
- func (m *Model) EqNumber(fieldName string) int
- func (m *Model) GetDenum(fieldNo int, freq Frequency, t float64) []complex128
- func (m *Model) GetRHS(fieldNo int, freq Frequency, t float64) []complex128
- func (m *Model) Init()
- func (m *Model) IsBrickName(name string) bool
- func (m *Model) IsExplicitTerm(desc string) bool
- func (m *Model) IsFieldName(name string) bool
- func (m *Model) IsImplicitTerm(desc string) bool
- func (m *Model) IsMixedTerm(desc string) bool
- func (m *Model) IsUserDefinedTerm(desc string) bool
- func (m *Model) NumNodes() int
- func (m *Model) RegisterDerivedField(d DerivedField)
- func (m *Model) RegisterExplicitTerm(name string, t PureTerm, dFields []DerivedField)
- func (m *Model) RegisterFunction(name string, F GenericFunction)
- func (m *Model) RegisterImplicitTerm(name string, t PureTerm, dFields []DerivedField)
- func (m *Model) RegisterMixedTerm(name string, t MixedTerm, dFields []DerivedField)
- func (m *Model) RegisterRHSModifier(eqNumber int, modifier func(data []complex128))
- func (m *Model) Summarize()
- func (m *Model) SyncDerivedFields()
- func (m *Model) UpdateDerivedFields(eq string)
- type Monitor
- type NegativeValuePenalty
- type PairCorrlationTerm
- type PointMonitor
- type PureTerm
- type RHS
- type RHSModifier
- type RK4
- type RealAmplitudeIterator
- type SDD
- func (sdd *SDD) DimerLength(t float64) float64
- func (sdd *SDD) FieldNorm(fields []Field) float64
- func (sdd *SDD) GetTime() float64
- func (sdd *SDD) Init(init []Field, final []Field)
- func (sdd *SDD) RequiredDimerLengthTime(l float64) float64
- func (sdd *SDD) SaveOrientation(fname string)
- func (sdd *SDD) SetFilter(filter ModalFilter)
- func (sdd *SDD) SetInitialOrientation(orient []float64)
- func (sdd *SDD) ShiftFieldsAlongDimer(fields []Field, scale float64)
- func (sdd *SDD) Step(m *Model)
- type SDDMonitor
- type SDDTimeConstants
- type Scalar
- type Solver
- type SolverCB
- type Source
- type Sources
- type SpectralViscosity
- type SquaredGradient
- type SubStringDelimiter
- type TensorialHessian
- type Term
- type TimeDepSource
- type TimeStepper
- type Uint8IO
- type UniqueFreqIterator
- type Vandeven
- type VolumeConservingLP
- type WeightedLaplacian
- type WhiteNoise
- type XDMF
- type XDMFAttribute
- type XDMFDataItem
- type XDMFDomain
- type XDMFGeometry
- type XDMFGrid
- type XDMFTime
- type XDMFTopology
Examples ¶
Constants ¶
This section is empty.
Variables ¶
This section is empty.
Functions ¶
func ApplyModalFilter ¶
func ApplyModalFilter(filter ModalFilter, freq Frequency, data []complex128)
ApplyModalFilter applies the filter f in-place to data
func DefaultNonLinSolver ¶
func DefaultNonLinSolver() nonlin.NewtonKrylov
DefaultNonLinSolver returns the default non-linear solver used in Implicit Euler. Internally, the GMRES method is used. With the default settings, this method can consume a lot of memory depending on the problem. If the program uses too much memory, try to use a restarted version of GMRES (e.g. InnerMethod: &linsolve.GMRES{Restart: 50}) See GMRES description at https://godoc.org/github.com/gonum/exp/linsolve for further details
func GetFieldName ¶
GetFieldName returns the field name of a term.
func GetNonLinearFieldExpressions ¶
GetNonLinearFieldExpressions returns constructions that must be FFT separateley. Note that if the expression is bilinear in field, field is omitted from the returned expression
func IndicatorDeriv ¶
IndicatorDeriv is the derivative of the indicator function
func LoadFloat64 ¶
LoadFloat64 loads an array of float64 encoded as binary data it is assumed that the it is stored with BigEndian
func RealPartAsUint8 ¶
func RealPartAsUint8(data []complex128, min float64, max float64) []uint8
RealPartAsUint8 return the real part of the field as uint8. The data is scaled such that min --> 0 and max --> 255
func SaveCsv ¶
SaveCsv stores data to a csv file. The format is X, Y, Z, field1, field2, field3 DomainSize gives the shape of the domain
func SaveFloat64 ¶
SaveFloat64 writes a float sice to a binary file. BigEndian is used. Files stored with this function can be read using LoadFloat64
func SortFactors ¶
SortFactors sorts the factors in an expression in alphabetical order The passed string is splitted on * and then the resulting slice is sorted. Finally, it is join together. Example: solute*temperature*conc is converted to conc*solute*temperature
Types ¶
type Advection ¶
Advection is a type that is used to build. This represents a term of the form v dot grad field, where v is a velocity vector and field is the name of a field Since terms in GOPF are assumed to enter on the right hand side of equations of the form dy/dt = ..., this term return -v dot grad field. Field is the name of the field and VelocifyFields is a slice with the names of the velocity components. The number of velocity components has to be exactly equal to the number of dimensions (e.g. if it is 2D model then the length of this slice should be two)
func (*Advection) AllFieldsExist ¶
AllFieldsExist returns true if all the fields are registered in the model
func (*Advection) OnStepFinished ¶
OnStepFinished does nothing as we don't need any updates in between steps
func (*Advection) PrepareModel ¶
func (ad *Advection) PrepareModel(N int, m *Model, FT FourierTransform)
PrepareModel adds the nessecary fields to the model in order to be able to use the Advection term
type BoundTransform ¶
BoundTransform is a type that maps variables on the real-line (-inf, inf) into the interval (Min, Max). It does so by applying the following transformation: y = Min + (Max-Min)*(arctan(x)/pi + 1) x is a number on the real line.
func (*BoundTransform) Backward ¶
func (bt *BoundTransform) Backward(y float64) float64
Backward maps the coordinate y into x
func (*BoundTransform) Forward ¶
func (bt *BoundTransform) Forward(x float64) float64
Forward returns the mapped coordinate
type Brick ¶
type Brick interface {
Get(i int) complex128
}
Brick is a generic interface to terms in a PDE
type BrickPlaceholder ¶
type BrickPlaceholder struct{}
BrickPlaceholder is struct that is used as a brick, in case no brick is specified
type ChargeTransport ¶
type ChargeTransport struct { // Conductivity return the conductivity tensor at node i. If 3D // the order should be s_xx, s_yy, s_zz, s_xz, s_yz, s_xy and if 2D // it should be s_xx, s_yy, s_xy. The vacuum permittivity should be embeded // in the conductivity. Thus, if the "normal" conductiviy in Ohm's law is // labeled sigma, the current density is given by j = sigma*E, this function // should return sigma/eps_0, where eps_0 is the vacumm permittivity Conductivity func(i int) []float64 // ExternalField represents the external electric field ExternalField []float64 // Field is the name of the field that corresponds to the charge density Field string // FT is a fourier transformer FT FourierTransform }
ChargeTransport implements a term corresponding to the rate of change of local charge density in an electric field. The implementation follows closely the approach in
Jin, Y.M., 2013. Phase field modeling of current density distribution and effective electrical conductivity in complex microstructures. Applied Physics Letters, 103(2), p.021906. https://doi.org/10.1063/1.4813392
The evolution of the charge density is described by ¶
dp/dt = -div j, where p is the charge density and j is the current density given by j_i = sigma_ij E_j, where sigma is the conductivity tensor and E is the local electric field. The local electriv field is calculated by solving poisson equation LAP phi = -p/eps_0, eps_0 is the vacuum permittivity The net result of this is that the local electric field is given by
** i * d^3 k p(k)k
E = E_ext - ------- * --------- -------- e^(ik*r)
eps_0 * (2pi)^3 k^2 **
where p(k) is the fourier transform of the density, k is the reciprocal wave vector and r is the real space position vector. E_ext is the external electric field
func (*ChargeTransport) Construct ¶
func (ct *ChargeTransport) Construct(bricks map[string]Brick) Term
Construct returns the fourier transformed value of
func (*ChargeTransport) Current ¶
func (ct *ChargeTransport) Current(density Brick, N int, realspace bool) [][]float64
Current returns the current in real space. Density should give the fourier transformed or real space density. If it represents the realspace density, the realspace flag should be set to true (in which case a fourier transform is performed internally) otherwise it should be set to false charge density
func (*ChargeTransport) OnStepFinished ¶
func (ct *ChargeTransport) OnStepFinished(t float64, b map[string]Brick)
OnStepFinished does nothing for this term
type ConservativeNoise ¶
ConservativeNoise adds noise in a such a way that the field it is added to is conserved
func NewConservativeNoise ¶
func NewConservativeNoise(strength float64, dim int) ConservativeNoise
NewConservativeNoise returns an instance of ConservativeNoise with a correctly initialized UniquePrefix which is used to identify derived fields associated with the conservative noise
func (*ConservativeNoise) Construct ¶
func (cn *ConservativeNoise) Construct(bricks map[string]Brick) Term
Construct builds the right hand side term
func (*ConservativeNoise) CurrentFieldsAreRegistered ¶
func (cn *ConservativeNoise) CurrentFieldsAreRegistered(bricks map[string]Brick) bool
CurrentFieldsAreRegistered checks that all the current fields are available among the bricks
func (*ConservativeNoise) GetCurrentName ¶
func (cn *ConservativeNoise) GetCurrentName(comp int) string
GetCurrentName returns the name of the current field
func (*ConservativeNoise) OnStepFinished ¶
func (cn *ConservativeNoise) OnStepFinished(t float64, bricks map[string]Brick)
OnStepFinished does nothing
func (*ConservativeNoise) RequiredDerivedFields ¶
func (cn *ConservativeNoise) RequiredDerivedFields(numNodes int) []DerivedField
RequiredDerivedFields returns the a set of derived fields that is nessecary to register in order to use conservative noise. The number of nodes in the grid is specified via the numNodes argument
type CsvData ¶
type CsvData struct { Name string Data pfutil.ImmutableSlice }
CsvData is a type that can be used to store data to a CSV file Name will be placed as header
type CsvIO ¶
CsvIO writes data to text csv text files
func (*CsvIO) SaveFields ¶
SaveFields stores the results in CSV files. The format X, Y, Z, field1, field2, field3 etc.
type DerivedField ¶
type DerivedField struct { Data []complex128 Name string Calc DerivedFieldCalc }
DerivedField is a type that is derived from Fields via multiplications and power operations
func (DerivedField) Get ¶
func (d DerivedField) Get(i int) complex128
Get returns the value at position i
func (*DerivedField) Update ¶
func (d *DerivedField) Update()
Update recalculates the derived fields and places the result in Data
type DerivedFieldCalc ¶
type DerivedFieldCalc func(data []complex128)
DerivedFieldCalc is a function that calculates the derived field
func DerivedFieldCalcFromDesc ¶
func DerivedFieldCalcFromDesc(desc string, fields []Field) DerivedFieldCalc
DerivedFieldCalcFromDesc returns a derived field calculator based on its description (e.g. conc^2*eta, if conc and eta are two fields)
type DiffOp ¶
type DiffOp func(freq Frequency, ft []complex128) []complex128
DiffOp is a differential operator. It takes an object that has the frequency method implemented and applies the operator in-place to the passed array
type DisplacementGetter ¶
type DisplacementGetter func(force [][]complex128, freq elasticity.Frequency, matProp elasticity.Rank4Tensor) [][]complex128
DisplacementGetter is a functon that type that returns the displacement, when given fourier transformed body forces (force), a corresponding frequency getter and a set elasticity tensor
type DivGrad ¶
type DivGrad struct { Field string F GenericFunction }
DivGrad is a type used to represent the term Div F(c)Grad <field>
func (*DivGrad) OnStepFinished ¶
OnStepFinished does nothing as we don't need any updates in between steps
func (*DivGrad) PrepareModel ¶
func (dg *DivGrad) PrepareModel(N int, m *Model, FT FourierTransform)
PrepareModel impose nessecary changes in the model in order to use the DivGrad term. N is the number of grid points in the simulation domain, dim is the dimension of the simulation domain. FT is a fourier transformer required for gradient evaluations
type Euler ¶
type Euler struct { Dt float64 FT FourierTransform Filter ModalFilter CurrentStep int }
Euler performs semi-implicit euler method
func (*Euler) SetFilter ¶
func (eu *Euler) SetFilter(filter ModalFilter)
SetFilter set a new modal filter
type ExplicitPairCorrelationTerm ¶
type ExplicitPairCorrelationTerm struct {
PairCorrlationTerm
}
ExplicitPairCorrelationTerm implement the pair correlation function, but the construct method returns the expression corresponding to an explicit treatment of the term in the PDE. The only difference between ExplicitPairCorrelationTerm and PairCorrelationTerm is that the Construct method returns the explicit and implicit variant, respectively.
type Field ¶
type Field struct { Data []complex128 Name string }
Field is a type that is used to represent a field in the context of phase field models
func NewField ¶
func NewField(name string, N int, data []complex128) Field
NewField initializes a new field
type FieldDB ¶
type FieldDB struct { DB *sql.DB // DomainSize of the simulation domain. This is needed in order to populate the // position table DomainSize []int // contains filtered or unexported fields }
FieldDB is a type used for storing field data in a SQl database. A field database consists of the following tables
- positions (id int, x int, y int, z int) Describes the positions in 3D space of all nodes. The first column is the node number in the simulation cell and the remaining columns represent the x, y and z index, respectively. If the simulation domain is 2D, the z column is always zero.
- fields (id int, name text, value real, positionId int, timestep int, simID int) Describes the value of the fields at a given timestep and point - id: Row identifier that is auto-incremented when new rows are added - name: Name of the field - value: Value of the field - positionId: ID pointing to the positions table. The position in 3D space of the node represented by this row, is given by the corresponding row in the positions table - timestep: Timestep of the record - simID: Unique ID identifying all entries in the database written by the current simulation
- simAttributes (key text, value float, simID int) Additional attributes that belongs to the simulations. Input arguments can for example be stored here. - key: Name of the attribute - value: Value of the attribute - simID: Unique ID which is common to all entries written by the current simulation
- comments (simID int, value text) Describes comments about the simulation. - simID: Unique ID which is common to all entries written by the current simulation - value: A text string describing the simulation
- simIDs (simID int, creationTime text) List of all simulation IDs in the database - simID: Unique ID - creationTime: Timestamp for when the simulation ID was created
- simTextAttributes (key TEXT, value TEXT, simID int) Same as simAttributes, apart from the value field is a string.
- timeseries (key TEXT, value float, timestep int, simID int) Describes time varying data typically derived from the fields in the calculations. Some examples: peak concentration in a diffusion calculation, peak stress in an elasiticy calculation, average domain size in a spinoidal calculation etc. - key: Name of the item - value: Value of the data point - timestep: Timestep - simID: Unique ID which is common to all entries written by the current simulation
func (*FieldDB) Load ¶
Load loads all the fields from a database and return a list of Field simID is the ID of the simulation that the field should be loaded from, and timestep is the timestep from which the fields should be initialized.
func (*FieldDB) LoadLast ¶
LoadLast loads the fields from the latest timestep available for the passed simulation ID
func (*FieldDB) SaveFields ¶
SaveFields stores all the field to the database. This function satisfies the SolverCB type, and can thus be attached as a callback to a solver
func (*FieldDB) SetTextAttr ¶
SetTextAttr sets text attributes associated with the current simulation
type Float64IO ¶
type Float64IO struct {
Prefix string
}
Float64IO stores the fields as raw binary files using BigEndian. The datatype is float64
func NewFloat64IO ¶
NewFloat64IO returns a new Float64IO. All files are prepended ay prefix
func (*Float64IO) SaveFields ¶
SaveFields stores all fields as raw binary files. It can be passed as a callback to the solver
type FourierTransform ¶
type FourierTransform interface { FFT(data []complex128) []complex128 IFFT(data []complex128) []complex128 Freq(i int) []float64 }
FourierTransform is a type used to represent fourier transforms
type GenericFunction ¶
type GenericFunction func(i int, bricks map[string]Brick) complex128
GenericFunction is a generic function that may depend on any of the fields
type GradientCalculator ¶
type GradientCalculator struct { FT FourierTransform Comp int KeepNyquist bool }
GradientCalculator calculates the gradient of a field
func (*GradientCalculator) Calculate ¶
func (g *GradientCalculator) Calculate(indata []complex128, data []complex128)
Calculate calculates the gradient of the data passed data contain the field in real-space
func (*GradientCalculator) ToDerivedField ¶
func (g *GradientCalculator) ToDerivedField(name string, N int, brick Brick) DerivedField
ToDerivedField constructs a derived field from the gradient calculator. N is the number of grid points and brick is the brick that should be differentiated
type HomogeneousModulusLinElast ¶
type HomogeneousModulusLinElast struct { FieldName string Field []float64 Misfit *mat.Dense EffForce elasticity.EffectiveForce MatProp elasticity.Rank4 Disps DisplacementGetter FT FourierTransform Dim int N int }
HomogeneousModulusLinElast is a type that is used in phase field models where the elastic constants are homogeneous throughout the domain. In other words it represents the energy term E = (1/2)*C_{jikl}(eps_{ij} - eps^*_{ij}H(x))(eps_{kl} - eps^*_{kl}H(x)) where C is the elastic tensor, eps^* is the misfit strain and x is a field that is one if we are inside the domain where misfit strains exists and zero elswhere. The specific name of the name of the field x is stored in FieldName. When this term is used in an equation it is assumed that the right hand side consists of dx/dt = -dE/dx + ..., i.e. the time evolution of the field should be such that it minimizes the strain energy. The strain eps_{ij} is determined enforcing mechanical equillibrium
func NewHomogeneousModolus ¶
func NewHomogeneousModolus(fieldName string, domainSize []int, matProp elasticity.Rank4, misfit *mat.Dense) *HomogeneousModulusLinElast
NewHomogeneousModolus initializes a new instance of the linear elasticity model
func (*HomogeneousModulusLinElast) Construct ¶
func (h *HomogeneousModulusLinElast) Construct(bricks map[string]Brick) Term
Construct returns the function needed to build the term on the right hand side
func (*HomogeneousModulusLinElast) Force ¶
func (h *HomogeneousModulusLinElast) Force(indicator []complex128) [][]complex128
Force returns the effective force
func (*HomogeneousModulusLinElast) Freq ¶
func (h *HomogeneousModulusLinElast) Freq(i int) []float64
Freq wraps the passed frequency method such that the length of the returned frequency is always 3
func (*HomogeneousModulusLinElast) OnStepFinished ¶
func (h *HomogeneousModulusLinElast) OnStepFinished(t float64, bricks map[string]Brick)
OnStepFinished update the real space version of the field
type IdealMixtureTerm ¶
type IdealMixtureTerm struct { IdealMix pfc.IdealMix Field string Prefactor float64 Laplacian bool }
IdealMixtureTerm implements the ideal mixture model used in the paper by Greenwood et al. To use this term in a model, register the function Eval as a function in the model. Prefactor is a constant factor that is multiplied with the energy. When used together with a model, this term should be registered as a mixed term. Laplacian indicates whether the laplacian should be applied to the term. If true then the laplacian will be applied internally
func (*IdealMixtureTerm) ConstructLinear ¶
func (idt *IdealMixtureTerm) ConstructLinear(bricks map[string]Brick) Term
ConstructLinear builds the linear part
func (*IdealMixtureTerm) ConstructNonLinear ¶
func (idt *IdealMixtureTerm) ConstructNonLinear(bricks map[string]Brick) Term
ConstructNonLinear builds the non linear part
func (*IdealMixtureTerm) DerivedField ¶
func (idt *IdealMixtureTerm) DerivedField(numNodes int, bricks map[string]Brick) DerivedField
DerivedField returns the required derived field that is nessecary in order to use this model
func (*IdealMixtureTerm) Eval ¶
func (idt *IdealMixtureTerm) Eval(i int, bricks map[string]Brick) complex128
Eval returns the derivative of the underlying ideal mixture term
func (*IdealMixtureTerm) GetEnergy ¶
func (idt *IdealMixtureTerm) GetEnergy(bricks map[string]Brick, nodes int) float64
GetEnergy evaluates the energy contribution from this term. The fields in bricks should be the real space representation. The number of nodes in the simulation cell is specified in the nodes argument
func (*IdealMixtureTerm) OnStepFinished ¶
func (idt *IdealMixtureTerm) OnStepFinished(t float64, bricks map[string]Brick)
OnStepFinished does nothing
type ImplicitEuler ¶
type ImplicitEuler struct { Dt float64 FT *pfutil.FFTWWrapper Filter ModalFilter CurrentStep int // each time step. If not given (or nil), a solver with sensible default values // will be used NonlinSolver *nonlin.NewtonKrylov }
ImplicitEuler performs a full implicit euler method. We have a problem of the form dy/dt = Ly + N(y), where L is a linear operator and N is a non-linear operator. The implicit update scheme is y_{n+1} = y_n + dt*Ly_{n+1} + N(y_{n+1}), this gives rise to a non-linear set of equations that must be solved on each time step. Provided that the non-linear solver is able to find the solution, this scheme is stable for all time steps dt. However, in practice the non-linear solver converges faster it dt is small. Thus, if the solution changes rapidly it might be wise to use a small time step.
func (*ImplicitEuler) GetTime ¶
func (ie *ImplicitEuler) GetTime() float64
GetTime returns the current time
func (*ImplicitEuler) SetFilter ¶
func (ie *ImplicitEuler) SetFilter(f ModalFilter)
SetFilter sets a new filter. Currently this has no effect in this timestepper
func (*ImplicitEuler) Step ¶
func (ie *ImplicitEuler) Step(m *Model)
Step evolves the equation one timestep
type LaplacianN ¶
type LaplacianN struct {
Power int
}
LaplacianN is a type used for the Laplacian operator raised to some power
func (LaplacianN) Eval ¶
func (l LaplacianN) Eval(freq Frequency, ft []complex128) []complex128
Eval implements the fourier transformed laplacian operator. freq is a function that returns the frequency at index i of the passed array. ft is the fourier transformed field
type MixedTerm ¶
type MixedTerm interface { // ConstructLinear builds the function to evaluate the linear part of the // term. The function returned should give the fourier transform of // f({otherfields}). The bricks parameter contains the fourier transform // of all bricks ConstructLinear(bricks map[string]Brick) Term // ConstructNonLinear returns a function that calculates the non-linear // part of the expression. The function being returned should calculate // the fourier transform of f({allfields}). The bricks parameter contains // the fourier transform of all known bricks ConstructNonLinear(bricks map[string]Brick) Term // OnStepFinished is called after each step is finished. If a term needs // be updated based on how the fields evolves, the update should happen // inside this method OnStepFinished(t float64, bricks map[string]Brick) }
MixedTerm is a type that can be used to represents terms that have both a linear part and a non-linear part. Mixed terms are on the form f({otherfields})*field + g({allfields}), where f is a function (or operator) that depends on all the other fields in the model and g is a function/operator that depends on all the fields in the model
type ModalFilter ¶
ModalFilter is a generic interface for modal filters
type Model ¶
type Model struct { Fields []Field DerivedFields []DerivedField Bricks map[string]Brick ImplicitTerms map[string]PureTerm ExplicitTerms map[string]PureTerm MixedTerms map[string]MixedTerm Equations []string RHS []RHS AllSources []Sources RHSModifiers []eqModifier }
Model is a type used to represent a general set of equations
Example ¶
N := 16 f1 := NewField("conc", N*N, nil) for i := range f1.Data { if i > 5 { f1.Data[i] = complex(0.1, 0.0) } } model := NewModel() model.AddField(f1) model.AddEquation("dconc/dt = conc^3 - conc + LAP conc") dt := 0.01 solver := NewSolver(&model, []int{N, N}, dt) return &model, solver
Output:
func (*Model) AddEquation ¶
AddEquation adds equations to the model
func (*Model) AllFieldNames ¶
AllFieldNames returns all field names (including derived fields)
func (*Model) AllVariableNames ¶
AllVariableNames return all known variable names
func (*Model) EqNumber ¶
EqNumber returns the equation number corresponding to the passed field name
func (*Model) GetDenum ¶
func (m *Model) GetDenum(fieldNo int, freq Frequency, t float64) []complex128
GetDenum evaluates the denuminator
func (*Model) GetRHS ¶
func (m *Model) GetRHS(fieldNo int, freq Frequency, t float64) []complex128
GetRHS evaluates the right hand side of one of the equations
func (*Model) IsBrickName ¶
IsBrickName returns true if a brick with the passed name exists
func (*Model) IsExplicitTerm ¶
IsExplicitTerm returns true if the passed string is a non-linear term
func (*Model) IsFieldName ¶
IsFieldName checks if the passed name is a field name
func (*Model) IsImplicitTerm ¶
IsImplicitTerm checks if the given term is a linear term
func (*Model) IsMixedTerm ¶
IsMixedTerm returns true if the passed string is a mixed term
func (*Model) IsUserDefinedTerm ¶
IsUserDefinedTerm returns true if desc matches either one of the linear terms, non-linear terms or mixed terms
func (*Model) NumNodes ¶
NumNodes returns the number of nodes in the simulation cell. It panics if no fields are added
func (*Model) RegisterDerivedField ¶
func (m *Model) RegisterDerivedField(d DerivedField)
RegisterDerivedField registers a new derived field
func (*Model) RegisterExplicitTerm ¶
func (m *Model) RegisterExplicitTerm(name string, t PureTerm, dFields []DerivedField)
RegisterExplicitTerm defines a new term. To add the term to an equation add the name as one of the terms.
Example: If there is a user defined term called LINEAR_ELASTICITY, and we have a PDE where the term is present on the right hand side, the equation would look like dx/dt = LINEAR_ELASTICITY where x is the name of the field. The additional derived fields (which are fields that are contructed from the original fields) is specified via dFields
func (*Model) RegisterFunction ¶
func (m *Model) RegisterFunction(name string, F GenericFunction)
RegisterFunction registers a function that may be used in the equations
func (*Model) RegisterImplicitTerm ¶
func (m *Model) RegisterImplicitTerm(name string, t PureTerm, dFields []DerivedField)
RegisterImplicitTerm can be used to register terms if the form f({otherfields])*field
func (*Model) RegisterMixedTerm ¶
func (m *Model) RegisterMixedTerm(name string, t MixedTerm, dFields []DerivedField)
RegisterMixedTerm is used to register terms that contains a linear part and a non-linear part. The linear part will be treated implicitly during time evolution, while the non-linear part is treated explicitly
func (*Model) RegisterRHSModifier ¶
func (m *Model) RegisterRHSModifier(eqNumber int, modifier func(data []complex128))
RegisterRHSModifier adds a modifier that is applied to the fourier transformed right hand side of the equation. The modifier will be applied to the right hand side of the equation passed. If several modifiers are added to the same equation, the one that is added first will be applied first
Simple example case: If we want to convert the units of the right from let's milli joules to joules, we should multiply the current right hand side by a factor 1e-3. We therefore supply the following modifier
func (data []complex128) { for i := range data { data[i] *= 1e-3 } }
The frequency corresponding to each node can be obtained from the frequency method of the fourier transformer used (see for example FFTWWrapper)
func (*Model) SyncDerivedFields ¶
func (m *Model) SyncDerivedFields()
SyncDerivedFields updates all derived fields
func (*Model) UpdateDerivedFields ¶
UpdateDerivedFields update fields that needs to be handle with FFT (required for non-linear equations)
type NegativeValuePenalty ¶
NegativeValuePenalty is a an explicit term that adds a penalty for negative values. The function represents f(x) = Prefactor*(|x|^Exponent - x^Exponent)
func NewDefaultNegativeValuePenalty ¶
func NewDefaultNegativeValuePenalty(field string) NegativeValuePenalty
NewDefaultNegativeValuePenalty returns a new constraint type using the parameters from
Chan, P.Y., Goldenfeld, N. and Dantzig, J., 2009. Molecular dynamics on diffusive time scales from the phase-field-crystal equation. Physical Review E, 79(3), p.035701.
func (*NegativeValuePenalty) Evaluate ¶
func (nvp *NegativeValuePenalty) Evaluate(i int, bricks map[string]Brick) complex128
Evaluate can be registered as a function in any model
func (*NegativeValuePenalty) Penalty ¶
func (nvp *NegativeValuePenalty) Penalty(x float64) float64
Penalty returns the derivative of the underlying function
type PairCorrlationTerm ¶
type PairCorrlationTerm struct { PairCorrFunc pfc.ReciprocalSpacePairCorrelation Field string Prefactor float64 Laplacian bool }
PairCorrlationTerm implements the functional deriviative with respect to Q of the functional
** ** A * * g[Q] = - --- * dr Q(r) * dr' C(|r - r'|)Q(r') 2 * * ** **
where C(|r-r'|) is a pair correlation function. PairCorrFunc is the fourier transform of the pair correlation function and Field is the name of the field (e.g. name of Q in the equation above). Prefactor is a constant factor that is multiplied with the integral (A in the equation above). The attribute Laplacian determines whether the the laplacian should be applied to the functional derivative or not. If true, then this term represents nabla^2 dg/dQ, otherwise it simply represents dg/dQ. When used in a model this term should be registered as an implicit term
func (*PairCorrlationTerm) Construct ¶
func (pct *PairCorrlationTerm) Construct(bricks map[string]Brick) Term
Construct builds the rhs required to represent the term
func (*PairCorrlationTerm) GetEnergy ¶
func (pct *PairCorrlationTerm) GetEnergy(bricks map[string]Brick, ft FourierTransform, domainSize []int) float64
GetEnergy evaluates the energy contribution from this term. The fields in bricks should be the real space representation. The number of nodes in the simulation cell is specified in the nodes argument
func (*PairCorrlationTerm) OnStepFinished ¶
func (pct *PairCorrlationTerm) OnStepFinished(t float64, bricks map[string]Brick)
OnStepFinished is simply included to satisfy the UserDefinedTerm interface
type PointMonitor ¶
A PointMonitor is a monitor that monitors the solution at a given pixel/voxel
func NewPointMonitor ¶
func NewPointMonitor(site int, field string) PointMonitor
NewPointMonitor returns a new instance of PointMonitor
func (*PointMonitor) Add ¶
func (p *PointMonitor) Add(bricks map[string]Brick)
Add adds a new value to the monitor
type PureTerm ¶
type PureTerm interface { // Construct creates the righ hand side function for the term. // The returned function should should give the fourier transformed // quantity. If this method is a linear term, it should return the // fourier transform of f({otherfields}) (excluding the multiplication // with the field in question). If it is a non-linear term, it should // return the fourier transform of f({allfields}). See the documentation // of the interface for a detailed definition of f({otherfields}) and // f({allfields}). When bricks is passed to this method, all fields // have already been fourier transformed. Construct(bricks map[string]Brick) Term // OnStepFinished is gets called after each time step. It can be used // to do nessecary updates after the fields have been updated. If no such // updates are nessecary, just implement an empty function OnStepFinished(t float64, bricks map[string]Brick) }
PureTerm is an interface that can be used to add special terms to a PDE that cannot easily be defined via a string representation. This term type is either linear or non-linear. Since, a linear term can be treated as a non-linear term, it is beneficial to treat the linear and non-linear part differently. Terms that exhibit both linear and non-linear parts should be treated via the MixedTerm interface. In this context a linear term is a term on the form f({otherfields})*field, where f is a function (or operator) that is independent of the field in question, but can depend on the other fields in the model. Terms of this form can easily be treated implicitly.
Non-linear terms, are terms on the form f({all fields}), where f is a function (or operator) that depends on all the fields. Terms on this form will be treated explicitly when evolving the fields
type RHSModifier ¶
type RHSModifier func(data []complex128)
RHSModifier is a function type used to modify the right hand side prior to adding it
type RK4 ¶
type RK4 struct { Dt float64 FT FourierTransform Filter ModalFilter CurrentStep int }
RK4 implements the fourth order Runge-Kutta scheme. Dt is the timestep FT is a fourier transform object used to translate back-and fourth between fourier domain.
func (*RK4) PrepareNextCorrection ¶
func (rk *RK4) PrepareNextCorrection(initial []Field, final []Field, kFactor []Field, m *Model, factor float64)
PrepareNextCorrection updates the final fields and rests the fields of the model to the original
func (*RK4) SetFilter ¶
func (rk *RK4) SetFilter(filter ModalFilter)
SetFilter sets a new modal filter
func (*RK4) Step ¶
Step performs one RK4 time step. If the equation is given by dy/dt = A*y + N(y) where N(y) is some non-linear function, the update consists of the following steps y_{n+1} = (y_n + dt*(k1 + 2*k2 + 2*k3 + k4)/6)/(1 - dt*A), where the k coefficients are given by
k1 = N(y_n) k2 = N( (y_n + dt*k1/2)/(1 - 0.5*dt*A) ) k3 = N( (y_n + dt*k2/2)/(1 - 0.5*dt*A) ) k4 = N( (y_n + dt*k3)/(1 - dt*A) )
This leads to a scheme that is first order accurate in dt, but with much better stability properties than the Euler scheme
type RealAmplitudeIterator ¶
type RealAmplitudeIterator struct { Freq Frequency End int // contains filtered or unexported fields }
RealAmplitudeIterator iterates over all frequencies that has a real fourier amplitude when the input signal has a real amplitude
func (*RealAmplitudeIterator) Next ¶
func (rai *RealAmplitudeIterator) Next() int
Next returns the next index. It returns -1 when the iterator is exhausted
type SDD ¶
type SDD struct { TimeConstants SDDTimeConstants // Alpha contains a weight for how to evaluate the force at the center of the dimer // the "end points" are given by // x_1 = x - alpha*l*v/2 and x_2 = x + alpha*l*v/2 // where l is the length of the dimer and v is a unit orientation vector // The force exerted at the center is given by // F(x) = alpha*F(x_1) + (1-alpha)*F(x_2). Default when the SDD is constructed // from NewSDD is alpha = 0.5. Alpha float64 // Dt is the timestep used to evlolve the equation Dt float64 // CurrentStep contains the current starting iteration CurrentStep int // Minimum dimer length (default 0). When the length of dimer reaches this value // it will stop to shrink and maintain the minimum length MinDimerLength float64 // Monitor holds data on the status of the SDD stepper Monitor SDDMonitor // InitDimerLength is the length of the dimer at the first timestep. // Default is to use the length of the orientation vector passed to // Init or SetOrientationVector InitDimerLength float64 // contains filtered or unexported fields }
SDD implements Shrinking-Dimer-Dynamics
func (*SDD) DimerLength ¶
DimerLength returns the the length of the timer at the given time
func (*SDD) RequiredDimerLengthTime ¶
RequiredDimerLengthTime returns the the time needed to reach the passed length
func (*SDD) SaveOrientation ¶
SaveOrientation stores the current orientation vector to a csv file
func (*SDD) SetFilter ¶
func (sdd *SDD) SetFilter(filter ModalFilter)
SetFilter is implemented to satisfy the TimeStepper interface. However, a call to this method will panic because model filters are currently not supported by this stepper.
func (*SDD) SetInitialOrientation ¶
SetInitialOrientation sets the starting orientation vector the initial dimer length is set to the length of the passed orientation vector Subsequently, the passed orientation vector will be normalized to unit length
func (*SDD) ShiftFieldsAlongDimer ¶
ShiftFieldsAlongDimer shifts the field along the dimer. The length is given by scale*v, where v is the dimer orientation vector
type SDDMonitor ¶
type SDDMonitor struct { // MaxForce is the maximum fourier transformed force on the center MaxForce float64 // ForcePowerSpectrum contains the power spectrum of the force // sum_{m=0}^{N-1} F(k_m)/N, where N is the number of nodes. This quantity // is the same as the real space integram // sum_ {m=0}^{N-1} F(k_m)/N = sum_{m=0}^{N-1} f(r_m), where F(k) = FFT(f(r)) ForcePowerSpectrum float64 // MaxTorque is the maximum torque exerted on the dimer MaxTorque float64 // FieldNorm is the L2 norm of the fields which is given by // sum_{f=0}^{M_f-1} sum_{n=0}^{N-1} | Field[f].Data[n] |^2, // where M_f is the number of fields and N is the number of nodes FieldNorm float64 // FieldNormChange contains the difference in L2 norm of between two // sucessive time steps FieldNormChange float64 // LogFile is a writable file where the values of the quantities // above will be logged at each timestep. If not given (or nil), // no information will be logged LogFile *os.File }
SDDMonitor is a type that is used to monitor the progress of the SDD stepper
func (*SDDMonitor) Log ¶
func (sddm *SDDMonitor) Log(s *Solver, epoch int)
Log prints the solver to a csv based format. This method can be attached as a callback to the solver
type SDDTimeConstants ¶
type SDDTimeConstants struct { // Orientation is the time constant used to evolve the orientation vector equation Orientation float64 // DimerLength is the time constant used to evolve the length of the dimer // (e.g. l(t) = l_0 exp(-t/tau)) DimerLength float64 }
SDDTimeConstants is a struct for storing time-constants for the additional equations used in SDD
type Scalar ¶
type Scalar struct { Value complex128 Name string }
Scalar represents a scalar value
func NewScalar ¶
func NewScalar(name string, value complex128) Scalar
NewScalar returns a new scalar value
type Solver ¶
type Solver struct { Model *Model Dt float64 FT FourierTransform Stepper TimeStepper Callbacks []SolverCB Monitors []Monitor StartEpoch int }
Solver is a type used to solve phase field equations
func (*Solver) AddCallback ¶
AddCallback appends a new callback function to the solver
func (*Solver) AddMonitor ¶
AddMonitor adds a new monitor to the solver
func (*Solver) JSONifyMonitors ¶
JSONifyMonitors return a JSON representation of all the monitors
func (*Solver) SetStepper ¶
SetStepper updates the stepper method based on a string. name has to be one of ["euler", "rk4"]
type SolverCB ¶
SolverCB is function type that can be added to the solver it is executed after each iteration
type Source ¶
type Source struct { Pos []float64 // contains filtered or unexported fields }
Source is structure used to add source terms to an equation
func NewSource ¶
func NewSource(pos []float64, f TimeDepSource) Source
NewSource returns a new source instance
type SpectralViscosity ¶
SpectralViscosity implements the spectral viscosity method proposed in Tadmor, E., 1989. Convergence of spectral methods for nonlinear conservation laws. SIAM Journal on Numerical Analysis, 26(1), pp.30-44. In the fourier transformed domain, the following the spectral viscosity term can be written as -eps*Q(k)*k^Power <field>, where field is an arbitrary field name. The function Q(k) is an interpolating function defined by. For simplicity m = DissipationThreshold is introduced and x = 3k/2m - 1/2
** * 0, if k < m/3
Q(k) = ** 2x^2 - 3x^3, if k >= m/3 and k <= m
- 1, if k > m **
In the paper by E. Tadmor it is found that Eps*m = 0.25 yields good results
func (*SpectralViscosity) Construct ¶
func (sv *SpectralViscosity) Construct(bricks map[string]Brick) Term
Construct return the denomonator that is required for an implicit treatment of the spectral viscosity term
func (*SpectralViscosity) OnStepFinished ¶
func (sv *SpectralViscosity) OnStepFinished(t float64, bricks map[string]Brick)
OnStepFinished empty function implemented to satisfy interface
type SquaredGradient ¶
type SquaredGradient struct { Field string FT FourierTransform Factor float64 }
SquaredGradient an be used for terms that has a gradient raised to second power (e.g (grad u)^2, where u is a field). The power has to be an even number Factor is multiplied with the gradient. Thus, the full term reads Factor*(grad Field)^2
func NewSquareGradient ¶
func NewSquareGradient(field string, domainSize []int) SquaredGradient
NewSquareGradient returns a new instance of square gradient
func (*SquaredGradient) Construct ¶
func (s *SquaredGradient) Construct(bricks map[string]Brick) Term
Construct return a function the calculates the squared of the gradient of
func (*SquaredGradient) OnStepFinished ¶
func (s *SquaredGradient) OnStepFinished(t float64, bricks map[string]Brick)
OnStepFinished does not need to perform any work
type SubStringDelimiter ¶
SubStringDelimiter is a type that represents a substring as well as the delimiter preceeding it in the string it was extracted from
func SplitOnMany ¶
func SplitOnMany(value string, delimiters []string) []SubStringDelimiter
SplitOnMany splits a string on several delimiters. The resulted split is returned together with the delimiter preceeding it
type TensorialHessian ¶
TensorialHessian is a type used to represent the term K_ij d^2c/dx_idx_j (sum of i and j). The tensial hessian brick can be added to an equation in two ways as it can be treated both explicitly and implicitly (in some cases). If the term can be treated implicitly, is recommended that one specify it such that the equation parser understands that it should be dealt with implicitly. Cases where it can be treated implicitly is if the hessian should be applied to the same field as located on the left hand side
d<field>/dt = K_ij d^2<field>/dx_idx_j
Note that this term can only be applied to the field that is also on the left hand side of the equation. The following will not work
d<fieldA>/dt = K_ij d^2<fieldB>/dx_idx_j
func (*TensorialHessian) Construct ¶
func (th *TensorialHessian) Construct(bricks map[string]Brick) Term
Construct builds the correct RHS term
func (*TensorialHessian) GetCoeff ¶
func (th *TensorialHessian) GetCoeff(i, j int) float64
GetCoeff return the i, j element of the coefficient tensor
func (*TensorialHessian) OnStepFinished ¶
func (th *TensorialHessian) OnStepFinished(t float64, bricks map[string]Brick)
OnStepFinished does nothing as there is no need for updates for this term
type Term ¶
type Term func(freq Frequency, t float64, field []complex128)
Term is generic function type that evaluates the right hand side of a set of ODE used to evolve the phase fields. The function should place the values into the field array
func ConcreteTerm ¶
func ConcreteTerm(termDelim SubStringDelimiter, m *Model) Term
ConcreteTerm returns a function representing the passed term
type TimeDepSource ¶
TimeDepSource is a generic interface for functions that can be used as sources
type TimeStepper ¶
type TimeStepper interface { Step(m *Model) SetFilter(filter ModalFilter) GetTime() float64 }
TimeStepper is a generic interface for a the time stepper types
type Uint8IO ¶
type Uint8IO struct {
Prefix string
}
Uint8IO is a struct used to store fields as uint8
func (*Uint8IO) SaveFields ¶
SaveFields can be passed as a callback to the solver. It stores each field in a raw binary file.
type UniqueFreqIterator ¶
UniqueFreqIterator is an iterator that can be used to iterate over all the unique frequencies of the fourier transform of a real-valued function. This iterator can be used in a for loop as follows for i := iterator.Next(); i != -1; i = iterator.Next()
func (*UniqueFreqIterator) Next ¶
func (ufi *UniqueFreqIterator) Next() int
Next returns the index of unique site. It returns -1 when exhausted
type Vandeven ¶
type Vandeven struct {
Data []float64
}
Vandeven implements the family of filters described in Vandeven, H., 1991. Family of spectral filters for discontinuous problems. Journal of Scientific Computing, 6(2), pp.159-192.
func NewVandeven ¶
NewVandeven constructs a new vanden filter of the passed order
type VolumeConservingLP ¶
type VolumeConservingLP struct { Multiplier float64 Indicator string Field string CurrentIntegral float64 Dt float64 NumNodes int IsFirstUpdate bool }
VolumeConservingLP is a term that can be added to a PDE such that the volume averaged integral of a field is zero. It does so by tuning a time dependent Lagrange multiplier.
Note that due the way the multiplier is updated, the field is only approximately conserved. This is most easily seen when the fields changes rapidly
func NewVolumeConservingLP ¶
func NewVolumeConservingLP(fieldName string, indicator string, dt float64, numNodes int) VolumeConservingLP
NewVolumeConservingLP returns a new instance of the volume conserving LP
func (*VolumeConservingLP) Construct ¶
func (v *VolumeConservingLP) Construct(bricks map[string]Brick) Term
Construct build the function needed to evaluate the term
func (*VolumeConservingLP) OnStepFinished ¶
func (v *VolumeConservingLP) OnStepFinished(t float64, bricks map[string]Brick)
OnStepFinished updates the Lagrange multiplier
type WeightedLaplacian ¶
type WeightedLaplacian struct { // Name of the field to the right of the laplacian (must be registered in the model) Field string // Name of the function in front of the laplacian (must be registered in the model) PreFactor string // FourierTransformer of the correct domain size. This is used to // go back and fourth between the real domain and the fourier domain // internally FT FourierTransform }
WeightedLaplacian is a type used to represent terms of the form F(c) LAP <field>, where F is a function of all the fields.
type WhiteNoise ¶
type WhiteNoise struct {
Strength float64
}
WhiteNoise is a type that can be used to add white noise to a model if Y(x, t) is a noise term, the correlation function is defined by <Y(x, t)Y(x', t')> = 2*Strength*delta(x-x')delta(t-t') (e.g. uncorrelated in time and space)
func (*WhiteNoise) Generate ¶
func (w *WhiteNoise) Generate(i int, bricks map[string]Brick) complex128
Generate returns random gaussian distributed noise
type XDMF ¶
type XDMF struct { XMLName xml.Name `xml:"Xdmf"` Domain XDMFDomain }
XDMF represents the domain attribute in paraview xdmf format
type XDMFAttribute ¶
type XDMFAttribute struct { XMLName xml.Name `xml:"Attribute"` DataItem XDMFDataItem Name string `xml:"Name,attr,omitempty"` Center string `xml:"Center,attr,omitempty"` }
XDMFAttribute is type used to represent the attribute item in paraview xdmf format
type XDMFDataItem ¶
type XDMFDataItem struct { XMLName xml.Name `xml:"DataItem"` Dim int `xml:"Dimension,attr,omitempty"` Format string `xml:"Format,attr,omitempty"` Value string `xml:",chardata"` NumberType string `xml:"NumberType,attr,omitempty"` Dimensions string `xml:"Dimensions,attr,omitempty"` Precision int `xml:"Precision,attr,omitempty"` DataType string `xml:"DataType,attr,omitempty"` Endian string `xml:"Endian,attr,omitempty"` }
XDMFDataItem is a type used to represent the data item in paraview xdmf format
type XDMFDomain ¶
type XDMFDomain struct { XMLName xml.Name `xml:"Domain"` Topology XDMFTopology Geometry XDMFGeometry Grid XDMFGrid }
XDMFDomain represent domain attribute in paraview xdmf format
type XDMFGeometry ¶
type XDMFGeometry struct { XMLName xml.Name `xml:"Geometry"` Name string `xml:"name,attr,omitempty"` DataItems []XDMFDataItem `xml:"DataItem"` Type string `xml:"Type,attr,omitempty"` Reference string `xml:"Reference,attr,omitempty"` }
XDMFGeometry is a structure used to represent eh Geometry item in paraview xdmf format
type XDMFGrid ¶
type XDMFGrid struct { XMLName xml.Name `xml:"Grid"` Name string `xml:"Name,attr,omitempty"` Type string `xml:"GridType,attr,omitempty"` CollectionType string `xml:"CollectionType,attr,omitempty"` Grids []XDMFGrid `xml:"Grid"` Attributes []XDMFAttribute `xml:"Attribute"` Topology XDMFTopology `xml:"Topology,omitempty"` Geometry XDMFGeometry `xml:"Geometry,omitempty"` Time *XDMFTime `xml:"Time"` }
XDMFGrid is used to represent the grid structure in paraview xdmf format
type XDMFTime ¶
type XDMFTime struct { XMLName xml.Name `xml:"Time"` Type string `xml:"TimeType,attr"` DataItem XDMFDataItem }
XDMFTime is used to represent the time attribute in paraview xfdmf format
type XDMFTopology ¶
type XDMFTopology struct { XMLName xml.Name `xml:"Topology"` Name string `xml:"name,attr,omitempty"` Type string `xml:"TopologyType,attr,omitempty"` Dimensions string `xml:"Dimensions,attr,omitempty"` Reference string `xml:"Reference,attr,omitempty"` }
XDMFTopology is a type used to represent the Topology item in paraview xdmf format
Source Files ¶
- advection.go
- boundTransform.go
- chargeTransport.go
- database.go
- diffOp.go
- euler.go
- fileIO.go
- gradientCalculator.go
- homoLinElast.go
- implicitEuler.go
- model.go
- monitor.go
- negative_value_penalty.go
- noise.go
- pairCorrelationTerm.go
- rhsBuilder.go
- rk4.go
- sdd.go
- solver.go
- sourceTerm.go
- spectralViscosity.go
- squareGradientTerm.go
- tensorialHessian.go
- userDefinedTerm.go
- util.go
- vandeven.go
- volumeConserving.go
- xdmf.go