esys.modellib.input Package

Classes

class esys.modellib.input.GaussianProfile(**kwargs)

Bases: esys.escriptcore.modelframe.ParameterSet

Generates a Gaussian profile at center x_c, width width and height A over a domain

Note:Instance variable domain - domain
Note:Instance variable x_c - center of the Gaussian profile (default [0.,0.,0.])
Note:Instance variable A - (in) height of the profile. A maybe a vector. (default 1.)
Note:Instance variable width - (in) width of the profile (default 0.1)
Note:Instance variable r - (in) radius of the circle (default = 0)

In the case that the spatial dimension is two, The third component of x_c is dropped.

__init__(**kwargs)

Creates a ParameterSet with given parameters.

checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

hasAttribute(name)

Returns True if self has attribute name.

out()

Generate the Gaussian profile

Link against this method to get the output of this model.

releaseParameters(name)

Removes parameter name from the parameters.

showParameters()

Returns a description of the parameters.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.InterpolateOverBox(**kwargs)

Bases: esys.escriptcore.modelframe.ParameterSet

Returns values at each time. The values are defined through given values at time node. For two dimensional domains back values are ignored.

Note:Instance variable domain - domain
Note:Instance variable value_left_bottom_front - (in) value at left,bottom,front corner
Note:Instance variable value_right_bottom_front - (in) value at right, bottom, front corner
Note:Instance variable value_left_top_front - (in) value at left,top,front corner
Note:Instance variable value_right_top_front - (in) value at right,top,front corner
Note:Instance variable value_left_bottom_back - (in) value at left,bottom,back corner
Note:Instance variable value_right_bottom_back - (in) value at right,bottom,back corner
Note:Instance variable value_left_top_back - (in) value at left,top,back corner
Note:Instance variable value_right_top_back - (in) value at right,top,back corner
__init__(**kwargs)

Creates a ParameterSet with given parameters.

checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

hasAttribute(name)

Returns True if self has attribute name.

out()

values at domain locations by bilinear interpolation of the given values.

Link against this method to get the output of this model.

releaseParameters(name)

Removes parameter name from the parameters.

showParameters()

Returns a description of the parameters.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.InterpolatedTimeProfile(**kwargs)

Bases: esys.escriptcore.modelframe.ParameterSet

Returns values at each time. The values are defined through given values at time node.

value[i] defines the value at time nodes[i]. Between nodes linear interpolation is used.

For time t<nodes[0], value[0] is used and for t>nodes[l], values[l] is used where l=len(nodes)-1.

Note:Instance variable t - (in) current time
Note:Instance variable node - (in) list of time nodes
Note:Instance variable values - (in) list of values at time nodes
__init__(**kwargs)

Creates a ParameterSet with given parameters.

checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

hasAttribute(name)

Returns True if self has attribute name.

out()

current value

Link against this method to get the output of this model.

releaseParameters(name)

Removes parameter name from the parameters.

showParameters()

Returns a description of the parameters.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.LinearCombination(**kwargs)

Bases: esys.escriptcore.modelframe.ParameterSet

Returns a linear combination of the f0*v0+f1*v1+f2*v2+f3*v3+f4*v4

Variables:
  • f0 – numerical object or None, default=None (in)
  • v0 – numerical object or None, default=None (in)
  • f1 – numerical object or None, default=None (in)
  • v1 – numerical object or None, default=None (in)
  • f2 – numerical object or None, default=None (in)
  • v2 – numerical object or None, default=None (in)
  • f3 – numerical object or None, default=None (in)
  • v3 – numerical object or None, default=None (in)
  • f4 – numerical object or None, default=None (in)
  • v4 – numerical object or None, default=None (in)
__init__(**kwargs)

Creates a ParameterSet with given parameters.

checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

hasAttribute(name)

Returns True if self has attribute name.

out()

returns f0*v0+f1*v1+f2*v2+f3*v3+f4*v4. Link against this method to get the output of this model.

releaseParameters(name)

Removes parameter name from the parameters.

showParameters()

Returns a description of the parameters.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.LinearPDE(domain, numEquations=None, numSolutions=None, isComplex=False, debug=False)

Bases: esys.escriptcore.linearPDEs.LinearProblem

This class is used to define a general linear, steady, second order PDE for an unknown function u on a given domain defined through a Domain object.

For a single PDE having a solution with a single component the linear PDE is defined in the following form:

-(grad(A[j,l]+A_reduced[j,l])*grad(u)[l]+(B[j]+B_reduced[j])u)[j]+(C[l]+C_reduced[l])*grad(u)[l]+(D+D_reduced)=-grad(X+X_reduced)[j,j]+(Y+Y_reduced)

