esys.dudley Package

A domain meshed with triangles or tetrahedra only. Imports submodules into its namespace

Classes

class esys.dudley.DudleyDomain

Bases: esys.escriptcore.escriptcpp.ContinuousDomain

A concrete class representing a dudley domain. For more details, please consult the C++ documentation.

__init__((object)arg1, (DudleyDomain)arg2) → None
MPIBarrier((DudleyDomain)arg1) → None :

Wait until all processes have reached this point

addPDEToLumpedSystem((DudleyDomain)arg1, (Data)arg2, (Data)arg3, (Data)mat, (Data)D, (bool)d) → None :

adds a PDE onto the lumped stiffness matrix

Parameters:
  • mat (Data) –
  • D (Data) –
  • d (Data) –
  • useHRZ (bool) –
addPDEToRHS((DudleyDomain)arg1, (Data)arg2, (Data)rhs, (Data)X, (Data)Y, (Data)y, (Data)y_contact) → None :

adds a PDE onto the stiffness matrix mat and a rhs

Parameters:
  • rhs (Data) –
  • X (Data) –
  • Y (Data) –
  • y (Data) –
  • y_contact (Data) –
addPDEToSystem((DudleyDomain)arg1, (Operator)arg2, (Data)arg3, (Data)mat, (Data)rhs, (Data)A, (Data)B, (Data)C, (Data)D, (Data)X, (Data)Y, (Data)d, (Data)y, (Data)d_contact, (Data)y_contact) → None :

adds a PDE onto the stiffness matrix mat and a rhs

Parameters:
  • mat (OperatorAdapter) –
  • rhs (Data) –
  • A (Data) –
  • B (Data) –
  • C (Data) –
  • D (Data) –
  • X (Data) –
  • Y (Data) –
  • d (Data) –
  • d_contact (Data) –
  • y_contact (Data) –
addPDEToTransportProblem((DudleyDomain)arg1, (TransportProblem)arg2, (Data)arg3, (Data)tp, (Data)source, (Data)M, (Data)A, (Data)B, (Data)C, (Data)D, (Data)X, (Data)Y, (Data)d, (Data)y, (Data)d_contact, (Data)y_contact) → None :
Parameters:
  • tp (AbstractTransportProblem) –
  • source (Data) –
  • M (Data) –
  • A (Data) –
  • B (Data) –
  • C (Data) –
  • D (Data) –
  • X (Data) –
  • Y (Data) –
  • d (Data) –
  • y (Data) –
  • d_contact (Data) –
  • y_contact (Data) –
dump((DudleyDomain)arg1, (str)fileName) → None :

dumps the mesh to a file with the given name.

getDataShape((DudleyDomain)arg1, (int)functionSpaceCode) → object :
Returns:a pair (dps, ns) where dps=the number of data points per sample, and ns=the number of samples
Return type:tuple
getDescription((DudleyDomain)arg1) → str :
Returns:a description for this domain
Return type:string
getDim((DudleyDomain)arg1) → int :
Return type:int
getMPIRank((DudleyDomain)arg1) → int :
Returns:the rank of this process
Return type:int
getMPISize((DudleyDomain)arg1) → int :
Returns:the number of processes used for this Domain
Return type:int
getNormal((DudleyDomain)arg1) → Data :
Returns:boundary normals at the quadrature point on the face elements
Return type:Data
getNumDataPointsGlobal((DudleyDomain)arg1) → int :
Returns:the number of data points summed across all MPI processes
Return type:int
getSize((DudleyDomain)arg1) → Data :
Returns:the element size
Return type:Data
getStatus((Domain)arg1) → int :

The status of a domain changes whenever the domain is modified

Return type:int
getSystemMatrixTypeId((DudleyDomain)arg1, (object)options) → int :
Returns:the identifier of the matrix type to be used for the global stiffness matrix when particular solver options are used.
Return type:int
Parameters:options (SolverBuddy) –
getTag((DudleyDomain)arg1, (str)name) → int :
Returns:tag id for name
Return type:string
getTransportTypeId((DudleyDomain)arg1, (int)solver, (int)preconditioner, (int)package, (bool)symmetry) → int :
Returns:

the identifier of the transport problem type to be used when a particular solver, preconditioner, package and symmetric matrix is used.

Return type:

int

Parameters:
  • solver (int) –
  • preconditioner (int) –
  • package (int) –
  • symmetry (int) –
