getfem (version 2534)
index
/home/pommier/getfem/getfem++/gcc-4.0/interface/src/python/getfem.py

Getfem-interface classes.
 
Provides access to the pseudo-objects exported by the getfem-python interface.

 
Modules
       
numarray
sys

 
Classes
       
CvStruct
Eltm
Fem
GeoTrans
Integ
LevelSet
MdBrick
MdState
Mesh
MeshFem
MeshIm
Poly
Precond
Slice
Spmat

 
class CvStruct
    Descriptor for a convex structure.
 
  Methods defined here:
__del__(self)
__init__(self, *args)
__repr__(self)
basic_structure(self)
Get the simplest structure.
 
Synopsis:  cs=CvStruct.basic_structure()
For example, the 'basic structure' of the 6-node triangle, is the canonical
3-noded triangle.
dim(self)
Get the dimension of the convex structure.
face(self, F)
Return the structure of the face F.
facepts(self, F)
Return the list of point indices for the face F.
get(self, *args)
nbpts(self)
Get the number of points of the convex structure.

 
class Eltm
    Descriptor for an elementary matrix type.
 
  Methods defined here:
__del__(self)
__init__(self, *args)
Generates a descriptor for an elementary matrix type.
 
  * Eltm('base', @tfem FEM)
 Integration of shape functions on elements, using the fem FEM.
 
  * Eltm('grad', @tfem FEM)
 Integration of gradient of shape functions on elements, using the fem FEM.
 
  * Eltm('hessian', @tfem FEM)
 Integration of hessian of shape functions on elements, using the fem FEM.
 
  * Eltm('normal')
 The unit normal to the current convex face
 
  * Eltm('grad_geotrans')
 The gradient of the geometric transformation.  * Eltm('grad_geotrans_inv')
 The inverse of the gradient of the geometric transformation.  *
Eltm('product', @eltm A, @eltm B)
 Integration of the tensorial product of elementary matrices A and B.
 
  In order to obtain a numerical value of theses matrices, see MeshIm.eltm().

 
class Fem
    FEM (Finite Element Method) objects.
 
  Methods defined here:
__del__(self)
__init__(self, fem_name)
Build a FEM object from a string description.
 
* FEM_PK(N,K)
   classical Lagrange element PK on a simplex of dimension N
 * FEM_PK_DISCONTINUOUS(N,K[,alpha])}
    discontinuous Lagrange element PK on a simplex of dim N
 * FEM_QK(N,K)}
   classical Lagrange element QK on quadrangles, hexahedrons etc
 * FEM_QK_DISCONTINUOUS(N,K[,alpha])}
   discontinuous Lagrange element QK on quadrangles, hexahedrons etc
 * FEM_Q2_INCOMPLETE
   incomplete 2D Q2 element with 8 dof (serendipity Quad 8 element)  *
FEM_PK_PRISM(N,K)
   classical Lagrange element PK on a prism
 * FEM_PK_PRISM_DISCONTINUOUS(N,K[,alpha])
   classical discontinuous Lagrange element PK on a prism.
  * FEM_PK_WITH_CUBIC_BUBBLE(N,K)
   classical Lagrange element PK on a simplex with an additional volumic
bubble function.
 * FEM_P1_NONCONFORMING
   non-conforming P1 method on a triangle.
 * FEM_P1_BUBBLE_FACE(N)
   P1 method on a simplex with an additional bubble function on face 0.
 * FEM_P1_BUBBLE_FACE_LAG
   P1 method on a simplex with an additional lagrange dof on face 0.
  * FEM_PK_HIERARCHICAL(N,K)
   PK element with a hierarchical basis
 * FEM_QK_HIERARCHICAL(N,K)
   QK element with a hierarchical basis
 * FEM_PK_PRISM_HIERARCHICAL(N,K)
   PK element on a prism with a hierarchical basis
 * FEM_STRUCTURED_COMPOSITE(FEM, K)
   Composite fem on a grid with K divisions
 * FEM_PK_HIERARCHICAL_COMPOSITE(N,K,S)
   PK composite element on a grid with S subdivisions and with a hierarchical
basis
 * FEM_PK_FULL_HIERARCHICAL_COMPOSITE(N,K,S)
   PK composite element with S subdivisions and a hierarchical basis on both
