Add A New Element: Difference between revisions

From OpenSeesWiki
Jump to navigation Jump to search
No edit summary
No edit summary
 
Line 10: Line 10:
# [[Add a New Element C]]
# [[Add a New Element C]]
# [[Add a New Element Fortran]]
# [[Add a New Element Fortran]]
To add a new Element into the interpreted applications, the developer must:
# provide a new C++ subclass of the Element class
# provide an interface function that will be used to parse the input and create the new element.
=== Element Class ===
The Element class itself is an abstract base class. It inherits from both the DomainComponent class, which is itself a subclass of TaggedObject class and the MovableObject class. The class has a large number of methods defined in the interface, not all these methods need to be included in a new Element class.
The following is the minimal interface that should be considered:
The Element Class:
<source lang="cpp">
class Element : public DomainComponent
{
  public:
    Element(int tag, int classTag);
    virtual ~Element();
    // methods dealing with nodes and number of external dof                         
    virtual int getNumExternalNodes(void) const =0;
    virtual const ID &getExternalNodes(void)  =0;
    virtual Node **getNodePtrs(void)  =0;
    virtual int getNumDOF(void) =0;
    // methods dealing with committed state and update                               
    virtual int commitState(void);
    virtual int revertToLastCommit(void) = 0;
    virtual int revertToStart(void);
    virtual int update(void);
    // methods dealing with element stiffness                                                     
    virtual const Matrix &getTangentStiff(void) =0;
    virtual const Matrix &getInitialStiff(void) =0;
 
