Salome HOME
Merge branch 'omu/insitu'
[modules/paravis.git] / src / Insitu / VisualizationComponent / VisualizationComponent.cxx
1
2 #include "VisualizationComponent.hxx"
3 #include <string>
4 #include <unistd.h>
5 #include <signal.h>
6 #include <SALOME_NamingService.hxx>
7 #include <Utils_SALOME_Exception.hxx>
8 #include "Utils_CorbaException.hxx"
9 #include <pthread.h>
10 #include <execinfo.h>
11
12 typedef struct
13 {
14   bool exception;
15   std::string msg;
16 } exception_st;
17
18 //DEFS
19
20 #include "visu.hxx"
21 #include "MEDCouplingFieldDoubleClient.hxx"
22 #include "ParaMEDCouplingFieldDoubleServant.hxx"
23 #include "omniORB4/poa.h"
24
25 //ENDDEF
26
27
28 using namespace std;
29
30 //! Constructor for component "VisualizationComponent" instance
31 /*!
32  *
33  */
34 VisualizationComponent_i::VisualizationComponent_i(CORBA::ORB_ptr orb,
35                      PortableServer::POA_ptr poa,
36                      PortableServer::ObjectId * contId,
37                      const char *instanceName,
38                      const char *interfaceName,
39                      bool regist)
40           : Engines_Component_i(orb, poa, contId, instanceName, interfaceName,
41                                 false, regist)
42 {
43   _thisObj = this ;
44   _id = _poa->activate_object(_thisObj);
45 }
46
47 VisualizationComponent_i::VisualizationComponent_i(CORBA::ORB_ptr orb,
48                      PortableServer::POA_ptr poa,
49                      Engines::Container_ptr container,
50                      const char *instanceName,
51                      const char *interfaceName,
52                      bool regist)
53           : Engines_Component_i(orb, poa, container, instanceName, interfaceName,
54                                 false, regist)
55 {
56   _thisObj = this ;
57   _id = _poa->activate_object(_thisObj);
58 }
59
60 //! Destructor for component "VisualizationComponent" instance
61 VisualizationComponent_i::~VisualizationComponent_i()
62 {
63 }
64
65
66 void * th_Visualize(void * s)
67 {
68   std::ostringstream msg;
69   exception_st *est = new exception_st;
70   est->exception = false;
71   
72   thread_Visualize_struct *st = (thread_Visualize_struct *)s;
73   
74   try
75   {
76     
77     PARAVIS_ORB::VisualizationComponent_var compo = PARAVIS_ORB::VisualizationComponent::_narrow((*(st->tior))[st->ip]);
78     compo->Visualize(st->field,st->path_python_file);
79   }
80   catch(const SALOME::SALOME_Exception &ex)
81   {
82     est->exception = true;
83     est->msg = ex.details.text;
84   }
85   catch(const CORBA::Exception &ex)
86   {
87     est->exception = true;
88     msg << "CORBA::Exception: " << ex;
89     est->msg = msg.str();
90   }
91   
92   delete st;
93   return ((void*)est);
94 }
95
96 void VisualizationComponent_i::Visualize(SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface_ptr field,const char* path_python_file)
97 {
98   beginService("VisualizationComponent_i::Visualize");
99   void *ret_th;
100   pthread_t *th;
101   exception_st *est;
102
103   try
104     {
105       // Run the service in every MPI process
106       if(_numproc == 0)
107       {
108         th = new pthread_t[_nbproc];
109         for(int ip=1;ip<_nbproc;ip++)
110         {
111             thread_Visualize_struct *st = new thread_Visualize_struct;
112             st->ip = ip;
113             st->tior = _tior;
114             st->field = field;
115 st->path_python_file = path_python_file;
116             pthread_create(&(th[ip]),NULL,th_Visualize,(void*)st);
117         }
118       }
119       
120 //BODY
121
122 const MEDCoupling::MEDCouplingFieldDouble * local_field(NULL);
123 int nb_fields = field->tior()->length();
124 if(nb_fields == _nbproc)
125 {
126   SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface_ptr local_corba_field = 
127                SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface::_narrow((*(field->tior()))[_numproc]);
128
129
130   PortableServer::ServantBase *ret;
131   try {
132     ret=PortableServer::POA::_the_root_poa()->reference_to_servant(local_corba_field);
133     std::cerr << "Servant succeeded!" << std::endl;
134     MEDCoupling::ParaMEDCouplingFieldDoubleServant* servant_field = 
135                dynamic_cast<MEDCoupling::ParaMEDCouplingFieldDoubleServant*>(ret);
136     if(servant_field != NULL)
137     {
138       std::cerr << "In-situ configuration!" << std::endl;
139       // same container, same mpi proc, use the pointer directly.
140       local_field = servant_field->getPointer();
141       ret->_remove_ref();
142     }
143   }
144   catch(...){
145     // different container - need to make a copy of the field.
146     ret = NULL;
147     std::cerr << "Co-processing configuration!" << std::endl;
148     local_field = MEDCoupling::MEDCouplingFieldDoubleClient::New(local_corba_field);
149   }
150 }
151 else if(nb_fields < _nbproc)
152 {
153   if(_numproc < nb_fields)
154   {
155     SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface_ptr local_corba_field = 
156                SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface::_narrow((*(field->tior()))[_numproc]);
157     local_field = MEDCoupling::MEDCouplingFieldDoubleClient::New(local_corba_field);
158   }
159 }
160 else //nb_fields > _nbproc
161 {
162   int q = nb_fields / _nbproc; // int division
163   int r = nb_fields - q * _nbproc;
164   int start, end;
165   
166   if(_numproc < r)
167   {
168     // get one more field to process
169     start = _numproc * (q + 1);
170     end = start + q + 1;
171   }
172   else
173   {
174     start = r * (q + 1) + (_numproc - r) * q;
175     end = start + q;
176   }
177   
178   std::cerr << "Proc n° " << _numproc << ". Merge fields from " << start << " to " << end << std::endl;
179   
180   std::vector<const MEDCoupling::MEDCouplingFieldDouble *> fieldsToProcess;
181   for(int i = start; i < end; i++)
182   {
183     SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface_ptr local_corba_field = 
184                SALOME_MED::ParaMEDCouplingFieldDoubleCorbaInterface::_narrow((*(field->tior()))[i]);
185     fieldsToProcess.push_back(MEDCoupling::MEDCouplingFieldDoubleClient::New(local_corba_field));
186   }
187   
188   local_field = MEDCoupling::MEDCouplingFieldDouble::MergeFields(fieldsToProcess);
189 }
190
191 Visualization v;
192 v.run(const_cast<MEDCoupling::MEDCouplingFieldDouble*>(local_field), path_python_file);
193
194 //ENDBODY
195       if(_numproc == 0)
196       {
197         for(int ip=1;ip<_nbproc;ip++)
198         {
199           pthread_join(th[ip],&ret_th);
200           est = (exception_st*)ret_th;
201           if(est->exception)
202           {
203               std::ostringstream msg;
204               msg << "[" << ip << "] " << est->msg;
205               delete est;
206               delete[] th;
207               THROW_SALOME_CORBA_EXCEPTION(msg.str().c_str(),SALOME::INTERNAL_ERROR);
208           }
209           delete est;
210         }
211         delete[] th;
212       }
213     }
214   catch ( const SALOME_Exception & ex)
215     {
216       THROW_SALOME_CORBA_EXCEPTION(CORBA::string_dup(ex.what()), SALOME::INTERNAL_ERROR);
217     }
218   catch ( const SALOME::SALOME_Exception & ex)
219     {
220       throw;
221     }
222   catch ( const std::exception& ex)
223     {
224       THROW_SALOME_CORBA_EXCEPTION(CORBA::string_dup(ex.what()), SALOME::INTERNAL_ERROR);
225     }
226   catch (...)
227     {
228       THROW_SALOME_CORBA_EXCEPTION("unknown exception", SALOME::INTERNAL_ERROR);
229     }
230   endService("VisualizationComponent_i::Visualize");
231 }
232
233
234
235 extern "C"
236 {
237   PortableServer::ObjectId * VisualizationComponentEngine_factory( CORBA::ORB_ptr orb,
238                                                     PortableServer::POA_ptr poa,
239                                                     PortableServer::ObjectId * contId,
240                                                     const char *instanceName,
241                                                     const char *interfaceName)
242   {
243     MESSAGE("PortableServer::ObjectId * VisualizationComponentEngine_factory()");
244     int is_mpi_container;
245     bool regist;
246     int numproc;
247     
248     MPI_Initialized(&is_mpi_container);
249     if (!is_mpi_container)
250     {
251       int argc = 0;
252       char ** argv = NULL;
253       MPI_Init(&argc, &argv);
254     }
255     
256     MPI_Comm_rank( MPI_COMM_WORLD, &numproc );
257     regist = ( numproc == 0 );
258     VisualizationComponent_i * myEngine = new VisualizationComponent_i(orb, poa, contId, instanceName, interfaceName, regist);
259     return myEngine->getId() ;
260   }
261 }