// Derived from VTK/Examples/Cxx/Medical2.cxx // The example reads a volume dataset, extracts two isosurfaces that // represent the skin and bone, and then displays them. // // Modified heavily by Peter Boyle to display lattice field theory data as movies and compare multiple files #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #define MPEG #ifdef MPEG #include #endif #include #include #include #include #include #include #include #define USE_FLYING_EDGES #ifdef USE_FLYING_EDGES #include typedef vtkFlyingEdges3D isosurface; #else #include typedef vtkMarchingCubes isosurface; #endif int mpeg = 0 ; int framerate = 10; template void readFile(T& out, std::string const fname){ Grid::emptyUserRecord record; Grid::ScidacReader RD; RD.open(fname); RD.readScidacFieldRecord(out,record); RD.close(); } using namespace Grid; class FrameUpdater : public vtkCallbackCommand { public: FrameUpdater() { ffile=0; TimerCount = 0; xoff = 0; t = 0; imageData = nullptr; timerId = 0; maxCount = -1; old_file=-1; } static FrameUpdater* New() { FrameUpdater* cb = new FrameUpdater; cb->TimerCount = 0; return cb; } // // Must map a x,y,z + frame index into // i) a d-dimensional site Coordinate // ii) a file name // Need a: // loop_ranges // sum_ranges // loop_vol -- map loop_idx -> loop_coor // sum_vol -- map sum_idx -> sum_coor with Lexicographic // /* * Just set this up */ int old_file ; // Cache, avoid reread Coordinate latt; Coordinate xyz_dims ; // List lattice dimensions corresponding to xyz_dims displayed Coordinate xyz_ranges ; // 3-vector Coordinate g_xyz_ranges; // Nd-vector uint64_t xyz_vol ; Coordinate loop_dims; // List lattice dimensions put into movie time Coordinate loop_ranges; // movie time ranges uint64_t loop_vol; Coordinate sum_dims; // List lattice dimensions summed Coordinate sum_ranges; // summation ranges uint64_t sum_vol; Coordinate slice_dims; // List slice dimensions Coordinate Slice; std::vector files; // file list that is looped over int Nd; GridBase *grid; Grid::LatticeComplexD *grid_data; void SetGrid(GridBase *_grid) { grid = _grid; Nd=grid->Nd(); latt = grid->GlobalDimensions(); grid_data = new Grid::LatticeComplexD(grid); } void SetFiles(std::vector list) { files = list; old_file = -1; } void SetSlice(Coordinate _Slice) { Slice = _Slice;} // Offset / skew for lattice coords void SetSumDimensions (Coordinate _SumDims ) { sum_ranges=Coordinate(Nd); sum_dims = _SumDims; // 1 hot for dimensions summed sum_vol = 1; for(int d=0;d1) is_slice = 0; if(loop_dims[d]) is_slice = 0; if(sum_dims[d] ) is_slice = 0; if(is_slice) _slice_dims.push_back(d); } slice_dims = _slice_dims; std::cout << " Setting Slice Dimensions to "<GlobalDimensions(); if ( vtkCommand::KeyPressEvent == eventId ) { vtkRenderWindowInteractor* iren = static_cast(caller); std::string key = iren->GetKeySym(); std::cout << "Pressed: " << key << std::endl; if (slice_dims.size()>0) { int vert = slice_dims[slice_dims.size()-1]; int horz = slice_dims[0]; if ( key == "Up" ) { Slice[vert] = (Slice[vert]+1)%latt[vert]; } if ( key == "Down" ) { Slice[vert] = (Slice[vert]+latt[vert]-1)%latt[vert]; } if ( key == "Right" ) { Slice[horz] = (Slice[horz]+1)%latt[horz]; } if ( key == "Left" ) { Slice[horz] = (Slice[horz]+latt[horz]-1)%latt[horz]; } } if ( key == "greater" ) { ffile = (ffile + 1) % files.size(); } if ( key == "less" ) { ffile = (ffile - 1 + files.size()) % files.size(); } std::cout <<"Slice " <TimerCount / loop_vol) + ffile )%files.size(); if ( file != old_file ) { readFile(*grid_data,files[file]); old_file = file; } RealD max, min, max_abs,min_abs; Coordinate max_site; Coordinate min_site; Coordinate max_abs_site; Coordinate min_abs_site; for(int idx=0;idxgSites();idx++){ Coordinate site; Lexicographic::CoorFromIndex (site,idx,latt); RealD val=real(peekSite(*grid_data,site)); if (idx==0){ max = min = val; max_abs = min_abs = fabs(val); max_site = site; min_site = site; min_abs_site = site; max_abs_site = site; } else { if ( val > max ) { max=val; max_site = site; } if ( fabs(val) > max_abs ) { max_abs=fabs(val); max_abs_site = site; } if ( val < min ) { min=val; min_site = site; } if ( fabs(val) < min_abs ) { min_abs=fabs(val); min_abs_site = site; } } } std::cout << " abs_max "<TimerCount % loop_vol; Coordinate loop_coor; Lexicographic::CoorFromIndex (loop_coor,loop_idx,loop_ranges); // Loop over xyz sites Coordinate xyz_coor(3); Coordinate g_xyz_coor(Nd); Coordinate sum_coor(Nd); for(uint64_t xyz = 0 ; xyz< xyz_vol; xyz++){ Lexicographic::CoorFromIndex (xyz_coor,xyz,xyz_ranges); Lexicographic::CoorFromIndex (g_xyz_coor,xyz,g_xyz_ranges); RealD sum_value = 0.0; for(uint64_t sum_idx = 0 ; sum_idx< sum_vol; sum_idx++){ Lexicographic::CoorFromIndex (sum_coor,sum_idx,sum_ranges); Coordinate site(Nd); for(int d=0;dSetScalarComponentFromDouble(xyz_coor[0],xyz_coor[1],xyz_coor[2],0,sum_value); } imageData->Modified(); std::stringstream ss; ss<< files[file] <<"\nSlice "<SetInput(ss.str().c_str()); vtkRenderWindowInteractor* iren = dynamic_cast(caller); iren->GetRenderWindow()->Render(); } if ( vtkCommand::TimerEvent == eventId ) { ++this->TimerCount; std::cout << " This was a timer event count "<TimerCount << std::endl; } if (this->TimerCount >= this->maxCount) { vtkRenderWindowInteractor* iren = dynamic_cast(caller); if (this->timerId > -1) { iren->DestroyTimer(this->timerId); } } } private: int TimerCount; int ffile; int xoff; int t; public: vtkImageData* imageData = nullptr; vtkTextActor* text = nullptr; vtkFFMPEGWriter *writer = nullptr; int timerId ; int maxCount ; double rms; isosurface * posExtractor; isosurface * negExtractor; }; class SliderCallback : public vtkCommand { public: static SliderCallback* New() { return new SliderCallback; } virtual void Execute(vtkObject* caller, unsigned long eventId, void* callData) { vtkSliderWidget *sliderWidget = vtkSliderWidget::SafeDownCast(caller); if (sliderWidget) { contour = ((vtkSliderRepresentation *)sliderWidget->GetRepresentation())->GetValue(); } fu->posExtractor->SetValue(0, SliderCallback::contour*fu->rms); fu->negExtractor->SetValue(0, -SliderCallback::contour*fu->rms); fu->posExtractor->Modified(); fu->negExtractor->Modified(); } public: static double contour; FrameUpdater * fu; }; FrameUpdater * KBfu; void KeypressCallbackFunction(vtkObject* caller, long unsigned int eventId, void* clientData, void* callData) { std::cout << "Keypress callback" << std::endl; vtkRenderWindowInteractor* iren = static_cast(caller); std::cout << "Pressed: " << iren->GetKeySym() << std::endl; // imageData->Modified(); } double SliderCallback::contour; int main(int argc, char* argv[]) { using namespace Grid; Grid_init(&argc, &argv); GridLogLayout(); auto latt_size = GridDefaultLatt(); auto simd_layout = GridDefaultSimd(latt_size.size(), vComplex::Nsimd()); auto mpi_layout = GridDefaultMpi(); GridCartesian Grid(latt_size, simd_layout, mpi_layout); double default_contour = 1.0; std::string arg; std::cout << argc << " command Line arguments "< file_list({ "file1", "file2", "file3", "file4", "file5", "file6", "file7", "file8" }); if( GridCmdOptionExists(argv,argv+argc,"--files") ){ arg=GridCmdOptionPayload(argv,argv+argc,"--files"); GridCmdOptionCSL(arg, file_list); } #ifdef MPEG if( GridCmdOptionExists(argv,argv+argc,"--mpeg") ){ mpeg = 1; } #endif if( GridCmdOptionExists(argv,argv+argc,"--fps") ){ arg=GridCmdOptionPayload(argv,argv+argc,"--fps"); GridCmdOptionInt(arg,framerate); } if( GridCmdOptionExists(argv,argv+argc,"--isosurface") ){ arg=GridCmdOptionPayload(argv,argv+argc,"--isosurface"); GridCmdOptionFloat(arg,default_contour); } for(int c=0;c colors; std::array posColor{{240, 184, 160, 255}}; colors->SetColor("posColor", posColor.data()); std::array bkg{{51, 77, 102, 255}}; colors->SetColor("BkgColor", bkg.data()); // Create the renderer, the render window, and the interactor. The renderer // draws into the render window, the interactor enables mouse- and // keyboard-based interaction with the data within the render window. // vtkNew renWin; vtkNew iren; iren->SetRenderWindow(renWin); // std::vector data(file_list.size(),&Grid); // FieldMetaData header; int frameCount = file_list.size(); for(int d=0;d aCamera; aCamera->SetViewUp(0, 0, -1); aCamera->SetPosition(0, -1000, 0); aCamera->SetFocalPoint(0, 0, 0); aCamera->ComputeViewPlaneNormal(); aCamera->Azimuth(30.0); aCamera->Elevation(30.0); vtkNew aRenderer; renWin->AddRenderer(aRenderer); double vol = Grid.gSites(); std::cout << "Reading "< imageData; imageData->SetDimensions(latt_size[0],latt_size[1],latt_size[2]); imageData->AllocateScalars(VTK_DOUBLE, 1); for(int xx=0;xxSetScalarComponentFromDouble(xx,yy,zz,0,value); }}} vtkNew posExtractor; posExtractor->SetInputData(imageData); posExtractor->SetValue(0, contour); vtkNew posStripper; posStripper->SetInputConnection(posExtractor->GetOutputPort()); vtkNew posMapper; posMapper->SetInputConnection(posStripper->GetOutputPort()); posMapper->ScalarVisibilityOff(); vtkNew pos; pos->SetMapper(posMapper); pos->GetProperty()->SetDiffuseColor(colors->GetColor3d("posColor").GetData()); pos->GetProperty()->SetSpecular(0.3); pos->GetProperty()->SetSpecularPower(20); pos->GetProperty()->SetOpacity(0.5); // An isosurface, or contour value is set // The triangle stripper is used to create triangle strips from the // isosurface; these render much faster on may systems. vtkNew negExtractor; negExtractor->SetInputData(imageData); negExtractor->SetValue(0, -contour); vtkNew negStripper; negStripper->SetInputConnection(negExtractor->GetOutputPort()); vtkNew negMapper; negMapper->SetInputConnection(negStripper->GetOutputPort()); negMapper->ScalarVisibilityOff(); vtkNew neg; neg->SetMapper(negMapper); neg->GetProperty()->SetDiffuseColor(colors->GetColor3d("Ivory").GetData()); // An outline provides context around the data. vtkNew outlineData; outlineData->SetInputData(imageData); vtkNew mapOutline; mapOutline->SetInputConnection(outlineData->GetOutputPort()); vtkNew outline; outline->SetMapper(mapOutline); outline->GetProperty()->SetColor(colors->GetColor3d("Black").GetData()); vtkNew Text; // Text->SetInput(file_list[f].c_str()); Text->SetPosition2(0,0); Text->GetTextProperty()->SetFontSize(24); Text->GetTextProperty()->SetColor(colors->GetColor3d("Gold").GetData()); vtkNew TextT; TextT->SetInput("T=0"); TextT->SetPosition(0,.7*1025); TextT->GetTextProperty()->SetFontSize(24); TextT->GetTextProperty()->SetColor(colors->GetColor3d("Gold").GetData()); // Actors are added to the renderer. An initial camera view is created. // The Dolly() method moves the camera towards the FocalPoint, // thereby enlarging the image. // aRenderer->AddActor(Text); aRenderer->AddActor(TextT); aRenderer->AddActor(outline); aRenderer->AddActor(pos); aRenderer->AddActor(neg); // Sign up to receive TimerEvent vtkNew fu; fu->SetGrid(&Grid); fu->SetFiles(file_list); fu->SetSlice(Slice); fu->SetSumDimensions (SumDims); fu->SetLoopDimensions(LoopDims); fu->SetDisplayDimensions(XYZDims); fu->SetSliceDimensions(); fu->imageData = imageData; // fu->grid_data = &data[f]; fu->text = TextT; fu->maxCount = frameCount; fu->posExtractor = posExtractor; fu->negExtractor = negExtractor; fu->rms = rms; iren->AddObserver(vtkCommand::TimerEvent, fu); iren->AddObserver(vtkCommand::KeyPressEvent, fu); aRenderer->SetActiveCamera(aCamera); aRenderer->ResetCamera(); aRenderer->SetBackground(colors->GetColor3d("BkgColor").GetData()); aCamera->Dolly(1.0); // double nf = file_list.size(); // std::cout << " Adding renderer " <SetViewport(0.0, 0.0,1.0 , 1.0); // Note that when camera movement occurs (as it does in the Dolly() // method), the clipping planes often need adjusting. Clipping planes // consist of two planes: near and far along the view direction. The // near plane clips out objects in front of the plane; the far plane // clips out objects behind the plane. This way only what is drawn // between the planes is actually rendered. aRenderer->ResetCameraClippingRange(); // Set a background color for the renderer and set the size of the // render window (expressed in pixels). // Initialize the event loop and then start it. renWin->SetSize(1024, 1024); renWin->SetWindowName("FieldDensity"); renWin->Render(); // Take a pointer to the FrameUpdater for keypress mgt. // KBfu = fu; // vtkNew keypressCallback; // keypressCallback->SetCallback(KeypressCallbackFunction); // iren->AddObserver(vtkCommand::KeyPressEvent,keypressCallback); iren->Initialize(); if ( mpeg ) { #ifdef MPEG vtkWindowToImageFilter *imageFilter = vtkWindowToImageFilter::New(); imageFilter->SetInput( renWin ); imageFilter->SetInputBufferTypeToRGB(); vtkFFMPEGWriter *writer = vtkFFMPEGWriter::New(); writer->SetFileName("movie.avi"); writer->SetRate(framerate); writer->SetInputConnection(imageFilter->GetOutputPort()); writer->Start(); for(int i=0;imaxCount;i++){ fu->Execute(iren,vtkCommand::TimerEvent,nullptr); imageFilter->Modified(); writer->Write(); } writer->End(); writer->Delete(); #else assert(-1 && "MPEG support not compiled"); #endif } else { // Add control of contour threshold // Create a slider widget vtkSmartPointer sliderRep = vtkSmartPointer::New(); sliderRep->SetMinimumValue(0.0); sliderRep->SetMaximumValue(10.0); sliderRep->SetValue(1.0); sliderRep->SetTitleText("Fraction RMS"); // Set color properties: // Change the color of the knob that slides // sliderRep->GetSliderProperty()->SetColor(colors->GetColor3d("Green").GetData()); sliderRep->GetTitleProperty()->SetColor(colors->GetColor3d("AliceBlue").GetData()); sliderRep->GetLabelProperty()->SetColor(colors->GetColor3d("AliceBlue").GetData()); sliderRep->GetSelectedProperty()->SetColor(colors->GetColor3d("DeepPink").GetData()); // Change the color of the bar sliderRep->GetTubeProperty()->SetColor(colors->GetColor3d("MistyRose").GetData()); sliderRep->GetCapProperty()->SetColor(colors->GetColor3d("Yellow").GetData()); sliderRep->SetSliderLength(0.05); sliderRep->SetSliderWidth(0.025); sliderRep->SetEndCapLength(0.02); sliderRep->GetPoint1Coordinate()->SetCoordinateSystemToNormalizedDisplay(); sliderRep->GetPoint1Coordinate()->SetValue(0.1, 0.1); sliderRep->GetPoint2Coordinate()->SetCoordinateSystemToNormalizedDisplay(); sliderRep->GetPoint2Coordinate()->SetValue(0.9, 0.1); vtkSmartPointer sliderWidget = vtkSmartPointer::New(); sliderWidget->SetInteractor(iren); sliderWidget->SetRepresentation(sliderRep); sliderWidget->SetAnimationModeToAnimate(); sliderWidget->EnabledOn(); // Create the slider callback vtkSmartPointer slidercallback = vtkSmartPointer::New(); slidercallback->fu = fu; sliderWidget->AddObserver(vtkCommand::InteractionEvent, slidercallback); if ( NoTime==0 ) { int timerId = iren->CreateRepeatingTimer(10000/framerate); std::cout << "timerId "<Start(); } Grid_finalize(); return EXIT_SUCCESS; }