35 #ifndef vtkxmlexportmodule_h 36 #define vtkxmlexportmodule_h 48 #include <vtkUnstructuredGrid.h> 49 #include <vtkSmartPointer.h> 57 #define _IFT_VTKXMLExportModule_Name "vtkxml" 58 #define _IFT_VTKXMLExportModule_cellvars "cellvars" 59 #define _IFT_VTKXMLExportModule_vars "vars" 60 #define _IFT_VTKXMLExportModule_primvars "primvars" 61 #define _IFT_VTKXMLExportModule_externalForces "externalforces" 62 #define _IFT_VTKXMLExportModule_ipvars "ipvars" 63 #define _IFT_VTKXMLExportModule_stype "stype" 64 #define _IFT_VTKXMLExportModule_particleexportflag "particleexportflag" 129 std :: vector< std :: vector< FloatArray > >
nodeVars;
133 std :: vector< std :: vector< FloatArray > >
elVars;
184 virtual void doOutput(
TimeStep *tStep,
bool forcedOutput =
false);
185 virtual void initialize();
186 virtual void terminate();
191 void exportPointDataHeader(FILE *fileStream,
TimeStep *tStep);
192 void giveDataHeaders(std :: string &pointHeader, std :: string &cellHeader);
200 vtkSmartPointer< vtkUnstructuredGrid >fileStream;
201 vtkSmartPointer< vtkPoints >nodes;
202 vtkSmartPointer< vtkIdList >elemNodeArray;
204 vtkSmartPointer< vtkDoubleArray >intVarArray;
205 vtkSmartPointer< vtkDoubleArray >primVarArray;
226 std :: string giveOutputFileName(
TimeStep *tStep);
229 FILE *giveOutputStream(
TimeStep *tStep);
238 int giveNumberOfNodesPerCell(
int cellType);
272 virtual void setupVTKPiece(
VTKPiece &vtkPiece,
TimeStep *tStep,
int region);
273 void writeIntVars(
VTKPiece &vtkPiece);
274 void writeXFEMVars(
VTKPiece &vtkPiece);
275 void writePrimaryVars(
VTKPiece &vtkPiece);
276 void writeCellVars(
VTKPiece &vtkPiece);
277 void writeExternalForces(
VTKPiece &vtkPiece);
323 void writeVTKCollection();
326 void writeGPVTKCollection();
329 void writeVTKPointData(
const char *name, vtkSmartPointer< vtkDoubleArray >varArray);
331 void writeVTKPointData(
FloatArray &valueArray);
335 void writeVTKCellData(
const char *name, vtkSmartPointer< vtkDoubleArray >varArray);
337 void writeVTKCellData(
FloatArray &valueArray);
342 bool isElementComposite(
Element *elem);
344 void exportCompositeElement(std::vector< VTKPiece > &vtkPieces,
Element *el,
TimeStep *tStep);
360 virtual const char *
giveClassName()
const {
return "VTKXMLExportModuleElementInterface"; }
365 #endif // vtkxmlexportmodule_h void setInternalVarInNode(int varNum, int nodeNum, FloatArray valueArray)
void setInternalXFEMVarInNode(int varNum, int eiNum, int nodeNum, FloatArray valueArray)
InternalStateType
Type representing the physical meaning of element or constitutive model internal variable.
int giveCellType(int cellNum)
void setCellVar(int varNum, int cellNum, FloatArray valueArray)
Abstract class representing entity, which is included in the FE model using one (or more) global func...
IntArray ipInternalVarsToExport
List of internal variables to export directly in Integration Points (no smoothing to nodes) ...
std::vector< std::vector< FloatArray > > nodeVars
NodalRecoveryModel * smoother
Smoother.
void setNumberOfNodes(int numNodes)
void setNumberOfInternalVarsToExport(int numVars, int numNodes)
Elements with geometry defined as EGT_Composite are exported using individual pieces.
Abstract base class for all finite elements.
InternalStateValueType
Determines the type of internal variable.
Base class for dof managers.
std::vector< IntArray > connectivity
int giveCellOffset(int cellNum)
Represents export output module - a base class for all output modules.
Class implementing an array of integers.
int & at(int i)
Coefficient access function.
void setNumberOfCells(int numCells)
FloatArray & giveInternalVarInNode(int varNum, int nodeNum)
std::vector< FloatArray > nodeCoords
Abstract base class representing integration rule.
void setNumberOfLoadsToExport(int numVars, int numNodes)
std::vector< std::vector< FloatArray > > nodeVarsFromIS
std::list< std::string > pvdBuffer
Buffer for earlier time steps exported to *.pvd file.
void setNodeCoords(int nodeNum, FloatArray &coords)
FloatArray & giveLoadInNode(int varNum, int nodeNum)
void setLoadInNode(int varNum, int nodeNum, FloatArray valueArray)
IntArray cellVarsToExport
List of cell data to export.
std::vector< std::vector< std::vector< FloatArray > > > nodeVarsFromXFEMIS
void setCellType(int cellNum, int type)
void setNumberOfInternalXFEMVarsToExport(int numVars, int numEnrichmentItems, int numNodes)
IntArray internalVarsToExport
List of InternalStateType values, identifying the selected vars for export.
void setNumberOfPrimaryVarsToExport(int numVars, int numNodes)
void setNumberOfCellVarsToExport(int numVars, int numCells)
FloatArray & giveInternalXFEMVarInNode(int varNum, int eiNum, int nodeNum)
FloatArray & givePrimaryVarInNode(int varNum, int nodeNum)
UnknownType
Type representing particular unknown (its physical meaning).
virtual void giveCompositeExportData(std::vector< VTKPiece > &vtkPieces, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep)
virtual const char * giveClassName() const
Returns class name of the receiver.
NodalRecoveryModel * primVarSmoother
Smoother for primary variables.
Represents VTK (Visualization Toolkit) export module.
bool particleExportFlag
particle export flag
IntArray externalForcesToExport
List of primary unknowns to export.
void setPrimaryVarInNode(int varNum, int nodeNum, FloatArray valueArray)
virtual void giveCompositeExportData(VTKPiece &vtkPiece, IntArray &primaryVarsToExport, IntArray &internalVarsToExport, IntArray cellVarsToExport, TimeStep *tStep)
Class representing vector of real numbers.
IRResultType
Type defining the return values of InputRecord reading operations.
virtual const char * giveClassName() const
IntArray primaryVarsToExport
List of primary unknowns to export.
void setOffset(int cellNum, int offset)
IntArray & giveCellConnectivity(int cellNum)
std::list< std::string > gpPvdBuffer
Buffer for earlier time steps with gauss points exported to *.gp.pvd file.
NodalRecoveryModel::NodalRecoveryModelType stype
Smoother type.
std::vector< std::vector< FloatArray > > elVars
static IntArray redToFull
Map from Voigt to full tensor.
std::vector< VTKPiece > defaultVTKPieces
FloatArray & giveCellVar(int field, int cellNum)
void setConnectivity(int cellNum, IntArray &nodes)
std::vector< std::vector< FloatArray > > nodeLoads
Abstract base class representing the "problem" under consideration.
the oofem namespace is to define a context or scope in which all oofem names are defined.
FloatArray & giveNodeCoords(int nodeNum)
Class implementing node in finite element mesh.
The base class for all recovery models, which perform nodal averaging or projection processes for int...
Class representing solution step.
VTKXMLExportModuleElementInterface()