1 // Copyright (C) 2021 CEA/DEN, EDF R&D
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // Author : Anthony Geay (EDF R&D)
21 #include "vtkZJFilter.h"
23 #include <vtkAdjacentVertexIterator.h>
24 #include <vtkAlgorithmOutput.h>
25 #include <vtkCallbackCommand.h>
27 #include <vtkCellData.h>
28 #include <vtkCellType.h>
29 #include <vtkCharArray.h>
30 #include <vtkCompositeDataToUnstructuredGridFilter.h>
31 #include <vtkDataArraySelection.h>
32 #include <vtkAOSDataArrayTemplate.h>
33 #include <vtkDataObjectTreeIterator.h>
34 #include <vtkDataSet.h>
35 #include <vtkDataSetAttributes.h>
36 #include <vtkDemandDrivenPipeline.h>
37 #include <vtkDoubleArray.h>
38 #include <vtkExecutive.h>
39 #include <vtkFloatArray.h>
40 #include <vtkInformation.h>
41 #include <vtkInformationDataObjectKey.h>
42 #include <vtkInformationDataObjectMetaDataKey.h>
43 #include <vtkInformationStringKey.h>
44 #include <vtkInformationVector.h>
45 #include <vtkLongArray.h>
46 #include <vtkMultiBlockDataSet.h>
47 #include <vtkMutableDirectedGraph.h>
48 #include <vtkObjectFactory.h>
49 #include <vtkPointData.h>
50 #include <vtkStreamingDemandDrivenPipeline.h>
51 #include <vtkStringArray.h>
52 #include <vtkThreshold.h>
53 #include <vtkUnstructuredGrid.h>
55 #include "InterpKernelException.hxx"
56 #include "MEDCouplingRefCountObject.hxx"
63 vtkStandardNewMacro(vtkZJFilter);
67 vtkInformationDataObjectMetaDataKey* GetMEDReaderMetaDataIfAny()
69 static const char ZE_KEY[] = "vtkMEDReader::META_DATA";
70 MEDCoupling::GlobalDict* gd(MEDCoupling::GlobalDict::GetInstance());
71 if (!gd->hasKey(ZE_KEY))
73 std::string ptSt(gd->value(ZE_KEY));
75 std::istringstream iss(ptSt);
77 return reinterpret_cast<vtkInformationDataObjectMetaDataKey*>(pt);
80 bool IsInformationOK(vtkInformation* info)
82 vtkInformationDataObjectMetaDataKey* key(GetMEDReaderMetaDataIfAny());
85 // Check the information contain meta data key
89 vtkMutableDirectedGraph* sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(key)));
93 vtkAbstractArray* verticesNames(sil->GetVertexData()->GetAbstractArray("Names", idNames));
94 vtkStringArray* verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
97 for (int i = 0; i < verticesNames2->GetNumberOfValues(); i++)
99 vtkStdString& st(verticesNames2->GetValue(i));
100 if (st == "MeshesFamsGrps")
109 Grp(const std::string& name)
113 void setFamilies(const std::vector<std::string>& fams) { _fams = fams; }
114 std::string getName() const { return _name; }
115 std::vector<std::string> getFamilies() const { return _fams; }
119 std::vector<std::string> _fams;
125 Fam(const std::string& name)
127 constexpr char ZE_SEP[] = "@@][@@";
128 std::size_t pos(name.find(ZE_SEP));
129 std::string name0(name.substr(0, pos)), name1(name.substr(pos + strlen(ZE_SEP)));
130 std::istringstream iss(name1);
134 std::string getName() const { return _name; }
135 int getID() const { return _id; }
142 void ExtractInfo(vtkInformationVector* inputVector, vtkUnstructuredGrid*& usgIn)
144 vtkInformation* inputInfo(inputVector->GetInformationObject(0));
145 vtkDataSet* input(0);
146 vtkDataSet* input0(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
147 vtkMultiBlockDataSet* input1(
148 vtkMultiBlockDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
154 throw INTERP_KERNEL::Exception(
155 "Input dataSet must be a DataSet or single elt multi block dataset expected !");
156 if (input1->GetNumberOfBlocks() != 1)
157 throw INTERP_KERNEL::Exception(
158 "Input dataSet is a multiblock dataset with not exactly one block ! Use MergeBlocks or "
159 "ExtractBlocks filter before calling this filter !");
160 vtkDataObject* input2(input1->GetBlock(0));
162 throw INTERP_KERNEL::Exception("Input dataSet is a multiblock dataset with exactly one block "
163 "but this single element is NULL !");
164 vtkDataSet* input2c(vtkDataSet::SafeDownCast(input2));
166 throw INTERP_KERNEL::Exception(
167 "Input dataSet is a multiblock dataset with exactly one block but this single element is "
168 "not a dataset ! Use MergeBlocks or ExtractBlocks filter before calling this filter !");
172 throw INTERP_KERNEL::Exception("Input data set is NULL !");
173 usgIn = vtkUnstructuredGrid::SafeDownCast(input);
175 throw INTERP_KERNEL::Exception("Input data set is not an unstructured mesh ! This filter works "
176 "only on unstructured meshes !");
179 void LoadFamGrpMapInfo(vtkMutableDirectedGraph* sil, std::string& meshName,
180 std::vector<Grp>& groups, std::vector<Fam>& fams)
183 throw INTERP_KERNEL::Exception("LoadFamGrpMapInfo : internal error !");
185 vtkAbstractArray* verticesNames(sil->GetVertexData()->GetAbstractArray("Names", idNames));
186 vtkStringArray* verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
189 for (int i = 0; i < verticesNames2->GetNumberOfValues(); i++)
191 vtkStdString& st(verticesNames2->GetValue(i));
192 if (st == "MeshesFamsGrps")
199 throw INTERP_KERNEL::Exception(
200 "There is an internal error ! The tree on server side has not the expected look !");
201 vtkAdjacentVertexIterator* it0(vtkAdjacentVertexIterator::New());
202 sil->GetAdjacentVertices(id0, it0);
204 while (it0->HasNext())
206 vtkIdType id1(it0->Next());
207 std::string mName((const char*)verticesNames2->GetValue(id1));
209 vtkAdjacentVertexIterator* it1(vtkAdjacentVertexIterator::New());
210 sil->GetAdjacentVertices(id1, it1);
211 vtkIdType idZeGrps(it1->Next()); // zeGroups
212 vtkAdjacentVertexIterator* itGrps(vtkAdjacentVertexIterator::New());
213 sil->GetAdjacentVertices(idZeGrps, itGrps);
214 while (itGrps->HasNext())
216 vtkIdType idg(itGrps->Next());
217 Grp grp((const char*)verticesNames2->GetValue(idg));
218 vtkAdjacentVertexIterator* itGrps2(vtkAdjacentVertexIterator::New());
219 sil->GetAdjacentVertices(idg, itGrps2);
220 std::vector<std::string> famsOnGroup;
221 while (itGrps2->HasNext())
223 vtkIdType idgf(itGrps2->Next());
224 famsOnGroup.push_back(std::string((const char*)verticesNames2->GetValue(idgf)));
226 grp.setFamilies(famsOnGroup);
228 groups.push_back(grp);
231 vtkIdType idZeFams(it1->Next()); // zeFams
233 vtkAdjacentVertexIterator* itFams(vtkAdjacentVertexIterator::New());
234 sil->GetAdjacentVertices(idZeFams, itFams);
235 while (itFams->HasNext())
237 vtkIdType idf(itFams->Next());
238 Fam fam((const char*)verticesNames2->GetValue(idf));
246 std::vector<std::string> FindConds(const std::vector<Grp>& grps)
248 constexpr char PAT[] = "COND_";
249 constexpr std::size_t SZ_PAT(sizeof(PAT) - 1);
250 std::vector<std::string> ret;
251 for (std::vector<Grp>::const_iterator it = grps.begin(); it != grps.end(); it++)
253 std::string name((*it).getName());
254 std::string part(name.substr(0, SZ_PAT));
256 ret.push_back(name.substr(SZ_PAT, std::string::npos));
261 constexpr char EPORT_PAT[] = "EPORT_";
263 std::vector<std::string> FindEports(
264 const std::string& condEntry, const std::vector<Grp>& grps, std::vector<std::string>& eportsZip)
266 std::vector<std::string> ret;
267 std::string commonPart(std::string(EPORT_PAT) + condEntry);
268 std::size_t commonPart_sz(commonPart.length());
269 for (std::vector<Grp>::const_iterator it = grps.begin(); it != grps.end(); it++)
271 std::string name((*it).getName());
272 std::string part(name.substr(0, commonPart_sz));
273 if (part == commonPart)
276 eportsZip.push_back(name.substr(commonPart_sz, std::string::npos));
282 std::string BigestCommonPart(const std::string& s1, const std::string& s2)
284 std::size_t ls1(s1.length()), ls2(s2.length()), lb(0), lt(0);
300 for (std::size_t l0 = lt; l0 > 0; l0--)
302 for (std::size_t l1 = 0; l1 < lt - l0 + 1; l1++)
304 std::string cand(t.substr(l1, l0));
305 if (b.find(cand) != std::string::npos)
309 return std::string();
312 std::vector<int> DeduceIdsFrom(const std::vector<std::string>& eportsZip)
314 if (eportsZip.empty())
315 return std::vector<int>();
316 std::string ref(eportsZip[0]);
317 std::size_t sz(eportsZip.size());
318 for (std::size_t i = 1; i < sz; i++)
319 ref = BigestCommonPart(ref, eportsZip[i]);
320 std::vector<int> ret(sz);
321 for (std::size_t i = 0; i < sz; i++)
323 std::size_t pos(eportsZip[i].find(ref));
324 if (pos == std::string::npos)
325 throw INTERP_KERNEL::Exception("DeduceIdsFrom : internal error !");
327 eportsZip[i].substr(0, pos) + eportsZip[i].substr(pos + ref.length(), std::string::npos));
328 std::istringstream iss(res);
336 std::set<int> FamiliesIdsFromGrp(
337 const std::vector<Grp>& grps, const std::vector<Fam>& fams, const std::string& grp)
340 std::vector<std::string> locFams;
341 for (std::vector<Grp>::const_iterator it = grps.begin(); it != grps.end(); it++)
343 if ((*it).getName() == grp)
345 locFams = (*it).getFamilies();
351 throw INTERP_KERNEL::Exception("FamiliesIdsFromGrp : internal error !");
353 for (std::vector<Fam>::const_iterator it = fams.begin(); it != fams.end(); it++)
355 if (std::find(locFams.begin(), locFams.end(), (*it).getName()) != locFams.end())
356 ret.insert((*it).getID());
361 vtkDataSet* FilterFamilies(vtkZJFilter* zeBoss, vtkDataSet* input, const std::set<int>& idsToKeep)
363 bool catchAll, catchSmth;
364 vtkNew<vtkThreshold> thres;
365 thres->AddObserver(vtkCommand::ProgressEvent, zeBoss->InternalProgressObserver);
366 constexpr int VTK_DATA_ARRAY_DELETE = vtkAOSDataArrayTemplate<double>::VTK_DATA_ARRAY_DELETE;
367 constexpr char ZE_SELECTION_ARR_NAME[] = "@@ZeSelection@@";
368 constexpr char arrNameOfFamilyField[] = "FamilyIdCell";
369 constexpr char associationForThreshold[] = "vtkDataObject::FIELD_ASSOCIATION_CELLS";
370 vtkDataSet* output(input->NewInstance());
371 output->ShallowCopy(input);
372 thres->SetInputData(output);
373 vtkDataSetAttributes *dscIn(input->GetCellData()), *dscIn2(input->GetPointData());
374 vtkDataSetAttributes *dscOut(output->GetCellData()), *dscOut2(output->GetPointData());
376 constexpr double vMin(1.), vMax(2.);
377 thres->ThresholdBetween(vMin, vMax);
380 vtkDataArray* da(input->GetCellData()->GetScalars(arrNameOfFamilyField));
383 std::string daName(da->GetName());
384 vtkLongArray* dai(vtkLongArray::SafeDownCast(da));
385 if (daName != arrNameOfFamilyField || !dai)
388 int nbOfTuples(dai->GetNumberOfTuples());
389 vtkCharArray* zeSelection(vtkCharArray::New());
390 zeSelection->SetName(ZE_SELECTION_ARR_NAME);
391 zeSelection->SetNumberOfComponents(1);
392 char* pt(new char[nbOfTuples]);
393 zeSelection->SetArray(pt, nbOfTuples, 0, VTK_DATA_ARRAY_DELETE);
394 const long* inPtr(dai->GetPointer(0));
395 std::fill(pt, pt + nbOfTuples, 0);
398 std::vector<bool> pt2(nbOfTuples, false);
399 for (std::set<int>::const_iterator it = idsToKeep.begin(); it != idsToKeep.end(); it++)
401 bool catchFid(false);
402 for (int i = 0; i < nbOfTuples; i++)
413 for (int ii = 0; ii < nbOfTuples; ii++)
416 int idx(output->GetCellData()->AddArray(zeSelection));
417 output->GetCellData()->SetActiveAttribute(idx, vtkDataSetAttributes::SCALARS);
418 output->GetCellData()->CopyScalarsOff();
419 zeSelection->Delete();
421 thres->SetInputArrayToProcess(idx, 0, 0, associationForThreshold, ZE_SELECTION_ARR_NAME);
423 vtkUnstructuredGrid* zeComputedOutput(thres->GetOutput());
424 zeComputedOutput->GetCellData()->RemoveArray(idx);
426 zeComputedOutput->Register(0);
427 thres->RemoveObserver(zeBoss->InternalProgressObserver);
428 return zeComputedOutput;
433 vtkZJFilter::vtkZJFilter()
434 : InternalProgressObserver(0)
436 this->InternalProgressObserver = vtkCallbackCommand::New();
437 this->InternalProgressObserver->SetCallback(&vtkZJFilter::InternalProgressCallbackFunction);
438 this->InternalProgressObserver->SetClientData(this);
441 vtkZJFilter::~vtkZJFilter()
443 this->InternalProgressObserver->Delete();
446 void vtkZJFilter::InternalProgressCallbackFunction(
447 vtkObject* arg, unsigned long, void* clientdata, void*)
449 reinterpret_cast<vtkZJFilter*>(clientdata)
450 ->InternalProgressCallback(static_cast<vtkAlgorithm*>(arg));
453 void vtkZJFilter::InternalProgressCallback(vtkAlgorithm* algorithm)
455 /*this->UpdateProgress(algorithm->GetProgress()); // To intercept progression of Threshold filters
456 if (this->AbortExecute)
458 algorithm->SetAbortExecute(1);
462 int vtkZJFilter::RequestInformation(
463 vtkInformation* request, vtkInformationVector** inputVector, vtkInformationVector* outputVector)
465 // std::cerr << "########################################## vtkZJFilter::RequestInformation
466 // ##########################################" << std::endl;
469 vtkUnstructuredGrid* usgIn = nullptr;
470 ExtractInfo(inputVector[0], usgIn);
472 catch (INTERP_KERNEL::Exception& e)
474 std::ostringstream oss;
475 oss << "Exception has been thrown in vtkZJFilter::RequestInformation : " << e.what()
477 if (this->HasObserver("ErrorEvent"))
479 this->InvokeEvent("ErrorEvent", const_cast<char*>(oss.str().c_str()));
483 vtkOutputWindowDisplayErrorText(const_cast<char*>(oss.str().c_str()));
485 vtkObject::BreakOnError();
491 int vtkZJFilter::RequestData(
492 vtkInformation* request, vtkInformationVector** inputVector, vtkInformationVector* outputVector)
494 // std::cerr << "########################################## vtkZJFilter::RequestData
495 // ##########################################" << std::endl;
498 vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
499 vtkInformation *outInfo(outputVector->GetInformationObject(0));
500 vtkUnstructuredGrid* usgIn(nullptr);
501 ExtractInfo(inputVector[0], usgIn);
502 std::string meshName;
503 std::vector<Grp> groups;
504 std::vector<Fam> fams;
505 if (IsInformationOK(inputInfo))
507 vtkMutableDirectedGraph* famGrpGraph(
508 vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(GetMEDReaderMetaDataIfAny())));
509 LoadFamGrpMapInfo(famGrpGraph, meshName, groups, fams);
511 std::vector<std::string> conds(FindConds(groups));
512 vtkUnstructuredGrid* output(
513 vtkUnstructuredGrid::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
514 vtkNew<vtkMultiBlockDataSet> mb2;
515 std::size_t iblock(0);
516 for (std::vector<std::string>::const_iterator it = conds.begin(); it != conds.end();
519 std::vector<std::string> eports2;
520 std::vector<std::string> eports(FindEports(*it, groups, eports2));
521 std::vector<int> ids(DeduceIdsFrom(eports2));
522 vtkNew<vtkMultiBlockDataSet> mb;
523 std::size_t sz(eports.size());
524 for (std::size_t i = 0; i < sz; i++)
526 std::set<int> zeIds(FamiliesIdsFromGrp(groups, fams, eports[i]));
528 vtkSmartPointer<vtkDataSet> ds(FilterFamilies(this, usgIn, zeIds));
530 vtkNew<vtkLongArray> arr;
531 arr->SetName((*it).c_str());
532 arr->SetNumberOfComponents(1);
533 int nbTuples(ds->GetNumberOfCells());
534 arr->SetNumberOfTuples(nbTuples);
535 long* pt(arr->GetPointer(0));
536 std::fill(pt, pt + nbTuples, ids[i]);
537 ds->GetCellData()->AddArray(arr);
539 this->UpdateProgress(double(i) / double(sz));
542 vtkNew<vtkCompositeDataToUnstructuredGridFilter> cd;
543 cd->SetInputData(mb);
544 cd->SetMergePoints(0);
546 mb2->SetBlock(iblock, cd->GetOutput());
549 vtkNew<vtkCompositeDataToUnstructuredGridFilter> cd2;
550 cd2->SetInputData(mb2);
551 cd2->SetMergePoints(0);
553 output->ShallowCopy(cd2->GetOutput());
556 catch (INTERP_KERNEL::Exception& e)
558 std::ostringstream oss;
559 oss << "Exception has been thrown in vtkZJFilter::RequestInformation : " << e.what()
561 if (this->HasObserver("ErrorEvent"))
562 this->InvokeEvent("ErrorEvent", const_cast<char*>(oss.str().c_str()));
564 vtkOutputWindowDisplayErrorText(const_cast<char*>(oss.str().c_str()));
565 vtkObject::BreakOnError();
571 void vtkZJFilter::PrintSelf(ostream& os, vtkIndent indent)
573 this->Superclass::PrintSelf(os, indent);