OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
surfacetensionbc.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 "surfacetensionbc.h"
36 #include "element.h"
37 #include "node.h"
38 #include "crosssection.h"
39 #include "error.h"
40 #include "feinterpol.h"
41 #include "feinterpol2d.h"
42 #include "feinterpol3d.h"
43 #include "classfactory.h"
44 #include "set.h"
45 #include "sparsemtrx.h"
46 
47 #include "integrationdomain.h"
48 #include "integrationrule.h"
49 #include "gausspoint.h"
50 #include "mathfem.h"
51 
52 #include <utility>
53 #include <list>
54 #include <memory>
55 
56 namespace oofem {
57 REGISTER_BoundaryCondition(SurfaceTensionBoundaryCondition);
58 
60 {
61  IRResultType result;
62 
64 
66 
68 }
69 
70 void SurfaceTensionBoundaryCondition :: giveLocationArrays(std :: vector< IntArray > &rows, std :: vector< IntArray > &cols, CharType type,
71  const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s)
72 {
73  if ( !this->useTangent || type != TangentStiffnessMatrix ) {
74  return;
75  }
76 
77  Set *set = this->giveDomain()->giveSet(this->set);
78  const IntArray &boundaries = set->giveBoundaryList();
79  IntArray bNodes;
80 
81  rows.resize(boundaries.giveSize() / 2);
82  cols.resize(boundaries.giveSize() / 2);
83 
84  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
85  Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
86  int boundary = boundaries.at(pos * 2);
87 
88  e->giveInterpolation()->boundaryGiveNodes(bNodes, boundary);
89 
90  e->giveBoundaryLocationArray(rows [ pos ], bNodes, this->dofs, r_s);
91  e->giveBoundaryLocationArray(cols [ pos ], bNodes, this->dofs, c_s);
92  }
93 }
94 
96  CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, double scale)
97 {
98  if ( !this->useTangent || type != TangentStiffnessMatrix ) {
99  return;
100  }
101 
102  OOFEM_ERROR("Not implemented yet.");
103 
104  FloatMatrix Ke;
105  IntArray r_loc, c_loc, bNodes;
106 
107  Set *set = this->giveDomain()->giveSet(this->set);
108  const IntArray &boundaries = set->giveBoundaryList();
109 
110  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
111  Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
112  int boundary = boundaries.at(pos * 2);
113 
114  e->giveInterpolation()->boundaryGiveNodes(bNodes, boundary);
115 
116  e->giveBoundaryLocationArray(r_loc, bNodes, this->dofs, r_s);
117  e->giveBoundaryLocationArray(c_loc, bNodes, this->dofs, c_s);
118  this->computeTangentFromElement(Ke, e, boundary, tStep);
119  Ke.times(scale);
120  answer.assemble(r_loc, c_loc, Ke);
121  }
122 }
123 
125  CharType type, ValueModeType mode,
126  const UnknownNumberingScheme &s, FloatArray *eNorms)
127 {
128  if ( type != ExternalForcesVector ) {
129  return;
130  }
131 
132  FloatArray fe;
133  IntArray loc, masterdofids, bNodes;
134 
135  Set *set = this->giveDomain()->giveSet(this->set);
136  const IntArray &boundaries = set->giveBoundaryList();
137 
138  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
139  Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
140  int boundary = boundaries.at(pos * 2);
141 
142  e->giveInterpolation()->boundaryGiveNodes(bNodes, boundary);
143 
144  e->giveBoundaryLocationArray(loc, bNodes, this->dofs, s, & masterdofids);
145  this->computeLoadVectorFromElement(fe, e, boundary, tStep);
146  answer.assemble(fe, loc);
147  if ( eNorms ) {
148  eNorms->assembleSquared(fe, masterdofids);
149  }
150  }
151 }
152 
154 {
156  if ( !fei ) {
157  OOFEM_ERROR("No interpolation available for element.");
158  }
159  std :: unique_ptr< IntegrationRule > iRule( fei->giveBoundaryIntegrationRule(fei->giveInterpolationOrder()-1, side) );
160 
161  int nsd = e->giveDomain()->giveNumberOfSpatialDimensions();
162  int nodes = e->giveNumberOfNodes();
163  if ( side == -1 ) {
164  side = 1;
165  }
166 
167  answer.clear();
168 
169  if ( nsd == 2 ) {
170  FEInterpolation2d *fei2d = static_cast< FEInterpolation2d * >(fei);
171 
173  FloatMatrix xy(2, nodes);
174  Node *node;
175  for ( int i = 1; i <= nodes; i++ ) {
176  node = e->giveNode(i);
177  xy.at(1, i) = node->giveCoordinate(1);
178  xy.at(2, i) = node->giveCoordinate(2);
179  }
180 
181  FloatArray tmpA(2 *nodes);
182  FloatArray es; // Tangent vector to curve
183  FloatArray dNds;
184  FloatMatrix B(2, 2 *nodes);
185  B.zero();
186 
187  if ( e->giveDomain()->isAxisymmetric() ) {
188  FloatArray N;
189  FloatArray gcoords;
190  FloatArray tmpB(2 *nodes);
191  for ( GaussPoint *gp: *iRule ) {
192  fei2d->edgeEvaldNds( dNds, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
193  fei->boundaryEvalN( N, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
194  double J = fei->boundaryGiveTransformationJacobian( side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
195  fei->boundaryLocal2Global( gcoords, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
196  double r = gcoords(0); // First coordinate is the radial coord.
197 
198  es.beProductOf(xy, dNds);
199 
200  // Construct the different matrices in the integrand;
201  for ( int i = 0; i < nodes; i++ ) {
202  tmpA(i * 2 + 0) = dNds(i) * es(0);
203  tmpA(i * 2 + 1) = dNds(i) * es(1);
204  tmpB(i * 2 + 0) = N(i);
205  tmpB(i * 2 + 1) = 0;
206  B(i * 2, 0) = B(i * 2 + 1, 1) = dNds(i);
207  }
208 
209  double dV = 2 *M_PI *gamma *J *gp->giveWeight();
210  answer.plusDyadUnsym(tmpA, tmpB, dV);
211  answer.plusDyadUnsym(tmpB, tmpA, dV);
212  answer.plusProductSymmUpper(B, B, r * dV);
213  answer.plusDyadUnsym(tmpA, tmpA, -r * dV);
214  }
215  } else {
216  for ( GaussPoint *gp: *iRule ) {
217  double t = e->giveCrossSection()->give(CS_Thickness, gp);
218  fei2d->edgeEvaldNds( dNds, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
219  double J = fei->boundaryGiveTransformationJacobian( side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
220 
221  es.beProductOf(xy, dNds);
222 
223  // Construct the different matrices in the integrand;
224  for ( int i = 0; i < nodes; i++ ) {
225  tmpA(i * 2 + 0) = dNds(i) * es(0);
226  tmpA(i * 2 + 1) = dNds(i) * es(1);
227  B(i * 2, 0) = B(i * 2 + 1, 1) = dNds(i);
228  }
229 
230  double dV = t * gamma * J * gp->giveWeight();
231  answer.plusProductSymmUpper(B, B, dV);
232  answer.plusDyadSymmUpper(tmpA, -dV);
233  }
234  }
235 
236  answer.symmetrized();
237  } else if ( nsd == 3 ) {
238 
239  FEInterpolation3d *fei3d = static_cast< FEInterpolation3d * >(fei);
240 
241  OOFEM_ERROR("3D tangents not implemented yet.");
242 
243  //FloatMatrix tmp(3 *nodes, 3 *nodes);
244  FloatMatrix dNdx;
245  FloatArray n;
246  for ( GaussPoint *gp: *iRule ) {
247  fei3d->surfaceEvaldNdx( dNdx, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
248  /*double J = */ fei->boundaryEvalNormal( n, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
249  //double dV = gamma * J * gp->giveWeight();
250 
251  for ( int i = 0; i < nodes; i++ ) {
252  //tmp(3*i+0) = dNdx(i,0) - (dNdx(i,0)*n(0)* + dNdx(i,1)*n(1) + dNdx(i,2)*n(2))*n(0);
253  //tmp(3*i+1) = dNdx(i,1) - (dNdx(i,0)*n(0)* + dNdx(i,1)*n(1) + dNdx(i,2)*n(2))*n(1);
254  //tmp(3*i+2) = dNdx(i,2) - (dNdx(i,0)*n(0)* + dNdx(i,1)*n(1) + dNdx(i,2)*n(2))*n(2);
255  }
256  //answer.plusProductSymmUpper(A,B, dV);
258  }
259  } else {
260  OOFEM_WARNING("Only 2D or 3D is possible!");
261  }
262 }
263 
265 {
267  if ( !fei ) {
268  OOFEM_ERROR("No interpolation or default integration available for element.");
269  }
270  std :: unique_ptr< IntegrationRule > iRule( fei->giveBoundaryIntegrationRule(fei->giveInterpolationOrder()-1, side) );
271 
272  int nsd = e->giveDomain()->giveNumberOfSpatialDimensions();
273 
274  if ( side == -1 ) {
275  side = 1;
276  }
277 
278  answer.clear();
279 
280  if ( nsd == 2 ) {
281 
282  FEInterpolation2d *fei2d = static_cast< FEInterpolation2d * >(fei);
283 
285  IntArray bnodes;
286  fei2d->boundaryGiveNodes(bnodes, side);
287  int nodes = bnodes.giveSize();
288  FloatMatrix xy(2, nodes);
289  for ( int i = 1; i <= nodes; i++ ) {
290  Node *node = e->giveNode(bnodes.at(i));
292  xy.at(1, i) = node->giveCoordinate(1);
293  xy.at(2, i) = node->giveCoordinate(2);
294  }
295 
296  FloatArray tmp; // Integrand
297  FloatArray es; // Tangent vector to curve
298  FloatArray dNds;
299 
300  if ( e->giveDomain()->isAxisymmetric() ) {
301  FloatArray N;
302  FloatArray gcoords;
303  for ( GaussPoint *gp: *iRule ) {
304  fei2d->edgeEvaldNds( dNds, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
305  fei->boundaryEvalN( N, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
306  double J = fei->boundaryGiveTransformationJacobian( side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
307  fei->boundaryLocal2Global( gcoords, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
308  double r = gcoords(0); // First coordinate is the radial coord.
309 
310  es.beProductOf(xy, dNds);
311 
312  tmp.resize( 2 * nodes);
313  for ( int i = 0; i < nodes; i++ ) {
314  tmp(2 * i) = dNds(i) * es(0) * r + N(i);
315  tmp(2 * i + 1) = dNds(i) * es(1) * r;
316  }
317 
318  answer.add(- 2 * M_PI * gamma * J * gp->giveWeight(), tmp);
319  }
320  } else {
321  for ( GaussPoint *gp: *iRule ) {
322  double t = e->giveCrossSection()->give(CS_Thickness, gp);
323  fei2d->edgeEvaldNds( dNds, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
324  double J = fei->boundaryGiveTransformationJacobian( side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
325  es.beProductOf(xy, dNds);
326 
327  tmp.resize( 2 * nodes);
328  for ( int i = 0; i < nodes; i++ ) {
329  tmp(2 * i) = dNds(i) * es(0);
330  tmp(2 * i + 1) = dNds(i) * es(1);
331  //B.at(1, 1+i*2) = B.at(2, 2+i*2) = dNds(i);
332  }
333  //tmp.beTProductOf(B, es);
334 
335  answer.add(- t * gamma * J * gp->giveWeight(), tmp);
336  }
337  }
338  } else if ( nsd == 3 ) {
339 
340  FEInterpolation3d *fei3d = static_cast< FEInterpolation3d * >(fei);
341  FloatArray n, surfProj;
342  FloatMatrix dNdx, B;
343  for ( GaussPoint *gp: *iRule ) {
344  fei3d->surfaceEvaldNdx( dNdx, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
345  double J = fei->boundaryEvalNormal( n, side, gp->giveNaturalCoordinates(), FEIElementGeometryWrapper(e) );
346 
347  // [I - n(x)n] in voigt form:
348  surfProj = {1. - n(0)*n(0), 1. - n(1)*n(1), 1. - n(2)*n(2),
349  - n(1)*n(2), - n(0)*n(2), - n(0)*n(1),
350  };
351 
352  // Construct B matrix of the surface nodes
353  B.resize(6, 3 * dNdx.giveNumberOfRows());
354  for ( int i = 1; i <= dNdx.giveNumberOfRows(); i++ ) {
355  B.at(1, 3 * i - 2) = dNdx.at(i, 1);
356  B.at(2, 3 * i - 1) = dNdx.at(i, 2);
357  B.at(3, 3 * i - 0) = dNdx.at(i, 3);
358 
359  B.at(5, 3 * i - 2) = B.at(4, 3 * i - 1) = dNdx.at(i, 3);
360  B.at(6, 3 * i - 2) = B.at(4, 3 * i - 0) = dNdx.at(i, 2);
361  B.at(6, 3 * i - 1) = B.at(5, 3 * i - 0) = dNdx.at(i, 1);
362  }
363 
364  answer.plusProduct(B, surfProj, -gamma * J * gp->giveWeight() );
365  }
366  }
367 }
368 } // end namespace oofem
CrossSection * giveCrossSection()
Definition: element.C:495
#define _IFT_SurfaceTensionBoundaryCondition_gamma
virtual void assembleVector(FloatArray &answer, TimeStep *tStep, CharType type, ValueModeType mode, const UnknownNumberingScheme &s, FloatArray *eNorms=NULL)
Assembles B.C.
REGISTER_BoundaryCondition(BoundaryCondition)
virtual int assemble(const IntArray &loc, const FloatMatrix &mat)=0
Assembles sparse matrix from contribution of local elements.
Base class for all matrices stored in sparse format.
Definition: sparsemtrx.h:60
IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Definition: activebc.h:75
virtual void boundaryGiveNodes(IntArray &answer, int boundary)
Gives the boundary nodes for requested boundary number.
Definition: feinterpol2d.C:139
virtual double boundaryEvalNormal(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the normal on the requested boundary.
ValueModeType
Type representing the mode of UnknownType or CharType, or similar types.
Definition: valuemodetype.h:78
Class representing a general abstraction for surface finite element interpolation class...
Definition: feinterpol2d.h:45
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
int giveInterpolationOrder()
Returns the interpolation order.
Definition: feinterpol.h:159
bool isAxisymmetric()
Returns true of axisymmetry is in effect.
Definition: domain.C:1074
virtual void assemble(SparseMtrx &answer, TimeStep *tStep, CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, double scale=1.0)
Assembles B.C.
virtual bool hasField(InputFieldType id)=0
Returns true if record contains field identified by idString keyword.
Abstract base class for all finite elements.
Definition: element.h:145
void plusDyadUnsym(const FloatArray &a, const FloatArray &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:811
virtual double giveCoordinate(int i)
Definition: node.C:82
int giveNumberOfSpatialDimensions()
Returns number of spatial dimensions.
Definition: domain.C:1067
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
virtual FEInterpolation * giveInterpolation() const
Definition: element.h:629
void computeLoadVectorFromElement(FloatArray &answer, Element *e, int side, TimeStep *tStep)
Helper function for computing the contributions to the load vector.
virtual int giveNumberOfNodes() const
Returns number of nodes of receiver.
Definition: element.h:662
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
Element * giveElement(int n)
Service for accessing particular domain fe element.
Definition: domain.C:160
#define M_PI
Definition: mathfem.h:52
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual void boundaryLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Maps the local boundary coordinates to global.
void times(double f)
Multiplies receiver by factor f.
Definition: floatmatrix.C:1594
Set of elements, boundaries, edges and/or nodes.
Definition: set.h:66
Abstract base class allowing to control the way, how equations are assigned to individual DOFs...
Wrapper around element definition to provide FEICellGeometry interface.
Definition: feinterpol.h:95
Set * giveSet(int n)
Service for accessing particular domain set.
Definition: domain.C:363
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
virtual double boundaryGiveTransformationJacobian(int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the determinant of the transformation Jacobian on the requested boundary.
#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
IntArray dofs
Dofs that b.c. is applied to (relevant for Dirichlet type b.c.s).
double at(int i, int j) const
Coefficient access function.
Definition: floatmatrix.h:176
void resize(int n)
Checks size of receiver towards requested bounds.
Definition: intarray.C:124
virtual void scale(double s)
Scales the receiver according to given value.
virtual void boundaryGiveNodes(IntArray &answer, int boundary)=0
Gives the boundary nodes for requested boundary number.
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
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 void surfaceEvaldNdx(FloatMatrix &answer, int isurf, const FloatArray &lcoords, const FEICellGeometry &cellgeo)
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point...
Definition: feinterpol3d.C:128
virtual void giveLocationArrays(std::vector< IntArray > &rows, std::vector< IntArray > &cols, CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s)
Gives a list of location arrays that will be assembled.
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
Class representing the general Input Record.
Definition: inputrecord.h:101
Class representing a general abstraction for surface finite element interpolation class...
Definition: feinterpol3d.h:44
virtual IntegrationRule * giveBoundaryIntegrationRule(int order, int boundary)
Sets up a suitable integration rule for integrating over the requested boundary.
Definition: feinterpol.C:63
void computeTangentFromElement(FloatMatrix &answer, Element *e, int side, TimeStep *tStep)
Helper function for computing the tangent ( )
void assembleSquared(const FloatArray &fe, const IntArray &loc)
Assembles the array fe with each component squared.
Definition: floatarray.C:572
virtual void giveBoundaryLocationArray(IntArray &locationArray, const IntArray &bNodes, const UnknownNumberingScheme &s, IntArray *dofIds=NULL)
Returns the location array for the boundary of the element.
Definition: element.C:446
virtual void edgeEvaldNds(FloatArray &answer, int iedge, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the matrix of derivatives of edge interpolation functions (shape functions) at given point...
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
Domain * giveDomain() const
Definition: femcmpnn.h:100
#define _IFT_SurfaceTensionBoundaryCondition_useTangent
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
Class implementing node in finite element mesh.
Definition: node.h:87
#define IR_GIVE_FIELD(__ir, __value, __id)
Macro facilitating the use of input record reading methods.
Definition: inputrecord.h:69
Node * giveNode(int i) const
Returns reference to the i-th node of element.
Definition: element.h:610
bool useTangent
Determines if tangent should be used.
int giveNumberOfRows() const
Returns number of rows of receiver.
Definition: floatmatrix.h:156
void symmetrized()
Initializes the lower half of the receiver according to the upper half.
Definition: floatmatrix.C:1576
virtual void boundaryEvalN(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the basis functions on the requested boundary.
Class representing integration point in finite element program.
Definition: gausspoint.h:93
#define OOFEM_WARNING(...)
Definition: error.h:62
Class representing solution step.
Definition: timestep.h:80
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
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