OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
mixedgradientpressureweakperiodic.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 "dofiditem.h"
37 #include "dofmanager.h"
38 #include "dof.h"
39 #include "valuemodetype.h"
40 #include "floatarray.h"
41 #include "floatmatrix.h"
42 #include "engngm.h"
43 #include "node.h"
44 #include "element.h"
45 #include "integrationrule.h"
46 #include "gaussintegrationrule.h"
47 #include "gausspoint.h"
48 #include "masterdof.h"
49 #include "classfactory.h" // For sparse matrix creation.
50 #include "sparsemtrxtype.h"
51 #include "mathfem.h"
52 #include "sparsemtrx.h"
53 #include "sparselinsystemnm.h"
54 #include "set.h"
55 #include "dynamicinputrecord.h"
56 #include "feinterpol.h"
57 #include "unknownnumberingscheme.h"
58 
59 namespace oofem {
60 REGISTER_BoundaryCondition(MixedGradientPressureWeakPeriodic);
61 
63  voldman( new Node(-1, d) ),
64  tractionsdman( new Node(-1, d) )
65 {
66  int dofid = this->domain->giveNextFreeDofID();
67  v_id.followedBy(dofid);
68  this->voldman->appendDof( new MasterDof( voldman.get(), ( DofIDItem ) dofid ) );
69 }
70 
71 
73 {
74 }
75 
76 
78 {
79  return 2;
80 }
81 
82 
84 {
85  if ( i == 1 ) {
86  return this->tractionsdman.get();
87  } else {
88  return this->voldman.get();
89  }
90 }
91 
92 
94 {
95  IRResultType result;
96 
97 
99  if ( this->order < 0 ) {
100  OOFEM_ERROR("order must be at least 0");
101  }
102 
103  int nsd = this->domain->giveNumberOfSpatialDimensions();
104  int total = nsd * nsd * ( int ) pow(double ( order + 1 ), nsd - 1);
105  this->tractionsdman->setNumberOfDofs(0);
106  t_id.clear();
107  for ( int i = 1; i <= total; ++i ) {
108  int dofid = this->domain->giveNextFreeDofID();
109  t_id.followedBy(dofid);
110  // Simply use t_i = S_i . n, where S_1 = [1,0,0;0,0,0;0,0,0], S_2 = [0,1,0;0,0,0;0,0,0], etc.
111  // then the linear terms, [x,0,0], [0,x,0] ... and so on
112  this->tractionsdman->appendDof( new MasterDof( tractionsdman.get(), ( DofIDItem ) dofid ) );
113  }
114 
116 }
117 
118 
120 {
121  this->constructFullMatrixForm(this->devGradient, t);
122  this->volGradient = 0.;
123  for ( int i = 1; i <= this->devGradient.giveNumberOfRows(); i++ ) {
124  this->volGradient += this->devGradient.at(i, i);
125  }
126 }
127 
128 
130 {
131  int n = d_voigt.giveSize();
132  if ( n == 6 ) { // Then 3D
133  d.resize(3, 3);
134  d.at(1, 1) = d_voigt.at(1);
135  d.at(2, 2) = d_voigt.at(2);
136  d.at(3, 3) = d_voigt.at(3);
137  d.at(2, 3) = d.at(3, 2) = d_voigt.at(4) * 0.5;
138  d.at(1, 3) = d.at(3, 1) = d_voigt.at(5) * 0.5;
139  d.at(1, 2) = d.at(2, 1) = d_voigt.at(6) * 0.5;
140  } else if ( n == 3 ) { // Then 2D
141  d.resize(2, 2);
142  d.at(1, 1) = d_voigt.at(1);
143  d.at(2, 2) = d_voigt.at(2);
144  d.at(1, 2) = d.at(2, 1) = d_voigt.at(3) * 0.5;
145  } else if ( n == 1 ) { // Then 1D
146  d.resize(1, 1);
147  d.at(1, 1) = d_voigt.at(1);
148  } else {
149  OOFEM_ERROR("Tensor is in strange voigt format.");
150  }
151 }
152 
153 
154 void MixedGradientPressureWeakPeriodic :: giveLocationArrays(std :: vector< IntArray > &rows, std :: vector< IntArray > &cols, CharType type,
155  const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s)
156 {
157  if ( type != TangentStiffnessMatrix ) {
158  return;
159  }
160 
161  IntArray loc_r, loc_c, t_loc_r, t_loc_c, e_loc_r, e_loc_c;
162 
163  // Fetch the columns/rows for the tractions;
164  this->tractionsdman->giveLocationArray(t_id, t_loc_r, r_s);
165  this->tractionsdman->giveLocationArray(t_id, t_loc_c, c_s);
166 
167  // Fetch the columns/rows for dvol;
168  this->voldman->giveLocationArray(v_id, e_loc_r, r_s);
169  this->voldman->giveLocationArray(v_id, e_loc_c, c_s);
170 
171  Set *set = this->giveDomain()->giveSet(this->set);
172  IntArray bNodes;
173  const IntArray &boundaries = set->giveBoundaryList();
174 
175  rows.resize(boundaries.giveSize() + 2);
176  cols.resize(boundaries.giveSize() + 2);
177  int i = 0;
178  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
179  Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
180  int boundary = boundaries.at(pos * 2);
181 
182  e->giveInterpolation()->boundaryGiveNodes(bNodes, boundary);
183  e->giveBoundaryLocationArray(loc_r, bNodes, this->dofs, r_s);
184  e->giveBoundaryLocationArray(loc_c, bNodes, this->dofs, c_s);
185  // For most uses, *loc_r == *loc_c
186  rows [ i ] = loc_r;
187  cols [ i ] = t_loc_c;
188  i++;
189  // and the symmetric part (usually the transpose of above)
190  rows [ i ] = t_loc_r;
191  cols [ i ] = loc_c;
192  i++;
193  }
194 
195  // The volumetric part:
196  rows [ i ] = t_loc_r;
197  cols [ i ] = e_loc_c;
198  i++;
199 
200  rows [ i ] = e_loc_r;
201  cols [ i ] = t_loc_c;
202  i++;
203 }
204 
205 
207 {
208  // This evaluates the function x^a * y^b, for a, b in [0,order]
209  int total = ( int ) pow(double ( order + 1 ), surfCoords.giveSize());
210  answer.resize(total);
211  int pos = 0;
212  double tx = 1.;
213  for ( int xi = 0; xi < ( this->order + 1 ); ++xi ) {
214  double ty = 1.;
215  for ( int yi = 0; yi < ( this->order + 1 ); ++yi ) {
216  answer(pos) = tx * ty;
217  pos++;
218  ty *= surfCoords(1);
219  }
220  tx *= surfCoords(0);
221  }
222 }
223 
224 
226 {
227  FloatArray t, surfCoords;
228  int nsd = this->giveDomain()->giveNumberOfSpatialDimensions();
229  int total = nsd * nsd * ( int ) pow(double ( this->order + 1 ), nsd - 1);
230 
231  surfCoords.resize(nsd - 1);
232  mMatrix.resize(nsd, total);
233  mMatrix.zero();
234  int pos = 0;
235  for ( int kj = 0; kj < nsd; ++kj ) {
236  // First we compute the surface coordinates. Its the global coordinates without the kj component.
237  for ( int c = 0, q = 0; c < nsd; ++c ) {
238  if ( c != kj ) {
239  surfCoords(q) = coords(c);
240  q++;
241  }
242  }
243  // Evaluate all the basis functions for a given column.
244  this->evaluateTractionBasisFunctions(t, surfCoords);
245  for ( int ti = 0; ti < t.giveSize(); ++ti ) {
246  for ( int ki = 0; ki < nsd; ++ki ) {
247  mMatrix(ki, pos) = t(ti) * normal(kj);
248  pos++;
249  }
250  }
251  }
252 }
253 
254 
256 {
257  // Computes the integral: int dt . dv dA
258  FloatArray normal, n, coords;
259  FloatMatrix nMatrix, mMatrix;
260 
261  FEInterpolation *interp = el->giveInterpolation(); // Geometry interpolation. The displacements or velocities must have the same interpolation scheme (on the boundary at least).
262 
263  int maxorder = this->order + interp->giveInterpolationOrder() * 3;
264  std :: unique_ptr< IntegrationRule >ir( interp->giveBoundaryIntegrationRule(maxorder, boundary) );
265  int nsd = this->giveDomain()->giveNumberOfSpatialDimensions();
266 
267  answer.clear();
268  for ( GaussPoint *gp: *ir ) {
269  const FloatArray &lcoords = gp->giveNaturalCoordinates();
270  FEIElementGeometryWrapper cellgeo(el);
271 
272  double detJ = interp->boundaryEvalNormal(normal, boundary, lcoords, cellgeo);
273  interp->boundaryEvalN(n, boundary, lcoords, cellgeo);
274  interp->boundaryLocal2Global(coords, boundary, lcoords, cellgeo);
275 
276  // Construct the basis functions for the tractions:
277  this->constructMMatrix(mMatrix, coords, normal);
278  nMatrix.beNMatrixOf(n, nsd);
279 
280  answer.plusProductUnsym( mMatrix, nMatrix, detJ * gp->giveWeight() );
281  }
282 }
283 
284 
286 {
287  // Computes the integral: int dt . dx_m dA
288  FloatMatrix mMatrix;
289  FloatArray normal, coords, vM_vol;
290 
291  FEInterpolation *interp = el->giveInterpolation(); // Geometry interpolation. The displacements or velocities must have the same interpolation scheme (on the boundary at least).
292 
293  int maxorder = this->order + interp->giveInterpolationOrder() * 3;
294  std :: unique_ptr< IntegrationRule >ir( interp->giveBoundaryIntegrationRule(maxorder, boundary) );
295 
296  FloatArray tmpAnswer;
297  for ( GaussPoint *gp: *ir ) {
298  const FloatArray &lcoords = gp->giveNaturalCoordinates();
299  FEIElementGeometryWrapper cellgeo(el);
300 
301  double detJ = interp->boundaryEvalNormal(normal, boundary, lcoords, cellgeo);
302  interp->boundaryLocal2Global(coords, boundary, lcoords, cellgeo);
303 
304  vM_vol.beScaled(1.0/3.0, coords);
305  this->constructMMatrix(mMatrix, coords, normal);
306 
307  tmpAnswer.plusProduct(mMatrix, vM_vol, detJ * gp->giveWeight());
308  }
309  answer.resize(tmpAnswer.giveSize(), 1);
310  answer.setColumn(tmpAnswer, 1);
311 }
312 
313 
315 {
316  // Computes the integral: int dt . dx dA
317  FloatMatrix mMatrix;
318  FloatArray normal, coords, vM_dev;
319 
320  FEInterpolation *interp = el->giveInterpolation(); // Geometry interpolation. The displacements or velocities must have the same interpolation scheme (on the boundary at least).
321 
322  int maxorder = this->order + interp->giveInterpolationOrder() * 3;
323  std :: unique_ptr< IntegrationRule >ir( interp->giveBoundaryIntegrationRule(maxorder, boundary) );
324  answer.clear();
325 
326  for ( GaussPoint *gp: *ir ) {
327  const FloatArray &lcoords = gp->giveNaturalCoordinates();
328  FEIElementGeometryWrapper cellgeo(el);
329 
330  double detJ = interp->boundaryEvalNormal(normal, boundary, lcoords, cellgeo);
331  // Compute v_m = d_dev . x
332  interp->boundaryLocal2Global(coords, boundary, lcoords, cellgeo);
333  vM_dev.beProductOf(ddev, coords);
334 
335  this->constructMMatrix(mMatrix, coords, normal);
336 
337  answer.plusProduct(mMatrix, vM_dev, detJ * gp->giveWeight());
338  }
339 }
340 
341 
343  CharType type, ValueModeType mode, const UnknownNumberingScheme &s, FloatArray *eNorms)
344 {
345  Set *set = this->giveDomain()->giveSet(this->set);
346  const IntArray &boundaries = set->giveBoundaryList();
347 
348  IntArray v_loc, t_loc, e_loc; // For the velocities and stress respectively
349  IntArray velocityDofIDs, bNodes;
350  this->tractionsdman->giveLocationArray(t_id, t_loc, s);
351  this->voldman->giveLocationArray(v_id, e_loc, s);
352 
353  if ( type == ExternalForcesVector ) {
354  // The external forces have two contributions. on the traction and on dvol.
355  double rve_size = this->domainSize();
356 
357  if ( e_loc.at(1) ) {
358  answer.at( e_loc.at(1) ) -= rve_size * pressure; // Note the negative sign (pressure as opposed to mean stress)
359  if ( eNorms ) {
360  eNorms->at( v_id.at(1) ) += rve_size * pressure * rve_size * pressure;
361  }
362  }
363 
364  // The second contribution is on the momentumbalance equation; int t . [[ d_dev . x ]] dA = int t . [[ d_dev . x ]] dA
365  FloatArray fe;
366  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
367  Element *e = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
368  int boundary = boundaries.at(pos * 2);
369 
370  this->integrateTractionDev(fe, e, boundary, this->devGradient);
371  fe.negated();
372 
373  answer.assemble(fe, t_loc);
374  if ( eNorms ) {
375  eNorms->assembleSquared(fe, velocityDofIDs);
376  }
377  }
378  } else if ( type == InternalForcesVector ) {
379  FloatMatrix Ke_v, Ke_e;
380  FloatArray fe_v, fe_t, fe_t2, fe_e(1);
381  FloatArray t, e, v;
382 
383  // Fetch the current values of internal dofs and their master dof ids;
384  this->tractionsdman->giveUnknownVector(t, t_id, mode, tStep);
385  this->voldman->giveUnknownVector(e, v_id, mode, tStep);
386 
387  // Assemble: -int t . [[ delta_v ]] dA
388  // int delta_t . [[ e.x - v ]] dA
389  // int t . [[ x ]] dA delta_e
390  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
391  Element *el = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
392  int boundary = boundaries.at(pos * 2);
393 
394  // Fetch the element information;
395  el->giveInterpolation()->boundaryGiveNodes(bNodes, boundary);
396  el->giveBoundaryLocationArray(v_loc, bNodes, this->dofs, s, & velocityDofIDs);
397  el->computeBoundaryVectorOf(bNodes, this->dofs, mode, tStep, v);
398 
399  // Integrate the tangents;
400  this->integrateTractionVelocityTangent(Ke_v, el, boundary);
401  this->integrateTractionXTangent(Ke_e, el, boundary);
402 
403  // We just use the tangent, less duplicated code
404  fe_v.beTProductOf(Ke_v, t);
405  fe_v.negated();
406  fe_t.beProductOf(Ke_v, v);
407  fe_t.negated();
408  fe_t2.beProductOf(Ke_e, e);
409  fe_t.add(fe_t2);
410  fe_e.beTProductOf(Ke_e, t);
411 
412  answer.assemble(fe_v, v_loc); // Contributions to delta_v equations
413  answer.assemble(fe_t, t_loc); // Contribution to delta_t equations
414  answer.assemble(fe_e, e_loc); // Contribution to delta_e equations
415  if ( eNorms ) {
416  eNorms->assembleSquared(fe_v, velocityDofIDs);
417  eNorms->assembleSquared(fe_t, t_id);
418  eNorms->assembleSquared(fe_e, v_id);
419  }
420  }
421  }
422 }
423 
424 
426  CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, double scale)
427 {
428  if ( type == TangentStiffnessMatrix || type == SecantStiffnessMatrix || type == ElasticStiffnessMatrix ) {
429  FloatMatrix Ke_v, Ke_vT, Ke_e, Ke_eT;
430  IntArray v_loc_r, v_loc_c, t_loc_r, t_loc_c, e_loc_r, e_loc_c;
431  IntArray bNodes;
432  Set *set = this->giveDomain()->giveSet(this->set);
433  const IntArray &boundaries = set->giveBoundaryList();
434 
435  // Fetch the columns/rows for the stress contributions;
436  this->tractionsdman->giveLocationArray(t_id, t_loc_r, r_s);
437  this->tractionsdman->giveLocationArray(t_id, t_loc_c, c_s);
438 
439  this->voldman->giveLocationArray(v_id, e_loc_r, r_s);
440  this->voldman->giveLocationArray(v_id, e_loc_c, c_s);
441 
442  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
443  Element *el = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
444  int boundary = boundaries.at(pos * 2);
445 
446  // Fetch the element information;
447  el->giveInterpolation()->boundaryGiveNodes(bNodes, boundary);
448  el->giveBoundaryLocationArray(v_loc_r, bNodes, this->dofs, r_s);
449  el->giveBoundaryLocationArray(v_loc_c, bNodes, this->dofs, c_s);
450 
451  this->integrateTractionVelocityTangent(Ke_v, el, boundary);
452  this->integrateTractionXTangent(Ke_e, el, boundary);
453  Ke_v.negated();
454  Ke_v.times(scale);
455  Ke_e.times(scale);
456  Ke_vT.beTranspositionOf(Ke_v);
457  Ke_eT.beTranspositionOf(Ke_e);
458 
459  answer.assemble(t_loc_r, v_loc_c, Ke_v);
460  answer.assemble(v_loc_r, t_loc_c, Ke_vT);
461 
462  answer.assemble(t_loc_r, e_loc_c, Ke_e);
463  answer.assemble(e_loc_r, t_loc_c, Ke_eT);
464  }
465  }
466 }
467 
468 
470 {
471  double rve_size = this->domainSize();
472  FloatArray tractions;
473  this->tractionsdman->giveUnknownVector(tractions, t_id, VM_Total, tStep);
474  this->computeStress(sigmaDev, tractions, rve_size);
475 
476  // And the volumetric part is much easier.
477  vol = (*this->voldman->begin())->giveUnknown(VM_Total, tStep);
478  vol /= rve_size;
479  vol -= volGradient; // This is needed for consistency; We return the volumetric "residual" if a gradient with volumetric contribution is set.
480 }
481 
482 
484 {
485  FloatMatrix mMatrix;
486  FloatArray normal, coords, t;
487 
489  Set *set = this->giveDomain()->giveSet(this->set);
490  const IntArray &boundaries = set->giveBoundaryList();
491  // Reminder: sigma = int t * n dA, where t = sum( C_i * n t_i ).
492  // This loop will construct sigma in matrix form.
493 
494  FloatMatrix sigma;
495 
496  for ( int pos = 1; pos <= boundaries.giveSize() / 2; ++pos ) {
497  Element *el = this->giveDomain()->giveElement( boundaries.at(pos * 2 - 1) );
498  int boundary = boundaries.at(pos * 2);
499 
500  FEInterpolation *interp = el->giveInterpolation(); // Geometry interpolation. The displacements or velocities must have the same interpolation scheme (on the boundary at least).
501 
502  int maxorder = this->order + interp->giveInterpolationOrder() * 3;
503  std :: unique_ptr< IntegrationRule >ir( interp->giveBoundaryIntegrationRule(maxorder, boundary) );
504 
505  for ( GaussPoint *gp: *ir ) {
506  const FloatArray &lcoords = gp->giveNaturalCoordinates();
507  FEIElementGeometryWrapper cellgeo(el);
508 
509  double detJ = interp->boundaryEvalNormal(normal, boundary, lcoords, cellgeo);
510  // Compute v_m = d_dev . x
511  interp->boundaryLocal2Global(coords, boundary, lcoords, cellgeo);
512 
513  this->constructMMatrix(mMatrix, coords, normal);
514  t.beProductOf(mMatrix, tractions);
515  sigma.plusDyadUnsym(t, coords, detJ * gp->giveWeight());
516  }
517  }
518  sigma.times(1. / rve_size);
519 
520  double pressure = 0.;
521  for ( int i = 1; i <= nsd; i++ ) {
522  pressure += sigma.at(i, i);
523  }
524  pressure /= 3; // Not 100% sure about this for 2D cases.
525  if ( nsd == 3 ) {
526  sigmaDev.resize(6);
527  sigmaDev.at(1) = sigma.at(1, 1) - pressure;
528  sigmaDev.at(2) = sigma.at(2, 2) - pressure;
529  sigmaDev.at(3) = sigma.at(3, 3) - pressure;
530  sigmaDev.at(4) = 0.5 * ( sigma.at(2, 3) + sigma.at(3, 2) );
531  sigmaDev.at(5) = 0.5 * ( sigma.at(1, 3) + sigma.at(3, 1) );
532  sigmaDev.at(6) = 0.5 * ( sigma.at(1, 2) + sigma.at(2, 1) );
533  } else if ( nsd == 2 ) {
534  sigmaDev.resize(3);
535  sigmaDev.at(1) = sigma.at(1, 1) - pressure;
536  sigmaDev.at(2) = sigma.at(2, 2) - pressure;
537  sigmaDev.at(3) = 0.5 * ( sigma.at(1, 2) + sigma.at(2, 1) );
538  } else {
539  sigmaDev.resize(1);
540  sigmaDev.at(1) = sigma.at(1, 1) - pressure;
541  }
542 }
543 
545 {
546  //double size = this->domainSize();
547  // Fetch some information from the engineering model
548  EngngModel *rve = this->giveDomain()->giveEngngModel();
550  std :: unique_ptr< SparseLinearSystemNM > solver(
551  classFactory.createSparseLinSolver( ST_Petsc, this->domain, this->domain->giveEngngModel() ) ); // = rve->giveLinearSolver();
552  SparseMtrxType stype = solver->giveRecommendedMatrix(true);
554  Set *set = this->giveDomain()->giveSet(this->set);
555  const IntArray &boundaries = set->giveBoundaryList();
556  double rve_size = this->domainSize();
557 
558  // Set up and assemble tangent FE-matrix which will make up the sensitivity analysis for the macroscopic material tangent.
559  std :: unique_ptr< SparseMtrx > Kff( classFactory.createSparseMtrx(stype) );
560  if ( !Kff ) {
561  OOFEM_ERROR("Couldn't create sparse matrix of type %d\n", stype);
562  }
563  Kff->buildInternalStructure(rve, this->domain->giveNumber(), fnum);
564  rve->assemble(*Kff, tStep, TangentAssembler(TangentStiffness), fnum, this->domain);
565 
566  // Setup up indices and locations
567  int neq = Kff->giveNumberOfRows();
568 
569  // Indices and such of internal dofs
570  int nsd = this->devGradient.giveNumberOfColumns();
571  int nd = nsd * ( 1 + nsd ) / 2;
572 
573  IntArray t_loc, e_loc;
574 
575  // Fetch the positions;
576  this->tractionsdman->giveLocationArray( t_id, t_loc, fnum );
577  this->voldman->giveLocationArray( v_id, e_loc, fnum );
578 
579  // Matrices and arrays for sensitivities
580  FloatMatrix ddev_pert(neq, nd);
581  FloatArray p_pert(neq);
582 
583  // Solutions for the pertubations:
584  FloatMatrix s_d(neq, nd);
585  FloatArray s_p(neq);
586 
587  // Unit pertubations for d_dev
588  for ( int i = 1; i <= nd; ++i ) {
589  FloatMatrix ddev_tmp;
590  FloatArray tmp(neq), fe, fetot, ddevVoigt(nd);
591  ddevVoigt.at(i) = 1.0;
592  this->constructFullMatrixForm(ddev_tmp, ddevVoigt);
593  for ( int k = 1; k <= boundaries.giveSize() / 2; ++k ) {
594  Element *e = this->giveDomain()->giveElement( boundaries.at(k * 2 - 1) );
595  int boundary = boundaries.at(k * 2);
596 
597  this->integrateTractionDev(fe, e, boundary, ddev_tmp);
598  fetot.subtract(fe);
599  }
600  tmp.assemble(fetot, t_loc);
601  ddev_pert.setColumn(tmp, i);
602  }
603 
604  // Unit pertubation for p
605  p_pert.zero();
606  p_pert.at( e_loc.at(1) ) = - 1.0 * rve_size;
607 
608  // Solve all sensitivities
609  solver->solve(*Kff, ddev_pert, s_d);
610  solver->solve(*Kff, p_pert, s_p);
611 
612  // Extract the tractions from the sensitivity solutions s_d and s_p:
613  FloatArray tractions_p( t_loc.giveSize() );
614  FloatMatrix tractions_d(t_loc.giveSize(), nd);
615  tractions_p.beSubArrayOf(s_p, t_loc);
616  for ( int i = 1; i <= nd; ++i ) {
617  for ( int k = 1; k <= t_loc.giveSize(); ++k ) {
618  tractions_d.at(k, i) = s_d.at(t_loc.at(k), i);
619  }
620  }
621 
622  // The de/dp tangent:
623  Cp = - s_p.at( e_loc.at(1) );
624  // The de/dd tangent:
625  Cd.resize(nd);
626  if ( nsd == 3 ) {
627  Cd.at(1) = s_d.at(e_loc.at(1), 1);
628  Cd.at(2) = s_d.at(e_loc.at(1), 2);
629  Cd.at(3) = s_d.at(e_loc.at(1), 3);
630  Cd.at(4) = s_d.at(e_loc.at(1), 4);
631  Cd.at(5) = s_d.at(e_loc.at(1), 5);
632  Cd.at(6) = s_d.at(e_loc.at(1), 6);
633  } else if ( nsd == 2 ) {
634  Cd.at(1) = s_d.at(e_loc.at(1), 1);
635  Cd.at(2) = s_d.at(e_loc.at(1), 2);
636  Cd.at(3) = s_d.at(e_loc.at(1), 3);
637  }
638  // The ds/de tangent:
639  Ep = Cd; // This is much simpler, but it's "only" correct if there is a potential (i.e. symmetric problems). This could be generalized if needed.
640  // The ds/dd tangent:
641  Ed.resize(nd, nd);
642  for ( int dpos = 1; dpos <= nd; ++dpos ) {
643  FloatArray tractions, sigmaDev;
644  tractions.beColumnOf(tractions_d, dpos);
645  this->computeStress(sigmaDev, tractions, rve_size);
646  Ed.setColumn(sigmaDev, dpos);
647  }
648 }
649 
651 {
654  OOFEM_ERROR("Not supported yet");
655  //FloatArray devGradientVoigt;
656  //input.setField(devGradientVoigt, _IFT_MixedGradientPressure_devGradient);
657 }
658 
660 {
661  devGradient.times(s);
662  pressure *= s;
663 }
664 } // end namespace oofem
void setField(int item, InputFieldType id)
int giveNumberOfColumns() const
Returns number of columns of receiver.
Definition: floatmatrix.h:158
The representation of EngngModel default unknown numbering.
void subtract(const FloatArray &src)
Subtracts array src to receiver.
Definition: floatarray.C:258
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.
#define _IFT_MixedGradientPressureWeakPeriodic_order
Order of global polynomial in the unknown tractions. Should be at least 1.
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
SparseMtrx * createSparseMtrx(SparseMtrxType type)
Creates new instance of sparse matrix corresponding to given keyword.
Definition: classfactory.C:105
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
double domainSize()
Computes the size (including pores) by surface integral over the domain.
void clear()
Clears receiver (zero size).
Definition: floatarray.h:206
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.
std::unique_ptr< Node > tractionsdman
DOF-manager containing the unknown tractions (Lagrange mult. for micro-periodic velocity) ...
void constructMMatrix(FloatMatrix &mMatrix, FloatArray &coords, FloatArray &normal)
Abstract base class for all finite elements.
Definition: element.h:145
Base class for dof managers.
Definition: dofmanager.h:113
MixedGradientPressureWeakPeriodic(int n, Domain *d)
Creates boundary condition with given number, belonging to given domain.
void negated()
Changes sign of receiver values.
Definition: floatmatrix.C:1605
int giveNumber()
Returns domain number.
Definition: domain.h:266
virtual double giveUnknown(PrimaryField &field, ValueModeType mode, TimeStep *tStep, ActiveDof *dof)
Computes the value of the dof.
Definition: activebc.h:200
void plusDyadUnsym(const FloatArray &a, const FloatArray &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:811
int giveNumberOfSpatialDimensions()
Returns number of spatial dimensions.
Definition: domain.C:1067
virtual void assemble(SparseMtrx &answer, TimeStep *tStep, CharType type, const UnknownNumberingScheme &r_s, const UnknownNumberingScheme &c_s, double scale=1.0)
Assembles B.C.
Class implementing an array of integers.
Definition: intarray.h:61
int & at(int i)
Coefficient access function.
Definition: intarray.h:103
virtual FEInterpolation * giveInterpolation() const
Definition: element.h:629
void beColumnOf(const FloatMatrix &mat, int col)
Reciever will be set to a given column in a matrix.
Definition: floatarray.C:1114
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
Class representing a general abstraction for finite element interpolation class.
Definition: feinterpol.h:132
void plusProductUnsym(const FloatMatrix &a, const FloatMatrix &b, double dV)
Adds to the receiver the product .
Definition: floatmatrix.C:780
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
void evaluateTractionBasisFunctions(FloatArray &answer, const FloatArray &coords)
Element * giveElement(int n)
Service for accessing particular domain fe element.
Definition: domain.C:160
double volGradient
The volumetric part of what was sent in (needed to return the difference).
void integrateTractionVelocityTangent(FloatMatrix &answer, Element *el, int boundary)
#define OOFEM_ERROR(...)
Definition: error.h:61
void computeStress(FloatArray &sigmaDev, FloatArray &tractions, double rve_size)
void integrateTractionDev(FloatArray &answer, Element *el, int boundary, const FloatMatrix &ddev)
void clear()
Clears the array (zero size).
Definition: intarray.h:177
DofIDItem
Type representing particular dof type.
Definition: dofiditem.h:86
SparseMtrxType
Enumerative type used to identify the sparse matrix type.
#define _IFT_MixedGradientPressure_pressure
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
IntArray dofs
Dofs that b.c. is applied to (relevant for Dirichlet type b.c.s).
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.
std::unique_ptr< Node > voldman
DOF-manager containing the unknown volumetric gradient (always exactly one dof).
void beNMatrixOf(const FloatArray &n, int nsd)
Assigns the receiver to be a repeated diagonal matrix.
Definition: floatmatrix.C:505
virtual void boundaryGiveNodes(IntArray &answer, int boundary)=0
Gives the boundary nodes for requested boundary number.
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 integrateTractionXTangent(FloatMatrix &answer, Element *el, int boundary)
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 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.
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
virtual void assembleVector(FloatArray &answer, TimeStep *tStep, CharType type, ValueModeType mode, const UnknownNumberingScheme &s, FloatArray *eNorm=NULL)
Assembles B.C.
void followedBy(const IntArray &b, int allocChunk=0)
Appends array b at the end of receiver.
Definition: intarray.C:145
void zero()
Zeroes all coefficients of receiver.
Definition: floatarray.C:658
void assembleSquared(const FloatArray &fe, const IntArray &loc)
Assembles the array fe with each component squared.
Definition: floatarray.C:572
virtual void scale(double s)
Scales the receiver according to given value.
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
void setColumn(const FloatArray &src, int c)
Sets the values of the matrix in specified column.
Definition: floatmatrix.C:648
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
#define new
virtual void computeTangents(FloatMatrix &Ed, FloatArray &Ep, FloatArray &Cd, double &Cp, TimeStep *tStep)
Computes the macroscopic tangents through sensitivity analysis.
void zero()
Zeroes all coefficient of receiver.
Definition: floatmatrix.C:1326
Domain * giveDomain() const
Definition: femcmpnn.h:100
virtual void giveInputRecord(DynamicInputRecord &input)
Setups the input record string of receiver.
virtual IRResultType initializeFrom(InputRecord *ir)
Initializes receiver according to object description stored in input record.
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.
virtual DofManager * giveInternalDofManager(int i)
Returns the volumetric DOF manager for i == 1, and the deviatoric manager for i == 2...
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
virtual int giveNumberOfInternalDofManagers()
Returns the number of internal DOF managers (=2).
void negated()
Switches the sign of every coefficient of receiver.
Definition: floatarray.C:739
virtual void computeFields(FloatArray &sigmaDev, double &vol, TimeStep *tStep)
Computes the homogenized fields through sensitivity analysis.
int giveNumberOfRows() const
Returns number of rows of receiver.
Definition: floatmatrix.h:156
General class for boundary condition that prolongates macroscopic fields to incompressible flow...
void computeBoundaryVectorOf(const IntArray &bNodes, const IntArray &dofIDMask, ValueModeType u, TimeStep *tStep, FloatArray &answer, bool padding=false)
Boundary version of computeVectorOf.
Definition: element.C:139
void constructFullMatrixForm(FloatMatrix &d, const FloatArray &d_voigt) const
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
virtual void setPrescribedDeviatoricGradientFromVoigt(const FloatArray &ddev)
Sets the prescribed tensor from the matrix from given Voigt notation.
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:30 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011