OOFEM 3.0
Loading...
Searching...
No Matches
planstrssxfem.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
40#include "xfem/enrichmentitem.h"
41#include "vtkxmlexportmodule.h"
42#include "dynamicinputrecord.h"
43#include "feinterpol.h"
44#include "classfactory.h"
45#include "mathfem.h"
46
47#ifdef __OOFEG
48 #include "oofeggraphiccontext.h"
49#endif
50
51namespace oofem {
53
54void PlaneStress2dXfem :: updateYourself(TimeStep *tStep)
55{
56 PlaneStress2d :: updateYourself(tStep);
57 XfemStructuralElementInterface :: updateYourselfCZ(tStep);
58}
59
60void PlaneStress2dXfem :: postInitialize()
61{
62 PlaneStress2d :: postInitialize();
63 XfemStructuralElementInterface :: initializeCZMaterial();
64}
65
67PlaneStress2dXfem :: giveInterface(InterfaceType it)
68{
69 if ( it == XfemElementInterfaceType ) {
70 return static_cast< XfemElementInterface * >(this);
71 } else if ( it == VTKXMLExportModuleElementInterfaceType ) {
72 return static_cast< VTKXMLExportModuleElementInterface * >(this);
73 } else {
74 return PlaneStress2d :: giveInterface(it);
75 }
76}
77
78
79void
80PlaneStress2dXfem :: computeGaussPoints()
81{
82 if ( this->giveDomain()->hasXfemManager() ) {
83 XfemManager *xMan = this->giveDomain()->giveXfemManager();
84
85 if ( xMan->isElementEnriched(this) ) {
87
88 // Blending element
89// increaseNumGP();
90 PlaneStress2d :: computeGaussPoints();
91 }
92 } else {
93 PlaneStress2d :: computeGaussPoints();
94 }
95 } else {
96 PlaneStress2d :: computeGaussPoints();
97 }
98}
99
100void PlaneStress2dXfem :: increaseNumGP()
101{
102 switch(numberOfGaussPoints)
103 {
104 case(4):
106 break;
107 case(9):
109 break;
110 case(16):
112 break;
113 }
114}
115
116void PlaneStress2dXfem :: computeBmatrixAt(GaussPoint *gp, FloatMatrix &answer, int li, int ui)
117{
118 XfemElementInterface_createEnrBmatrixAt(answer, * gp, * this);
119}
120
121void PlaneStress2dXfem :: computeBHmatrixAt(GaussPoint *gp, FloatMatrix &answer)
122{
123 XfemElementInterface_createEnrBHmatrixAt(answer, * gp, * this);
124}
125
126
127void PlaneStress2dXfem :: computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer)
128{
129 XfemElementInterface_createEnrNmatrixAt(answer, iLocCoord, * this, false);
130}
131
132
133
134int PlaneStress2dXfem :: computeNumberOfDofs()
135{
136 int ndofs = 0;
137 for ( int i = 1; i <= this->giveNumberOfDofManagers(); i++ ) {
138 ndofs += this->giveDofManager(i)->giveNumberOfDofs();
139 }
140
141 return ndofs;
142}
143
144
145void
146PlaneStress2dXfem :: giveDofManDofIDMask(int inode, IntArray &answer) const
147{
148 // Continuous part
149 PlaneStress2d :: giveDofManDofIDMask(inode, answer);
150
151 // Discontinuous part
152 if( this->giveDomain()->hasXfemManager() ) {
153 DofManager *dMan = giveDofManager(inode);
154 XfemManager *xMan = giveDomain()->giveXfemManager();
155
156 const std::vector<int> &nodeEiIndices = xMan->giveNodeEnrichmentItemIndices( dMan->giveGlobalNumber() );
157 for ( size_t i = 0; i < nodeEiIndices.size(); i++ ) {
158 EnrichmentItem *ei = xMan->giveEnrichmentItem(nodeEiIndices[i]);
159 if ( ei->isDofManEnriched(* dMan) ) {
160 IntArray eiDofIdArray;
161 ei->computeEnrichedDofManDofIdArray(eiDofIdArray, *dMan);
162 answer.followedBy(eiDofIdArray);
163 }
164 }
165 }
166}
167
168
169
170void PlaneStress2dXfem :: computeConstitutiveMatrixAt(FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
171{
172 XfemStructuralElementInterface :: XfemElementInterface_computeConstitutiveMatrixAt(answer, rMode, gp, tStep);
173}
174
175void
176PlaneStress2dXfem :: computeStressVector(FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
177{
178 XfemStructuralElementInterface :: XfemElementInterface_computeStressVector(answer, strain, gp, tStep);
179}
180
181void PlaneStress2dXfem :: computeStiffnessMatrix(FloatMatrix &answer, MatResponseMode rMode, TimeStep *tStep)
182{
183 PlaneStress2d :: computeStiffnessMatrix(answer, rMode, tStep);
184 XfemStructuralElementInterface :: computeCohesiveTangent(answer, tStep);
185
186 const double tol = 1.0e-6;
187 const double regularizationCoeff = 1.0e-6;
188 int numRows = answer.giveNumberOfRows();
189 for(int i = 0; i < numRows; i++) {
190 if( fabs(answer(i,i)) < tol ) {
191 answer(i,i) += regularizationCoeff;
192 //printf("Found zero on diagonal.\n");
193 }
194 }
195}
196
197void PlaneStress2dXfem :: computeDeformationGradientVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
198{
200}
201
202void
203PlaneStress2dXfem :: giveInternalForcesVector(FloatArray &answer, TimeStep *tStep, int useUpdatedGpRecord)
204{
205 PlaneStress2d :: giveInternalForcesVector(answer, tStep, useUpdatedGpRecord);
206 XfemStructuralElementInterface :: computeCohesiveForces(answer, tStep);
207}
208
210PlaneStress2dXfem :: giveGeometryType() const
211{
212 if ( this->giveDomain()->hasXfemManager() ) {
213 XfemManager *xMan = this->giveDomain()->giveXfemManager();
214 if ( xMan->isElementEnriched(this) ) {
215 return EGT_Composite;
216 //return EGT_quad_1;
217 } else {
218 return EGT_Composite;
219 //return EGT_quad_1;
220 }
221 } else {
222 return EGT_quad_1;
223 }
224}
225
226
227
228#ifdef __OOFEG
229void PlaneStress2dXfem :: drawRawGeometry(oofegGraphicContext &gc, TimeStep *tStep)
230{
231 if ( !gc.testElementGraphicActivity(this) ) {
232 return;
233 }
234#if 0
235 XfemManager *xf = this->giveDomain()->giveXfemManager();
236 if ( !xf->isElementEnriched(this) ) {
237 PlaneStress2d :: drawRawGeometry(gc);
238 } else {
239 if ( integrationRulesArray.size() > 1 ) {
240 // TODO: Implement visualization
241 for ( auto &iRule: integrationRulesArray ) {
242 PatchIntegrationRule *piRule = dynamic_cast< PatchIntegrationRule * >( ir );
243 if ( piRule ) {
244 piRule->givePatch()->draw(context);
245 }
246 }
247 } else {
248 PlaneStress2d :: drawRawGeometry(gc);
249 }
250 }
251#endif
252}
253
254void PlaneStress2dXfem :: drawScalar(oofegGraphicContext &gc, TimeStep *tStep)
255{
256 if ( !gc.testElementGraphicActivity(this) ) {
257 return;
258 }
259#if 0
260 XfemManager *xf = this->giveDomain()->giveXfemManager();
261 if ( !xf->isElementEnriched(this) ) {
262 PlaneStress2d :: drawScalar(context);
263 } else {
264 if ( gc.giveIntVarMode() == ISM_local ) {
265 int indx;
266 double val;
267 FloatArray s(3), v;
268
269 indx = gc.giveIntVarIndx();
270
271 for ( auto &ir: integrationRulesArray ) {
272 PatchIntegrationRule *iRule = dynamic_cast< PatchIntegrationRule * >(ir);
273
274 #if 0
275 val = iRule->giveMaterial();
276 #else
277 val = 0.0;
278 for ( GaussPoint *gp: *iRule ) {
279 giveIPValue(v, gp, gc.giveIntVarType(), tStep);
280 val += v.at(indx);
281 }
282
283 val /= iRule->giveNumberOfIntegrationPoints();
284 #endif
285 s.at(1) = s.at(2) = s.at(3) = val;
286 // TODO: Implement visualization
287 // iRule->givePatch()->drawWD(context, s);
288 }
289 } else {
290 PlaneStress2d :: drawScalar(context);
291 }
292 }
293#endif
294}
295#endif
296
297
298void
299PlaneStress2dXfem :: initializeFrom(InputRecord &ir, int priority)
300{
301 PlaneStress2d :: initializeFrom(ir, priority);
302 XfemStructuralElementInterface :: initializeCZFrom(ir, priority);
303}
304
305MaterialMode PlaneStress2dXfem :: giveMaterialMode()
306{
307 return XfemStructuralElementInterface :: giveMaterialMode();
308}
309
310void PlaneStress2dXfem :: giveInputRecord(DynamicInputRecord &input)
311{
312 PlaneStress2d :: giveInputRecord(input);
313 XfemStructuralElementInterface :: giveCZInputRecord(input);
314}
315
316void
317PlaneStress2dXfem :: giveCompositeExportData(std::vector< ExportRegion > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep)
318{
319 vtkPieces.resize(1);
320
321 const int numCells = mSubTri.size();
322
323 if(numCells == 0) {
324 // Enriched but uncut element
325 // Visualize as a quad
326 vtkPieces[0].setNumberOfCells(1);
327
328 int numTotalNodes = 4;
329 vtkPieces[0].setNumberOfNodes(numTotalNodes);
330
331 // Node coordinates
332 std :: vector< FloatArray >nodeCoords;
333 for(int i = 1; i <= 4; i++) {
334 const auto &x = giveDofManager(i)->giveCoordinates();
335 nodeCoords.push_back(x);
336
337 vtkPieces[0].setNodeCoords(i, x);
338 }
339
340 // Connectivity
341 IntArray nodes1 = {1, 2, 3, 4};
342 vtkPieces[0].setConnectivity(1, nodes1);
343
344 // Offset
345 int offset = 4;
346 vtkPieces[0].setOffset(1, offset);
347
348 // Cell types
349 vtkPieces[0].setCellType(1, 9); // Linear quad
350
351
352
353
354 // Export nodal variables from primary fields
355 vtkPieces[0].setNumberOfPrimaryVarsToExport(primaryVarsToExport, numTotalNodes);
356
357 for ( int fieldNum = 1; fieldNum <= primaryVarsToExport.giveSize(); fieldNum++ ) {
358 UnknownType type = ( UnknownType ) primaryVarsToExport.at(fieldNum);
359
360 for ( int nodeInd = 1; nodeInd <= numTotalNodes; nodeInd++ ) {
361
362 if ( type == DisplacementVector ) { // compute displacement
363
364 FloatArray u = {0.0, 0.0, 0.0};
365
366 // Fetch global coordinates (in undeformed configuration)
367 const FloatArray &x = nodeCoords[nodeInd-1];
368
369 // Compute local coordinates
370 FloatArray locCoord;
371 computeLocalCoordinates(locCoord, x);
372
373 // Compute displacement in point
374 FloatMatrix NMatrix;
375 computeNmatrixAt(locCoord, NMatrix);
376 FloatArray solVec;
377 computeVectorOf(VM_Total, tStep, solVec);
378 FloatArray uTemp;
379 uTemp.beProductOf(NMatrix, solVec);
380
381 if(uTemp.giveSize() == 3) {
382 u = uTemp;
383 }
384 else {
385 u = {uTemp[0], uTemp[1], 0.0};
386 }
387
388 vtkPieces[0].setPrimaryVarInNode(type, nodeInd, u);
389 } else {
390 printf("fieldNum: %d\n", fieldNum);
391 // TODO: Implement
392// ZZNodalRecoveryMI_recoverValues(values, layer, ( InternalStateType ) 1, tStep); // does not work well - fix
393// for ( int j = 1; j <= numCellNodes; j++ ) {
394// vtkPiece.setPrimaryVarInNode(fieldNum, nodeNum, values [ j - 1 ]);
395// nodeNum += 1;
396// }
397 }
398 }
399 }
400
401
402 // Export nodal variables from internal fields
403 vtkPieces[0].setNumberOfInternalVarsToExport(internalVarsToExport, numTotalNodes);
404
405
406 // Export cell variables
407 vtkPieces[0].setNumberOfCellVarsToExport(cellVarsToExport, 1);
408 for ( int i = 1; i <= cellVarsToExport.giveSize(); i++ ) {
409 InternalStateType type = ( InternalStateType ) cellVarsToExport.at(i);
410 FloatArray average;
411 std :: unique_ptr< IntegrationRule >&iRule = integrationRulesArray [ 0 ];
412 VTKXMLExportModule :: computeIPAverage(average, iRule.get(), this, type, tStep);
413
414 FloatArray averageVoigt;
415
416 if( average.giveSize() == 6 ) {
417
418 averageVoigt.resize(9);
419
420 averageVoigt.at(1) = average.at(1);
421 averageVoigt.at(5) = average.at(2);
422 averageVoigt.at(9) = average.at(3);
423 averageVoigt.at(6) = averageVoigt.at(8) = average.at(4);
424 averageVoigt.at(3) = averageVoigt.at(7) = average.at(5);
425 averageVoigt.at(2) = averageVoigt.at(4) = average.at(6);
426 }
427 else {
428 if(average.giveSize() == 1) {
429 averageVoigt.resize(1);
430 averageVoigt.at(1) = average.at(1);
431 }
432 }
433
434 vtkPieces[0].setCellVar( type, 1, averageVoigt );
435 }
436
437
438 // Export of XFEM related quantities
439 if ( domain->hasXfemManager() ) {
440 XfemManager *xMan = domain->giveXfemManager();
441
442 int nEnrIt = xMan->giveNumberOfEnrichmentItems();
443 vtkPieces[0].setNumberOfInternalXFEMVarsToExport(xMan->vtkExportFields.giveSize(), nEnrIt, numTotalNodes);
444
445 const int nDofMan = giveNumberOfDofManagers();
446
447
448 for ( int field = 1; field <= xMan->vtkExportFields.giveSize(); field++ ) {
449 XFEMStateType xfemstype = ( XFEMStateType ) xMan->vtkExportFields [ field - 1 ];
450
451 for ( int enrItIndex = 1; enrItIndex <= nEnrIt; enrItIndex++ ) {
452 EnrichmentItem *ei = xMan->giveEnrichmentItem(enrItIndex);
453 for ( int nodeInd = 1; nodeInd <= numTotalNodes; nodeInd++ ) {
454
455 const FloatArray &x = nodeCoords[nodeInd-1];
456 FloatArray locCoord;
457 computeLocalCoordinates(locCoord, x);
458
461 interp->evalN( N, locCoord, FEIElementGeometryWrapper(this) );
462
463
464 if ( xfemstype == XFEMST_LevelSetPhi ) {
465 double levelSet = 0.0, levelSetInNode = 0.0;
466
467 for(int elNodeInd = 1; elNodeInd <= nDofMan; elNodeInd++) {
468 DofManager *dMan = giveDofManager(elNodeInd);
469 ei->evalLevelSetNormalInNode(levelSetInNode, dMan->giveGlobalNumber(), dMan->giveCoordinates() );
470
471 levelSet += N.at(elNodeInd)*levelSetInNode;
472 }
473
474
475 FloatArray valueArray = {levelSet};
476 vtkPieces[0].setInternalXFEMVarInNode(field, enrItIndex, nodeInd, valueArray);
477
478 } else if ( xfemstype == XFEMST_LevelSetGamma ) {
479 double levelSet = 0.0, levelSetInNode = 0.0;
480
481 for(int elNodeInd = 1; elNodeInd <= nDofMan; elNodeInd++) {
482 DofManager *dMan = giveDofManager(elNodeInd);
483 ei->evalLevelSetTangInNode(levelSetInNode, dMan->giveGlobalNumber(), dMan->giveCoordinates() );
484
485 levelSet += N.at(elNodeInd)*levelSetInNode;
486 }
487
488
489 FloatArray valueArray = {levelSet};
490 vtkPieces[0].setInternalXFEMVarInNode(field, enrItIndex, nodeInd, valueArray);
491
492 } else if ( xfemstype == XFEMST_NodeEnrMarker ) {
493 double nodeEnrMarker = 0.0, nodeEnrMarkerInNode = 0.0;
494
495 for(int elNodeInd = 1; elNodeInd <= nDofMan; elNodeInd++) {
496 DofManager *dMan = giveDofManager(elNodeInd);
497 ei->evalNodeEnrMarkerInNode(nodeEnrMarkerInNode, dMan->giveGlobalNumber() );
498
499 nodeEnrMarker += N.at(elNodeInd)*nodeEnrMarkerInNode;
500 }
501
502
503 FloatArray valueArray = {nodeEnrMarker};
504 vtkPieces[0].setInternalXFEMVarInNode(field, enrItIndex, nodeInd, valueArray);
505 }
506
507 }
508 }
509 }
510 }
511
512 }
513 else {
514 // Enriched and cut element
515
516 XfemStructuralElementInterface::giveSubtriangulationCompositeExportData(vtkPieces, primaryVarsToExport, internalVarsToExport, cellVarsToExport, tStep);
517
518
519 }
520}
521
522} // end namespace oofem
#define N(a, b)
#define REGISTER_Element(class)
int giveGlobalNumber() const
Definition dofmanager.h:515
const FloatArray & giveCoordinates() const
Definition dofmanager.h:390
void computeVectorOf(ValueModeType u, TimeStep *tStep, FloatArray &answer)
Definition element.C:103
std::vector< std ::unique_ptr< IntegrationRule > > integrationRulesArray
Definition element.h:157
int numberOfGaussPoints
Definition element.h:175
virtual int giveNumberOfDofManagers() const
Definition element.h:695
DofManager * giveDofManager(int i) const
Definition element.C:553
virtual bool computeLocalCoordinates(FloatArray &answer, const FloatArray &gcoords)
Definition element.C:1240
bool evalNodeEnrMarkerInNode(double &oNodeEnrMarker, int iNodeInd) const
bool evalLevelSetTangInNode(double &oLevelSet, int iNodeInd, const FloatArray &iGlobalCoord) const
bool isDofManEnriched(const DofManager &iDMan) const
virtual void computeEnrichedDofManDofIdArray(IntArray &oDofIdArray, DofManager &iDMan)
bool evalLevelSetNormalInNode(double &oLevelSet, int iNodeInd, const FloatArray &iGlobalCoord) const
virtual void evalN(FloatArray &answer, const FloatArray &lcoords, const FEICellGeometry &cellgeo) const =0
Domain * giveDomain() const
Definition femcmpnn.h:97
Domain * domain
Link to domain object, useful for communicating with other FEM components.
Definition femcmpnn.h:79
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
Index giveSize() const
Returns the size of receiver.
Definition floatarray.h:261
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
int giveNumberOfRows() const
Returns number of rows of receiver.
void followedBy(const IntArray &b, int allocChunk=0)
Definition intarray.C:94
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
int giveNumberOfIntegrationPoints() const
void computeNmatrixAt(const FloatArray &iLocCoord, FloatMatrix &answer) override
FEInterpolation * giveInterpolation() const override
Definition planstrss.C:76
int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
void XfemElementInterface_createEnrBmatrixAt(FloatMatrix &oAnswer, GaussPoint &iGP, Element &iEl)
Creates enriched B-matrix.
void XfemElementInterface_createEnrBHmatrixAt(FloatMatrix &oAnswer, GaussPoint &iGP, Element &iEl)
Creates enriched BH-matrix.
void XfemElementInterface_createEnrNmatrixAt(FloatMatrix &oAnswer, const FloatArray &iLocCoord, Element &iEl, bool iSetDiscontContribToZero)
Creates enriched N-matrix.
bool isElementEnriched(const Element *elem)
IntArray vtkExportFields
const std ::vector< int > & giveNodeEnrichmentItemIndices(int iNodeIndex) const
EnrichmentItem * giveEnrichmentItem(int n)
int giveNumberOfEnrichmentItems() const
virtual void XfemElementInterface_computeDeformationGradientVector(FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
void giveSubtriangulationCompositeExportData(std ::vector< ExportRegion > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep)
VTK Interface.
bool XfemElementInterface_updateIntegrationRule() override
Updates integration rule based on the triangulation.
XFEMStateType
Definition xfemmanager.h:92
@ XfemElementInterfaceType
@ VTKXMLExportModuleElementInterfaceType
oofem::oofegGraphicContext gc[OOFEG_LAST_LAYER]

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