OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
prescribedgradientbcneumann.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 
36 #include "classfactory.h"
37 #include "node.h"
38 #include "masterdof.h"
39 #include "element.h"
40 #include "feinterpol.h"
41 #include "feinterpol2d.h"
42 #include "gausspoint.h"
43 #include "sparsemtrx.h"
46 #include "timestep.h"
47 #include "function.h"
48 #include "sparselinsystemnm.h"
49 #include "unknownnumberingscheme.h"
50 #include "engngm.h"
51 #include "mathfem.h"
52 
53 namespace oofem {
54 REGISTER_BoundaryCondition(PrescribedGradientBCNeumann);
55 
59  mpSigmaHom( new Node(0, d) )
60 {
61  int nsd = d->giveNumberOfSpatialDimensions();
62  for ( int i = 0; i < nsd * nsd; i++ ) {
63  // Just putting in X_i id-items since they don't matter.
64  int dofId = d->giveNextFreeDofID();
65  mSigmaIds.followedBy(dofId);
66  mpSigmaHom->appendDof( new MasterDof( mpSigmaHom.get(), ( DofIDItem ) ( dofId ) ) );
67  }
68 }
69 
71 {
72 }
73 
74 
76 {
79 }
80 
81 
83 {
86 }
87 
88 
90 {
91  return mpSigmaHom.get();
92 }
93 
95 {
96  this->mGradient.times(s);
97 }
98 
100  CharType type, ValueModeType mode,
101  const UnknownNumberingScheme &s, FloatArray *eNorm)
102 {
103  IntArray sigma_loc; // For the displacements and stress respectively
104  mpSigmaHom->giveLocationArray(mSigmaIds, sigma_loc, s);
105 
106  if ( type == ExternalForcesVector ) {
107  // The external forces have two contributions. On the additional equations for sigma, the load is simply the prescribed gradient.
108  double rve_size = this->domainSize(this->giveDomain(), this->giveSetNumber());
109  FloatArray stressLoad;
110  FloatArray gradVoigt;
111  giveGradientVoigt(gradVoigt);
112 
113  double loadLevel = this->giveTimeFunction()->evaluateAtTime(tStep->giveTargetTime());
114  stressLoad.beScaled(-rve_size*loadLevel, gradVoigt);
115 
116  answer.assemble(stressLoad, sigma_loc);
117  } else if ( type == InternalForcesVector ) {
118  FloatMatrix Ke;
119  FloatArray fe_v, fe_s;
120  FloatArray sigmaHom, e_u;
121  IntArray loc, masterDofIDs, sigmaMasterDofIDs;
122 
123  // Fetch the current values of the stress;
124  mpSigmaHom->giveUnknownVector(sigmaHom, mSigmaIds, mode, tStep);
125  // and the master dof ids for sigmadev used for the internal norms
126  mpSigmaHom->giveMasterDofIDArray(mSigmaIds, sigmaMasterDofIDs);
127 
128  // Assemble
129  Set *setPointer = this->giveDomain()->giveSet(this->set);
130  const IntArray &boundaries = setPointer->giveBoundaryList();
131  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
132  Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
133  int boundary = boundaries.at(pos * 2);
134 
135  // Fetch the element information;
136  e->giveLocationArray(loc, s, & masterDofIDs);
137  // Here, we could use only the nodes actually located on the boundary, but we don't.
138  // Instead, we use all nodes belonging to the element, which is allowed because the
139  // basis functions related to the interior nodes will be zero on the boundary.
140  // Obviously, this is less efficient, so why do we want to do it this way?
141  // Because it is easier when XFEM enrichments are present. /ES
142  e->computeVectorOf(mode, tStep, e_u);
143  this->integrateTangent(Ke, e, boundary);
144 
145  // We just use the tangent, less duplicated code (the addition of sigmaDev is linear).
146  fe_v.beProductOf(Ke, e_u);
147  fe_s.beTProductOf(Ke, sigmaHom);
148 
149  // Note: The terms appear negative in the equations:
150  fe_v.negated();
151  fe_s.negated();
152 
153  answer.assemble(fe_s, loc); // Contributions to delta_v equations
154  answer.assemble(fe_v, sigma_loc); // Contribution to delta_s_i equations
155  if ( eNorm != NULL ) {
156  eNorm->assembleSquared(fe_s, masterDofIDs);
157  eNorm->assembleSquared(fe_v, sigmaMasterDofIDs);
158  }
159  }
160  }
161 }
162 
164  CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, double scale)
165 {
166  if ( type == TangentStiffnessMatrix || type == SecantStiffnessMatrix || type == ElasticStiffnessMatrix ) {
167  FloatMatrix Ke, KeT;
168  IntArray loc_r, loc_c, sigma_loc_r, sigma_loc_c;
169  Set *set = this->giveDomain()->giveSet(this->set);
170  const IntArray &boundaries = set->giveBoundaryList();
171 
172  // Fetch the columns/rows for the stress contributions;
173  mpSigmaHom->giveLocationArray(mSigmaIds, sigma_loc_r, r_s);
174  mpSigmaHom->giveLocationArray(mSigmaIds, sigma_loc_c, c_s);
175 
176  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
177  Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
178  int boundary = boundaries.at(pos * 2);
179 
180  // Here, we could use only the nodes actually located on the boundary, but we don't.
181  // Instead, we use all nodes belonging to the element, which is allowed because the
182  // basis functions related to the interior nodes will be zero on the boundary.
183  // Obviously, this is less efficient, so why do we want to do it this way?
184  // Because it is easier when XFEM enrichments are present. /ES
185  e->giveLocationArray(loc_r, r_s);
186  e->giveLocationArray(loc_c, c_s);
187 
188  this->integrateTangent(Ke, e, boundary);
189  Ke.negated();
190  Ke.times(scale);
191  KeT.beTranspositionOf(Ke);
192 
193  answer.assemble(sigma_loc_r, loc_c, Ke); // Contribution to delta_s_i equations
194  answer.assemble(loc_r, sigma_loc_c, KeT); // Contributions to delta_v equations
195  }
196  } else {
197  OOFEM_LOG_DEBUG("Skipping assembly in PrescribedGradientBCNeumann::assemble().");
198  }
199 }
200 
201 void PrescribedGradientBCNeumann :: giveLocationArrays(std :: vector< IntArray > &rows, std :: vector< IntArray > &cols, CharType type,
202  const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s)
203 {
204  IntArray loc_r, loc_c, sigma_loc_r, sigma_loc_c;
205 
206  // Fetch the columns/rows for the stress contributions;
207  mpSigmaHom->giveLocationArray(mSigmaIds, sigma_loc_r, r_s);
208  mpSigmaHom->giveLocationArray(mSigmaIds, sigma_loc_c, c_s);
209 
210  Set *set = this->giveDomain()->giveSet(this->set);
211  const IntArray &boundaries = set->giveBoundaryList();
212 
213  rows.resize( boundaries.giveSize() );
214  cols.resize( boundaries.giveSize() );
215  int i = 0;
216  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
217  Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
218 
219  // Here, we could use only the nodes actually located on the boundary, but we don't.
220  // Instead, we use all nodes belonging to the element, which is allowed because the
221  // basis functions related to the interior nodes will be zero on the boundary.
222  // Obviously, this is less efficient, so why do we want to do it this way?
223  // Because it is easier when XFEM enrichments are present. /ES
224  e->giveLocationArray(loc_r, r_s);
225  e->giveLocationArray(loc_c, c_s);
226 
227  // For most uses, loc_r == loc_c, and sigma_loc_r == sigma_loc_c.
228  rows [ i ] = loc_r;
229  cols [ i ] = sigma_loc_c;
230  i++;
231  // and the symmetric part (usually the transpose of above)
232  rows [ i ] = sigma_loc_r;
233  cols [ i ] = loc_c;
234  i++;
235  }
236 }
237 
239 {
240  mpSigmaHom->giveUnknownVector(sigma, mSigmaIds, VM_Total, tStep);
241 }
242 
243 
245 {
246  EngngModel *rve = this->giveDomain()->giveEngngModel();
248  std :: unique_ptr< SparseLinearSystemNM > solver(
249  classFactory.createSparseLinSolver( ST_Petsc, this->domain, this->domain->giveEngngModel() ) ); // = rve->giveLinearSolver();
250  SparseMtrxType stype = solver->giveRecommendedMatrix(true);
251  double rve_size = this->domainSize(this->giveDomain(), this->giveSetNumber());
252 
253  // 1. Kuu*us = -Kus*s => us = -Kuu\Ku where u = us*s
254  // 2. Ks = Kus'*us
255  // 3. Ks*lambda = I
256 
257  // 1.
258  // This is not very good. We have to keep Kuu and Kff in memory at the same time. Not optimal
259  // Consider changing this approach.
261  std :: unique_ptr< SparseMtrx > Kff( classFactory.createSparseMtrx(stype) );
262  if ( !Kff ) {
263  OOFEM_ERROR("Couldn't create sparse matrix of type %d\n", stype);
264  }
265  Kff->buildInternalStructure(rve, this->domain->giveNumber(), fnum);
266 
267  rve->assemble(*Kff, tStep, TangentAssembler(TangentStiffness), fnum, this->domain);
268 
269  IntArray loc_u, loc_s;
270  this->mpSigmaHom->giveLocationArray(this->mSigmaIds, loc_s, fnum);
271  int neq = Kff->giveNumberOfRows();
272  loc_u.resize(neq - loc_s.giveSize());
273  int k = 0;
274  for ( int i = 1; i <= neq; ++i ) {
275  if ( !loc_s.contains(i) ) {
276  loc_u.at(k++) = i;
277  }
278  }
279 
280  std :: unique_ptr< SparseMtrx > Kuu(Kff->giveSubMatrix(loc_u, loc_u));
281  // NOTE: Kus is actually a dense matrix, but we have to make it a dense matrix first
282  std :: unique_ptr< SparseMtrx > Kus(Kff->giveSubMatrix(loc_u, loc_s));
283  FloatMatrix eye(Kus->giveNumberOfColumns(), Kus->giveNumberOfColumns());
284  eye.beUnitMatrix();
285  FloatMatrix KusD;
286  Kus->times(KusD, eye);
287 
288  // Release a large chunk of redundant memory early.
289  Kus.reset();
290  Kff.reset();
291 
292  // 1.
293  FloatMatrix us;
294  solver->solve(*Kuu, KusD, us);
295  us.negated();
296 
297  // 2.
298  FloatMatrix Ks;
299  Ks.beTProductOf(KusD, us);
300  Kus->times(Ks, us);
301 
302  // 3.
303  tangent.beInverseOf(Ks);
304  tangent.times(rve_size);
305 }
306 
307 
309 {
310  mpSigmaHom->giveLocationArray(mSigmaIds, oCols, r_s);
311 }
312 
314 {
315  FloatArray normal, n;
316  FloatMatrix nMatrix, E_n;
317  FloatMatrix contrib;
318 
319  Domain *domain = e->giveDomain();
320 
321  FEInterpolation *interp = e->giveInterpolation(); // Geometry interpolation
322 
323  int nsd = e->giveDomain()->giveNumberOfSpatialDimensions();
324 
325  // Interpolation order
326  int order = interp->giveInterpolationOrder();
327  std :: unique_ptr< IntegrationRule > ir;
328 
329  XfemElementInterface *xfemElInt = dynamic_cast< XfemElementInterface * >( e );
330  if ( xfemElInt != NULL && domain->hasXfemManager() ) {
331  IntArray edgeNodes;
332  FEInterpolation2d *interp2d = dynamic_cast< FEInterpolation2d * >( interp );
333  if ( interp2d == NULL ) {
334  OOFEM_ERROR("failed to cast to FEInterpolation2d.")
335  }
336  interp2d->computeLocalEdgeMapping(edgeNodes, iBndIndex);
337 
338 // const FloatArray &xS = * ( e->giveDofManager( edgeNodes.at(1) )->giveCoordinates() );
339 // const FloatArray &xE = * ( e->giveDofManager( edgeNodes.at( edgeNodes.giveSize() ) )->giveCoordinates() );
340 
341  std :: vector< Line >segments;
342  std :: vector< FloatArray >intersecPoints;
343  xfemElInt->partitionEdgeSegment(iBndIndex, segments, intersecPoints);
344  MaterialMode matMode = e->giveMaterialMode();
345  ir.reset( new DiscontinuousSegmentIntegrationRule(1, e, segments) );
346  int numPointsPerSeg = 1;
347  ir->SetUpPointsOnLine(numPointsPerSeg, matMode);
348  } else {
349  ir.reset( interp->giveBoundaryIntegrationRule(order, iBndIndex) );
350  }
351 
352  oTangent.clear();
353 
354  for ( auto &gp: *ir ) {
355  const FloatArray &lcoords = gp->giveNaturalCoordinates();
356  FEIElementGeometryWrapper cellgeo(e);
357 
358  // Evaluate the normal;
359  double detJ = interp->boundaryEvalNormal(normal, iBndIndex, lcoords, cellgeo);
360 
361  interp->boundaryEvalN(n, iBndIndex, lcoords, cellgeo);
362  // If cracks cross the edge, special treatment is necessary.
363  // Exploit the XfemElementInterface to minimize duplication of code.
364  if ( xfemElInt != NULL && domain->hasXfemManager() ) {
365  // Compute global coordinates of Gauss point
366  FloatArray globalCoord;
367 
368  interp->boundaryLocal2Global(globalCoord, iBndIndex, lcoords, cellgeo);
369 
370  // Compute local coordinates on the element
371  FloatArray locCoord;
372  e->computeLocalCoordinates(locCoord, globalCoord);
373 
374  xfemElInt->XfemElementInterface_createEnrNmatrixAt(nMatrix, locCoord, * e, false);
375  } else {
376  // Evaluate the velocity/displacement coefficients
377  nMatrix.beNMatrixOf(n, nsd);
378  }
379 
380  E_n.resize(nsd*nsd, nsd);
381  E_n.zero();
382 
383  if ( nsd == 3 ) {
384  E_n.at(1, 1) = normal.at(1);
385  E_n.at(2, 2) = normal.at(2);
386  E_n.at(3, 3) = normal.at(3);
387  E_n.at(4, 1) = normal.at(2);
388  E_n.at(5, 1) = normal.at(3);
389  E_n.at(6, 2) = normal.at(3);
390  E_n.at(7, 2) = normal.at(1);
391  E_n.at(8, 3) = normal.at(1);
392  E_n.at(9, 3) = normal.at(2);
393  } else if ( nsd == 2 ) {
394  E_n.at(1, 1) = normal.at(1);
395  E_n.at(2, 2) = normal.at(2);
396  E_n.at(3, 1) = normal.at(2);
397  E_n.at(4, 2) = normal.at(1);
398  } else {
399  E_n.at(1, 1) = normal.at(1);
400  }
401 
402  contrib.beProductOf(E_n, nMatrix);
403 
404  oTangent.add(detJ * gp->giveWeight(), contrib);
405  }
406 }
407 } /* namespace oofem */
bool contains(int value) const
Definition: intarray.h:283
The representation of EngngModel default unknown numbering.
virtual void computeField(FloatArray &sigma, TimeStep *tStep)
Computes the homogenized, macroscopic, field (stress).
void integrateTangent(FloatMatrix &oTangent, Element *e, int iBndIndex)
Help function that integrates the tangent contribution from a single element boundary.
REGISTER_BoundaryCondition(BoundaryCondition)
Class and object Domain.
Definition: domain.h:115
virtual int assemble(const IntArray &loc, const FloatMatrix &mat)=0
Assembles sparse matrix from contribution of local elements.
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Returns local vector of unknowns.
Definition: element.C:86
virtual void computeTangent(FloatMatrix &tangent, TimeStep *tStep)
Computes the macroscopic tangent for homogenization problems through sensitivity analysis.
virtual void giveInputRecord(DynamicInputRecord &input)
virtual DofManager * giveInternalDofManager(int i)
Gives an internal dof manager from receiver.
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition: femcmpnn.h:82
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
SparseMtrx * createSparseMtrx(SparseMtrxType type)
Creates new instance of sparse matrix corresponding to given keyword.
Definition: classfactory.C:105
const IntArray & giveBoundaryList()
Returns list of element boundaries within set.
Definition: set.C:140
Class for homogenization of applied gradients.
Provides Xfem interface for an element.
virtual double boundaryEvalNormal(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the normal on the requested boundary.
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 representing a general abstraction for surface finite element interpolation class...
Definition: feinterpol2d.h:45
int giveInterpolationOrder()
Returns the interpolation order.
Definition: feinterpol.h:159
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
Implementation for assembling tangent matrices in standard monolithic FE-problems.
DiscontinuousSegmentIntegrationRule provides integration over a discontinuous boundary segment...
double giveTargetTime()
Returns target time.
Definition: timestep.h:146
Abstract base class for all finite elements.
Definition: element.h:145
void giveStressLocationArray(IntArray &oCols, const UnknownNumberingScheme &r_s)
Base class for dof managers.
Definition: dofmanager.h:113
virtual MaterialMode giveMaterialMode()
Returns material mode for receiver integration points.
Definition: element.h:691
void negated()
Changes sign of receiver values.
Definition: floatmatrix.C:1605
int giveNumber()
Returns domain number.
Definition: domain.h:266
int giveNumberOfSpatialDimensions()
Returns number of spatial dimensions.
Definition: domain.C:1067
MaterialMode
Type representing material mode of integration point.
Definition: materialmode.h:89
#define OOFEM_LOG_DEBUG(...)
Definition: logger.h:128
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual void computeLocalEdgeMapping(IntArray &edgeNodes, int iedge)=0
virtual FEInterpolation * giveInterpolation() const
Definition: element.h:629
virtual void giveInputRecord(DynamicInputRecord &input)
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
virtual void assemble(SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain)
Assembles characteristic matrix of required type into given sparse matrix.
Definition: engngm.C:803
Class representing "master" degree of freedom.
Definition: masterdof.h:92
void beScaled(double s, const FloatArray &b)
Sets receiver to be .
Definition: floatarray.C:146
int giveSetNumber()
Gives the set number which boundary condition is applied to.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
Element * giveElement(int n)
Service for accessing particular domain fe element.
Definition: domain.C:160
#define OOFEM_ERROR(...)
Definition: error.h:61
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
SparseMtrxType
Enumerative type used to identify the sparse matrix type.
virtual void boundaryLocal2Global(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Maps the local boundary coordinates to global.
void giveLocationArray(IntArray &locationArray, const UnknownNumberingScheme &s, IntArray *dofIds=NULL) const
Returns the location array (array of code numbers) of receiver for given numbering scheme...
Definition: element.C:390
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...
void XfemElementInterface_createEnrNmatrixAt(FloatMatrix &oAnswer, const FloatArray &iLocCoord, Element &iEl, bool iSetDiscontContribToZero)
Creates enriched N-matrix.
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
Abstract base class for all active boundary conditions.
Definition: activebc.h:63
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
Definition: floatarray.C:708
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 IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
void partitionEdgeSegment(int iBndIndex, std::vector< Line > &oSegments, std::vector< FloatArray > &oIntersectionPoints, const double &iTangDistPadding=0.0)
Partition a boundary segment to account for cracks cutting the boundary.
void beNMatrixOf(const FloatArray &n, int nsd)
Assigns the receiver to be a repeated diagonal matrix.
Definition: floatmatrix.C:505
bool hasXfemManager()
Definition: domain.C:386
int giveNextFreeDofID(int increment=1)
Gives the next free dof ID.
Definition: domain.C:1411
Class representing vector of real numbers.
Definition: floatarray.h:82
SparseLinearSystemNM * createSparseLinSolver(LinSystSolverType st, Domain *d, EngngModel *m)
Creates new instance of SparseLinearSystemNM corresponding to given type.
Definition: classfactory.C:120
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
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
virtual IntegrationRule * giveBoundaryIntegrationRule(int order, int boundary)
Sets up a suitable integration rule for integrating over the requested boundary.
Definition: feinterpol.C:63
void add(const FloatMatrix &a)
Adds matrix to the receiver.
Definition: floatmatrix.C:1023
virtual void scale(double s)
Scales the receiver according to given value.
void followedBy(const IntArray &b, int allocChunk=0)
Appends array b at the end of receiver.
Definition: intarray.C:145
virtual void assemble(SparseMtrx &answer, TimeStep *tStep, CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, double scale=1.0)
Assembles B.C.
void beTProductOf(const FloatMatrix &a, const FloatMatrix &b)
Assigns to the receiver product of .
Definition: floatmatrix.C:367
void assembleSquared(const FloatArray &fe, const IntArray &loc)
Assembles the array fe with each component squared.
Definition: floatarray.C:572
Class representing the a dynamic Input Record.
void beTranspositionOf(const FloatMatrix &src)
Assigns to the receiver the transposition of parameter.
Definition: floatmatrix.C:323
ClassFactory & classFactory
Definition: classfactory.C:59
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Computes the element local coordinates from given global coordinates.
Definition: element.C:1222
#define new
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
void giveGradientVoigt(FloatArray &oGradient) const
Gives back the applied gradient in Voigt form.
std::unique_ptr< Node > mpSigmaHom
DOF-manager containing the unknown homogenized stress.
Domain * giveDomain() const
Definition: femcmpnn.h:100
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
void beUnitMatrix()
Sets receiver to unity matrix.
Definition: floatmatrix.C:1332
Abstract base class representing the "problem" under consideration.
Definition: engngm.h:181
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
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
void beInverseOf(const FloatMatrix &src)
Modifies receiver to become inverse of given parameter.
Definition: floatmatrix.C:835
void negated()
Switches the sign of every coefficient of receiver.
Definition: floatarray.C:739
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.
virtual double evaluateAtTime(double t)
Returns the value of the function at given time.
Definition: function.C:76
virtual void boundaryEvalN(FloatArray &answer, int boundary, const FloatArray &lcoords, const FEICellGeometry &cellgeo)=0
Evaluates the basis functions on the requested boundary.
Class representing solution step.
Definition: timestep.h:80
virtual void assembleVector(FloatArray &answer, TimeStep *tStep, CharType type, ValueModeType mode, const UnknownNumberingScheme &s, FloatArray *eNorm=NULL)
Assembles B.C.

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