OOFEM 3.0
Loading...
Searching...
No Matches
vtkmemoryexportmodule.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
36#include "element.h"
37#include "timestep.h"
38#include "engngm.h"
39#include "node.h"
42#include "classfactory.h"
43
44namespace oofem {
49
50void
61
62void
63VTKMemoryExportModule::doOutput(TimeStep *tStep, bool forcedOutput)
64{
65 if ( !( testTimeStepOutput(tStep) || forcedOutput ) ) {
66 return;
67 }
68
69 int nPiecesToExport = this->giveNumberOfRegions(); //old name: region, meaning: sets
70 ZZNodalRecoveryModel smoother(emodel->giveDomain(1));
71 NodalAveragingRecoveryModel primVarSmoother(emodel->giveDomain(1));
72
73 this->vtkPieces.resize(nPiecesToExport);
74
75 // loop over regular pieces only (no support for composite elements at present)
76 for ( int pieceNum = 1; pieceNum <= nPiecesToExport; pieceNum++ ) {
77 ExportRegion& p = this->vtkPieces[pieceNum-1];
78 p.clear();
79 // Fills a data struct (VTKPiece) with all the necessary data.
80 Set* region = this->giveRegionSet(pieceNum);
81 this->setupVTKPiece(p, tStep, *region);
82 // Export primary, internal and XFEM variables as nodal quantities
83 this->exportPrimaryVars(p, *region, primaryVarsToExport, primVarSmoother, tStep);
84 this->exportIntVars(p, *region, internalVarsToExport, smoother, tStep);
85 this->exportExternalForces(p, *region, externalForcesToExport, tStep);
86 this->exportCellVars(p, *region, cellVarsToExport, tStep);
87
88 }
89}
90
91std::vector< ExportRegion>&
95
96
97} // end namespace oofem
#define REGISTER_ExportModule(class)
virtual void initializeFrom(InputRecord &ir)
Initializes receiver according to object description stored in input record.
Set * giveRegionSet(int i)
Returns element set.
int giveNumberOfRegions()
Returns number of regions (aka regionSets).
EngngModel * emodel
Problem pointer.
bool testTimeStepOutput(TimeStep *tStep)
Stores all neccessary data (of a region) in a VTKPiece so it can be exported later.
void exportCellVars(ExportRegion &piece, Set &region, IntArray &cellVarsToExport, TimeStep *tStep)
Exports cell variables (typically internal variables).
virtual void exportIntVars(ExportRegion &piece, Set &region, IntArray &internalVarsToExport, NodalRecoveryModel &smoother, TimeStep *tStep)
void exportExternalForces(ExportRegion &piece, int region, TimeStep *tStep)
virtual void setupVTKPiece(ExportRegion &vtkPiece, TimeStep *tStep, Set &region)
virtual void exportPrimaryVars(ExportRegion &piece, Set &region, IntArray &primaryVarsToExport, NodalRecoveryModel &smoother, TimeStep *tStep)
std::vector< ExportRegion > & getExportRegions()
void doOutput(TimeStep *tStep, bool forcedOutput=false) override
VTKMemoryExportModule(int n, EngngModel *e)
Constructor. Creates empty Output Manager. By default all components are selected.
IntArray internalVarsToExport
List of InternalStateType values, identifying the selected vars for export.
IntArray cellVarsToExport
List of cell data to export.
virtual ~VTKMemoryExportModule()
Destructor.
void initializeFrom(InputRecord &ir) override
Initializes receiver according to object description stored in input record.
IntArray externalForcesToExport
List of primary unknowns to export.
std::vector< ExportRegion > vtkPieces
IntArray primaryVarsToExport
List of primary unknowns to export.
IntArray ipInternalVarsToExport
List of internal variables to export directly in Integration Points (no smoothing to nodes).
#define IR_GIVE_OPTIONAL_FIELD(__ir, __value, __id)
Definition inputrecord.h:75
#define _IFT_VTKXMLExportModule_ipvars
#define _IFT_VTKXMLExportModule_vars
#define _IFT_VTKXMLExportModule_externalForces
#define _IFT_VTKXMLExportModule_primvars
#define _IFT_VTKXMLExportModule_cellvars

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