degree and subdivision
 * FEM_PRODUCT(FEM1,FEM2)
   tensorial product of two polynomial elements
 * FEM_HERMITE(N)
   Hermite element P3 on a simplex of dimension N=1,2,3
 * FEM_ARGYRIS
   Argyris element P5 on the triangle.
 * FEM_HCT_TRIANGLE
    Hsieh-Clough-Tocher element on the triangle (composite P3   element which
is C^1), should be used with IM_HCT_COMPOSITE()   integration method.
 * FEM_QUADC1_COMPOSITE
   Quadrilateral element, composite P3 element and C^1 (16 dof).
 * FEM_REDUCED_QUADC1_COMPOSITE
   Quadrilateral element, composite P3 element and C^1 (12 dof).
 * FEM(:_(BRT0(N)
   Raviart-Thomas element of order 0 on a simplex of dimension N.
 * FEM(:_(BNEDELEC(N)
   Nedelec edge element of order 0 on a simplex of dimension N.
 
 * f=Fem('interpolated_fem', @mf MF1, @mim MIM2, [ivec blocked_dof])
  Build a special fem which is interpolated from another mesh_fem.  Using this
special finite element, it is possible to interpolate a given mesh_fem MF1 on
another mesh, given the integration method that will be used on this mesh.
Note that this finite element may be quite slow, and eats much memory.
__repr__(self)
__str__(self)
base_value(self, X)
Evaluate all base functions at the given point.
char(self)
Ouput a (unique) string representation of the FEM.
 
Synopsis: Fem.char()
This can be used to perform comparisons between two different FEM objects.
dim(self)
Return the dimension (dimension of the reference convex) of the FEM.
estimated_degree(self)
Return an estimation of the polynomial degree of the FEM.
get(self, *args)
grad_base_value(self, X)
Evaluate the gradient of all base functions at the given point.
hess_base_value(self, X)
Evaluate the Hessian of all base functions at the given point.
is_equivalent(self)
Return 0 if the FEM is not equivalent.
 
Synopsis: Fem.is_equivalent()
Equivalent FEM are evaluated on the reference convex. This is the case of most
classical FEMs.
is_lagrange(self)
Return 0 if the FEM is not of Lagrange type.
is_polynomial(self)
Return 0 if the basis functions are not polynomials.
nbdof(self, CV=None)
Return the number of DOF for the FEM.
 
Synopsis: Fem.nbdof([, int CV])
Some specific FEM (for example 'interpolated_fem') may require a convex number
CV to give their result (interpolated fems for example). In most of the case,
you can omit this convex number.
poly_str(self)
Return (if possible) the polynomial expressions of the base functions in the
reference element.
 
Synopsis: Fem.poly_str()
This function will of course fail for non-polynomial FEMs.
pts(self, CV=None)
Get the location of the degrees of freedom on the reference element.
 
Synopsis: Fem.pts([, int CV])
Some specific FEM may require a convex number CV to give their result
(interpolated fems for example). In most of the case, you can omit this convex
number.
target_dim(self)
Return the dimension of the target space.
 
Synopsis: Fem.target_dim()
The target space dimension is usually 1, except for vector FEMs (none of them
has been implemented in getfem++ for now)

 
class GeoTrans
     Methods defined here:
__del__(self)
__init__(self, *args)
General function for building descriptors to geometric transformations.  Name
can be:  * 'GT_PK(N,K)'   : Transformation on simplexes, dim N, degree K *
'GT_QK(N,K)'   : Transformation on parallelepipeds, dim N, degree K *
'GT_PRISM(N,K)'          : Transformation on prisms, dim N, degree K *
'GT_PRODUCT(a,b)'        : tensorial product of two transformations *
'GT_LINEAR_PRODUCT(a,b)' : Linear tensorial product of two transformations
__repr__(self)
__str__(self)
char(self)
Output a (unique) string representation of the @tgeotrans.
 
Synopsis: CvStruct.char()
This can be used to perform comparisons between two different @tgeotrans
objects.
dim(self)
Get the dimension of the geometric transformation.
 
Synopsis: CvStruct.dim()
This is the dimension of the source space, i.e. the dimension of the reference
convex.
get(self, *args)
is_linear(self)
Return 0 if the geometric transformation is not linear.
nbpts(self)
return the number of points of the geometric transformation.
normals(self)
Get the normals for each face of the reference convex
pts(self)
Return the reference convex points.
 
Synopsis: CvStruct.pts()
The points are stored in the columns of the output matrix.
transform(self, G, pts)
Apply the geometric transformation to a set of points.
 
Synopsis: CvStruct.transform(G,pts)
G is the set of vertices of the real convex, pts is the set of points (in the
reference convex) that are to be transformed. The corresponding set of points
in the real convex is returned.

 
class Integ
    Integration Method Objects.
 
  Methods defined here:
__del__(self)
__init__(self, *args)
Integ(stringing method) Return a FEM Integration Method from a string
description.
 
  Possible descriptions are (non exhaustive list):
 * IM_EXACT_SIMPLEX(n)
   Exact integration on simplices (works only with linear geometric
transformations and PK FEMs). Note that 'exact integration'   should be
avoided in general, since they only apply to linear   geometric
transformations, are quite slow, and subject to   numerical stability problems
for high degree FEMs.
 * IM_PRODUCT(a, b)   Product of two integration methods.
 
  * IM_EXACT_PARALLELEPIPED(n)
   Exact integration on parallelepipeds.
 
  * IM_EXACT_PRISM(n)
   Exact integration on prisms.
 
  * IM_GAUSS1D(K)
   Gauss method on the segment, order K (K=1,3,...99).
 
  * IM_NC(N,K)
   Newton-Cotes approximative integration on simplexes, order K.
 
  * IM_NC_PARALLELEPIPED(N,K)
   Product of Newton-Cotes integration on parallelepipeds.
 
  * IM_NC_PRISM(N,K)
   Product of Newton-Cotes integration on prisms.
 
  * IM_GAUSS_PARALLELEPIPED(N,K)
   Product of Gauss1D integration on parallelepipeds.
 
  * IM_TRIANGLE(K)
   Gauss methods on triangles (K=1,3,5,6,7,8,9,10,13,17,19).
 
  * IM_QUAD(K)
   Gauss methods on quadrilaterons (K=2, 3, 5, .. 17). Note that
IM_GAUSS_PARALLELEPIPED should be prefered for QK FEMs.
 
  * IM_TETRAHEDRON(K)
   Gauss methods on tetrahedrons (K=1, 2, 3, 5, 6 or 8).
 
  * IM_SIMPLEX4D(3)
   Gauss method on a 4-dimensional simplex.
 
  * IM_STRUCTURED_COMPOSITE(IM1, K)
    Composite method on a grid with K divisions.
 
  * IM_HCT_COMPOSITE(IM1)
   Composite integration suited to the HCT composite finite element.
 
  example:
   Integ('IM_PRODUCT(IM_GAUSS1D(5),IM_GAUSS1D(5))')
   is the same as:
   Integ('IM_GAUSS_PARALLELEPIPED(2,5)')
__repr__(self)
__str__(self)
char(self)
Ouput a (unique) string representation of the integration method.
 
Synopsis: Integ.char()
This can be used to  comparisons between two different @tinteg objects.
coeffs(self)
Returns the coefficients associated to each integration point.
dim(self)
Return the dimension of the ref. convex of the method.
face_coeffs(self, F)
Returns the coefficients associated to each integration of a face.
face_pts(self, F)
Return the list of integration points for a face.
get(self, *args)
is_exact(self)
Return 0 if the integration is an approximate one.
nbpts(self)
Return the total number of integration points.
 
Synopsis: Integ.nbpts()
Count the points for the volume integration, and points for surface
integration on each face of the reference convex. Raises an error for exact
integration methods.
pts(self)
Return the list of integration points (only for approximate methods).

 
class LevelSet
    Getfem Level-set.
 
  Methods defined here:
__del__(self)
__init__(self, *args)
General constructor for LevelSet objects.
get(self, *args)
set(self, *args)

 
class MdBrick
     Methods defined here:
__del__(self)
__init__(self, *args)
General constructor for MdBrick objects.
 
 * B=MdBrick('constraint', mdbrick parent, string CTYPE [, int numfem])
  Build a generic constraint brick.   It may be useful in some situations,
such as the Stokes problem where the pressure in defined modulo a constant. In
such a situation, this brick can be used to add an additional constraint on
the pressure value.       CTYPE has to be chosen among 'augmented',
'penalized', and 'eliminated'.  The constraint can be specified with
MdBrick.set_constraints(). Note that Dirichlet bricks (except the 'generalized
Dirichlet' one) are also specializations of the 'constraint' brick.
 * B=MdBrick('dirichlet', mdbrick parent, int BNUM, meshfem MFMULT, string
CTYPE [, int numfem])
  Build a Dirichlet condition brick which impose the value of a field along a
mesh boundary.  The BNUM parameter selects on which mesh region the Dirichlet
condition is imposed. CTYPE has to be chosen among 'augmented', 'penalized',
and 'eliminated'. The MFMULT may generally be taken as the meshfem of the
unknown, but for 'augmented' Dirichlet conditions, you may have to respect the
Inf-Sup condition and choose an adequate meshfem.
 * B=MdBrick('dirichlet_on_normal_component', mdbrick parent, int BNUM,
meshfem MFMULT, string CTYPE [, int numfem])
  Build a Dirichlet condition brick which imposes the value of the normal
component of a vector field.
 * B=MdBrick('dirichlet_on_normal_derivative', mdbrick parent, int BNUM,
meshfem MFMULT, string CTYPE [, int numfem])
  Build a Dirichlet condition brick which imposes the value of the normal
derivative of the unknown.
 * B=MdBrick('generalized_dirichlet', mdbrick parent, int BNUM [, int numfem])
  This is the "old" Dirichlet brick of getfem.  This brick can be used to
impose general Dirichlet conditions 'h(x)u(x) = r(x)' , however it may have
some issues with elaborated FEM (such as Argyris, etc). It should be avoided
when possible.
 * B=MdBrick('source_term', mdbrick parent, [, int BNUM=-1[, int numfem]])
  Add a boundary or volumic source term ( \int B.v ).  If BNUM is omitted (or
set to -1) , the brick adds a volumic source term on the whole mesh. For BNUM
>= 0, the source term is imposed on the mesh region BNUM. Use
MdBrick.set_param('source term',mf,B) to set the source term field. The source
term is expected as a vector field of size Q (with Q = qdim).
 * B=MdBrick('normal_source_term', mdbrick parent, int BNUM [, int numfem])
  Add a boundary source term ( \int (Bn).v ).  The source term is imposed on
the mesh region BNUM (which of course is not allowed to be a volumic region,
only boundary regions are allowed). Use MdBrick.set_param('source term',mf,B)
to set the source term field. The source term B is expected as tensor field of
size QxN (with Q = qdim, N = mesh dim). For example, if you consider an
elasticity problem, this brick may be used to impose a force on the boundary
with B as the stress tensor.
 * B=MdBrick('normal_derivative_source_term', mdbrick parent, int BNUM [, int
numfem])
  Add a boundary source term ( \int (\partial_n B).v ).  The source term is
imposed on the mesh region BNUM. Use MdBrick.set_param('source term',mf,B) to
set the source term field, which is expected as a vector field  of size Q
(with Q = qdim).
 * B=MdBrick('neumann KirchhoffLove source term', mdbrick parent, int BNUM [,
int numfem])
   Add a boundary source term for neumann Kirchhoff-Love plate problems
(should be used with the Kirchhoff-Love flavour of the bilaplacian brick).
 * B=MdBrick('qu_term', mdbrick parent, [, int BNUM [, int numfem]])
  Update the tangent matrix with a \int (Qu).v term.  The Q(x) parameter is a
matrix field of size qdim x qdim. An example of use is for the "iku" part of
Robin boundary conditions \partial_n u + iku = ...
 * B=MdBrick('mass_matrix', meshim mim, meshfem mf_u [,'real'|'complex'])
  Build a mass-matrix brick.
 * B=MdBrick('generic_elliptic', meshim MIM, meshfem mfu
[,'scalar'|'matrix'|'tensor'][,'real'|'complex'])
  Setup a generic elliptic problem ( (A*grad(U)).grad(V) ).  The brick
parameter 'A' may be a scalar field, a matrix field, or a tensor field
(default is scalar).
 * B=MdBrick('helmholtz', meshim MIM, meshfem mfu [,'real'|'complex'])
  Setup a Helmholtz problem.  The brick has one parameter, 'wave_number'.
 * B=MdBrick('isotropic_linearized_elasticity', meshim MIM, meshfem mfu)
  Setup a linear elasticity problem.  The brick has two scalar parameter,
'lambda' and 'mu' (the Lame coefficients).
 * B=MdBrick('linear_incompressibility_term', mdbrick parent, meshfem mf_p [,
int numfem])
  Add an incompressibily constraint (div u = 0).
 * B=MdBrick('nonlinear_elasticity', meshim MIM, meshfem mfu, string lawname)
  Setup a nonlinear elasticity (large deformations) problem.  The material law
can be chosen among - 'SaintVenant Kirchhoff' (linearized material law) -
'Mooney Rivlin' (to be used with the nonlinear incompressibily term) -
'Ciarlet Geymonat'
 * B=MdBrick('nonlinear_elasticity_incompressibility_term', mdbrick parent,
meshfem mf_p [, int numfem])
  Add an incompressibily constraint to a large strain elasticity problem.
 * B=MdBrick('small_deformations_plasticity', meshim MIM, meshfem mfu, scalar
THRESHOLD)
  Setup a plasticity problem (with small deformations).  The THRESHOLD
parameter is the maximum value of the Von Mises stress before 'plastification'
of the material.
 * B=MdBrick('dynamic', mdbrick parent, scalar rho [, int numfem])
  Dynamic brick. This brick is not ready.
 * B=MdBrick('navier_stokes', meshim MIM, meshfem mfu, meshfem mfp)
  Setup a Navier-Stokes problem (this brick is not ready, do not use it).
 * B=MdBrick('bilaplacian', meshim MIM, meshfem mfu, ['Kirchhoff-Love'])
  Setup a bilaplacian problem.  If the Kirchhoff-Love option is specified, the
Kirchhoff-Love plate model is used.
 * B=MdBrick('isotropic_linearized_plate', meshim MIM, meshim MIMSUB, meshfem
MF_UT, meshfem MF_U3, meshfem MF_THETA, scalar EPSILON)
   Setup a linear plate model brick (for moderately thick plates, using the
Reissner-Mindlin model). EPSILON is the plate thinkness, the @mf MF_UT and
MF_U3 are used respectively for the membrane displacement and the transverse
displacement of the plate. The @mf MF_THETA is the rotation of the normal
("section rotations").  The second integration method MIMSUB can be chosen
equal to MIM, or different if you want to perform sub-integration on the
transverse shear term (mitc4 projection).  This brick has two parameters
"lambda" and "mu" (the Lamé coefficients)
 * B=MdBrick('mixed_isotropic_linearized_plate', meshim MIM, meshfem MF_UT,
meshfem MF_U3, meshfem MF_THETA, scalar EPSILON)
   Setup a mixed linear plate model brick (for thin plates, using Kirchhoff-
Love model).  For a non-mixed version, use the bilaplacian brick.
 * B=MdBrick('plate_source_term', mdbrick parent, [, int BNUM=-1[, int
numfem]])
   Add a boundary or a volumic source term to a plate problem. This brick has
two parameters: "B" is the displacement (ut and u3) source term, "M" is the
moment source term (i.e. the source term on the rotation of the normal).
 * B=MdBrick('plate_simple_support', mdbrick parent, int BNUM, string CTYPE [,
int numfem])
   Add a "simple support" boundary condition to a plate problem (homogeneous
Dirichlet condition on the displacement, free rotation). CTYPE specifies how
the constraint is enforced ('penalized', 'augmented' or 'eliminated').
 * B=MdBrick('plate_clamped_support', mdbrick parent, int BNUM, string CTYPE
[, int numfem])
   Add a "clamped support" boundary condition to a plate problem (homogeneous
Dirichlet condition on the displacement and on the rotation). CTYPE specifies
how the constraint is enforced ('penalized', 'augmented' or 'eliminated').
 * B=MdBrick('plate_closing', mdbrick parent [, int numfem])
   Add a free edges condition for the mixed plate model brick.  This brick is
required when the mixed linearized plate brick is used. It must be inserted
after all other boundary conditions (the reason is that the brick has to
inspect all other boundary conditions to determine the number of disconnected
boundary parts which are free edges).
__repr__(self)
__str__(self)
dim(self)
Get the dimension of the main mesh (2 for a 2D mesh, etc).
get(self, *args)
is_coercive(self)
Return true if the problem is coercive.
is_complex(self)
Return true if the problem uses complex numbers.
is_linear(self)
Return true if the problem is linear.
is_symmetric(self)
Return true if the problem is symmetric.
memsize(self)
Return the amount of memory (in bytes) used by the model brick.
mixed_variables(self)
Identify the indices of mixed variables (typically the pressure, etc.) in the
tangent matrix.
nbdof(self)
Get the total number of dof of the current problem.
 
Synopsis: MdBrick.nbdof()
This is the sum of the brick specific dof plus the dof of the parent bricks.
param(self, parameter_name)
Get the parameter value.
 
Synopsis: MdBrick.param( string parameter_name)
When the parameter has been assigned a specific mesh_fem, it is returned  as a
large array (the last dimension being the mesh_fem dof). When no mesh_fem has
been assigned, the parameter is considered to be constant over the mesh.
param_list(self)
Get the list of parameters names.
 
Synopsis: MdBrick.param_list()
Each brick embeds a number of parameters (the Lamé coefficients for the
linearized elasticity brick, the wave number for the Helmholtz brick,...),
described as a (scalar, or vector, tensor etc) field on a mesh_fem. You can
read/change the parameter values with MdBrick.param() and MdBrick.set_param().
penalization_epsilon(self, eps)
Change the penalization coefficient of a constraint brick.
 
Synopsis: MdBrick.penalization_epsilon( eps)
This is only applicable to the bricks which inherit from the constraint brick,
such as the Dirichlet ones. And of course it is not effective when the
constraint is enforced via direct elimination or via Lagrange multipliers. The
default value of eps is 1e-9.
set(self, *args)
set_constraints(self, H, R)
Set the constraints imposed by a constraint brick.
 
Synopsis: MdBrick.set_constraints( mat H, vec R)
This is only applicable to the bricks which inherit from the constraint brick,
such as the Dirichlet ones. Imposes HU=R.
set_constraints_rhs(self, H, R)
Set the right hand side of the constraints imposed by a constraint brick.
 
Synopsis: MdBrick.set_constraints_rhs( mat H, vec R)
This is only applicable to the bricks which inherit from the constraint brick,
such as the Dirichlet ones.
set_param(self, name, *args)
Change the value of a brick parameter.
 
Synopsis: MdBrick.set_param( string name, {meshfem MF,V | V})
V should contain the new parameter value. If a meshfem is given , V should
hold the field values over that meshfem (i.e. its last dimension should be
MeshFem.nbdof()).
solve(self, mds, *args)
Run the standard getfem solver.
 
Synopsis: MdBrick.solve( @mdstate mds [,...])
Note that you should be able to use your own solver if you want (it is
possible to obtain the tangent matrix and its right hand side with the
MdState.tangent_matrix() etc.).   Various options can be specified:
 
  - 'noisy' or 'very noisy' : the solver will display some
 information showing the progress (residual values etc.).
 - 'max_iter', NIT         : set the maximum iterations numbers.
 - 'max_res', RES          : set the target residual value.
 - 'lsolver', SOLVERNAME   : select explicitely the solver used for< the
linear systems (the default value is 'auto', which lets getfem choose itself).
Possible values are 'superlu', 'mumps' (if supported), 'cg/ildlt', 'gmres/ilu'
and 'gmres/ilut'.
subclass(self)
Get the typename of the brick.
tresca(self, mds, MFVM)
Compute the Tresca stress criterion on the mesh_fem MFVM.
 
Synopsis: VM=MdBrick.tresca( @mdstate mds, meshfem MFVM)
Only available on bricks where it has a meaning: linearized elasticity,
plasticity, nonlinear elasticity..
von_mises(self, mds, MFVM)
Compute the Von Mises stress on the mesh_fem MFVM.
 
Synopsis: VM=MdBrick.von_mises( @mdstate mds, meshfem MFVM)
Only available on bricks where it has a meaning: linearized elasticity,
plasticity, nonlinear elasticity.. Note that in 2D it is not the "real" Von
Mises (which should take into account the 'plane stress' or 'plane strain'
aspect), but a pure 2D Von Mises.

 
class MdState
     Methods defined here:
__del__(self)
__init__(self, *args)
General constructor for MdState objects.
 
       These objects hold the global model data of a chain of
       MdBricks, such as the right hand side, the tangent matrix and
       the constraints.
 
       * MDS=MdState(mdbrick B)
       Build a modelstate for the brick B (selects the real or
       complex state from the complexity of B).
 
* MDS=MdState('real')
 Build a model state for real unknowns.
 
* MDS=MdState('complex')
 Build a model state for complex unknowns.
__repr__(self)
__str__(self)
clear(self)
Clear the model state.
compute_reduced_residual(self)
Compute the reduced residual from the residual and constraints.
compute_reduced_system(self)
Compute the reduced system from the tangent matrix and constraints.
compute_residual(self, B)
Compute the residual for the brick B.
compute_tangent_matrix(self, B)
Update the tangent matrix from the brick B.
constraints_matrix(self)
Return the constraints matrix stored in the model state.
constraints_nullspace(self)
Return the nullspace of the constraints matrix.
get(self, *args)
is_complex(self)
Return 0 is the model state is real, 1 if it is complex.
memsize(self)
Return the amount of memory (in bytes) used by the model state.
reduced_residual(self)
Return the residual on the reduced system.
reduced_tangent_matrix(self)
Return the reduced tangent matrix (i.e. the tangent matrix after elimination
of the constraints).
residual(self)
Return the residual.
set(self, *args)
state(self)
Return the vector of unknowns, which contains the solution after
MdBrick.solve().
tangent_matrix(self)
Return the tangent matrix stored in the model state.
unreduce(self, U)
Reinsert the constraint eliminated from the system.

 
class Mesh
    Class for getfem mesh objects.
 
  Methods defined here:
__del__(self)
__init__(self, *args)
General constructor for Mesh objects.
 
 *  M=Mesh('empty', int dim)
  Create a new empty mesh.
 *  M=Mesh('cartesian', vec X[, vec Y[, vec Z,..]])
  Build quickly a regular mesh of quadrangles, cubes, etc.
 * M=Mesh('regular_simplices', vec X[, vec Y[, vec Z,.., ]]['degree', int
K]['noised'])
  Mesh a n-dimensionnal parallelepipeded with simplices (triangles,
tetrahedrons etc) .  The optional degree may be used to build meshes with non
linear geometric transformations.
 *  M=Mesh('triangles_grid', vec X, vec Y)
  Build quickly a regular mesh of triangles.   This is a very limited and
somehow deprecated function (See also Mesh('ptND'), Mesh('regular_simplices')
and Mesh('cartesian')).
 *  M=Mesh('curved', M0, vec F)
  Build a curved (n+1)-dimensions mesh from a n-dimensions mesh M0.   The
points of the new mesh have one additional coordinate, given by the vector F.
This can be used to obtain meshes for shells. M0 may be a meshfem object, in
that case its linked mesh will be used.
 *  M=Mesh('prismatic', M0, int NLAY)
  Extrude a prismatic mesh M from a mesh M0.   In the additional dimension
there are NLAY layers of elements built from 0 to 1.
 *  M=Mesh('pt2D', mat P, ivec T[, int N])
  Build a mesh from a 2D triangulation.  Each column of P contains a point
coordinate, and each column of T contains the point indices of a triangle. N
is optional and is a zone number. If N is specified then only the zone number
'N' is converted (in that case, T is expected to have 4 rows, the fourth
containing these zone numbers).
 *  M=Mesh('ptND', mat P, imat T)
  Build a mesh from a N-dimensional "triangulation".  Similar function to
'pt2D', for building simplexes meshes from a triangulation given in T, and a
list of points given in P. The dimension of the mesh will be the number of
rows of P, and the dimension of the simplexes will be the number of rows of T.
 *  M=Mesh('load', string FILENAME)
  Load a mesh from a GETFEM++ ascii mesh file. See also Mesh.save(FILENAME).
 *  M=Mesh('from_string', string S)
  Load a mesh from a string description. For example, a string returned by
Mesh.char().
 *  M=Mesh('import', string FORMAT, string FILENAME)
  Import a mesh, FORMAT may be:
 
 - 'gmsh'   for a mesh created with gmsh ( http://www.geuz.org/gmsh )
 - 'gid'    for a mesh created with GiD  ( http://gid.cimne.upc.es )
 - 'am_fmt' for a mesh created with emc2 (
http://pauillac.inria.fr/cdrom/www/emc2/fra.htm )
 *  M=Mesh('clone', mesh M2)
  Create a copy of a mesh.
__repr__(self)
__str__(self)
add_convex(self, CVTR, CVPTS)
Add a new convex into the mesh.
 
Synopsis: IDX = Mesh.add_convex( @geotrans CVTR, mat CVPTS)
The convex structure (triangle, prism,...) is given by CVTR (obtained with
GeoTrans('...')), and its points are given by the columns of CVPTS. On return,
IDX contains the convex number. CVPTS might be a 3-dimensional array in order
to insert more than one convex (or a two dimensional array correctly shaped
according to Fortran ordering).
add_point(self, PT)
Insert new points in the mesh and return their #id.
 
Synopsis: [IDX]=Mesh.add_point( mat PT)
PT should be an [n x m] matrix , where n is the mesh dimension, and m is the
number of points that will be added to the mesh. On output, IDX contains the
indices of these new points.
 
 Remark: if some points are already part of the mesh (with a small tolerance
of approximately 1e-8), they won't be inserted again, and IDX will contain the
previously assigned indices of the points.
char(self)
Output a string description of the mesh.
curved_edges(self, N, CVLST=None)
[OBSOLETE FUNCTION! will be removed in a future release]\\
 
Synopsis: [E,C]=Mesh.curved_edges( int N [, CVLST])
More sophisticated version of Mesh.edges() designed for curved elements. This
one will return N (N>=2) points of the (curved) edges. With N==2, this is
equivalent to Mesh.edges(). Since the points are no more always part of the
mesh, their coordinates are returned instead of points number, in the array E
which is a [ mesh_dim x 2 x nb_edges ] array.  If the optional output argument
C is specified, it will contain the convex number associated with each edge.
cvid(self)
Return the list of all convex #id.\\
 
Synopsis: [CVID] = Mesh.cvid()
Note that their numbering is not supposed to be contiguous from 0 to
Mesh.nbcvs()-1,  especially if some points have been removed from the mesh.
You can use Mesh.optimize_structure() to enforce a contiguous numbering.
cvid_from_pid(self, PIDLST)
Returns the convex #ids that share the point #ids given in PIDLST.
cvstruct(self, CVLST=None)
Return an array of the convex structures.
 
Synopsis: [CVS, ivec CV2STRUC]=Mesh.cvstruct([ivec CVLST])
If CVLST is not given, all convexes are considered. Each convex structure is
listed once in CVS, and CV2STRUC maps the convexes indice in CVLST to the
indice of its structure in CVS.
del_convex(self, IDX)
Remove one or more convexes from the mesh.
 
Synopsis: Mesh.del_convex( IDX)
IDX should contain the convexes #ids, such as the ones returned bu
Mesh.add_convex().
del_convex_of_dim(self, DIM)
Remove all convexes of dimension listed in DIM.
 
Synopsis: Mesh.del_convex_of_dim( ivec DIM)
For example Mesh.del_convex_of_dim( [1,2]) remove all line segments, triangles
and quadrangles.
del_point(self, PIDLST)
Removes one or more points from the mesh.
 
Synopsis: Mesh.del_point( ivec PIDLST)
PIDLST should contain the  point #id, such as the one returned by the 'add
point' command.
delete_region(self, RLST)
Remove the regions whose #ids are listed in RLST
dim(self)
Get the dimension of the mesh (2 for a 2D mesh, etc).
edges(self, *args)
[OBSOLETE FUNCTION! will be removed in a future release]\\
 
Synopsis: [E,C]=Mesh.edges([, CVLST][,'merge'])
Return the list of edges of mesh M for the convexes listed in the row vector
CVLST. E is a 2 x nb_edges matrix containing point indices. If CVLST is
omitted, then the edges of all convexes are returned. If CVLST has two rows
then the first row is supposed to contain convex numbers, and the second face
numbers, of which the edges will be returned.  If 'merge' is indicated, all
common edges of convexes are merged in a single edge.  If the optional output
argument C is specified, it will contain the convex number associated with
each edge.
export_to_dx(self, filename, *args)
Exports a mesh to an OpenDX file.
 
Synopsis: Mesh.export_to_dx( string filename, ...
[,'ascii'][,'append'][,'as',string name,[,'serie',string
serie_name]][,'edges'])
See also MeshFem.export_to_dx() , Slice.export_to_dx().
export_to_vtk(self, filename, *args)
Exports a mesh to a VTK file .
 
Synopsis: Mesh.export_to_vtk( string filename, ... [,'ascii'][,'quality'])
If 'quality' is specified, an estimation of the quality of each convex will be
written to the file.  See also MeshFem.export_to_vtk() ,
Slice.export_to_vtk().
faces_from_cvid(self, CVLST)
Return a list of convexes faces from a list of convex #id.\\
 
Synopsis: [imat CVFLST]=Mesh.faces_from_cvid( ivec CVLST,[ 'merge'])
CVFLST is a two-rows matrix, the first row lists convex #ids, and the second
lists face numbers. The optional argument 'merge' merges faces shared by two
convexes of CVLST.
faces_from_pid(self, PIDLST)
Return the convex faces whose vertex #ids are in PIDLST.\\
 
Synopsis: [imat CVFLST]=Mesh.faces_from_pid( ivec PIDLST)
For a convex face to be returned, EACH of its points have to be listed in
PIDLST. On output, the first row of CVFLST contains the convex number, and the
second row contains the face number (local number in the convex).
geotrans(self, CVLST=None)
Returns an array of the geometric transformations.
 
Synopsis: [GT,GT2STRUCT]=Mesh.geotrans([CVLST])
See also Mesh.cvstruct().
get(self, *args)
max_cvid(self)
Return the maximum #id of all convexes in the mesh (see 'max pid').
max_pid(self)
Return the maximum #id of all points in the mesh (see 'max cvid').
memsize(self)
Return the amount of memory (in bytes) used by the mesh.
merge(self, M2)
Merge with the mesh M2.
 
Synopsis: Mesh.merge( @mesh M2)
Overlapping points won't be duplicated. If M2 is a mesh_fem object, its linked
mesh will be used.
nbcvs(self)
Get the number of convexes of the mesh.
nbpts(self)
Get the number of points of the mesh.
normal_of_face(self, CV, F, FPTNUM=None)
Evaluates the normal of convex CV, face F at the FPTNUMth point of the face.
If FPTNUM is not specified, then the normal is evaluated at each geometrical
node of the face.
normal_of_faces(self, CVFLST)
Evaluates (at face centers) the normals of convexes.
 
Synopsis: [mat N] = Mesh.normal_of_faces( imat CVFLST)
CVFLST is supposed to contain convex numbers in its first row and convex face
number in its second row.
optimize_structure(self)
Reset point and convex numbering.
 
Synopsis: Mesh.optimize_structure()
After optimisation, the points (resp. convexes) will be consecutively numbered
from 0 to Mesh.max_pid()-1 (resp. Mesh.max_cvid()-1).
outer_faces(self, CVLST=None)
Return the faces which are not shared by two convexes.
 
Synopsis: [CVFLST]=Mesh.outer_faces([, CVLST])
If CVLST is not given, it basically returns the mesh boundary. If CVLST is
given, it returns the boundary of the convex set whose #ids are listed in
CVLST.
pid(self)
Return the list of points #id of the mesh.\\
 
Synopsis: [ivec PID]=Mesh.pid()
Note that their numbering is not supposed to be contiguous from 0 to
Mesh.nbpts()-1,  especially if you destroyed some convexes. You can use
Mesh.optimize_structure() to enforce a contiguous numbering.
pid_from_coords(self, PT)
Search point #id whose coordinates are listed in PT.\\
 
Synopsis: [ivec PIDLST]=Mesh.pid_from_coords( mat PT)
PT is an array containing a list of point coordinates. On return, PIDLST is a
vector containing points #id for each point found, and -1  -1 for those which
where not found in the mesh.
pid_from_cvid(self, CVLST=None)
Return the points attached to each convex of the mesh.\\
 
Synopsis:  [PID,IDX]=Mesh.pid_from_cvid([,CVLST])
If CVLST is omitted, all the convexes will be considered (equivalent to CVLST
= @Mesh.max_cvid()). IDX is a vector, length(IDX) = length(CVLIST)+1. PID is a
vector containing the concatenated list of points of each convex in cvlst.
Each entry of IDX is the position of the corresponding convex point list in
PID. Hence, for example, the list of points of the second convex is
PID[IDX(2):IDX(3)].\  If CVLST contains convex #id which do not exist in the
mesh, their point list will be empty.
pts(self, PIDLST=None)
Return the list of point coordinates of the mesh.\\
 
Synopsis: [mat PT]=Mesh.pts([, ivec PIDLST])
Each column of the returned matrix contains the coordinates of one point.  If
the optional argument PIDLST was given, only the points whose #id is listed in
this vector are returned. Otherwise, the returned matrix will have
Mesh.max_pid() columns, which might be greater than Mesh.nbpts() (if some
points of the mesh have been destroyed and no call to
Mesh.optimize_structure() have been issued).  The columns corresponding to
deleted points will be filled with NaN. You can use Mesh.pid() to filter such
invalid points.
quality(self, CVLST=None)
Return an estimation of the quality of each convex (0 <= Q <= 1).
refine(self, CVLST=None)
Use a Bank strategy for mesh refinement.
 
Synopsis: Mesh.refine([, ivec CVLST])
If CVLST is not given, the whole mesh is refined. Note that the regions, and
the finite element methods and integration methods of the mesh_fem and mesh_im
objects linked to this mesh will be automagically refined.
region(self, rnum)
Return the list of convexes/faces on the region 'rnum'.
 
Synopsis: CVFLST = Mesh.region( int rnum)
On output, the first row of I contains the convex numbers, and the second row
contains the face numbers (and -1 when the whole convex is in the region).
region_intersect(self, R1, R2)
Replace the region number R1 with its intersection with region number R2.
region_merge(self, R1, R2)
Merge region number R2 into region number R1.
region_substract(self, R1, R2)
Replace the region number R1 with its difference with region number R2.
regions(self)
Return the list of valid regions stored in the mesh.
save(self, FILENAME)
Save the mesh object to an ascii file.
 
Synopsis: Mesh.save( stringing FILENAME)
This mesh can be restored with Mesh('load', FILENAME).
set(self, *args)
set_pts(self, P)
Replace the coordinates of the mesh points with those given in P.
set_region(self, rnum, CVFLST)
Assigns the region number rnum to the convex faces stored in each column of
the matrix CVFLST.   The first row of CVFLST contains a convex number, and the
second row contains a face number in the convex (or -1 for the whole convex --
regions are usually used to store a list of convex faces, but you may also use
them to store a list of convexes).
transform(self, T)
Applies the matrix T to each point of the mesh.
 
Synopsis: Mesh.transform( mat T)
Note that T is not required to be a NxN matrix (with N=Mesh.dim()). Hence it
is possible to transform a 2D mesh into a 3D one (and reciprocally).
translate(self, V)
Translates each point of the mesh from V.
triangulated_surface(self, Nrefine, CVLIST=None)
[OBSOLETE FUNCTION! will be removed in a future release]
 
Synopsis: [mat T]=Mesh.triangulated_surface( int Nrefine [,CVLIST])
Similar function to Mesh.curved_edges() : split (if necessary, i.e. if the
geometric transformation if non-linear) each face into sub-triangles and
return their coordinates in T (see also compute_eval on P1 tri mesh(mf, U, ))

 
class MeshFem
     Methods defined here:
__del__(self)
__init__(self, *args)
General constructor for MeshFem objects.
 
MeshFem(mesh M [, int Qdim=1])
Build a new MeshFem object. The Qdim parameter is optional.
 * MeshFem('load', fname[, mesh M])
  Load a meshfem from a file.   If the mesh M is not supplied (this kind of
file does not store the mesh), then it is read from the file and its
descriptor is returned as the second output argument.
 * MeshFem('from_string', str[, mesh M])
  Create a mesh_fem object from its string description.  See also
MeshFem.char()
 * MeshFem('clone', meshfem MF2)
  Create a copy of a mesh_fem.
__repr__(self)
__str__(self)
char(self, opt=None)
Output a string description of the mesh_fem.
 
Synopsis: MeshFem.char([, opt])
By default, it does not include the description of the linked mesh object,
except if opt is 'with_mesh'
dof_from_cv(self, CVLST)
Return the DoF of the convexes listed in CVLST.
 
Synopsis: I = MeshFem.dof_from_cv( CVLST)
WARNING: the Degree of Freedom might be returned in ANY order, do not use this
function in your assembly routines. Use 'dof from cvid' instead, if you want
to be able to map a convex number with its associated degrees of freedom.  One
can also get the list of dof on a set on convex faces, by indicating on the
second row of CVLST the faces numbers (with respect to the convex number on
the first row).
dof_from_cvid(self, CVLST=None)
Return the degrees of freedom attached to each convex of the mesh.\\
 
Synopsis:  [DOFS,IDX]=MeshFem.dof_from_cvid([,CVLST])
If CVLST is omitted, all the convexes will be considered (equivalent to CVLST
= 1 ... @Mesh.max_cvid()).  IDX is a vector, length(IDX) = length(CVLIST)+1.
DOFS is a vector containing the concatenated list of dof of each convex in
CVLST. Each entry of IDX is the position of the corresponding convex point
list in DOFS. Hence, for example, the list of points of the second convex is
DOFS[IDX(2):IDX(3)].\  If CVLST contains convex #id which do not exist in the
mesh, their point list will be empty.
dof_nodes(self, DOFLST=None)
Get location of Degrees of Freedom.
 
Synopsis: [DOF_XY] = MeshFem.dof_nodes([, DOFLST])
Return the list of interpolation points for the specified dof #IDs in DOFLST
(if DOFLST is omitted, all DoF are considered).
dof_on_region(self, RLIST)
Return the list of dof lying on one of the mesh regions listed in RLIST.
 
Synopsis: DOFLST = MeshFem.dof_on_region( RLIST)
More precisely, this function returns the DoF whose support is non-null on one
of regions whose #ids are listed in RLIST (note that for boundary regions,
some dof nodes may not lie exactly on the boundary, for example the dof of
PK(n,0) lies on the center of the convex, but the base function in not null on
the convex border).
dof_partition(self, DOFP)
Change the dof_partition array.
 
Synopsis: MeshFem.dof_partition( ivec DOFP)
DOFP is a vector holding a integer value for each convex of the mesh_fem. See
MeshFem.dof_partition() for a description of "dof partition".
eval(self, expression)
interpolate an expression on the (lagrangian) MeshFem.
 
Examples:
  mf.eval('x[0]*x[1]') interpolates the function 'x*y'
  mf.eval('[x[0],x[1]]') interpolates the vector field '[x,y]'
export_to_dx(self, FILENAME, *args)
Export a mesh_fem and some fields to an OpenDX file.
 
Synopsis: MeshFem.export_to_dx( string FILENAME, ... ['as', string
mesh_name][,'edges']['serie',string serie_name][,'ascii'][,'append'], U,
'name'...)
This function will fail if the mesh_fem mixes different convex types (i.e.
quads and triangles), or if OpenDX does not handle a specific element type
(i.e. prism connections are not known by OpenDX).  The FEM will be mapped to
order 1 PK (or QK) FEMs. If you need to represent high-order FEMs or high-
order geometric transformations, you should consider @Slice.export_to_dx().
export_to_vtk(self, FILENAME, *args)
Export a mesh_fem and some fields to a vtk file.
 
Synopsis: MeshFem.export_to_vtk( string FILENAME, ... ['ascii'], U, 'name'...)
The FEM and geometric transformations will be mapped to order 1 or 2
isoparametric PK (or QK) FEMs (as VTK does not handle higher order elements).
If you need to represent high-order FEMs or high-order geometric
transformations, you should consider @Slice.export_to_vtk().
fem(self, CVLST=None)
Return a list of FEM used by the meshfem.
 
Synopsis: [FEMLST, CV2F] = MeshFem.fem([, CVLST])
FEMLST is an array of all @tfem objects found in the convexes given in CVLST.
If CV2F was supplied as an output argument, it contains, for each convex
listed in CVLST, the index of its correspounding FEM in FEMLST.   Convexes
which are not part of the mesh, or convexes which do not have any FEM have
their correspounding entry in CV2F set to -1.
get(self, *args)
interpolate_convex_data(self, Ucv)
Interpolate data given on each convex of the mesh to the meshfem dof.  The
meshfem has to be lagrangian, and should be discontinuous (typically a
FEM_PK(N,0) or FEM_QK(N,0) should be used).  The last dimension of the input
vector Ucv should have Mesh.max_cvid() elements.  Example of use:
MeshFem.interpolate_convex_dataMesh.quality())
is_equivalent(self, CVLST=None)
Test if the meshfem is equivalent.
 
Synopsis: I = MeshFem.is_equivalent([,CVLST])
See MeshFem.is_lagrangian()
is_lagrangian(self, CVLST=None)
Test if the meshfem is Lagrangian.
 
Synopsis: I = MeshFem.is_lagrangian([,CVLST])
Lagrangian means that each base function Phi[i] is such that
    Phi[i](P[j]) = delta(i,j),
 where P[j] is the DoF location of the jth base function, and delta(i,j) = 1
if i==j, else 0.
  If CVLST is omitted, it returns 1 if all convexes in the mesh are
Lagrangian. If CVLST is used, it returns the convex indices (with respect to
CVLST) which are Lagrangian.
is_polynomial(self, CVLST=None)
Test if all base functions are polynomials.
 
Synopsis: I = MeshFem.is_polynomial([,CVLST])
See MeshFem.is_lagrangian()
linked_mesh(self)
Return a reference to the mesh object linked to MF.
memsize(self)
Return the amount of memory (in bytes) used by the mesh_fem object.
 
Synopsis: MeshFem.memsize()
The result does not take into account the linked mesh object.
nbdof(self)
Return the number of Degrees of Freedom (DoF) of the meshfem MF.
non_conformal_dof(self, CVLST=None)
Return partially linked degrees of freedom.
 
Synopsis: MeshFem.non_conformal_dof([,CVLST])
Return the dof located on the border of a convex and which belong to only one
convex, except the ones which are located on the border of the mesh.  For
example, if the convex 'a' and 'b' share a common face, 'a' has a P1 FEM, and
'b' has a P2 FEM, then the dof on the middle of the face will be returned by
this function (this can be useful when searching the interfaces between
classical FEM and hierarchical FEM).
qdim(self)
Return the dimension Q of the field interpolated by the mesh_fem.
 
Synopsis: MeshFem.qdim()
By default, Q=1 (scalar field). This has an impact on the DOF numbering.
save(self, filename, opt=None)
Save a meshfem in a text file (and optionaly its linked mesh object if opt is
the string 'with_mesh').
set(self, *args)
set_classical_discontinuous_fem(self, *args)
Assigns a classical (Lagrange polynomial) discontinuous fem or order K.
 
Synopsis: MeshFem.set_classical_discontinuous_fem( int K, [int IM_DEGREE
[,ivec CVIDX]])
Similar to MeshFem.set_classical_fem() except that FEM_PK_DISCONTINUOUS is
used.
set_classical_fem(self, K, CVIDX=None)
Assign a classical (Lagrange polynomial) fem of order K to the meshfem.
 
Synopsis: MeshFem.set_classical_fem( int K, [,ivec CVIDX])
Uses FEM_PK for simplexes, FEM_QK for parallelepipeds etc.
set_fem(self, FEM, CVIDX=None)
Set the Finite Element Method.
 
Synopsis: MeshFem.set_fem( @tfem FEM [, ivec CVIDX])
Assign a FEM to all convexes whose #ids are listed in CVIDX. If CVIDX is not
given, the integration is assigned to all convexes.  See the help of Fem to
obtain a list of available FEM methods.
set_qdim(self, Q)
Change the Q dimension of the field that is interpolated by the meshfem.
 
Synopsis: MeshFem.set_qdim( int Q)
Q=1 means that the meshfem describes a scalar field, Q=N means that the
meshfem describes a vector field of dimension N.

 
class MeshIm
     Methods defined here:
__del__(self)
__init__(self, *args)
General constructor for MeshIm objects.
 
MeshIm(mesh M, [{Integ Im|int IM_DEGREE}])
Build a new MeshIm object. For convenience, optional arguments (IM or IM_DEGREE) can be provided, in that case a call to MESHIM:SET('integ') is issued with these arguments.
 
 * MeshIm('load', fname[, mesh M])
  Load a meshim from a file.   If the mesh M is not supplied (this kind of
file does not store the mesh), then it is read from the file and its
descriptor is returned as the second output argument.
 * MeshIm('from_string', str[, mesh M])
  Create a meshim object from its string description.  See also MeshIm.char()
 * MeshIm('clone', meshim MIM2)
  Create a copy of a meshim.
__repr__(self)
__str__(self)
char(self)
Output a string description of the mesh_im.
 
Synopsis: MeshIm.char([,'with mesh'])
By default, it does not include the description of the linked mesh object.
eltm(self, MET, CV, F=None)
Return the elementary matrix (or tensor) integrated on the convex CV.
 
Synopsis: M = MeshIm.eltm( @eltm MET, int CV [int F])
!!WARNING!! Be sure that the fem used for the construction of MET is
compatible with the fem assigned to element CV ! This is not checked by the
function ! If the argument F is given, then the elementary tensor is
integrated on the face F of CV instead of the whole convex.
get(self, *args)
integ(self, CVLIST=None)
Return a list of integration methods used by the meshim.
 
Synopsis: [INTEG, CV2I] = MeshIm.integ([, CVLIST])
INTEG is an array of all @tinteg objects found in the convexes given in CVLST.
If CV2F was supplied as an output argument, it contains, for each convex
listed in CVLST, the index of its correspounding integration method in INTEG.
Convexes which are not part of the mesh, or convexes which do not have any
integration method have their correspounding entry in CV2I set to -1.
linked_mesh(self)
Returns a reference to the mesh object linked to MIM.
memsize(self)
Return the amount of memory (in bytes) used by the mesh_im object.
 
Synopsis: MeshIm.memsize()
The result does not take into account the linked mesh object.
save(self, filename)
Saves a meshim in a text file (and optionaly its linked mesh object).
set(self, *args)
set_integ(self, *args)
Set the integration method.
 
Synopsis: MeshIm.set_integ( {integ IM|int IM_DEGREE} [, ivec CVIDX])
Assign an integration method to all convexes whose #ids are listed in CVIDX.
If CVIDX is not given, the integration is assigned to all convexes. It is
possible to assign a specific integration method with an integration method
handle IM obtained via Integ('IM_SOMETHING'), or to let getfem choose a
suitable integration method with IM_DEGREE (choosen such that polynomials of
degree <= IM_DEGREE are exactly integrated. If IM_DEGREE=-1, then the dummy
integration method IM_NONE will be used.)

 
class Poly
    

 
class Precond
    Getfem preconditioner.
 
  Methods defined here:
__del__(self)
__init__(self, *args)
General constructor for getfem preconditioners.
 
  The preconditioners may store REAL or COMPLEX values. They accept getfem
sparse matrices and Matlab sparse matrices.
 
  * Precond('identity')
 Create a REAL identity precondioner.
 
  * Precond('cidentity')
 Create a COMPLEX identity precondioner.
 
  * Precond('diagonal', @dcvec D)
 Create a diagonal precondioner.
 
  * Precond('ildlt', @spmat M)
 Create an ILDLT (Cholesky) preconditioner for the (symmetric) sparse matrix
M.
 This preconditioner has the same sparsity pattern than M (no fill-in).
 
  * Precond('ilu', @spmat M)
 Create an ILU (Incomplete LU) preconditioner for the sparse matrix M.
 This preconditioner has the same sparsity pattern than M (no fill-in).
 
  * Precond('ildltt', @spmat M [, int fillin [, scalar threshold]])
 Create an ILDLT (Cholesky with filling) preconditioner for the (symmetric)
sparse matrix M.
  The preconditioner may add at most 'fillin' additional non-zero entries on
each line. The default value for 'fillin' is 10, and the default threshold is
1e-7.
 
  * Precond('ilut', @spmat M [, int fillin [, scalar threshold]])
 Create an ILUT (Incomplete LU with filling) preconditioner for the sparse
matrix M.
  The preconditioner may add at most 'fillin' additional non-zero entries on
each line. The default value for 'fillin' is 10, and the default threshold is
1e-7.
 
  * Precond('superlu', @spmat M)
 Uses SuperLU to build an exact factorization of the sparse matrix M.  This
preconditioner is only available if the getfem-interface was built with
SuperLU support. Note that LU factorization is likely to eat all your memory
for 3D problems.
get(self, *args)
info(self)
Return a short informative string about the preconditioner.
is_complex(self)
Return 1 if the preconditioner stores complex values.
mult(self, V)
Apply the preconditioner to the supplied vector.
set(self, *args)
size(self)
Return the dimensions of the preconditioner.
tmult(self, V)
Apply the transposed preconditioner to the supplied vector.
type(self)
Return a string describing the type of the preconditioner ('ilu', 'ildlt',..).

 
class Slice
    Mesh slices.
 
The slices may be considered as a (non-conformal) mesh of
simplexes which provides fast interpolation on a P1-discontinuous
MeshFem.
 
It is used mainly for post-processing purposes.
 
  Methods defined here:
__del__(self)
__init__(self, *args)
General constructor for Slice objects.
 
sl = Slice(sliceop, mesh M, int REFINE [, CVFLST])
 sl = Slice(sliceop, mesh_fem MF, vec U, int REFINE [, CVFLST])
 sl = Slice(sliceop, slice SL)
 sl = Slice('streamlines', mesh_fem MF, vec U, mat SEEDS)
 sl = Slice('points', mesh M, mat PTS)
 sl = Slice('load', mesh M)
 
  Creation of a mesh slice. Mesh slices are very similar to a P1-discontinuous
mesh_fem on which interpolation is very fast. The slice is built from a mesh
object, and a description of the slicing operation, for example,
 
  sl = Slice(('planar',+1,[0;0],[1;0]), m, 5);
 
  cuts the original mesh with the half space {y>0}. Each convex of the
original mesh m is simplexified (for example a quadrangle is splitted into 2
triangles), and each simplex is refined 5 times.
 
  Slicing operations can be:
  - cutting with a plane, a sphere or a cylinder
  - intersection or union of slices
  - isovalues surfaces/volumes
  - "points", "streamlines" (see below)
 
  If the first argument is a mesh_fem mf instead of a mesh, and if it is
followed by a field U , then the deformation U will be applied to the mesh
before the slicing operation.
 
  The first argument can also be a slice.
 
  Slicing operations:
===============
 They are specified with TUPLES, do not forget the extra parentheses!. The
first element is the name of the operation, followed the slicing options.
 
   * ('none')
  Does not cut the mesh.
 
  * ('planar', orient, p, n)
   Planar cut. p and n define a half-space, p being a point belong to the
boundary of the half-space, and n being its normal. If orient is equal to -1
(resp. 0, +1), then the slicing operation will cut the mesh with the
"interior" (resp. "boundary", "exterior") of the half-space. Orient may also
be set to +2 which means that the mesh will be sliced, but both the outer and
inner parts will be kept.
 
  * ('ball', orient, c, r)
  Cut with a ball of center c and radius r.
 
  * ('cylinder', orient, p1, p2, r)
  Cut with a cylinder whose axis is the line (p1,p2) and whose radius is r.
 
  * ('isovalues',orient, mesh_fem MF, vec U, scalar V)
  Cut using the isosurface of the field U (defined on the mesh_fem MF). The
result is the set {x such that U(x) <= V} or {x such that U(x) == V} or {x
such that U(x) <= V} depending on the value of ORIENT.
 
  * ('boundary'[, SLICEOP])
  Return the boundary of the result of SLICEOP, where SLICEOP is any slicing
operation. If SLICEOP is not specified, then the whole mesh is considered
(i.e. it is equivalent to ('boundary',{'none'})).
 
  * ('explode', coef)
  Build an 'exploded' view of the mesh: each convex is shrinked (0 < coef <=
1). In the case of 3D convexes, only their faces are kept.
 
  * ('union', SLICEOP1, SLICEOP2)
 * ('intersection', SLICEOP1, SLICEOP2)
 * ('comp', SLICEOP)
 * ('diff', SLICEOP1, SLICEOP2)
  Boolean operations: returns the union,intersection,complementary or
difference of slicing operations.
 
  * ('mesh', MESH)
  Build a slice which is the intersection of the sliced mesh with another
mesh. The slice is such that all of its simplexes are stricl