OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr1darcy.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 
35 
36 #include "tr1darcy.h"
37 #include "fei2dtrlin.h"
38 #include "gaussintegrationrule.h"
39 #include "gausspoint.h"
40 #include "bcgeomtype.h"
42 #include "transportmaterial.h"
43 #include "load.h"
44 #include "boundaryload.h"
45 #include "mathfem.h"
46 #include "crosssection.h"
47 #include "classfactory.h"
48 
49 namespace oofem {
50 REGISTER_Element(Tr1Darcy);
51 
52 FEI2dTrLin Tr1Darcy :: interpolation_lin(1, 2);
53 
54 Tr1Darcy :: Tr1Darcy(int n, Domain *aDomain) : TransportElement(n, aDomain)
55 {
56  numberOfDofMans = 3;
57 }
58 
60 { }
61 
63 {
64  this->numberOfGaussPoints = 1;
66 }
67 
70 {
71  return & interpolation_lin;
72 }
73 
75 {
76  if ( integrationRulesArray.size() == 0 ) {
77  integrationRulesArray.resize( 1 );
78  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
80  }
81 }
82 
84 {
85  /*
86  * Return Ke = integrate(B^T K B)
87  */
88 
89  FloatMatrix B, BT, K, KB;
90 
91  TransportMaterial *mat = static_cast< TransportMaterial * >( this->giveMaterial() );
92 
93  answer.clear();
94 
95  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
96  const FloatArray &lcoords = gp->giveNaturalCoordinates();
98  double detJ = fabs( this->interpolation_lin.evaldNdx( BT, lcoords, FEIElementGeometryWrapper(this) ) );
99 
100  mat->giveCharacteristicMatrix(K, mode, gp, tStep);
101 
102  B.beTranspositionOf(BT);
103  KB.beProductOf(K, B);
104  answer.plusProductUnsym( B, KB, detJ * gp->giveWeight() ); // Symmetric part is just a single value, not worth it.
105  }
106 }
107 
109 {
110  if ( mtrx == ExternalForcesVector ) {
111  this->computeExternalForcesVector(answer, tStep, mode);
112  } else if ( mtrx == InternalForcesVector ) {
113  this->computeInternalForcesVector(answer, tStep);
114  } else {
115  OOFEM_ERROR("Unknown Type of characteristic mtrx.");
116  }
117 }
118 
120 {
121  FloatArray w, a, gradP, P(1), n;
122  FloatMatrix B, BT;
123 
124  TransportMaterial *mat = static_cast< TransportMaterial * >( this->giveMaterial() );
125 
126  this->computeVectorOf(VM_Total, tStep, a);
127 
128  answer.resize(3);
129  answer.zero();
130 
131  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
132  const FloatArray &lcoords = gp->giveNaturalCoordinates();
133 
134  double detJ = fabs( this->interpolation_lin.giveTransformationJacobian( lcoords, FEIElementGeometryWrapper(this) ) );
135  this->interpolation_lin.evaldNdx( BT, lcoords, FEIElementGeometryWrapper(this) );
136  this->interpolation_lin.evalN( n, lcoords, FEIElementGeometryWrapper(this) );
137  B.beTranspositionOf(BT);
138  P.at(1) = n.dotProduct(a); // Evaluates the field at this point.
139 
140  gradP.beProductOf(B, a);
141 
142  mat->giveFluxVector(w, gp, gradP, P, tStep);
143 
144  answer.plusProduct(B, w, -gp->giveWeight() * detJ);
145  }
146 }
147 
149 {
150  // TODO: Implement support for body forces
151 
152  FloatArray vec;
153 
154  answer.resize(3);
155  answer.zero();
156 
157  // Compute characteristic vector for Neumann boundary conditions.
158  int load_number, load_id;
159  Load *load;
160  bcGeomType ltype;
161 
162  int nLoads = boundaryLoadArray.giveSize() / 2;
163 
164  for ( int i = 1; i <= nLoads; i++ ) { // For each Neumann boundary condition ....
165  load_number = boundaryLoadArray.at(2 * i - 1);
166  load_id = boundaryLoadArray.at(2 * i);
167  load = domain->giveLoad(load_number);
168  ltype = load->giveBCGeoType();
169 
170  if ( ltype == EdgeLoadBGT ) {
171  this->computeEdgeBCSubVectorAt(vec, load, load_id, tStep, mode, 0);
172  }
173 
174  answer.add(vec);
175  }
176 
177  answer.negated();
178 }
179 
180 void Tr1Darcy :: computeEdgeBCSubVectorAt(FloatArray &answer, Load *load, int iEdge, TimeStep *tStep, ValueModeType mode, int indx)
181 {
182  /*
183  * Given the load *load, return it's contribution.
184  *
185  */
186 
187  answer.resize(3);
188  answer.zero();
189 
190  if ( load->giveType() == TransmissionBC ) { // Neumann boundary conditions (traction)
191  BoundaryLoad *boundaryLoad;
192  boundaryLoad = static_cast< BoundaryLoad * >(load);
193 
194  int numberOfEdgeIPs;
195  numberOfEdgeIPs = ( int ) ceil( ( boundaryLoad->giveApproxOrder() + 1. ) / 2. ) * 2;
196 
197  GaussIntegrationRule iRule(1, this, 1, 1);
198  FloatArray N, loadValue, reducedAnswer;
199  reducedAnswer.resize(3);
200  reducedAnswer.zero();
201  IntArray mask;
202 
203  iRule.SetUpPointsOnLine(numberOfEdgeIPs, _Unknown);
204 
205  for ( GaussPoint *gp: iRule ) {
206  const FloatArray &lcoords = gp->giveNaturalCoordinates();
207  this->interpolation_lin.edgeEvalN( N, iEdge, lcoords, FEIElementGeometryWrapper(this) );
208  double dV = this->computeEdgeVolumeAround(gp, iEdge);
209 
210  if ( boundaryLoad->giveFormulationType() == Load :: FT_Entity ) { // Edge load in xi-eta system
211  boundaryLoad->computeValueAt(loadValue, tStep, lcoords, mode);
212  } else { // Edge load in x-y system
213  FloatArray gcoords;
214  this->interpolation_lin.edgeLocal2global( gcoords, iEdge, lcoords, FEIElementGeometryWrapper(this) );
215  boundaryLoad->computeValueAt(loadValue, tStep, gcoords, mode);
216  }
217 
218  reducedAnswer.add(loadValue.at(1) * dV, N);
219  }
220 
221  this->interpolation_lin.computeLocalEdgeMapping(mask, iEdge);
222  answer.assemble(reducedAnswer, mask);
223  }
224 }
225 
227 {
228  return this->giveCrossSection()->give(CS_Thickness, gcoords, this, false);
229 }
230 
231 
233 {
234  double thickness = 1;
236  return detJ *thickness *gp->giveWeight();
237 }
238 
240 {
241  /*
242  * Compute characteristic matrix for this element. The only option is the stiffness matrix...
243  */
244  if ( mtrx == ConductivityMatrix || mtrx == TangentStiffnessMatrix ) {
245  this->computeStiffnessMatrix(answer, TangentStiffness, tStep);
246  } else {
247  OOFEM_ERROR("Unknown Type of characteristic mtrx %d", mtrx);
248  }
249 }
250 
251 void Tr1Darcy :: giveDofManDofIDMask(int inode, IntArray &answer) const
252 {
253  answer = {P_f};
254 }
255 
256 void
258  InternalStateType type, TimeStep *tStep)
259 {
260  CrossSection *cs = this->giveCrossSection();
262  GaussPoint *gp = integrationRulesArray [ 0 ]->getIntegrationPoint(0);
263  cs->giveIPValue(answer, gp, type, tStep);
264 }
265 
266 Interface *
268 {
269  if ( interface == NodalAveragingRecoveryModelInterfaceType ) {
270  return static_cast< NodalAveragingRecoveryModelInterface * >(this);
271  }
272 
273  return NULL;
274 }
275 
276 int
278 {
279  return 3;
280 }
281 }
virtual void edgeLocal2global(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates edge global coordinates from given local ones.
Definition: fei2dtrlin.C:213
CrossSection * giveCrossSection()
Definition: element.C:495
virtual void edgeEvalN(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of edge interpolation functions (shape functions) at given point.
Definition: fei2dtrlin.C:184
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: tr1darcy.C:277
void computeEdgeBCSubVectorAt(FloatArray &answer, Load *load, int iEdge, TimeStep *tStep, ValueModeType mode, int indx)
Definition: tr1darcy.C:180
virtual double evaldNdx(FloatMatrix &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of interpolation functions (shape functions) at given point...
Definition: fei2dtrlin.C:53
Class and object Domain.
Definition: domain.h:115
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dtrlin.C:43
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
virtual void computeLocalEdgeMapping(IntArray &edgeNodes, int iedge)
Definition: fei2dtrlin.C:230
const FloatArray & giveSubPatchCoordinates()
Returns local sub-patch coordinates of the receiver.
Definition: gausspoint.h:144
bcGeomType
Type representing the geometric character of loading.
Definition: bcgeomtype.h:40
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
virtual void computeExternalForcesVector(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Definition: tr1darcy.C:148
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
virtual double giveThicknessAt(const FloatArray &gcoords)
Gives the thickness at some global coordinate.
Definition: tr1darcy.C:226
virtual double edgeGiveTransformationJacobian(int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the edge Jacobian of transformation between local and global coordinates.
Definition: feinterpol2d.C:175
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
MatResponseMode
Describes the character of characteristic material matrix.
Base abstract class representing cross section in finite element mesh.
Definition: crosssection.h:107
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: tr1darcy.C:62
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:780
Distributed edge load.
Definition: bcgeomtype.h:44
Abstract base class representing a boundary load (force, momentum, ...) that acts directly on a bound...
Definition: boundaryload.h:110
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: element.C:638
virtual void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
Definition: tr1darcy.C:83
double dotProduct(const FloatArray &x) const
Computes the dot product (or inner product) of receiver and argument.
Definition: floatarray.C:463
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
virtual double computeEdgeVolumeAround(GaussPoint *gp, int iEdge)
Computes the length around a integration point on a edge.
Definition: tr1darcy.C:232
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
virtual FEInterpolation * giveInterpolation() const
Definition: tr1darcy.C:69
virtual FormulationType giveFormulationType()
Specifies is load should take local or global coordinates.
Definition: load.h:155
This abstract class represent a general base element class for transport problems.
Neumann type (prescribed flux).
Definition: bctype.h:43
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: tr1darcy.C:74
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
#define N(p, q)
Definition: mdm.C:367
virtual int giveApproxOrder()=0
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Definition: tr1darcy.C:251
virtual ~Tr1Darcy()
Definition: tr1darcy.C:59
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
int numberOfGaussPoints
Number of integration points as specified by nip.
Definition: element.h:188
virtual Interface * giveInterface(InterfaceType interface)
Interface requesting service.
Definition: tr1darcy.C:267
Class representing vector of real numbers.
Definition: floatarray.h:82
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType mtrx, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
Definition: tr1darcy.C:239
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 int giveIPValue(FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep)
Returns the integration point corresponding value in Reduced form.
Definition: crosssection.C:89
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei2dtrlin.C:147
void assemble(const FloatArray &fe, const IntArray &loc)
Assembles the array fe (typically, the load vector of a finite element) into the receiver, using loc as location array.
Definition: floatarray.C:551
CharType
Definition: chartype.h:87
Class representing the general Input Record.
Definition: inputrecord.h:101
Class Interface.
Definition: interface.h:82
virtual void computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
Definition: tr1darcy.C:119
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
Abstract base class for all constitutive models for transport problems.
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
Definition: floatmatrix.C:323
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
virtual void giveCharacteristicMatrix(FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
Computes characteristic matrix of receiver in given integration point.
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
Load is base abstract class for all loads.
Definition: load.h:61
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:337
virtual void giveFluxVector(FloatArray &answer, GaussPoint *gp, const FloatArray &grad, const FloatArray &field, TimeStep *tStep)=0
Returns the flux for the field and its gradient.
int giveSize() const
Definition: intarray.h:203
the oofem namespace is to define a context or scope in which all oofem names are defined.
void clear()
Sets size of receiver to be an empty matrix. It will have zero rows and zero columns size...
Definition: floatmatrix.h:516
void negated()
Switches the sign of every coefficient of receiver.
Definition: floatarray.C:739
static FEI2dTrLin interpolation_lin
Definition: tr1darcy.h:54
virtual Material * giveMaterial()
virtual void computeValueAt(FloatArray &answer, TimeStep *tStep, const FloatArray &coords, ValueModeType mode)
Computes components values of load at given point - global coordinates (coordinates given)...
Definition: boundaryload.C:58
Load * giveLoad(int n)
Service for accessing particular domain load.
Definition: domain.C:222
virtual int SetUpPointsOnLine(int nPoints, MaterialMode mode)
Sets up receiver&#39;s integration points on unit line integration domain.
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: tr1darcy.C:257
Class representing integration point in finite element program.
Definition: gausspoint.h:93
IntArray boundaryLoadArray
Definition: element.h:160
Class representing solution step.
Definition: timestep.h:80
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
virtual void giveCharacteristicVector(FloatArray &answer, CharType mtrx, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
Definition: tr1darcy.C:108
Class representing Gaussian-quadrature integration rule.
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
Tr1Darcy(int, Domain *)
Definition: tr1darcy.C:54
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .
Definition: floatarray.C:226

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:32 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011