]> SALOME platform Git repositories - modules/paravis.git/blob - src/Plugins/ParaMEDCorba/vtkParaMEDCorbaSource.cxx
Salome HOME
Add missing MEDFile dependency in DevelopedSurface filter CMakeLists.txt
[modules/paravis.git] / src / Plugins / ParaMEDCorba / vtkParaMEDCorbaSource.cxx
1 // Copyright (C) 2010-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
20 #include "vtkParaMEDCorbaSource.h"
21
22 #include "vtkPoints.h"
23 #include "vtkIntArray.h"
24 #include "vtkCellData.h"
25 #include "vtkCellTypes.h"
26 #include "vtkCharArray.h"
27 #include "vtkPointData.h"
28 #include "vtkDoubleArray.h"
29 #include "vtkMultiBlockDataSet.h"
30 #include "vtkUnstructuredGrid.h"
31 //
32 #include "vtkStreamingDemandDrivenPipeline.h"
33 #include "vtkGenericAttributeCollection.h"
34 #include "vtkInformationVector.h"
35 #include "vtkObjectFactory.h"
36 #include "vtkInformation.h"
37 //
38 #include "string"
39 #include "fstream"
40 #include "algorithm"
41
42 #include "VTKMEDCouplingMeshClient.hxx"
43 #include "VTKMEDCouplingFieldClient.hxx"
44 #include "VTKMEDCouplingMultiFieldsClient.hxx"
45 #include "VTKParaMEDFieldClient.hxx"
46
47 //Work with IOR.
48 #include "ParaMEDCouplingCorbaServant.hh"
49 //
50
51 vtkStandardNewMacro(vtkParaMEDCorbaSource);
52 //vtkCxxRevisionMacro(vtkParaMEDCorbaSource,"$Revision$");
53
54 void *vtkParaMEDCorbaSource::Orb=0;
55
56 vtkParaMEDCorbaSource::vtkParaMEDCorbaSource():mfieldsFetcher(0)
57 {
58   this->MyDataSet=0;
59   if(!Orb)
60   {
61     CORBA::ORB_var *OrbC=new CORBA::ORB_var;
62     int argc=0;
63     *OrbC=CORBA::ORB_init(argc,0);
64     this->Orb=OrbC;
65   }
66   this->SetNumberOfInputPorts(0);
67   this->SetNumberOfOutputPorts(1);
68 }
69
70 vtkParaMEDCorbaSource::~vtkParaMEDCorbaSource()
71 {
72   delete mfieldsFetcher;
73 }
74
75 const char *vtkParaMEDCorbaSource::GetIORCorba()
76 {
77   return &IOR[0];
78 }
79
80 void vtkParaMEDCorbaSource::SetIORCorba(char *ior)
81 {
82   if(!ior)
83     return;
84   if(ior[0]=='\0')
85     return;
86   int length=strlen(ior);
87   IOR.resize(length+1);
88   std::copy(ior,ior+length+1,&IOR[0]);
89   this->Modified();
90 }
91
92 void vtkParaMEDCorbaSource::SetBufferingPolicy(int pol)
93 {
94   BufferingPolicy=pol;
95 }
96
97 int vtkParaMEDCorbaSource::GetBufferingPolicy()
98 {
99   return BufferingPolicy;
100 }
101
102 //int vtkParaMEDCorbaSource::RequestUpdateExtent( vtkInformation* request, vtkInformationVector** inInfo, vtkInformationVector* outInfo )
103 //{
104 //return this->Superclass::RequestUpdateExtent(request,inInfo,outInfo);
105
106 /*vtkParaMEDCorbaDataSet* output = vtkParaMEDCorbaDataSet::SafeDownCast( info->Get( vtkDataObject::DATA_OBJECT() ) );
107   if ( ! output )
108     {
109     output = vtkParaMEDCorbaDataSet::New();
110     output->SetPipelineInformation( info );
111     output->Delete();
112     this->GetOutputPortInformation( 0 )->Set( vtkDataObject::DATA_EXTENT_TYPE(), output->GetExtentType() );
113     }*/
114
115 // return 1;
116 //}
117
118 int vtkParaMEDCorbaSource::ProcessRequest(vtkInformation* request,
119     vtkInformationVector** inputVector,
120     vtkInformationVector* outputVector)
121 {
122   // generate the data
123   if(request->Has(vtkDemandDrivenPipeline::REQUEST_DATA()))
124   {
125     return this->RequestData(request, inputVector, outputVector);
126   }
127   if(request->Has(vtkDemandDrivenPipeline::REQUEST_INFORMATION()))
128   {
129     return this->RequestInformation(request, inputVector, outputVector);
130   }
131   return this->Superclass::ProcessRequest(request, inputVector, outputVector);
132 }
133
134 int vtkParaMEDCorbaSource::FillOutputPortInformation(int vtkNotUsed(port), vtkInformation* info)
135 {
136   info->Set(vtkDataObject::DATA_TYPE_NAME(), "vtkMultiBlockDataSet");
137   return 1;
138 }
139
140 int vtkParaMEDCorbaSource::RequestInformation(vtkInformation* request, vtkInformationVector** inInfo, vtkInformationVector* outInfo)
141 {
142   vtkInformation* myInfo=outInfo->GetInformationObject(0);
143   //myInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkUnstructuredGrid");
144   if(!IOR.empty())
145   {
146     //myInfo->Remove(vtkDataObject::DATA_TYPE_NAME());
147     //myInfo->Remove(PORT_REQUIREMENTS_FILLED());
148     //myInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkUnstructuredGrid");
149     //myInfo->Set(PORT_REQUIREMENTS_FILLED(),1);
150     //vtkUnstructuredGrid *tony=vtkUnstructuredGrid::New();
151     //tony->SetInformation(myInfo);
152     //myInfo->Set(vtkDataObject::DATA_OBJECT(),tony);
153     //
154     try {
155       CORBA::ORB_var *OrbC=(CORBA::ORB_var *)this->Orb;
156       CORBA::Object_var obj=(*OrbC)->string_to_object(&IOR[0]);
157       //
158       Engines::MPIObject_ptr objPara=Engines::MPIObject::_narrow(obj);
159       if(CORBA::is_nil(objPara))
160       {//sequential
161         this->TotalNumberOfPieces=1;
162         SALOME_MED::MEDCouplingMultiFieldsCorbaInterface_var multiPtr=SALOME_MED::MEDCouplingMultiFieldsCorbaInterface::_narrow(obj);
163         if(!CORBA::is_nil(multiPtr))
164         {//Request for multiFields
165           delete mfieldsFetcher;
166           mfieldsFetcher=new ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher(BufferingPolicy,multiPtr);
167           std::vector<double> tsteps=mfieldsFetcher->getTimeStepsForPV();
168           double timeRange[2];
169           timeRange[0]=tsteps.front();
170           timeRange[1]=tsteps.back();
171           myInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&tsteps[0],tsteps.size());
172           myInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
173         }
174       }
175       else
176       {
177         Engines::IORTab *iorTab=objPara->tior();
178         this->TotalNumberOfPieces=iorTab->length();
179         delete iorTab;
180         CORBA::release(objPara);
181       }
182       myInfo->Set(CAN_HANDLE_PIECE_REQUEST(), 1);
183     }
184     catch(CORBA::Exception&) {
185       vtkErrorMacro("On fetching object error occurs");
186     }
187   }
188   return 1;
189 }
190
191 int vtkParaMEDCorbaSource::RequestData(vtkInformation* request, vtkInformationVector** inInfo, vtkInformationVector* outputVector)
192 {
193   vtkInformation *outInfo=outputVector->GetInformationObject(0);
194   //
195   this->UpdatePiece = vtkStreamingDemandDrivenPipeline::GetUpdatePiece(outInfo);
196   this->NumberOfPieces = vtkStreamingDemandDrivenPipeline::GetUpdateNumberOfPieces(outInfo);
197   this->GhostLevel = vtkStreamingDemandDrivenPipeline::GetUpdateGhostLevel(outInfo);
198   this->StartPiece=((this->UpdatePiece*this->TotalNumberOfPieces)/this->NumberOfPieces);
199   this->EndPiece=(((this->UpdatePiece+1)*this->TotalNumberOfPieces)/this->NumberOfPieces);
200   vtkMultiBlockDataSet *ret0=vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
201   double reqTS = 0;
202   if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
203     reqTS = outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
204   try {
205     //Client request on ORB.
206     CORBA::ORB_var *OrbC=(CORBA::ORB_var *)this->Orb;
207     CORBA::Object_var obj=(*OrbC)->string_to_object(&IOR[0]);
208     //
209     Engines::MPIObject_var objPara=Engines::MPIObject::_narrow(obj);
210     if(CORBA::is_nil(objPara))
211     {//sequential
212       SALOME_MED::MEDCouplingMeshCorbaInterface_var meshPtr=SALOME_MED::MEDCouplingMeshCorbaInterface::_narrow(obj);
213       if(!CORBA::is_nil(meshPtr))
214       {
215         bool dummy;//bug VTK
216         vtkDataSet *ret=ParaMEDMEM2VTK::BuildFromMEDCouplingMeshInstance(meshPtr,dummy);//bug VTK
217         if(!ret)
218           return 0;
219         ret0->SetBlock(0,ret);
220         ret->Delete();
221         return 1;
222       }
223       SALOME_MED::MEDCouplingFieldDoubleCorbaInterface_var fieldPtr=SALOME_MED::MEDCouplingFieldDoubleCorbaInterface::_narrow(obj);
224       if(!CORBA::is_nil(fieldPtr))
225       {
226         std::vector<double> ret2;
227         vtkDataSet *ret=ParaMEDMEM2VTK::BuildFullyFilledFromMEDCouplingFieldDoubleInstance(fieldPtr,ret2);
228         if(!ret)
229         {
230           vtkErrorMacro("On single field CORBA fetching an error occurs !");
231           return 0;
232         }
233         ret0->SetBlock(0,ret);
234         ret->Delete();
235         //
236         double timeRange[2];
237         timeRange[0]=ret2[0];
238         timeRange[1]=ret2[0];
239         outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&ret2[0],1);
240         outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
241         ret0->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),ret2[0]);
242         return 1;
243       }
244       SALOME_MED::MEDCouplingMultiFieldsCorbaInterface_var multiPtr=SALOME_MED::MEDCouplingMultiFieldsCorbaInterface::_narrow(obj);
245       if(!CORBA::is_nil(multiPtr))
246       {
247         vtkDataSet *ret=mfieldsFetcher->buildDataSetOnTime(reqTS);
248         if(!ret)
249         {
250           vtkErrorMacro("On multi fields CORBA fetching an error occurs !");
251           return 0;
252         }
253         ret0->SetBlock(0,ret);
254         ret->Delete();
255         ret0->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),reqTS);
256         return 1;
257       }
258       vtkErrorMacro("Unrecognized sequential CORBA reference !");
259       return 0;
260     }
261     else
262     {
263       SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface_var paraFieldCorba=SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface::_narrow(obj);
264       if(!CORBA::is_nil(paraFieldCorba))
265       {
266         ParaMEDMEM2VTK::FillMEDCouplingParaFieldDoubleInstanceFrom(paraFieldCorba,this->StartPiece,this->EndPiece,ret0);
267         return 1;
268       }
269       vtkErrorMacro("Unrecognized parallel CORBA reference !");
270       return 0;
271     }
272   }
273   catch(CORBA::Exception&) {
274     vtkErrorMacro("On fetching object error occurs");
275   }
276 }
277
278 void vtkParaMEDCorbaSource::PrintSelf(ostream& os, vtkIndent indent)
279 {
280   this->Superclass::PrintSelf( os, indent );
281   os << "Data: " << this->MyDataSet << "\n";
282 }
283