    // methods dealing with element forces
    virtual void zeroLoad(void);
    virtual int addLoad(ElementalLoad *theLoad, double loadFactor);       
    virtual const Vector &getResistingForce(void) =0;
    // public methods for output                                                     
    int sendSelf(int commitTag, Channel &theChannel);
    int recvSelf(int commitTag, Channel &theChannel, FEM_ObjectBroker &theBroker);
    void Print(OPS_Stream &s, int flag =0);
    // method for obtaining information specific to an element                 
    virtual Response *setResponse(const char **argv, int argc, OPS_Stream &theHandler);
    virtual int getResponse(int responseID, Information &eleInformation);
}
</source>
=== Example - Truss2D ===
In the following section we will provide all necessary code to add a new 2d planar truss element into an OpenSees interpreter. To demonstrate the power of object-oriented programming, the stress-strain relationship will be provided by a UniaxialMaterial object.
==== Header ====
The header for thew new class, which we will call Truss2D is as follows:
<source lang="cpp">
// include directives
#include <Element.h>
#include <Matrix.h>
#include <Vector.h>
// forward declarations
class Node;
class UniaxialMaterial;
class Truss2D : public Element
{
  public:
    // constructors                                                                 
    Truss2D(int tag,
            int Nd1, int Nd2,
            UniaxialMaterial &theMaterial,
            double A);
    Truss2D();
    // destructor                                                                   
    ~Truss2D();
    // public methods to obtain inforrmation about dof & connectivity               
    int getNumExternalNodes(void) const;
    const ID &getExternalNodes(void);
    Node **getNodePtrs(void);
    int getNumDOF(void);
    void setDomain(Domain *theDomain);
    // public methods to set the state of the element                               
    int commitState(void);
    int revertToLastCommit(void);
    int revertToStart(void);
    int update(void);
  // public methods to obtain stiffness, mass, damping and residual information   
    const Matrix &getTangentStiff(void);
    const Matrix &getInitialStiff(void);
    const Vector &getResistingForce(void);
    // public methods for output                                                     
    int sendSelf(int commitTag, Channel &theChannel);
    int recvSelf(int commitTag, Channel &theChannel, FEM_ObjectBroker &theBroker);
    void Print(OPS_Stream &s, int flag =0);
    // method for obtaining information specific to an element
    Response *setResponse(const char **argv, int argc, OPS_Stream &s);
    int getResponse(int responseID, Information &eleInformation);
  protected:
  private:
    // private member functions - only available to objects of the class             
    double computeCurrentStrain(void) const;
    // private attributes - a copy for each object of the class                     
    UniaxialMaterial *theMaterial;      // pointer to a material                   
    ID  externalNodes;                  // contains the id's of end nodes           
    Matrix trans;      // hold the transformation matrix                     
    double L;          // length of truss (undeformed configuration)                                                                             
    double A;          // area of truss                                             
    Node *theNodes[2];  // node pointers                                             
    // static data - single copy for all objects of the class                       
    static Matrix trussK;  // class wide matrix for returning stiffness                           
    static Vector trussR;  // class wide vector for returning residual             
};
#endif
</source>
The header file defines the interface and variables for the class Truss2D. It defines the new class to be a sublass of the Element class. In the public interface are 2 constructors and 1 destructor in addition to minimal set of methods we showed for the Element class.
There are no protected data or methods as we do not expect this class to be subclassed. In the private section we define one private method, computeCurrentStrain(), we define a number of private variables, and we also define a number of static variables.
The header has a number of #include directives, one is needed for the base class and every class used as a variable in the list of data (except those that are used as pointers). For those classes that only appear as pointers in the header file (Node, UniaxialMaterial) a forward declaration is all that is needed (the include could also have been used, but using the forward declaration simplifies dependencies and reduces the amount of code that ha to be recompiled later if changes are made).
==== Implementation ====
It another file, Truss2D.cpp, we place the code that details what the constructors, destructor and methods do. In addition we provide one additional procedure OPS_Truss2D() (NOTE it has the same name as the class with an OPS_ prefix). We will go through each part of the file.
===== Include Directives =====
The first part of the file contains the list of includes. It is necessary to have an #include directive for each class and api file that is used within the .cpp file and is not included in the header.
<source lang="cpp">
#include "Truss2D.h"
#include <elementAPI.h>
#include <G3Globals.h>
#include <Matrix.h>
#include <Vector.h>
#include <Information.h>
#include <Domain.h>
#include <Node.h>
#include <Channel.h>
#include <Message.h>
#include <FEM_ObjectBroker.h>
#include <UniaxialMaterial.h>
#include <Renderer.h>
#include <ElementResponse.h>
#include <math.h>
#include <stdlib.h>
#include <string.h>
</source>
===== Static Variables =====
Next, we initialize the static variables. For this example we are using 2 static-variables (objects shared by each Truss2D object that is created), one to return the tangent matrix and the other to return the resisting force.
<source lang="cpp">
// initialise the class wide variables                                               
Matrix Truss2D::trussK(4,4);
Vector Truss2D::trussR(4);
</source>
===== Constructors =====
After the list of includes, we provide the 2 constructors. The constructors are rather simple. They just initialize all the data variables defined in the header. Note it is very important to set all pointer values to 0.
The first constructor is the one most typically used. The arguments provide the elements tag, the tags of the two end nodes, the element's area and a copy of the element's material. The code in the constructor does the following:
# The elements tag and a 0 are passed to the Element constructor.
# The matreial pointer, theMaterial, is set to a copy of the material obtained from the material that is passed in the arguments.
# The externalNodes array is set to be an array of size 2 and it's values are set to the nodal tags of the 2 nodes.
# The theNodes array components are set to be 0.
<source lang="cpp">                                                             
TrussCPP::TrussCPP(int tag,
                int Nd1, int Nd2,
                UniaxialMaterial &theMat,
                double a)
