xc
EightNodeBrick_u_p_U.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 //----------------------------------------------------------------------------
29 //
30 // COPYRIGHT (C): :-))
31 // PROJECT: Object Oriented Finite Element Program
32 // FILE: EightNodeBrick_u_p_U.h
33 // CLASS: EightNodeBrick_u_p_U
34 // MEMBER FUNCTIONS:
35 //
36 // MEMBER VARIABLES
37 //
38 // PURPOSE: Finite Element Class for coupled system
39 // "Coupled system": Solid and fluid coexist.
40 // u-- Solid displacement
41 // p-- Pore pressure
42 // U-- Absolute fluid displacement
43 // RETURN:
44 // VERSION:
45 // LANGUAGE: C++.ver >= 3.0
46 // TARGET OS: DOS || UNIX || . . .
47 // DESIGNER: Boris Jeremic, Xiaoyan Wu, Chao Cheng (main for last revision)
48 // PROGRAMMER: Boris Jeremic, Xiaoyan Wu, Zhaohui Yang, Zhao Cheng (main for last revision)
49 // DATE: Sept. 2001
50 // UPDATE HISTORY: Modified from EightNodeBrick.h. Reorganized a lot by Xiaoyan
51 // 01/24/2002 Xiaoyan
52 // Add the permeability tensor and ks, kf to the constructor Xiaoyan
53 //
54 //
55 //
56 // Clean-up and re-write by Zhao Cheng, 10/20/2004
57 //
58 //
60 
61 
62 #ifndef EIGHTNODEBRICK_U_P_U_H
63 #define EIGHTNODEBRICK_U_P_U_H
64 
65 #include <domain/mesh/element/volumetric/BrickBase.h>
66 #include "domain/mesh/element/utils/body_forces/BodyForces3D.h"
67 
68 namespace XC {
69 class BJtensor;
70 class NDMaterial;
71 
81  {
82  private:
83  static Matrix K;
84  static Matrix C;
85  static Matrix M;
86  static Vector P;
87 
88  static const int Num_IntegrationPts;
89  static const int Num_TotalGaussPts;
90  static const int Num_Nodes;
91  static const int Num_Dim;
92  static const int Num_Dof;
93  static const int Num_ElemDof;
94  static const double pts[2];
95  static const double wts[2];
96  static BJtensor perm;
97 
98  BodyForces3D bf;
99  double poro;
100  double alpha;
101  double rho_s;
102  double rho_f;
103  double ks;
104  double kf;
105  double pressure;
106 
107  Vector *eleQ;
108  mutable Matrix *Ki;
109  public:
110  EightNodeBrick_u_p_U(int element_number,
111  int node_numb_1, int node_numb_2, int node_numb_3, int node_numb_4,
112  int node_numb_5, int node_numb_6, int node_numb_7, int node_numb_8,
113  NDMaterial *Globalmmodel, const BodyForces3D &,
114  double nn, double alf, double rs,double rf,
115  double permb_x, double permb_y, double permb_z,
116  double kks, double kkf, double pp);
117  EightNodeBrick_u_p_U(void);
118  Element *getCopy(void) const;
119  ~EightNodeBrick_u_p_U(void);
120  // public methods to obtain information about dof & connectivity
121  int getNumDOF(void) const;
122  void setDomain(Domain *theDomain);
123 
124  // public methods to set the state of the element
125  int update(void);
126 
127  // public methods to obtain stiffness, mass, damping and residual information
128  const Matrix &getTangentStiff(void) const;
129  const Matrix &getInitialStiff(void) const;
130  const Matrix &getDamp(void) const;
131  const Matrix &getMass(void) const;
132 
133  void alive(void);
134  void zeroLoad(void);
135  int addLoad(ElementalLoad *theLoad, double loadFactor);
136  int addInertiaLoadToUnbalance(const Vector &accel);
137  const Vector &getResistingForce(void) const;
138  const Vector &getResistingForceIncInertia(void) const;
139 
140  int sendSelf(Communicator &);
141  int recvSelf(const Communicator &);
142 
143  void Print(std::ostream &s, int flag =0) const;
144 
145  Response *setResponse(const std::vector<std::string> &argv, Information &eleInfo);
146  int getResponse(int responseID, Information &eleInformation);
147 
148  //int setParameter(const std::vector<std::string> &argv, Parameter &param);
149  //int updateParameter (int parameterID, Information &info);
150 
151  private:
152  BJtensor shapeFunction(double, double, double) const;
153  BJtensor shapeFunctionDerivative(double, double, double) const;
154  BJtensor getGaussPts(void);
155  BJtensor getNodesCrds(void) const;
156  BJtensor getNodesDisp(void) const;
157  BJtensor Jacobian_3D(const BJtensor &dh) const;
158  BJtensor Jacobian_3Dinv(const BJtensor &dh) const;
159  BJtensor dh_Global(const BJtensor &dh) const;
160 
161  BJtensor getStiffnessTensorKep() const;
162  BJtensor getStiffnessTensorG12() const;
163  BJtensor getStiffnessTensorP() const;
164  BJtensor getMassTensorMsf() const;
165  BJtensor getDampTensorC123(void) const;
166  const Matrix &getStiff(int Ki_flag) const;
167  double getPorePressure(double, double, double);
168  Vector getExForceS();
169  Vector getExForceF();
170  };
171 } // end of XC namespace
172 
173 
174 #endif
175 
const Matrix & getTangentStiff(void) const
Return the tangent stiffness matrix.
Definition: EightNodeBrick_u_p_U.cpp:159
Float vector abstraction.
Definition: Vector.h:94
int getResponse(int responseID, Information &eleInformation)
Obtain information from an analysis.
Definition: EightNodeBrick_u_p_U.cpp:440
int sendSelf(Communicator &)
Send the object.
Definition: EightNodeBrick_u_p_U.cpp:399
int update(void)
Updates the element state.
Definition: EightNodeBrick_u_p_U.cpp:523
Element * getCopy(void) const
Virtual constructor.
Definition: EightNodeBrick_u_p_U.cpp:136
Information about an element.
Definition: Information.h:81
Communication parameters between processes.
Definition: Communicator.h:66
const Vector & getResistingForce(void) const
Returns the resisting force vector for the element.
Definition: EightNodeBrick_u_p_U.cpp:324
Base class for hexahedral elements.
Definition: BrickBase.h:69
Base class response objects.
Definition: Response.h:81
void alive(void)
Reactivates the element.
Definition: EightNodeBrick_u_p_U.cpp:252
const Matrix & getDamp(void) const
Returns the damping matrix.
Definition: EightNodeBrick_u_p_U.cpp:179
Boris Jeremic tensor class.
Definition: BJtensor.h:112
Base class for the finite elements.
Definition: Element.h:112
void zeroLoad(void)
Zeroes loads on element.
Definition: EightNodeBrick_u_p_U.cpp:264
Body forces over an element.
Definition: BodyForces3D.h:40
Response * setResponse(const std::vector< std::string > &argv, Information &eleInfo)
setResponse() is a method invoked to determine if the element will respond to a request for a certain...
Definition: EightNodeBrick_u_p_U.cpp:413
void Print(std::ostream &s, int flag=0) const
Print stuff.
Definition: EightNodeBrick_u_p_U.cpp:491
int recvSelf(const Communicator &)
Receive the object.
Definition: EightNodeBrick_u_p_U.cpp:406
Base class for 2D and 3D materials.
Definition: NDMaterial.h:101
const Matrix & getMass(void) const
Returns the mass matrix.
Definition: EightNodeBrick_u_p_U.cpp:223
int addLoad(ElementalLoad *theLoad, double loadFactor)
Adds to the element the load being passed as parameter.
Definition: EightNodeBrick_u_p_U.cpp:272
Eight node hexahedral element for three-dimensional coupled problems.
Definition: EightNodeBrick_u_p_U.h:80
Open source finite element program for structural analysis.
Definition: ContinuaReprComponent.h:35
Matrix of floats.
Definition: Matrix.h:111
void setDomain(Domain *theDomain)
Sets the domain for the element.
Definition: EightNodeBrick_u_p_U.cpp:151
Base class for loads over elements.
Definition: ElementalLoad.h:79
Domain (mesh and boundary conditions) of the finite element model.
Definition: Domain.h:117
const Vector & getResistingForceIncInertia(void) const
Returns the resisting force vector including inertia forces.
Definition: EightNodeBrick_u_p_U.cpp:358
int getNumDOF(void) const
return the number of DOF associated with the element.
Definition: EightNodeBrick_u_p_U.cpp:147