To store and manipulate finite element meshes. More...
Public Member Functions | |
Mesh () | |
Default constructor (Empty mesh) | |
Mesh (const string &file, bool bc=false, int opt=NODE_DOF, int nb_dof=1) | |
Constructor using a mesh file. More... | |
Mesh (real_t xmin, real_t xmax, size_t nb_el, bool bc=false, size_t p=1, size_t nb_dof=1, int c1=0, int c2=0) | |
Constructor for a 1-D mesh. The domain is the interval [xmin,xmax]. More... | |
Mesh (const Grid &g, int opt=QUADRILATERAL) | |
Constructor for a uniform finite difference grid given by and instance of class Grid. More... | |
Mesh (const Grid &g, int shape, int opt) | |
Constructor of dual mesh for a uniform finite difference grid given by and instance of class Grid. More... | |
Mesh (real_t xmin, real_t xmax, size_t ne, int c1, int c2, int p=1, size_t nb_dof=1) | |
Constructor for a uniform 1-D finite element mesh. More... | |
Mesh (real_t xmin, real_t xmax, real_t ymin, real_t ymax, size_t nx, size_t ny, int cx0, int cxN, int cy0, int cyN, int opt=0, size_t nb_dof=1) | |
Constructor for a uniform 2-D structured finite element mesh. More... | |
Mesh (real_t xmin, real_t xmax, real_t ymin, real_t ymax, real_t zmin, real_t zmax, size_t nx, size_t ny, size_t nz, int cx0, int cxN, int cy0, int cyN, int cz0, int czN, int opt=0, size_t nb_dof=1) | |
Constructor for a uniform 3-D structured finite element mesh. More... | |
Mesh (const Mesh &m, const Point< real_t > &x_bl, const Point< real_t > &x_tr) | |
Constructor that extracts the mesh of a rectangular region from an initial mesh. More... | |
Mesh (const Mesh &mesh, int opt, size_t dof1, size_t dof2, bool bc=false) | |
Constructor that copies the input mesh and selects given degrees of freedom. More... | |
Mesh (const Mesh &ms) | |
Copy Constructor. More... | |
~Mesh () | |
Destructor. | |
void | setDim (size_t dim) |
Define space dimension. Normally, between 1 and 3. More... | |
void | Add (Node *nd) |
Add a node to mesh. More... | |
void | Add (Element *el) |
Add an element to mesh. More... | |
void | Add (Side *sd) |
Add a side to mesh. More... | |
void | Add (Edge *ed) |
Add an edge to mesh. More... | |
Mesh & | operator*= (real_t a) |
Operator *= More... | |
void | get (const string &mesh_file) |
Read mesh data in file. More... | |
void | get (const string &mesh_file, int ff, int nb_dof=1) |
Read mesh data in file with giving its format. More... | |
void | setDOFSupport (int opt, int nb_nodes=1) |
Define supports of degrees of freedom. More... | |
void | setNbDOFPerNode (size_t nb_dof=1) |
Define number of degrees of freedom for each node. More... | |
void | setPointInDomain (Point< real_t > x) |
Define a point in the domain. This function makes sense only if boundary mesh is given without internal mesh (Case of Boundary Elements) More... | |
void | removeImposedDOF () |
Eliminate equations corresponding to imposed DOF. | |
size_t | NumberEquations (size_t dof=0) |
Renumber Equations. More... | |
size_t | NumberEquations (size_t dof, int c) |
Renumber Equations. More... | |
int | getAllSides (int opt=0) |
Determine all mesh sides. More... | |
size_t | getNbSideNodes () const |
Return the number of nodes on each side. | |
size_t | getNbElementNodes () const |
Return the number of nodes in each element. | |
int | getBoundarySides () |
Determine all boundary sides. More... | |
int | createBoundarySideList () |
Create list of boundary sides. More... | |
int | getBoundaryNodes () |
Determine all boundary nodes. More... | |
int | createInternalSideList () |
Create list of internal sides (not on the boundary). More... | |
int | getAllEdges () |
Determine all edges. More... | |
void | getNodeNeighborElements () |
Create node neighboring elements. More... | |
void | getElementNeighborElements () |
Create element neighboring elements. More... | |
void | setMaterial (int code, const string &mname) |
Associate material to code of element. More... | |
void | Reorder (size_t m=GRAPH_MEMORY) |
Renumber mesh nodes according to reverse Cuthill Mc Kee algorithm. More... | |
void | Add (size_t num, real_t *x) |
Add a node by giving its label and an array containing its coordinates. More... | |
void | DeleteNode (size_t label) |
Remove a node given by its label. More... | |
void | DeleteElement (size_t label) |
Remove an element given by its label. More... | |
void | DeleteSide (size_t label) |
Remove a side given by its label. More... | |
void | Delete (Node *nd) |
Remove a node given by its pointer. More... | |
void | Delete (Element *el) |
Remove a node given by its pointer. More... | |
void | Delete (Side *sd) |
Remove a side given by its pointer. More... | |
void | Delete (Edge *ed) |
Remove an edge given by its pointer. More... | |
void | RenumberNode (size_t n1, size_t n2) |
Renumber a node. More... | |
void | RenumberElement (size_t n1, size_t n2) |
Renumber an element. More... | |
void | RenumberSide (size_t n1, size_t n2) |
Renumber a side. More... | |
void | RenumberEdge (size_t n1, size_t n2) |
Renumber an edge. More... | |
void | setList (const vector< Node *> &nl) |
Initialize list of mesh nodes using the input vector. More... | |
void | setList (const vector< Element *> &el) |
Initialize list of mesh elements using the input vector. More... | |
void | setList (const vector< Side *> &sl) |
Initialize list of mesh sides using the input vector. More... | |
void | Rescale (real_t sx, real_t sy=0., real_t sz=0.) |
Rescale mesh by multiplying node coordinates by constants. More... | |
size_t | getDim () const |
Return space dimension. | |
size_t | getNbNodes () const |
Return number of nodes. | |
size_t | getNbMarkedNodes () const |
Return number of marked nodes. | |
size_t | getNbVertices () const |
Return number of vertices. | |
size_t | getNbDOF () const |
Return total number of degrees of freedom (DOF) | |
size_t | getNbEq () const |
Return number of equations. | |
size_t | getNbEq (int i) const |
Return number of equations for the i-th set of degrees of freedom. | |
size_t | getNbElements () const |
Return number of elements. | |
size_t | getNbSides () const |
Return number of sides. | |
size_t | getNbEdges () const |
Return number of sides. | |
size_t | getNbBoundarySides () const |
Return number of boundary sides. More... | |
size_t | getNbInternalSides () const |
Return number of internal sides. More... | |
size_t | getNbMat () const |
Return number of materials. | |
void | AddMidNodes (int g=0) |
Add mid-side nodes. More... | |
Point< real_t > | getMaxCoord () const |
Return maximum coordinates of nodes. | |
Point< real_t > | getMinCoord () const |
Return minimum coordinates of nodes. | |
void | set (Node *nd) |
Replace node in the mesh. More... | |
void | set (Element *el) |
Replace element in the mesh. More... | |
void | set (Side *sd) |
Choose side in the mesh. More... | |
bool | NodesAreDOF () const |
Return information about DOF type. More... | |
bool | SidesAreDOF () const |
Return information about DOF type. More... | |
bool | EdgesAreDOF () const |
Return information about DOF type. More... | |
bool | ElementsAreDOF () const |
Return information about DOF type. More... | |
int | getDOFSupport () const |
Return information on dof support Return an integer according to enumerated values: NODE_DOF, ELEMENT_DOF SIDE_DOF. | |
void | Deform (const Vect< real_t > &u, real_t rate=0.2) |
Deform mesh according to a displacement vector. More... | |
void | put (const string &mesh_file) const |
Write mesh data on file. More... | |
void | save (const string &mesh_file) const |
Write mesh data on file in various formats. More... | |
bool | withImposedDOF () const |
Return true if imposed DOF count in equations, false if not. | |
bool | isStructured () const |
Return true is mesh is structured, false if not. | |
size_t | getNodeNewLabel (size_t n) const |
Return new label of node of a renumbered node. | |
void | getList (vector< Node *> &nl) const |
Fill vector nl with list of pointers to nodes. More... | |
void | getList (vector< Element *> &el) const |
Fill vector el with list of pointers to elements. More... | |
void | getList (vector< Side *> &sl) const |
Fill vector sl with list of pointers to sides. More... | |
Node * | getPtrNode (size_t i) const |
Return pointer to node with label i . | |
Node & | getNode (size_t i) const |
Return refenrece to node with label i | |
Element * | getPtrElement (size_t i) const |
Return pointer to element with label i | |
Side * | getPtrSide (size_t i) const |
Return pointer to side with label i | |
Side & | getSide (size_t i) const |
Return reference to side with label i | |
Edge * | getPtrEdge (size_t i) const |
Return pointer to edge with label i | |
size_t | getNodeLabel (size_t i) const |
Return label of i -th node. More... | |
size_t | getElementLabel (size_t i) const |
Return label of i -th element. More... | |
size_t | getSideLabel (size_t i) const |
Return label of i -th side. More... | |
size_t | getEdgeLabel (size_t i) const |
Return label of i -th edge. More... | |
int | getShape () const |
Determine shape of elements Return Shape index (see enum ElementShape) if all elements have the same shape, 0 if not. | |
Element * | operator() (size_t i) const |
Operator () : Return pointer to i-th element. | |
Node * | operator[] (size_t i) const |
Operator [] : Return pointer to i-th node. | |
size_t | operator() (size_t i, size_t n) const |
Operator () : Return pointer to i-th node of n-th element. | |
Mesh & | operator= (Mesh &ms) |
Operator = : Assign a Mesh instance. | |
Friends | |
void | Refine (Mesh &in_mesh, Mesh &out_mesh) |
Refine mesh. Subdivide each triangle into 4 subtriangles. This member function is valid for 2-D triangular meshes only. More... | |
To store and manipulate finite element meshes.
Class Mesh enables defining as an object a finite element mesh. A finite element mesh is characterized by its nodes, elements and sides. Each of these types of data constitutes a class in the OFELI library.
The standard procedure to introduce the finite element mesh is to provide an input file containing its data. For this, we have defined our own mesh data file (following the XML syntax). Of course, a developer can write his own function to read his finite element mesh file using the methods in Mesh.
Constructor using a mesh file.
[in] | file | File containing mesh data. The extension of the file yields the file format: The extension .m implies OFELI file format and .msh implies GMSH msh file. |
[in] | bc | Flag to remove (true) or not (false) imposed Degrees of Freedom [default: false] |
[in] | opt | Type of DOF support: To choose among enumerated values NODE_DOF , SIDE_DOF or ELEMENT_DOF .Say if degrees of freedom (unknowns) are supported by nodes, sides or elements. |
[in] | nb_dof | Number of degrees of freedom per node [Default: 1 ]. |
Mesh | ( | real_t | xmin, |
real_t | xmax, | ||
size_t | nb_el, | ||
bool | bc = false , |
||
size_t | p = 1 , |
||
size_t | nb_dof = 1 , |
||
int | c1 = 0 , |
||
int | c2 = 0 |
||
) |
Constructor for a 1-D mesh. The domain is the interval [xmin,xmax].
[in] | xmin | Value of xmin |
[in] | xmax | Value of xmax |
[in] | nb_el | Number of elements to generate |
[in] | bc | Flag to remove (true) or not (false) imposed Degrees of Freedom [default: false] |
[in] | p | Degree of finite element polynomial [default: 1] |
[in] | nb_dof | Number of degrees of freedom for each node [default: 1] |
[in] | c1 | Code to assign to first node (at x=xmin ) [default: 0] |
[in] | c2 | Code to assign to last node (at x=xmax ) [default: 0] |
Mesh | ( | const Grid & | g, |
int | opt = QUADRILATERAL |
||
) |
Constructor of dual mesh for a uniform finite difference grid given by and instance of class Grid.
[in] | g | Grid instance |
[in] | shape | Value to say which type of elements to generate
|
[in] | opt | This argument can take any value. It is here only to distinguish from the other constructor using Grid instance. |
Constructor for a uniform 1-D finite element mesh.
The domain is the line (xmin,xmax)
[in] | xmin | Minimal coordinate |
[in] | xmax | Maximal coordinate |
[in] | ne | Number of elements |
[in] | c1 | Code for the first node (x=xmin) |
[in] | c2 | Code for the last node (x=xmax) |
[in] | p | Degree of approximation polynomial [Default: 1 ]. |
[in] | nb_dof | Number of degrees of freedom per node [Default: 1 ]. |
Mesh | ( | real_t | xmin, |
real_t | xmax, | ||
real_t | ymin, | ||
real_t | ymax, | ||
size_t | nx, | ||
size_t | ny, | ||
int | cx0, | ||
int | cxN, | ||
int | cy0, | ||
int | cyN, | ||
int | opt = 0 , |
||
size_t | nb_dof = 1 |
||
) |
Constructor for a uniform 2-D structured finite element mesh.
The domain is the rectangle (xmin,xmax)x(ymin,ymax)
[in] | xmin | Minimal x-coordinate |
[in] | xmax | Maximal x-coordinate |
[in] | ymin | Minimal y-coordinate |
[in] | ymax | Maximal y-coordinate |
[in] | nx | Number of subintervals on the x-axis |
[in] | ny | Number of subintervals on the y-axis |
[in] | cx0 | Code for nodes generated on the line x=x0 if >0, for sides on this line if <0 |
[in] | cxN | Code for nodes generated on the line x=xN if >0, for sides on this line if <0 |
[in] | cy0 | Code for nodes generated on the line y=y0 if >0, for sides on this line if <0 |
[in] | cyN | Code for nodes generated on the line y=yN if >0, for sides on this line if <0 |
[in] | opt | Flag to generate elements as well (if not zero) [Default: 0 ]. If the flag is not 0, it can take one of the enumerated values: TRIANGLE or QUADRILATERAL, with obvious meaning. |
[in] | nb_dof | Number of degrees of freedom per node [Default: 1 ]. |
Mesh | ( | real_t | xmin, |
real_t | xmax, | ||
real_t | ymin, | ||
real_t | ymax, | ||
real_t | zmin, | ||
real_t | zmax, | ||
size_t | nx, | ||
size_t | ny, | ||
size_t | nz, | ||
int | cx0, | ||
int | cxN, | ||
int | cy0, | ||
int | cyN, | ||
int | cz0, | ||
int | czN, | ||
int | opt = 0 , |
||
size_t | nb_dof = 1 |
||
) |
Constructor for a uniform 3-D structured finite element mesh.
The domain is the parallepiped (xmin,xmax)x(ymin,ymax)x(zmin,zmax)
[in] | xmin | Minimal x-coordinate |
[in] | xmax | Maximal x-coordinate |
[in] | ymin | Minimal y-coordinate |
[in] | ymax | Maximal y-coordinate |
[in] | zmin | Minimal z-coordinate |
[in] | zmax | Maximal z-coordinate |
[in] | nx | Number of subintervals on the x-axis |
[in] | ny | Number of subintervals on the y-axis |
[in] | nz | Number of subintervals on the z-axis |
[in] | cx0 | Code for nodes generated on the line x=xmin if >0, for sides on this line if <0 |
[in] | cxN | Code for nodes generated on the line x=xmax if >0, for sides on this line if <0 |
[in] | cy0 | Code for nodes generated on the line y=ymin if >0, for sides on this line if <0 |
[in] | cyN | Code for nodes generated on the line y=ymax if >0, for sides on this line if <0 |
[in] | cz0 | Code for nodes generated on the line z=zmin if >0, for sides on this line if <0 |
[in] | czN | Code for nodes generated on the line z=zmax if >0, for sides on this line if <0 |
[in] | opt | Flag to generate elements as well (if not zero) [Default: 0]. If the flag is not 0, it can take one of the enumerated values: HEXAHEDRON or TETRAHEDRON, with obvious meaning. |
[in] | nb_dof | Number of degrees of freedom per node [Default: 1 ]. |
Constructor that extracts the mesh of a rectangular region from an initial mesh.
This constructor is useful for zooming purposes for instance.
[in] | m | Initial mesh from which the submesh is extracted |
[in] | x_bl | Coordinate of bottom left vertex of the rectangle |
[in] | x_tr | Coordinate of top right vertex of the rectangle |
Constructor that copies the input mesh and selects given degrees of freedom.
This constructor is to be used for coupled problems where each subproblem uses a choice of degrees of freedom.
[in] | mesh | Initial mesh from which the submesh is extracted |
[in] | opt | Type of DOF support: To choose among enumerated values NODE_DOF , SIDE_DOF or ELEMENT_DOF . |
[in] | dof1 | Label of first degree of freedom to select to the output mesh |
[in] | dof2 | Label of last degree of freedom to select to the output mesh |
[in] | bc | Flag to remove (true ) or not (false ) imposed Degrees of Freedom [Default: false ] |
void Add | ( | size_t | num, |
real_t * | x | ||
) |
Add a node by giving its label and an array containing its coordinates.
[in] | num | Label of node to add |
[in] | x | C-array of node coordinates |
void AddMidNodes | ( | int | g = 0 | ) |
Add mid-side nodes.
This is function is valid for triangles only
[in] | g | Option to say of barycentre node is to be added (>0) or not (=0) |
int createBoundarySideList | ( | ) |
Create list of boundary sides.
This function is useful to loop over boundary sides without testing Once this one is called, the function getNbBoundarySides() is available. Moreover, looping over boundary sides is available via the member functions topBoundarySide() and getBoundarySide()
int createInternalSideList | ( | ) |
Create list of internal sides (not on the boundary).
This function is useful to loop over internal sides without testing Once this one is called, the function getNbInternalSides() is available. Moreover, looping over internal sides is available via the member functions topInternalSide() and getInternalSide()
Deform mesh according to a displacement vector.
This function modifies node coordinates according to given displacement vector and given rate
[in] | u | Displacement vector |
[in] | rate | Maximal rate of deformation of resulting mesh. Its default value is 0.2, i.e. The resulting mesh has a maximum of deformation rate of 20% |
void Delete | ( | Node * | nd | ) |
Remove a node given by its pointer.
This function does not release the space previously occupied
[in] | nd | Pointer to node to delete |
void Delete | ( | Element * | el | ) |
Remove a node given by its pointer.
This function does not release the space previously occupied
[in] | el | Pointer to element to delete |
void Delete | ( | Side * | sd | ) |
Remove a side given by its pointer.
This function does not release the space previously occupied
[in] | sd | Pointer to side to delete |
void Delete | ( | Edge * | ed | ) |
Remove an edge given by its pointer.
This function does not release the space previously occupied
[in] | ed | Pointer to edge to delete |
void DeleteElement | ( | size_t | label | ) |
Remove an element given by its label.
This function does not release the space previously occupied
[in] | label | Label of element to delete |
void DeleteNode | ( | size_t | label | ) |
Remove a node given by its label.
This function does not release the space previously occupied
[in] | label | Label of node to delete |
void DeleteSide | ( | size_t | label | ) |
Remove a side given by its label.
This function does not release the space previously occupied
[in] | label | Label of side to delete |
bool EdgesAreDOF | ( | ) | const |
Return information about DOF type.
false
otherwise bool ElementsAreDOF | ( | ) | const |
Return information about DOF type.
false
otherwise void get | ( | const string & | mesh_file | ) |
void get | ( | const string & | mesh_file, |
int | ff, | ||
int | nb_dof = 1 |
||
) |
Read mesh data in file with giving its format.
File format can be chosen among a variety of choices. See "File Formats" page
[in] | mesh_file | Mesh file name |
[in] | ff | File format: Integer to chose among enumerated values: OFELI_FF , GMSH , MATLAB , EASYMESH , GAMBIT , BAMG , NETGEN , TRIANGLE_FF |
[in] | nb_dof | Number of degrees of freedom per node (Default value: 1) |
int getAllEdges | ( | ) |
Determine all edges.
int getAllSides | ( | int | opt = 0 | ) |
Determine all mesh sides.
int getBoundaryNodes | ( | ) |
Determine all boundary nodes.
int getBoundarySides | ( | ) |
Determine all boundary sides.
size_t getEdgeLabel | ( | size_t | i | ) | const |
Return label of i
-th edge.
[in] | i | Edge index |
size_t getElementLabel | ( | size_t | i | ) | const |
Return label of i
-th element.
[in] | i | Element index |
void getElementNeighborElements | ( | ) |
Create element neighboring elements.
This function creates for each element the list of elements that share a side with it. Once this function is invoked, one can retrieve the list of neighboring elements of any element (Element::getNeigborElement)
void getList | ( | vector< Node *> & | nl | ) | const |
Fill vector nl
with list of pointers to nodes.
[out] | nl | Instance of class vector that contain on output the list |
void getList | ( | vector< Element *> & | el | ) | const |
Fill vector el
with list of pointers to elements.
[out] | el | Instance of class vector that contain on output the list |
void getList | ( | vector< Side *> & | sl | ) | const |
Fill vector sl
with list of pointers to sides.
[out] | sl | Instance of class vector that contain on output the list |
size_t getNbBoundarySides | ( | ) | const |
Return number of boundary sides.
This function is valid if member function getAllSides or getBoundarySides has been invoked before
size_t getNbInternalSides | ( | ) | const |
Return number of internal sides.
This function is valid if member functions getAllSides and createInternalSideList have been invoked before
size_t getNodeLabel | ( | size_t | i | ) | const |
Return label of i
-th node.
[in] | i | Node index |
void getNodeNeighborElements | ( | ) |
Create node neighboring elements.
This function is generally useful when, for a numerical method, one looks for a given node to the list of elements that share this node. Once this function is invoked, one can retrieve the list of neighboring elements of any node (Node::getNeigEl)
size_t getSideLabel | ( | size_t | i | ) | const |
Return label of i
-th side.
[in] | i | Side index |
bool NodesAreDOF | ( | ) | const |
Return information about DOF type.
false
otherwise size_t NumberEquations | ( | size_t | dof = 0 | ) |
Renumber Equations.
[in] | dof | Label of degree of freedom for which numbering is performed. Default value (0 ) means that all degrees of freedom are taken into account |
size_t NumberEquations | ( | size_t | dof, |
int | c | ||
) |
Renumber Equations.
[in] | dof | Label of degree of freedom for which numbering is performed. |
[in] | c | code for which degrees of freedom are enforced. |
Operator *=
Rescale mesh coordinates by myltiplying by a factor
[in] | a | Value to multiply by |
void put | ( | const string & | mesh_file | ) | const |
Write mesh data on file.
[in] | mesh_file | Mesh file name |
void RenumberEdge | ( | size_t | n1, |
size_t | n2 | ||
) |
Renumber an edge.
[in] | n1 | Old label |
[in] | n2 | New label |
void RenumberElement | ( | size_t | n1, |
size_t | n2 | ||
) |
Renumber an element.
[in] | n1 | Old label |
[in] | n2 | New label |
void RenumberNode | ( | size_t | n1, |
size_t | n2 | ||
) |
Renumber a node.
[in] | n1 | Old label |
[in] | n2 | New label |
void RenumberSide | ( | size_t | n1, |
size_t | n2 | ||
) |
Renumber a side.
[in] | n1 | Old label |
[in] | n2 | New label |
void Reorder | ( | size_t | m = GRAPH_MEMORY | ) |
Renumber mesh nodes according to reverse Cuthill Mc Kee algorithm.
[in] | m | Memory size needed for matrix graph (default value is GRAPH_MEMORY, see OFELI_Config.h) |
Rescale mesh by multiplying node coordinates by constants.
This function can be used e.g. for changing coordinate units
[in] | sx | Factor to multiply by x coordinates |
[in] | sy | Factor to multiply by y coordinates [Default: sx ] |
[in] | sz | Factor to multiply by z coordinates [Default: sx ] |
void save | ( | const string & | mesh_file | ) | const |
Write mesh data on file in various formats.
File format depends on the extension in file name
[in] | mesh_file | Mesh file name If the extension is '.m', the output file is an OFELI file If the extension is '.gpl', the output file is a Gnuplot file If the extension is '.msh' or '.geo', the output file is a Gmsh file If the extension is '.vtk', the output file is a VTK file |
void set | ( | Node * | nd | ) |
Replace node in the mesh.
If the node label exists already, the existing node pointer will be replaced by the current one. If not, an error message is displayed.
[in] | nd | Pointer to node |
void set | ( | Element * | el | ) |
Replace element in the mesh.
If the element label exists already, the existing element pointer will be replaced by the current one. If not, an error message is displayed.
[in] | el | Pointer to element |
void set | ( | Side * | sd | ) |
Choose side in the mesh.
If the side label exists already, the existing side pointer will be replaced by the current one. If not, an error message is displayed.
[in] | sd | Pointer to side |
void setDim | ( | size_t | dim | ) |
Define space dimension. Normally, between 1 and 3.
[in] | dim | Space dimension to set (must be between 1 and 3) |
void setDOFSupport | ( | int | opt, |
int | nb_nodes = 1 |
||
) |
Define supports of degrees of freedom.
[in] | opt | DOF type:
|
[in] | nb_nodes | Number of nodes on sides or elements (default=1). This parameter is useful only if dofs are supported by sides or elements |
ELEMENT_DOF
or SIDE_DOF
is selected. So it not necessary to call getAllSides() after void setList | ( | const vector< Node *> & | nl | ) |
Initialize list of mesh nodes using the input vector.
[in] | nl | vector instance that contains the list of pointers to nodes |
void setList | ( | const vector< Element *> & | el | ) |
Initialize list of mesh elements using the input vector.
[in] | el | vector instance that contains the list of pointers to elements |
void setList | ( | const vector< Side *> & | sl | ) |
Initialize list of mesh sides using the input vector.
[in] | sl | vector instance that contains the list of pointers to sides |
void setMaterial | ( | int | code, |
const string & | mname | ||
) |
Associate material to code of element.
[in] | code | Element code for which material is assigned |
[in] | mname | Name of material |
void setNbDOFPerNode | ( | size_t | nb_dof = 1 | ) |
Define number of degrees of freedom for each node.
[in] | nb_dof | Number of degrees of freedom (unknowns) for each mesh node (Default value is 1 ) |
Define a point in the domain. This function makes sense only if boundary mesh is given without internal mesh (Case of Boundary Elements)
[in] | x | Coordinates of point to define |
bool SidesAreDOF | ( | ) | const |
Return information about DOF type.
false
otherwise