:Element(tag, 0),
externalNodes(2),
trans(1,4), L(0.0), A(a)
{
  // get a copy of the material object for our own use                               
  theMaterial = theMat.getCopy();
  if (theMaterial == 0) {
  opserr << "FATAL TrussCPP::TrussCPP() - out of memory, could not get a copy of the Material\n";
  exit(-1);
  }
  // fill in the ID containing external node info with node id's                     
  if (externalNodes.Size() != 2) {
    opserr << "FATAL TrussCPP::TrussCPP() - out of memory, could not create an ID of size 2\n";
    exit(-1);
  }
  externalNodes(0) = Nd1;
  externalNodes(1) = Nd2;
  theNodes[0] = 0;
  theNodes[1] = 0;
}
</source>
The second constructor is called when paralell processing or the database feature of the OpenSees application is used. It's pupose is to create a blank Truss2D object, that will be filled in when the recvSelf() method is invoked on the object.
<source lang="cpp">
TrussCPP::TrussCPP()
:Element(0, 0),
theMaterial(0),
externalNodes(2),
trans(1,4), L(0.0), A(0.0)
{
  theNodes[0] = 0;
  theNodes[1] = 0;
}
</source>
===== Destructor =====
The we provide the destructor. In the destructor all memory that the Truss2D created or was passed to it in the constructor must be destroyed.
For our example, we need to invoke the destructor on the copy of the material object.
<source lang="cpp"> 
TrussCPP::~TrussCPP()
{                                     
    if (theMaterial != 0)
        delete theMaterial;
}
</source>
===== Methods Dealing With Nodes =====
Next comes 4 rather simple methods that return basic information about the elements nodes. These are one line methods that should not need any explanation!
<source lang="cpp">
int
TrussCPP::getNumExternalNodes(void) const
{
    return 2;
}
const ID &
TrussCPP::getExternalNodes(void)
{
  return externalNodes;
}
Node **
TrussCPP::getNodePtrs(void)
{
  return theNodes;
}
int
TrussCPP::getNumDOF(void) {
    return 4;
}
</source>
===== Methods Dealing With Current State =====
<source lang="cpp">
int
TrussCPP::commitState()
{
    return theMaterial->commitState();
}
int
TrussCPP::revertToLastCommit()
{
    return theMaterial->revertToLastCommit();
}
int
TrussCPP::revertToStart()
{
    return theMaterial->revertToStart();
}
int
TrussCPP::update()
{
  // determine the current strain given trial displacements at nodes                 
  double strain = this->computeCurrentStrain();
  // set the strain in the materials                                                 
  theMaterial->setTrialStrain(strain);
  return 0;
}
</source>
===== Methods To Return Tangent Matrix =====
In both methods, we obtain the appropriate tangent from the material and use this to return the transformed matrix.
<source lang="cpp">
const Matrix &
TrussCPP::getTangentStiff(void)
{
    if (L == 0.0) { // length = zero - problem in setDomain() warning message already printed                                         
        trussK.Zero();
        return trussK;
    }
    // get the current E from the material for the last updated strain                                                               
    double E = theMaterial->getTangent();
    // form the tangent stiffness matrix                                                                                             
    trussK = trans^trans;
    trussK *= A*E/L;
    // return the matrix                                                                                                             
    return trussK;
}
const Matrix &
TrussCPP::getInitialStiff(void)
{
    if (L == 0.0) { // length = zero - problem in setDomain() warning message already printed                                         
        trussK.Zero();
        return trussK;
    }
    // get the current E from the material for the last updated strain                                                               
    double E = theMaterial->getInitialTangent();
    // form the tangent stiffness matrix                                                                                             
    trussK = trans^trans;
    trussK *= A*E/L;
    // return the matrix                                                                                                             
    return trussK;
}
</source>
===== Methods To Return Resisting Force =====
In this method we obtain the stress from the material and use this to return the transformed force vector.
<source lang="cpp">
const Vector &
TrussCPP::getResistingForce()
{
    if (L == 0.0) { // if length == 0, problem in setDomain()                                                                         
        trussR.Zero();
        return trussR;
    }
    // want: R = Ku - Pext                                                                                                           
    // force = F * transformation                                                                                                     
    double force = A*theMaterial->getStress();
    for (int i=0; i<4; i++)
        trussR(i) = trans(0,i)*force;
    return trussR;
}
</source>

Latest revision as of 23:04, 13 July 2010

The OpenSees applications allows developers to use their own element modules to the application. Unlike most other programs, the elements are added at run-time and not at compile time. The advantage of this is the the developers:

  1. Do not need the OpenSees source files or libraries to compile and link the application.
  2. Can share their modules with others without having to provide the source code.

The element modules can be written using either C++, C, or Fortran. Whatever the language the developer wishes to use, the element modules make use of the OpenSees API to for example find nodal coordinates and displacements, and make use of material modules that exist in OpenSees or are also added as external routines.

  1. Add a New Element C++
  2. Add a New Element C
  3. Add a New Element Fortran