1 /*=========================================================================
3 Program: Visualization Toolkit
4 Module: $RCSfile: vtkPV3FoamReader.cxx,v $
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
14 =========================================================================*/
16 #include "vtkPV3FoamReader.h"
18 #include "pqApplicationCore.h"
19 #include "pqRenderView.h"
20 #include "pqServerManagerModel.h"
23 #include "vtkCallbackCommand.h"
24 #include "vtkDataArraySelection.h"
25 #include "vtkInformation.h"
26 #include "vtkInformationVector.h"
27 #include "vtkMultiBlockDataSet.h"
28 #include "vtkObjectFactory.h"
29 #include "vtkSMRenderViewProxy.h"
30 #include "vtkStreamingDemandDrivenPipeline.h"
31 #include "vtkStringArray.h"
34 #include "vtkPV3Foam.H"
36 vtkCxxRevisionMacro(vtkPV3FoamReader
, "$Revision: 1.5$");
37 vtkStandardNewMacro(vtkPV3FoamReader
);
39 #undef EXPERIMENTAL_TIME_CACHING
41 vtkPV3FoamReader::vtkPV3FoamReader()
44 vtkDebugMacro(<<"Constructor");
46 SetNumberOfInputPorts(0);
53 #ifdef VTKPV3FOAM_DUALPORT
54 // Add second output for the Lagrangian
55 this->SetNumberOfOutputPorts(2);
56 vtkMultiBlockDataSet
*lagrangian
= vtkMultiBlockDataSet::New();
57 lagrangian
->ReleaseData();
59 this->GetExecutive()->SetOutputData(1, lagrangian
);
68 ExtrapolatePatches
= 0;
75 PartSelection
= vtkDataArraySelection::New();
76 VolFieldSelection
= vtkDataArraySelection::New();
77 PointFieldSelection
= vtkDataArraySelection::New();
78 LagrangianFieldSelection
= vtkDataArraySelection::New();
80 // Setup the selection callback to modify this object when an array
81 // selection is changed.
82 SelectionObserver
= vtkCallbackCommand::New();
83 SelectionObserver
->SetCallback
85 &vtkPV3FoamReader::SelectionModifiedCallback
87 SelectionObserver
->SetClientData(this);
89 PartSelection
->AddObserver
91 vtkCommand::ModifiedEvent
,
92 this->SelectionObserver
94 VolFieldSelection
->AddObserver
96 vtkCommand::ModifiedEvent
,
97 this->SelectionObserver
99 PointFieldSelection
->AddObserver
101 vtkCommand::ModifiedEvent
,
102 this->SelectionObserver
104 LagrangianFieldSelection
->AddObserver
106 vtkCommand::ModifiedEvent
,
107 this->SelectionObserver
112 vtkPV3FoamReader::~vtkPV3FoamReader()
114 vtkDebugMacro(<<"Deconstructor");
129 PartSelection
->RemoveObserver(this->SelectionObserver
);
130 VolFieldSelection
->RemoveObserver(this->SelectionObserver
);
131 PointFieldSelection
->RemoveObserver(this->SelectionObserver
);
132 LagrangianFieldSelection
->RemoveObserver(this->SelectionObserver
);
134 SelectionObserver
->Delete();
136 PartSelection
->Delete();
137 VolFieldSelection
->Delete();
138 PointFieldSelection
->Delete();
139 LagrangianFieldSelection
->Delete();
143 // Do everything except set the output info
144 int vtkPV3FoamReader::RequestInformation
146 vtkInformation
* vtkNotUsed(request
),
147 vtkInformationVector
** vtkNotUsed(inputVector
),
148 vtkInformationVector
* outputVector
151 vtkDebugMacro(<<"RequestInformation");
153 if (Foam::vtkPV3Foam::debug
)
155 cout
<<"REQUEST_INFORMATION\n";
160 vtkErrorMacro("FileName has to be specified!");
164 int nInfo
= outputVector
->GetNumberOfInformationObjects();
166 if (Foam::vtkPV3Foam::debug
)
168 cout
<<"RequestInformation with " << nInfo
<< " item(s)\n";
169 for (int infoI
= 0; infoI
< nInfo
; ++infoI
)
171 outputVector
->GetInformationObject(infoI
)->Print(cout
);
177 foamData_
= new Foam::vtkPV3Foam(FileName
, this);
181 foamData_
->updateInfo();
185 double* timeSteps
= foamData_
->findTimes(nTimeSteps
);
189 vtkErrorMacro("could not find valid OpenFOAM mesh");
191 // delete foamData and flag it as fatal error
197 // set identical time steps for all ports
198 for (int infoI
= 0; infoI
< nInfo
; ++infoI
)
200 outputVector
->GetInformationObject(infoI
)->Set
202 vtkStreamingDemandDrivenPipeline::TIME_STEPS(),
211 timeRange
[0] = timeSteps
[0];
212 timeRange
[1] = timeSteps
[nTimeSteps
-1];
214 if (Foam::vtkPV3Foam::debug
> 1)
216 cout
<<"nTimeSteps " << nTimeSteps
<< "\n"
217 <<"timeRange " << timeRange
[0] << " to " << timeRange
[1]
220 for (int timeI
= 0; timeI
< nTimeSteps
; ++timeI
)
222 cout
<< "step[" << timeI
<< "] = " << timeSteps
[timeI
] << "\n";
226 for (int infoI
= 0; infoI
< nInfo
; ++infoI
)
228 outputVector
->GetInformationObject(infoI
)->Set
230 vtkStreamingDemandDrivenPipeline::TIME_RANGE(),
243 // Set the output info
244 int vtkPV3FoamReader::RequestData
246 vtkInformation
* vtkNotUsed(request
),
247 vtkInformationVector
** vtkNotUsed(inputVector
),
248 vtkInformationVector
* outputVector
251 vtkDebugMacro(<<"RequestData");
255 vtkErrorMacro("FileName has to be specified!");
259 // catch previous error
262 vtkErrorMacro("Reader failed - perhaps no mesh?");
266 int nInfo
= outputVector
->GetNumberOfInformationObjects();
268 if (Foam::vtkPV3Foam::debug
)
270 cout
<<"RequestData with " << nInfo
<< " item(s)\n";
271 for (int infoI
= 0; infoI
< nInfo
; ++infoI
)
273 outputVector
->GetInformationObject(infoI
)->Print(cout
);
277 // Get the requested time step.
278 // We only support requests for a single time step
280 int nRequestTime
= 0;
281 double requestTime
[nInfo
];
283 // taking port0 as the lead for other outputs would be nice, but fails when
284 // a filter is added - we need to check everything
285 // but since PREVIOUS_UPDATE_TIME_STEPS() is protected, relay the logic
286 // to the vtkPV3Foam::setTime() method
287 for (int infoI
= 0; infoI
< nInfo
; ++infoI
)
289 vtkInformation
*outInfo
= outputVector
->GetInformationObject(infoI
);
295 vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEPS()
299 vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEPS()
303 requestTime
[nRequestTime
++] = outInfo
->Get
305 vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEPS()
312 foamData_
->setTime(nRequestTime
, requestTime
);
316 vtkMultiBlockDataSet
* output
= vtkMultiBlockDataSet::SafeDownCast
318 outputVector
->GetInformationObject(0)->Get
320 vtkMultiBlockDataSet::DATA_OBJECT()
324 if (Foam::vtkPV3Foam::debug
)
326 cout
<< "update output with "
327 << output
->GetNumberOfBlocks() << " blocks\n";
331 #ifdef EXPERIMENTAL_TIME_CACHING
332 bool needsUpdate
= false;
336 output0_
= vtkMultiBlockDataSet::New();
340 // This experimental bit of code seems to work for the geometry,
341 // but trashes the fields and still triggers the GeometryFilter
344 foamData_
->Update(output
);
345 output0_
->ShallowCopy(output
);
349 output
->ShallowCopy(output0_
);
352 if (Foam::vtkPV3Foam::debug
)
356 cout
<< "full UPDATE ---------\n";
360 cout
<< "cached UPDATE ---------\n";
363 cout
<< "UPDATED output: ";
366 cout
<< "UPDATED output0_: ";
367 output0_
->Print(cout
);
372 #ifdef VTKPV3FOAM_DUALPORT
376 vtkMultiBlockDataSet::SafeDownCast
378 outputVector
->GetInformationObject(1)->Get
380 vtkMultiBlockDataSet::DATA_OBJECT()
385 foamData_
->Update(output
, output
);
390 addPatchNamesToView();
394 removePatchNamesFromView();
399 // Do any cleanup on the Foam side
400 foamData_
->CleanUp();
406 void vtkPV3FoamReader::addPatchNamesToView()
408 pqApplicationCore
* appCore
= pqApplicationCore::instance();
410 // Server manager model for querying items in the server manager
411 pqServerManagerModel
* smModel
= appCore
->getServerManagerModel();
413 // Get all the pqRenderView instances
414 QList
<pqRenderView
*> renderViews
= smModel
->findItems
<pqRenderView
*>();
416 for (int viewI
=0; viewI
<renderViews
.size(); viewI
++)
418 foamData_
->addPatchNames
420 renderViews
[viewI
]->getRenderViewProxy()->GetRenderer()
426 void vtkPV3FoamReader::removePatchNamesFromView()
428 pqApplicationCore
* appCore
= pqApplicationCore::instance();
430 // Server manager model for querying items in the server manager
431 pqServerManagerModel
* smModel
= appCore
->getServerManagerModel();
433 // Get all the pqRenderView instances
434 QList
<pqRenderView
*> renderViews
= smModel
->findItems
<pqRenderView
*>();
436 for (int viewI
=0; viewI
<renderViews
.size(); viewI
++)
438 foamData_
->removePatchNames
440 renderViews
[viewI
]->getRenderViewProxy()->GetRenderer()
446 void vtkPV3FoamReader::PrintSelf(ostream
& os
, vtkIndent indent
)
448 vtkDebugMacro(<<"PrintSelf");
450 this->Superclass::PrintSelf(os
,indent
);
451 os
<< indent
<< "File name: "
452 << (this->FileName
? this->FileName
: "(none)") << "\n";
454 foamData_
->PrintSelf(os
, indent
);
456 os
<< indent
<< "Time step range: "
457 << this->TimeStepRange
[0] << " - " << this->TimeStepRange
[1]
459 os
<< indent
<< "Time step: " << this->GetTimeStep() << endl
;
463 int vtkPV3FoamReader::GetTimeStep()
465 return foamData_
? foamData_
->timeIndex() : -1;
469 // ----------------------------------------------------------------------
470 // Parts selection list control
472 vtkDataArraySelection
* vtkPV3FoamReader::GetPartSelection()
474 vtkDebugMacro(<<"GetPartSelection");
475 return PartSelection
;
479 int vtkPV3FoamReader::GetNumberOfPartArrays()
481 vtkDebugMacro(<<"GetNumberOfPartArrays");
482 return PartSelection
->GetNumberOfArrays();
486 const char* vtkPV3FoamReader::GetPartArrayName(int index
)
488 vtkDebugMacro(<<"GetPartArrayName");
489 return PartSelection
->GetArrayName(index
);
493 int vtkPV3FoamReader::GetPartArrayStatus(const char* name
)
495 vtkDebugMacro(<<"GetPartArrayStatus");
496 return PartSelection
->ArrayIsEnabled(name
);
500 void vtkPV3FoamReader::SetPartArrayStatus(const char* name
, int status
)
502 vtkDebugMacro(<<"SetPartArrayStatus");
505 PartSelection
->EnableArray(name
);
509 PartSelection
->DisableArray(name
);
514 // ----------------------------------------------------------------------
515 // volField selection list control
517 vtkDataArraySelection
* vtkPV3FoamReader::GetVolFieldSelection()
519 vtkDebugMacro(<<"GetVolFieldSelection");
520 return VolFieldSelection
;
524 int vtkPV3FoamReader::GetNumberOfVolFieldArrays()
526 vtkDebugMacro(<<"GetNumberOfVolFieldArrays");
527 return VolFieldSelection
->GetNumberOfArrays();
531 const char* vtkPV3FoamReader::GetVolFieldArrayName(int index
)
533 vtkDebugMacro(<<"GetVolFieldArrayName");
534 return VolFieldSelection
->GetArrayName(index
);
538 int vtkPV3FoamReader::GetVolFieldArrayStatus(const char* name
)
540 vtkDebugMacro(<<"GetVolFieldArrayStatus");
541 return VolFieldSelection
->ArrayIsEnabled(name
);
545 void vtkPV3FoamReader::SetVolFieldArrayStatus(const char* name
, int status
)
547 vtkDebugMacro(<<"SetVolFieldArrayStatus");
550 VolFieldSelection
->EnableArray(name
);
554 VolFieldSelection
->DisableArray(name
);
559 // ----------------------------------------------------------------------
560 // pointField selection list control
562 vtkDataArraySelection
* vtkPV3FoamReader::GetPointFieldSelection()
564 vtkDebugMacro(<<"GetPointFieldSelection");
565 return PointFieldSelection
;
569 int vtkPV3FoamReader::GetNumberOfPointFieldArrays()
571 vtkDebugMacro(<<"GetNumberOfPointFieldArrays");
572 return PointFieldSelection
->GetNumberOfArrays();
576 const char* vtkPV3FoamReader::GetPointFieldArrayName(int index
)
578 vtkDebugMacro(<<"GetPointFieldArrayName");
579 return PointFieldSelection
->GetArrayName(index
);
583 int vtkPV3FoamReader::GetPointFieldArrayStatus(const char* name
)
585 vtkDebugMacro(<<"GetPointFieldArrayStatus");
586 return PointFieldSelection
->ArrayIsEnabled(name
);
590 void vtkPV3FoamReader::SetPointFieldArrayStatus(const char* name
, int status
)
592 vtkDebugMacro(<<"SetPointFieldArrayStatus");
595 PointFieldSelection
->EnableArray(name
);
599 PointFieldSelection
->DisableArray(name
);
604 // ----------------------------------------------------------------------
605 // lagrangianField selection list control
607 vtkDataArraySelection
* vtkPV3FoamReader::GetLagrangianFieldSelection()
609 vtkDebugMacro(<<"GetLagrangianFieldSelection");
610 return LagrangianFieldSelection
;
614 int vtkPV3FoamReader::GetNumberOfLagrangianFieldArrays()
616 vtkDebugMacro(<<"GetNumberOfLagrangianFieldArrays");
617 return LagrangianFieldSelection
->GetNumberOfArrays();
621 const char* vtkPV3FoamReader::GetLagrangianFieldArrayName(int index
)
623 vtkDebugMacro(<<"GetLagrangianFieldArrayName");
624 return LagrangianFieldSelection
->GetArrayName(index
);
628 int vtkPV3FoamReader::GetLagrangianFieldArrayStatus(const char* name
)
630 vtkDebugMacro(<<"GetLagrangianFieldArrayStatus");
631 return LagrangianFieldSelection
->ArrayIsEnabled(name
);
635 void vtkPV3FoamReader::SetLagrangianFieldArrayStatus
641 vtkDebugMacro(<<"SetLagrangianFieldArrayStatus");
644 LagrangianFieldSelection
->EnableArray(name
);
648 LagrangianFieldSelection
->DisableArray(name
);
653 // ----------------------------------------------------------------------
655 void vtkPV3FoamReader::SelectionModifiedCallback
663 static_cast<vtkPV3FoamReader
*>(clientdata
)->SelectionModified();
667 void vtkPV3FoamReader::SelectionModified()
669 vtkDebugMacro(<<"SelectionModified");
674 int vtkPV3FoamReader::FillOutputPortInformation
682 return this->Superclass::FillOutputPortInformation(port
, info
);
684 info
->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
689 // ************************************************************************* //