XC Open source finite element analysis program
TwentyNodeBrick_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: TwentyNodeBrick_u_p_U.h
32 // CLASS: TwentyNodeBrick_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
47 // PROGRAMMER: Boris Jeremic, Xiaoyan Wu
48 // DATE: Aug. 2001
49 // UPDATE HISTORY: Modified from EightNodeBrick_u_p_U.h. Aug. 2001
50 // 01/16/2002 Xiaoyan
51 // Add the permeability tensor and ks, kf to the constructor Xiaoyan
52 //
53 // Clean-up and re-write by Zhao Cheng, 10/20/2004
54 //
56 
57 #ifndef TWENTYNODEBRICK_U_P_U_H
58 #define TWENTYNODEBRICK_U_P_U_H
59 
60 
61 
62 #include <domain/mesh/element/ElemWithMaterial.h>
63 #include "domain/mesh/element/utils/physical_properties/NDMaterialPhysicalProperties.h"
64 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
65 
66 namespace XC {
67  class BJtensor;
68  class NDMaterial;
69 
71 //
73 class TwentyNodeBrick_u_p_U: public ElemWithMaterial<20,NDMaterialPhysicalProperties>
74  {
75  private:
76  static Matrix K; // Stiffness
77  static Matrix C; // Damping
78  static Matrix M; // Mass
79  static Vector P; // Residual
80 
81  static const int Num_IntegrationPts;
82  static const int Num_TotalGaussPts;
83  static const int Num_Nodes;
84  static const int Num_Dim;
85  static const int Num_Dof;
86  static const int Num_ElemDof;
87  static const double pts[3]; // Stores quadrature points
88  static const double wts[3]; // Stores quadrature weights
89  static BJtensor perm; // Permeability = k/(rho_f*g)
90 
91  BodyForces3D bf;
92  double poro; // Porosity
93  double alpha; // Coefficient for soil (approximate equal 1)
94  double rho_s; // Solid density
95  double rho_f; // Fluid density
96  double ks; // Bulk modulus of solid
97  double kf; // Bulk modulus of fluid
98  double pressure; // Normal surface traction (pressure) over entire element //?
99 
100  Vector *eleQ;
101  mutable Matrix *Ki;
102  private:
103  static BJtensor shapeFunction(double, double, double);
104  static BJtensor shapeFunctionDerivative(double, double, double);
105  BJtensor getGaussPts(void);
106  BJtensor getNodesCrds(void) const;
107  BJtensor getNodesDisp(void);
108  BJtensor Jacobian_3D(BJtensor dh) const;
109  BJtensor Jacobian_3Dinv(BJtensor dh) const;
110  BJtensor dh_Global(BJtensor dh) const;
111 
112  BJtensor getStiffnessTensorKep(void) const;
113  BJtensor getStiffnessTensorG12(void) const;
114  BJtensor getStiffnessTensorP(void) const;
115  BJtensor getMassTensorMsf(void) const;
116  BJtensor getDampTensorC123(void) const;
117  const Matrix &getStiff(int Ki_flag) const;
118  double getPorePressure(double, double, double) const;
119  Vector getExForceS();
120  Vector getExForceF();
121 
122  public:
123  TwentyNodeBrick_u_p_U(int element_number,
124  int node_numb_1, int node_numb_2, int node_numb_3, int node_numb_4,
125  int node_numb_5, int node_numb_6, int node_numb_7, int node_numb_8,
126  int node_numb_9, int node_numb_10, int node_numb_11, int node_numb_12,
127  int node_numb_13, int node_numb_14, int node_numb_15, int node_numb_16,
128  int node_numb_17, int node_numb_18, int node_numb_19, int node_numb_20,
129  NDMaterial *Globalmmodel, const BodyForces3D &,
130  double nn, double alf, double rs,double rf,
131  double permb_x, double permb_y, double permb_z,
132  double kks, double kkf, double pp);
133  TwentyNodeBrick_u_p_U(void);
134  Element *getCopy(void) const;
135  ~TwentyNodeBrick_u_p_U(void);
136 
137  // public methods to obtain information about dof & connectivity
138  int getNumDOF(void) const;
139  void setDomain(Domain *theDomain);
140 
141  // public methods to set the state of the element
142  int update(void);
143 
144  // public methods to obtain stiffness, mass, damping and residual information
145  const Matrix &getTangentStiff(void) const;
146  const Matrix &getInitialStiff(void) const;
147  const Matrix &getDamp(void) const;
148  const Matrix &getMass(void) const;
149 
150  void zeroLoad(void);
151  int addLoad(ElementalLoad *theLoad, double loadFactor);
152  int addInertiaLoadToUnbalance(const Vector &accel);
153  const Vector &getResistingForce(void) const;
154  const Vector &getResistingForceIncInertia(void) const;
155 
156  int sendSelf(CommParameters &);
157  int recvSelf(const CommParameters &);
158 
159  void Print(std::ostream &s, int flag =0);
160 
161  Response *setResponse(const std::vector<std::string> &argv, Information &eleInfo);
162  int getResponse(int responseID, Information &eleInformation);
163 
164  //int setParameter(const std::vector<std::string> &argv, Parameter &param);
165  //int updateParameter (int parameterID, Information &info);
166 };
167 } // end of XC namespace
168 
169 
170 #endif
171 
Element with material.
Definition: ElemWithMaterial.h:40
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
int update(void)
Actualiza el estado of the element.
Definition: TwentyNodeBrick_u_p_U.cpp:535
Information about an element.
Definition: Information.h:80
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: TwentyNodeBrick_u_p_U.cpp:165
Base calass for the finite elements.
Definition: Element.h:104
Base class for loads over elements.
Definition: ElementalLoad.h:73
Hexaedro de veinte nodos.
Definition: TwentyNodeBrick_u_p_U.h:73
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: TwentyNodeBrick_u_p_U.cpp:238
const Matrix & getDamp(void) const
Returns the matriz de amortiguamiento.
Definition: TwentyNodeBrick_u_p_U.cpp:193
Definition: Matrix.h:82
Body forces over an element.
Definition: BodyForces3D.h:39
const Vector & getResistingForceIncInertia(void) const
Returns the action of the element over its attached nodes. Computes damping matrix.
Definition: TwentyNodeBrick_u_p_U.cpp:358
void Print(std::ostream &s, int flag=0)
Imprime el objeto.
Definition: TwentyNodeBrick_u_p_U.cpp:491
Communication parameters between processes.
Definition: CommParameters.h:65
================================================================================
Definition: ContinuaReprComponent.h:34
void zeroLoad(void)
Anula el load vector aplicadas of the element.
Definition: TwentyNodeBrick_u_p_U.cpp:265
Definition: Response.h:71
Element * getCopy(void) const
Virtual constructor.
Definition: TwentyNodeBrick_u_p_U.cpp:150