OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tr21stokes.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 #include "tr21stokes.h"
36 #include "node.h"
37 #include "dof.h"
38 #include "domain.h"
39 #include "gaussintegrationrule.h"
40 #include "gausspoint.h"
41 #include "bcgeomtype.h"
42 #include "load.h"
43 #include "bodyload.h"
44 #include "boundaryload.h"
45 #include "mathfem.h"
46 #include "fluiddynamicmaterial.h"
47 #include "fei2dtrlin.h"
48 #include "fei2dtrquad.h"
49 #include "fluidcrosssection.h"
50 #include "assemblercallback.h"
51 #include "classfactory.h"
52 #include "engngm.h"
53 
54 namespace oofem {
55 REGISTER_Element(Tr21Stokes);
56 
57 // Set up interpolation coordinates
58 FEI2dTrLin Tr21Stokes :: interpolation_lin(1, 2);
59 FEI2dTrQuad Tr21Stokes :: interpolation_quad(1, 2);
60 // Set up ordering vectors (for assembling)
61 IntArray Tr21Stokes :: momentum_ordering = {1, 2, 4, 5, 7, 8, 10, 11, 12, 13, 14, 15};
62 IntArray Tr21Stokes :: conservation_ordering = {3, 6, 9};
63 IntArray Tr21Stokes :: edge_ordering [ 3 ] = {
64  {1, 2, 4, 5, 10, 11},
65  {4, 5, 7, 8, 12, 13},
66  {7, 8, 1, 2, 14, 15}
67 };
68 
70 {
71  this->numberOfDofMans = 6;
72  if ( aDomain->giveEngngModel()->giveProblemScale() == macroScale ) { // Pick the lowest default value for multiscale computations.
73  this->numberOfGaussPoints = 3;
74  } else {
75  this->numberOfGaussPoints = 4;
76  }
77 }
78 
80 { }
81 
83 {
84  if ( integrationRulesArray.size() == 0 ) {
85  integrationRulesArray.resize(1);
86  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
88  }
89 }
90 
92 {
93  return 15;
94 }
95 
96 void Tr21Stokes :: giveDofManDofIDMask(int inode, IntArray &answer) const
97 {
98  if ( inode <= 3 ) {
99  answer = {V_u, V_v, P_f};
100  } else {
101  answer = {V_u, V_v};
102  }
103 }
104 
106 {
108  return detJ *gp->giveWeight();
109 }
110 
112  TimeStep *tStep)
113 {
114  // Compute characteristic vector for this element. I.e the load vector(s)
115  if ( mtrx == ExternalForcesVector ) {
116  this->computeExternalForcesVector(answer, tStep);
117  } else if ( mtrx == InternalForcesVector ) {
118  this->computeInternalForcesVector(answer, tStep);
119  } else {
120  OOFEM_ERROR("Unknown Type of characteristic mtrx.");
121  }
122 }
123 
125  CharType mtrx, TimeStep *tStep)
126 {
127  // Compute characteristic matrix for this element. The only option is the stiffness matrix...
128  if ( mtrx == TangentStiffnessMatrix ) {
129  this->computeStiffnessMatrix(answer, TangentStiffness, tStep);
130  } else {
131  OOFEM_ERROR("Unknown Type of characteristic mtrx.");
132  }
133 }
134 
136 {
137  FluidDynamicMaterial *mat = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial();
138  FloatArray a_pressure, a_velocity, devStress, epsp, Nh, dNv(12);
139  double r_vol, pressure;
140  FloatMatrix dN, B(3, 12);
141  B.zero();
142 
143  this->computeVectorOfVelocities(VM_Total, tStep, a_velocity);
144  this->computeVectorOfPressures(VM_Total, tStep, a_pressure);
145 
146  FloatArray momentum, conservation;
147 
148  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
149  const FloatArray &lcoords = gp->giveNaturalCoordinates();
150 
151  double detJ = fabs( this->interpolation_quad.evaldNdx( dN, lcoords, FEIElementGeometryWrapper(this) ) );
152  this->interpolation_lin.evalN( Nh, lcoords, FEIElementGeometryWrapper(this) );
153  double dA = detJ * gp->giveWeight();
154 
155  for ( int j = 0, k = 0; j < dN.giveNumberOfRows(); j++, k += 2 ) {
156  dNv(k) = B(0, k) = B(2, k + 1) = dN(j, 0);
157  dNv(k + 1) = B(1, k + 1) = B(2, k) = dN(j, 1);
158  }
159 
160  pressure = Nh.dotProduct(a_pressure);
161  epsp.beProductOf(B, a_velocity);
162 
163  mat->computeDeviatoricStressVector(devStress, r_vol, gp, epsp, pressure, tStep);
164 
165  momentum.plusProduct(B, devStress, dA);
166  momentum.add(-pressure * dA, dNv);
167  conservation.add(r_vol * dA, Nh);
168  }
169 
170  answer.resize(15);
171  answer.zero();
172  answer.assemble(momentum, this->momentum_ordering);
173  answer.assemble(conservation, this->conservation_ordering);
174 }
175 
177 {
178  FloatArray vec;
179 
180  answer.clear();
181 
182  int nLoads = this->boundaryLoadArray.giveSize() / 2;
183  for ( int i = 1; i <= nLoads; i++ ) { // For each Neumann boundary condition
184  int load_number = this->boundaryLoadArray.at(2 * i - 1);
185  int load_id = this->boundaryLoadArray.at(2 * i);
186  Load *load = this->domain->giveLoad(load_number);
187  bcGeomType ltype = load->giveBCGeoType();
188 
189  if ( ltype == EdgeLoadBGT ) {
190  this->computeBoundarySurfaceLoadVector(vec, static_cast< BoundaryLoad * >(load), load_id, ExternalForcesVector, VM_Total, tStep);
191  answer.add(vec);
192  }
193  }
194 
195  BodyLoad *bload;
196  nLoads = this->giveBodyLoadArray()->giveSize();
197  for ( int i = 1; i <= nLoads; i++ ) {
198  Load *load = domain->giveLoad( bodyLoadArray.at(i) );
199  if ((bload = dynamic_cast<BodyLoad*>(load))) {
200  bcGeomType ltype = load->giveBCGeoType();
201  if ( ltype == BodyLoadBGT && load->giveBCValType() == ForceLoadBVT ) {
202  this->computeLoadVector(vec, bload, ExternalForcesVector, VM_Total, tStep);
203  answer.add(vec);
204  }
205  }
206  }
207 }
208 
210 {
211  if ( type != ExternalForcesVector ) {
212  answer.clear();
213  return;
214  }
215 
216  FluidDynamicMaterial *mat = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial();
217  FloatArray N, gVector, temparray(12);
218 
219  load->computeComponentArrayAt(gVector, tStep, VM_Total);
220  temparray.zero();
221  if ( gVector.giveSize() ) {
222  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
223  const FloatArray &lcoords = gp->giveNaturalCoordinates();
224 
225  double rho = mat->give('d', gp);
226  double detJ = fabs( this->interpolation_quad.giveTransformationJacobian( lcoords, FEIElementGeometryWrapper(this) ) );
227  double dA = detJ * gp->giveWeight();
228 
229  this->interpolation_quad.evalN( N, lcoords, FEIElementGeometryWrapper(this) );
230  for ( int j = 0; j < 6; j++ ) {
231  temparray(2 * j) += N(j) * rho * gVector(0) * dA;
232  temparray(2 * j + 1) += N(j) * rho * gVector(1) * dA;
233  }
234  }
235  }
236 
237  answer.resize(15);
238  answer.zero();
239  answer.assemble(temparray, this->momentum_ordering);
240 }
241 
242  void Tr21Stokes :: computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
243 {
244  if ( type != ExternalForcesVector ) {
245  answer.clear();
246  return;
247  }
248 
249  if ( load->giveType() == TransmissionBC ) { // Neumann boundary conditions (traction)
250 
251  int numberOfEdgeIPs = ( int ) ceil( ( load->giveApproxOrder() + 1. ) / 2. ) * 2;
252 
253  GaussIntegrationRule iRule(1, this, 1, 1);
254  FloatArray N, t, f(6);
255 
256  f.zero();
257  iRule.SetUpPointsOnLine(numberOfEdgeIPs, _Unknown);
258 
259  for ( GaussPoint *gp: iRule ) {
260  const FloatArray &lcoords = gp->giveNaturalCoordinates();
261 
262  this->interpolation_quad.edgeEvalN( N, boundary, lcoords, FEIElementGeometryWrapper(this) );
263  double detJ = fabs( this->interpolation_quad.boundaryGiveTransformationJacobian( boundary, lcoords, FEIElementGeometryWrapper(this) ) );
264  double dS = gp->giveWeight() * detJ;
265 
266  if ( load->giveFormulationType() == Load :: FT_Entity ) { // Edge load in xi-eta system
267  load->computeValueAt(t, tStep, lcoords, VM_Total);
268  } else { // Edge load in x-y system
269  FloatArray gcoords;
270  this->interpolation_quad.boundaryLocal2Global( gcoords, boundary, lcoords, FEIElementGeometryWrapper(this) );
271  load->computeValueAt(t, tStep, gcoords, VM_Total);
272  }
273 
274  // Reshape the vector
275  for ( int j = 0; j < 3; j++ ) {
276  f(2 * j) += N(j) * t(0) * dS;
277  f(2 * j + 1) += N(j) * t(1) * dS;
278  }
279  }
280 
281  answer.resize(15);
282  answer.zero();
283  answer.assemble(f, this->edge_ordering [ boundary - 1 ]);
284  } else {
285  OOFEM_ERROR("Strange boundary condition type");
286  }
287 }
288 
290 {
291  // Note: Working with the components; [K, G+Dp; G^T+Dv^T, C] . [v,p]
292  FluidDynamicMaterial *mat = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial();
293  FloatMatrix B(3, 12), EdB, K, G, Dp, DvT, C, Ed, dN;
294  FloatArray dN_V(12), Nlin, Ep, Cd, tmpA, tmpB;
295  double Cp;
296 
297  B.zero();
298 
299  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
300  // Compute Gauss point and determinant at current element
301  const FloatArray &lcoords = gp->giveNaturalCoordinates();
302 
303  this->interpolation_lin.evalN( Nlin, lcoords, FEIElementGeometryWrapper(this) );
304  double detJ = fabs( this->interpolation_quad.evaldNdx( dN, lcoords, FEIElementGeometryWrapper(this) ) );
305  double dA = detJ * gp->giveWeight();
306 
307  for ( int j = 0, k = 0; j < 6; j++, k += 2 ) {
308  dN_V(k) = B(0, k) = B(2, k + 1) = dN(j, 0);
309  dN_V(k + 1) = B(1, k + 1) = B(2, k) = dN(j, 1);
310  }
311 
312  // Computing the internal forces should have been done first.
313  // dsigma_dev/deps_dev dsigma_dev/dp deps_vol/deps_dev deps_vol/dp
314  mat->giveStiffnessMatrices(Ed, Ep, Cd, Cp, mode, gp, tStep);
315 
316  EdB.beProductOf(Ed, B);
317  K.plusProductSymmUpper(B, EdB, dA);
318  G.plusDyadUnsym(dN_V, Nlin, -dA);
319  C.plusDyadSymmUpper(Nlin, Cp * dA);
320 
321  tmpA.beTProductOf(B, Ep);
322  Dp.plusDyadUnsym(tmpA, Nlin, dA);
323 
324  tmpB.beTProductOf(B, Cd);
325  DvT.plusDyadUnsym(Nlin, tmpB, dA);
326  }
327 
328  K.symmetrized();
329  C.symmetrized();
330  FloatMatrix GTDvT, GDp;
331  GTDvT.beTranspositionOf(G);
332  GTDvT.add(DvT);
333  GDp = G;
334  GDp.add(Dp);
335 
336  answer.resize(15, 15);
337  answer.zero();
338  answer.assemble(K, this->momentum_ordering);
339  answer.assemble(GTDvT, this->conservation_ordering, this->momentum_ordering);
340  answer.assemble(GDp, this->momentum_ordering, this->conservation_ordering);
341  answer.assemble(C, this->conservation_ordering);
342 }
343 
345 {
346  return & interpolation_quad;
347 }
348 
350 {
351  if ( id == P_f ) {
352  return & interpolation_lin;
353  } else {
354  return & interpolation_quad;
355  }
356 }
357 
359 {
361 }
362 
363 // Some extension Interfaces to follow:
364 
366 {
367  switch ( it ) {
369  return static_cast< NodalAveragingRecoveryModelInterface * >(this);
370 
372  return static_cast< ZZNodalRecoveryModelInterface * >(this);
373 
375  return static_cast< SpatialLocalizerInterface * >(this);
376 
377  default:
378  return FMElement :: giveInterface(it);
379  }
380 }
381 
382 void Tr21Stokes :: computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
383 {
384  FloatArray n, n_lin, pressures, velocities;
385  this->interpolation_quad.evalN( n, lcoords, FEIElementGeometryWrapper(this) );
386  this->interpolation_lin.evalN( n_lin, lcoords, FEIElementGeometryWrapper(this) );
387  this->computeVectorOf({P_f}, mode, tStep, pressures);
388  this->computeVectorOf({V_u, V_v, V_w}, mode, tStep, velocities);
389 
390  answer.resize(4);
391  answer.zero();
392  for ( int i = 1; i <= n.giveSize(); i++ ) {
393  answer(0) += n.at(i) * velocities.at(i*2-1);
394  answer(1) += n.at(i) * velocities.at(i*2);
395  }
396 
397  for ( int i = 1; i <= n_lin.giveSize(); i++ ) {
398  answer(2) += n_lin.at(i) * pressures.at(i);
399  }
400 }
401 
402 
404 {
405  if ( type == IST_Pressure ) {
406  answer.resize(1);
407  if ( node == 1 || node == 2 || node == 3 ) {
408  answer.at(1) = this->giveNode(node)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep);
409  } else {
410  double a, b;
411  if ( node == 4 ) {
412  a = this->giveNode(1)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep);
413  b = this->giveNode(2)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep);
414  } else if ( node == 5 ) {
415  a = this->giveNode(2)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep);
416  b = this->giveNode(3)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep);
417  } else { /*if ( node == 6 )*/
418  a = this->giveNode(3)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep);
419  b = this->giveNode(1)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep);
420  }
421 
422  answer.at(1) = ( a + b ) / 2;
423  }
424  } else {
425  answer.clear();
426  }
427 }
428 
429 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
Tr21Stokes(int n, Domain *d)
Definition: tr21stokes.C:69
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
The element interface required by NodalAvergagingRecoveryModel.
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
Fluid cross-section.
problemScale giveProblemScale()
Returns scale in multiscale simulation.
Definition: engngm.h:418
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: tr21stokes.C:365
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: feinterpol.C:43
Abstract base class for all fluid materials.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Gives the dof ID mask for the element.
Definition: tr21stokes.C:96
static FEI2dTrQuad interpolation_quad
Interpolation for geometry and velocity.
Definition: tr21stokes.h:66
The element interface required by ZZNodalRecoveryModel.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: tr21stokes.C:358
bcGeomType
Type representing the geometric character of loading.
Definition: bcgeomtype.h:40
double & at(int i)
Coefficient access function.
Definition: floatarray.h:131
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
Class implementing element body load, acting over whole element volume (e.g., the dead weight)...
Definition: bodyload.h:49
virtual ~Tr21Stokes()
Definition: tr21stokes.C:79
void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
Definition: tr21stokes.C:289
void computeExternalForcesVector(FloatArray &answer, TimeStep *tStep)
Definition: tr21stokes.C:176
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual void computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int boundary, CharType type, ValueModeType mode, TimeStep *tStep, bool global=true)
Computes the contribution of the given load at the given boundary surface in global coordinate system...
Definition: tr21stokes.C:242
virtual double giveUnknown(ValueModeType mode, TimeStep *tStep)=0
The key method of class Dof.
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: element.C:779
virtual void computeComponentArrayAt(FloatArray &answer, TimeStep *tStep, ValueModeType mode)
Computes boundary condition value - its components values at given time.
Definition: load.C:82
void plusDyadUnsym(const FloatArray &a, const FloatArray &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:811
virtual void giveStiffnessMatrices(FloatMatrix &dsdd, FloatArray &dsdp, FloatArray &dedd, double &dedp, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
Computes the 4 tangents of the compressible material response in 3D.
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.
Distributed body load.
Definition: bcgeomtype.h:43
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: tr21stokes.C:82
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: tr21stokes.C:105
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
virtual void computeVectorOfPressures(ValueModeType mode, TimeStep *tStep, FloatArray &pressures)
Definition: fmelement.C:55
virtual double give(int aProperty, GaussPoint *gp)
Returns the value of material property &#39;aProperty&#39;.
Definition: material.C:52
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei2dtrquad.C:43
static IntArray momentum_ordering
Ordering of dofs in element. Used to assemble the element stiffness.
Definition: tr21stokes.h:68
virtual void computeDeviatoricStressVector(FloatArray &stress_dev, double &epsp_vol, GaussPoint *gp, const FloatArray &eps, double pressure, TimeStep *tStep)
Computes the deviatoric stress vector and volumetric strain rate from given deviatoric strain and pre...
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 void computeVectorOfVelocities(ValueModeType mode, TimeStep *tStep, FloatArray &velocities)
Definition: fmelement.C:48
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 void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
Definition: tr21stokes.C:111
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
virtual FormulationType giveFormulationType()
Specifies is load should take local or global coordinates.
Definition: load.h:155
Neumann type (prescribed flux).
Definition: bctype.h:43
static FEI2dTrLin interpolation_lin
Interpolation for pressure.
Definition: tr21stokes.h:64
virtual double boundaryGiveTransformationJacobian(int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation Jacobian on the requested boundary.
Definition: feinterpol2d.C:154
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: fei2dtrquad.C:60
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
IntArray * giveBodyLoadArray()
Returns array containing load numbers of loads acting on element.
Definition: element.C:372
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
IntArray bodyLoadArray
Array containing indexes of loads (body loads and boundary loads are kept separately), that apply on receiver.
Definition: element.h:160
#define N(p, q)
Definition: mdm.C:367
void plusProductSymmUpper(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:698
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
Definition: tr21stokes.C:124
virtual int giveApproxOrder()=0
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
Definition: floatarray.C:708
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: tr21stokes.C:403
virtual bcType giveType() const
Returns receiver load type.
Definition: boundaryload.h:166
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
Class representing vector of real numbers.
Definition: floatarray.h:82
void plusDyadSymmUpper(const FloatArray &a, double dV)
Adds to the receiver the dyadic product .
Definition: floatmatrix.C:756
This abstract class represent a general base element class for fluid dynamic problems.
Definition: fmelement.h:54
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
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
void resize(int rows, int cols)
Checks size of receiver towards requested bounds.
Definition: floatmatrix.C:1358
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
Dof * giveDofWithID(int dofID) const
Returns DOF with given dofID; issues error if not present.
Definition: dofmanager.C:119
Class Interface.
Definition: interface.h:82
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
virtual bcGeomType giveBCGeoType() const
Returns geometry character of boundary condition.
virtual Interface * giveInterface(InterfaceType t)
Interface requesting service.
Definition: femcmpnn.h:179
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
Definition: floatmatrix.C:323
The spatial localizer element interface associated to spatial localizer.
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
InterfaceType
Enumerative type, used to identify interface type.
Definition: interfacetype.h:43
void computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
Definition: tr21stokes.C:135
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
int giveSize() const
Definition: intarray.h:203
virtual FEInterpolation * giveInterpolation() const
Definition: tr21stokes.C:344
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.
void assemble(const FloatMatrix &src, const IntArray &loc)
Assembles the contribution using localization array into receiver.
Definition: floatmatrix.C:251
virtual bcValType giveBCValType() const
Returns receiver load type.
virtual int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
Definition: tr21stokes.C:91
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
virtual void boundaryLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Maps the local boundary coordinates to global.
Definition: feinterpol2d.C:159
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
int giveNumberOfRows() const
Returns number of rows of receiver.
Definition: floatmatrix.h:156
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.
void symmetrized()
Initializes the lower half of the receiver according to the upper half.
Definition: floatmatrix.C:1576
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: fei2dtrquad.C:146
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
virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
Computes the contribution of the given body load (volumetric).
Definition: tr21stokes.C:209
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 computeField(ValueModeType u, TimeStep *tStep, const FloatArray &coords, FloatArray &answer)
Computes the unknown vector interpolated at the specified local coordinates.
Definition: tr21stokes.C:382
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
Class representing Gaussian-quadrature integration rule.
static IntArray edge_ordering[3]
Ordering of dofs on edges. Used to assemble edge loads.
Definition: tr21stokes.h:70
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631
void plusProduct(const FloatMatrix &b, const FloatArray &s, double dV)
Adds the product .
Definition: floatarray.C:226
static IntArray conservation_ordering
Definition: tr21stokes.h:68

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