XC Open source finite element analysis program
EightNodeBrick_u_p_U.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 //----------------------------------------------------------------------------
28 //
29 // COPYRIGHT (C): :-))
30 // PROJECT: Object Oriented Finite Element Program
31 // FILE: EightNodeBrick_u_p_U.h
32 // CLASS: EightNodeBrick_u_p_U
33 // MEMBER FUNCTIONS:
34 //
35 // MEMBER VARIABLES
36 //
37 // PURPOSE: Finite Element Class for coupled system
38 // "Coupled system": Solid and fluid coexist.
39 // u-- Solid displacement
40 // p-- Pore pressure
41 // U-- Absolute fluid displacement
42 // RETURN:
43 // VERSION:
44 // LANGUAGE: C++.ver >= 3.0
45 // TARGET OS: DOS || UNIX || . . .
46 // DESIGNER: Boris Jeremic, Xiaoyan Wu, Chao Cheng (main for last revision)
47 // PROGRAMMER: Boris Jeremic, Xiaoyan Wu, Zhaohui Yang, Zhao Cheng (main for last revision)
48 // DATE: Sept. 2001
49 // UPDATE HISTORY: Modified from EightNodeBrick.h. Reorganized a lot by Xiaoyan
50 // 01/24/2002 Xiaoyan
51 // Add the permeability tensor and ks, kf to the constructor Xiaoyan
52 //
53 //
54 //
55 // Clean-up and re-write by Zhao Cheng, 10/20/2004
56 //
57 //
59 
60 
61 #ifndef EIGHTNODEBRICK_U_P_U_H
62 #define EIGHTNODEBRICK_U_P_U_H
63 
64 #include <domain/mesh/element/volumen/BrickBase.h>
65 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
66 
67 namespace XC {
68  class BJtensor;
69  class NDMaterial;
70 
72 //
75  {
76  private:
77  static Matrix K;
78  static Matrix C;
79  static Matrix M;
80  static Vector P;
81 
82  static const int Num_IntegrationPts;
83  static const int Num_TotalGaussPts;
84  static const int Num_Nodes;
85  static const int Num_Dim;
86  static const int Num_Dof;
87  static const int Num_ElemDof;
88  static const double pts[2];
89  static const double wts[2];
90  static BJtensor perm;
91 
92  BodyForces3D bf;
93  double poro;
94  double alpha;
95  double rho_s;
96  double rho_f;
97  double ks;
98  double kf;
99  double pressure;
100 
101  Vector *eleQ;
102  mutable Matrix *Ki;
103  public:
104  EightNodeBrick_u_p_U(int element_number,
105  int node_numb_1, int node_numb_2, int node_numb_3, int node_numb_4,
106  int node_numb_5, int node_numb_6, int node_numb_7, int node_numb_8,
107  NDMaterial *Globalmmodel, const BodyForces3D &,
108  double nn, double alf, double rs,double rf,
109  double permb_x, double permb_y, double permb_z,
110  double kks, double kkf, double pp);
111  EightNodeBrick_u_p_U(void);
112  Element *getCopy(void) const;
113  ~EightNodeBrick_u_p_U(void);
114  // public methods to obtain information about dof & connectivity
115  int getNumDOF(void) const;
116  void setDomain(Domain *theDomain);
117 
118  // public methods to set the state of the element
119  int update(void);
120 
121  // public methods to obtain stiffness, mass, damping and residual information
122  const Matrix &getTangentStiff(void) const;
123  const Matrix &getInitialStiff(void) const;
124  const Matrix &getDamp(void) const;
125  const Matrix &getMass(void) const;
126 
127  void zeroLoad(void);
128  int addLoad(ElementalLoad *theLoad, double loadFactor);
129  int addInertiaLoadToUnbalance(const Vector &accel);
130  const Vector &getResistingForce(void) const;
131  const Vector &getResistingForceIncInertia(void) const;
132 
133  int sendSelf(CommParameters &);
134  int recvSelf(const CommParameters &);
135 
136  void Print(std::ostream &s, int flag =0);
137 
138  Response *setResponse(const std::vector<std::string> &argv, Information &eleInfo);
139  int getResponse(int responseID, Information &eleInformation);
140 
141  //int setParameter(const std::vector<std::string> &argv, Parameter &param);
142  //int updateParameter (int parameterID, Information &info);
143 
144  private:
145  BJtensor shapeFunction(double, double, double) const;
146  BJtensor shapeFunctionDerivative(double, double, double) const;
147  BJtensor getGaussPts(void);
148  BJtensor getNodesCrds(void) const;
149  BJtensor getNodesDisp(void) const;
150  BJtensor Jacobian_3D(BJtensor dh) const;
151  BJtensor Jacobian_3Dinv(BJtensor dh) const;
152  BJtensor dh_Global(BJtensor dh) const;
153 
154  BJtensor getStiffnessTensorKep() const;
155  BJtensor getStiffnessTensorG12() const;
156  BJtensor getStiffnessTensorP() const;
157  BJtensor getMassTensorMsf() const;
158  BJtensor getDampTensorC123(void) const;
159  const Matrix &getStiff(int Ki_flag) const;
160  double getPorePressure(double, double, double);
161  Vector getExForceS();
162  Vector getExForceF();
163  };
164 } // end of XC namespace
165 
166 
167 #endif
168 
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: EightNodeBrick_u_p_U.cpp:346
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
Element * getCopy(void) const
Virtual constructor.
Definition: EightNodeBrick_u_p_U.cpp:136
Hexaedro de ocho nodos.
Definition: EightNodeBrick_u_p_U.h:74
Information about an element.
Definition: Information.h:80
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: EightNodeBrick_u_p_U.cpp:151
const Matrix & getDamp(void) const
Returns the matriz de amortiguamiento.
Definition: EightNodeBrick_u_p_U.cpp:179
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
int update(void)
Actualiza el estado of the element.
Definition: EightNodeBrick_u_p_U.cpp:511
void zeroLoad(void)
Anula el load vector aplicadas of the element.
Definition: EightNodeBrick_u_p_U.cpp:252
void Print(std::ostream &s, int flag=0)
Imprime el objeto.
Definition: EightNodeBrick_u_p_U.cpp:479
Definition: Matrix.h:82
Body forces over an element.
Definition: BodyForces3D.h:39
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: EightNodeBrick_u_p_U.cpp:223
Communication parameters between processes.
Definition: CommParameters.h:65
Base class for hexahedra.
Definition: BrickBase.h:44
================================================================================
Definition: ContinuaReprComponent.h:34
Definition: Response.h:71