OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
intelline1phf.C
Go to the documentation of this file.
1 /*
2  *
3  * ##### ##### ###### ###### ### ###
4  * ## ## ## ## ## ## ## ### ##
5  * ## ## ## ## #### #### ## # ##
6  * ## ## ## ## ## ## ## ##
7  * ## ## ## ## ## ## ## ##
8  * ##### ##### ## ###### ## ##
9  *
10  *
11  * OOFEM : Object Oriented Finite Element Code
12  *
13  * Copyright (C) 1993 - 2013 Borek Patzak
14  *
15  *
16  *
17  * Czech Technical University, Faculty of Civil Engineering,
18  * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19  *
20  * This library is free software; you can redistribute it and/or
21  * modify it under the terms of the GNU Lesser General Public
22  * License as published by the Free Software Foundation; either
23  * version 2.1 of the License, or (at your option) any later version.
24  *
25  * This program is distributed in the hope that it will be useful,
26  * but WITHOUT ANY WARRANTY; without even the implied warranty of
27  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28  * Lesser General Public License for more details.
29  *
30  * You should have received a copy of the GNU Lesser General Public
31  * License along with this library; if not, write to the Free Software
32  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33  */
34 
38 #include "node.h"
39 #include "gausspoint.h"
40 #include "gaussintegrationrule.h"
41 #include "floatmatrix.h"
42 #include "floatarray.h"
43 #include "intarray.h"
44 #include "mathfem.h"
45 #include "fei2dlinelin.h"
46 #include "classfactory.h"
47 
48 
49 namespace oofem {
50 REGISTER_Element(IntElLine1PhF);
51 
52 FEI2dLineLin IntElLine1PhF :: interp(1, 1);
53 
54 
56 {
57  numberOfDofMans = 4;
58  axisymmode = false;
59 
61 }
62 
63 
64 void
66 {
67  // Returns the modified N-matrix which multiplied with u give the spatial jump.
68 
69  FloatArray N;
71  interp->evalN( N, ip->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
72 
73  answer.resize(2, 8);
74  answer.zero();
75  answer.at(1, 1) = answer.at(2, 2) = -N.at(1);
76  answer.at(1, 3) = answer.at(2, 4) = -N.at(2);
77 
78  answer.at(1, 5) = answer.at(2, 6) = N.at(1);
79  answer.at(1, 7) = answer.at(2, 8) = N.at(2);
80 }
81 
82 
83 void
85 // Sets up the array of Gauss Points of the receiver.
86 {
87  if ( integrationRulesArray.size() == 0 ) {
88  integrationRulesArray.resize( 1 );
89  //integrationRulesArray[ 0 ].reset( new LobattoIntegrationRule (1,domain, 1, 2) );
90  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 2) );
91  integrationRulesArray [ 0 ]->SetUpPointsOnLine(this->numberOfGaussPoints, _2dInterface);
92  }
93 }
94 
95 void
97 {
98  FloatMatrix dNdxi;
100  interp->evaldNdxi( dNdxi, ip->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) );
101  G.resize(2);
102  G.zero();
103  int numNodes = this->giveNumberOfNodes();
104  for ( int i = 1; i <= dNdxi.giveNumberOfRows(); i++ ) {
105  double X1_i = 0.5 * ( this->giveNode(i)->giveCoordinate(1) + this->giveNode(i + numNodes / 2)->giveCoordinate(1) ); // (mean) point on the fictious mid surface
106  double X2_i = 0.5 * ( this->giveNode(i)->giveCoordinate(2) + this->giveNode(i + numNodes / 2)->giveCoordinate(2) );
107  G.at(1) += dNdxi.at(i, 1) * X1_i;
108  G.at(2) += dNdxi.at(i, 1) * X2_i;
109  }
110 }
111 
112 double
114 {
115  FloatArray G;
116  this->computeCovarBaseVectorAt(ip, G);
117 
118  double weight = ip->giveWeight();
119  double ds = sqrt( G.dotProduct(G) ) * weight;
120  if ( this->axisymmode ) {
121  int numNodes = this->giveNumberOfNodes();
122  FloatArray N;
124  // interpolate radius
125  double r = 0.0;
126  for ( int i = 1; i <= N.giveSize(); i++ ) {
127  double X_i = 0.5 * ( this->giveNode(i)->giveCoordinate(1) + this->giveNode(i + numNodes / 2)->giveCoordinate(1) ); // X-coord of the fictious mid surface
128  r += N.at(i) * X_i;
129  }
130  return ds * r;
131 
132  } else { // regular 2d
133  double thickness = this->giveCrossSection()->give(CS_Thickness, ip);
134  return ds * thickness;
135  }
136 }
137 
138 
141 {
142  this->axisymmode = false;
145 
146 }
147 
148 
149 void
151 {
152  answer = {D_u, D_v, T_f};
153 }
154 
155 void
157 {
158  // Transformation matrix to the local coordinate system
159  FloatArray G;
160  this->computeCovarBaseVectorAt(gp, G);
161  G.normalize();
162 
163  answer.resize(2, 2);
164  answer.at(1, 1) = G.at(1);
165  answer.at(2, 1) = -G.at(2);
166  answer.at(1, 2) = G.at(2);
167  answer.at(2, 2) = G.at(1);
168 
169 }
170 
173 {
174  return & interp;
175 }
176 
177 
178 void
180 {
181  answer = {D_u, D_v};
182 }
183 
184 void
186 {
187  // use temperature for now
188  answer = {T_f};
189 }
190 
191 void
193 {
194  answer = {1, 2, 4, 5, 7, 8, 10, 11};
195 }
196 void
198 {
199  answer = {3, 6, 9, 12};
200 }
201 
202 
203 void
204 IntElLine1PhF :: giveEngTraction(FloatArray &answer, GaussPoint *gp, const FloatArray &jump, const double damage, TimeStep *tStep)
205 {
206 
208  //double damage = this->computeDamageAt(gp, VM_Total, tStep);
209 
210  //StructuralInterfaceMaterial *mat = this->giveInterfaceCrossSection()->giveInterfaceMaterial();
211  mat->giveEngTraction_2d(answer, gp, jump, damage, tStep);
212  //this->giveInterfaceCrossSection()->giveEngTraction_2d(answer, gp, jump, tStep);
213 }
214 
215 
216 void
218 {
219 
220  FloatArray G1;
221  this->computeCovarBaseVectorAt(gp, G1);
222 
223  G.resize(2, 1);
224  G.at(1, 1) = G1.at(1);
225  G.at(2, 1) = G1.at(2);
226 
227 
228 }
229 
230 
231 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
StructuralInterfaceMaterial * giveInterfaceMaterial()
virtual void evaldNdxi(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
Definition: feinterpol.h:193
virtual void giveEngTraction(FloatArray &answer, GaussPoint *gp, const FloatArray &jump, const double damage, TimeStep *tStep)
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the array of interpolation functions (shape functions) at given point.
Class and object Domain.
Definition: domain.h:115
virtual void computeCovarBaseVectorsAt(GaussPoint *gp, FloatMatrix &G)
virtual void getLocationArray_d(IntArray &answer)
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
virtual FEInterpolation * giveInterpolation() const
static FEI2dLineLin interp
Definition: intelline1phf.h:53
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dlinelin.C:42
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
virtual void giveDofManDofIDMask_d(IntArray &answer)
virtual void computeTransformationMatrixAt(GaussPoint *gp, FloatMatrix &answer)
virtual double giveCoordinate(int i)
Definition: node.C:82
virtual void giveDofManDofIDMask_u(IntArray &answer)
Class implementing an array of integers.
Definition: intarray.h:61
virtual void computeNmatrixAt(GaussPoint *gp, FloatMatrix &answer)
Computes modified interpolation matrix (N) for the element which multiplied with the unknowns vector ...
Definition: intelline1phf.C:65
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: intelline1phf.C:84
virtual int giveNumberOfNodes() const
Returns number of nodes of receiver.
Definition: element.h:662
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
double dotProduct(const FloatArray &x) const
Computes the dot product (or inner product) of receiver and argument.
Definition: floatarray.C:463
REGISTER_Element(LSpace)
virtual void getLocationArray_u(IntArray &answer)
virtual void computeCovarBaseVectorAt(GaussPoint *gp, FloatArray &G)
Definition: intelline1phf.C:96
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
#define N(p, q)
Definition: mdm.C:367
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
Class representing vector of real numbers.
Definition: floatarray.h:82
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
IRResultType
Type defining the return values of InputRecord reading operations.
Definition: irresulttype.h:47
virtual double give(CrossSectionProperty a, GaussPoint *gp)
Returns the value of cross section property at given point.
Definition: crosssection.C:151
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
IntElLine1PhF(int n, Domain *d)
Definition: intelline1phf.C:55
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
Class representing the general Input Record.
Definition: inputrecord.h:101
#define _IFT_IntElLine1PhF_axisymmode
Definition: intelline1phf.h:41
bool axisymmode
Flag controlling axisymmetric mode (integration over unit circumferential angle)
Definition: intelline1phf.h:55
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual double computeAreaAround(GaussPoint *gp)
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
std::vector< std::unique_ptr< IntegrationRule > > integrationRulesArray
List of integration rules of receiver (each integration rule contains associated integration points a...
Definition: element.h:170
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
the oofem namespace is to define a context or scope in which all oofem names are defined.
virtual void giveEngTraction_2d(FloatArray &answer, GaussPoint *gp, const FloatArray &jump, const double damage, TimeStep *tStep)
double normalize()
Normalizes receiver.
Definition: floatarray.C:828
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
int giveNumberOfRows() const
Returns number of rows of receiver.
Definition: floatmatrix.h:156
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Interface element class with phase field (PhF) modeling of damage.
Class representing solution step.
Definition: timestep.h:80
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
StructuralInterfaceCrossSection * giveInterfaceCrossSection()
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:29 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011