OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
tet1bubblestokes.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 "tet1bubblestokes.h"
36 #include "node.h"
37 #include "domain.h"
38 #include "gaussintegrationrule.h"
39 #include "gausspoint.h"
40 #include "bcgeomtype.h"
41 #include "load.h"
42 #include "bodyload.h"
43 #include "boundaryload.h"
44 #include "mathfem.h"
45 #include "fluiddynamicmaterial.h"
46 #include "fei3dtetlin.h"
47 #include "masterdof.h"
48 #include "fluidcrosssection.h"
49 #include "assemblercallback.h"
50 #include "classfactory.h"
51 
52 namespace oofem {
53 REGISTER_Element(Tet1BubbleStokes);
54 
55 FEI3dTetLin Tet1BubbleStokes :: interp;
56 // Set up ordering vectors (for assembling)
57 IntArray Tet1BubbleStokes :: momentum_ordering = {1, 2, 3, 5, 6, 7, 9, 10, 11, 13, 14, 15, 17, 18, 19};
58 IntArray Tet1BubbleStokes :: conservation_ordering = {4, 8, 12, 16};
59 IntArray Tet1BubbleStokes :: surf_ordering [ 4 ] = {
60  {1, 2, 3, 9, 10, 11, 5, 6, 7},
61  {1, 2, 3, 5, 6, 7, 13, 14, 15},
62  {5, 6, 7, 9, 10, 11, 13, 14, 15},
63  {1, 2, 3, 13, 14, 15, 9, 10, 11}
64 };
65 
67 {
68  this->numberOfDofMans = 4;
69  this->numberOfGaussPoints = 24;
70 
71  this->bubble.reset( new ElementDofManager(1, aDomain, this) );
72  this->bubble->appendDof( new MasterDof(this->bubble.get(), V_u) );
73  this->bubble->appendDof( new MasterDof(this->bubble.get(), V_v) );
74  this->bubble->appendDof( new MasterDof(this->bubble.get(), V_w) );
75 }
76 
78 {
79 }
80 
82 {
83  if ( integrationRulesArray.size() == 0 ) {
84  integrationRulesArray.resize(1);
85  integrationRulesArray [ 0 ].reset( new GaussIntegrationRule(1, this, 1, 3) );
87  }
88 }
89 
91 {
92  return 19;
93 }
94 
95 void Tet1BubbleStokes :: giveDofManDofIDMask(int inode, IntArray &answer) const
96 {
97  answer = {V_u, V_v, V_w, P_f};
98 }
99 
101 {
102  answer = {V_u, V_v, V_w};
103 }
104 
106 {
107  double detJ = fabs( this->interp.giveTransformationJacobian( gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(this) ) );
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, N, dNv(15);
139  double r_vol, pressure;
140  FloatMatrix dN, B(6, 15);
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->interp.evaldNdx( dN, lcoords, FEIElementGeometryWrapper(this) ) );
152  this->interp.evalN( N, lcoords, FEIElementGeometryWrapper(this) );
153  double dV = detJ * gp->giveWeight();
154 
155  for ( int j = 0, k = 0; j < dN.giveNumberOfRows(); j++, k += 3 ) {
156  dNv(k + 0) = B(0, k + 0) = B(5, k + 1) = B(4, k + 2) = dN(j, 0);
157  dNv(k + 1) = B(1, k + 1) = B(5, k + 0) = B(3, k + 2) = dN(j, 1);
158  dNv(k + 2) = B(2, k + 2) = B(4, k + 0) = B(3, k + 1) = dN(j, 2);
159  }
160 
161  // Bubble contribution;
162  dNv(12) = B(0, 12) = B(5, 13) = B(4, 14) = dN(0, 0) * N(1) * N(2) * N(3) + N(0) * dN(1, 0) * N(2) * N(3) + N(0) * N(1) * dN(2, 0) * N(3) + N(0) * N(1) * N(2) * dN(3, 0);
163  dNv(13) = B(1, 13) = B(5, 12) = B(3, 14) = dN(0, 1) * N(1) * N(2) * N(3) + N(0) * dN(1, 1) * N(2) * N(3) + N(0) * N(1) * dN(2, 1) * N(3) + N(0) * N(1) * N(2) * dN(3, 1);
164  dNv(14) = B(2, 14) = B(4, 12) = B(3, 13) = dN(0, 2) * N(1) * N(2) * N(3) + N(0) * dN(1, 2) * N(2) * N(3) + N(0) * N(1) * dN(2, 2) * N(3) + N(0) * N(1) * N(2) * dN(3, 2);
165 
166  pressure = N.dotProduct(a_pressure);
167  epsp.beProductOf(B, a_velocity);
168 
169  mat->computeDeviatoricStressVector(devStress, r_vol, gp, epsp, pressure, tStep);
170 
171  momentum.plusProduct(B, devStress, dV);
172  momentum.add(-pressure * dV, dNv);
173  conservation.add(r_vol * dV, N);
174  }
175 
176  answer.resize(19);
177  answer.zero();
178  answer.assemble(momentum, this->momentum_ordering);
179  answer.assemble(conservation, this->conservation_ordering);
180 }
181 
183 {
184  FloatArray vec;
185 
186  int nLoads = this->boundaryLoadArray.giveSize() / 2;
187  answer.clear();
188  for ( int i = 1; i <= nLoads; i++ ) { // For each Neumann boundary condition
189  int load_number = this->boundaryLoadArray.at(2 * i - 1);
190  int load_id = this->boundaryLoadArray.at(2 * i);
191  Load *load = this->domain->giveLoad(load_number);
192  bcGeomType ltype = load->giveBCGeoType();
193 
194  if ( ltype == SurfaceLoadBGT ) {
195  this->computeBoundarySurfaceLoadVector(vec, static_cast< BoundaryLoad * >(load), load_id, ExternalForcesVector, VM_Total, tStep);
196  } else {
197  OOFEM_ERROR("Unsupported boundary condition: %d", load_id);
198  }
199 
200  answer.add(vec);
201  }
202 
203  BodyLoad *bload;
204  nLoads = this->giveBodyLoadArray()->giveSize();
205  for ( int i = 1; i <= nLoads; i++ ) {
206  Load *load = domain->giveLoad( bodyLoadArray.at(i) );
207  if ((bload = dynamic_cast<BodyLoad*>(load))) {
208  bcGeomType ltype = load->giveBCGeoType();
209  if ( ltype == BodyLoadBGT && load->giveBCValType() == ForceLoadBVT ) {
210  this->computeLoadVector(vec, bload, ExternalForcesVector, VM_Total, tStep);
211  answer.add(vec);
212  }
213  } else {
214  OOFEM_ERROR("Unsupported body load: %d", load);
215  }
216  }
217 }
218 
219 
221 {
222  if ( type != ExternalForcesVector ) {
223  answer.clear();
224  return;
225  }
226 
227  FloatArray N, gVector, temparray(15);
228  double dV, detJ, rho;
229 
230  // This is assumed to be the dead weight (thus multiplied by rho)
231  load->computeComponentArrayAt(gVector, tStep, VM_Total);
232  temparray.zero();
233  if ( gVector.giveSize() ) {
234  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
235  const FloatArray &lcoords = gp->giveNaturalCoordinates();
236 
237  rho = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveDensity(gp);
238  detJ = fabs( this->interp.giveTransformationJacobian( lcoords, FEIElementGeometryWrapper(this) ) );
239  dV = detJ * gp->giveWeight() * rho;
240 
241  this->interp.evalN( N, lcoords, FEIElementGeometryWrapper(this) );
242  for ( int j = 0; j < N.giveSize(); j++ ) {
243  temparray(3 * j + 0) += N(j) * rho * gVector(0) * dV;
244  temparray(3 * j + 1) += N(j) * rho * gVector(1) * dV;
245  temparray(3 * j + 2) += N(j) * rho * gVector(2) * dV;
246  }
247 
248  temparray(12) += N(0) * N(1) * N(2) * N(3) * rho * gVector(0) * dV;
249  temparray(13) += N(0) * N(1) * N(2) * N(3) * rho * gVector(1) * dV;
250  temparray(14) += N(0) * N(1) * N(2) * N(3) * rho * gVector(1) * dV;
251  }
252  }
253 
254  answer.resize(19);
255  answer.zero();
256  answer.assemble(temparray, this->momentum_ordering);
257 }
258 
259 
260  void Tet1BubbleStokes :: computeBoundarySurfaceLoadVector(FloatArray &answer, BoundaryLoad *load, int iSurf, CharType type, ValueModeType mode, TimeStep *tStep, bool global)
261 {
262  if ( type != ExternalForcesVector ) {
263  answer.clear();
264  return;
265  }
266 
267  if ( load->giveType() == TransmissionBC ) { // Neumann boundary conditions (traction)
268  int numberOfIPs = ( int ) ceil( ( load->giveApproxOrder() + 2. ) / 2. );
269 
270  GaussIntegrationRule iRule(1, this, 1, 1);
271  FloatArray N, t, f(9);
272 
273  f.zero();
274  iRule.SetUpPointsOnTriangle(numberOfIPs, _Unknown);
275 
276  for ( GaussPoint *gp: iRule ) {
277  const FloatArray &lcoords = gp->giveNaturalCoordinates();
278 
279  this->interp.surfaceEvalN( N, iSurf, lcoords, FEIElementGeometryWrapper(this) );
280  double detJ = fabs( this->interp.surfaceGiveTransformationJacobian( iSurf, lcoords, FEIElementGeometryWrapper(this) ) );
281  double dA = gp->giveWeight() * detJ;
282 
283  if ( load->giveFormulationType() == Load :: FT_Entity ) { // Edge load in xi-eta system
284  load->computeValueAt(t, tStep, lcoords, VM_Total);
285  } else { // Edge load in x-y system
286  FloatArray gcoords;
287  this->interp.boundaryLocal2Global( gcoords, iSurf, lcoords, FEIElementGeometryWrapper(this) );
288  load->computeValueAt(t, tStep, gcoords, VM_Total);
289  }
290 
291  // Reshape the vector
292  for ( int j = 0; j < N.giveSize(); j++ ) {
293  f(3 * j + 0) += N(j) * t(0) * dA;
294  f(3 * j + 1) += N(j) * t(1) * dA;
295  f(3 * j + 2) += N(j) * t(2) * dA;
296  }
297 
298  this->interp.surfaceEvalNormal( N, iSurf, lcoords, FEIElementGeometryWrapper(this) );
299  }
300 
301  answer.resize(19);
302  answer.zero();
303  answer.assemble(f, this->surf_ordering [ iSurf - 1 ]);
304  } else {
305  OOFEM_ERROR("Strange boundary condition type");
306  }
307 }
308 
310 {
311  // Note: Working with the components; [K, G+Dp; G^T+Dv^T, C] . [v,p]
312  FluidDynamicMaterial *mat = static_cast< FluidCrossSection * >( this->giveCrossSection() )->giveFluidMaterial();
313  FloatMatrix B(6, 15), EdB, K, G, Dp, DvT, C, Ed, dN;
314  FloatArray dNv(15), N, Ep, Cd, tmpA, tmpB;
315  double Cp;
316 
317  B.zero();
318 
319  for ( GaussPoint *gp: *integrationRulesArray [ 0 ] ) {
320  // Compute Gauss point and determinant at current element
321  const FloatArray &lcoords = gp->giveNaturalCoordinates();
322 
323  double detJ = fabs( this->interp.evaldNdx( dN, lcoords, FEIElementGeometryWrapper(this) ) );
324  double dV = detJ * gp->giveWeight();
325  this->interp.evalN( N, lcoords, FEIElementGeometryWrapper(this) );
326 
327  for ( int j = 0, k = 0; j < 4; j++, k += 3 ) {
328  dNv(k + 0) = B(0, k + 0) = B(5, k + 1) = B(4, k + 2) = dN(j, 0);
329  dNv(k + 1) = B(1, k + 1) = B(5, k + 0) = B(3, k + 2) = dN(j, 1);
330  dNv(k + 2) = B(2, k + 2) = B(4, k + 0) = B(3, k + 1) = dN(j, 2);
331  }
332 
333  // Bubble contribution;
334  dNv(12) = B(0, 12) = B(5, 13) = B(4, 14) = dN(0, 0) * N(1) * N(2) * N(3) + N(0) * dN(1, 0) * N(2) * N(3) + N(0) * N(1) * dN(2, 0) * N(3) + N(0) * N(1) * N(2) * dN(3, 0);
335  dNv(13) = B(1, 13) = B(5, 12) = B(3, 14) = dN(0, 1) * N(1) * N(2) * N(3) + N(0) * dN(1, 1) * N(2) * N(3) + N(0) * N(1) * dN(2, 1) * N(3) + N(0) * N(1) * N(2) * dN(3, 1);
336  dNv(14) = B(2, 14) = B(4, 12) = B(3, 13) = dN(0, 2) * N(1) * N(2) * N(3) + N(0) * dN(1, 2) * N(2) * N(3) + N(0) * N(1) * dN(2, 2) * N(3) + N(0) * N(1) * N(2) * dN(3, 2);
337 
338  // Computing the internal forces should have been done first.
339  // dsigma_dev/deps_dev dsigma_dev/dp deps_vol/deps_dev deps_vol/dp
340  mat->giveStiffnessMatrices(Ed, Ep, Cd, Cp, mode, gp, tStep);
341 
342  EdB.beProductOf(Ed, B);
343  K.plusProductSymmUpper(B, EdB, dV);
344  G.plusDyadUnsym(dNv, N, -dV);
345  C.plusDyadSymmUpper(N, Cp * dV);
346 
347  tmpA.beTProductOf(B, Ep);
348  Dp.plusDyadUnsym(tmpA, N, dV);
349 
350  tmpB.beTProductOf(B, Cd);
351  DvT.plusDyadUnsym(N, tmpB, dV);
352  }
353 
354  K.symmetrized();
355  C.symmetrized();
356  FloatMatrix GTDvT, GDp;
357  GTDvT.beTranspositionOf(G);
358  GTDvT.add(DvT);
359  GDp = G;
360  GDp.add(Dp);
361 
362  answer.resize(19, 19);
363  answer.zero();
364  answer.assemble(K, this->momentum_ordering);
365  answer.assemble(GTDvT, this->conservation_ordering, this->momentum_ordering);
366  answer.assemble(GDp, this->momentum_ordering, this->conservation_ordering);
367  answer.assemble(C, this->conservation_ordering);
368 }
369 
371 {
372  return & interp;
373 }
374 
376 {
377  return & interp;
378 }
379 
381 {
383 }
384 
385 // Some extension Interfaces to follow:
386 
388 {
389  switch ( it ) {
391  return static_cast< ZZNodalRecoveryModelInterface * >(this);
392 
394  return static_cast< SpatialLocalizerInterface * >(this);
395 
396  default:
397  return FMElement :: giveInterface(it);
398  }
399 }
400 
401 
403 {
404  FloatArray n, n_lin, pressures, velocities;
405  this->interp.evalN( n, lcoords, FEIElementGeometryWrapper(this) );
406  this->interp.evalN( n_lin, lcoords, FEIElementGeometryWrapper(this) );
407  this->computeVectorOf({P_f}, mode, tStep, pressures);
408  this->computeVectorOf({V_u, V_v, V_w}, mode, tStep, velocities);
409 
410  answer.resize(4);
411  answer.zero();
412  for ( int i = 1; i <= n.giveSize(); i++ ) {
413  answer(0) += n.at(i) * velocities.at(i*3-2);
414  answer(1) += n.at(i) * velocities.at(i*3-1);
415  answer(2) += n.at(i) * velocities.at(i*3);
416  }
417  answer(0) += n.at(1) * n.at(2) * n.at(3) * n.at(4) * velocities.at(13);
418  answer(1) += n.at(1) * n.at(2) * n.at(3) * n.at(4) * velocities.at(14);
419  answer(2) += n.at(1) * n.at(2) * n.at(3) * n.at(4) * velocities.at(15);
420 
421  for ( int i = 1; i <= n_lin.giveSize(); i++ ) {
422  answer(3) += n_lin.at(i) * pressures.at(i);
423  }
424 }
425 
426 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
virtual Interface * giveInterface(InterfaceType it)
Interface requesting service.
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.
static IntArray momentum_ordering
Ordering of dofs in element. Used to assemble the element stiffness.
Abstract base class for all fluid materials.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
void computeInternalForcesVector(FloatArray &answer, TimeStep *tStep)
virtual void computeField(ValueModeType u, TimeStep *tStep, const FloatArray &coords, FloatArray &answer)
Computes the unknown vector interpolated at the specified local coordinates.
Class implementing internal element dof manager having some DOFs.
The element interface required by ZZNodalRecoveryModel.
void computeExternalForcesVector(FloatArray &answer, TimeStep *tStep)
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
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Definition: element.C:779
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the array of interpolation functions (shape functions) at given point.
Definition: fei3dtetlin.C:43
virtual void giveDofManDofIDMask(int inode, IntArray &answer) const
Returns dofmanager dof mask for node.
static IntArray surf_ordering[4]
Ordering of dofs on surfaces. Used to assemble surface loads.
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 void giveInternalDofManDofIDMask(int i, IntArray &answer) const
Returns internal dofmanager dof mask for node.
virtual double giveTransformationJacobian(const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the determinant of the transformation.
Definition: fei3dtetlin.C:178
Distributed body load.
Definition: bcgeomtype.h:43
virtual void boundaryLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Maps the local boundary coordinates to global.
Definition: feinterpol3d.C:86
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
Class representing "master" degree of freedom.
Definition: masterdof.h:92
virtual void computeVectorOfPressures(ValueModeType mode, TimeStep *tStep, FloatArray &pressures)
Definition: fmelement.C:55
virtual void updateYourself(TimeStep *tStep)
Updates element state after equilibrium in time step has been reached.
Tet1BubbleStokes(int n, Domain *d)
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 int computeNumberOfDofs()
Computes or simply returns total number of element&#39;s local DOFs.
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
virtual double surfaceEvalNormal(FloatArray &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the normal out of the surface at given point.
Definition: fei3dtetlin.C:366
#define OOFEM_ERROR(...)
Definition: error.h:61
REGISTER_Element(LSpace)
static FEI3dTetLin interp
Interpolation for pressure.
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
virtual double giveWeight()
Returns integration weight of receiver.
Definition: gausspoint.h:181
virtual void giveCharacteristicVector(FloatArray &answer, CharType type, ValueModeType mode, TimeStep *tStep)
Computes characteristic vector of receiver of requested type in given time step.
virtual FormulationType giveFormulationType()
Specifies is load should take local or global coordinates.
Definition: load.h:155
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Neumann type (prescribed flux).
Definition: bctype.h:43
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
virtual void giveCharacteristicMatrix(FloatMatrix &answer, CharType type, TimeStep *tStep)
Computes characteristic matrix of receiver of requested type in given time step.
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
static IntArray conservation_ordering
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: fei3dtetlin.C:54
virtual double surfaceGiveTransformationJacobian(int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the edge jacobian of transformation between local and global coordinates.
Definition: fei3dtetlin.C:380
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 FEInterpolation * giveInterpolation() const
virtual bcType giveType() const
Returns receiver load type.
Definition: boundaryload.h:166
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
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 computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode mode, TimeStep *tStep)
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
virtual void computeLoadVector(FloatArray &answer, BodyLoad *load, CharType type, ValueModeType mode, TimeStep *tStep)
Computes the contribution of the given body load (volumetric).
Implementation of matrix containing floating point numbers.
Definition: floatmatrix.h:94
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.
Definition: fei3dtetlin.C:309
std::unique_ptr< ElementDofManager > bubble
The extra dofs from the bubble.
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
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...
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
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
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
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 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
void symmetrized()
Initializes the lower half of the receiver according to the upper half.
Definition: floatmatrix.C:1576
virtual void computeGaussPoints()
Initializes the array of integration rules member variable.
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
const FloatArray & giveNaturalCoordinates()
Returns coordinate array of receiver.
Definition: gausspoint.h:138
virtual int SetUpPointsOnTriangle(int nPoints, MaterialMode mode)
Sets up receiver&#39;s integration points on triangular (area coords) integration domain.
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:31 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011