OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
primvarmapper.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 "primvarmapper.h"
36 #include "domain.h"
37 #include "dofmanager.h"
38 #include "linsystsolvertype.h"
39 #include "../sm/Elements/structuralelement.h"
40 #include "engngm.h"
41 #include "gausspoint.h"
42 #include "feinterpol.h"
43 #include "spatiallocalizer.h"
44 #include "sparsemtrx.h"
45 #include "sparselinsystemnm.h"
46 #include "classfactory.h"
47 #include "timestep.h"
48 #include "activebc.h"
51 #include "unknownnumberingscheme.h"
52 #include "mathfem.h"
53 
54 #include <fstream>
55 
56 namespace oofem {
58 
60 
62 { }
63 
65 { }
66 
68 {
69  EngngModel *engngMod = iNewDom.giveEngngModel();
71 
72 
73  const int dim = iNewDom.giveNumberOfSpatialDimensions();
74 
75  int numElNew = iNewDom.giveNumberOfElements();
76 
77  // Count dofs
78  int numDofsNew = engngMod->giveNumberOfDomainEquations( 1, num );
79 
80 
81  oU.resize(numDofsNew);
82  oU.zero();
83 
84  FloatArray du(numDofsNew);
85  du.zero();
86 
87  FloatArray res(numDofsNew);
88 
89  std :: unique_ptr< SparseMtrx > K;
90  std :: unique_ptr< SparseLinearSystemNM > solver;
91 
92  solver.reset( classFactory.createSparseLinSolver(ST_Petsc, & iOldDom, engngMod) );
93  if (!solver) {
94  solver.reset( classFactory.createSparseLinSolver(ST_Direct, & iOldDom, engngMod) );
95  }
96  K.reset( classFactory.createSparseMtrx(solver->giveRecommendedMatrix(true)) );
97 
98  K->buildInternalStructure( engngMod, iNewDom.giveNumber(), num );
99 
100  int maxIter = 1;
101 
102  for ( int iter = 0; iter < maxIter; iter++ ) {
103  K->zero();
104  res.zero();
105 
106 
107  // Contribution from elements
108  for ( int elIndex = 1; elIndex <= numElNew; elIndex++ ) {
109  StructuralElement *elNew = dynamic_cast< StructuralElement * >( iNewDom.giveElement(elIndex) );
110  if ( elNew == NULL ) {
111  OOFEM_ERROR("Failed to cast Element new to StructuralElement.");
112  }
113 
115  // Compute residual
116 
117  // Count element dofs
118  int numElNodes = elNew->giveNumberOfDofManagers();
119  int numElDofs = 0;
120  for ( int i = 1; i <= numElNodes; i++ ) {
121  numElDofs += elNew->giveDofManager(i)->giveNumberOfDofs();
122  }
123 
124  FloatArray elRes(numElDofs);
125  elRes.zero();
126 
127  IntArray elDofsGlob;
128  elNew->giveLocationArray( elDofsGlob, num );
129 
130 
131  // Loop over Gauss points
132  for ( int intRuleInd = 0; intRuleInd < elNew->giveNumberOfIntegrationRules(); intRuleInd++ ) {
133 
134  for ( GaussPoint *gp: *elNew->giveIntegrationRule(intRuleInd) ) {
135 
136  // New N-matrix
137  FloatMatrix NNew;
138  elNew->computeNmatrixAt(gp->giveNaturalCoordinates(), NNew);
139 
140 
142  // Global coordinates of GP
143  const FloatArray &localCoord = gp->giveNaturalCoordinates();
144  FloatArray globalCoord;
145  elNew->computeGlobalCoordinates(globalCoord, localCoord);
147 
148 
149  // Localize element and point in the old domain
150  FloatArray localCoordOld(dim), pointCoordOld(dim);
151  StructuralElement *elOld = dynamic_cast< StructuralElement * >( iOldDom.giveSpatialLocalizer()->giveElementClosestToPoint(localCoordOld, pointCoordOld, globalCoord, 0) );
152  if ( elOld == NULL ) {
153  OOFEM_ERROR("Failed to cast Element old to StructuralElement.");
154  }
155 
156 
157  // Compute N-Matrix for the old element
158  FloatMatrix NOld;
159  elOld->computeNmatrixAt(localCoordOld, NOld);
160 
161  // Fetch nodal displacements for the new element
162  FloatArray nodeDispNew( elDofsGlob.giveSize() );
163 
164 
165  int dofsPassed = 1;
166  for ( int elNode: elNew->giveDofManArray() ) {
167 // DofManager *dMan = iNewDom.giveNode(elNode);
168  DofManager *dMan = iNewDom.giveDofManager(elNode);
169 
170  for ( Dof *dof: *dMan ) {
171  if ( elDofsGlob.at(dofsPassed) != 0 ) {
172  nodeDispNew.at(dofsPassed) = oU.at( elDofsGlob.at(dofsPassed) );
173  } else {
174  if ( dof->hasBc(& iTStep) ) {
175  nodeDispNew.at(dofsPassed) = dof->giveBcValue(iMode, & iTStep);
176  }
177  }
178 
179  dofsPassed++;
180  }
181  }
182 
183 
184  FloatArray newDisp;
185  newDisp.beProductOf(NNew, nodeDispNew);
186 
187 
188  // Fetch nodal displacements for the old element
189  FloatArray nodeDispOld;
190  dofsPassed = 1;
191  IntArray elDofsGlobOld;
192  elOld->giveLocationArray( elDofsGlobOld, num );
193 
194 // elOld->computeVectorOf(iMode, &(iTStep), nodeDisp);
195  int numElNodesOld = elOld->giveNumberOfDofManagers();
196  for(int nodeIndOld = 1; nodeIndOld <= numElNodesOld; nodeIndOld++) {
197  DofManager *dManOld = elOld->giveDofManager(nodeIndOld);
198 
199  for ( Dof *dof: *dManOld ) {
200  if ( elDofsGlobOld.at(dofsPassed) != 0 ) {
201  FloatArray dofUnknowns;
202 
203  if(dof->giveEqn() > 0) {
204  dof->giveUnknowns(dofUnknowns, iMode, &iTStep);
205 
206 #ifdef DEBUG
207  if(!dofUnknowns.isFinite()) {
208  OOFEM_ERROR("!dofUnknowns.isFinite()")
209  }
210 
211  if(dofUnknowns.giveSize() < 1) {
212  OOFEM_ERROR("dofUnknowns.giveSize() < 1")
213  }
214 #endif
215  nodeDispOld.push_back(dofUnknowns.at(1));
216  }
217  else {
218  // TODO: Why does this case occur?
219  nodeDispOld.push_back(0.0);
220  }
221  } else {
222  if ( dof->hasBc(& iTStep) ) {
223 // printf("hasBC.\n");
224 #ifdef DEBUG
225  if(!std::isfinite(dof->giveBcValue(iMode, & iTStep))) {
226  OOFEM_ERROR("!std::isfinite(dof->giveBcValue(iMode, & iTStep))")
227  }
228 #endif
229  nodeDispOld.push_back( dof->giveBcValue(iMode, & iTStep) );
230  }
231  else {
232 // printf("Unhandled case in LSPrimaryVariableMapper :: mapPrimaryVariables().\n");
233  nodeDispOld.push_back( 0.0 );
234  }
235  }
236 
237  dofsPassed++;
238  }
239 
240  }
241 
242 
243  FloatArray oldDisp;
244  oldDisp.beProductOf(NOld, nodeDispOld);
245 
246  FloatArray temp, duu;
247 
248 #ifdef DEBUG
249  if(!oldDisp.isFinite()) {
250  OOFEM_ERROR("!oldDisp.isFinite()")
251  }
252 
253  if(!newDisp.isFinite()) {
254  OOFEM_ERROR("!newDisp.isFinite()")
255  }
256 #endif
257 
258  duu.beDifferenceOf(oldDisp, newDisp);
259  temp.beTProductOf(NNew, duu);
260  double dV = elNew->computeVolumeAround(gp);
261  elRes.add(dV, temp);
262  }
263  }
264 
265 
267  // Compute matrix
268  FloatMatrix me;
269 
270  double mass = 0.0;
271  double density = 1.0;
272  elNew->computeConsistentMassMatrix(me, & iTStep, mass, & density);
273 
274 #ifdef DEBUG
275  if(!me.isFinite()) {
276  OOFEM_ERROR("!me.isFinite()")
277  }
278 #endif
279 
280  K->assemble(elDofsGlob, me);
281  res.assemble(elRes, elDofsGlob);
282  }
283 
284 
285  // Contribution from active boundary conditions
286  int numBC = iNewDom.giveNumberOfBoundaryConditions();
287  for(int bcInd = 1; bcInd <= numBC; bcInd++) {
288 
289  PrescribedGradientBCWeak *activeBC = dynamic_cast<PrescribedGradientBCWeak*> ( iNewDom.giveBc(bcInd) );
290 
291  if(activeBC != NULL) {
292  IntArray tractionRows;
293  activeBC->giveTractionLocationArray(tractionRows, num);
294 
295  // TODO: Compute correct mass matrix and add residual contribution.
296 
297  FloatMatrix mNode(1,1);
298  mNode.beUnitMatrix();
299 
300  for(int tracDofInd : tractionRows) {
301  const IntArray tracDofArray = {tracDofInd};
302  K->assemble(tracDofArray, tracDofArray, mNode);
303  }
304 
305  }
306 
307 
308  PrescribedGradientBCNeumann *neumannBC = dynamic_cast<PrescribedGradientBCNeumann*> ( iNewDom.giveBc(bcInd) );
309 
310  if(neumannBC != NULL) {
311  IntArray stressRows;
312  neumannBC->giveStressLocationArray(stressRows, num);
313 
314  FloatMatrix massMtrxBc( stressRows.giveSize(), stressRows.giveSize() );
315  massMtrxBc.beUnitMatrix(); // TODO: Compute correct mass matrix and add residual contribution.
316  K->assemble(stressRows, massMtrxBc);
317  }
318 
319  }
320 
321 #ifdef DEBUG
322  if(!res.isFinite()) {
323  OOFEM_ERROR("!res.isFinite()")
324  }
325 #endif
326 
327 // printf("iter: %d res norm: %e\n", iter, res.computeNorm() );
328  K->assembleBegin();
329  K->assembleEnd();
330 // K->writeToFile("Kmapping.txt");
331 
332  // Solve
333  solver->solve(*K, res, du);
334 
335 #ifdef DEBUG
336  if(!du.isFinite()) {
337  OOFEM_ERROR("!du.isFinite()")
338  }
339 #endif
340 
341  oU.add(du);
342  }
343 }
344 } /* namespace oofem */
The representation of EngngModel default unknown numbering.
Class and object Domain.
Definition: domain.h:115
virtual int giveNumberOfDomainEquations(int di, const UnknownNumberingScheme &num)
Returns number of equations for given domain in active (current time step) time step.
Definition: engngm.C:391
virtual IntegrationRule * giveIntegrationRule(int i)
Definition: element.h:835
SparseMtrx * createSparseMtrx(SparseMtrxType type)
Creates new instance of sparse matrix corresponding to given keyword.
Definition: classfactory.C:105
int giveNumberOfBoundaryConditions() const
Returns number of boundary conditions in domain.
Definition: domain.h:440
Imposes a prescribed gradient weakly on the boundary with a Neumann boundary condition.
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
bool isFinite() const
Returns true if no element is NAN or infinite.
Definition: floatmatrix.C:209
EngngModel * giveEngngModel()
Returns engineering model to which receiver is associated.
Definition: domain.C:433
void giveStressLocationArray(IntArray &oCols, const UnknownNumberingScheme &r_s)
Base class for dof managers.
Definition: dofmanager.h:113
int giveNumber()
Returns domain number.
Definition: domain.h:266
int giveNumberOfElements() const
Returns number of elements in domain.
Definition: domain.h:434
int giveNumberOfSpatialDimensions()
Returns number of spatial dimensions.
Definition: domain.C:1067
virtual void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
Computes interpolation matrix for element unknowns.
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual int giveNumberOfDofManagers() const
Definition: element.h:656
bool isFinite() const
Returns true if no element is NAN or infinite.
Definition: floatarray.C:72
void beDifferenceOf(const FloatArray &a, const FloatArray &b)
Sets receiver to be a - b.
Definition: floatarray.C:341
int giveNumberOfIntegrationRules()
Definition: element.h:830
GeneralBoundaryCondition * giveBc(int n)
Service for accessing particular domain bc.
Definition: domain.C:243
Abstract base class for all "structural" finite elements.
virtual double computeVolumeAround(GaussPoint *gp)
Returns volume related to given integration point.
Definition: element.h:518
const IntArray & giveDofManArray() const
Definition: element.h:592
int giveNumberOfDofs() const
Definition: dofmanager.C:279
Element * giveElement(int n)
Service for accessing particular domain fe element.
Definition: domain.C:160
#define OOFEM_ERROR(...)
Definition: error.h:61
virtual void giveTractionLocationArray(IntArray &rows, const UnknownNumberingScheme &s)
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
SpatialLocalizer * giveSpatialLocalizer()
Returns receiver&#39;s associated spatial localizer.
Definition: domain.C:1184
virtual void mapPrimaryVariables(FloatArray &oU, Domain &iOldDom, Domain &iNewDom, ValueModeType iMode, TimeStep &iTStep)
Definition: primvarmapper.C:67
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix and anArray.
Definition: floatarray.C:676
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Receiver becomes the result of the product of aMatrix^T and anArray.
Definition: floatarray.C:708
Imposes a prescribed gradient weakly on the boundary with an independent traction discretization...
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
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 Element * giveElementClosestToPoint(FloatArray &lcoords, FloatArray &closest, const FloatArray &coords, int region=0)=0
Returns the element closest to a given point.
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
ClassFactory & classFactory
Definition: classfactory.C:59
void push_back(const double &iVal)
Add one element.
Definition: floatarray.h:118
void beUnitMatrix()
Sets receiver to unity matrix.
Definition: floatmatrix.C:1332
Abstract base class representing the "problem" under consideration.
Definition: engngm.h:181
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.
DofManager * giveDofManager(int i) const
Definition: element.C:514
Abstract class Dof represents Degree Of Freedom in finite element mesh.
Definition: dof.h:93
DofManager * giveDofManager(int n)
Service for accessing particular domain dof manager.
Definition: domain.C:314
Class representing integration point in finite element program.
Definition: gausspoint.h:93
Class representing solution step.
Definition: timestep.h:80
virtual int computeGlobalCoordinates(FloatArray &answer, const FloatArray &lcoords)
Computes the global coordinates from given element&#39;s local coordinates.
Definition: element.C:1207
void add(const FloatArray &src)
Adds array src to receiver.
Definition: floatarray.C:156
virtual void computeConsistentMassMatrix(FloatMatrix &answer, TimeStep *tStep, double &mass, const double *ipDensity=NULL)
Computes consistent mass matrix of receiver using numerical integration over element volume...
void resize(int s)
Resizes receiver towards requested size.
Definition: floatarray.C:631

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