OOFEM 3.0
Loading...
Searching...
No Matches
vtkexportmodule.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/* COMPONENTS in TENSORS like stress or strain
36 * x y z
37 * x 1 6 5
38 * y 6 2 4
39 * z 5 4 3
40 *
41 * PARAVIEW - stresses and strains in global c.s., damage tensor in local c.s.
42 * x y z
43 * x 0 1 2
44 * y 3 4 5
45 * z 6 7 8
46 *
47 */
48
49
50#include "vtkexportmodule.h"
51#include "timestep.h"
52#include "gausspoint.h"
53#include "engngm.h"
54#include "node.h"
55#include "materialinterface.h"
56#include "mathfem.h"
57#include "cltypes.h"
58#include "element.h"
59#include "material.h"
60#include "classfactory.h"
61#include "crosssection.h"
62#include "dof.h"
63
64#include <vector>
65
66namespace oofem {
68
69VTKExportModule :: VTKExportModule(int n, EngngModel *e) : ExportModule(n, e), internalVarsToExport(), primaryVarsToExport() { }
70
71
72VTKExportModule :: ~VTKExportModule() { }
73
74
75void
76VTKExportModule :: initializeFrom(InputRecord &ir)
77{
78 ExportModule :: initializeFrom(ir);
79
80 int val;
81
85
86 val = NodalRecoveryModel :: NRM_ZienkiewiczZhu;
88 stype = ( NodalRecoveryModel :: NodalRecoveryModelType ) val;
89}
90
91
92void
93VTKExportModule :: doOutput(TimeStep *tStep, bool forcedOutput)
94{
95 if ( !( testTimeStepOutput(tStep) || forcedOutput ) ) {
96 return;
97 }
98
99 FILE *stream = this->giveOutputStream(tStep);
100
101 fprintf(stream, "# vtk DataFile Version 2.0\n");
102 fprintf( stream, "Output for time %f\n", tStep->giveTargetTime() );
103 fprintf(stream, "ASCII\n");
104
105 fprintf(stream, "DATASET UNSTRUCTURED_GRID\n");
106
107
108 Domain *d = emodel->giveDomain(1);
109 int inode, nnodes = d->giveNumberOfDofManagers();
110 this->giveSmoother(); // make sure smoother is created
111
112 // output points
113
114 fprintf(stream, "POINTS %d double\n", nnodes);
115 int ireg = -1;
116 int regionDofMans;
118
119 // asemble local->global region map
120 this->initRegionNodeNumbering(map, regionDofMans, 0, d, ireg, 1);
121
122 OOFEM_LOG_DEBUG("vktexportModule: %d %d\n", nnodes, regionDofMans);
123 for ( inode = 1; inode <= regionDofMans; inode++ ) {
124 const auto &coords = d->giveNode( map.at(inode) )->giveCoordinates();
125 for ( double s : coords ) {
126 fprintf( stream, "%e ", s );
127 }
128
129 for ( int i = coords.giveSize() + 1; i <= 3; i++ ) {
130 fprintf(stream, "%e ", 0.0);
131 }
132
133 fprintf(stream, "\n");
134 }
135
136 int elemToProcess = 0;
137 int ncells, celllistsize = 0;
138 for ( auto &elem : d->giveElements() ) {
139 if ( elem->giveParallelMode() != Element_local ) {
140 continue;
141 }
142
143 elemToProcess++;
144 // element composed from same-type cells asumed
145 ncells = this->giveNumberOfElementCells( elem.get() );
146 celllistsize += ncells + ncells *this->giveNumberOfNodesPerCell( this->giveCellType ( elem.get() ) );
147 }
148
149 int nelemNodes;
150 int vtkCellType;
151 IntArray cellNodes;
152 // output cells
153 fprintf(stream, "\nCELLS %d %d\n", elemToProcess, celllistsize);
154
155 IntArray regionNodalNumbers(nnodes);
156
157 // assemble global->local map
158 this->initRegionNodeNumbering(regionNodalNumbers, regionDofMans, 0, d, ireg, 0);
159 for ( auto &elem : d->giveElements() ) {
160 if ( elem->giveParallelMode() != Element_local ) {
161 continue;
162 }
163
164 vtkCellType = this->giveCellType(elem.get());
165
166 nelemNodes = this->giveNumberOfNodesPerCell(vtkCellType); //elem->giveNumberOfNodes(); // It HAS to be the same size as giveNumberOfNodesPerCell, otherwise the file will be incorrect.
167 this->giveElementCell(cellNodes, elem.get(), 0);
168 fprintf(stream, "%d ", nelemNodes);
169 for ( int i = 1; i <= nelemNodes; i++ ) {
170 fprintf(stream, "%d ", regionNodalNumbers.at( cellNodes.at(i) ) - 1);
171 }
172
173 fprintf(stream, "\n");
174 }
175
176 // output cell types
177 fprintf(stream, "\nCELL_TYPES %d\n", elemToProcess);
178 for ( auto &elem : d->giveElements() ) {
179 if ( elem->giveParallelMode() != Element_local ) {
180 continue;
181 }
182
183 vtkCellType = this->giveCellType(elem.get());
184 fprintf(stream, "%d\n", vtkCellType);
185 }
186
187 // output cell data (Material ID ...)
188 if ( cellVarsToExport.giveSize() ) {
189 exportCellVars(stream, elemToProcess, tStep);
190 }
191
192 if ( primaryVarsToExport.giveSize() || internalVarsToExport.giveSize() ) {
193 fprintf( stream, "\n\nPOINT_DATA %d\n", this->giveTotalRBRNumberOfNodes(d) );
194 }
195
196 this->exportPrimaryVars(stream, tStep);
197 this->exportIntVars(stream, tStep);
198
199 fclose(stream);
200}
201
202void
203VTKExportModule :: initialize()
204{
205 this->smoother.reset();
206 ExportModule :: initialize();
207}
208
209
210void
211VTKExportModule :: terminate()
212{ }
213
214
215FILE *
216VTKExportModule :: giveOutputStream(TimeStep *tStep)
217{
218 std :: string fileName;
219 FILE *answer;
220 fileName = this->giveOutputBaseFileName(tStep) + ".vtk";
221 if ( ( answer = fopen(fileName.c_str(), "w") ) == NULL ) {
222 OOFEM_ERROR("failed to open file %s", fileName.c_str());
223 }
224 return answer;
225}
226
227int
228VTKExportModule :: giveCellType(Element *elem)
229{
231 int vtkCellType = 0;
232 if ( elemGT == EGT_point ) {
233 vtkCellType = 1;
234 } else if ( elemGT == EGT_line_1 ) {
235 vtkCellType = 3;
236 } else if ( elemGT == EGT_line_2 ) {
237 vtkCellType = 21;
238 } else if ( elemGT == EGT_triangle_1 ) {
239 vtkCellType = 5;
240 } else if ( elemGT == EGT_triangle_2 ) {
241 vtkCellType = 22;
242 } else if ( elemGT == EGT_tetra_1 ) {
243 vtkCellType = 10;
244 } else if ( elemGT == EGT_tetra_2 ) {
245 vtkCellType = 24;
246 } else if ( elemGT == EGT_quad_1 ) {
247 vtkCellType = 9;
248 } else if ( elemGT == EGT_quad_2 ) {
249 vtkCellType = 23;
250 } else if ( elemGT == EGT_hexa_1 ) {
251 vtkCellType = 12;
252 } else if ( elemGT == EGT_hexa_2 ) {
253 vtkCellType = 25;
254 } else if ( elemGT == EGT_wedge_1 ) {
255 vtkCellType = 13;
256 } else if ( elemGT == EGT_wedge_2 ) {
257 vtkCellType = 26;
258 } else {
259 OOFEM_ERROR("unsupported element gemetry type");
260 }
261
262 return vtkCellType;
263}
264
265int
266VTKExportModule :: giveNumberOfNodesPerCell(int cellType)
267{
268 switch ( cellType ) {
269 case 1:
270 return 1;
271
272 case 3:
273 return 2;
274
275 case 5:
276 case 21:
277 return 3;
278
279 case 9:
280 case 10:
281 return 4;
282
283 case 14:
284 return 5;
285
286 case 13:
287 case 22:
288 return 6;
289
290 case 12:
291 case 23:
292 return 8;
293
294 case 24:
295 return 10;
296
297 case 25:
298 return 20;
299
300 default:
301 OOFEM_ERROR("unsupported cell type ID");
302 }
303
304 // return 0; // to make compiler happy
305}
306
307
308void
309VTKExportModule :: giveElementCell(IntArray &answer, Element *elem, int cell)
310{
312 int i, nelemNodes;
313
314 if ( ( elemGT == EGT_point ) ||
315 ( elemGT == EGT_line_1 ) || ( elemGT == EGT_line_2 ) ||
316 ( elemGT == EGT_triangle_1 ) || ( elemGT == EGT_triangle_2 ) ||
317 ( elemGT == EGT_tetra_1 ) || ( elemGT == EGT_tetra_2 ) ||
318 ( elemGT == EGT_quad_1 ) || ( elemGT == EGT_quad_2 ) ||
319 ( elemGT == EGT_hexa_1 ) || ( elemGT == EGT_wedge_1 ) ) {
320 nelemNodes = elem->giveNumberOfNodes();
321 answer.resize(nelemNodes);
322 for ( i = 1; i <= nelemNodes; i++ ) {
323 answer.at(i) = elem->giveNode(i)->giveNumber();
324 }
325 } else if ( elemGT == EGT_hexa_2 ) {
326 int HexaQuadNodeMapping [] = {
327 5, 8, 7, 6, 1, 4, 3, 2, 16, 15, 14, 13, 12, 11, 10, 9, 17, 20, 19, 18
328 };
329 nelemNodes = elem->giveNumberOfNodes();
330 answer.resize(nelemNodes);
331 for ( i = 1; i <= nelemNodes; i++ ) {
332 answer.at(i) = elem->giveNode(HexaQuadNodeMapping [ i - 1 ])->giveNumber();
333 }
334 } else if ( elemGT == EGT_wedge_2 ) {
335 int WedgeQuadNodeMapping [] = {
336 4, 6, 5, 1, 3, 2, 12, 11, 10, 9, 8, 7, 13, 15, 14
337 };
338 nelemNodes = elem->giveNumberOfNodes();
339 answer.resize(nelemNodes);
340 for ( i = 1; i <= nelemNodes; i++ ) {
341 answer.at(i) = elem->giveNode(WedgeQuadNodeMapping [ i - 1 ])->giveNumber();
342 }
343 } else {
344 OOFEM_ERROR("unsupported element geometry type");
345 }
346}
347
348
349int
350VTKExportModule :: giveNumberOfElementCells(Element *elem)
351{
353
354 if ( ( elemGT == EGT_point ) ||
355 ( elemGT == EGT_line_1 ) || ( elemGT == EGT_line_2 ) ||
356 ( elemGT == EGT_triangle_1 ) || ( elemGT == EGT_triangle_2 ) ||
357 ( elemGT == EGT_tetra_1 ) || ( elemGT == EGT_tetra_2 ) ||
358 ( elemGT == EGT_quad_1 ) || ( elemGT == EGT_quad_2 ) ||
359 ( elemGT == EGT_hexa_1 ) || ( elemGT == EGT_hexa_2 ) ||
360 ( elemGT == EGT_wedge_1 ) || ( elemGT == EGT_wedge_2 ) ) {
361 return 1;
362 } else {
363 OOFEM_ERROR("unsupported element geometry type");
364 }
365
366 // return 0;
367}
368
369//keyword "vars" in OOFEM input file
370void
371VTKExportModule :: exportIntVars(FILE *stream, TimeStep *tStep)
372{
373 // should be performed over regions
374
375 int i, n = internalVarsToExport.giveSize();
376 //int nnodes;
377 //Domain *d = emodel->giveDomain(1);
379
380 if ( n == 0 ) {
381 return;
382 }
383
384 for ( i = 1; i <= n; i++ ) {
387 this->exportIntVarAs(type, iType, stream, tStep);
388 }
389}
390
391void
392VTKExportModule :: exportCellVars(FILE *stream, int elemToProcess, TimeStep *tStep)
393{
394 int pos;
396 Domain *d = emodel->giveDomain(1);
397 FloatMatrix mtrx(3, 3);
398 FloatArray temp, vec;
399 double gptot;
400
401 fprintf(stream, "\nCELL_DATA %d\n", elemToProcess);
402 for ( int i = 1; i <= cellVarsToExport.giveSize(); i++ ) {
403 type = ( InternalStateType ) cellVarsToExport.at(i);
404 switch ( type ) {
405 case IST_MaterialNumber:
406 case IST_ElementNumber:
407 fprintf( stream, "SCALARS %s int\nLOOKUP_TABLE default\n", __InternalStateTypeToString(type) );
408 for ( auto &elem : d->giveElements() ) {
409 if ( elem->giveParallelMode() != Element_local ) {
410 continue;
411 }
412
413 if ( type == IST_MaterialNumber || type == IST_CrossSectionNumber ) {
414 OOFEM_WARNING("Material numbers are deprecated, outputing cross section number instead...");
415 fprintf( stream, "%d\n", elem->giveCrossSection()->giveNumber() );
416 } else if ( type == IST_ElementNumber ) {
417 fprintf( stream, "%d\n", elem->giveNumber() );
418 } else {
419 OOFEM_ERROR("Unsupported Cell variable %s\n", __InternalStateTypeToString(type));
420 }
421 }
422
423 break;
424
425 case IST_MaterialOrientation_x:
426 case IST_MaterialOrientation_y:
427 case IST_MaterialOrientation_z:
428 if ( type == IST_MaterialOrientation_x ) {
429 pos = 1;
430 } else if ( type == IST_MaterialOrientation_y ) {
431 pos = 2;
432 } else {
433 pos = 3;
434 }
435
436 fprintf( stream, "VECTORS %s double\n", __InternalStateTypeToString(type) );
437 for ( auto &elem : d->giveElements() ) {
438 if ( !elem->giveLocalCoordinateSystem(mtrx) ) {
439 mtrx.resize(3, 3);
440 mtrx.zero();
441 }
442
443 fprintf( stream, "%f %f %f\n", mtrx.at(1, pos), mtrx.at(2, pos), mtrx.at(3, pos) );
444 }
445
446 break;
447 default:
449 switch ( vt ) {
450 case ISVT_TENSOR_S3: // These could be written as tensors as well, if one wants to.
451 case ISVT_TENSOR_S3E:
452 case ISVT_VECTOR:
453 case ISVT_SCALAR:
454 if ( vt == ISVT_SCALAR ) {
455 fprintf( stream, "SCALARS %s double\nLOOKUP_TABLE default\n", __InternalStateTypeToString(type) );
456 } else {
457 fprintf( stream, "VECTORS %s double\nLOOKUP_TABLE default\n", __InternalStateTypeToString(type) );
458 }
459 for ( auto &elem : d->giveElements() ) {
460 if ( elem->giveParallelMode() != Element_local ) {
461 continue;
462 }
463
464 gptot = 0;
465 vec.clear();
466 for ( GaussPoint *gp: *elem->giveDefaultIntegrationRulePtr() ) {
467 elem->giveIPValue(temp, gp, type, tStep);
468 gptot += gp->giveWeight();
469 vec.add(gp->giveWeight(), temp);
470 }
471 vec.times(1 / gptot);
472 for ( int j = 1; j <= vec.giveSize(); ++j ) {
473 fprintf( stream, "%e ", vec.at(j) );
474 }
475 fprintf(stream, "\n");
476 }
477 break;
478#if 0 // Hardly even worth the effort...
479 case ISVT_TENSOR_G:
480 for ( int indx = 1; indx < 9; ++indx ) {
481 fprintf(stream, "SCALARS %s_%d double 1\n", __InternalStateTypeToString(valID), indx);
482
483 for ( ielem = 1; ielem <= nelem; ielem++ ) {
484 elem = d->giveElement(ielem);
485 elem->giveIPValue(vec, elem->giveDefaultIntegrationRulePtr()->getIntegrationPoint(1), type, tStep);
486 fprintf( stream, "%e ", vec.at(indx) );
487 }
488 }
489#endif
490 default:
491 OOFEM_ERROR("Quantity %s not handled yet.", __InternalStateTypeToString(type));
492 }
493 }
494
495 fprintf(stream, "\n\n");
496 }
497}
498
499
500int
501VTKExportModule :: giveTotalRBRNumberOfNodes(Domain *d)
502//
503// Returns total number of nodes with region multiplicity
504// (RBR = Region by Region)
505// The nodes on region boundaries are for each region.
506// We need unique node numbers for each region, to prevent
507// vtk from smoothing the nodal values at region boundaries.
508//
509{
510 int rbrnodes = 0, nnodes = d->giveNumberOfDofManagers();
511 std :: vector< char > map(nnodes);
512 //char map[nnodes];
513 int elemnodes;
514
515 for ( int j = 0; j < nnodes; j++ ) {
516 map [ j ] = 0;
517 }
518
519 for ( auto &elem : d->giveElements() ) {
520 if ( elem->giveParallelMode() != Element_local ) {
521 continue;
522 }
523
524 elemnodes = elem->giveNumberOfNodes();
525 for ( int ielemnode = 1; ielemnode <= elemnodes; ielemnode++ ) {
526 map [ elem->giveNode(ielemnode)->giveNumber() - 1 ] = 1;
527 }
528 }
529
530 for ( int j = 0; j < nnodes; j++ ) {
531 rbrnodes += map [ j ];
532 }
533
534 return rbrnodes;
535}
536
537
538int
539VTKExportModule :: initRegionNodeNumbering(IntArray &regionNodalNumbers, int &regionDofMans,
540 int offset, Domain *domain, int reg, int mode)
541{
542 // if mode == 0 then regionNodalNumbers is array with mapping from global numbering to local region numbering.
543 // The i-th value contains the corresponding local region number (or zero, if global number is not in region).
544
545 // if mode == 1 then regionNodalNumbers is array with mapping from local to global numbering.
546 // The i-th value contains the corresponding global node number.
547
548
549 int nnodes = domain->giveNumberOfDofManagers();
550 int elemNodes;
551 int currOffset = offset + 1;
552
553 regionNodalNumbers.resize(nnodes);
554 regionNodalNumbers.zero();
555 regionDofMans = 0;
556
557 for ( auto &element : domain->giveElements() ) {
558
559 if ( element->giveParallelMode() != Element_local ) {
560 continue;
561 }
562
563 elemNodes = element->giveNumberOfNodes();
564 // elemSides = element->giveNumberOfSides();
565
566 // determine local region node numbering
567 for ( int elementNode = 1; elementNode <= elemNodes; elementNode++ ) {
568 int node = element->giveNode(elementNode)->giveNumber();
569 if ( regionNodalNumbers.at(node) == 0 ) { // assign new number
570 /* mark for assignement. This is done later, as it allows to preserve
571 * natural node numbering.
572 */
573 regionNodalNumbers.at(node) = 1;
574 regionDofMans++;
575 }
576 }
577 }
578
579 if ( mode == 1 ) {
580 IntArray answer(nnodes);
581 for ( int i = 1; i <= nnodes; i++ ) {
582 if ( regionNodalNumbers.at(i) ) {
583 regionNodalNumbers.at(i) = currOffset++;
584 answer.at( regionNodalNumbers.at(i) ) = i;
585 }
586 }
587
588 regionNodalNumbers = answer;
589 } else {
590 for ( int i = 1; i <= nnodes; i++ ) {
591 if ( regionNodalNumbers.at(i) ) {
592 regionNodalNumbers.at(i) = currOffset++;
593 }
594 }
595 }
596
597 return 1;
598}
599
600//keyword "vars" in OOFEM input file
601void
602VTKExportModule :: exportIntVarAs(InternalStateType valID, InternalStateValueType type, FILE *stream, TimeStep *tStep)
603{
604 Domain *d = emodel->giveDomain(1);
605 int ireg;
606 int nnodes = d->giveNumberOfDofManagers(), inode;
607 int j, jsize;
608 FloatArray iVal(3);
609 FloatMatrix t(3, 3);
610 const FloatArray *val = NULL;
611
612
613 // create a new set containing all elements
614 Set elemSet(0, d);
615 elemSet.addAllElements();
616
617 this->giveSmoother();
618
619 int nindx = giveInternalStateTypeSize(type);
620
621 for ( int indx = 1; indx <= nindx; indx++ ) {
622 // print header
623 if ( type == ISVT_SCALAR ) {
624 fprintf( stream, "SCALARS %s double 1\n", __InternalStateTypeToString(valID) );
625 } else if ( type == ISVT_VECTOR ) {
626 fprintf( stream, "VECTORS %s double\n", __InternalStateTypeToString(valID) );
627 } else if ( ( type == ISVT_TENSOR_S3 ) || ( type == ISVT_TENSOR_S3E ) ) {
628 fprintf( stream, "TENSORS %s double\n", __InternalStateTypeToString(valID) );
629 } else if ( type == ISVT_TENSOR_G ) {
630 fprintf(stream, "SCALARS %s_%d double 1\n", __InternalStateTypeToString(valID), indx);
631 } else {
632 fprintf( stderr, "VTKExportModule :: exportIntVarAs: unsupported variable type %s\n", __InternalStateTypeToString(valID) );
633 }
634
635 if ( ( type == ISVT_SCALAR ) || ( type == ISVT_TENSOR_G ) ) {
636 fprintf(stream, "LOOKUP_TABLE default\n");
637 }
638
639 if ( !( ( valID == IST_DisplacementVector ) || ( valID == IST_MaterialInterfaceVal ) ) ) {
640 this->smoother->recoverValues(elemSet, valID, tStep);
641 }
642
643 IntArray regionNodalNumbers(nnodes);
644 int regionDofMans = 0, offset = 0;
645 ireg = -1;
646 int defaultSize = 0;
647
648 this->initRegionNodeNumbering(regionNodalNumbers, regionDofMans, offset, d, ireg, 1);
649 if ( !( ( valID == IST_DisplacementVector ) || ( valID == IST_MaterialInterfaceVal ) ) ) {
650 // assemble local->global map
651 defaultSize = giveInternalStateTypeSize(type);
652 } else {
653 regionDofMans = nnodes;
654 }
655
656 for ( inode = 1; inode <= regionDofMans; inode++ ) {
657 if ( valID == IST_DisplacementVector ) {
658 iVal.resize(3);
659 val = & iVal;
660 for ( j = 1; j <= 3; j++ ) {
661 iVal.at(j) = d->giveNode( regionNodalNumbers.at(inode) )->giveUpdatedCoordinate(j, tStep, 1.0) -
662 d->giveNode( regionNodalNumbers.at(inode) )->giveCoordinate(j);
663 }
664 } else if ( valID == IST_MaterialInterfaceVal ) {
665 MaterialInterface *mi = emodel->giveMaterialInterface(1);
666 if ( mi ) {
667 iVal.resize(1);
668 val = & iVal;
669 iVal.at(1) = mi->giveNodalScalarRepresentation( regionNodalNumbers.at(inode) );
670 }
671 } else {
672 this->smoother->giveNodalVector( val, regionNodalNumbers.at(inode) );
673 }
674
675 if ( val == NULL ) {
676 iVal.resize(defaultSize);
677 iVal.zero();
678 val = & iVal;
679 //OOFEM_ERROR("internal error: invalid dofman data");
680 }
681 }
682
683
684 if ( type == ISVT_SCALAR ) {
685 if ( val->giveSize() ) {
686 fprintf( stream, "%e ", val->at(1) );
687 } else {
688 fprintf(stream, "%e ", 0.0);
689 }
690 } else if ( type == ISVT_VECTOR ) {
691 jsize = min( 3, (int)val->giveSize() );
692 for ( j = 1; j <= jsize; j++ ) {
693 fprintf( stream, "%e ", val->at(j) );
694 }
695
696 for ( j = jsize + 1; j <= 3; j++ ) {
697 fprintf(stream, "0.0 ");
698 }
699 } else if ( type == ISVT_TENSOR_S3 ) {
700 t.zero();
701 for ( int ii = 1; ii <= 6; ii++ ) {
702 if ( ii == 1 ) {
703 t.at(1, 1) = val->at(ii);
704 } else if ( ii == 2 ) {
705 t.at(2, 2) = val->at(ii);
706 } else if ( ii == 3 ) {
707 t.at(3, 3) = val->at(ii);
708 } else if ( ii == 4 ) {
709 t.at(2, 3) = val->at(ii);
710 t.at(3, 2) = val->at(ii);
711 } else if ( ii == 5 ) {
712 t.at(1, 3) = val->at(ii);
713 t.at(3, 1) = val->at(ii);
714 } else if ( ii == 6 ) {
715 t.at(1, 2) = val->at(ii);
716 t.at(2, 1) = val->at(ii);
717 }
718 }
719
720 for ( int ii = 1; ii <= 3; ii++ ) {
721 for ( int jj = 1; jj <= 3; jj++ ) {
722 fprintf( stream, "%e ", t.at(ii, jj) );
723 }
724
725 fprintf(stream, "\n");
726 }
727 } else if ( type == ISVT_TENSOR_G ) { // export general tensor values as scalars
728 fprintf( stream, "%e ", val->at(indx) );
729 }
730
731 fprintf(stream, "\n");
732 }
733}
734
735
737VTKExportModule :: giveSmoother()
738{
739 Domain *d = emodel->giveDomain(1);
740
741 if ( !this->smoother ) {
742 this->smoother = classFactory.createNodalRecoveryModel(this->stype, d);
743 }
744
745 return this->smoother.get();
746}
747
748
749//keyword "primvars" in OOFEM input file
750void
751VTKExportModule :: exportPrimaryVars(FILE *stream, TimeStep *tStep)
752{
753 for ( auto &primvar : primaryVarsToExport ) {
754 UnknownType type = ( UnknownType ) primvar;
755 this->exportPrimVarAs(type, stream, tStep);
756 }
757}
758
759//keyword "primvars" in OOFEM input file
760void
761VTKExportModule :: exportPrimVarAs(UnknownType valID, FILE *stream, TimeStep *tStep)
762{
763 Domain *d = emodel->giveDomain(1);
764 int ireg;
765 int nnodes = d->giveNumberOfDofManagers();
766 int j, jsize;
767 FloatArray iVal, iValLCS;
768 IntArray dofIDMask;
770 int nScalarComp = 1;
771
772 if ( ( valID == DisplacementVector ) || ( valID == EigenVector ) || ( valID == VelocityVector ) ) {
773 type = ISVT_VECTOR;
774 } else if ( ( valID == FluxVector ) || ( valID == PressureVector ) || ( valID == Temperature ) || ( valID == Humidity )) {
775 type = ISVT_SCALAR;
776 //nScalarComp = d->giveNumberOfDefaultNodeDofs();
777 } else {
778 OOFEM_ERROR("unsupported UnknownType (%s)", __UnknownTypeToString(valID) );
779 }
780
781 // print header
782 if ( type == ISVT_SCALAR ) {
783 fprintf(stream, "SCALARS %s double %d\n", __UnknownTypeToString(valID), nScalarComp);
784 } else if ( type == ISVT_VECTOR ) {
785 fprintf( stream, "VECTORS %s double\n", __UnknownTypeToString(valID) );
786 } else {
787 fprintf(stderr, "exportPrimVarAs: unsupported variable type\n");
788 }
789
790 if ( type == ISVT_SCALAR ) {
791 fprintf(stream, "LOOKUP_TABLE default\n");
792 }
793
794 DofIDItem id;
795
796 IntArray regionNodalNumbers(nnodes);
797 int regionDofMans = 0, offset = 0;
798 ireg = -1;
799
800
801 // assemble local->global map
802 this->initRegionNodeNumbering(regionNodalNumbers, regionDofMans, offset, d, ireg, 1);
803
804 // Used for special nodes/elements with mixed number of dofs.
805 this->giveSmoother();
806
807 for ( int inode = 1; inode <= regionDofMans; inode++ ) {
808 DofManager *dman = d->giveNode( regionNodalNumbers.at(inode) );
809
810 if ( ( valID == DisplacementVector ) || ( valID == EigenVector ) || ( valID == VelocityVector ) ) {
811 iVal.resize(3);
812 iVal.zero();
813
814 for ( Dof *dof: *dman ) {
815 id = dof->giveDofID();
816 if ( ( id == V_u ) || ( id == D_u ) ) {
817 iVal.at(1) = dof->giveUnknown(VM_Total, tStep);
818 } else if ( ( id == V_v ) || ( id == D_v ) ) {
819 iVal.at(2) = dof->giveUnknown(VM_Total, tStep);
820 } else if ( ( id == V_w ) || ( id == D_w ) ) {
821 iVal.at(3) = dof->giveUnknown(VM_Total, tStep);
822 }
823 }
824 } else if ( (valID == FluxVector) || (valID == Humidity) ) {
825 iVal.resize(1);
826
827 for ( Dof *dof: *dman ) {
828 id = dof->giveDofID();
829 if ( id == C_1 ) {
830 iVal.at(1) = dof->giveUnknown(VM_Total, tStep);
831 }
832 }
833 } else if ( valID == Temperature ) {
834 iVal.resize(1);
835
836 for ( Dof *dof: *dman ) {
837 id = dof->giveDofID();
838 if ( id == T_f ) {
839 iVal.at(1) = dof->giveUnknown(VM_Total, tStep);
840 }
841 }
842 } else if ( valID == PressureVector ) {
843 dofIDMask = {P_f};
844 this->getDofManPrimaryVariable(iVal, dman, dofIDMask, VM_Total, tStep, IST_Pressure);
845 } else {
846 OOFEM_ERROR("unsupported unknownType (%s)", __UnknownTypeToString(valID) );
847 }
848
849 if ( type == ISVT_SCALAR ) {
850 if ( iVal.giveSize() ) {
851 for ( j = 1; j <= nScalarComp; j++ ) {
852 fprintf( stream, "%e ", iVal.at(j) );
853 }
854 } else {
855 fprintf(stream, "%e ", 0.0);
856 }
857 } else if ( type == ISVT_VECTOR ) {
858 jsize = min( 3, (int)iVal.giveSize() );
859 //rotate back from nodal CS to global CS if applies
860 if ( d->giveNode( dman->giveNumber() )->hasLocalCS() ) {
861 iVal.resize(3);
862 iValLCS = iVal;
863 iVal.beTProductOf(* d->giveNode( dman->giveNumber() )->giveLocalCoordinateTriplet(), iValLCS);
864 }
865
866 for ( j = 1; j <= jsize; j++ ) {
867 fprintf( stream, "%e ", iVal.at(j) );
868 }
869
870 for ( j = jsize + 1; j <= 3; j++ ) {
871 fprintf(stream, "0.0 ");
872 }
873 }
874
875 fprintf(stream, "\n");
876 }
877
878 fprintf(stream, "\n");
879}
880
881void
882VTKExportModule :: getDofManPrimaryVariable(FloatArray &answer, DofManager *dman, IntArray &dofIDMask,
883 ValueModeType mode, TimeStep *tStep, InternalStateType iType)
884{
885 int size = dofIDMask.giveSize();
886 const FloatArray *recoveredVal;
887 answer.resize(size);
888 // all values zero by default
889 answer.zero();
890
891
892
893 for ( int j = 1; j <= size; j++ ) {
894 if ( dman->hasDofID( ( DofIDItem ) dofIDMask.at(j) ) ) {
895 // primary variable available directly in dof manager
896 answer.at(j) = dman->giveDofWithID( dofIDMask.at(j) )->giveUnknown(VM_Total, tStep);
897 } else if ( iType != IST_Undefined ) {
898 // ok primary variable not directly available
899 // but equivalent InternalStateType provided
900 // in this case use smoother to recover nodal value
901
902 // create a new set containing all elements
903 Set elemSet( 0, dman->giveDomain() );
904 elemSet.addAllElements();
905
906 this->giveSmoother()->recoverValues(elemSet, iType, tStep);
907 this->giveSmoother()->giveNodalVector( recoveredVal, dman->giveNumber() );
908 // here we have a lack of information about how to convert recoveredVal to response
909 // if the size is compatible we accept it, otherwise an error is thrown
910 if ( size == recoveredVal->giveSize() ) {
911 answer.at(j) = recoveredVal->at(j);
912 } else {
913 OOFEM_WARNING("recovered variable size mismatch for %d", iType);
914 answer.at(j) = 0.0;
915 }
916 }
917 }
918}
919} // end namespace oofem
#define REGISTER_ExportModule(class)
bool hasDofID(DofIDItem id) const
Definition dofmanager.C:174
double giveCoordinate(int i) const
Definition dofmanager.h:383
const FloatArray & giveCoordinates() const
Definition dofmanager.h:390
Dof * giveDofWithID(int dofID) const
Definition dofmanager.C:127
virtual double giveUnknown(ValueModeType mode, TimeStep *tStep)=0
int giveNumberOfDofManagers() const
Returns number of dof managers in domain.
Definition domain.h:461
Element * giveElement(int n)
Definition domain.C:165
Node * giveNode(int n)
Definition domain.h:398
std ::vector< std ::unique_ptr< Element > > & giveElements()
Definition domain.h:294
Node * giveNode(int i) const
Definition element.h:629
virtual int giveNumberOfNodes() const
Definition element.h:703
virtual int giveIPValue(FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
Definition element.C:1298
virtual Element_Geometry_Type giveGeometryType() const =0
std::string giveOutputBaseFileName(TimeStep *tStep)
EngngModel * emodel
Problem pointer.
bool testTimeStepOutput(TimeStep *tStep)
Domain * giveDomain() const
Definition femcmpnn.h:97
int giveNumber() const
Definition femcmpnn.h:104
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 zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beTProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:721
void add(const FloatArray &src)
Definition floatarray.C:218
void times(double s)
Definition floatarray.C:834
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void zero()
Zeroes all coefficient of receiver.
double at(std::size_t i, std::size_t j) const
void resize(int n)
Definition intarray.C:73
void zero()
Sets all component to zero.
Definition intarray.C:52
int & at(std::size_t i)
Definition intarray.h:104
int giveSize() const
Definition intarray.h:211
virtual double giveNodalScalarRepresentation(int)=0
FloatMatrix * giveLocalCoordinateTriplet()
Definition node.h:149
bool hasLocalCS()
Returns nonzero if node has prescribed local coordinate system.
Definition node.h:141
virtual double giveUpdatedCoordinate(int ic, TimeStep *tStep, double scale=1.)
Definition node.C:234
void addAllElements()
Definition set.C:239
double giveTargetTime()
Returns target time.
Definition timestep.h:164
NodalRecoveryModel::NodalRecoveryModelType stype
Smoother type.
int giveTotalRBRNumberOfNodes(Domain *d)
FILE * giveOutputStream(TimeStep *tStep)
Returns the output stream for given solution step.
NodalRecoveryModel * giveSmoother()
Returns the internal smoother.
int giveCellType(Element *tStep)
std::unique_ptr< NodalRecoveryModel > smoother
Smoother.
void exportPrimVarAs(UnknownType valID, FILE *stream, TimeStep *tStep)
void exportIntVarAs(InternalStateType valID, InternalStateValueType type, FILE *stream, TimeStep *tStep)
void giveElementCell(IntArray &answer, Element *elem, int cell)
void getDofManPrimaryVariable(FloatArray &answer, DofManager *dman, IntArray &dofIDMask, ValueModeType mode, TimeStep *tStep, InternalStateType iType)
IntArray primaryVarsToExport
List of primary unknowns to export.
IntArray internalVarsToExport
List of InternalStateType values, identifying the selected vars for export.
IntArray cellVarsToExport
List of cell data to export.
void exportCellVars(FILE *stream, int elemToProcess, TimeStep *tStep)
void exportIntVars(FILE *stream, TimeStep *tStep)
void exportPrimaryVars(FILE *stream, TimeStep *tStep)
int giveNumberOfNodesPerCell(int cellType)
int initRegionNodeNumbering(IntArray &regionNodalNumbers, int &regionDofMans, int offset, Domain *domain, int reg, int mode)
int giveNumberOfElementCells(Element *)
#define OOFEM_WARNING(...)
Definition error.h:80
#define OOFEM_ERROR(...)
Definition error.h:79
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#define OOFEM_LOG_DEBUG(...)
Definition logger.h:144
const char * __InternalStateTypeToString(InternalStateType _value)
Definition cltypes.C:309
FloatArrayF< N > min(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
@ Element_local
Element is local, there are no contributions from other domains to this element.
Definition element.h:88
int giveInternalStateTypeSize(InternalStateValueType valType)
Definition cltypes.C:229
InternalStateValueType
Determines the type of internal variable.
@ ISVT_TENSOR_S3E
symmetric 3x3 tensor, packed with off diagonal components multiplied by 2 (engineering strain vector,...
@ ISVT_SCALAR
Scalar.
@ ISVT_UNDEFINED
Undefined.
@ ISVT_TENSOR_S3
Symmetric 3x3 tensor.
@ ISVT_VECTOR
Vector.
@ ISVT_TENSOR_G
General tensor.
const char * __UnknownTypeToString(UnknownType _value)
Definition cltypes.C:313
ClassFactory & classFactory
InternalStateValueType giveInternalStateValueType(InternalStateType type)
Definition cltypes.C:75
#define _IFT_VTKExportModule_stype
#define _IFT_VTKExportModule_cellvars
#define _IFT_VTKExportModule_vars
#define _IFT_VTKExportModule_primvars

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