ISAAC
Overview :: Library Doc :: Server Doc :: JSON Commands

In Situ Animation of Accelerated Computations

isaac.hpp
Go to the documentation of this file.
1 /* This file is part of ISAAC.
2  *
3  * ISAAC is free software: you can redistribute it and/or modify
4  * it under the terms of the GNU Lesser General Public License as
5  * published by the Free Software Foundation, either version 3 of the
6  * License, or (at your option) any later version.
7  *
8  * ISAAC 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
11  * GNU 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 ISAAC. If not, see <www.gnu.org/licenses/>. */
15 
16 #pragma once
17 
18 //Hack for a bug, which occurs only in Cuda 7.0
19 #if (__CUDACC_VER__ < 70500) && !defined(BOOST_RESULT_OF_USE_TR1)
20  #define BOOST_RESULT_OF_USE_TR1
21 #endif
22 
23 #include <boost/config/select_compiler_config.hpp>
24 
25 #include <string>
26 #include <string.h>
27 #include <jansson.h>
28 #include <pthread.h>
29 #include <list>
30 #include <vector>
31 #include <memory>
32 #include <mpi.h>
33 //Against annoying C++11 warning in mpi.h
34 #include <IceT.h>
35 #pragma GCC diagnostic push
36 #pragma GCC diagnostic error "-Wall"
37 #include <IceTMPI.h>
38 #pragma GCC diagnostic pop
39 #include <math.h>
40 #include <stdio.h>
41 #include <assert.h>
42 #include <map>
43 
44 #if ISAAC_ALPAKA == 1
45  #include <alpaka/alpaka.hpp>
46  #include <boost/type_traits.hpp>
47  #include <boost/mpl/not.hpp>
48 #else
49  #include <boost/mpl/int.hpp>
50 #endif
51 
52 #include "isaac/isaac_kernel.hpp"
54 #include "isaac/isaac_helper.hpp"
58 #include "isaac/isaac_version.hpp"
59 
60 namespace isaac
61 {
62 
63 template <
64 #if ISAAC_ALPAKA == 1
65  typename THost,
66  typename TAcc,
67  typename TStream,
68  typename TAccDim,
69 #endif
70  typename TSimDim,
71  typename TSourceList,
72  typename TDomainSize,
73  size_t TTransfer_size,
74  typename TScale,
75  typename TController,
76  typename TCompositor
77 >
79 {
80  public:
81  #if ISAAC_ALPAKA == 1
82  using TDevAcc = alpaka::dev::Dev<TAcc>;
83  using TFraDim = alpaka::dim::DimInt<1>;
84  using TTexDim = alpaka::dim::DimInt<1>;
85  #endif
87  {
88  template
89  <
90  typename TSource,
91  typename TJsonRoot
92  >
93  ISAAC_HOST_INLINE void operator()( const int I,const TSource& s, TJsonRoot& jsonRoot) const
94  {
95  json_t *content = json_object();
96  json_array_append_new( jsonRoot, content );
97  json_object_set_new( content, "name", json_string ( TSource::getName().c_str() ) );
98  json_object_set_new( content, "feature dimension", json_integer ( s.feature_dim ) );
99  }
100  };
101 
103  {
104  template
105  <
106  typename TFunctor,
107  typename TJsonRoot
108  >
109  ISAAC_HOST_INLINE void operator()( const int I,const TFunctor& f, TJsonRoot& jsonRoot) const
110  {
111  json_t *content = json_object();
112  json_array_append_new( jsonRoot, content );
113  json_object_set_new( content, "name", json_string ( TFunctor::getName().c_str() ) );
114  json_object_set_new( content, "description", json_string ( TFunctor::getDescription().c_str() ) );
115  json_object_set_new( content, "uses parameter", json_boolean ( TFunctor::uses_parameter ) );
116  }
117  };
118 
120  {
121  template
122  <
123  typename TFunctor,
124  typename TName,
125  typename TValue,
126  typename TFound
127  >
128  ISAAC_HOST_INLINE void operator()( const int I,TFunctor& f,const TName& name, TValue& value, TFound& found) const
129  {
130  if (!found && name == TFunctor::getName())
131  {
132  value = I;
133  found = true;
134  }
135  }
136  };
137 
139  {
140  template
141  <
142  typename TSource,
143  typename TFunctions,
144  typename TDest
145  >
147  const int I,
148  const TSource& source,
149  const TFunctions& functions,
150  TDest& dest
151  ) const
152  {
153  isaac_int chain_nr = 0;
154  for (int i = 0; i < ISAAC_MAX_FUNCTORS; i++)
155  {
156  chain_nr *= ISAAC_FUNCTOR_COUNT;
157  chain_nr += functions[I].bytecode[i];
158  }
159  dest.nr[I] = chain_nr * 4 + TSource::feature_dim - 1;
160  }
161  };
162 
164  {
165  template
166  <
167  typename TSource,
168  typename TArray,
169  typename TLocalSize
170  #if ISAAC_ALPAKA == 1
171  ,typename TVector
172  ,typename TDevAcc__
173  #endif
174  >
176  const int I,
177  const TSource& source,
178  TArray& pointer_array,
179  const TLocalSize& local_size
180  #if ISAAC_ALPAKA == 1
181  ,TVector& alpaka_vector
182  ,const TDevAcc__& acc
183  #endif
184  ) const
185  {
186  if (TSource::persistent)
187  pointer_array.pointer[I] = NULL;
188  else
189  {
190  #if ISAAC_ALPAKA == 1
191  alpaka_vector.push_back(
192  alpaka::mem::buf::Buf< TDevAcc, isaac_float, TFraDim, size_t> (
193  alpaka::mem::buf::alloc<isaac_float, size_t> (
194  acc,
195  alpaka::vec::Vec<TFraDim, size_t> (
196  TSource::feature_dim * (local_size[0] + 2 * ISAAC_GUARD_SIZE) * (local_size[1] + 2 * ISAAC_GUARD_SIZE) * (local_size[2] + 2 * ISAAC_GUARD_SIZE)
197  )
198  )
199  )
200  );
201  pointer_array.pointer[I] = alpaka::mem::view::getPtrNative( alpaka_vector.back() );
202  #else
203  ISAAC_CUDA_CHECK(cudaMalloc((void**)&(pointer_array.pointer[I]), sizeof(isaac_float_dim< TSource::feature_dim >) * (local_size[0] + 2 * ISAAC_GUARD_SIZE) * (local_size[1] + 2 * ISAAC_GUARD_SIZE) * (local_size[2] + 2 * ISAAC_GUARD_SIZE) ) );
204  #endif
205  }
206  }
207  };
208 
210  {
211  template
212  <
213  typename TSource,
214  typename TArray,
215  typename TLocalSize,
216  typename TWeight,
217  typename TPointer
218  #if ISAAC_ALPAKA == 1
219  ,typename TStream__
220  #endif
221  >
223  const int I,
224  TSource& source,
225  TArray& pointer_array,
226  const TLocalSize& local_size,
227  const TWeight& weight,
228  const TPointer& pointer
229  #if ISAAC_ALPAKA == 1
230  ,TStream__& stream
231  #endif
232  ) const
233  {
234  bool enabled = weight.value[ I ] != isaac_float(0);
235  source.update( enabled, pointer );
236  if (!TSource::persistent && enabled)
237  {
238  isaac_size2 grid_size=
239  {
240  size_t((local_size[0] + ISAAC_GUARD_SIZE * 2 + 15)/16),
241  size_t((local_size[1] + ISAAC_GUARD_SIZE * 2 + 15)/16),
242  };
243  isaac_size2 block_size=
244  {
245  size_t(16),
246  size_t(16),
247  };
248  isaac_int3 local_size_array = { isaac_int(local_size[0]), isaac_int(local_size[1]), isaac_int(local_size[2]) };
249  #if ISAAC_ALPAKA == 1
250 #if ALPAKA_ACC_GPU_CUDA_ENABLED == 1
251  if ( mpl::not_<boost::is_same<TAcc, alpaka::acc::AccGpuCudaRt<TAccDim, size_t> > >::value )
252 #endif
253  {
254  grid_size.x = size_t(local_size[0] + ISAAC_GUARD_SIZE * 2);
255  grid_size.y = size_t(local_size[0] + ISAAC_GUARD_SIZE * 2);
256  block_size.x = size_t(1);
257  block_size.y = size_t(1);
258  }
259  const alpaka::vec::Vec<TAccDim, size_t> threads (size_t(1), size_t(1), size_t(1));
260  const alpaka::vec::Vec<TAccDim, size_t> blocks (size_t(1), block_size.x, block_size.y);
261  const alpaka::vec::Vec<TAccDim, size_t> grid (size_t(1), grid_size.x, grid_size.y);
262  auto const workdiv(alpaka::workdiv::WorkDivMembers<TAccDim, size_t>(grid,blocks,threads));
263  updateBufferKernel<TSource> kernel;
264  auto const instance
265  (
266  alpaka::exec::create<TAcc>
267  (
268  workdiv,
269  kernel,
270  source,
271  pointer_array.pointer[ I ],
272  local_size_array
273  )
274  );
275  alpaka::stream::enqueue(stream, instance);
276  alpaka::wait::wait(stream);
277  #else
278  dim3 block (block_size.x, block_size.y);
279  dim3 grid (grid_size.x, grid_size.y);
280  updateBufferKernel<<<grid,block>>>( source, pointer_array.pointer[ I ], local_size_array );
281  ISAAC_CUDA_CHECK(cudaDeviceSynchronize());
282  #endif
283  }
284  }
285  };
286 
288  {
289  template
290  <
291  typename TSource,
292  typename TArray,
293  typename TMinmax,
294  typename TLocalMinmax,
295  typename TLocalSize
296  #if ISAAC_ALPAKA == 1
297  ,typename TStream__
298  ,typename THost__
299  #endif
300  >
302  const int I,
303  const TSource& source,
304  TArray& pointer_array,
305  TMinmax& minmax,
306  TLocalMinmax& local_minmax,
307  TLocalSize& local_size
308  #if ISAAC_ALPAKA == 1
309  ,TStream__& stream
310  ,const THost__& host
311  #endif
312  ) const
313  {
314  isaac_size2 grid_size=
315  {
316  size_t((local_size[0]+15)/16),
317  size_t((local_size[1]+15)/16),
318  };
319  isaac_size2 block_size=
320  {
321  size_t(16),
322  size_t(16),
323  };
324  isaac_int3 local_size_array = { isaac_int(local_size[0]), isaac_int(local_size[1]), isaac_int(local_size[2]) };
325  minmax_struct local_minmax_array_h[ local_size_array.x * local_size_array.y ];
326  #if ISAAC_ALPAKA == 1
327 #if ALPAKA_ACC_GPU_CUDA_ENABLED == 1
328  if ( mpl::not_<boost::is_same<TAcc, alpaka::acc::AccGpuCudaRt<TAccDim, size_t> > >::value )
329 #endif
330  {
331  grid_size.x = size_t(local_size[0]);
332  grid_size.y = size_t(local_size[0]);
333  block_size.x = size_t(1);
334  block_size.y = size_t(1);
335  }
336  const alpaka::vec::Vec<TAccDim, size_t> threads (size_t(1), size_t(1), size_t(1));
337  const alpaka::vec::Vec<TAccDim, size_t> blocks (size_t(1), block_size.x, block_size.y);
338  const alpaka::vec::Vec<TAccDim, size_t> grid (size_t(1), grid_size.x, grid_size.y);
339  auto const workdiv(alpaka::workdiv::WorkDivMembers<TAccDim, size_t>(grid,blocks,threads));
340  minMaxKernel<TSource> kernel;
341  auto const instance
342  (
343  alpaka::exec::create<TAcc>
344  (
345  workdiv,
346  kernel,
347  source,
348  I,
349  alpaka::mem::view::getPtrNative(local_minmax),
350  local_size_array,
351  pointer_array.pointer[ I ]
352  )
353  );
354  alpaka::stream::enqueue(stream, instance);
355  alpaka::wait::wait(stream);
356  alpaka::mem::view::ViewPlainPtr<THost, minmax_struct, TFraDim, size_t> minmax_buffer(local_minmax_array_h, host, alpaka::vec::Vec<TFraDim, size_t>(size_t(local_size_array.x * local_size_array.y)));
357  alpaka::mem::view::copy( stream, minmax_buffer, local_minmax, alpaka::vec::Vec<TFraDim, size_t>(size_t(local_size_array.x * local_size_array.y)));
358  #else
359  dim3 block (block_size.x, block_size.y);
360  dim3 grid (grid_size.x, grid_size.y);
361  minMaxKernel<<<grid,block>>>( source, I, local_minmax, local_size_array, pointer_array.pointer[ I ]);
362  ISAAC_CUDA_CHECK(cudaMemcpy( local_minmax_array_h, local_minmax, sizeof(minmax_struct)*local_size_array.x * local_size_array.y, cudaMemcpyDeviceToHost));
363  #endif
364  minmax.min[ I ] = FLT_MAX;
365  minmax.max[ I ] = -FLT_MAX;
366  for (int i = 0; i < local_size_array.x * local_size_array.y; i++)
367  {
368  if ( local_minmax_array_h[i].min < minmax.min[ I ])
369  minmax.min[ I ] = local_minmax_array_h[i].min;
370  if ( local_minmax_array_h[i].max > minmax.max[ I ])
371  minmax.max[ I ] = local_minmax_array_h[i].max;
372  }
373  }
374  };
375 
377  #if ISAAC_ALPAKA == 1
378  THost host,
379  TDevAcc acc,
380  TStream stream,
381  #endif
382  const std::string name,
383  const isaac_int master,
384  const std::string server_url,
385  const isaac_uint server_port,
386  const isaac_size2 framebuffer_size,
387  const TDomainSize global_size,
388  const TDomainSize local_size,
389  const TDomainSize position,
390  TSourceList& sources,
391  TScale scale
392  ) :
393  #if ISAAC_ALPAKA == 1
394  host(host),
395  acc(acc),
396  stream(stream),
397  #endif
398  global_size(global_size),
399  local_size(local_size),
400  position(position),
401  name(name),
402  master(master),
403  server_url(server_url),
404  server_port(server_port),
405  framebuffer_size(framebuffer_size),
406  compbuffer_size(TCompositor::getCompositedbufferSize(framebuffer_size)),
407  compositor(framebuffer_size),
408  metaNr(0),
409  visualizationThread(0),
410  kernel_time(0),
411  merge_time(0),
412  video_send_time(0),
413  copy_time(0),
414  sorting_time(0),
415  buffer_time(0),
416  interpolation(false),
417  iso_surface(false),
419  framebuffer_prod(size_t(framebuffer_size.x) * size_t(framebuffer_size.y)),
420  sources( sources ),
421  scale( scale ),
422  icet_bounding_box( true )
423  #if ISAAC_ALPAKA == 1
424  ,framebuffer(alpaka::mem::buf::alloc<uint32_t, size_t>(acc, framebuffer_prod))
425  ,functor_chain_d(alpaka::mem::buf::alloc<isaac_functor_chain_pointer_N, size_t>(acc, size_t( ISAAC_FUNCTOR_COMPLEX * 4)))
426  ,functor_chain_choose_d(alpaka::mem::buf::alloc<isaac_functor_chain_pointer_N, size_t>(acc, size_t( boost::mpl::size< TSourceList >::type::value )))
427  ,local_minmax_array_d(alpaka::mem::buf::alloc<minmax_struct, size_t>(acc, size_t( local_size[0] * local_size[1] )))
428  {
429  #else
430  {
431  ISAAC_CUDA_CHECK(cudaMalloc((uint32_t**)&framebuffer, sizeof(uint32_t)*framebuffer_prod));
433  ISAAC_CUDA_CHECK(cudaMalloc((isaac_functor_chain_pointer_N**)&functor_chain_choose_d, sizeof(isaac_functor_chain_pointer_N) * boost::mpl::size< TSourceList >::type::value));
434  ISAAC_CUDA_CHECK(cudaMalloc((minmax_struct**)&local_minmax_array_d, sizeof(minmax_struct) * local_size[0] * local_size[1]));
435  #endif
436  #if ISAAC_VALGRIND_TWEAKS == 1
437  //Jansson has some optimizations for 2 and 4 byte aligned
438  //memory. However valgrind complains then sometimes about reads
439  //in not allocated memory. Valgrind is right, but nevertheless
440  //reads will never crash and be much faster. But for
441  //debugging reasons let's alloc 4 extra bytes for valgrind:
442  json_set_alloc_funcs(extra_malloc, extra_free);
443  #endif
444  json_object_seed(0);
445  for (int i = 0; i < 3; i++)
446  {
447  global_size_scaled.push_back( isaac_int( (isaac_float)global_size[i] * (isaac_float)scale[i] ) );
448  local_size_scaled.push_back( isaac_int( (isaac_float) local_size[i] * (isaac_float)scale[i] ) );
449  position_scaled.push_back( isaac_int( (isaac_float) position[i] * (isaac_float)scale[i] ) );
450  }
451 
452  background_color[0] = 0;
453  background_color[1] = 0;
454  background_color[2] = 0;
455  background_color[3] = 1;
456 
457  //INIT
458  MPI_Comm_dup(MPI_COMM_WORLD, &mpi_world);
459  MPI_Comm_rank(mpi_world, &rank);
460  MPI_Comm_size(mpi_world, &numProc);
461  if (rank == master)
462  {
463  this->communicator = new IsaacCommunicator(server_url,server_port);
464  }
465  else
466  {
467  this->communicator = NULL;
468  }
469  recreateJSON();
470  controller.updateProjection( projection, framebuffer_size, NULL, true );
471  look_at[0] = 0.0f;
472  look_at[1] = 0.0f;
473  look_at[2] = 0.0f;
474  ISAAC_SET_IDENTITY(3,rotation)
475  distance = -4.5f;
476  updateModelview();
477 
478  //Create functor chain pointer lookup table
479  #if ISAAC_ALPAKA == 1
480  const alpaka::vec::Vec<TAccDim, size_t> threads (size_t(1), size_t(1), size_t(1));
481  const alpaka::vec::Vec<TAccDim, size_t> blocks (size_t(1), size_t(1), size_t(1));
482  const alpaka::vec::Vec<TAccDim, size_t> grid (size_t(1), size_t(1), size_t(1));
483  auto const workdiv(alpaka::workdiv::WorkDivMembers<TAccDim, size_t>(grid,blocks,threads));
485  auto const instance
486  (
487  alpaka::exec::create<TAcc>
488  (
489  workdiv,
490  kernel,
491  alpaka::mem::view::getPtrNative(functor_chain_d)
492  )
493  );
494  alpaka::stream::enqueue(stream, instance);
495  alpaka::wait::wait(stream);
496  #else
497  dim3 grid(1);
498  dim3 block(1);
499  fillFunctorChainPointerKernel<<<grid,block>>>( functor_chain_d );
500  ISAAC_CUDA_CHECK(cudaDeviceSynchronize());
501  #endif
502  //Init functions:
503  for (int i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
504  functions[i].source = std::string("idem");
505  updateFunctions();
506 
507  //non persistent buffer memory
508  isaac_for_each_params(sources,allocate_pointer_array_iterator(),pointer_array,local_size
509  #if ISAAC_ALPAKA == 1
510  ,pointer_array_alpaka
511  ,acc
512  #endif
513  );
514 
515  //Transfer func memory:
516  for (int i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
517  {
518  source_weight.value[i] = ISAAC_DEFAULT_WEIGHT;
519  #if ISAAC_ALPAKA == 1
520  transfer_d_buf.push_back( alpaka::mem::buf::Buf<TDevAcc, isaac_float4, TTexDim, size_t> ( alpaka::mem::buf::alloc<isaac_float4, size_t>( acc, alpaka::vec::Vec<TTexDim, size_t> ( TTransfer_size ) ) ) );
521  transfer_h_buf.push_back( alpaka::mem::buf::Buf< THost, isaac_float4, TTexDim, size_t> ( alpaka::mem::buf::alloc<isaac_float4, size_t>(host, alpaka::vec::Vec<TTexDim, size_t> ( TTransfer_size ) ) ) );
522  transfer_d.pointer[i] = alpaka::mem::view::getPtrNative( transfer_d_buf[i] );
523  transfer_h.pointer[i] = alpaka::mem::view::getPtrNative( transfer_h_buf[i] );
524  #else
525  ISAAC_CUDA_CHECK(cudaMalloc((isaac_float4**)&(transfer_d.pointer[i]), sizeof(isaac_float4)*TTransfer_size));
526  transfer_h.pointer[i] = (isaac_float4*)malloc( sizeof(isaac_float4)*TTransfer_size );
527  #endif
528  transfer_h.description[i].insert( std::pair< isaac_uint, isaac_float4> (0 , getHSVA(isaac_float(2*i)*M_PI/isaac_float(boost::mpl::size< TSourceList >::type::value),1,1,0) ));
529  transfer_h.description[i].insert( std::pair< isaac_uint, isaac_float4> (TTransfer_size, getHSVA(isaac_float(2*i)*M_PI/isaac_float(boost::mpl::size< TSourceList >::type::value),1,1,1) ));
530  }
531  updateTransfer();
532 
533  max_size = ISAAC_MAX(uint32_t(global_size[0]),uint32_t(global_size[1]));
534  if (TSimDim::value > 2)
535  max_size = ISAAC_MAX(uint32_t(global_size[2]),uint32_t(max_size));
536  max_size_scaled = ISAAC_MAX(uint32_t(global_size_scaled[0]),uint32_t(global_size_scaled[1]));
537  if (TSimDim::value > 2)
538  max_size_scaled = ISAAC_MAX(uint32_t(global_size_scaled[2]),uint32_t(max_size_scaled));
539 
540  //ICET:
541  IceTCommunicator icetComm;
542  icetComm = icetCreateMPICommunicator(mpi_world);
543  for (int pass = 0; pass < TController::pass_count; pass++)
544  {
545  icetContext[pass] = icetCreateContext(icetComm);
546  icetResetTiles();
547  icetAddTile(0, 0, framebuffer_size.x, framebuffer_size.y, master);
548  //icetStrategy(ICET_STRATEGY_DIRECT);
549  icetStrategy(ICET_STRATEGY_SEQUENTIAL);
550  //icetStrategy(ICET_STRATEGY_REDUCE);
551 
552  //icetSingleImageStrategy( ICET_SINGLE_IMAGE_STRATEGY_AUTOMATIC );
553  icetSingleImageStrategy( ICET_SINGLE_IMAGE_STRATEGY_BSWAP );
554  //icetSingleImageStrategy( ICET_SINGLE_IMAGE_STRATEGY_RADIXK );
555  //icetSingleImageStrategy( ICET_SINGLE_IMAGE_STRATEGY_TREE );
556 
557  /*IceTBoolean supports;
558  icetGetBooleanv( ICET_STRATEGY_SUPPORTS_ORDERING, &supports );
559  if (supports)
560  printf("yes\n");
561  else
562  printf("no\n");*/
563 
564  icetSetColorFormat(ICET_IMAGE_COLOR_RGBA_UBYTE);
565  icetSetDepthFormat(ICET_IMAGE_DEPTH_NONE);
566  icetCompositeMode(ICET_COMPOSITE_MODE_BLEND);
567  icetEnable(ICET_ORDERED_COMPOSITE);
568  icetPhysicalRenderSize(framebuffer_size.x, framebuffer_size.y);
569  icetDrawCallback( drawCallBack );
570  }
571  icetDestroyMPICommunicator(icetComm);
572  updateBounding( );
573 
574  //JSON
575  if (rank == master)
576  {
577  json_object_set_new( json_root, "type", json_string( "register" ) );
578  json_object_set_new( json_root, "name", json_string( name.c_str() ) );
579  json_object_set_new( json_root, "nodes", json_integer( numProc ) );
580  json_object_set_new( json_root, "framebuffer width", json_integer ( compbuffer_size.x ) );
581  json_object_set_new( json_root, "framebuffer height", json_integer ( compbuffer_size.y ) );
582 
583  json_object_set_new( json_root, "max functors", json_integer( ISAAC_MAX_FUNCTORS ) );
584  json_t *json_functors_array = json_array();
585  json_object_set_new( json_root, "functors", json_functors_array );
586  IsaacFunctorPool functors;
587  isaac_for_each_params( functors, functor_2_json_iterator(), json_functors_array );
588 
589  json_t *matrix;
590  json_object_set_new( json_root, "projection count", json_integer ( TController::pass_count ) );
591  json_object_set_new( json_root, "projection", matrix = json_array() );
592  ISAAC_JSON_ADD_MATRIX(matrix,projection,16 * TController::pass_count)
593  json_object_set_new( json_root, "position", matrix = json_array() );
594  ISAAC_JSON_ADD_MATRIX(matrix,look_at,3)
595  json_object_set_new( json_root, "rotation", matrix = json_array() );
596  ISAAC_JSON_ADD_MATRIX(matrix,rotation,9)
597  json_object_set_new( json_root, "distance", json_real( distance ) );
598 
599  json_t *json_sources_array = json_array();
600  json_object_set_new( json_root, "sources", json_sources_array );
601 
602  isaac_for_each_params( sources, source_2_json_iterator(), json_sources_array );
603 
604  json_object_set_new( json_root, "interpolation", json_boolean( interpolation ) );
605  json_object_set_new( json_root, "iso surface", json_boolean( iso_surface ) );
606  json_object_set_new( json_root, "step", json_real( step ) );
607 
608  json_object_set_new( json_root, "dimension", json_integer ( TSimDim::value ) );
609  json_object_set_new( json_root, "width", json_integer ( global_size_scaled[0] ) );
610  if (TSimDim::value > 1)
611  json_object_set_new( json_root, "height", json_integer ( global_size_scaled[1] ) );
612  if (TSimDim::value > 2)
613  json_object_set_new( json_root, "depth", json_integer ( global_size_scaled[2] ) );
614  json_t *json_version_array = json_array();
615  json_array_append_new( json_version_array, json_integer( ISAAC_PROTOCOL_VERSION_MAJOR ) );
616  json_array_append_new( json_version_array, json_integer( ISAAC_PROTOCOL_VERSION_MINOR ) );
617  json_object_set_new( json_root, "protocol", json_version_array);
618  }
619  }
620  void setJpegQuality(isaac_uint jpeg_quality)
621  {
623  if (communicator)
624  communicator->setJpegQuality(jpeg_quality);
625  }
627  {
629  if (nr >= ISAAC_MAX_CLIPPING)
630  return false;
631  isaac_float nx_s = nx * scale[0];
632  isaac_float ny_s = ny * scale[1];
633  isaac_float nz_s = nz * scale[2];
634  isaac_float l = sqrt(nx_s*nx_s+ny_s*ny_s+nz_s*nz_s);
635  if (l == 0.0f)
636  return false;
637  nx_s /= l;
638  ny_s /= l;
639  nz_s /= l;
640  clipping.elem[ nr ].position.x = px;
641  clipping.elem[ nr ].position.y = py;
642  clipping.elem[ nr ].position.z = pz;
643  clipping.elem[ nr ].normal.x = nx_s;
644  clipping.elem[ nr ].normal.y = ny_s;
645  clipping.elem[ nr ].normal.z = nz_s;
646  clipping_saved_normals[ nr ].x = nx;
647  clipping_saved_normals[ nr ].y = ny;
648  clipping_saved_normals[ nr ].z = nz;
649  return true;
650  }
652  {
653  if ( editClipping( clipping.count, px, py, pz, nx, ny, nz ) )
654  clipping.count++;
655  }
657  {
659  if (nr >= clipping.count)
660  return;
661  clipping.count--;
662  for (isaac_uint i = nr; i < clipping.count; i++)
663  {
664  clipping.elem[i] = clipping.elem[i+1];
665  clipping_saved_normals[i] = clipping_saved_normals[i+1];
666  }
667  }
669  {
671  for (int pass = 0; pass < TController::pass_count; pass++)
672  {
673  icetSetContext( icetContext[pass] );
674  if (icet_bounding_box)
675  {
676  isaac_float f_l_width = (isaac_float)local_size_scaled[0]/(isaac_float)max_size_scaled * 2.0f;
677  isaac_float f_l_height = (isaac_float)local_size_scaled[1]/(isaac_float)max_size_scaled * 2.0f;
678  isaac_float f_l_depth = 0.0f;
679  if (TSimDim::value > 2)
680  f_l_depth = (isaac_float)local_size_scaled[2]/(isaac_float)max_size_scaled * 2.0f;
681  isaac_float f_x = (isaac_float)position_scaled[0]/(isaac_float)max_size_scaled * 2.0f - (isaac_float)global_size_scaled[0]/(isaac_float)max_size_scaled;
682  isaac_float f_y = (isaac_float)position_scaled[1]/(isaac_float)max_size_scaled * 2.0f - (isaac_float)global_size_scaled[1]/(isaac_float)max_size_scaled;
683  isaac_float f_z = 0.0f;
684  if (TSimDim::value > 2)
685  f_z = (isaac_float)position_scaled[2]/(isaac_float)max_size_scaled * isaac_float(2) - (isaac_float)global_size_scaled[2]/(isaac_float)max_size_scaled;
686  icetBoundingBoxf( f_x, f_x + f_l_width, f_y, f_y + f_l_height, f_z, f_z + f_l_depth);
687  }
688  else
689  icetBoundingVertices(0,0,0,0,NULL);
690  }
691  }
692  void updatePosition( const TDomainSize position )
693  {
695  this->position = position;
696  for (int i = 0; i < 3; i++)
697  position_scaled[i] = isaac_int( (isaac_float) position[i] * (isaac_float)scale[i] );
698  }
699  void updateLocalSize( const TDomainSize local_size )
700  {
702  this->local_size = local_size;
703  for (int i = 0; i < 3; i++)
704  local_size_scaled[i] = isaac_int( (isaac_float) local_size[i] * (isaac_float)scale[i] );
705  }
707  {
709  IsaacFunctorPool functors;
710  isaac_float4 isaac_parameter_h[boost::mpl::size< TSourceList >::type::value * ISAAC_MAX_FUNCTORS];
711  for (int i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
712  {
713  functions[i].error_code = 0;
714  //Going from | to |...
715  std::string source = functions[i].source;
716  size_t pos = 0;
717  bool again = true;
718  int elem = 0;
719  while (again && ((pos = source.find("|")) != std::string::npos || ((again = false) == false)) )
720  {
721  if (elem >= ISAAC_MAX_FUNCTORS)
722  {
723  functions[i].error_code = 1;
724  break;
725  }
726  std::string token = source.substr(0, pos);
727  source.erase(0, pos + 1);
728  //ignore " " in token
729  token.erase(remove_if(token.begin(), token.end(), isspace), token.end());
730  //search "(" and parse parameters
731  size_t t_begin = token.find("(");
732  if (t_begin == std::string::npos)
733  memset(&(isaac_parameter_h[i * ISAAC_MAX_FUNCTORS + elem]), 0, sizeof(isaac_float4));
734  else
735  {
736  size_t t_end = token.find(")");
737  if (t_end == std::string::npos)
738  {
739  functions[i].error_code = -1;
740  break;
741  }
742  if (t_end - t_begin == 1) //()
743  memset(&(isaac_parameter_h[i * ISAAC_MAX_FUNCTORS + elem]), 0, sizeof(isaac_float4));
744  else
745  {
746  std::string parameters = token.substr(t_begin+1, t_end-t_begin-1);
747  size_t p_pos = 0;
748  bool p_again = true;
749  int p_elem = 0;
750  isaac_float* parameter_array = (float*)&(isaac_parameter_h[i * ISAAC_MAX_FUNCTORS + elem]);
751  while (p_again && ((p_pos = parameters.find(",")) != std::string::npos || ((p_again = false) == false)) )
752  {
753  if (p_elem >= 4)
754  {
755  functions[i].error_code = 2;
756  break;
757  }
758  std::string par = parameters.substr(0, p_pos);
759  parameters.erase(0, p_pos + 1);
760  parameter_array[p_elem] = std::stof( par );
761  p_elem++;
762  }
763  for (;p_elem < 4; p_elem++)
764  parameter_array[p_elem] = parameter_array[p_elem - 1]; //last one repeated
765  }
766  }
767  //parse token itself
768  if (t_begin != std::string::npos)
769  token = token.substr(0, t_begin);
770  bool found = false;
771  isaac_for_each_params( functors, parse_functor_iterator(), token, functions[i].bytecode[elem], found );
772  if (!found)
773  {
774  functions[i].error_code = -2;
775  break;
776  }
777 
778  elem++;
779  }
780  for (;elem < ISAAC_MAX_FUNCTORS; elem++)
781  {
782  functions[i].bytecode[elem] = 0; //last one idem
783  memset(&(isaac_parameter_h[i * ISAAC_MAX_FUNCTORS + elem]), 0, sizeof(isaac_float4));
784  }
785  }
786 
787  //Calculate functor chain nr per source
789  isaac_for_each_params( sources, update_functor_chain_iterator(), functions, dest);
790  #if ISAAC_ALPAKA == 1
791  alpaka::mem::view::ViewPlainPtr<THost, isaac_float4, TFraDim, size_t> parameter_buffer(isaac_parameter_h, host, alpaka::vec::Vec<TFraDim, size_t>(size_t(ISAAC_MAX_FUNCTORS * boost::mpl::size< TSourceList >::type::value)));
792 
793  alpaka::vec::Vec<alpaka::dim::DimInt<1u>, size_t> const parameter_d_extent(size_t(16));
794  auto parameter_d_view(alpaka::mem::view::createStaticDevMemView(&isaac_parameter_d[0u],acc,parameter_d_extent));
795  alpaka::mem::view::copy( stream, parameter_d_view, parameter_buffer, alpaka::vec::Vec<TFraDim, size_t>(size_t(ISAAC_MAX_FUNCTORS * boost::mpl::size< TSourceList >::type::value)) );
796 
797  const alpaka::vec::Vec<TAccDim, size_t> threads (size_t(1), size_t(1), size_t(1));
798  const alpaka::vec::Vec<TAccDim, size_t> blocks (size_t(1), size_t(1), size_t(1));
799  const alpaka::vec::Vec<TAccDim, size_t> grid (size_t(1), size_t(1), size_t(1));
800  auto const workdiv(alpaka::workdiv::WorkDivMembers<TAccDim, size_t>(grid,blocks,threads));
802  <
803  boost::mpl::size< TSourceList >::type::value,
805  > kernel;
806  auto const instance
807  (
808  alpaka::exec::create<TAcc>
809  (
810  workdiv,
811  kernel,
812  alpaka::mem::view::getPtrNative(functor_chain_choose_d),
813  alpaka::mem::view::getPtrNative(functor_chain_d),
814  dest
815  )
816  );
817  alpaka::stream::enqueue(stream, instance);
818  alpaka::wait::wait(stream);
819 
820  alpaka::vec::Vec<alpaka::dim::DimInt<1u>, size_t> const function_chain_d_extent( size_t(ISAAC_MAX_SOURCES) );
821  auto function_chain_d_view(alpaka::mem::view::createStaticDevMemView(&isaac_function_chain_d[0u],acc,function_chain_d_extent));
822  alpaka::mem::view::copy( stream, function_chain_d_view, functor_chain_choose_d, size_t( boost::mpl::size< TSourceList >::type::value ) );
823  #else
824  ISAAC_CUDA_CHECK(cudaMemcpyToSymbol( isaac_parameter_d, isaac_parameter_h, sizeof( isaac_parameter_h )));
825  dim3 grid(1);
826  dim3 block(1);
827  updateFunctorChainPointerKernel< boost::mpl::size< TSourceList >::type::value > <<<grid,block>>>(functor_chain_choose_d, functor_chain_d, dest);
828  ISAAC_CUDA_CHECK(cudaDeviceSynchronize());
829  isaac_functor_chain_pointer_N* constant_ptr;
830  ISAAC_CUDA_CHECK(cudaGetSymbolAddress((void **)&constant_ptr, isaac_function_chain_d));
831  ISAAC_CUDA_CHECK(cudaMemcpy(constant_ptr, functor_chain_choose_d, boost::mpl::size< TSourceList >::type::value * sizeof( isaac_functor_chain_pointer_N ), cudaMemcpyDeviceToDevice));
832  #endif
833  }
835  {
837  for (int i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
838  {
839  auto next = transfer_h.description[i].begin();
840  auto before = next;
841  for(next++; next != transfer_h.description[i].end(); next++)
842  {
843  isaac_uint width = next->first - before->first;
844  for (size_t j = 0; j < width && j + before->first < TTransfer_size; j++)
845  {
846  transfer_h.pointer[i][before->first + j] = (
847  before->second * isaac_float(width-j) +
848  next->second * isaac_float(j)
849  ) / isaac_float( width );
850  }
851  before = next;
852  }
853  #if ISAAC_ALPAKA == 1
854  alpaka::mem::view::copy(stream, transfer_d_buf[i], transfer_h_buf[i], TTransfer_size );
855  #else
856  ISAAC_CUDA_CHECK(cudaMemcpy(transfer_d.pointer[i], transfer_h.pointer[i], sizeof(isaac_float4)*TTransfer_size, cudaMemcpyHostToDevice));
857  #endif
858  }
859  }
860  json_t* getJsonMetaRoot()
861  {
863  return json_meta_root;
864  }
865  int init(CommunicatorSetting communicatorBehaviour = ReturnAtError)
866  {
867  int failed = 0;
868  if (communicator && (communicator->serverConnect(communicatorBehaviour) < 0))
869  failed = 1;
870  MPI_Bcast(&failed,1, MPI_INT, master, mpi_world);
871  if (failed)
872  return -1;
873  if (rank == master)
874  {
875  json_init_root = json_root;
876  communicator->serverSendRegister(&json_init_root);
877  recreateJSON();
878  }
879  return 0;
880  }
881  json_t* doVisualization( const IsaacVisualizationMetaEnum metaTargets = META_MASTER, void* pointer = NULL, bool redraw = true)
882  {
883  if (redraw)
884  {
886  isaac_for_each_params( sources, update_pointer_array_iterator(), pointer_array, local_size, source_weight, pointer
887  #if ISAAC_ALPAKA == 1
888  ,stream
889  #endif
890  );
892  }
893  //if (rank == master)
894  // printf("-----\n");
896 
897  myself = this;
898 
899  send_distance = false;
900  send_look_at = false;
901  send_projection = false;
902  send_rotation = false;
903  send_transfer = false;
904  send_interpolation = false;
905  send_step = false;
906  send_iso_surface = false;
907  send_functions = false;
908  send_weight = false;
909  send_minmax = false;
910  send_background_color = false;
911  send_clipping = false;
912  send_controller = false;
913  send_init_json = false;
914 
915  //Handle messages
916  json_t* message;
917  char message_buffer[ISAAC_MAX_RECEIVE] = "{}";
918  //Master merges all messages and broadcasts it.
919 
920  if (rank == master)
921  {
922  message = json_object();
923  bool add_modelview = false;
924  while (json_t* last = communicator->getLastMessage())
925  {
926  json_t* js;
927  size_t index;
928  json_t *value;
929  //search for update requests
930  if ( js = json_object_get(last, "request") )
931  {
932  const char* target = json_string_value( js );
933  if ( strcmp( target, "rotation" ) == 0 )
934  send_rotation = true;
935  if ( strcmp( target, "position" ) == 0 )
936  send_look_at = true;
937  if ( strcmp( target, "distance" ) == 0 )
938  send_distance = true;
939  if ( strcmp( target, "projection" ) == 0 )
940  send_projection = true;
941  if ( strcmp( target, "transfer" ) == 0 )
942  send_transfer = true;
943  if ( strcmp( target, "interpolation" ) == 0 )
944  send_interpolation = true;
945  if ( strcmp( target, "step" ) == 0 )
946  send_step = true;
947  if ( strcmp( target, "iso surface" ) == 0 )
948  send_iso_surface = true;
949  if ( strcmp( target, "functions" ) == 0 )
950  send_functions = true;
951  if ( strcmp( target, "weight" ) == 0 )
952  send_weight = true;
953  if ( strcmp( target, "background color" ) == 0 )
954  send_background_color = true;
955  if ( strcmp( target, "clipping" ) == 0 )
956  send_clipping = true;
957  if ( strcmp( target, "controller" ) == 0 )
958  send_controller = true;
959  if ( strcmp( target, "init" ) == 0 )
960  send_init_json = true;
961  }
962  //Search for scene changes
963  if (json_array_size( js = json_object_get(last, "rotation absolute") ) == 9)
964  {
965  add_modelview = true;
966  send_rotation = true;
967  json_array_foreach(js, index, value)
968  rotation[index] = json_number_value( value );
969  json_object_del( last, "rotation absolute" );
970  }
971  if (json_array_size( js = json_object_get(last, "rotation relative") ) == 9)
972  {
973  add_modelview = true;
974  send_rotation = true;
975  IceTDouble relative[9];
976  IceTDouble new_rotation[9];
977  json_array_foreach(js, index, value)
978  relative[index] = json_number_value( value );
979  for (isaac_int x = 0; x < 3; x++)
980  for (isaac_int y = 0; y < 3; y++)
981  new_rotation[y+x*3] = relative[y+0*3] * rotation[0+x*3]
982  + relative[y+1*3] * rotation[1+x*3]
983  + relative[y+2*3] * rotation[2+x*3];
984  memcpy(rotation, new_rotation, 9 * sizeof(IceTDouble) );
985  json_object_del( last, "rotation relative" );
986  }
987  if (json_array_size( js = json_object_get(last, "rotation axis") ) == 4)
988  {
989  IceTDouble relative[9];
990  IceTDouble x = json_number_value( json_array_get( js, 0 ) );
991  IceTDouble y = json_number_value( json_array_get( js, 1 ) );
992  IceTDouble z = json_number_value( json_array_get( js, 2 ) );
993  IceTDouble rad = json_number_value( json_array_get( js, 3 ) );
994  IceTDouble s = sin( rad * M_PI / 180.0);
995  IceTDouble c = cos( rad * M_PI / 180.0);
996  IceTDouble l = sqrt( x * x + y * y + z * z);
997  if ( l != 0.0 )
998  {
999  add_modelview = true;
1000  send_rotation = true;
1001  x /= l;
1002  y /= l;
1003  z /= l;
1004  relative[0] = c + x * x * ( 1 - c );
1005  relative[3] = x * y * ( 1 - c ) - z * s;
1006  relative[6] = x * z * ( 1 - c ) + y * s;
1007  relative[1] = y * x * ( 1 - c ) + z * s;
1008  relative[4] = c + y * y * ( 1 - c );
1009  relative[7] = y * z * ( 1 - c ) - x * s;
1010  relative[2] = z * x * ( 1 - c ) - y * s;
1011  relative[5] = z * y * ( 1 - c ) + x * s;
1012  relative[8] = c + z * z * ( 1 - c );
1013  IceTDouble new_rotation[9];
1014  for (isaac_int x = 0; x < 3; x++)
1015  for (isaac_int y = 0; y < 3; y++)
1016  new_rotation[y+x*3] = relative[y+0*3] * rotation[0+x*3]
1017  + relative[y+1*3] * rotation[1+x*3]
1018  + relative[y+2*3] * rotation[2+x*3];
1019  memcpy(rotation, new_rotation, 9 * sizeof(IceTDouble) );
1020  }
1021  json_object_del( last, "rotation axis" );
1022  }
1023  if (json_array_size( js = json_object_get(last, "position absolute") ) == 3)
1024  {
1025  add_modelview = true;
1026  send_look_at = true;
1027  json_array_foreach(js, index, value)
1028  look_at[index] = json_number_value( value );
1029  json_object_del( last, "position absolute" );
1030  }
1031  if (json_array_size( js = json_object_get(last, "position relative") ) == 3)
1032  {
1033  add_modelview = true;
1034  send_look_at = true;
1035  IceTDouble add[3];
1036  json_array_foreach(js, index, value)
1037  add[index] = json_number_value( value );
1038  IceTDouble add_p[3] =
1039  {
1040  rotation[0] * add[0] + rotation[1] * add[1] + rotation[2] * add[2],
1041  rotation[3] * add[0] + rotation[4] * add[1] + rotation[5] * add[2],
1042  rotation[6] * add[0] + rotation[7] * add[1] + rotation[8] * add[2]
1043  };
1044  look_at[0] += add_p[0];
1045  look_at[1] += add_p[1];
1046  look_at[2] += add_p[2];
1047  json_object_del( last, "position relative" );
1048  }
1049  if ( js = json_object_get(last, "distance absolute") )
1050  {
1051  add_modelview = true;
1052  send_distance = true;
1053  distance = json_number_value( js );
1054  json_object_del( last, "distance absolute" );
1055  }
1056  if ( js = json_object_get(last, "distance relative") )
1057  {
1058  add_modelview = true;
1059  send_distance = true;
1060  distance += json_number_value( js );
1061  json_object_del( last, "distance relative" );
1062  }
1063  //Giving the Controller the chance to grep for controller specific messages:
1064  if ( controller.updateProjection( projection, framebuffer_size, last ) )
1065  {
1066  redraw = true;
1067  send_projection = true;
1068  json_t *matrix;
1069  json_object_set_new( message, "projection", matrix = json_array() );
1070  ISAAC_JSON_ADD_MATRIX(matrix,projection,16 * TController::pass_count)
1071  }
1072  mergeJSON(message,last);
1073  json_decref( last );
1074  }
1075  if (add_modelview)
1076  {
1077  redraw = true;
1078  updateModelview();
1079  json_t *matrix;
1080  json_object_set_new( message, "modelview", matrix = json_array() );
1081  ISAAC_JSON_ADD_MATRIX(matrix,modelview,16)
1082  }
1083  char* buffer = json_dumps( message, 0 );
1084  strncpy( message_buffer, buffer, ISAAC_MAX_RECEIVE-1);
1085  message_buffer[ ISAAC_MAX_RECEIVE-1 ] = 0;
1086  free(buffer);
1087  int l = strlen( message_buffer );
1088  MPI_Bcast( &l, 1, MPI_INT, master, mpi_world);
1089  MPI_Bcast( message_buffer, l, MPI_CHAR, master, mpi_world);
1090  }
1091  else //The others just get the message
1092  {
1093  int l;
1094  MPI_Bcast( &l, 1, MPI_INT, master, mpi_world);
1095  MPI_Bcast( message_buffer, l, MPI_CHAR, master, mpi_world);
1096  message_buffer[l] = 0;
1097  message = json_loads(message_buffer, 0, NULL);
1098  }
1099 
1100  json_t* js;
1101  size_t index;
1102  json_t *value;
1103 
1104  //search for requests for all ranks
1105  if ( js = json_object_get(message, "request") )
1106  {
1107  const char* target = json_string_value( js );
1108  if ( strcmp( target, "redraw" ) == 0 )
1109  redraw = true;
1110  if ( strcmp( target, "minmax" ) == 0 )
1111  send_minmax = true;
1112  }
1113 
1114  //Scene set?
1115  if (json_array_size( js = json_object_get(message, "projection") ) == 16 * TController::pass_count)
1116  {
1117  redraw = true;
1118  send_projection = true;
1119  json_array_foreach(js, index, value)
1120  projection[index] = json_number_value( value );
1121  }
1122  if (rank!= master && json_array_size( js = json_object_get(message, "modelview") ) == 16)
1123  {
1124  redraw = true;
1125  json_array_foreach(js, index, value)
1126  modelview[index] = json_number_value( value );
1127  }
1128  if (json_array_size( js = json_object_get(message, "transfer points") ) )
1129  {
1130  redraw = true;
1131  json_array_foreach(js, index, value)
1132  {
1133  transfer_h.description[index].clear();
1134  size_t index_2;
1135  json_t *element;
1136  json_array_foreach(value, index_2, element)
1137  {
1138  transfer_h.description[index].insert( std::pair< isaac_uint, isaac_float4> (
1139  isaac_uint( json_number_value( json_object_get( element, "value" ) ) ), {
1140  isaac_float( json_number_value( json_object_get( element, "r" ) ) ),
1141  isaac_float( json_number_value( json_object_get( element, "g" ) ) ),
1142  isaac_float( json_number_value( json_object_get( element, "b" ) ) ),
1143  isaac_float( json_number_value( json_object_get( element, "a" ) ) ) } ) );
1144  }
1145  }
1146  updateTransfer();
1147  send_transfer = true;
1148  }
1149  if ( js = json_object_get(message, "interpolation") )
1150  {
1151  redraw = true;
1152  interpolation = json_boolean_value ( js );
1153  send_interpolation = true;
1154  }
1155  if ( js = json_object_get(message, "step") )
1156  {
1157  redraw = true;
1158  step = json_number_value ( js );
1159  if (step < isaac_float(0.01))
1160  step = isaac_float(0.01);
1161  send_step = true;
1162  }
1163  if ( js = json_object_get(message, "iso surface") )
1164  {
1165  redraw = true;
1166  iso_surface = json_boolean_value ( js );
1167  send_iso_surface = true;
1168  }
1169  if (json_array_size( js = json_object_get(message, "functions") ) )
1170  {
1171  redraw = true;
1172  json_array_foreach(js, index, value)
1173  functions[index].source = std::string(json_string_value(value));
1174  updateFunctions();
1175  send_functions = true;
1176  }
1177  if (json_array_size( js = json_object_get(message, "weight") ) )
1178  {
1179  redraw = true;
1180  json_array_foreach(js, index, value)
1181  source_weight.value[index] = json_number_value(value);
1182  send_weight = true;
1183  }
1184  if (js = json_object_get(message, "bounding box") )
1185  {
1186  redraw = true;
1187  icet_bounding_box = !icet_bounding_box;
1188  updateBounding( );
1189  }
1190  if (json_array_size( js = json_object_get(message, "background color") ) == 3 )
1191  {
1192  redraw = true;
1193  json_array_foreach(js, index, value)
1194  background_color[index] = json_number_value( value );
1195  for (int pass = 0; pass < TController::pass_count; pass++)
1196  {
1197  icetSetContext( icetContext[pass] );
1198  if (background_color[0] == 0.0f &&
1199  background_color[1] == 0.0f &&
1200  background_color[2] == 0.0f )
1201  icetDisable(ICET_CORRECT_COLORED_BACKGROUND);
1202  else
1203  icetEnable(ICET_CORRECT_COLORED_BACKGROUND);
1204  }
1205  send_background_color = true;
1206  }
1207  if (js = json_object_get(message, "clipping add") )
1208  {
1209  redraw = true;
1210  send_clipping = true;
1211  json_t* position = json_object_get(js, "position");
1212  json_t* normal = json_object_get(js, "normal");
1213  addClipping(
1214  json_number_value( json_array_get( position, 0 ) ),
1215  json_number_value( json_array_get( position, 1 ) ),
1216  json_number_value( json_array_get( position, 2 ) ),
1217  json_number_value( json_array_get( normal, 0 ) ),
1218  json_number_value( json_array_get( normal, 1 ) ),
1219  json_number_value( json_array_get( normal, 2 ) )
1220  );
1221  }
1222  if (js = json_object_get(message, "clipping remove") )
1223  {
1224  redraw = true;
1225  send_clipping = true;
1226  removeClipping( json_integer_value( js ) );
1227  }
1228  if (js = json_object_get(message, "clipping edit") )
1229  {
1230  redraw = true;
1231  send_clipping = true;
1232  json_t* nr = json_object_get(js, "nr");
1233  json_t* position = json_object_get(js, "position");
1234  json_t* normal = json_object_get(js, "normal");
1235  editClipping(
1236  json_integer_value( nr ),
1237  json_number_value( json_array_get( position, 0 ) ),
1238  json_number_value( json_array_get( position, 1 ) ),
1239  json_number_value( json_array_get( position, 2 ) ),
1240  json_number_value( json_array_get( normal, 0 ) ),
1241  json_number_value( json_array_get( normal, 1 ) ),
1242  json_number_value( json_array_get( normal, 2 ) )
1243  );
1244  }
1245 
1246  json_t* metadata = json_object_get( message, "metadata" );
1247  if (metadata)
1248  json_incref(metadata);
1249  json_decref(message);
1250  thr_metaTargets = metaTargets;
1251 
1252  if (send_minmax)
1253  {
1254  isaac_for_each_params( sources, calc_minmax_iterator(), pointer_array, minmax_array, local_minmax_array_d, local_size
1255  #if ISAAC_ALPAKA == 1
1256  ,stream
1257  ,host
1258  #endif
1259  );
1260  if (rank == master)
1261  {
1262  MPI_Reduce( MPI_IN_PLACE, minmax_array.min, boost::mpl::size< TSourceList >::type::value, MPI_FLOAT, MPI_MIN, master, mpi_world);
1263  MPI_Reduce( MPI_IN_PLACE, minmax_array.max, boost::mpl::size< TSourceList >::type::value, MPI_FLOAT, MPI_MAX, master, mpi_world);
1264  }
1265  else
1266  {
1267  MPI_Reduce( minmax_array.min, NULL, boost::mpl::size< TSourceList >::type::value, MPI_FLOAT, MPI_MIN, master, mpi_world);
1268  MPI_Reduce( minmax_array.max, NULL, boost::mpl::size< TSourceList >::type::value, MPI_FLOAT, MPI_MAX, master, mpi_world);
1269  }
1270  }
1271 
1272  for (int pass = 0; pass < TController::pass_count; pass++)
1273  image[pass].opaque_internals = NULL;
1274 
1275  if (redraw)
1276  {
1277  for (int pass = 0; pass < TController::pass_count; pass++)
1278  {
1279  icetSetContext( icetContext[pass] );
1280  //Calc order
1282  //Every rank calculates it's distance to the camera
1283  IceTDouble point[4] =
1284  {
1285  (IceTDouble(position_scaled[0]) + (IceTDouble(local_size_scaled[0]) - IceTDouble(global_size_scaled[0])) / 2.0) / IceTDouble(max_size_scaled/2),
1286  (IceTDouble(position_scaled[1]) + (IceTDouble(local_size_scaled[1]) - IceTDouble(global_size_scaled[1])) / 2.0) / IceTDouble(max_size_scaled/2),
1287  (IceTDouble(position_scaled[2]) + (IceTDouble(local_size_scaled[2]) - IceTDouble(global_size_scaled[2])) / 2.0) / IceTDouble(max_size_scaled/2),
1288  1.0
1289  };
1290  IceTDouble result[4];
1291  mulMatrixVector(result,modelview,point);
1292  float point_distance = sqrt(result[0] * result[0] + result[1] * result[1] + result[2] * result[2]);
1293  //Allgather of the distances
1294  float receive_buffer[numProc];
1295  MPI_Allgather( &point_distance, 1, MPI_FLOAT, receive_buffer, 1, MPI_FLOAT, mpi_world);
1296  //Putting to a std::multimap of {rank, distance}
1297  std::multimap<float, isaac_int, std::less< float > > distance_map;
1298  for (isaac_int i = 0; i < numProc; i++)
1299  distance_map.insert( std::pair<float, isaac_int>( receive_buffer[i], i ) );
1300  //Putting in an array for IceT
1301  IceTInt icet_order_array[numProc];
1302  {
1303  isaac_int i = 0;
1304  for (auto it = distance_map.begin(); it != distance_map.end(); it++)
1305  {
1306  icet_order_array[i] = it->second;
1307  i++;
1308  }
1309  }
1310  icetCompositeOrder( icet_order_array );
1312 
1313  //Drawing
1315  image[pass] = icetDrawFrame(&(projection[pass * 16]),modelview,background_color);
1317  }
1318  }
1319  else
1320  usleep(10000);
1321 
1322  //Message merging
1323  char* buffer = json_dumps( json_root, 0 );
1324  strcpy( message_buffer, buffer );
1325  free(buffer);
1326  if ( metaTargets == META_MERGE )
1327  {
1328  if (rank == master)
1329  {
1330  char receive_buffer[numProc][ISAAC_MAX_RECEIVE];
1331  MPI_Gather( message_buffer, ISAAC_MAX_RECEIVE, MPI_CHAR, receive_buffer, ISAAC_MAX_RECEIVE, MPI_CHAR, master, mpi_world);
1332  for (isaac_int i = 0; i < numProc; i++)
1333  {
1334  if (i == master)
1335  continue;
1336  json_t* js = json_loads(receive_buffer[i], 0, NULL);
1337  mergeJSON( json_root, js );
1338  json_decref( js );
1339  }
1340  }
1341  else
1342  MPI_Gather( message_buffer, ISAAC_MAX_RECEIVE, MPI_CHAR, NULL, 0, MPI_CHAR, master, mpi_world);
1343  }
1344 
1345  #ifdef ISAAC_THREADING
1346  pthread_create(&visualizationThread,NULL,visualizationFunction,NULL);
1347  #else
1348  visualizationFunction(NULL);
1349  #endif
1350  return metadata;
1351  }
1353  {
1355  json_decref( json_root );
1356  if (rank == master)
1357  {
1358  json_root = json_object();
1359  json_object_set_new( json_root, "type", json_string( "exit" ) );
1360  char* buffer = json_dumps( json_root, 0 );
1361  communicator->serverSend(buffer,true,true);
1362  free(buffer);
1363  json_decref( json_root );
1364  }
1365  for (int pass = 0; pass < TController::pass_count; pass++)
1366  icetDestroyContext(icetContext[pass]);
1367  #if ISAAC_ALPAKA == 0
1368  for (int i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
1369  {
1370  if ( pointer_array.pointer[i] )
1371  ISAAC_CUDA_CHECK(cudaFree( pointer_array.pointer[i] ) );
1372  ISAAC_CUDA_CHECK( cudaFree( transfer_d.pointer[i] ) );
1373  free( transfer_h.pointer[i] );
1374  }
1375  ISAAC_CUDA_CHECK(cudaFree( framebuffer ) );
1376  ISAAC_CUDA_CHECK(cudaFree( functor_chain_d ) );
1377  ISAAC_CUDA_CHECK(cudaFree( functor_chain_choose_d ) );
1378  ISAAC_CUDA_CHECK(cudaFree( local_minmax_array_d ) );
1379  #endif
1380  delete communicator;
1381  json_decref(json_init_root);
1382  }
1383  uint64_t getTicksUs()
1384  {
1385  struct timespec ts;
1386  if (clock_gettime(CLOCK_MONOTONIC_RAW,&ts) == 0)
1387  return ts.tv_sec*1000000 + ts.tv_nsec/1000;
1388  return 0;
1389  }
1390  uint64_t kernel_time;
1391  uint64_t merge_time;
1393  uint64_t copy_time;
1394  uint64_t sorting_time;
1395  uint64_t buffer_time;
1396  private:
1397  static void drawCallBack(
1398  const IceTDouble * projection_matrix,
1399  const IceTDouble * modelview_matrix,
1400  const IceTFloat * background_color,
1401  const IceTInt * readback_viewport,
1402  IceTImage result)
1403  {
1404  #if ISAAC_ALPAKA == 1
1405  alpaka::mem::buf::Buf<THost, isaac_float, TFraDim, size_t> inverse_h_buf ( alpaka::mem::buf::alloc<isaac_float, size_t>(myself->host, size_t(16)));
1406  alpaka::mem::buf::Buf<THost, isaac_size_struct< TSimDim::value >, TFraDim, size_t> size_h_buf ( alpaka::mem::buf::alloc<isaac_size_struct< TSimDim::value >, size_t>(myself->host, size_t(1)));
1407  isaac_float* inverse_h = reinterpret_cast<float*>(alpaka::mem::view::getPtrNative(inverse_h_buf));
1408  isaac_size_struct< TSimDim::value >* size_h = reinterpret_cast<isaac_size_struct< TSimDim::value >*>(alpaka::mem::view::getPtrNative(size_h_buf));
1409  #else
1410  isaac_float inverse_h[16];
1412  #endif
1413  IceTDouble inverse[16];
1414  calcInverse(inverse,projection_matrix,modelview_matrix);
1415  for (int i = 0; i < 16; i++)
1416  inverse_h[i] = static_cast<float>(inverse[i]);
1417 
1418  size_h[0].global_size.value.x = myself->global_size[0];
1419  size_h[0].global_size.value.y = myself->global_size[1];
1420  if (TSimDim::value > 2)
1421  size_h[0].global_size.value.z = myself->global_size[2];
1422  size_h[0].position.value.x = myself->position[0];
1423  size_h[0].position.value.y = myself->position[1];
1424  if (TSimDim::value > 2)
1425  size_h[0].position.value.z = myself->position[2];
1426  size_h[0].local_size.value.x = myself->local_size[0];
1427  size_h[0].local_size.value.y = myself->local_size[1];
1428  if (TSimDim::value > 2)
1429  size_h[0].local_size.value.z = myself->local_size[2];
1430  size_h[0].max_global_size = static_cast<float>(ISAAC_MAX(ISAAC_MAX(uint32_t(myself->global_size[0]),uint32_t(myself->global_size[1])),uint32_t(myself->global_size[2])));
1431 
1432  size_h[0].global_size_scaled.value.x = myself->global_size_scaled[0];
1433  size_h[0].global_size_scaled.value.y = myself->global_size_scaled[1];
1434  if (TSimDim::value > 2)
1435  size_h[0].global_size_scaled.value.z = myself->global_size_scaled[2];
1436  size_h[0].position_scaled.value.x = myself->position_scaled[0];
1437  size_h[0].position_scaled.value.y = myself->position_scaled[1];
1438  if (TSimDim::value > 2)
1439  size_h[0].position_scaled.value.z = myself->position_scaled[2];
1440  size_h[0].local_size_scaled.value.x = myself->local_size_scaled[0];
1441  size_h[0].local_size_scaled.value.y = myself->local_size_scaled[1];
1442  if (TSimDim::value > 2)
1443  size_h[0].local_size_scaled.value.z = myself->local_size_scaled[2];
1444  size_h[0].max_global_size_scaled = static_cast<float>(ISAAC_MAX(ISAAC_MAX(uint32_t(myself->global_size_scaled[0]),uint32_t(myself->global_size_scaled[1])),uint32_t(myself->global_size_scaled[2])));
1445 
1446  isaac_float3 isaac_scale =
1447  {
1448  myself->scale[0],
1449  myself->scale[1],
1450  myself->scale[2]
1451  };
1452 
1453  #if ISAAC_ALPAKA == 1
1454  alpaka::vec::Vec<alpaka::dim::DimInt<1u>, size_t> const inverse_d_extent(size_t(16));
1455  auto inverse_d_view(alpaka::mem::view::createStaticDevMemView(&isaac_inverse_d[0u],myself->acc,inverse_d_extent));
1456  alpaka::mem::view::copy(myself->stream, inverse_d_view, inverse_h_buf, size_t(16));
1457 
1458  alpaka::vec::Vec<alpaka::dim::DimInt<1u>, size_t> const size_d_extent(size_t(1));
1459  auto size_d_view(alpaka::mem::view::createStaticDevMemView(&isaac_size_d[0u],myself->acc,size_d_extent));
1460  alpaka::mem::view::copy(myself->stream, size_d_view, size_h_buf, size_t(1));
1461  #else
1462  ISAAC_CUDA_CHECK(cudaMemcpyToSymbol( isaac_inverse_d, inverse_h, 16 * sizeof(float)));
1463  ISAAC_CUDA_CHECK(cudaMemcpyToSymbol( isaac_size_d, size_h, sizeof(isaac_size_struct< TSimDim::value >)));
1464  #endif
1465  IceTUByte* pixels = icetImageGetColorub(result);
1466  ISAAC_START_TIME_MEASUREMENT( kernel, myself->getTicksUs() )
1467  isaac_float4 bg_color =
1468  {
1469  isaac_float(background_color[3]),
1470  isaac_float(background_color[2]),
1471  isaac_float(background_color[1]),
1472  isaac_float(background_color[0])
1473  };
1474  isaac_uint2 framebuffer_start =
1475  {
1476  isaac_uint( readback_viewport[0] ),
1477  isaac_uint( readback_viewport[1] )
1478  };
1479  #if ISAAC_ALPAKA == 1
1481  <
1482  TSimDim,
1483  TSourceList,
1487  mpl::vector<>,
1488  alpaka::mem::buf::Buf<TDevAcc, uint32_t, TFraDim, size_t>,
1489  TTransfer_size,
1490  isaac_float3,
1491  TAccDim,
1492  TAcc,
1493  TStream,
1494  alpaka::mem::buf::Buf<TDevAcc, isaac_functor_chain_pointer_N, TFraDim, size_t>,
1495  boost::mpl::size< TSourceList >::type::value
1496  >
1497  ::call(
1498  myself->stream,
1499  myself->framebuffer,
1500  myself->framebuffer_size,
1501  framebuffer_start,
1502  myself->sources,
1503  myself->step,
1504  bg_color,
1505  myself->transfer_d,
1506  myself->source_weight,
1507  myself->pointer_array,
1508  readback_viewport,
1509  myself->interpolation,
1510  myself->iso_surface,
1511  isaac_scale,
1512  myself->clipping
1513  );
1514  alpaka::wait::wait(myself->stream);
1515  ISAAC_STOP_TIME_MEASUREMENT( myself->kernel_time, +=, kernel, myself->getTicksUs() )
1516  ISAAC_START_TIME_MEASUREMENT( copy, myself->getTicksUs() )
1517  alpaka::mem::view::ViewPlainPtr<THost, uint32_t, TFraDim, size_t> result_buffer((uint32_t*)(pixels), myself->host, alpaka::vec::Vec<TFraDim, size_t>(myself->framebuffer_prod));
1518  alpaka::mem::view::copy(myself->stream, result_buffer, myself->framebuffer, alpaka::vec::Vec<TFraDim, size_t>(myself->framebuffer_prod));
1519  #else
1521  <
1522  TSimDim,
1523  TSourceList,
1524  transfer_d_struct< boost::mpl::size< TSourceList >::type::value >,
1525  source_weight_struct< boost::mpl::size< TSourceList >::type::value >,
1526  pointer_array_struct< boost::mpl::size< TSourceList >::type::value >,
1527  mpl::vector<>,
1528  uint32_t*,
1529  TTransfer_size,
1530  isaac_float3,
1531  boost::mpl::size< TSourceList >::type::value
1532  >
1533  ::call(
1534  myself->framebuffer,
1535  myself->framebuffer_size,
1536  framebuffer_start,
1537  myself->sources,
1538  myself->step,
1539  bg_color,
1540  myself->transfer_d,
1541  myself->source_weight,
1542  myself->pointer_array,
1543  readback_viewport,
1544  myself->interpolation,
1545  myself->iso_surface,
1546  isaac_scale,
1547  myself->clipping
1548  );
1549  ISAAC_CUDA_CHECK(cudaDeviceSynchronize());
1550  ISAAC_STOP_TIME_MEASUREMENT( myself->kernel_time, +=, kernel, myself->getTicksUs() )
1551  ISAAC_START_TIME_MEASUREMENT( copy, myself->getTicksUs() )
1552  ISAAC_CUDA_CHECK(cudaMemcpy((uint32_t*)(pixels), myself->framebuffer, sizeof(uint32_t)*myself->framebuffer_prod, cudaMemcpyDeviceToHost));
1553  #endif
1554  ISAAC_STOP_TIME_MEASUREMENT( myself->copy_time, +=, copy, myself->getTicksUs() )
1555  }
1556 
1557  static void* visualizationFunction(void* dummy)
1558  {
1559  //Message sending
1560  if (myself->rank == myself->master)
1561  {
1562  json_object_set_new( myself->json_root, "type", json_string( "period" ) );
1563  json_object_set_new( myself->json_root, "meta nr", json_integer( myself->metaNr ) );
1564 
1565  json_t *matrix;
1566  if ( myself->send_projection )
1567  {
1568  json_object_set_new( myself->json_root, "projection", matrix = json_array() );
1569  ISAAC_JSON_ADD_MATRIX(matrix,myself->projection,16 * TController::pass_count)
1570  json_object_set( myself->json_init_root, "projection", matrix );
1571  myself->send_init_json = true;
1572  }
1573  if ( myself->send_look_at )
1574  {
1575  json_object_set_new( myself->json_root, "position", matrix = json_array() );
1576  ISAAC_JSON_ADD_MATRIX(matrix,myself->look_at,3)
1577  json_object_set( myself->json_init_root, "position", matrix );
1578  myself->send_init_json = true;
1579  }
1580  if ( myself->send_rotation )
1581  {
1582  json_object_set_new( myself->json_root, "rotation", matrix = json_array() );
1583  ISAAC_JSON_ADD_MATRIX(matrix, myself->rotation,9)
1584  json_object_set( myself->json_init_root, "rotation", matrix );
1585  myself->send_init_json = true;
1586  }
1587  if ( myself->send_distance )
1588  {
1589  json_object_set_new( myself->json_root, "distance", json_real( myself->distance ) );
1590  json_object_set_new( myself->json_init_root, "distance", json_real( myself->distance ) );
1591  myself->send_init_json = true;
1592  }
1593  if ( myself->send_transfer )
1594  {
1595  json_object_set_new( myself->json_root, "transfer array", matrix = json_array() );
1596  for (size_t i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
1597  {
1598  json_t* transfer = json_array();
1599  json_array_append_new( matrix, transfer );
1600  for (size_t j = 0; j < TTransfer_size; j++)
1601  {
1602  json_t* color = json_array();
1603  json_array_append_new( transfer, color );
1604  json_array_append_new( color, json_integer( isaac_uint( myself->transfer_h.pointer[i][j].x * isaac_float(255) ) ) );
1605  json_array_append_new( color, json_integer( isaac_uint( myself->transfer_h.pointer[i][j].y * isaac_float(255) ) ) );
1606  json_array_append_new( color, json_integer( isaac_uint( myself->transfer_h.pointer[i][j].z * isaac_float(255) ) ) );
1607  json_array_append_new( color, json_integer( isaac_uint( myself->transfer_h.pointer[i][j].w * isaac_float(255) ) ) );
1608  }
1609  }
1610  json_object_set_new( myself->json_root, "transfer points", matrix = json_array() );
1611  for (size_t i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
1612  {
1613  json_t* points = json_array();
1614  json_array_append_new( matrix, points );
1615  for(auto it = myself->transfer_h.description[i].begin(); it != myself->transfer_h.description[i].end(); it++)
1616  {
1617  json_t* p = json_object();
1618  json_array_append_new( points, p);
1619  json_object_set_new(p, "value", json_integer( it->first ) );
1620  json_object_set_new(p, "r", json_real( it->second.x ) );
1621  json_object_set_new(p, "g", json_real( it->second.y ) );
1622  json_object_set_new(p, "b", json_real( it->second.z ) );
1623  json_object_set_new(p, "a", json_real( it->second.w ) );
1624  }
1625  }
1626  }
1627  if ( myself->send_functions )
1628  {
1629  json_object_set_new( myself->json_root, "functions", matrix = json_array() );
1630  for (size_t i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
1631  {
1632  json_t* f = json_object();
1633  json_array_append_new( matrix, f );
1634  json_object_set_new(f, "source", json_string( myself->functions[i].source.c_str() ) );
1635  json_object_set_new(f, "error", json_integer( myself->functions[i].error_code ) );
1636  }
1637  }
1638  if ( myself->send_weight )
1639  {
1640  json_object_set_new( myself->json_root, "weight", matrix = json_array() );
1641  for (size_t i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
1642  json_array_append_new( matrix, json_real( myself->source_weight.value[i] ) );
1643  }
1644  if ( myself->send_interpolation )
1645  {
1646  json_object_set_new( myself->json_root, "interpolation", json_boolean( myself->interpolation ) );
1647  json_object_set_new( myself->json_init_root, "interpolation", json_boolean( myself->interpolation ) );
1648  myself->send_init_json = true;
1649  }
1650  if ( myself->send_step )
1651  {
1652  json_object_set_new( myself->json_root, "step", json_real( myself->step ) );
1653  json_object_set_new( myself->json_init_root, "step", json_boolean( myself->step ) );
1654  myself->send_init_json = true;
1655  }
1656  if ( myself->send_iso_surface )
1657  {
1658  json_object_set_new( myself->json_root, "iso surface", json_boolean( myself->iso_surface ) );
1659  json_object_set_new( myself->json_init_root, "iso surface", json_boolean( myself->iso_surface ) );
1660  myself->send_init_json = true;
1661  }
1662  if ( myself->send_minmax )
1663  {
1664  json_object_set_new( myself->json_root, "minmax", matrix = json_array() );
1665  for (size_t i = 0; i < boost::mpl::size< TSourceList >::type::value; i++)
1666  {
1667  json_t* v = json_object();
1668  json_array_append_new( matrix, v );
1669  json_object_set_new(v, "min", json_real( myself->minmax_array.min[i] ) );
1670  json_object_set_new(v, "max", json_real( myself->minmax_array.max[i] ) );
1671  }
1672  }
1673  if ( myself->send_background_color )
1674  {
1675  json_object_set_new( myself->json_root, "background color", matrix = json_array() );
1676  for (size_t i = 0; i < 3; i++)
1677  json_array_append_new( matrix, json_real( myself->background_color[i] ) );
1678  json_object_set( myself->json_init_root, "background color", matrix );
1679  myself->send_init_json = true;
1680  }
1681  if ( myself->send_clipping )
1682  {
1683  json_object_set_new( myself->json_root, "clipping", matrix = json_array() );
1684  for (size_t i = 0; i < myself->clipping.count; i++)
1685  {
1686  json_t* f = json_object();
1687  json_array_append_new( matrix, f );
1688  json_t* inner = json_array();
1689  json_object_set_new(f, "position", inner );
1690  json_array_append_new( inner, json_real( myself->clipping.elem[i].position.x ) );
1691  json_array_append_new( inner, json_real( myself->clipping.elem[i].position.y ) );
1692  json_array_append_new( inner, json_real( myself->clipping.elem[i].position.z ) );
1693  inner = json_array();
1694  json_object_set_new(f, "normal", inner );
1695  json_array_append_new( inner, json_real( myself->clipping_saved_normals[i].x ) );
1696  json_array_append_new( inner, json_real( myself->clipping_saved_normals[i].y ) );
1697  json_array_append_new( inner, json_real( myself->clipping_saved_normals[i].z ) );
1698  }
1699  }
1700  myself->controller.sendFeedback( myself->json_root, myself->send_controller );
1701  if ( myself->send_init_json )
1702  json_object_set( myself->json_root,"init",myself->json_init_root );
1703  char* buffer = json_dumps( myself->json_root, 0 );
1704  myself->communicator->serverSend(buffer);
1705  free(buffer);
1706  }
1707  json_decref( myself->json_root );
1708  myself->recreateJSON();
1709 
1710  //Sending video
1711  ISAAC_START_TIME_MEASUREMENT( video_send, myself->getTicksUs() )
1712  if (myself->communicator)
1713  {
1714  if (myself->image[0].opaque_internals)
1715  myself->communicator->serverSendFrame(myself->compositor.doCompositing(myself->image),myself->compbuffer_size.x,myself->compbuffer_size.y,4);
1716  else
1717  {
1718  myself->communicator->serverSend(NULL,false,true);
1719  }
1720  }
1721  ISAAC_STOP_TIME_MEASUREMENT( myself->video_send_time, +=, video_send, myself->getTicksUs() )
1722  myself->metaNr++;
1723  return 0;
1724  }
1725  void recreateJSON()
1726  {
1727  json_root = json_object();
1728  json_meta_root = json_object();
1729  json_object_set_new( json_root, "metadata", json_meta_root );
1730  }
1731  void updateModelview()
1732  {
1733  IceTDouble look_at_m[16];
1734  ISAAC_SET_IDENTITY(4,look_at_m)
1735  look_at_m[12] = look_at[0];
1736  look_at_m[13] = look_at[1];
1737  look_at_m[14] = look_at[2];
1738 
1739  IceTDouble rotation_m[16];
1740  for (isaac_int x = 0; x < 4; x++)
1741  for (isaac_int y = 0; y < 4; y++)
1742  {
1743  if (x < 3 && y < 3)
1744  rotation_m[x+y*4] = rotation[x+y*3];
1745  else
1746  if (x!=3 || y!=3)
1747  rotation_m[x+y*4] = 0.0;
1748  else
1749  rotation_m[x+y*4] = 1.0;
1750  }
1751 
1752  IceTDouble distance_m[16];
1753  ISAAC_SET_IDENTITY(4,distance_m)
1754  distance_m[14] = distance;
1755 
1756  IceTDouble temp[16];
1757 
1758  mulMatrixMatrix( temp, rotation_m, look_at_m );
1759  mulMatrixMatrix( modelview, distance_m, temp );
1760  }
1761  #if ISAAC_ALPAKA == 1
1762  THost host;
1763  TDevAcc acc;
1764  TStream stream;
1765  #endif
1766  std::string name;
1767  std::string server_url;
1768  isaac_uint server_port;
1769  isaac_size2 framebuffer_size;
1770  isaac_size2 compbuffer_size;
1771  #if ISAAC_ALPAKA == 1
1772  alpaka::vec::Vec<TFraDim, size_t> framebuffer_prod;
1773  alpaka::mem::buf::Buf<TDevAcc, uint32_t, TFraDim, size_t> framebuffer;
1774  alpaka::mem::buf::Buf<TDevAcc, isaac_functor_chain_pointer_N, TFraDim, size_t> functor_chain_d;
1775  alpaka::mem::buf::Buf<TDevAcc, isaac_functor_chain_pointer_N, TFraDim, size_t> functor_chain_choose_d;
1776  alpaka::mem::buf::Buf<TDevAcc, minmax_struct, TFraDim, size_t> local_minmax_array_d;
1777  #else
1778  size_t framebuffer_prod;
1779  isaac_uint* framebuffer;
1780  isaac_functor_chain_pointer_N* functor_chain_d;
1781  isaac_functor_chain_pointer_N* functor_chain_choose_d;
1782  minmax_struct* local_minmax_array_d;
1783  #endif
1784  TDomainSize global_size;
1785  TDomainSize local_size;
1786  TDomainSize position;
1787  std::vector<size_t> global_size_scaled;
1788  std::vector<size_t> local_size_scaled;
1789  std::vector<size_t> position_scaled;
1790  MPI_Comm mpi_world;
1791  IceTDouble projection[16 * TController::pass_count];
1792  IceTDouble look_at[3];
1793  IceTDouble rotation[9];
1794  IceTDouble distance;
1795  bool send_look_at;
1796  bool send_rotation;
1797  bool send_distance;
1798  bool send_projection;
1799  bool send_transfer;
1800  bool send_interpolation;
1801  bool send_step;
1802  bool send_iso_surface;
1803  bool send_functions;
1804  bool send_weight;
1805  bool send_minmax;
1806  bool send_background_color;
1807  bool send_clipping;
1808  bool send_controller;
1809  bool send_init_json;
1810  bool interpolation;
1811  bool iso_surface;
1812  bool icet_bounding_box;
1813  isaac_float step;
1814  IceTDouble modelview[16];
1815  IsaacCommunicator* communicator;
1816  json_t *json_root;
1817  json_t *json_init_root;
1818  json_t *json_meta_root;
1819  isaac_int rank;
1820  isaac_int master;
1821  isaac_int numProc;
1822  isaac_uint metaNr;
1823  TSourceList& sources;
1824  IceTContext icetContext[TController::pass_count];
1825  IsaacVisualizationMetaEnum thr_metaTargets;
1826  pthread_t visualizationThread;
1827  #if ISAAC_ALPAKA == 1
1828  std::vector< alpaka::mem::buf::Buf<TDevAcc, isaac_float4, TTexDim, size_t> > transfer_d_buf;
1829  std::vector< alpaka::mem::buf::Buf< THost, isaac_float4, TTexDim, size_t> > transfer_h_buf;
1830  std::vector< alpaka::mem::buf::Buf< TDevAcc, isaac_float, TFraDim, size_t> > pointer_array_alpaka;
1831  #endif
1837  const static size_t transfer_size = TTransfer_size;
1838  functions_struct functions[boost::mpl::size< TSourceList >::type::value];
1839  size_t max_size;
1840  size_t max_size_scaled;
1841  IceTFloat background_color[4];
1842  static IsaacVisualization *myself;
1843  TScale scale;
1844  clipping_struct clipping;
1845  isaac_float3 clipping_saved_normals[ISAAC_MAX_CLIPPING];
1846  TController controller;
1847  TCompositor compositor;
1848  IceTImage image[TController::pass_count];
1849 };
1850 
1851 #if ISAAC_ALPAKA == 1
1852  template <typename THost,typename TAcc,typename TStream,typename TAccDim,typename TSimDim, typename TSourceList, typename TDomainSize, size_t TTransfer_size,typename TScale,typename TController,typename TCompositor>
1854 #else
1855  template <typename TSimDim, typename TSourceList, typename TDomainSize, size_t TTransfer_size,typename TScale,typename TController,typename TCompositor>
1857 #endif
1858 
1859 } //namespace isaac;
float isaac_float
Definition: isaac_types.hpp:25
#define ISAAC_HOST_INLINE
void calcInverse(IceTDouble *inv, const IceTDouble *projection, const IceTDouble *modelview)
json_t * getJsonMetaRoot()
Definition: isaac.hpp:860
isaac_float4 getHSVA(isaac_float h, isaac_float s, isaac_float v, isaac_float a)
isaac_size_dim< simdim > local_size
void addClipping(isaac_float px, isaac_float py, isaac_float pz, isaac_float nx, isaac_float ny, isaac_float nz)
Definition: isaac.hpp:651
ISAAC_CONSTANT isaac_functor_chain_pointer_N isaac_function_chain_d[ISAAC_MAX_SOURCES]
IsaacVisualizationMetaEnum
ISAAC_HOST_INLINE void operator()(const int I, const TSource &source, TArray &pointer_array, const TLocalSize &local_size) const
Definition: isaac.hpp:175
Definition: isaac.hpp:60
void mergeJSON(json_t *result, json_t *candidate)
isaac_size_dim< simdim > position_scaled
isaac_size_dim< simdim > local_size_scaled
isaac_size_dim< simdim > global_size_scaled
ISAAC_HOST_INLINE void operator()(const int I, const TSource &source, TArray &pointer_array, TMinmax &minmax, TLocalMinmax &local_minmax, TLocalSize &local_size) const
Definition: isaac.hpp:301
#define ISAAC_MAX_SOURCES
int32_t isaac_int
Definition: isaac_types.hpp:26
#define ISAAC_STOP_TIME_MEASUREMENT(result, operand, unique_name, time_function)
void updatePosition(const TDomainSize position)
Definition: isaac.hpp:692
#define ISAAC_FUNCTOR_COUNT
ISAAC_HOST_INLINE void operator()(const int I, const TFunctor &f, TJsonRoot &jsonRoot) const
Definition: isaac.hpp:109
#define ISAAC_FUNCTOR_COMPLEX
fus::list< IsaacFunctorIdem,IsaacFunctorAdd,IsaacFunctorMul,IsaacFunctorLength,IsaacFunctorPow,IsaacFunctorSum > IsaacFunctorPool
#define ISAAC_JSON_ADD_MATRIX(array, matrix, count)
ISAAC_CONSTANT isaac_float4 isaac_parameter_d[ISAAC_MAX_SOURCES *ISAAC_MAX_FUNCTORS]
uint32_t isaac_uint
Definition: isaac_types.hpp:27
#define ISAAC_MAX_FUNCTORS
void updateLocalSize(const TDomainSize local_size)
Definition: isaac.hpp:699
#define ISAAC_PROTOCOL_VERSION_MAJOR
void mulMatrixMatrix(IceTDouble *result, const IceTDouble *matrix1, const IceTDouble *matrix2)
__global__ void fillFunctorChainPointerKernel(isaac_functor_chain_pointer_N *const functor_chain_d)
ISAAC_HOST_INLINE void operator()(const int I, const TSource &source, const TFunctions &functions, TDest &dest) const
Definition: isaac.hpp:146
json_t * doVisualization(const IsaacVisualizationMetaEnum metaTargets=META_MASTER, void *pointer=NULL, bool redraw=true)
Definition: isaac.hpp:881
ISAAC_NO_HOST_DEVICE_WARNING ISAAC_HOST_DEVICE_INLINE void isaac_for_each_params(Sequence &seq, F const &f, P &... p)
#define ISAAC_SET_IDENTITY(size, matrix)
#define ISAAC_DEFAULT_WEIGHT
#define ISAAC_MAX
isaac_float(* isaac_functor_chain_pointer_N)(void *, isaac_int)
ISAAC_CONSTANT isaac_size_struct< 3 > isaac_size_d[1]
ISAAC_CONSTANT isaac_float isaac_inverse_d[16]
void mulMatrixVector(IceTDouble *result, const IceTDouble *matrix, const IceTDouble *vector)
ISAAC_HOST_INLINE void operator()(const int I, TSource &source, TArray &pointer_array, const TLocalSize &local_size, const TWeight &weight, const TPointer &pointer) const
Definition: isaac.hpp:222
#define ISAAC_CUDA_CHECK(call)
bool editClipping(isaac_uint nr, isaac_float px, isaac_float py, isaac_float pz, isaac_float nx, isaac_float ny, isaac_float nz)
Definition: isaac.hpp:626
#define ISAAC_MAX_RECEIVE
#define ISAAC_START_TIME_MEASUREMENT(unique_name, time_function)
#define ISAAC_MAX_CLIPPING
void removeClipping(isaac_uint nr)
Definition: isaac.hpp:656
#define ISAAC_DEFAULT_STEP
__global__ void updateFunctorChainPointerKernel(isaac_functor_chain_pointer_N *const functor_chain_choose_d, isaac_functor_chain_pointer_N const *const functor_chain_d, TDest dest)
int init(CommunicatorSetting communicatorBehaviour=ReturnAtError)
Definition: isaac.hpp:865
isaac_size_dim< simdim > global_size
#define ISAAC_PROTOCOL_VERSION_MINOR
ISAAC_HOST_INLINE void operator()(const int I, const TSource &s, TJsonRoot &jsonRoot) const
Definition: isaac.hpp:93
#define ISAAC_WAIT_VISUALIZATION
void setJpegQuality(isaac_uint jpeg_quality)
Definition: isaac.hpp:620
ISAAC_HOST_INLINE void operator()(const int I, TFunctor &f, const TName &name, TValue &value, TFound &found) const
Definition: isaac.hpp:128
#define ISAAC_GUARD_SIZE
isaac_size_dim< simdim > position
IsaacVisualization(const std::string name, const isaac_int master, const std::string server_url, const isaac_uint server_port, const isaac_size2 framebuffer_size, const TDomainSize global_size, const TDomainSize local_size, const TDomainSize position, TSourceList &sources, TScale scale)
Definition: isaac.hpp:376