Salome HOME
Porting to ParaView 5.8
[modules/paravis.git] / src / Plugins / MEDWriter / plugin / MEDWriterIO / vtkMEDWriter.cxx
1 // Copyright (C) 2016  CEA/DEN, EDF R&D
2 //
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.
7 //
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.
12 //
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
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19 // Author : Anthony Geay (EDF R&D)
20
21 #include "vtkMEDWriter.h"
22 #include "VTKToMEDMem.h"
23
24 #include <vtkAdjacentVertexIterator.h>
25 #include <vtkAlgorithmOutput.h>
26 #include <vtkCellArray.h>
27 #include <vtkCellData.h>
28 #include <vtkCharArray.h>
29 #include <vtkDataArraySelection.h>
30 #include <vtkDataObjectTreeIterator.h>
31 #include <vtkDataSet.h>
32 #include <vtkDataSetAttributes.h>
33 #include <vtkDemandDrivenPipeline.h>
34 #include <vtkDoubleArray.h>
35 #include <vtkExecutive.h>
36 #include <vtkFloatArray.h>
37 #include <vtkInEdgeIterator.h>
38 #include <vtkInformation.h>
39 #include <vtkInformationDataObjectKey.h>
40 #include <vtkInformationDataObjectMetaDataKey.h>
41 #include <vtkInformationStringKey.h>
42 #include <vtkInformationVector.h>
43 #include <vtkIntArray.h>
44 #include <vtkLongArray.h>
45 #include <vtkMultiBlockDataSet.h>
46 #include <vtkMutableDirectedGraph.h>
47 #include <vtkObjectFactory.h>
48 #include <vtkPointData.h>
49 #include <vtkPolyData.h>
50 #include <vtkRectilinearGrid.h>
51 #include <vtkStreamingDemandDrivenPipeline.h>
52 #include <vtkStringArray.h>
53 #include <vtkTimeStamp.h>
54 #include <vtkUnsignedCharArray.h>
55 #include <vtkUnstructuredGrid.h>
56 #include <vtkVariantArray.h>
57 #include <vtkWarpScalar.h>
58
59 #include "MEDCouplingFieldDouble.hxx"
60 #include "MEDCouplingFieldFloat.hxx"
61 #include "MEDCouplingFieldInt.hxx"
62 #include "MEDCouplingMemArray.hxx"
63 #include "MEDCouplingRefCountObject.hxx"
64 #include "MEDFileData.hxx"
65 #include "MEDFileField.hxx"
66 #include "MEDFileMesh.hxx"
67 #include "MEDLoaderTraits.hxx"
68
69 #include <deque>
70 #include <map>
71 #include <sstream>
72
73 vtkStandardNewMacro(vtkMEDWriter);
74
75 using MEDCoupling::MCAuto;
76 using MEDCoupling::MEDFileData;
77 using VTKToMEDMem::Fam;
78 using VTKToMEDMem::Grp;
79
80 vtkInformationDataObjectMetaDataKey* GetMEDReaderMetaDataIfAny()
81 {
82   static const char ZE_KEY[] = "vtkMEDReader::META_DATA";
83   MEDCoupling::GlobalDict* gd(MEDCoupling::GlobalDict::GetInstance());
84   if (!gd->hasKey(ZE_KEY))
85     return 0;
86   std::string ptSt(gd->value(ZE_KEY));
87   void* pt(0);
88   std::istringstream iss(ptSt);
89   iss >> pt;
90   return reinterpret_cast<vtkInformationDataObjectMetaDataKey*>(pt);
91 }
92
93 bool IsInformationOK(vtkInformation* info)
94 {
95   vtkInformationDataObjectMetaDataKey* key(GetMEDReaderMetaDataIfAny());
96   if (!key)
97     return false;
98   // Check the information contain meta data key
99   if (!info->Has(key))
100     return false;
101   // Recover Meta Data
102   vtkMutableDirectedGraph* sil(vtkMutableDirectedGraph::SafeDownCast(info->Get(key)));
103   if (!sil)
104     return false;
105   int idNames(0);
106   vtkAbstractArray* verticesNames(sil->GetVertexData()->GetAbstractArray("Names", idNames));
107   vtkStringArray* verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
108   if (!verticesNames2)
109     return false;
110   for (int i = 0; i < verticesNames2->GetNumberOfValues(); i++)
111   {
112     vtkStdString& st(verticesNames2->GetValue(i));
113     if (st == "MeshesFamsGrps")
114       return true;
115   }
116   return false;
117 }
118
119 void LoadFamGrpMapInfo(vtkMutableDirectedGraph* sil, std::string& meshName,
120   std::vector<Grp>& groups, std::vector<Fam>& fams)
121 {
122   if (!sil)
123     throw MZCException("LoadFamGrpMapInfo : internal error !");
124   int idNames(0);
125   vtkAbstractArray* verticesNames(sil->GetVertexData()->GetAbstractArray("Names", idNames));
126   vtkStringArray* verticesNames2(vtkStringArray::SafeDownCast(verticesNames));
127   vtkIdType id0;
128   bool found(false);
129   for (int i = 0; i < verticesNames2->GetNumberOfValues(); i++)
130   {
131     vtkStdString& st(verticesNames2->GetValue(i));
132     if (st == "MeshesFamsGrps")
133     {
134       id0 = i;
135       found = true;
136     }
137   }
138   if (!found)
139     throw INTERP_KERNEL::Exception(
140       "There is an internal error ! The tree on server side has not the expected look !");
141   vtkAdjacentVertexIterator* it0(vtkAdjacentVertexIterator::New());
142   sil->GetAdjacentVertices(id0, it0);
143   int kk(0), ll(0);
144   while (it0->HasNext())
145   {
146     vtkIdType id1(it0->Next());
147     std::string mName((const char*)verticesNames2->GetValue(id1));
148     meshName = mName;
149     vtkAdjacentVertexIterator* it1(vtkAdjacentVertexIterator::New());
150     sil->GetAdjacentVertices(id1, it1);
151     vtkIdType idZeGrps(it1->Next()); // zeGroups
152     vtkAdjacentVertexIterator* itGrps(vtkAdjacentVertexIterator::New());
153     sil->GetAdjacentVertices(idZeGrps, itGrps);
154     while (itGrps->HasNext())
155     {
156       vtkIdType idg(itGrps->Next());
157       Grp grp((const char*)verticesNames2->GetValue(idg));
158       vtkAdjacentVertexIterator* itGrps2(vtkAdjacentVertexIterator::New());
159       sil->GetAdjacentVertices(idg, itGrps2);
160       std::vector<std::string> famsOnGroup;
161       while (itGrps2->HasNext())
162       {
163         vtkIdType idgf(itGrps2->Next());
164         famsOnGroup.push_back(std::string((const char*)verticesNames2->GetValue(idgf)));
165       }
166       grp.setFamilies(famsOnGroup);
167       itGrps2->Delete();
168       groups.push_back(grp);
169     }
170     itGrps->Delete();
171     vtkIdType idZeFams(it1->Next()); // zeFams
172     it1->Delete();
173     vtkAdjacentVertexIterator* itFams(vtkAdjacentVertexIterator::New());
174     sil->GetAdjacentVertices(idZeFams, itFams);
175     while (itFams->HasNext())
176     {
177       vtkIdType idf(itFams->Next());
178       Fam fam((const char*)verticesNames2->GetValue(idf));
179       fams.push_back(fam);
180     }
181     itFams->Delete();
182   }
183   it0->Delete();
184 }
185
186 vtkMEDWriter::vtkMEDWriter()
187   : WriteAllTimeSteps(0)
188   , NumberOfTimeSteps(0)
189   , CurrentTimeIndex(0)
190   , FileName(0)
191 {
192 }
193
194 vtkMEDWriter::~vtkMEDWriter() {}
195
196 int vtkMEDWriter::RequestInformation(
197   vtkInformation* request, vtkInformationVector** inputVector, vtkInformationVector* outputVector)
198 {
199   // std::cerr << "########################################## vtkMEDWriter::RequestInformation
200   // ########################################## " << (const char *) this->FileName << std::endl;
201   vtkInformation* inInfo(inputVector[0]->GetInformationObject(0));
202   if (inInfo->Has(vtkStreamingDemandDrivenPipeline::TIME_STEPS()))
203     this->NumberOfTimeSteps = inInfo->Length(vtkStreamingDemandDrivenPipeline::TIME_STEPS());
204   else
205     this->NumberOfTimeSteps = 0;
206   return 1;
207 }
208
209 int vtkMEDWriter::RequestUpdateExtent(vtkInformation* vtkNotUsed(request),
210   vtkInformationVector** inputVector, vtkInformationVector* vtkNotUsed(outputVector))
211 {
212   double* inTimes(
213     inputVector[0]->GetInformationObject(0)->Get(vtkStreamingDemandDrivenPipeline::TIME_STEPS()));
214   if (inTimes && this->WriteAllTimeSteps)
215   {
216     double timeReq(inTimes[this->CurrentTimeIndex]);
217     inputVector[0]->GetInformationObject(0)->Set(
218       vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP(), timeReq);
219   }
220   return 1;
221 }
222
223 void ExceptionDisplayer(vtkMEDWriter* self, const std::string& fileName, std::exception& e)
224 {
225   std::ostringstream oss;
226   oss << "Exception has been thrown in vtkMEDWriter::RequestData : During writing of \"" << fileName
227       << "\", the following exception has been thrown : " << e.what() << std::endl;
228   if (self->HasObserver("ErrorEvent"))
229     self->InvokeEvent("ErrorEvent", const_cast<char*>(oss.str().c_str()));
230   else
231     vtkOutputWindowDisplayErrorText(const_cast<char*>(oss.str().c_str()));
232   vtkObject::BreakOnError();
233 }
234
235 int vtkMEDWriter::RequestData(
236   vtkInformation* request, vtkInformationVector** inputVector, vtkInformationVector* outputVector)
237 {
238   // std::cerr << "########################################## vtkMEDWriter::RequestData
239   // ########################################## " << (const char *) this->FileName << std::endl;
240   try
241   {
242     bool writeAll((this->WriteAllTimeSteps != 0 && this->NumberOfTimeSteps > 0));
243     if (writeAll)
244     {
245       if (this->CurrentTimeIndex == 0)
246         request->Set(vtkStreamingDemandDrivenPipeline::CONTINUE_EXECUTING(), 1);
247     }
248     else
249     {
250       request->Remove(vtkStreamingDemandDrivenPipeline::CONTINUE_EXECUTING());
251       this->CurrentTimeIndex = 0;
252     }
253     //////////////
254     vtkInformation* inputInfo(inputVector[0]->GetInformationObject(0));
255     std::string meshName;
256     std::vector<Grp> groups;
257     std::vector<Fam> fams;
258     if (IsInformationOK(inputInfo))
259     {
260       vtkMutableDirectedGraph* famGrpGraph(
261         vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(GetMEDReaderMetaDataIfAny())));
262       LoadFamGrpMapInfo(famGrpGraph, meshName, groups, fams);
263     }
264     vtkInformation* outInfo(outputVector->GetInformationObject(0));
265     vtkDataObject* input(vtkDataObject::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
266     if (!input)
267       throw MZCException(
268         "Not recognized data object in input of the MEDWriter ! Maybe not implemented yet !");
269     double timeStep;
270     {
271       vtkInformation* inInfo(inputVector[0]->GetInformationObject(0));
272       vtkDataObject* input(vtkDataObject::GetData(inInfo));
273       timeStep = input->GetInformation()->Get(vtkDataObject::DATA_TIME_STEP());
274     }
275     ////////////
276     MCAuto<MEDFileData> mfd(MEDFileData::New());
277     WriteMEDFileFromVTKGDS(mfd, input, timeStep, this->CurrentTimeIndex);
278     PutFamGrpInfoIfAny(mfd, meshName, groups, fams);
279     if (this->WriteAllTimeSteps == 0 ||
280       (this->WriteAllTimeSteps != 0 && this->CurrentTimeIndex == 0))
281       mfd->write(this->FileName, 2);
282     else
283     {
284       mfd->getFields()->write(this->FileName, 0);
285     }
286     outInfo->Set(vtkDataObject::DATA_OBJECT(), input);
287     ///////////
288     if (writeAll)
289     {
290       this->CurrentTimeIndex++;
291       if (this->CurrentTimeIndex >= this->NumberOfTimeSteps)
292       {
293         request->Remove(vtkStreamingDemandDrivenPipeline::CONTINUE_EXECUTING());
294         this->CurrentTimeIndex = 0;
295       }
296     }
297   }
298   catch (INTERP_KERNEL::Exception& e)
299   {
300     ExceptionDisplayer(this, (const char*)this->FileName, e);
301     return 0;
302   }
303   catch (MZCException& e)
304   {
305     ExceptionDisplayer(this, (const char*)this->FileName, e);
306     return 0;
307   }
308   return 1;
309 }
310
311 void vtkMEDWriter::PrintSelf(ostream& os, vtkIndent indent)
312 {
313   this->Superclass::PrintSelf(os, indent);
314 }
315
316 int vtkMEDWriter::Write()
317 {
318   this->Update();
319   return 0;
320 }