where grad(F) denotes the spatial derivative of F. Einstein’s summation convention, ie. summation over indexes appearing twice in a term of a sum performed, is used. The coefficients A, B, C, D, X and Y have to be specified through Data objects in Function and the coefficients A_reduced, B_reduced, C_reduced, D_reduced, X_reduced and Y_reduced have to be specified through Data objects in ReducedFunction. It is also allowed to use objects that can be converted into such Data objects. A and A_reduced are rank two, B, C, X, B_reduced, C_reduced and X_reduced are rank one and D, D_reduced, Y and Y_reduced are scalar.

The following natural boundary conditions are considered:

n[j]*((A[i,j]+A_reduced[i,j])*grad(u)[l]+(B+B_reduced)[j]*u)+(d+d_reduced)*u=n[j]*(X[j]+X_reduced[j])+y

where n is the outer normal field. Notice that the coefficients A, A_reduced, B, B_reduced, X and X_reduced are defined in the PDE. The coefficients d and y are each a scalar in FunctionOnBoundary and the coefficients d_reduced and y_reduced are each a scalar in ReducedFunctionOnBoundary.

Constraints for the solution prescribe the value of the solution at certain locations in the domain. They have the form

u=r where q>0

r and q are each scalar where q is the characteristic function defining where the constraint is applied. The constraints override any other condition set by the PDE or the boundary condition.

The PDE is symmetrical if

A[i,j]=A[j,i] and B[j]=C[j] and A_reduced[i,j]=A_reduced[j,i] and B_reduced[j]=C_reduced[j]

For a system of PDEs and a solution with several components the PDE has the form

