OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
hexa21stokes.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 "hexa21stokes.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"
43 #include "load.h"
44 #include "bodyload.h"
45 #include "boundaryload.h"
46 #include "mathfem.h"
47 #include "fluiddynamicmaterial.h"
48 #include "fei3dhexalin.h"
49 #include "fei3dhexatriquad.h"
50 #include "fluidcrosssection.h"
51 #include "assemblercallback.h"
52 #include "classfactory.h"
53 
54 namespace oofem {
55 REGISTER_Element(Hexa21Stokes);
56 
57 // Set up interpolation coordinates
59 FEI3dHexaTriQuad Hexa21Stokes :: interpolation_quad;
60 // Set up ordering vectors (for assembling)
62  1, 2, 3, 5, 6, 7, 9, 10, 11, 13, 14, 15, 17, 18, 19, 21, 22, 23, 25, 26, 27, 29, 30, 31, 33, 34, 35,
63  37, 38, 39, 41, 42, 43, 45, 46, 47, 49, 50, 51, 53, 54, 55, 57, 58, 59, 61, 62, 63, 65, 66, 67, 69, 70, 71,
64  73, 74, 75, 77, 78, 79, 81, 82, 83, 85, 86, 87, 89, 90, 91, 93, 94, 95, 97, 98, 99, 101, 102, 103, 105, 106, 107};
65 IntArray Hexa21Stokes :: conservation_ordering = {4, 8, 12, 16, 20, 24, 28, 32};
66 IntArray Hexa21Stokes :: surf_ordering [ 6 ] = {
67  { 5, 6, 7, 1, 2, 3, 13, 14, 15, 9, 10, 11, 33, 34, 35, 42, 43, 44, 39, 40, 41, 36, 37, 38, 69, 70, 71},
68  {17, 18, 19, 21, 22, 23, 25, 26, 27, 29, 30, 31, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 72, 73, 74},
69  { 1, 2, 3, 17, 18, 19, 21, 22, 23, 5, 6, 7, 57, 58, 59, 45, 46, 47, 60, 61, 62, 33, 34, 35, 75, 76, 77},
70  { 5, 6, 7, 9, 10, 11, 25, 26, 27, 21, 22, 23, 36, 37, 38, 63, 64, 65, 48, 49, 50, 60, 61, 62, 78, 79, 80},
71  { 9, 10, 11, 13, 14, 15, 29, 30, 31, 25, 26, 27, 39, 40, 41, 66, 67, 68, 51, 52, 53, 63, 64, 65, 81, 82, 83},
72  {13, 14, 15, 1, 2, 3, 17, 18, 19, 29, 30, 31, 42, 43, 44, 57, 58, 59, 54, 55, 56, 66, 67, 68, 84, 85, 86}
73 };
74 
76 {
77  this->numberOfDofMans = 27;
78  this->numberOfGaussPoints = 27;
79 }
80 
82 { }
83 
85 {
86  if ( integrationRulesArray.size() == 0 ) {
87  integrationRulesArray.resize(1);
88  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
90  }
91 }
92 
94 {
95  return 89;
96 }
97 
98 void Hexa21Stokes :: giveDofManDofIDMask(int inode, IntArray &answer) const
99 {
100  if ( inode <= 8 ) {
101  answer = {V_u, V_v, V_w, P_f};
102  } else {
103  answer = {V_u, V_v, V_w};
104  }
105 }
106 
108  TimeStep *tStep)
109 {
110  // Compute characteristic vector for this element. I.e the load vector(s)
111  if ( mtrx == ExternalForcesVector ) {
112  this->computeExternalForcesVector(answer, tStep);
113  } else if ( mtrx == InternalForcesVector ) {
114  this->computeInternalForcesVector(answer, tStep);
115  } else {
116  OOFEM_ERROR("Unknown Type of characteristic mtrx.");
117  }
118 }
119 
121  CharType mtrx, TimeStep *tStep)
122 {
123  // Compute characteristic matrix for this element. The only option is the stiffness matrix...
124  if ( mtrx == TangentStiffnessMatrix ) {
125  this->computeStiffnessMatrix(answer, TangentStiffness, tStep);
126  } else {
127  OOFEM_ERROR("Unknown Type of characteristic mtrx.");
128  }
129 }
130 
132 {
133  FluidDynamicMaterial *mat = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial();
134  FloatArray a_pressure, a_velocity, devStress, epsp, Nh, dN_V(81);
135  FloatMatrix dN, B(6, 81);
136  double r_vol, pressure;
137  this->computeVectorOfVelocities(VM_Total, tStep, a_velocity);
138  this->computeVectorOfPressures(VM_Total, tStep, a_pressure);
139  FloatArray momentum, conservation;
140 
141  B.zero();
142  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
143  const FloatArray &lcoords = gp->giveNaturalCoordinates();
144 
145  double detJ = fabs( this->interpolation_quad.evaldNdx( dN, lcoords, FEIElementGeometryWrapper(this) ) );
146  this->interpolation_lin.evalN( Nh, lcoords, FEIElementGeometryWrapper(this) );
147  double dV = detJ * gp->giveWeight();
148 
149  for ( int j = 0, k = 0; j < dN.giveNumberOfRows(); j++, k += 3 ) {
150  dN_V(k + 0) = B(0, k + 0) = B(5, k + 1) = B(4, k + 2) = dN(j, 0);
151  dN_V(k + 1) = B(1, k + 1) = B(5, k + 0) = B(3, k + 2) = dN(j, 1);
152  dN_V(k + 2) = B(2, k + 2) = B(4, k + 0) = B(3, k + 1) = dN(j, 2);
153  }
154 
155  epsp.beProductOf(B, a_velocity);
156  pressure = Nh.dotProduct(a_pressure);
157  mat->computeDeviatoricStressVector(devStress, r_vol, gp, epsp, pressure, tStep);
158 
159  momentum.plusProduct(B, devStress, dV);
160  momentum.add(-pressure * dV, dN_V);
161  conservation.add(r_vol * dV, Nh);
162  }
163 
164  answer.resize(89);
165  answer.zero();
166  answer.assemble(momentum, this->momentum_ordering);
167  answer.assemble(conservation, this->conservation_ordering);
168 }
169 
171 {
172  FloatArray vec;
173 
174  int nLoads = this->boundaryLoadArray.giveSize() / 2;
175  answer.clear();
176 
177  for ( int i = 1; i <= nLoads; i++ ) { // For each Neumann boundary condition
178  int load_number = this->boundaryLoadArray.at(2 * i - 1);
179  int load_id = this->boundaryLoadArray.at(2 * i);
180  Load *load = this->domain->giveLoad(load_number);
181 
182  if ( load->giveBCGeoType() == SurfaceLoadBGT ) {
183  this->computeBoundarySurfaceLoadVector(vec, static_cast< BoundaryLoad * >(load), load_id, ExternalForcesVector, VM_Total, tStep);
184  answer.add(vec);
185  }
186  }
187 
188  nLoads = this->giveBodyLoadArray()->giveSize();
189  for ( int i = 1; i <= nLoads; i++ ) {
190  BodyLoad *load;
191  if ((load = dynamic_cast<BodyLoad*>(domain->giveLoad( bodyLoadArray.at(i))))) {
192  if ( load->giveBCGeoType() == BodyLoadBGT && load->giveBCValType() == ForceLoadBVT ) {
193  this->computeLoadVector(vec, load, ExternalForcesVector, VM_Total, tStep);
194  answer.add(vec);
195  }
196  }
197  }
198 }
199 
201 {
202  if ( type != ExternalForcesVector ) {
203  answer.clear();
204  return;
205  }
206 
207  FloatArray N, gVector, temparray(81);
208 
209  load->computeComponentArrayAt(gVector, tStep, VM_Total);
210  temparray.zero();
211  if ( gVector.giveSize() ) {
212  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
213  const FloatArray &lcoords = gp->giveNaturalCoordinates();
214 
215  double rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
216  double detJ = fabs( this->interpolation_quad.giveTransformationJacobian( lcoords, FEIElementGeometryWrapper(this) ) );
217  double dA = detJ * gp->giveWeight();
218 
219  this->interpolation_quad.evalN( N, lcoords, FEIElementGeometryWrapper(this) );
220  for ( int j = 0; j < N.giveSize(); j++ ) {
221  temparray(3 * j + 0) += N(j) * rho * gVector(0) * dA;
222  temparray(3 * j + 1) += N(j) * rho * gVector(1) * dA;
223  temparray(3 * j + 2) += N(j) * rho * gVector(2) * dA;
224  }
225  }
226  }
227 
228  answer.resize(89);
229  answer.zero();
230  answer.assemble(temparray, this->momentum_ordering);
231 }
232 
233  void Hexa21Stokes :: computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int iSurf, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
234 {
235  if ( type != ExternalForcesVector ) {
236  answer.clear();
237  return;
238  }
239 
240  if ( load->giveType() == TransmissionBC ) { // Neumann boundary conditions (traction)
241 
242  int numberOfSurfaceIPs = ( int ) ceil( ( load->giveApproxOrder() + 1. ) / 2. ) * 2;
243 
244  GaussIntegrationRule iRule(1, this, 1, 1);
245  FloatArray N, t, f(27);
246 
247  f.zero();
248  iRule.SetUpPointsOnTriangle(numberOfSurfaceIPs, _Unknown);
249 
250  for ( GaussPoint *gp: iRule ) {
251  const FloatArray &lcoords = gp->giveNaturalCoordinates();
252 
253  this->interpolation_quad.surfaceEvalN( N, iSurf, lcoords, FEIElementGeometryWrapper(this) );
254  double dA = gp->giveWeight() * this->interpolation_quad.surfaceGiveTransformationJacobian( iSurf, lcoords, FEIElementGeometryWrapper(this) );
255 
256  if ( load->giveFormulationType() == Load :: FT_Entity ) { // load in xi-eta system
257  load->computeValueAt(t, tStep, lcoords, VM_Total);
258  } else { // Edge load in x-y system
259  FloatArray gcoords;
260  this->interpolation_quad.surfaceLocal2global( gcoords, iSurf, lcoords, FEIElementGeometryWrapper(this) );
261  load->computeValueAt(t, tStep, gcoords, VM_Total);
262  }
263 
264  // Reshape the vector
265  for ( int j = 0; j < N.giveSize(); j++ ) {
266  f(3 * j + 0) += N(j) * t(0) * dA;
267  f(3 * j + 1) += N(j) * t(1) * dA;
268  f(3 * j + 2) += N(j) * t(2) * dA;
269  }
270  }
271 
272  answer.resize(89);
273  answer.zero();
274  answer.assemble(f, this->surf_ordering [ iSurf - 1 ]);
275  } else {
276  OOFEM_ERROR("Strange boundary condition type");
277  }
278 }
279 
281 {
282  FluidDynamicMaterial *mat = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial();
283  FloatMatrix B(6, 81), EdB, K, G, Dp, DvT, C, Ed, dN;
284  FloatArray dN_V(81), Nlin, Ep, Cd, tmpA, tmpB;
285  double Cp;
286 
287  B.zero();
288 
289  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
290  // Compute Gauss point and determinant at current element
291  const FloatArray &lcoords = gp->giveNaturalCoordinates();
292 
293  double detJ = fabs( this->interpolation_quad.evaldNdx( dN, lcoords, FEIElementGeometryWrapper(this) ) );
294  double dV = detJ * gp->giveWeight();
295  this->interpolation_lin.evalN( Nlin, lcoords, FEIElementGeometryWrapper(this) );
296 
297  for ( int j = 0, k = 0; j < dN.giveNumberOfRows(); j++, k += 3 ) {
298  dN_V(k + 0) = B(0, k + 0) = B(5, k + 1) = B(4, k + 2) = dN(j, 0);
299  dN_V(k + 1) = B(1, k + 1) = B(5, k + 0) = B(3, k + 2) = dN(j, 1);
300  dN_V(k + 2) = B(2, k + 2) = B(4, k + 0) = B(3, k + 1) = dN(j, 2);
301  }
302 
303  // Computing the internal forces should have been done first.
304  // dsigma_dev/deps_dev dsigma_dev/dp deps_vol/deps_dev deps_vol/dp
305  mat->giveStiffnessMatrices(Ed, Ep, Cd, Cp, mode, gp, tStep);
306 
307  EdB.beProductOf(Ed, B);
308  K.plusProductSymmUpper(B, EdB, dV);
309  G.plusDyadUnsym(dN_V, Nlin, -dV);
310  C.plusDyadSymmUpper(Nlin, Cp * dV);
311 
312  tmpA.beTProductOf(B, Ep);
313  Dp.plusDyadUnsym(tmpA, Nlin, dV);
314 
315  tmpB.beTProductOf(B, Cd);
316  DvT.plusDyadUnsym(Nlin, tmpB, dV);
317  }
318 
319  K.symmetrized();
320  C.symmetrized();
321  FloatMatrix GTDvT, GDp;
322  GTDvT.beTranspositionOf(G);
323  GTDvT.add(DvT);
324  GDp = G;
325  GDp.add(Dp);
326 
327  answer.resize(89, 89);
328  answer.zero();
329  answer.assemble(K, this->momentum_ordering);
330  answer.assemble(GDp, this->momentum_ordering, this->conservation_ordering);
331  answer.assemble(GTDvT, this->conservation_ordering, this->momentum_ordering);
332  answer.assemble(C, this->conservation_ordering);
333 }
334 
336 {
337  return & interpolation_quad;
338 }
339 
341 {
342  if ( id == P_f ) {
343  return & interpolation_lin;
344  } else {
345  return & interpolation_quad;
346  }
347 }
348 
350 {
352 }
353 
354 // Some extension Interfaces to follow:
355 
357 {
358  switch ( it ) {
360  return static_cast< NodalAveragingRecoveryModelInterface * >(this);
361 
363  return static_cast< SpatialLocalizerInterface * >(this);
364 
365  default:
366  return FMElement :: giveInterface(it);
367  }
368 }
369 
370 
371 void Hexa21Stokes :: computeField(ValueModeType mode, TimeStep *tStep, const FloatArray &lcoords, FloatArray &answer)
372 {
373  FloatArray n, n_lin, pressures, velocities;
374  this->interpolation_quad.evalN( n, lcoords, FEIElementGeometryWrapper(this) );
375  this->interpolation_lin.evalN( n_lin, lcoords, FEIElementGeometryWrapper(this) );
376  this->computeVectorOf({P_f}, mode, tStep, pressures);
377  this->computeVectorOf({V_u, V_v, V_w}, mode, tStep, velocities);
378 
379  answer.resize(4);
380  answer.zero();
381  for ( int i = 1; i <= n.giveSize(); i++ ) {
382  answer(0) += n.at(i) * velocities.at(i*3-2);
383  answer(1) += n.at(i) * velocities.at(i*3-1);
384  answer(2) += n.at(i) * velocities.at(i*3);
385  }
386 
387  for ( int i = 1; i <= n_lin.giveSize(); i++ ) {
388  answer(3) += n_lin.at(i) * pressures.at(i);
389  }
390 }
391 
392 
394 {
395  answer.resize(1);
396  if ( type == IST_Pressure ) {
397  answer.resize(1);
398  if ( node <= 8 ) { // Corner nodes
399  answer.at(1) = this->giveNode(node)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep);
400  } else if ( node <= 20 ) { // Edge nodes
401  // Edges are numbered consistently with edge nodes, so node number - 8 = edge number
402  IntArray eNodes;
403  this->interpolation_quad.computeLocalEdgeMapping(eNodes, node - 8);
404  answer.at(1) = 0.5 * (
405  this->giveNode( eNodes.at(1) )->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
406  this->giveNode( eNodes.at(2) )->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) );
407  } else if ( node <= 26 ) { // Face nodes
408  // Faces are numbered consistently with edge nodes, so node number - 12 = face number
409  IntArray fNodes;
410  this->interpolation_quad.computeLocalSurfaceMapping(fNodes, node - 20);
411  answer.at(1) = 0.25 * (
412  this->giveNode( fNodes.at(1) )->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
413  this->giveNode( fNodes.at(2) )->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
414  this->giveNode( fNodes.at(3) )->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
415  this->giveNode( fNodes.at(4) )->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) );
416  } else { // Middle node
417  answer.at(1) = 0.125 * (
418  this->giveNode(1)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
419  this->giveNode(2)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
420  this->giveNode(3)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
421  this->giveNode(4)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
422  this->giveNode(5)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
423  this->giveNode(6)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
424  this->giveNode(7)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) +
425  this->giveNode(8)->giveDofWithID(P_f)->giveUnknown(VM_Total, tStep) );
426  }
427  } else {
428  answer.clear();
429  }
430 }
431 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
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: fei3dhexaquad.C:78
void computeExternalForcesVector(FloatArray &answer, TimeStep *tStep)
Definition: hexa21stokes.C:170
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
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
Fluid cross-section.
virtual void computeLocalSurfaceMapping(IntArray &nodes, int iSurf)
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
static IntArray momentum_ordering
Ordering of momentum balance dofs in element. Used to assemble the element stiffness.
Definition: hexa21stokes.h:64
static FEI3dHexaLin interpolation_lin
Interpolation for pressure.
Definition: hexa21stokes.h:60
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
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
static IntArray surf_ordering[6]
Ordering of dofs on surfaces. Used to assemble edge loads (only momentum balance) ...
Definition: hexa21stokes.h:68
Hexa21Stokes(int n, Domain *d)
Definition: hexa21stokes.C:75
virtual double giveUnknown(ValueModeType mode, TimeStep *tStep)=0
The key method of class Dof.
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: element.C:779
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: hexa21stokes.C:233
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.
virtual FEInterpolation * giveInterpolation() const
Definition: hexa21stokes.C:335
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei3dhexalin.C:43
Distributed body load.
Definition: bcgeomtype.h:43
virtual void NodalAveragingRecoveryMI_computeNodalValue(FloatArray &answer, int node, InternalStateType type, TimeStep *tStep)
Computes the element value in given node.
Definition: hexa21stokes.C:393
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 void computeGaussPoints()
Initializes the array of integration rules member variable.
Definition: hexa21stokes.C:84
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
Definition: hexa21stokes.C:107
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...
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
static FEI3dHexaTriQuad interpolation_quad
Interpolation for geometry and velocity.
Definition: hexa21stokes.h:62
virtual void computeField(ValueModeType u, TimeStep *tStep, const FloatArray &coords, FloatArray &answer)
Computes the unknown vector interpolated at the specified local coordinates.
Definition: hexa21stokes.C:371
REGISTER_Element(LSpace)
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
virtual FormulationType giveFormulationType()
Specifies is load should take local or global coordinates.
Definition: load.h:155
Neumann type (prescribed flux).
Definition: bctype.h:43
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
Definition: hexa21stokes.C:356
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 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 bcType giveType() const
Returns receiver load type.
Definition: boundaryload.h:166
void computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
Definition: hexa21stokes.C:131
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
Definition: hexa21stokes.C:98
Distributed surface load.
Definition: bcgeomtype.h:45
virtual int setupIntegrationPoints(IntegrationRule &irule, int npoints, Element *element)
Sets up integration rule for the given element.
Definition: crosssection.C:54
virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
Computes the contribution of the given body load (volumetric).
Definition: hexa21stokes.C:200
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
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: hexa21stokes.C:349
virtual void surfaceLocal2global(FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates edge global coordinates from given local ones.
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
void computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
Definition: hexa21stokes.C:280
int giveSize() const
Returns the size of receiver.
Definition: floatarray.h:218
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
Definition: hexa21stokes.C:120
virtual void computeLocalEdgeMapping(IntArray &edgeNodes, int iedge)
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: hexa21stokes.C:93
virtual void surfaceEvalN(FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of edge interpolation functions (shape functions) at given point.
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
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 double surfaceGiveTransformationJacobian(int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the edge jacobian of transformation between local and global coordinates.
void symmetrized()
Initializes the lower half of the receiver according to the upper half.
Definition: floatmatrix.C:1576
Class representing integration point in finite element program.
Definition: gausspoint.h:93
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
IntArray boundaryLoadArray
Definition: element.h:160
Class representing solution step.
Definition: timestep.h:80
static IntArray conservation_ordering
Ordering of conservation dofs in element. Used to assemble the element stiffness. ...
Definition: hexa21stokes.h:66
int numberOfDofMans
Number of dofmanagers.
Definition: element.h:149
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
virtual bcGeomType giveBCGeoType() const
Returns receiver&#39;s load geometry type.
Definition: bodyload.h:63
virtual int SetUpPointsOnTriangle(int nPoints, MaterialMode mode)
Sets up receiver&#39;s integration points on triangular (area coords) integration domain.
virtual ~Hexa21Stokes()
Definition: hexa21stokes.C:81
Class representing Gaussian-quadrature integration rule.
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

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