/*========================================================================= Program: Visualization Toolkit Module: $RCSfile: vtkPVFoamReader.cxx,v $ Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen All rights reserved. See Copyright.txt or http://www.kitware.com/Copyright.htm for details. This software is distributed WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the above copyright notice for more information. =========================================================================*/ #include "vtkPVFoamReader.h" #include "pqApplicationCore.h" #include "pqRenderView.h" #include "pqServerManagerModel.h" // VTK includes #include "vtkCallbackCommand.h" #include "vtkDataArraySelection.h" #include "vtkInformation.h" #include "vtkInformationVector.h" #include "vtkMultiBlockDataSet.h" #include "vtkObjectFactory.h" #include "vtkSMRenderViewProxy.h" #include "vtkStreamingDemandDrivenPipeline.h" #include "vtkStringArray.h" // Foam includes #include "vtkPVFoam.H" vtkStandardNewMacro(vtkPVFoamReader); #undef EXPERIMENTAL_TIME_CACHING vtkPVFoamReader::vtkPVFoamReader() { Debug = 0; vtkDebugMacro(<<"Constructor"); SetNumberOfInputPorts(0); FileName = NULL; foamData_ = NULL; output0_ = NULL; #ifdef VTKPVFOAM_DUALPORT // Add second output for the Lagrangian this->SetNumberOfOutputPorts(2); vtkMultiBlockDataSet *lagrangian = vtkMultiBlockDataSet::New(); lagrangian->ReleaseData(); this->GetExecutive()->SetOutputData(1, lagrangian); lagrangian->Delete(); #endif TimeStepRange[0] = 0; TimeStepRange[1] = 0; CacheMesh = 1; ExtrapolatePatches = 0; UseVTKPolyhedron = 0; IncludeSets = 0; IncludeZones = 0; ShowPatchNames = 0; UpdateGUI = 0; PartSelection = vtkDataArraySelection::New(); VolFieldSelection = vtkDataArraySelection::New(); PointFieldSelection = vtkDataArraySelection::New(); LagrangianFieldSelection = vtkDataArraySelection::New(); // Setup the selection callback to modify this object when an array // selection is changed. SelectionObserver = vtkCallbackCommand::New(); SelectionObserver->SetCallback ( &vtkPVFoamReader::SelectionModifiedCallback ); SelectionObserver->SetClientData(this); PartSelection->AddObserver ( vtkCommand::ModifiedEvent, this->SelectionObserver ); VolFieldSelection->AddObserver ( vtkCommand::ModifiedEvent, this->SelectionObserver ); PointFieldSelection->AddObserver ( vtkCommand::ModifiedEvent, this->SelectionObserver ); LagrangianFieldSelection->AddObserver ( vtkCommand::ModifiedEvent, this->SelectionObserver ); } vtkPVFoamReader::~vtkPVFoamReader() { vtkDebugMacro(<<"Deconstructor"); delete foamData_; if (FileName) { delete [] FileName; } if (output0_) { output0_->Delete(); } PartSelection->RemoveObserver(this->SelectionObserver); VolFieldSelection->RemoveObserver(this->SelectionObserver); PointFieldSelection->RemoveObserver(this->SelectionObserver); LagrangianFieldSelection->RemoveObserver(this->SelectionObserver); SelectionObserver->Delete(); PartSelection->Delete(); VolFieldSelection->Delete(); PointFieldSelection->Delete(); LagrangianFieldSelection->Delete(); } // Do everything except set the output info int vtkPVFoamReader::RequestInformation ( vtkInformation* vtkNotUsed(request), vtkInformationVector** vtkNotUsed(inputVector), vtkInformationVector* outputVector ) { vtkDebugMacro(<<"RequestInformation"); if (Foam::vtkPVFoam::debug) { cout<<"REQUEST_INFORMATION\n"; } if (!FileName) { vtkErrorMacro("FileName has to be specified!"); return 0; } int nInfo = outputVector->GetNumberOfInformationObjects(); if (Foam::vtkPVFoam::debug) { cout<<"RequestInformation with " << nInfo << " item(s)\n"; for (int infoI = 0; infoI < nInfo; ++infoI) { outputVector->GetInformationObject(infoI)->Print(cout); } } if (!foamData_) { foamData_ = new Foam::vtkPVFoam(FileName, this); } else { foamData_->updateInfo(); } int nTimeSteps = 0; double* timeSteps = foamData_->findTimes(nTimeSteps); if (!nTimeSteps) { vtkErrorMacro("could not find valid OpenFOAM mesh"); // delete foamData and flag it as fatal error delete foamData_; foamData_ = NULL; return 0; } // set identical time steps for all ports for (int infoI = 0; infoI < nInfo; ++infoI) { outputVector->GetInformationObject(infoI)->Set ( vtkStreamingDemandDrivenPipeline::TIME_STEPS(), timeSteps, nTimeSteps ); } if (nTimeSteps) { double timeRange[2]; timeRange[0] = timeSteps[0]; timeRange[1] = timeSteps[nTimeSteps-1]; if (Foam::vtkPVFoam::debug > 1) { cout<<"nTimeSteps " << nTimeSteps << "\n" <<"timeRange " << timeRange[0] << " to " << timeRange[1] << "\n"; for (int timeI = 0; timeI < nTimeSteps; ++timeI) { cout<< "step[" << timeI << "] = " << timeSteps[timeI] << "\n"; } } for (int infoI = 0; infoI < nInfo; ++infoI) { outputVector->GetInformationObject(infoI)->Set ( vtkStreamingDemandDrivenPipeline::TIME_RANGE(), timeRange, 2 ); } } delete timeSteps; return 1; } // Set the output info int vtkPVFoamReader::RequestData ( vtkInformation* vtkNotUsed(request), vtkInformationVector** vtkNotUsed(inputVector), vtkInformationVector* outputVector ) { vtkDebugMacro(<<"RequestData"); if (!FileName) { vtkErrorMacro("FileName has to be specified!"); return 0; } // catch previous error if (!foamData_) { vtkErrorMacro("Reader failed - perhaps no mesh?"); return 0; } int nInfo = outputVector->GetNumberOfInformationObjects(); if (Foam::vtkPVFoam::debug) { cout<<"RequestData with " << nInfo << " item(s)\n"; for (int infoI = 0; infoI < nInfo; ++infoI) { outputVector->GetInformationObject(infoI)->Print(cout); } } // Get the requested time step. // We only support requests for a single time step int nRequestTime = 0; double requestTime[nInfo]; // taking port0 as the lead for other outputs would be nice, but fails when // a filter is added - we need to check everything // but since PREVIOUS_UPDATE_TIME_STEPS() is protected, relay the logic // to the vtkPVFoam::setTime() method for (int infoI = 0; infoI < nInfo; ++infoI) { vtkInformation *outInfo = outputVector->GetInformationObject(infoI); if ( outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()) && outInfo->Length(vtkStreamingDemandDrivenPipeline::TIME_STEPS()) >= 1 ) { requestTime[nRequestTime++] = outInfo->Get ( vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP() ); } } if (nRequestTime) { foamData_->setTime(nRequestTime, requestTime); } vtkMultiBlockDataSet* output = vtkMultiBlockDataSet::SafeDownCast ( outputVector->GetInformationObject(0)->Get ( vtkMultiBlockDataSet::DATA_OBJECT() ) ); if (Foam::vtkPVFoam::debug) { cout<< "update output with " << output->GetNumberOfBlocks() << " blocks\n"; } #ifdef EXPERIMENTAL_TIME_CACHING bool needsUpdate = false; if (!output0_) { output0_ = vtkMultiBlockDataSet::New(); needsUpdate = true; } // This experimental bit of code seems to work for the geometry, // but trashes the fields and still triggers the GeometryFilter if (needsUpdate) { foamData_->Update(output); output0_->ShallowCopy(output); } else { output->ShallowCopy(output0_); } if (Foam::vtkPVFoam::debug) { if (needsUpdate) { cout<< "full UPDATE ---------\n"; } else { cout<< "cached UPDATE ---------\n"; } cout<< "UPDATED output: "; output->Print(cout); cout<< "UPDATED output0_: "; output0_->Print(cout); } #else #ifdef VTKPVFOAM_DUALPORT foamData_->Update ( output, vtkMultiBlockDataSet::SafeDownCast ( outputVector->GetInformationObject(1)->Get ( vtkMultiBlockDataSet::DATA_OBJECT() ) ); ); #else foamData_->Update(output, output); #endif if (ShowPatchNames) { addPatchNamesToView(); } else { removePatchNamesFromView(); } #endif // Do any cleanup on the Foam side foamData_->CleanUp(); return 1; } void vtkPVFoamReader::addPatchNamesToView() { pqApplicationCore* appCore = pqApplicationCore::instance(); // Server manager model for querying items in the server manager pqServerManagerModel* smModel = appCore->getServerManagerModel(); // Get all the pqRenderView instances QList renderViews = smModel->findItems(); for (int viewI=0; viewIaddPatchNames ( renderViews[viewI]->getRenderViewProxy()->GetRenderer() ); } } void vtkPVFoamReader::removePatchNamesFromView() { pqApplicationCore* appCore = pqApplicationCore::instance(); // Server manager model for querying items in the server manager pqServerManagerModel* smModel = appCore->getServerManagerModel(); // Get all the pqRenderView instances QList renderViews = smModel->findItems(); for (int viewI=0; viewIremovePatchNames ( renderViews[viewI]->getRenderViewProxy()->GetRenderer() ); } } void vtkPVFoamReader::PrintSelf(ostream& os, vtkIndent indent) { vtkDebugMacro(<<"PrintSelf"); this->Superclass::PrintSelf(os,indent); os<< indent << "File name: " << (this->FileName ? this->FileName : "(none)") << "\n"; foamData_->PrintSelf(os, indent); os<< indent << "Time step range: " << this->TimeStepRange[0] << " - " << this->TimeStepRange[1] << "\n"; os<< indent << "Time step: " << this->GetTimeStep() << endl; } int vtkPVFoamReader::GetTimeStep() { return foamData_ ? foamData_->timeIndex() : -1; } // ---------------------------------------------------------------------- // Parts selection list control vtkDataArraySelection* vtkPVFoamReader::GetPartSelection() { vtkDebugMacro(<<"GetPartSelection"); return PartSelection; } int vtkPVFoamReader::GetNumberOfPartArrays() { vtkDebugMacro(<<"GetNumberOfPartArrays"); return PartSelection->GetNumberOfArrays(); } const char* vtkPVFoamReader::GetPartArrayName(int index) { vtkDebugMacro(<<"GetPartArrayName"); return PartSelection->GetArrayName(index); } int vtkPVFoamReader::GetPartArrayStatus(const char* name) { vtkDebugMacro(<<"GetPartArrayStatus"); return PartSelection->ArrayIsEnabled(name); } void vtkPVFoamReader::SetPartArrayStatus(const char* name, int status) { vtkDebugMacro(<<"SetPartArrayStatus"); if (status) { PartSelection->EnableArray(name); } else { PartSelection->DisableArray(name); } } // ---------------------------------------------------------------------- // volField selection list control vtkDataArraySelection* vtkPVFoamReader::GetVolFieldSelection() { vtkDebugMacro(<<"GetVolFieldSelection"); return VolFieldSelection; } int vtkPVFoamReader::GetNumberOfVolFieldArrays() { vtkDebugMacro(<<"GetNumberOfVolFieldArrays"); return VolFieldSelection->GetNumberOfArrays(); } const char* vtkPVFoamReader::GetVolFieldArrayName(int index) { vtkDebugMacro(<<"GetVolFieldArrayName"); return VolFieldSelection->GetArrayName(index); } int vtkPVFoamReader::GetVolFieldArrayStatus(const char* name) { vtkDebugMacro(<<"GetVolFieldArrayStatus"); return VolFieldSelection->ArrayIsEnabled(name); } void vtkPVFoamReader::SetVolFieldArrayStatus(const char* name, int status) { vtkDebugMacro(<<"SetVolFieldArrayStatus"); if (status) { VolFieldSelection->EnableArray(name); } else { VolFieldSelection->DisableArray(name); } } // ---------------------------------------------------------------------- // pointField selection list control vtkDataArraySelection* vtkPVFoamReader::GetPointFieldSelection() { vtkDebugMacro(<<"GetPointFieldSelection"); return PointFieldSelection; } int vtkPVFoamReader::GetNumberOfPointFieldArrays() { vtkDebugMacro(<<"GetNumberOfPointFieldArrays"); return PointFieldSelection->GetNumberOfArrays(); } const char* vtkPVFoamReader::GetPointFieldArrayName(int index) { vtkDebugMacro(<<"GetPointFieldArrayName"); return PointFieldSelection->GetArrayName(index); } int vtkPVFoamReader::GetPointFieldArrayStatus(const char* name) { vtkDebugMacro(<<"GetPointFieldArrayStatus"); return PointFieldSelection->ArrayIsEnabled(name); } void vtkPVFoamReader::SetPointFieldArrayStatus(const char* name, int status) { vtkDebugMacro(<<"SetPointFieldArrayStatus"); if (status) { PointFieldSelection->EnableArray(name); } else { PointFieldSelection->DisableArray(name); } } // ---------------------------------------------------------------------- // lagrangianField selection list control vtkDataArraySelection* vtkPVFoamReader::GetLagrangianFieldSelection() { vtkDebugMacro(<<"GetLagrangianFieldSelection"); return LagrangianFieldSelection; } int vtkPVFoamReader::GetNumberOfLagrangianFieldArrays() { vtkDebugMacro(<<"GetNumberOfLagrangianFieldArrays"); return LagrangianFieldSelection->GetNumberOfArrays(); } const char* vtkPVFoamReader::GetLagrangianFieldArrayName(int index) { vtkDebugMacro(<<"GetLagrangianFieldArrayName"); return LagrangianFieldSelection->GetArrayName(index); } int vtkPVFoamReader::GetLagrangianFieldArrayStatus(const char* name) { vtkDebugMacro(<<"GetLagrangianFieldArrayStatus"); return LagrangianFieldSelection->ArrayIsEnabled(name); } void vtkPVFoamReader::SetLagrangianFieldArrayStatus ( const char* name, int status ) { vtkDebugMacro(<<"SetLagrangianFieldArrayStatus"); if (status) { LagrangianFieldSelection->EnableArray(name); } else { LagrangianFieldSelection->DisableArray(name); } } // ---------------------------------------------------------------------- void vtkPVFoamReader::SelectionModifiedCallback ( vtkObject*, unsigned long, void* clientdata, void* ) { static_cast(clientdata)->SelectionModified(); } void vtkPVFoamReader::SelectionModified() { vtkDebugMacro(<<"SelectionModified"); Modified(); } int vtkPVFoamReader::FillOutputPortInformation ( int port, vtkInformation* info ) { if (port == 0) { return this->Superclass::FillOutputPortInformation(port, info); } info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet"); return 1; } // ************************************************************************* //