-grad((A[i,j,k,l]+A_reduced[i,j,k,l])*grad(u[k])[l]+(B[i,j,k]+B_reduced[i,j,k])*u[k])[j]+(C[i,k,l]+C_reduced[i,k,l])*grad(u[k])[l]+(D[i,k]+D_reduced[i,k]*u[k] =-grad(X[i,j]+X_reduced[i,j])[j]+Y[i]+Y_reduced[i]

A and A_reduced are of rank four, B, B_reduced, C and C_reduced are each of rank three, D, D_reduced, X_reduced and X are each of rank two and Y and Y_reduced are of rank one. The natural boundary conditions take the form:

n[j]*((A[i,j,k,l]+A_reduced[i,j,k,l])*grad(u[k])[l]+(B[i,j,k]+B_reduced[i,j,k])*u[k])+(d[i,k]+d_reduced[i,k])*u[k]=n[j]*(X[i,j]+X_reduced[i,j])+y[i]+y_reduced[i]

The coefficient d is of rank two and y is of rank one both in FunctionOnBoundary. The coefficients d_reduced is of rank two and y_reduced is of rank one both in ReducedFunctionOnBoundary.

Constraints take the form

u[i]=r[i] where q[i]>0

r and q are each rank one. Notice that at some locations not necessarily all components must have a constraint.

The system of PDEs is symmetrical if

  • A[i,j,k,l]=A[k,l,i,j]
  • A_reduced[i,j,k,l]=A_reduced[k,l,i,j]
  • B[i,j,k]=C[k,i,j]
  • B_reduced[i,j,k]=C_reduced[k,i,j]
  • D[i,k]=D[i,k]
  • D_reduced[i,k]=D_reduced[i,k]
  • d[i,k]=d[k,i]
  • d_reduced[i,k]=d_reduced[k,i]

LinearPDE also supports solution discontinuities over a contact region in the domain. To specify the conditions across the discontinuity we are using the generalised flux J which, in the case of a system of PDEs and several components of the solution, is defined as

J[i,j]=(A[i,j,k,l]+A_reduced[[i,j,k,l])*grad(u[k])[l]+(B[i,j,k]+B_reduced[i,j,k])*u[k]-X[i,j]-X_reduced[i,j]

For the case of single solution component and single PDE J is defined as

J[j]=(A[i,j]+A_reduced[i,j])*grad(u)[j]+(B[i]+B_reduced[i])*u-X[i]-X_reduced[i]

In the context of discontinuities n denotes the normal on the discontinuity pointing from side 0 towards side 1 calculated from FunctionSpace.getNormal of FunctionOnContactZero. For a system of PDEs the contact condition takes the form

n[j]*J0[i,j]=n[j]*J1[i,j]=(y_contact[i]+y_contact_reduced[i])- (d_contact[i,k]+d_contact_reduced[i,k])*jump(u)[k]

where J0 and J1 are the fluxes on side 0 and side 1 of the discontinuity, respectively. jump(u), which is the difference of the solution at side 1 and at side 0, denotes the jump of u across discontinuity along the normal calculated by jump. The coefficient d_contact is of rank two and y_contact is of rank one both in FunctionOnContactZero or FunctionOnContactOne. The coefficient d_contact_reduced is of rank two and y_contact_reduced is of rank one both in ReducedFunctionOnContactZero or ReducedFunctionOnContactOne. In case of a single PDE and a single component solution the contact condition takes the form

n[j]*J0_{j}=n[j]*J1_{j}=(y_contact+y_contact_reduced)-(d_contact+y_contact_reduced)*jump(u)

In this case the coefficient d_contact and y_contact are each scalar both in FunctionOnContactZero or FunctionOnContactOne and the coefficient d_contact_reduced and y_contact_reduced are each scalar both in ReducedFunctionOnContactZero or ReducedFunctionOnContactOne.

Typical usage:

p = LinearPDE(dom)
p.setValue(A=kronecker(dom), D=1, Y=0.5)
u = p.getSolution()
__init__(domain, numEquations=None, numSolutions=None, isComplex=False, debug=False)

Initializes a new linear PDE.

Parameters:
  • domain (Domain) – domain of the PDE
  • numEquations – number of equations. If None the number of equations is extracted from the PDE coefficients.
  • numSolutions – number of solution components. If None the number of solution components is extracted from the PDE coefficients.
  • debug – if True debug information is printed
addPDEToLumpedSystem(operator, a, b, c, hrz_lumping)

adds a PDE to the lumped system, results depend on domain

Parameters:
  • mat (OperatorAdapter) –
  • rhs (Data) –
  • a (Data) –
  • b (Data) –
  • c (Data) –
  • hrz_lumping (bool) –
addPDEToRHS(righthandside, X, Y, y, y_contact, y_dirac)

adds a PDE to the right hand side, results depend on domain

Parameters:
  • mat (OperatorAdapter) –
  • righthandside (Data) –
  • X (Data) –
  • Y (Data) –
  • y (Data) –
  • y_contact (Data) –
  • y_dirac (Data) –
addPDEToSystem(operator, righthandside, A, B, C, D, X, Y, d, y, d_contact, y_contact, d_dirac, y_dirac)

adds a PDE to the system, results depend on domain

Parameters:
  • mat (OperatorAdapter) –
  • rhs (Data) –
  • A (Data) –
  • B (Data) –
  • C (Data) –
  • D (Data) –
  • X (Data) –
  • Y (Data) –
  • d (Data) –
  • y (Data) –
  • d_contact (Data) –
  • y_contact (Data) –
  • d_dirac (Data) –
  • y_dirac (Data) –
addToRHS(rhs, data)

adds a PDE to the right hand side, results depend on domain

Parameters:
  • mat (OperatorAdapter) –
  • righthandside (Data) –
  • data (list) –
addToSystem(op, rhs, data)

adds a PDE to the system, results depend on domain

Parameters:
  • mat (OperatorAdapter) –
  • rhs (Data) –
  • data (list) –
alteredCoefficient(name)

Announces that coefficient name has been changed.

Parameters:name (string) – name of the coefficient affected
Raises:IllegalCoefficient – if name is not a coefficient of the PDE
Note:if name is q or r, the method will not trigger a rebuild of the system as constraints are applied to the solved system.
checkReciprocalSymmetry(name0, name1, verbose=True)

Tests two coefficients for reciprocal symmetry.

Parameters:
  • name0 (str) – name of the first coefficient
  • name1 (str) – name of the second coefficient
  • verbose (bool) – if set to True or not present a report on coefficients which break the symmetry is printed
Returns:

True if coefficients name0 and name1 are reciprocally symmetric.

Return type:

bool

checkSymmetricTensor(name, verbose=True)

Tests a coefficient for symmetry.

Parameters:
  • name (str) – name of the coefficient
  • verbose (bool) – if set to True or not present a report on coefficients which break the symmetry is printed.
Returns:

True if coefficient name is symmetric

Return type:

bool

checkSymmetry(verbose=True)

Tests the PDE for symmetry.

Parameters:verbose (bool) – if set to True or not present a report on coefficients which break the symmetry is printed.
Returns:True if the PDE is symmetric
Return type:bool
Note:This is a very expensive operation. It should be used for degugging only! The symmetry flag is not altered.
createCoefficient(name)

Creates a Data object corresponding to coefficient name.

Returns:the coefficient name initialized to 0
Return type:Data
Raises:IllegalCoefficient – if name is not a coefficient of the PDE
createOperator()

Returns an instance of a new operator.

createRightHandSide()

Returns an instance of a new right hand side.

createSolution()

Returns an instance of a new solution.

getCoefficient(name)

Returns the value of the coefficient name.

Parameters:name (string) – name of the coefficient requested
Returns:the value of the coefficient
Return type:Data
Raises:IllegalCoefficient – if name is not a coefficient of the PDE
getCurrentOperator()

Returns the operator in its current state.

getCurrentRightHandSide()

Returns the right hand side in its current state.

getCurrentSolution()

Returns the solution in its current state.

getDim()

Returns the spatial dimension of the PDE.

Returns:the spatial dimension of the PDE domain
Return type:int
getDomain()

Returns the domain of the PDE.

Returns:the domain of the PDE
Return type:Domain
getDomainStatus()

Return the status indicator of the domain

getFlux(u=None)

Returns the flux J for a given u.

J[i,j]=(A[i,j,k,l]+A_reduced[A[i,j,k,l]]*grad(u[k])[l]+(B[i,j,k]+B_reduced[i,j,k])u[k]-X[i,j]-X_reduced[i,j]

or

J[j]=(A[i,j]+A_reduced[i,j])*grad(u)[l]+(B[j]+B_reduced[j])u-X[j]-X_reduced[j]

Parameters:u (Data or None) – argument in the flux. If u is not present or equals None the current solution is used.
Returns:flux
Return type:Data
getFunctionSpaceForCoefficient(name)

Returns the FunctionSpace to be used for coefficient name.

Parameters:name (string) – name of the coefficient enquired
Returns:the function space to be used for coefficient name
Return type:FunctionSpace
Raises:IllegalCoefficient – if name is not a coefficient of the PDE
getFunctionSpaceForEquation()

Returns the FunctionSpace used to discretize the equation.

Returns:representation space of equation
Return type:FunctionSpace
getFunctionSpaceForSolution()

Returns the FunctionSpace used to represent the solution.

Returns:representation space of solution
Return type:FunctionSpace
getNumEquations()

Returns the number of equations.

Returns:the number of equations
Return type:int
Raises:UndefinedPDEError – if the number of equations is not specified yet
getNumSolutions()

Returns the number of unknowns.

Returns:the number of unknowns
Return type:int
Raises:UndefinedPDEError – if the number of unknowns is not specified yet
getOperator()

Returns the operator of the linear problem.

Returns:the operator of the problem
getOperatorType()

Returns the current system type.

getRequiredOperatorType()

Returns the system type which needs to be used by the current set up.

getResidual(u=None)

Returns the residual of u or the current solution if u is not present.

Parameters:u (Data or None) – argument in the residual calculation. It must be representable in self.getFunctionSpaceForSolution(). If u is not present or equals None the current solution is used.
Returns:residual of u
Return type:Data
getRightHandSide()

Returns the right hand side of the linear problem.

Returns:the right hand side of the problem
Return type:Data
getShapeOfCoefficient(name)

Returns the shape of the coefficient name.

Parameters:name (string) – name of the coefficient enquired
Returns:the shape of the coefficient name
Return type:tuple of int
Raises:IllegalCoefficient – if name is not a coefficient of the PDE
getSolution()

Returns the solution of the PDE.

Returns:the solution
Return type:Data
getSolverOptions()

Returns the solver options

Return type:SolverOptions
getSystem()

Returns the operator and right hand side of the PDE.

Returns:the discrete version of the PDE
Return type:tuple of Operator and Data
getSystemStatus()

Return the domain status used to build the current system

hasCoefficient(name)

Returns True if name is the name of a coefficient.

Parameters:name (string) – name of the coefficient enquired
Returns:True if name is the name of a coefficient of the general PDE, False otherwise
Return type:bool
initializeSystem()

Resets the system clearing the operator, right hand side and solution.

insertConstraint(rhs_only=False)

Applies the constraints defined by q and r to the PDE.

Parameters:rhs_only (bool) – if True only the right hand side is altered by the constraint
introduceCoefficients(**coeff)

Introduces new coefficients into the problem.

Use:

p.introduceCoefficients(A=PDECoef(…), B=PDECoef(…))

to introduce the coefficients A and B.

invalidateOperator()

Indicates the operator has to be rebuilt next time it is used.

invalidateRightHandSide()

Indicates the right hand side has to be rebuilt next time it is used.

invalidateSolution()

Indicates the PDE has to be resolved if the solution is requested.

invalidateSystem()

Announces that everything has to be rebuilt.

isComplex()

Returns true if this is a complex-valued LinearProblem, false if real-valued.

Return type:bool
isHermitian()

Checks if the pde is indicated to be Hermitian.

Returns:True if a Hermitian PDE is indicated, False otherwise
Return type:bool
Note:the method is equivalent to use getSolverOptions().isHermitian()
isOperatorValid()

Returns True if the operator is still valid.

isRightHandSideValid()

Returns True if the operator is still valid.

isSolutionValid()

Returns True if the solution is still valid.

isSymmetric()

Checks if symmetry is indicated.

Returns:True if a symmetric PDE is indicated, False otherwise
Return type:bool
Note:the method is equivalent to use getSolverOptions().isSymmetric()
isSystemValid()

Returns True if the system (including solution) is still vaild.

isUsingLumping()

Checks if matrix lumping is the current solver method.

Returns:True if the current solver method is lumping
Return type:bool
preservePreconditioner(preserve=True)

Notifies the PDE that the preconditioner should not be reset when making changes to the operator.

Building the preconditioner data can be quite expensive (e.g. for multigrid methods) so if it is known that changes to the operator are going to be minor calling this method can speed up successive PDE solves.

Note:Not all operator types support this.
Parameters:preserve (bool) – if True, preconditioner will be preserved, otherwise it will be reset when making changes to the operator, which is the default behaviour.
reduceEquationOrder()

Returns the status of order reduction for the equation.

Returns:True if reduced interpolation order is used for the representation of the equation, False otherwise
Return type:bool
reduceSolutionOrder()

Returns the status of order reduction for the solution.

Returns:True if reduced interpolation order is used for the representation of the solution, False otherwise
Return type:bool
resetAllCoefficients()

Resets all coefficients to their default values.

resetOperator()

Makes sure that the operator is instantiated and returns it initialized with zeros.

resetRightHandSide()

Sets the right hand side to zero.

resetRightHandSideCoefficients()

Resets all coefficients defining the right hand side

resetSolution()

Sets the solution to zero.

setDebug(flag)

Switches debug output on if flag is True otherwise it is switched off.

Parameters:flag (bool) – desired debug status
setDebugOff()

Switches debug output off.

setDebugOn()

Switches debug output on.

setHermitian(flag=False)

Sets the Hermitian flag to flag.

Parameters:flag (bool) – If True, the Hermitian flag is set otherwise reset.
Note:The method overwrites the Hermitian flag set by the solver options
setHermitianOff()

Clears the Hermitian flag. :note: The method overwrites the Hermitian flag set by the solver options

setHermitianOn()

Sets the Hermitian flag. :note: The method overwrites the Hermitian flag set by the solver options

setReducedOrderForEquationOff()

Switches reduced order off for equation representation.

Raises:RuntimeError – if order reduction is altered after a coefficient has been set
setReducedOrderForEquationOn()

Switches reduced order on for equation representation.

Raises:RuntimeError – if order reduction is altered after a coefficient has been set
setReducedOrderForEquationTo(flag=False)

Sets order reduction state for equation representation according to flag.

Parameters:flag (bool) – if flag is True, the order reduction is switched on for equation representation, otherwise or if flag is not present order reduction is switched off
Raises:RuntimeError – if order reduction is altered after a coefficient has been set
setReducedOrderForSolutionOff()

Switches reduced order off for solution representation

Raises:RuntimeError – if order reduction is altered after a coefficient has been set.
setReducedOrderForSolutionOn()

Switches reduced order on for solution representation.

Raises:RuntimeError – if order reduction is altered after a coefficient has been set
setReducedOrderForSolutionTo(flag=False)

Sets order reduction state for solution representation according to flag.

Parameters:flag (bool) – if flag is True, the order reduction is switched on for solution representation, otherwise or if flag is not present order reduction is switched off
Raises:RuntimeError – if order reduction is altered after a coefficient has been set
setReducedOrderOff()

Switches reduced order off for solution and equation representation

Raises:RuntimeError – if order reduction is altered after a coefficient has been set
setReducedOrderOn()

Switches reduced order on for solution and equation representation.

Raises:RuntimeError – if order reduction is altered after a coefficient has been set
setReducedOrderTo(flag=False)

Sets order reduction state for both solution and equation representation according to flag.

Parameters:flag (bool) – if True, the order reduction is switched on for both solution and equation representation, otherwise or if flag is not present order reduction is switched off
Raises:RuntimeError – if order reduction is altered after a coefficient has been set
setSolution(u, validate=True)

Sets the solution assuming that makes the system valid with the tolrance defined by the solver options

setSolverOptions(options=None)

Sets the solver options.

Parameters:options (SolverOptions or None) – the new solver options. If equal None, the solver options are set to the default.
Note:The symmetry flag of options is overwritten by the symmetry flag of the LinearProblem.
setSymmetry(flag=False)

Sets the symmetry flag to flag.

Parameters:flag (bool) – If True, the symmetry flag is set otherwise reset.
Note:The method overwrites the symmetry flag set by the solver options
setSymmetryOff()

Clears the symmetry flag. :note: The method overwrites the symmetry flag set by the solver options

setSymmetryOn()

Sets the symmetry flag. :note: The method overwrites the symmetry flag set by the solver options

setSystemStatus(status=None)

Sets the system status to status if status is not present the current status of the domain is used.

setValue(**coefficients)

Sets new values to coefficients.

Parameters:
  • coefficients – new values assigned to coefficients
  • A (any type that can be cast to a Data object on Function) – value for coefficient A
  • A_reduced (any type that can be cast to a Data object on ReducedFunction) – value for coefficient A_reduced
  • B (any type that can be cast to a Data object on Function) – value for coefficient B
  • B_reduced (any type that can be cast to a Data object on ReducedFunction) – value for coefficient B_reduced
  • C (any type that can be cast to a Data object on Function) – value for coefficient C
  • C_reduced (any type that can be cast to a Data object on ReducedFunction) – value for coefficient C_reduced
  • D (any type that can be cast to a Data object on Function) – value for coefficient D
  • D_reduced (any type that can be cast to a Data object on ReducedFunction) – value for coefficient D_reduced
  • X (any type that can be cast to a Data object on Function) – value for coefficient X
  • X_reduced (any type that can be cast to a Data object on ReducedFunction) – value for coefficient X_reduced
  • Y (any type that can be cast to a Data object on Function) – value for coefficient Y
  • Y_reduced (any type that can be cast to a Data object on ReducedFunction) – value for coefficient Y_reduced
  • d (any type that can be cast to a Data object on FunctionOnBoundary) – value for coefficient d
  • d_reduced (any type that can be cast to a Data object on ReducedFunctionOnBoundary) – value for coefficient d_reduced
  • y (any type that can be cast to a Data object on FunctionOnBoundary) – value for coefficient y
  • d_contact (any type that can be cast to a Data object on FunctionOnContactOne or FunctionOnContactZero) – value for coefficient d_contact
  • d_contact_reduced (any type that can be cast to a Data object on ReducedFunctionOnContactOne or ReducedFunctionOnContactZero) – value for coefficient d_contact_reduced
  • y_contact (any type that can be cast to a Data object on FunctionOnContactOne or FunctionOnContactZero) – value for coefficient y_contact
  • y_contact_reduced (any type that can be cast to a Data object on ReducedFunctionOnContactOne or ReducedFunctionOnContactZero) – value for coefficient y_contact_reduced
  • d_dirac (any type that can be cast to a Data object on DiracDeltaFunctions) – value for coefficient d_dirac
  • y_dirac (any type that can be cast to a Data object on DiracDeltaFunctions) – value for coefficient y_dirac
  • r (any type that can be cast to a Data object on Solution or ReducedSolution depending on whether reduced order is used for the solution) – values prescribed to the solution at the locations of constraints
  • q (any type that can be cast to a Data object on Solution or ReducedSolution depending on whether reduced order is used for the representation of the equation) – mask for location of constraints
Raises:

IllegalCoefficient – if an unknown coefficient keyword is used

shouldPreservePreconditioner()

Returns true if the preconditioner / factorisation should be kept even when resetting the operator.

Return type:bool
trace(text)

Prints the text message if debug mode is switched on.

Parameters:text (string) – message to be printed
validOperator()

Marks the operator as valid.

validRightHandSide()

Marks the right hand side as valid.

validSolution()

Marks the solution as valid.

class esys.modellib.input.MergeConstraints(**kwargs)

Bases: esys.escriptcore.modelframe.ParameterSet

Returns a linear combination of the f0*v0+f1*v1+f2*v2+f3*v3+f4*v4

__init__(**kwargs)

Creates a ParameterSet with given parameters.

checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

hasAttribute(name)

Returns True if self has attribute name.

location_of_constraint()

return the values used to constrain a solution

Returns:the mask marking the locations of the constraints
Return type:escript.Scalar
releaseParameters(name)

Removes parameter name from the parameters.

showParameters()

Returns a description of the parameters.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

value_of_constraint()

return the values used to constrain a solution

Returns:values to be used at the locations of the constraints. If value is not given None is rerturned.
Return type:escript.Scalar
writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.Model(parameters=[], **kwargs)

Bases: esys.escriptcore.modelframe.ParameterSet

A Model object represents a process marching over time until a finalizing condition is fulfilled. At each time step an iterative process can be performed and the time step size can be controlled. A Model has the following work flow:

doInitialization()
while not terminateInitialIteration(): doInitialStep()
doInitialPostprocessing()
while not finalize():
    dt=getSafeTimeStepSize(dt)
    doStepPreprocessing(dt)
    while not terminateIteration(): doStep(dt)
    doStepPostprocessing(dt)
doFinalization()

where doInitialization, finalize, getSafeTimeStepSize, doStepPreprocessing, terminateIteration, doStepPostprocessing, doFinalization are methods of the particular instance of a Model. The default implementations of these methods have to be overwritten by the subclass implementing a Model.

__init__(parameters=[], **kwargs)

Creates a model.

Just calls the parent constructor.

UNDEF_DT = 1e+300
checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

doFinalization()

Finalizes the time stepping.

This function may be overwritten.

doInitialPostprocessing()

Finalises the initialization iteration process. This method is not called in case of a restart.

This function may be overwritten.

doInitialStep()

Performs an iteration step in the initialization phase. This method is not called in case of a restart.

This function may be overwritten.

doInitialization()

Initializes the time stepping scheme. This method is not called in case of a restart.

This function may be overwritten.

doStep(dt)

Executes an iteration step at a time step.

dt is the currently used time step size.

This function may be overwritten.

doStepPostprocessing(dt)

Finalises the time step.

dt is the currently used time step size.

This function may be overwritten.

doStepPreprocessing(dt)

Sets up a time step of step size dt.

This function may be overwritten.

finalize()

Returns False if the time stepping is finalized.

This function may be overwritten.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

getSafeTimeStepSize(dt)

Returns a time step size which can be safely used.

dt gives the previously used step size.

This function may be overwritten.

hasAttribute(name)

Returns True if self has attribute name.

releaseParameters(name)

Removes parameter name from the parameters.

setUp()

Sets up the model.

This function may be overwritten.

showParameters()

Returns a description of the parameters.

terminateInitialIteration()

Returns True if iteration at the inital phase is terminated.

terminateIteration()

Returns True if iteration on a time step is terminated.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.ParameterSet(parameters=[], **kwargs)

Bases: esys.escriptcore.modelframe.LinkableObject

A class which allows to emphasize attributes to be written and read to XML.

Leaves of an ESySParameters object can be:

  • a real number
  • an integer number
  • a string
  • a boolean value
  • a ParameterSet object
  • a Simulation object
  • a Model object
  • a numpy object
  • a list of booleans
  • any other object (not considered by writeESySXML and writeXML)

Example for how to create an ESySParameters object:

p11=ParameterSet(gamma1=1.,gamma2=2.,gamma3=3.)
p1=ParameterSet(dim=2,tol_v=0.001,output_file="/tmp/u.%3.3d.dx",runFlag=True,parm11=p11)
parm=ParameterSet(parm1=p1,parm2=ParameterSet(alpha=Link(p11,"gamma1")))

This can be accessed as:

parm.parm1.gamma=0.
parm.parm1.dim=2
parm.parm1.tol_v=0.001
parm.parm1.output_file="/tmp/u.%3.3d.dx"
parm.parm1.runFlag=True
parm.parm1.parm11.gamma1=1.
parm.parm1.parm11.gamma2=2.
parm.parm1.parm11.gamma3=3.
parm.parm2.alpha=1. (value of parm.parm1.parm11.gamma1)
__init__(parameters=[], **kwargs)

Creates a ParameterSet with given parameters.

checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

hasAttribute(name)

Returns True if self has attribute name.

releaseParameters(name)

Removes parameter name from the parameters.

showParameters()

Returns a description of the parameters.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.ScalarDistributionFromTags(**kwargs)

Bases: esys.escriptcore.modelframe.ParameterSet

creates a scalar distribution on a domain from tags, If tag_map is given the tags can be given a names and tag_map is used to map it into domain tags.

Variables:
  • domain – domain
  • default – default value
  • tag0 – tag 0
  • value0 – value for tag 0
  • tag1 – tag 1
  • value1 – value for tag 1
  • tag2 – tag 2
  • value2 – value for tag 2
  • tag3 – tag 3
  • value3 – value for tag 3
  • tag4 – tag 4
  • value4 – value for tag 4
  • tag5 – tag 5
  • value5 – value for tag 5
  • tag6 – tag 6
  • value6 – value for tag 6
  • tag7 – tag 7
  • value7 – value for tag 7
  • tag8 – tag 8
  • value8 – value for tag 8
  • tag9 – tag 9
  • value9 – value for tag 9
__init__(**kwargs)

Creates a ParameterSet with given parameters.

checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

hasAttribute(name)

Returns True if self has attribute name.

out()

returns a esys.escript.Data object Link against this method to get the output of this model.

releaseParameters(name)

Removes parameter name from the parameters.

showParameters()

Returns a description of the parameters.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.Sequencer(**kwargs)

Bases: esys.escriptcore.modelframe.Model

Runs through time until t_end is reached.

Variables:
  • t_end – model is terminated when t_end is passed, default 1 (in).
  • dt_max – maximum time step size, default Model.UNDEF_DT (in)
  • t – current time stamp (in/out). By default it is initialized with zero.
__init__(**kwargs)
UNDEF_DT = 1e+300
checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

doFinalization()

Finalizes the time stepping.

This function may be overwritten.

doInitialPostprocessing()

Finalises the initialization iteration process. This method is not called in case of a restart.

This function may be overwritten.

doInitialStep()

Performs an iteration step in the initialization phase. This method is not called in case of a restart.

This function may be overwritten.

doInitialization()

initialize time integration

doStep(dt)

Executes an iteration step at a time step.

dt is the currently used time step size.

This function may be overwritten.

doStepPostprocessing(dt)

Finalises the time step.

dt is the currently used time step size.

This function may be overwritten.

doStepPreprocessing(dt)

Sets up a time step of step size dt.

This function may be overwritten.

finalize()

returns true when t has reached t_end

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

getSafeTimeStepSize(dt)

returns dt_max

hasAttribute(name)

Returns True if self has attribute name.

releaseParameters(name)

Removes parameter name from the parameters.

setUp()

Sets up the model.

This function may be overwritten.

showParameters()

Returns a description of the parameters.

terminateInitialIteration()

Returns True if iteration at the inital phase is terminated.

terminateIteration()

Returns True if iteration on a time step is terminated.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

class esys.modellib.input.SmoothScalarDistributionFromTags(**kwargs)

Bases: esys.escriptcore.modelframe.ParameterSet

creates a smooth scalar distribution on a domain from region tags

Variables:
  • domain – domain
  • default – default value
  • tag0 – tag 0
  • value0 – value for tag 0
  • tag1 – tag 1
  • value1 – value for tag 1
  • tag2 – tag 2
  • value2 – value for tag 2
  • tag3 – tag 3
  • value3 – value for tag 3
  • tag4 – tag 4
  • value4 – value for tag 4
  • tag5 – tag 5
  • value5 – value for tag 5
  • tag6 – tag 6
  • value6 – value for tag 6
  • tag7 – tag 7
  • value7 – value for tag 7
  • tag8 – tag 8
  • value8 – value for tag 8
  • tag9 – tag 9
  • value9 – value for tag 9
__init__(**kwargs)

Creates a ParameterSet with given parameters.

checkLinkTargets(models, hash)

Returns a set of tuples (“<self>(<name>)”, <target model>) if the parameter <name> is linked to model <target model> but <target model> is not in the list of models. If a parameter is linked to another parameter set which is not in the hash list the parameter set is checked for its models. hash gives the call history.

declareParameter(**parameters)

Declares one or more new parameters and their initial value.

declareParameters(parameters)

Declares a set of parameters. parameters can be a list, a dictionary or a ParameterSet.

classmethod fromDom(esysxml, node)
getAttributeObject(name)

Returns the object stored for attribute name.

hasAttribute(name)

Returns True if self has attribute name.

out()

returns a esys.escript.Data object Link against this method to get the output of this model.

releaseParameters(name)

Removes parameter name from the parameters.

showParameters()

Returns a description of the parameters.

toDom(esysxml, node)

toDom method of Model class.

trace(msg)

If debugging is on, prints the message, otherwise does nothing.

writeXML(ostream=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='UTF-8'>)

Writes the object as an XML object into an output stream.

Functions

esys.modellib.input.exp(arg)

Returns e to the power of argument arg.

Parameters:arg (float, escript.Data, Symbol, numpy.ndarray.) – argument
Return type:float, escript.Data, Symbol, numpy.ndarray depending on the type of arg
Raises:TypeError – if the type of the argument is not expected
esys.modellib.input.inf(arg)

Returns the minimum value over all data points.

Parameters:arg (float, int, escript.Data, numpy.ndarray) – argument
Returns:minimum value of arg over all components and all data points
Return type:float
Raises:TypeError – if type of arg cannot be processed
esys.modellib.input.length(arg)

Returns the length (Euclidean norm) of argument arg at each data point.

Parameters:arg (float, escript.Data, Symbol, numpy.ndarray) – argument
Return type:float, escript.Data, Symbol depending on the type of arg
esys.modellib.input.sup(arg)

Returns the maximum value over all data points.

Parameters:arg (float, int, escript.Data, numpy.ndarray) – argument
Returns:maximum value of arg over all components and all data points
Return type:float
Raises:TypeError – if type of arg cannot be processed
esys.modellib.input.whereNegative(arg)

Returns mask of negative values of argument arg.

Parameters:arg (float, escript.Data, Symbol, numpy.ndarray) – argument
Return type:float, escript.Data, Symbol, numpy.ndarray depending on the type of arg
Raises:TypeError – if the type of the argument is not expected
esys.modellib.input.wherePositive(arg)

Returns mask of positive values of argument arg.

Parameters:arg (float, escript.Data, Symbol, numpy.ndarray.) – argument
Return type:float, escript.Data, Symbol, numpy.ndarray depending on the type of arg
Raises:TypeError – if the type of the argument is not expected

Others

  • log

Packages