OOFEM 3.0
Loading...
Searching...
No Matches
micromaterial.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 - 2025 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 "micromaterial.h"
38#include "domain.h"
39#include "dofmanager.h"
40#include "classfactory.h"
41#include "oofemtxtdatareader.h"
42#include "util.h"
43#include "classfactory.h"
44#include "node.h"
45#include "engngm.h"
46
47namespace oofem {
48//valgrind --leak-check=full --show-reachable=no -v --log-file=valgr.txt ./oofem -f Macrolspace_1.in
49
50// FloatArray A;
51// FloatMatrix K, B;
52// stiffnessMatrix->toFloatMatrix(K);
53// K.symmetrized();
54// K.printYourself();
55// displacementVector.printYourself();
56// A.beProductOf(K, displacementVector);
57// A.printYourself();
58//
59// B.beInverseOf(K);
60// A.beProductOf(B, displacementVector);
61// A.printYourself();
62
64
65// constructor
66//strainVector, tempStrainVector, stressVector, tempStressVector are defined on StructuralMaterialStatus
67MicroMaterialStatus :: MicroMaterialStatus(GaussPoint *gp) : StructuralMaterialStatus(gp) { }
68
69
70void MicroMaterialStatus :: initTempStatus()
71{
72 StructuralMaterialStatus :: initTempStatus();
73}
74
75void MicroMaterialStatus :: updateYourself(TimeStep *tStep)
76{
77 StructuralMaterialStatus :: updateYourself(tStep);
78}
79
80void MicroMaterialStatus :: printOutputAt(FILE *file, TimeStep *tStep) const
81{ }
82
83void MicroMaterialStatus :: saveContext(DataStream &stream, ContextMode mode)
84{
85}
86
87void MicroMaterialStatus :: restoreContext(DataStream &stream, ContextMode mode)
88{
89}
90
91
92MicroMaterial :: MicroMaterial(int n, Domain *d) : StructuralMaterial(n, d), UnknownNumberingScheme()
93{}
94
95
96void MicroMaterial :: initializeFrom(InputRecord &ir)
97{
99
100 OOFEM_LOG_INFO( "** Instanciating microproblem with BC from file %s\n", inputFileNameMicro.c_str() );
101 OOFEMTXTDataReader drMicro( inputFileNameMicro.c_str() );
102 this->problemMicro = InstanciateProblem(drMicro, _processor, 0); //0=contextFlag-store/resore
103 drMicro.finish();
104 OOFEM_LOG_INFO("Microproblem instanciated\n");
105}
106
107
108//original pure virtual function has to be declared here
109//this function should not be used, internal forces are calculated based on reactions not stresses in GPs
110FloatArrayF<6> MicroMaterial :: giveRealStressVector_3d(const FloatArrayF<6> &totalStrain, GaussPoint *gp, TimeStep *tStep) const
111{
112 //perform average over microproblem
113 // int index;
114 // double dV, VolTot = 0.;
115 // double scale = 1.;
116 // FloatArray VecStrain, VecStress, SumStrain(6), SumStress(6);
117 // GaussPoint *gpL;
118 // Domain *microDomain = problemMicro->giveDomain(1); //from engngm.h
119 // EngngModel *microEngngModel = microDomain->giveEngngModel();
120 // StructuralMaterialStatus *status = static_cast< StructuralMaterialStatus * >( this->giveStatus(gp) );
121
122 OOFEM_ERROR("Should not be called, use giveInternalForcesVector instead");
123
124 // int nelem = microDomain->giveNumberOfElements();
125 // //int nnodes = microDomain->giveNumberOfDofManagers();
126 //
127 // //need to update stress so change boundary conditions and recalculate
128 // OOFEM_LOG_INFO( "\n*** giveRealStress microproblem %p on macroElement %d GP %d, step %d, time %f\n", this, this->macroLSpaceElement->giveNumber(), gp->giveNumber(), microEngngModel->giveCurrentStep()->giveNumber(), microEngngModel->giveCurrentStep()->giveTime() );
129 //
130 // this->macroLSpaceElement->changeMicroBoundaryConditions(tStep);
131 // microEngngModel->solveYourselfAt( microEngngModel->giveCurrentStep() );
132 // microEngngModel->terminate( microEngngModel->giveCurrentStep() );
133 // OOFEM_LOG_INFO("\n*** giveRealStress microproblem %p done\n", this);
134 //
135 // for ( int ielem = 1; ielem <= nelem; ielem++ ) { //return stress as average through all elements of the same MicroMaterial
136 // Element *elem = microDomain->giveElement(ielem);
137 // for ( GaussPoint *gpL: *elem->giveDefaultIntegrationRulePtr() ) {
138 // dV = elem->computeVolumeAround(gpL);
139 // VolTot += dV;
140 // //OOFEM_LOG_INFO("Element %d GP %d Vol %f\n", elem->giveNumber(), gp->giveNumber(), dV);
141 // //fprintf(this->stream, "Element %d GP %d stress %f\n", elem->giveNumber(), gp->giveNumber(), 0.0);
142 // //((StructuralCrossSection*) gp->giveCrossSection())->giveFullCharacteristicVector(helpVec, gp, strainVector);
143 // elem->giveIPValue(VecStrain, gpL, IST_StrainTensor, tStep);
144 // elem->giveIPValue(VecStress, gpL, IST_StressTensor, tStep);
145 // elem->giveIntVarCompFullIndx(Mask, IST_StrainTensor);
146 //
147 // VecStrain.times(dV);
148 // VecStress.times(dV);
149 //
150 // for ( int j = 1; j <= 6; j++ ) {
151 // SumStrain.at(j) += VecStrain.at(j);
152 // SumStress.at(j) += VecStress.at(j);
153 // }
154 //
155 // //VecStrain.printYourself();
156 // //SumStrain.printYourself();
157 // }
158 // }
159 //
160 // //average
161 // SumStrain.times(1. / VolTot * scale);
162 // SumStress.times(1. / VolTot * scale);
163 // //SumStrain.printYourself();
164 // //SumStress.printYourself();
165 // answer.resize(6);
166 // answer = SumStress;
167 // //answer.printYourself();
168
169 // update gp
170 // status->letTempStrainVectorBe(totalStrain);
171 // status->letTempStressVectorBe(answer);
172}
173
174std::unique_ptr<MaterialStatus>
175MicroMaterial :: CreateStatus(GaussPoint *gp) const
176{
177 return std::make_unique<MicroMaterialStatus>(gp);
178}
179
180
181//answer must be of size 24x24 (linear brick 3*8=24)
182void MicroMaterial :: giveMacroStiffnessMatrix(FloatMatrix &answer, TimeStep *tStep, MatResponseMode rMode, const IntArray &microMasterNodes, const IntArray &microBoundaryNodes)
183{
184 int row, col;
185 double tmpDouble;
186 Domain *microDomain = this->problemMicro->giveDomain(1);
187 EngngModel *microEngngModel = microDomain->giveEngngModel();
188 DofManager *DofMan;
189 Dof *dof;
190
191 FloatMatrix Kbb; //contains reduced problem with boundary nodes and without interior nodes
192 FloatMatrix Kbi; //can be zero size if no internal DOFs
193 FloatMatrix Kii1KbiT; //can be zero size if no internal DOFs
194 FloatMatrix Kbb_1; //help matrix
195
196 FloatMatrix slaveMasterOnBoundary; //transformation matrix representing displacements on boundary tied to master nodes
197
198 SparseMtrxType sparseMtrxType = ( SparseMtrxType ) 0; //SMT_Skyline symmetric skyline
199 std::unique_ptr<SparseMtrx> stiffnessMatrixMicro; //full stiffness matrix without any constraint
200 std::unique_ptr<SparseMtrx> Kii; //submatrix of internal DOFs
201
203 Kbb.zero();
204
205 this->isDefaultNumbering = false; //total number of equations can be now set arbitrarily
207
208 //assemble sparse matrix K_ii of DOFS of internal nodes to be condensed out
209 //K_ii is generally large, inversion of FloatMatrix consumes a lot of time and memory, sparse matrix is used
211 printf( "Internal DOFs %d\n", this->giveRequiredNumberOfDomainEquation() );
213 if ( totalInternalDofs ) {
215 Kbi.zero();
217 Kii1KbiT.zero();
218 Kii = classFactory.createSparseMtrx(sparseMtrxType);
219 Kii->buildInternalStructure(microEngngModel, 1, * this);
220 microEngngModel->assemble(*Kii, tStep, TangentAssembler(rMode), * this, microDomain);
221 }
222
223
224
225 //build a full stiffness matrix for the extraction of submatrices
228
229 stiffnessMatrixMicro = classFactory.createSparseMtrx(sparseMtrxType);
230 stiffnessMatrixMicro->buildInternalStructure(microEngngModel, 1, * this);
231 microEngngModel->assemble(*stiffnessMatrixMicro, tStep, TangentAssembler(rMode), * this, microDomain);
232
233
234 for ( int i = 1; i <= totalBoundaryDofs; i++ ) {
235 for ( int j = 1; j <= totalBoundaryDofs; j++ ) { //Kbb
236 row = microBoundaryDofsArr.at(i);
237 col = microBoundaryDofsArr.at(j);
238 if ( stiffnessMatrixMicro->isAllocatedAt(row, col) ) {
239 //printf("%d %d ", row, col);
240 Kbb.at(i, j) = stiffnessMatrixMicro->at(row, col);
241 }
242 }
243
244 for ( int j = 1; j <= totalInternalDofs; j++ ) { //Kbi
245 row = microBoundaryDofsArr.at(i);
246 col = microInternalDofsArr.at(j);
247 if ( stiffnessMatrixMicro->isAllocatedAt(row, col) ) {
248 Kbi.at(i, j) = stiffnessMatrixMicro->at(row, col);
249 }
250 }
251 }
252
253 //Kbi.printYourself();
254
255 if ( totalInternalDofs ) {
256 FloatArray xVector( Kii->giveNumberOfColumns() );
257 //Kii->printYourself();
258 Kii->factorized();
259
260 for ( int i = 1; i <= totalInternalDofs; i++ ) {
261 xVector.zero();
262 xVector.at(i) = 1.;
263 Kii->backSubstitutionWith(xVector);
264 //xVector.printYourself();
265 for ( int b = 1; b <= totalBoundaryDofs; b++ ) { //do not store Kii^-1, it is a dense matrix, compute multiplication directly
266 tmpDouble = 0.;
267 for ( int j = 1; j <= totalInternalDofs; j++ ) {
268 tmpDouble += xVector.at(j) * Kbi.at(b, j);
269 }
270
271 Kii1KbiT.at(i, b) = tmpDouble;
272 }
273
274 //OOFEM_LOG_INFO("%d ", i);
275 }
276
277 OOFEM_LOG_INFO("\n");
278
279 Kbb_1.beProductOf( Kbi, Kii1KbiT );
280 for ( int i = 1; i <= totalBoundaryDofs; i++ ) {
281 for ( int j = 1; j <= totalBoundaryDofs; j++ ) {
282 Kbb.at(i, j) -= Kbb_1.at(i, j);
283 }
284 }
285 }
286
287 //Kbb_1-printYourself();
288 //OOFEM_ERROR("Stop");
289 this->isDefaultNumbering = true; //switch back to default numbering
291
292 // IntArray interiorDofNode; //equation numbers to be condensed out, sorted
293 // IntArray boundaryDofNode;//equation numbers in rows (or columns) of stiffness matrix without interior nodes, sorted
294 // interiorDofNode.clear();
295 // boundaryDofNode.clear();
296 // for ( i = 1; i <= microDomain->giveNumberOfDofManagers(); i++){
297 // DofMan = microDomain->giveDofManager(i);
298 // for ( Dof *dof: *DofMan ){
299 // eqNumber = giveDofEquationNumber(dof);
300 // if(microBoundaryNodes.contains( DofMan->giveGlobalNumber())){
301 // boundaryDofNode.followedBy(eqNumber);
302 // }
303 // else {
304 // interiorDofNode.followedBy(eqNumber);
305 // }
306 // }
307 // }
308 // //Tmp.beInverseOf(stiffnessMatrixMicroFloat);
309 // /*Perform static condensation of internal nodes, leave microBoundaryNodes which also contains microMasterNodes
310 // algorithm based on structuralelement.C, Rayleigh-Ritz method
311 // zeroes on particular rows and columns will appear in the FloatMatrix
312 // */
313 // int ii, k;
314 // int size = stiffnessMatrixMicroFloat.giveNumberOfRows();
315 // int ndofs = maxNumberOfDomainEquation;
316 // long double coeff, dii;
317 //
318 // for ( i = 1; i <= interiorDofNode.giveSize(); i++ ) {//how many DOFs will be condensed out
319 // ii = interiorDofNode.at(i);
320 // if ( ( ii > ndofs ) || ( ii <= 0 ) ) {
321 // OOFEM_ERROR("wrong DOF number");
322 // }
323 //
324 // dii = stiffnessMatrixMicroFloat.at(ii, ii);
325 //
326 // for ( j = 1; j <= size; j++ ) {
327 // coeff = -stiffnessMatrixMicroFloat.at(j, ii) / dii;
328 // if ( ii != j ) {
329 // for ( k = 1; k <= size; k++ ) {
330 // stiffnessMatrixMicroFloat.at(j, k) += stiffnessMatrixMicroFloat.at(ii, k) * coeff;
331 // }
332 // }
333 // }
334 //
335 // for ( k = 1; k <= size; k++ ) {
336 // stiffnessMatrixMicroFloat.at(ii, k) = 0.;
337 // stiffnessMatrixMicroFloat.at(k, ii) = 0.;
338 // }
339 // }
340 // //stiffnessMatrixMicroFloat.printYourself();
341 //
342 // //copy non-zeroed rows and columns to reduced stiffness matrix
343 // stiffnessMatrixMicroReducedFloat.resize(boundaryDofNode.giveSize(), boundaryDofNode.giveSize());
344 // stiffnessMatrixMicroReducedFloat.zero();
345 // for ( i = 1; i <= boundaryDofNode.giveSize(); i++ ) {
346 // for ( j = 1; j <= boundaryDofNode.giveSize(); j++ ) {
347 // stiffnessMatrixMicroReducedFloat.at(i,j) = stiffnessMatrixMicroFloat.at( boundaryDofNode.at(i), boundaryDofNode.at(j) );
348 // }
349 // }
350 //stiffnessMatrixMicroReducedFloat.printYourself();
351
352 //build the transformation matrix of size (x,24) relating slave nodes on the boundary to master nodes on the boundary
353 slaveMasterOnBoundary.resize(Kbb.giveNumberOfColumns(), 24);
354 slaveMasterOnBoundary.zero();
355
356 FloatArray n(8);
357 int node, nodePos;
358
359 //master nodes
360 // for ( i = 1; i <= microMasterNodes.giveSize(); i++ ) {//8 nodes
361 // node = microMasterNodes.at(i);
362 // DofMan = microDomain->giveDofManager(node);
363 // dof = DofMan->giveDofWithID(D_u);
364 // j = boundaryDofNode.findFirstIndexOf(giveDofEquationNumber(dof));
365 // if(!j)
366 // OOFEM_ERROR("Not found equation number %d in reduced stiffness matrix of node %d\n", giveDofEquationNumber(dof), DofMan->giveGlobalNumber());
367 // microMasterNodesLoc.followedBy(j);
368 // }
369
370
371 //boundary nodes
372 for ( int i = 1; i <= microBoundaryNodes.giveSize(); i++ ) {
373 node = microBoundaryNodes.at(i);
374 DofMan = microDomain->giveDofManager(node);
375 dof = DofMan->giveDofWithID(D_u);
376 nodePos = microBoundaryDofsArr.findFirstIndexOf( giveDofEquationNumber(dof) ); //row(column) of reduced stiffness matrix
377 if ( !nodePos ) {
378 OOFEM_ERROR("Not found equation number %d in reduced stiffness matrix of node %d\n", giveDofEquationNumber(dof), DofMan->giveGlobalNumber() );
379 }
380
381 this->macroLSpaceElement->evalInterpolation( n, this->microMasterCoords, DofMan->giveCoordinates() );
382
383 //construct transformation matrix relating displacement of slave boundary nodes to macroelement nodes
384 //the answer is returned to macroelement so the columns correspond to x,y,z DOFs of each macroelement node
385
386 for ( int j = 1; j <= this->macroLSpaceElement->giveNumberOfNodes(); j++ ) { //linhex - 8 nodes
387 for ( int k = 0; k < 3; k++ ) { //the same interpolation for x,y,z
388 row = nodePos + k;
389 col = 3 * j + k - 2;
390 if ( row > slaveMasterOnBoundary.giveNumberOfRows() ) {
391 OOFEM_ERROR("Row is outside the reduced stiffness matrix ");
392 }
393
394 if ( col > slaveMasterOnBoundary.giveNumberOfColumns() ) {
395 OOFEM_ERROR("Column is outside the reduced stiffness matrix ");
396 }
397
398 slaveMasterOnBoundary.at(row, col) = n.at(j);
399 }
400 }
401 }
402
403 //slaveMasterOnBoundary.printYourself();
404# if 0
405 //check of the transformation matrix - the sum of each third column must be either zero or one
406 double sum;
407 for ( int i = 1; i <= slaveMasterOnBoundary.giveNumberOfRows(); i++ ) {
408 sum = 0;
409 for ( int j = 1; j <= slaveMasterOnBoundary.giveNumberOfColumns(); j++ ) {
410 if ( j % 3 == 0 ) {
411 sum += slaveMasterOnBoundary.at(i, j);
412 }
413 }
414
415 OOFEM_LOG_INFO("Sum of %i row of transformation matrix row %f\n", i, sum);
416 }
417
418# endif
419
420 //slaveMasterOnBoundary.printYourself();
421
422 //perform K(24,24) = T^T * K * T
423 FloatMatrix A;
424 A.beProductOf(Kbb, slaveMasterOnBoundary);
425 //A.printYourself();
426 //slaveMasterOnBoundary.printYourself();
427 //answer.resize(24, 24);
428 //answer.zero();
429 //OOFEM_ERROR("Stop");
430 answer.beTProductOf(slaveMasterOnBoundary, A);
431}
432
433//should be called before the calculation of micromaterial
434void MicroMaterial :: setMacroProperties(Domain *macroDomain, MacroLSpace *macroLSpaceElement, const IntArray &microMasterNodes, const IntArray &microBoundaryNodes)
435{
436 Domain *microDomain = this->problemMicro->giveDomain(1);
437 EngngModel *microEngngModel = microDomain->giveEngngModel();
438 MetaStep *mstep;
439 DofManager *DofMan;
440
441 int numDofs, numDofMan;
442 int counterDefault = 1, counterBoundary = 1, counterInternal = 1;
443
444 this->macroDomain = macroDomain;
445 this->macroLSpaceElement = macroLSpaceElement;
446
447 this->microMasterCoords.clear();
448 this->microMasterCoords.reserve(microMasterNodes.giveSize());
449 for ( int nodeNum: microMasterNodes ) { //8 nodes
450 this->microMasterCoords.push_back( microDomain->giveNode( nodeNum )->giveCoordinates() );
451 }
452
453 microEngngModel->giveNextStep(); //set the first time step
454 mstep = microEngngModel->giveCurrentMetaStep();
455 mstep->setNumberOfSteps(this->macroDomain->giveEngngModel()->giveMetaStep(1)->giveNumberOfSteps() + 1);
456
457 //separate DOFs into boundary and internal
458 this->NumberOfDofManagers = microDomain->giveNumberOfDofManagers();
462 for ( int i = 0; i < this->NumberOfDofManagers; i++ ) {
463 microBoundaryDofs [ i ].resize(3);
464 microInternalDofs [ i ].resize(3);
465 microDefaultDofs [ i ].resize(3);
466 microBoundaryDofs [ i ].zero();
467 microInternalDofs [ i ].zero();
468 microDefaultDofs [ i ].zero();
469 }
470
471 for ( int i = 0; i < this->NumberOfDofManagers; i++ ) {
472 DofMan = microDomain->giveDofManager(i + 1);
473 numDofMan = DofMan->giveGlobalNumber();
474 numDofs = DofMan->giveNumberOfDofs();
475 if ( numDofs != 3 ) {
476 OOFEM_ERROR("Node %d does not have three degrees of freedom", DofMan->giveGlobalNumber() );
477 }
478
479 for ( int j = 0; j < numDofs; j++ ) {
480 microDefaultDofs [ i ] [ j ] = counterDefault; //equation number
481 if ( microBoundaryNodes.contains(numDofMan) ) { //boundary (and master) nodes
482 microBoundaryDofs [ i ] [ j ] = counterBoundary;
483 microBoundaryDofsArr.followedBy(counterDefault);
484 counterBoundary++;
485 } else { //internal nodes to be condensed out
486 microInternalDofs [ i ] [ j ] = counterInternal;
487 microInternalDofsArr.followedBy(counterDefault);
488 counterInternal++;
489 }
490
491 counterDefault++;
492 }
493 }
494
495 this->totalBoundaryDofs = counterBoundary - 1;
496 this->totalInternalDofs = counterInternal - 1;
497
498 // //the pointer to underlying problem is problemMicro
499 // DofManager *DofMan;
500 // int i, j, counter = 0;
501 // IntArray ut, dofIDArry = {D_u, D_v, D_w};1
502 //
503 // for ( i = 1; i <= problemMicro->giveDomain(1)->giveNumberOfDofManagers(); i++ ) { //for each node
504 // DofMan = problemMicro->giveDomain(1)->giveDofManager(i);
505 // //printf("%d\n",DofMan->giveNumberOfPrimaryMasterDofs(dofIDArry));
506 // //DofMan->giveLocationArray(dofIDArry, *this );
507 // DofMan->giveCompleteLocationArray(* this);
508 // for ( j = 1; j <= DofMan->giveNumberOfDofs(); j++ ) {
509 // counter++;
510 // }
511 // }
512
513 this->maxNumberOfDomainEquation = counterBoundary - 1 + counterInternal - 1;
514}
515
516
517
518//from class UnknownNumberingScheme
519void MicroMaterial :: init(void)
520{
521 if ( this->problemMicro->giveNumberOfDomains() != 1 ) {
522 OOFEM_ERROR("Number of domains on microproblem is greater than 1");
523 }
524
525 microMatIsUsed = true;
526}
527
528//from class UnknownNumberingScheme
529//Each node has three dofs (x,y,z direction)
530int MicroMaterial :: giveDofEquationNumber(Dof *dof) const
531{
532 int numDofMan = dof->giveDofManNumber();
533 int numDof = dof->giveDofID(); // Note: Relieas on D_u, D_v, D_w being 1, 2, 3
534
535 //depending on the assembly of submatrix, swith to equation numbering
536 switch ( DofEquationNumbering ) {
537 case AllNodes: //the default
538 return microDefaultDofs [ numDofMan - 1 ] [ numDof - 1 ];
539 break;
540 case BoundaryNodes:
541 return microBoundaryDofs [ numDofMan - 1 ] [ numDof - 1 ];
542 break;
543 case InteriorNodes:
544 return microInternalDofs [ numDofMan - 1 ] [ numDof - 1 ];
545 break;
546 default:
547 OOFEM_ERROR("Node numbering undefined");
548 }
549}
550
551//from class UnknownNumberingScheme
552int MicroMaterial :: giveRequiredNumberOfDomainEquation() const
553{
554 return this->reqNumberOfDomainEquation;
555}
556} // end namespace oofem
#define REGISTER_Material(class)
int giveGlobalNumber() const
Definition dofmanager.h:515
int giveNumberOfDofs() const
Definition dofmanager.C:287
const FloatArray & giveCoordinates() const
Definition dofmanager.h:390
Dof * giveDofWithID(int dofID) const
Definition dofmanager.C:127
DofIDItem giveDofID() const
Definition dof.h:276
int giveDofManNumber() const
Definition dof.C:72
int giveNumberOfDofManagers() const
Returns number of dof managers in domain.
Definition domain.h:461
DofManager * giveDofManager(int n)
Definition domain.C:317
Node * giveNode(int n)
Definition domain.h:398
EngngModel * giveEngngModel()
Definition domain.C:419
virtual TimeStep * giveNextStep()
Returns next time step (next to current step) of receiver.
Definition engngm.h:747
MetaStep * giveCurrentMetaStep()
Returns current meta step.
Definition engngm.C:1900
MetaStep * giveMetaStep(int i)
Returns the i-th meta step.
Definition engngm.h:773
virtual void assemble(SparseMtrx &answer, TimeStep *tStep, const MatrixAssembler &ma, const UnknownNumberingScheme &s, Domain *domain)
Definition engngm.C:889
double & at(Index i)
Definition floatarray.h:202
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
int giveNumberOfColumns() const
Returns number of columns of receiver.
void beProductOf(const FloatMatrix &a, const FloatMatrix &b)
void zero()
Zeroes all coefficient of receiver.
int giveNumberOfRows() const
Returns number of rows of receiver.
double at(std::size_t i, std::size_t j) const
void beTProductOf(const FloatMatrix &a, const FloatMatrix &b)
bool contains(int value) const
Definition intarray.h:292
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
GaussPoint * gp
Associated integration point.
void setNumberOfSteps(int newNumberOfSteps)
Sets the number of steps within the metastep.
Definition metastep.C:105
int giveNumberOfSteps()
Returns number of Steps it represent.
Definition metastep.h:118
int totalBoundaryDofs
Number of equations associated with boundary nodes.
IntArray microBoundaryDofsArr
Array of equation numbers associated to boundary nodes.
int NumberOfDofManagers
Number of DOF Managers.
std::vector< IntArray > microDefaultDofs
Array containing default equation numbers for all nodes [DofManagerNumber][DOF].
MacroLSpace * macroLSpaceElement
Pointer to the macroscale element.
std::unique_ptr< EngngModel > problemMicro
Pointer to the underlying micro problem.
int totalInternalDofs
Number of equations associated with boundary nodes.
Domain * macroDomain
Pointer to the macroscale domain.
std::vector< IntArray > microBoundaryDofs
Array containing equation numbers for boundary nodes [DofManagerNumber][DOF].
std::vector< IntArray > microInternalDofs
Array containing equation numbers for internal nodes to be condensed out [DofManagerNumber][DOF].
int reqNumberOfDomainEquation
Required number of domain equations.
std::string inputFileNameMicro
int giveDofEquationNumber(Dof *dof) const override
bool microMatIsUsed
Flag signalizing whether micromaterial is used by other element.
EquationNumbering DofEquationNumbering
int maxNumberOfDomainEquation
The maximum DOFs corresponding to released all of the boundary conditions.
int giveRequiredNumberOfDomainEquation() const override
IntArray microInternalDofsArr
Array of equation numbers associated to internal nodes.
std::vector< FloatArray > microMasterCoords
Array containing coordinates of 8 master nodes of microproblem.
StructuralMaterialStatus(GaussPoint *g)
Constructor. Creates new StructuralMaterialStatus with IntegrationPoint g.
StructuralMaterial(int n, Domain *d)
#define OOFEM_ERROR(...)
Definition error.h:79
#define IR_GIVE_FIELD(__ir, __value, __id)
Definition inputrecord.h:67
#define OOFEM_LOG_INFO(...)
Definition logger.h:143
#define _IFT_MicroMaterial_fileName
long ContextMode
Definition contextmode.h:43
std::unique_ptr< EngngModel > InstanciateProblem(DataReader &dr, problemMode mode, int contextFlag, EngngModel *_master, bool parallelFlag)
Definition util.C:153
double sum(const FloatArray &x)
ClassFactory & classFactory
@ _processor
Definition problemmode.h:40

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011