xc
TotalLagrangianFD20NodeBrick.h
1 // -*-c++-*-
2 //----------------------------------------------------------------------------
3 // XC program; finite element analysis code
4 // for structural analysis and design.
5 //
6 // Copyright (C) Luis C. Pérez Tato
7 //
8 // This program derives from OpenSees <http://opensees.berkeley.edu>
9 // developed by the «Pacific earthquake engineering research center».
10 //
11 // Except for the restrictions that may arise from the copyright
12 // of the original program (see copyright_opensees.txt)
13 // XC is free software: you can redistribute it and/or modify
14 // it under the terms of the GNU General Public License as published by
15 // the Free Software Foundation, either version 3 of the License, or
16 // (at your option) any later version.
17 //
18 // This software is distributed in the hope that it will be useful, but
19 // WITHOUT ANY WARRANTY; without even the implied warranty of
20 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 // GNU General Public License for more details.
22 //
23 //
24 // You should have received a copy of the GNU General Public License
25 // along with this program.
26 // If not, see <http://www.gnu.org/licenses/>.
27 //----------------------------------------------------------------------------
28 //===============================================================================
29 //# COPYRIGHT (C): Woody's license (by BJ):
30 // ``This source code is Copyrighted in
31 // U.S., for an indefinite period, and anybody
32 // caught using it without our permission, will be
33 // mighty good friends of ourn, cause we don't give
34 // a darn. Hack it. Compile it. Debug it. Run it.
35 // Yodel it. Enjoy it. We wrote it, that's all we
36 // wanted to do.''
37 //
38 //# PROJECT: Object Oriented Finite Element Program
39 //# PURPOSE: Finite Deformation Hyper-Elastic classes
40 //# CLASS:
41 //#
42 //# VERSION: 0.6_(1803398874989) (golden section)
43 //# LANGUAGE: C++
44 //# TARGET OS: all...
45 //# DESIGN: Zhao Cheng, Boris Jeremic (jeremic@ucdavis.edu)
46 //# PROGRAMMER(S): Zhao Cheng, Boris Jeremic
47 //#
48 //#
49 //# DATE: Sept2003
50 //# UPDATE HISTORY: 28May2004, Zhao Optimized the Stiffness Tensor
51 //#
52 //#
53 //===============================================================================
54 #ifndef TOTALLAGRANGIANFD20BRICK_H
55 #define TOTALLAGRANGIANFD20BRICK_H
56 
57 #include <domain/mesh/element/ElemWithMaterial.h>
58 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h"
59 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
60 
61 namespace XC {
62 class BJtensor;
63 class NDMaterial;
64 
71 class TotalLagrangianFD20NodeBrick: public ElemWithMaterial<20,NDMaterialPhysicalProperties>
72  {
73  private:
74  static Matrix K;
75 // static Matrix C; //!< Element damping matrix
76  static Matrix M;
77  static Vector P;
78  static const double pts[3];
79  static const double wts[3];
80  BodyForces3D bf;
81 
82  double det_of_Jacobian;
83 
84  mutable Matrix *Ki;
85 
86  static const int NumIntegrationPts;
87  static const int NumTotalGaussPts;
88  static const int NumNodes;
89  static const int NumDof;
90  static const int NumElemDof;
91 
92  static BJtensor shapeFunction(double , double , double );
93  static BJtensor shapeFunctionDerivative(double , double , double );
94 
95  BJtensor Jacobian_3D(const BJtensor &dh) const;
96  BJtensor Jacobian_3Dinv(const BJtensor &dh) const;
97  BJtensor dh_Global(const BJtensor &dh) const;
98  BJtensor getNodesCrds(void) const;
99  BJtensor getNodesDisp(void) const;
100 
101  BJtensor getStiffnessTensor(void) const;
102  BJtensor getRtensor(void) const;
103  BJtensor getBodyForce(void) const;
104  BJtensor getSurfaceForce(void) const;
105  BJtensor getForces(void) const;
106  public:
108  int node_numb_1, int node_numb_2, int node_numb_3, int node_numb_4,
109  int node_numb_5, int node_numb_6, int node_numb_7, int node_numb_8,
110  int node_numb_9, int node_numb_10, int node_numb_11, int node_numb_12,
111  int node_numb_13, int node_numb_14, int node_numb_15, int node_numb_16,
112  int node_numb_17, int node_numb_18, int node_numb_19, int node_numb_20,
113  NDMaterial &m, const BodyForces3D &bForces= BodyForces3D());
115  Element *getCopy(void) const;
117 
118  int getNumDOF(void) const;
119  void setDomain(Domain *theDomain);
120 
121  int update();
122 
123  const Matrix &getTangentStiff(void) const;
124  const Matrix &getInitialStiff(void) const;
125  const Matrix &getMass(void) const;
126 
127  void alive(void);
128  int addLoad(ElementalLoad *theLoad, double loadFactor);
129  int addInertiaLoadToUnbalance(const Vector &accel);
130 
131  const Vector &getResistingForce(void) const;
132  const Vector &getResistingForceIncInertia(void) const;
133 
134  virtual int sendSelf(Communicator &);
135  virtual int recvSelf(const Communicator &);
136 
137  void Print(std::ostream &s, int flag =0) const;
138 
139  Response *setResponse(const std::vector<std::string> &argv, Information &eleInformation);
140  int getResponse(int responseID, Information &eleInformation);
141 
142 // int setParameter(const std::vector<std::string> &argv, Parameter &param);
143 // int updateParameter(int parameterID, Information &info);
144  };
145 } // end of XC namespace
146 
147 
148 #endif
149 
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TotalLagrangianFD20NodeBrick.cpp:470
Float vector abstraction.
Definition: Vector.h:94
Total lagrangian formulation twenty node hexahedral element for three-dimensional problems...
Definition: TotalLagrangianFD20NodeBrick.h:71
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
Base class response objects.
Definition: Response.h:81
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: TotalLagrangianFD20NodeBrick.cpp:429
Boris Jeremic tensor class.
Definition: BJtensor.h:112
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: TotalLagrangianFD20NodeBrick.cpp:123
Response * setResponse(const std::vector< std::string > &argv, Information &eleInformation)
setResponse() is a method invoked to determine if the element will respond to a request for a certain...
Definition: TotalLagrangianFD20NodeBrick.cpp:659
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: TotalLagrangianFD20NodeBrick.cpp:564
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: TotalLagrangianFD20NodeBrick.cpp:127
Base class for the finite elements.
Definition: Element.h:112
virtual int recvSelf(const Communicator &)
Receive the object.
Definition: TotalLagrangianFD20NodeBrick.cpp:622
Body forces over an element.
Definition: BodyForces3D.h:40
Element * getCopy(void) const
Virtual constructor.
Definition: TotalLagrangianFD20NodeBrick.cpp:111
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: TotalLagrangianFD20NodeBrick.cpp:630
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
void alive(void)
Reactivates the element.
Definition: TotalLagrangianFD20NodeBrick.cpp:512
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: TotalLagrangianFD20NodeBrick.cpp:584
int update()
Updates the element state.
Definition: TotalLagrangianFD20NodeBrick.cpp:135
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Element with material.
Definition: ElemWithMaterial.h:45
Matrix of floats.
Definition: Matrix.h:111
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: TotalLagrangianFD20NodeBrick.cpp:685
virtual int sendSelf(Communicator &)
Send the object.
Definition: TotalLagrangianFD20NodeBrick.cpp:615