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