getX((DudleyDomain)arg1) → Data :
Returns:locations in the FEM nodes
Return type:Data
isValidTagName((DudleyDomain)arg1, (str)name) → bool :
Returns:True is name corresponds to a tag
Return type:bool
newOperator((DudleyDomain)arg1, (int)row_blocksize, (FunctionSpace)row_functionspace, (int)column_blocksize, (FunctionSpace)column_functionspace, (int)type) → Operator :

creates a stiffness matrix and initializes it with zeros

Parameters:
  • row_blocksize (int) –
  • row_functionspace (FunctionSpace) –
  • column_blocksize (int) –
  • column_functionspace (FunctionSpace) –
  • type (int) –
newTransportProblem((DudleyDomain)theta, (int)blocksize, (FunctionSpace)functionspace, (int)type) → TransportProblem :

creates a TransportProblem

Parameters:
  • theta (float) –
  • blocksize (int) –
  • functionspace (FunctionSpace) –
  • type (int) –
onMasterProcessor((DudleyDomain)arg1) → bool :
Returns:True if this code is executing on the master process
Return type:bool
print_mesh_info((DudleyDomain)arg1[, (bool)full=False]) → None :
Parameters:full (bool) –
setTagMap((DudleyDomain)arg1, (str)name, (int)tag) → None :

Give a tag number a name.

Parameters:
  • name (string) – Name for the tag
  • tag (int) – numeric id
Note:

Tag names must be unique within a domain

setX((DudleyDomain)arg1, (Data)arg) → None :

assigns new location to the domain

Parameters:arg (Data) –
showTagNames((DudleyDomain)arg1) → str :
Returns:A space separated list of tag names
Return type:string
supportsContactElements((Domain)arg1) → bool :

Does this domain support contact elements.

write((DudleyDomain)arg1, (str)filename) → None :

Write the current mesh to a file with the given name.

esys.dudley.GMSHDesign

alias of esys.pycad.gmsh.Design

Functions

esys.dudley.Brick(n0=1, n1=1, n2=1, order=1, l0=1.0, l1=1.0, l2=1.0, periodic0=False, periodic1=False, periodic2=False, integrationOrder=-1, reducedIntegrationOrder=-1, useElementsOnFace=False, useFullElementOrder=False, optimize=False, **kwargs)

__Brick_driver( (list)args) -> Domain

esys.dudley.LoadMesh([(str)fileName='file.nc']) → Domain :
Return type:DudleyDomain
esys.dudley.MakeDomain(design, integrationOrder=-1, reducedIntegrationOrder=-1, optimizeLabeling=True, useMacroElements=False)

Creates a Dudley Domain from a esys.pycad.design.Design object. Currently only gmsh is supported.

Parameters:
  • design (esys.pycad.design.Design) – the geometry
  • integrationOrder (int) – integration order. If -1 the default is used.
  • reducedIntegrationOrder (int) – reduced integration order. If -1 the default is used.
  • optimizeLabeling (bool) – if set the labeling of the mesh nodes is optimized
  • useMacroElements (bool) – for compatibility with finley. Must be False
Returns:

the Dudley domain defined by the design

Return type:

Domain

esys.dudley.ReadGmsh((str)fileName='file.msh', (int)numDim[, (int)integrationOrder=-1[, (int)reducedIntegrationOrder=-1[, (bool)optimize=True]]]) → Domain :

Read a gmsh mesh file

Return type:

Domain

Parameters:
  • fileName (string) –
  • integrationOrder (int) – order of the quadrature scheme. Always 2.
  • reducedIntegrationOrder – order of reduced quadrature scheme. Always 0.
  • optimize (bool) – Enable optimisation of node labels
esys.dudley.ReadMesh([(str)fileName='file.fly'[, (int)integrationOrder=-1[, (int)reducedIntegrationOrder=-1[, (bool)optimize=True]]]]) → Domain :

Read a mesh from a fly file. For MPI parallel runs fan out the mesh to multiple processes.

Return type:

Domain

Parameters:
  • fileName (string) –
  • integrationOrder (int) – order of the quadrature scheme. Ignored.
  • reducedIntegrationOrder – order of reduced quadrature scheme. Ignored.
  • optimize (bool) – Enable optimisation of node labels
esys.dudley.Rectangle(n0=1, n1=1, order=1, l0=1.0, l1=1.0, periodic0=False, periodic1=False, integrationOrder=-1, reducedIntegrationOrder=-1, useElementsOnFace=False, useFullElementOrder=False, optimize=False, **kwargs)

__Rectangle_driver( (list)args) -> Domain

esys.dudley.getTVTK(domain)

Returns the point and connectivity information of a finleymesh for tVTK :param domain: a finley domainsation of node labels :return connectivityinfo: numpy ndarray :return pointinfo: numpy ndarray

Others

Packages