1 //
2 // Created by Sergei Shudler on 2020-08-21.
3 //
4 
5 #include <vector>
6 #include <sstream>
7 #include <iomanip>
8 #include <iostream>
9 #include <algorithm>
10 #include <fstream>
11 #include <string>
12 
13 #include <conduit.hpp>
14 #include <conduit_relay.hpp>
15 #include <conduit_blueprint.hpp>
16 #include <ascent_data_object.hpp>
17 #include <ascent_logging.hpp>
18 #include <ascent_runtime_param_check.hpp>
19 #include <flow_workspace.hpp>
20 
21 #ifdef ASCENT_VTKM_ENABLED
22 #include <vtkh/vtkh.hpp>
23 #include <vtkh/DataSet.hpp>
24 #include <vtkh/filters/MarchingCubes.hpp>
25 #include <vtkh/rendering/RayTracer.hpp>
26 #include <vtkh/rendering/Scene.hpp>
27 #include <vtkm/rendering/Canvas.h>
28 
29 #include <vtkm/rendering/CanvasRayTracer.h>
30 #include <vtkm/rendering/MapperRayTracer.h>
31 
32 #include <ascent_vtkh_data_adapter.hpp>
33 #endif
34 
35 #ifdef ASCENT_MPI_ENABLED
36 #include <mpi.h>
37 #else
38 #include <mpidummy.h>
39 #define _NOMPI
40 #endif
41 
42 #include "BabelFlow/TypeDefinitions.h"
43 #include "ascent_runtime_babelflow_filters.hpp"
44 #include "ascent_runtime_babelflow_comp_utils.hpp"
45 
46 
47 // #define BFLOW_ISO_DEBUG
48 
49 #define DEF_IMG_WIDTH   1024
50 #define DEF_IMG_HEIGHT  1024
51 
52 
53 namespace BabelFlow
54 {
55 extern int relay_message(std::vector<Payload>& inputs, std::vector<Payload>& outputs, TaskId task);
56 }
57 
58 //-----------------------------------------------------------------------------
59 // -- begin ascent:: --
60 //-----------------------------------------------------------------------------
61 namespace ascent
62 {
63 
64 //-----------------------------------------------------------------------------
65 // -- begin bflow_comp:: --
66 //-----------------------------------------------------------------------------
67 namespace bflow_iso
68 {
69 
70 struct BoundsData
71 {
72   double m_rangeVec[6];   // Pairs of min, max for x, y, z axes
73 
operator =ascent::bflow_iso::BoundsData74   BoundsData& operator=( const vtkm::Bounds& bnd )
75   {
76     m_rangeVec[0] = bnd.X.Min;
77     m_rangeVec[1] = bnd.X.Max;
78     m_rangeVec[2] = bnd.Y.Min;
79     m_rangeVec[3] = bnd.Y.Max;
80     m_rangeVec[4] = bnd.Z.Min;
81     m_rangeVec[5] = bnd.Z.Max;
82 
83     return *this;
84   }
85 
getVtkmBoundsascent::bflow_iso::BoundsData86   vtkm::Bounds getVtkmBounds()
87   {
88     vtkm::Bounds bnd;
89 
90     bnd.X.Min = m_rangeVec[0];
91     bnd.X.Max = m_rangeVec[1];
92     bnd.Y.Min = m_rangeVec[2];
93     bnd.Y.Max = m_rangeVec[3];
94     bnd.Z.Min = m_rangeVec[4];
95     bnd.Z.Max = m_rangeVec[5];
96 
97     return bnd;
98   }
99 
serializeascent::bflow_iso::BoundsData100   BabelFlow::Payload serialize() const
101   {
102     uint32_t payl_size = sizeof(m_rangeVec);
103     char* out_buffer = new char[payl_size];
104     memcpy( out_buffer, (const char*)m_rangeVec, payl_size );
105 
106     return BabelFlow::Payload( payl_size, out_buffer );
107   }
108 
deserializeascent::bflow_iso::BoundsData109   void deserialize( BabelFlow::Payload payload )
110   {
111     memcpy( (char*)m_rangeVec, payload.buffer(), sizeof(m_rangeVec) );
112   }
113 };
114 
115 
116 struct IsoSurfaceData
117 {
118   vtkh::DataSet*      m_DataSet;
119   std::string         m_FieldName;
120   std::vector<double> m_IsoVals;
121   int                 m_Width;
122   int                 m_Height;
123   std::string         m_FileName;
124   BoundsData          m_bounds;
125 
126   IsoSurfaceData() = default;
127 
128   // IsoSurfaceData( vtkh::DataSet* data_set, const std::string& fld_name, const double* iso_vals, int num_iso_vals )
129   // : m_DataSet( data_set ), m_FieldName( fld_name ), m_IsoVals( iso_vals, iso_vals + num_iso_vals ) {}
130 
serializeascent::bflow_iso::IsoSurfaceData131   BabelFlow::Payload serialize() const
132   {
133     BabelFlow::Payload payl_bounds = m_bounds.serialize();
134 
135     uint32_t iso_vals_buff_size = sizeof(double) * m_IsoVals.size();
136     uint32_t payl_size =
137       sizeof(vtkh::DataSet*) +
138       sizeof(size_t) + m_FieldName.size() +
139       sizeof(size_t) + iso_vals_buff_size +
140       2 * sizeof(int) +
141       sizeof(size_t) + m_FileName.size() +
142       payl_bounds.size();
143 
144     char* out_buffer = new char[payl_size];
145     uint32_t offset = 0;
146     memcpy( out_buffer + offset, (const char*)&m_DataSet, sizeof(vtkh::DataSet*) );
147     offset += sizeof(vtkh::DataSet*);
148     size_t fld_name_size = m_FieldName.size();
149     memcpy( out_buffer + offset, (const char*)&fld_name_size, sizeof(size_t) );
150     offset += sizeof(size_t);
151     memcpy( out_buffer + offset, m_FieldName.c_str(), m_FieldName.size() );
152     offset += m_FieldName.size();
153     size_t num_iso_vals = m_IsoVals.size();
154     memcpy( out_buffer + offset, (const char*)&num_iso_vals, sizeof(size_t) );
155     offset += sizeof(size_t);
156     memcpy( out_buffer + offset, (const char*)m_IsoVals.data(), iso_vals_buff_size );
157     offset += iso_vals_buff_size;
158     memcpy( out_buffer + offset, (const char*)&m_Width, sizeof(int) );
159     offset += sizeof(int);
160     memcpy( out_buffer + offset, (const char*)&m_Height, sizeof(int) );
161     offset += sizeof(int);
162     fld_name_size = m_FileName.size();
163     memcpy( out_buffer + offset, (const char*)&fld_name_size, sizeof(size_t) );
164     offset += sizeof(size_t);
165     memcpy( out_buffer + offset, m_FileName.c_str(), m_FileName.size() );
166     offset += m_FileName.size();
167 
168     memcpy( out_buffer + offset, payl_bounds.buffer(), payl_bounds.size() );
169     payl_bounds.reset();
170 
171     return BabelFlow::Payload( payl_size, out_buffer );
172   }
173 
deserializeascent::bflow_iso::IsoSurfaceData174   void deserialize( BabelFlow::Payload payload )
175   {
176     uint32_t offset = 0;
177     size_t fld_name_len, num_iso_vals, file_name_len;
178 
179     memcpy( (char*)&m_DataSet, payload.buffer() + offset, sizeof(vtkh::DataSet*) );
180     offset += sizeof(vtkh::DataSet*);
181 
182     memcpy( (char*)(&fld_name_len), payload.buffer() + offset, sizeof(size_t) );
183     offset += sizeof(size_t);
184 
185     m_FieldName.assign( payload.buffer() + offset, fld_name_len );
186     offset += m_FieldName.size();
187 
188     memcpy( (char*)&num_iso_vals, payload.buffer() + offset, sizeof(size_t) );
189     offset += sizeof(size_t);
190 
191     uint32_t iso_vals_buff_size = sizeof(double) * num_iso_vals;
192     m_IsoVals.resize( num_iso_vals );
193     memcpy( (char*)m_IsoVals.data(), payload.buffer() + offset, iso_vals_buff_size );
194     offset += iso_vals_buff_size;
195 
196     memcpy( (char*)&m_Width, payload.buffer() + offset, sizeof(int) );
197     offset += sizeof(int);
198 
199     memcpy( (char*)&m_Height, payload.buffer() + offset, sizeof(int) );
200     offset += sizeof(int);
201 
202     memcpy( (char*)&file_name_len, payload.buffer() + offset, sizeof(size_t) );
203     offset += sizeof(size_t);
204 
205     m_FileName.assign( payload.buffer() + offset, file_name_len );
206     offset += m_FileName.size();
207 
208     m_bounds.deserialize( BabelFlow::Payload( payload.size() - offset, payload.buffer() + offset ) );
209   }
210 };
211 
212 
213 
marching_cubes(std::vector<BabelFlow::Payload> & inputs,std::vector<BabelFlow::Payload> & outputs,BabelFlow::TaskId task_id)214 int marching_cubes(std::vector<BabelFlow::Payload>& inputs,
215                    std::vector<BabelFlow::Payload>& outputs,
216                    BabelFlow::TaskId task_id)
217 {
218   assert( inputs.size() == 2 );
219 
220   IsoSurfaceData iso_surf_data;
221   iso_surf_data.deserialize( inputs[0] );
222 
223   BoundsData bounds_data;
224   bounds_data.deserialize( inputs[1] );
225 
226   vtkh::MarchingCubes marcher;
227 
228   marcher.SetInput( iso_surf_data.m_DataSet );
229   marcher.SetField( iso_surf_data.m_FieldName ) ;
230   marcher.SetIsoValues( iso_surf_data.m_IsoVals.data(), iso_surf_data.m_IsoVals.size() );
231   marcher.Update();
232 
233   iso_surf_data.m_DataSet = marcher.GetOutput();
234   iso_surf_data.m_bounds = bounds_data;
235 
236   outputs[0] = iso_surf_data.serialize();
237 
238   inputs[0].reset();
239 
240   return 1;
241 }
242 
243 
vtkm_rendering(std::vector<BabelFlow::Payload> & inputs,std::vector<BabelFlow::Payload> & outputs,BabelFlow::TaskId task_id)244 int vtkm_rendering(std::vector<BabelFlow::Payload>& inputs,
245                    std::vector<BabelFlow::Payload>& outputs,
246                    BabelFlow::TaskId task_id)
247 {
248   assert( inputs.size() == 1 );
249 
250   IsoSurfaceData iso_surf_data;
251   iso_surf_data.deserialize( inputs[0] );
252 
253   vtkh::Render render = vtkh::MakeRender( iso_surf_data.m_Width,
254                                           iso_surf_data.m_Height,
255                                           iso_surf_data.m_bounds.getVtkmBounds(),
256                                           iso_surf_data.m_FileName );
257 
258   vtkh::Scene scene;
259   vtkh::RayTracer renderer;
260 
261   renderer.SetInput( iso_surf_data.m_DataSet );
262   renderer.SetField( iso_surf_data.m_FieldName );
263 
264   scene.AddRenderer( &renderer );
265   scene.AddRender( render );
266 
267   vtkm::Range range;
268   range.Min = iso_surf_data.m_IsoVals.front();
269   range.Max = iso_surf_data.m_IsoVals.back();
270 
271   // {
272   //   auto ranges = iso_surf_data.m_DataSet->GetGlobalRange( iso_surf_data.m_FieldName );
273   //   int num_components = ranges.GetNumberOfValues();
274   //   std::cout << "vtkm_rendering -- range num components = " << num_components << std::endl;
275   //   vtkm::Range global_range = ranges.ReadPortal().Get(0);
276   //   // a min or max may be been set by the user, check to see
277   //   if(range.Min == vtkm::Infinity64())
278   //   {
279   //     range.Min = global_range.Min;
280   //   }
281   //   if(range.Max == vtkm::NegativeInfinity64())
282   //   {
283   //     range.Max = global_range.Max;
284   //   }
285   // }
286 
287 #ifdef BFLOW_ISO_DEBUG
288   {
289     std::cout << "vtkm_rendering --" << std::endl;
290     std::cout << "Num domains: " << iso_surf_data.m_DataSet->GetNumberOfDomains() << std::endl;
291     std::cout << "Color table: " << renderer.GetColorTable().GetName() << std::endl;
292     std::cout << "Range: " << range << std::endl;
293   }
294 #endif
295 
296   render.GetCanvas().Clear();
297 
298   // renderer.SetDoComposite( true );
299   renderer.AddRender ( render );
300 
301   // renderer.Update();
302 
303   vtkm::cont::DataSet data_set;
304   vtkm::Id domain_id;
305   iso_surf_data.m_DataSet->GetDomain( 0, data_set, domain_id );
306   const vtkm::cont::DynamicCellSet &cellset = data_set.GetCellSet();
307   const vtkm::cont::Field &field = data_set.GetField( renderer.GetFieldName() );
308   const vtkm::cont::CoordinateSystem &coords = data_set.GetCoordinateSystem();
309 
310   auto mapper = std::make_shared<vtkm::rendering::MapperRayTracer>();
311 
312   mapper->SetActiveColorTable( renderer.GetColorTable() );
313 
314   auto& canvas = render.GetCanvas();
315   auto& cam = render.GetCamera();
316   mapper->SetCanvas( &canvas );
317   mapper->RenderCells( cellset, coords, field, renderer.GetColorTable(), cam, range );
318 
319 #ifdef BFLOW_ISO_DEBUG
320   {
321     std::cout << "vtkm_rendering -- finished" << std::endl;
322   }
323 #endif
324 
325   auto color_portal = render.GetCanvas().GetColorBuffer().ReadPortal();
326   auto depth_portal = render.GetCanvas().GetDepthBuffer().ReadPortal();
327 
328   assert( color_portal.GetNumberOfValues() == iso_surf_data.m_Width*iso_surf_data.m_Height );
329 
330   bflow_comp::ImageData input_img;
331   input_img.image = new bflow_comp::ImageData::PixelType[iso_surf_data.m_Width*iso_surf_data.m_Height*bflow_comp::ImageData::sNUM_CHANNELS];
332   input_img.zbuf = new bflow_comp::ImageData::PixelType[iso_surf_data.m_Width*iso_surf_data.m_Height];
333   input_img.bounds = new uint32_t[4];
334   input_img.rend_bounds = new uint32_t[4];
335   input_img.bounds[0] = input_img.rend_bounds[0] = 0;
336   input_img.bounds[1] = input_img.rend_bounds[1] = iso_surf_data.m_Width - 1;
337   input_img.bounds[2] = input_img.rend_bounds[2] = 0;
338   input_img.bounds[3] = input_img.rend_bounds[3] = iso_surf_data.m_Height - 1;
339 
340   uint32_t img_offset = 0;
341 
342   for( vtkm::Id index = 0; index < color_portal.GetNumberOfValues(); ++index )
343   {
344     vtkm::Vec4f_32 cur_color = color_portal.Get( index );
345     vtkm::Float32 cur_z = depth_portal.Get( index );
346 
347     // input_img.image[img_offset + 0] = (unsigned char)(cur_color[0] * 255.f);
348     // input_img.image[img_offset + 1] = (unsigned char)(cur_color[1] * 255.f);
349     // input_img.image[img_offset + 2] = (unsigned char)(cur_color[2] * 255.f);
350     // input_img.image[img_offset + 3] = (unsigned char)(cur_color[3] * 255.f);
351     // input_img.zbuf[img_offset] = (unsigned char)(cur_z * 255.f);
352 
353     input_img.image[img_offset + 0] = cur_color[0];
354     input_img.image[img_offset + 1] = cur_color[1];
355     input_img.image[img_offset + 2] = cur_color[2];
356     input_img.image[img_offset + 3] = cur_color[3];
357 
358     input_img.zbuf[index] = cur_z;
359 
360     img_offset += bflow_comp::ImageData::sNUM_CHANNELS;
361   }
362 
363   outputs[0] = input_img.serialize();
364   input_img.delBuffers();
365 
366   return 1;
367 }
368 
allreduce(std::vector<BabelFlow::Payload> & inputs,std::vector<BabelFlow::Payload> & outputs,BabelFlow::TaskId task_id)369 int allreduce(std::vector<BabelFlow::Payload>& inputs,
370               std::vector<BabelFlow::Payload>& outputs,
371               BabelFlow::TaskId task_id)
372 {
373   std::vector<BoundsData> bounds_v( inputs.size() );
374 
375   for( uint32_t i = 0; i < inputs.size(); ++i )
376     bounds_v[i].deserialize( inputs[i] );
377 
378   BoundsData bounds_red = bounds_v[0];
379   for( uint32_t i = 1; i < inputs.size(); ++i )
380   {
381     // X axis
382     bounds_red.m_rangeVec[0] = std::min( bounds_red.m_rangeVec[0], bounds_v[i].m_rangeVec[0] );
383     bounds_red.m_rangeVec[1] = std::max( bounds_red.m_rangeVec[1], bounds_v[i].m_rangeVec[1] );
384     // Y axis
385     bounds_red.m_rangeVec[2] = std::min( bounds_red.m_rangeVec[2], bounds_v[i].m_rangeVec[2] );
386     bounds_red.m_rangeVec[3] = std::max( bounds_red.m_rangeVec[3], bounds_v[i].m_rangeVec[3] );
387     // Z axis
388     bounds_red.m_rangeVec[4] = std::min( bounds_red.m_rangeVec[4], bounds_v[i].m_rangeVec[4] );
389     bounds_red.m_rangeVec[5] = std::max( bounds_red.m_rangeVec[5], bounds_v[i].m_rangeVec[5] );
390   }
391 
392 #ifdef BFLOW_ISO_DEBUG
393   {
394     std::cout << "allreduce -- task_id = " << task_id << ", inputs = " << inputs.size() << ", outputs = " << outputs.size() << std::endl;
395   }
396 #endif
397 
398   for( uint32_t i = 0; i < outputs.size(); ++i )
399     outputs[i] = bounds_red.serialize();
400 
401   for( BabelFlow::Payload& payl : inputs )
402     payl.reset();
403 
404   return 1;
405 }
406 
407 class BabelIsoGraph : public bflow_comp::BabelCompRadixK
408 {
409 public:
BabelIsoGraph(const IsoSurfaceData & iso_surf_data,const std::string & img_name,int32_t rank_id,int32_t n_blocks,MPI_Comm mpi_comm,const std::vector<uint32_t> & radix_v)410   BabelIsoGraph( const IsoSurfaceData& iso_surf_data,
411                  const std::string& img_name,
412                  int32_t rank_id,
413                  int32_t n_blocks,
414                  MPI_Comm mpi_comm,
415                  const std::vector<uint32_t>& radix_v )
416   : BabelCompRadixK( bflow_comp::ImageData(), img_name, rank_id, n_blocks, 2, mpi_comm, radix_v ),
417     m_isoSurfData( iso_surf_data ) {}
418 
~BabelIsoGraph()419   virtual ~BabelIsoGraph() {}
420 
Initialize()421   virtual void Initialize() override
422   {
423     InitRadixKGraph();
424     InitGatherGraph();
425 
426     // Iso surf calc graph
427     m_isoCalcTaskGr = BabelFlow::SingleTaskGraph( m_nRanks, 2, 1, 1 );
428     m_isoCalcTaskMp = BabelFlow::ModuloMap( m_nRanks, m_nRanks );
429 
430     // Iso surface rendering graph
431     m_isoRenderTaskGr = BabelFlow::SingleTaskGraph( m_nRanks );
432     m_isoRenderTaskMp = BabelFlow::ModuloMap( m_nRanks, m_nRanks );
433 
434     m_redAllGr = BabelFlow::RadixKExchange( m_nRanks, m_Radices );
435     m_redAllMp = BabelFlow::RadixKExchangeTaskMap( m_nRanks, &m_radixkGr );
436 
437     m_redAllGr.setGraphId( 0 );
438     m_isoCalcTaskGr.setGraphId( 1 );
439     m_isoRenderTaskGr.setGraphId( 2 );
440     m_radixkGr.setGraphId( 3 );
441     m_gatherTaskGr.setGraphId( 4 );
442 
443 
444     BabelFlow::TaskGraph::registerCallback( 0, BabelFlow::RadixKExchange::LEAF_TASK_CB, allreduce );
445     BabelFlow::TaskGraph::registerCallback( 0, BabelFlow::RadixKExchange::MID_TASK_CB, allreduce );
446     BabelFlow::TaskGraph::registerCallback( 0, BabelFlow::RadixKExchange::ROOT_TASK_CB, allreduce );
447 
448     BabelFlow::TaskGraph::registerCallback( 1, BabelFlow::SingleTaskGraph::SINGLE_TASK_CB, marching_cubes );
449 
450     BabelFlow::TaskGraph::registerCallback( 2, BabelFlow::SingleTaskGraph::SINGLE_TASK_CB, vtkm_rendering );
451 
452     BabelFlow::TaskGraph::registerCallback( 3, BabelFlow::RadixKExchange::LEAF_TASK_CB, bflow_comp::volume_render_radixk );
453     BabelFlow::TaskGraph::registerCallback( 3, BabelFlow::RadixKExchange::MID_TASK_CB, bflow_comp::composite_radixk );
454     BabelFlow::TaskGraph::registerCallback( 3, BabelFlow::RadixKExchange::ROOT_TASK_CB, bflow_comp::composite_radixk );
455 
456     BabelFlow::TaskGraph::registerCallback( 4, BabelFlow::KWayReduction::LEAF_TASK_CB, BabelFlow::relay_message );
457     BabelFlow::TaskGraph::registerCallback( 4, BabelFlow::KWayReduction::MID_TASK_CB, bflow_comp::gather_results_radixk) ;
458     BabelFlow::TaskGraph::registerCallback( 4, BabelFlow::KWayReduction::ROOT_TASK_CB, bflow_comp::write_results_radixk );
459 
460     m_isoGrConnector_1 = BabelFlow::DefGraphConnector( &m_redAllGr, 0, &m_isoCalcTaskGr, 1 );
461     m_isoGrConnector_2 = BabelFlow::DefGraphConnector( &m_isoCalcTaskGr, 1, &m_isoRenderTaskGr, 2 );
462     m_isoGrConnector_3 = BabelFlow::DefGraphConnector( &m_isoRenderTaskGr, 2, &m_radixkGr, 3 );
463     m_defGraphConnector = BabelFlow::DefGraphConnector( &m_radixkGr, 3, &m_gatherTaskGr, 4 );
464 
465     std::vector<BabelFlow::TaskGraphConnector*> gr_connectors{ &m_isoGrConnector_1,
466                                                                &m_isoGrConnector_2,
467                                                                &m_isoGrConnector_3,
468                                                                &m_defGraphConnector };
469     std::vector<BabelFlow::TaskGraph*> gr_vec{ &m_redAllGr, &m_isoCalcTaskGr, &m_isoRenderTaskGr, &m_radixkGr, &m_gatherTaskGr };
470     std::vector<BabelFlow::TaskMap*> task_maps{ &m_redAllMp, &m_isoCalcTaskMp, &m_isoRenderTaskMp, &m_radixkMp, &m_gatherTaskMp };
471 
472     m_radGatherGraph = BabelFlow::ComposableTaskGraph( gr_vec, gr_connectors );
473     m_radGatherTaskMap = BabelFlow::ComposableTaskMap( task_maps );
474 
475 #ifdef BFLOW_ISO_DEBUG
476     if( m_rankId == 0 )
477     {
478       m_isoCalcTaskGr.outputGraphHtml( m_nRanks, &m_isoCalcTaskMp, "iso-gr.html" );
479       m_isoRenderTaskGr.outputGraphHtml( m_nRanks, &m_isoRenderTaskMp, "gather-task.html" );
480       m_radixkGr.outputGraphHtml( m_nRanks, &m_radixkMp, "radixk.html" );
481       // m_gatherTaskGr.outputGraphHtml( m_nRanks, &m_gatherTaskMp, "gather-task.html" );
482       m_radGatherGraph.outputGraphHtml( m_nRanks, &m_radGatherTaskMap, "bflow-iso.html" );
483     }
484 #endif
485 
486     ///
487     // MPI_Barrier(m_comm);
488     ///
489 
490     m_master.initialize( m_radGatherGraph, &m_radGatherTaskMap, m_comm, &m_contMap );
491 
492     m_inputs[ BabelFlow::TaskId(m_rankId, 0) ] = m_isoSurfData.m_bounds.serialize();
493     m_inputs[ BabelFlow::TaskId(m_rankId, 1) ] = m_isoSurfData.serialize();
494   }
495 
496 protected:
497   IsoSurfaceData m_isoSurfData;
498 
499   BabelFlow::SingleTaskGraph m_isoCalcTaskGr;
500   BabelFlow::ModuloMap m_isoCalcTaskMp;
501 
502   BabelFlow::SingleTaskGraph m_isoRenderTaskGr;
503   BabelFlow::ModuloMap m_isoRenderTaskMp;
504 
505   BabelFlow::DefGraphConnector m_isoGrConnector_1;
506   BabelFlow::DefGraphConnector m_isoGrConnector_2;
507   BabelFlow::DefGraphConnector m_isoGrConnector_3;
508 
509   BabelFlow::RadixKExchange m_redAllGr;
510   BabelFlow::RadixKExchangeTaskMap m_redAllMp;
511 };
512 
513 
514 //-----------------------------------------------------------------------------
515 }
516 //-----------------------------------------------------------------------------
517 // -- end bflow_iso --
518 //-----------------------------------------------------------------------------
519 
520 //-----------------------------------------------------------------------------
521 }
522 //-----------------------------------------------------------------------------
523 // -- end ascent --
524 //-----------------------------------------------------------------------------
525 
526 
527 //-----------------------------------------------------------------------------
528 ///
529 /// BFlowIso extract
530 ///
531 //-----------------------------------------------------------------------------
532 
declare_interface(conduit::Node & i)533 void ascent::runtime::filters::BFlowIso::declare_interface(conduit::Node &i)
534 {
535   i["type_name"] = "bflow_iso";
536   i["port_names"].append() = "in";
537   i["output_port"] = "false";  // true -- means filter, false -- means extract
538 }
539 
540 //-----------------------------------------------------------------------------
541 
verify_params(const conduit::Node & params,conduit::Node & info)542 bool ascent::runtime::filters::BFlowIso::verify_params(const conduit::Node &params, conduit::Node &info)
543 {
544   info.reset();
545 
546   bool res = true;
547 
548   res &= check_string("field", params, info, true);
549   res &= check_numeric("iso_values", params, info, true);
550   res &= check_string("image_name", params, info, true);
551   res &= check_numeric("radices", params, info, false);
552   res &= check_numeric("width", params, info, false);
553   res &= check_numeric("height", params, info, false);
554   // res &= check_string("col_field", params, info, true);
555 
556   return res;
557 }
558 
559 //-----------------------------------------------------------------------------
560 
execute()561 void ascent::runtime::filters::BFlowIso::execute()
562 {
563   if(!input(0).check_type<DataObject>())
564   {
565     ASCENT_ERROR("BFlowIso filter requires a DataObject");
566   }
567 
568   conduit::Node& p = params();
569   DataObject *data_object = input<DataObject>(0);
570   std::shared_ptr<VTKHCollection> collection = data_object->as_vtkh_collection();
571 
572   std::string image_name = p["image_name"].as_string();
573   std::string field_name = p["field"].as_string();
574   if( !collection->has_field(field_name) )
575   {
576     ASCENT_ERROR("BFlowIso cannot find the given field");
577   }
578   std::string topo_name = collection->field_topology(field_name);
579   vtkh::DataSet &data = collection->dataset_by_topology(topo_name);
580 
581   bflow_iso::IsoSurfaceData iso_surf_data;
582   iso_surf_data.m_DataSet = &data;
583   iso_surf_data.m_FieldName = field_name;
584   iso_surf_data.m_Width = DEF_IMG_WIDTH;
585   iso_surf_data.m_Height = DEF_IMG_HEIGHT;
586   iso_surf_data.m_FileName = image_name;
587   iso_surf_data.m_bounds = data.GetBounds();
588 
589   if( p.has_path("width") )
590   {
591     iso_surf_data.m_Width = p["width"].as_int64();
592   }
593 
594   if( p.has_path("height") )
595   {
596     iso_surf_data.m_Height = p["height"].as_int64();
597   }
598 
599   const conduit::Node &n_iso_vals = p["iso_values"];
600   conduit::Node n_iso_vals_dbls;
601   n_iso_vals.to_float64_array(n_iso_vals_dbls);
602   iso_surf_data.m_IsoVals.assign( n_iso_vals_dbls.as_double_ptr(), n_iso_vals_dbls.as_double_ptr() + n_iso_vals_dbls.dtype().number_of_elements() );
603 
604   MPI_Comm mpi_comm;
605   int my_rank = 0, n_ranks = 1;
606 #ifdef ASCENT_MPI_ENABLED
607   mpi_comm = MPI_Comm_f2c(flow::Workspace::default_mpi_comm());
608   MPI_Comm_rank(mpi_comm, &my_rank);
609   MPI_Comm_size(mpi_comm, &n_ranks);
610 #endif
611 
612   std::vector<uint32_t> radix_v(1);
613   radix_v[0] = n_ranks;
614   if( p.has_path("radices") )
615   {
616     conduit::DataArray<int64_t> radices_arr = p["radices"].as_int64_array();
617     radix_v.resize(radices_arr.number_of_elements());
618     for( uint32_t i = 0; i < radix_v.size(); ++i ) radix_v[i] = (uint32_t)radices_arr[i];
619   }
620 
621 #ifdef BFLOW_ISO_DEBUG
622   {
623     if( my_rank == 0 )
624     {
625       std::cout << ">> Iso surf data: <<" << std::endl;
626       std::cout << "Field name = " << iso_surf_data.m_FieldName << std::endl;
627       std::cout << "Iso vals = ";
628       for( uint32_t i = 0; i < iso_surf_data.m_IsoVals.size(); ++i )
629         std::cout << iso_surf_data.m_IsoVals[i] << "  ";
630       std::cout << std::endl;
631       std::cout << "Width = " << iso_surf_data.m_Width << std::endl;
632       std::cout << "Height = " << iso_surf_data.m_Height << std::endl;
633       std::cout << "File name = " << iso_surf_data.m_FileName << std::endl;
634       auto ranges = iso_surf_data.m_DataSet->GetRange( iso_surf_data.m_FileName );
635       int num_components = ranges.GetNumberOfValues();
636       std::cout << "Num components = " << num_components << std::endl;
637     }
638   }
639 #endif
640 
641   bflow_iso::BabelIsoGraph bflow_iso_gr( iso_surf_data, image_name, my_rank, n_ranks, mpi_comm, radix_v );
642   bflow_iso_gr.Initialize();
643   bflow_iso_gr.Execute();
644 }
645 
646