![]() |
MechSys
1.0
Computing library for simulations in continuum and discrete mechanics
|
00001 /************************************************************************ 00002 * MechSys - Open Library for Mechanical Systems * 00003 * Copyright (C) 2005 Dorival M. Pedroso, Raúl D. D. Farfan * 00004 * * 00005 * This program is free software: you can redistribute it and/or modify * 00006 * it under the terms of the GNU General Public License as published by * 00007 * the Free Software Foundation, either version 3 of the License, or * 00008 * any later version. * 00009 * * 00010 * This program is distributed in the hope that it will be useful, * 00011 * but WITHOUT ANY WARRANTY; without even the implied warranty of * 00012 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the * 00013 * GNU General Public License for more details. * 00014 * * 00015 * You should have received a copy of the GNU General Public License * 00016 * along with this program. If not, see <http://www.gnu.org/licenses/> * 00017 ************************************************************************/ 00018 00019 #ifndef MECHSYS_SGRID_H 00020 #define MECHSYS_SGRID_H 00021 00022 // VTK 00023 #include <vtkPoints.h> 00024 #include <vtkDoubleArray.h> 00025 #include <vtkStructuredGrid.h> 00026 #include <vtkDataSetMapper.h> 00027 #include <vtkActor.h> 00028 #include <vtkTextActor3D.h> 00029 #include <vtkProperty.h> 00030 #include <vtkTextProperty.h> 00031 #include <vtkPointData.h> 00032 #include <vtkStructuredGridWriter.h> 00033 #include <vtkColorTransferFunction.h> 00034 00035 // MechSys 00036 #include <mechsys/vtk/win.h> 00037 #include <mechsys/util/array.h> 00038 #include <mechsys/util/colors.h> 00039 #include <mechsys/util/string.h> 00040 #include <mechsys/linalg/matvec.h> 00041 00042 namespace VTK 00043 { 00044 00045 typedef void (*GridCallBack) (Vec3_t const & X, double & F, Vec3_t & V, void * UserData); 00046 00047 class SGrid 00048 { 00049 public: 00050 // Constructor & Destructor 00051 SGrid () : _func(NULL), _udat(NULL) { _create(); } 00052 ~SGrid (); 00053 00054 // Alternative constructors 00055 SGrid (int N[3], double L[6], GridCallBack Func=NULL, void * UserData=NULL); 00056 SGrid (Array<int> const & N, Array<double> const & L, GridCallBack Func=NULL, void * UserData=NULL); 00057 SGrid (int Nx, double Xmin, double Xmax, 00058 int Ny, double Ymin, double Ymax, 00059 int Nz, double Zmin, double Zmax, GridCallBack Func=NULL, void * UserData=NULL); 00060 00061 // Set methods 00062 SGrid & Resize (int N[3], double L[6]) { return Resize(N[0],L[0],L[1], N[1],L[2],L[3], N[2],L[4],L[5]); } 00063 SGrid & Resize (Array<int> const & N, Array<double> const & L) { return Resize(N[0],L[0],L[1], N[1],L[2],L[3], N[2],L[4],L[5]); } 00064 SGrid & Resize (int Nx, double Xmin, double Xmax, 00065 int Ny, double Ymin, double Ymax, 00066 int Nz, double Zmin, double Zmax); 00067 SGrid & SetColor (char const * Name="black", double Opacity=1.0); 00068 SGrid & SetFunc (GridCallBack Func, void * UserData=NULL) { _func=Func; _udat=UserData; _calc_f(); return (*this); } 00069 00070 // Additional methods 00071 double GetF (int i, int j, int k) const { return _scalars->GetTuple1 (i+j*_Nx+k*_Nx*_Ny); } 00072 void SetF (int i, int j, int k, double F) { _scalars->SetTuple1 (i+j*_Nx+k*_Nx*_Ny, F); } 00073 void SetCMap (double Fmin, double Fmax, char const * Name="Diverging"); 00074 void SetCMap ( char const * Name="Diverging") { SetCMap (_Fmin, _Fmax, Name); } 00075 void RescaleCMap (); 00076 00077 // Access methods 00078 int Size () const { return _points->GetNumberOfPoints(); } 00079 void GetPoint (int i, Vec3_t & x) const { _points->GetPoint (i, x.data()); } 00080 void SetPoint (int i, Vec3_t & x) { _points->SetPoint (i, x.data()); } 00081 vtkStructuredGrid * GetGrid () { return _sgrid; } 00082 vtkPoints * GetPoints () { return _points; } 00083 vtkDoubleArray * GetScalars () { return _scalars; } 00084 vtkDoubleArray * GetVectors () { return _vectors; } 00085 00086 // Methods 00087 void ShowWire () { _sgrid_actor->GetProperty()->SetRepresentationToWireframe(); } 00088 void ShowSurface () { _sgrid_actor->GetProperty()->SetRepresentationToSurface(); } 00089 void ShowPoints (int PtSize=4) { _sgrid_actor->GetProperty()->SetRepresentationToPoints(); _sgrid_actor->GetProperty()->SetPointSize(PtSize); } 00090 void ShowIds (double OriX=90, double OriY=90, double OriZ=45, double Scale=0.003, int SizePt=14, bool Shadow=true, char const * Color="blue"); 00091 void AddTo (VTK::Win & win); 00092 void WriteVTK (char const * Filekey); 00093 void FilterV (double F=0.0, double Tol=1.0e-3, bool Normalize=false); 00094 00095 private: 00096 GridCallBack _func; 00097 void * _udat; 00098 vtkPoints * _points; 00099 vtkDoubleArray * _scalars; 00100 vtkDoubleArray * _vectors; 00101 vtkStructuredGrid * _sgrid; 00102 vtkDataSetMapper * _sgrid_mapper; 00103 vtkActor * _sgrid_actor; 00104 vtkColorTransferFunction * _color_func; 00105 Array<vtkTextActor3D*> _text; 00106 double _Fmin; 00107 double _Fmax; 00108 int _Nx, _Ny, _Nz; 00109 double _Xmin,_Xmax, _Ymin,_Ymax, _Zmin,_Zmax; 00110 String _cmap_name; 00111 void _create (); 00112 void _calc_f (); 00113 }; 00114 00115 00117 00118 00119 inline SGrid::SGrid (int N[3], double L[6], GridCallBack Func, void * UserData) 00120 : _func(Func), _udat(UserData) 00121 { 00122 _create (); 00123 Resize (N,L); // also calculates F and set _Fmin and _Fmax 00124 SetCMap (_Fmin, _Fmax); 00125 } 00126 00127 inline SGrid::SGrid (Array<int> const & N, Array<double> const & L, GridCallBack Func, void * UserData) 00128 : _func(Func), _udat(UserData) 00129 { 00130 _create (); 00131 Resize (N,L); // also calculates F and set _Fmin and _Fmax 00132 SetCMap (_Fmin, _Fmax); 00133 } 00134 00135 inline SGrid::SGrid (int Nx, double Xmin, double Xmax, int Ny, double Ymin, double Ymax, int Nz, double Zmin, double Zmax, GridCallBack Func, void * UserData) 00136 : _func(Func), _udat(UserData) 00137 { 00138 _create (); 00139 Resize (Nx,Xmin,Xmax, Ny,Ymin,Ymax, Nz,Zmin,Zmax); // also calculates F and set _Fmin and _Fmax 00140 SetCMap (_Fmin, _Fmax); 00141 } 00142 00143 inline SGrid & SGrid::Resize (int Nx, double Xmin, double Xmax, int Ny, double Ymin, double Ymax, int Nz, double Zmin, double Zmax) 00144 { 00145 if (Nx<2) throw new Fatal("SGrid::Resize: Nx==N[0]=%d must be greater than 1",Nx); 00146 if (Ny<2) throw new Fatal("SGrid::Resize: Ny==N[1]=%d must be greater than 1",Ny); 00147 if (Nz<1) throw new Fatal("SGrid::Resize: Nz==N[2]=%d must be greater than 1",Nz); 00148 _Nx = Nx; _Xmin = Xmin; _Xmax = Xmax; 00149 _Ny = Ny; _Ymin = Ymin; _Ymax = Ymax; 00150 _Nz = Nz; _Zmin = Zmin; _Zmax = Zmax; 00151 _points -> Reset (); 00152 _points -> Allocate (Nx*Ny*Nz); 00153 _scalars -> Reset (); 00154 _scalars -> Allocate (Nx*Ny*Nz); 00155 _vectors -> Reset (); 00156 _vectors -> SetNumberOfComponents (3); 00157 _vectors -> SetNumberOfTuples (Nx*Ny*Nz); 00158 _sgrid -> SetDimensions (Nx,Ny,Nz); 00159 double dx = (Xmax-Xmin)/(Nx-1.0); 00160 double dy = (Ymax-Ymin)/(Ny-1.0); 00161 double dz = (Nz>1 ? (Zmax-Zmin)/(Nz-1.0) : 0.0); 00162 double f = 0.0; 00163 Vec3_t x, v; 00164 if (_func==NULL) 00165 { 00166 _Fmin = 0.0; 00167 _Fmax = 1.0; 00168 } 00169 else 00170 { 00171 x = Xmin, Ymin, Zmin; 00172 (*_func) (x, _Fmin, v, _udat); 00173 (*_func) (x, _Fmax, v, _udat); 00174 } 00175 for (int k=0; k<Nz; ++k) 00176 for (int j=0; j<Ny; ++j) 00177 for (int i=0; i<Nx; ++i) 00178 { 00179 int idx = i + j*_Nx + k*_Nx*_Ny; 00180 x = Xmin+i*dx, Ymin+j*dy, Zmin+k*dz; 00181 _points -> InsertPoint (idx, x.data()); 00182 if (_func==NULL) 00183 { 00184 _scalars -> InsertTuple1 (idx, 0); 00185 _vectors -> InsertTuple3 (idx, 0,0,0); 00186 } 00187 else 00188 { 00189 (*_func) (x, f, v, _udat); 00190 _scalars -> InsertTuple1 (idx, f); 00191 _vectors -> InsertTuple3 (idx, v(0), v(1), v(2)); 00192 if (f<_Fmin) _Fmin = f; 00193 if (f>_Fmax) _Fmax = f; 00194 } 00195 } 00196 _sgrid -> GetPointData() -> SetScalars (_scalars); 00197 _sgrid -> GetPointData() -> SetVectors (_vectors); 00198 return (*this); 00199 } 00200 00201 inline SGrid::~SGrid () 00202 { 00203 _points -> Delete(); 00204 _scalars -> Delete(); 00205 _vectors -> Delete(); 00206 _sgrid -> Delete(); 00207 _sgrid_mapper -> Delete(); 00208 _sgrid_actor -> Delete(); 00209 _color_func -> Delete(); 00210 for (size_t i=0; i<_text.Size(); ++i) _text[i] -> Delete(); 00211 } 00212 00213 inline SGrid & SGrid::SetColor (char const * Name, double Opacity) 00214 { 00215 Vec3_t c = Colors::Get(Name); 00216 _sgrid_actor->GetProperty()->SetColor (c.data()); 00217 _sgrid_actor->GetProperty()->SetOpacity (Opacity); 00218 return (*this); 00219 } 00220 00221 inline void SGrid::SetCMap (double Fmin, double Fmax, char const * Name) 00222 { 00223 _cmap_name = Name; 00224 double midpoint = 0.5; // halfway between the control points 00225 double sharpness = 0.0; // linear 00226 if (_color_func->GetSize()==2) _color_func->RemoveAllPoints(); // existent 00227 if (_cmap_name=="Rainbow") 00228 { 00229 _color_func -> SetColorSpaceToHSV (); 00230 _color_func -> HSVWrapOff (); 00231 _color_func -> AddHSVPoint (Fmin, 2.0/3.0, 1.0, 1.0, midpoint, sharpness); 00232 _color_func -> AddHSVPoint (Fmax, 0.0, 1.0, 1.0, midpoint, sharpness); 00233 } 00234 else 00235 { 00236 _color_func -> SetColorSpaceToDiverging (); 00237 _color_func -> HSVWrapOn (); 00238 _color_func -> AddRGBPoint (Fmin, 0.230, 0.299, 0.754, midpoint, sharpness); 00239 _color_func -> AddRGBPoint (Fmax, 0.706, 0.016, 0.150, midpoint, sharpness); 00240 } 00241 } 00242 00243 inline void SGrid::RescaleCMap () 00244 { 00245 _Fmin = _scalars->GetTuple1 (0); 00246 _Fmax = _Fmin; 00247 for (int i=0; i<_scalars->GetNumberOfTuples(); ++i) 00248 { 00249 double f = _scalars->GetTuple1(i); 00250 if (f<_Fmin) _Fmin = f; 00251 if (f>_Fmax) _Fmax = f; 00252 } 00253 SetCMap (_Fmin, _Fmax, _cmap_name.CStr()); 00254 } 00255 00256 inline void SGrid::ShowIds (double OriX, double OriY, double OriZ, double Scale, int SizePt, bool Shadow, char const * Color) 00257 { 00258 Vec3_t c(Colors::Get(Color)); 00259 for (size_t i=0; i<_text.Size(); ++i) _text[i] -> Delete(); 00260 String buf; 00261 _text.Resize (_points->GetNumberOfPoints()); 00262 for (int i=0; i<_points->GetNumberOfPoints(); ++i) 00263 { 00264 buf.Printf ("%d",i); 00265 _text[i] = vtkTextActor3D ::New(); 00266 _text[i] -> SetInput (buf.CStr()); 00267 _text[i] -> SetPosition (_points->GetPoint(i)); 00268 _text[i] -> SetOrientation (OriX, OriY, OriZ); 00269 _text[i] -> SetScale (Scale); 00270 _text[i] -> GetTextProperty()-> SetFontSize (SizePt); 00271 _text[i] -> GetTextProperty()-> SetShadow (Shadow); 00272 _text[i] -> GetTextProperty()-> SetColor (c.data()); 00273 } 00274 } 00275 00276 inline void SGrid::AddTo (VTK::Win & win) 00277 { 00278 win.AddActor (_sgrid_actor); 00279 for (size_t i=0; i<_text.Size(); ++i) win.AddActor (reinterpret_cast<vtkActor*>(_text[i])); 00280 } 00281 00282 inline void SGrid::WriteVTK (char const * Filekey) 00283 { 00284 String buf(Filekey); 00285 buf.append(".vtk"); 00286 vtkStructuredGridWriter * writer = vtkStructuredGridWriter::New(); 00287 writer -> SetInput (_sgrid); 00288 writer -> SetFileName (buf.CStr()); 00289 writer -> Write (); 00290 writer -> Delete (); 00291 printf("File <%s%s.vtk%s> written\n", TERM_CLR_BLUE_H, Filekey, TERM_RST); 00292 } 00293 00294 inline void SGrid::FilterV (double F, double Tol, bool Normalize) 00295 { 00296 for (int i=0; i<_scalars->GetNumberOfTuples(); ++i) 00297 { 00298 double f = _scalars->GetTuple1(i); 00299 if (Normalize) 00300 { 00301 double * v = _vectors->GetTuple3(i); 00302 double norm = sqrt(v[0]*v[0] + v[1]*v[1] + v[2]*v[2]); 00303 if (norm>0.0) _vectors->SetTuple3 (i, v[0]/norm, v[1]/norm, v[2]/norm); 00304 } 00305 if (fabs(f-F)>Tol) _vectors->SetTuple3 (i, 0.0, 0.0, 0.0); 00306 } 00307 } 00308 00309 inline void SGrid::_create () 00310 { 00311 _points = vtkPoints ::New(); 00312 _scalars = vtkDoubleArray ::New(); 00313 _vectors = vtkDoubleArray ::New(); 00314 _sgrid = vtkStructuredGrid ::New(); 00315 _sgrid_mapper = vtkDataSetMapper ::New(); 00316 _sgrid_actor = vtkActor ::New(); 00317 _color_func = vtkColorTransferFunction ::New(); 00318 _sgrid -> SetPoints (_points); 00319 _sgrid_mapper -> SetInput (_sgrid); 00320 _sgrid_mapper -> SetLookupTable (_color_func); 00321 _sgrid_actor -> SetMapper (_sgrid_mapper); 00322 _sgrid_actor -> GetProperty() -> SetPointSize (4); 00323 ShowWire (); 00324 SetColor (); 00325 } 00326 00327 inline void SGrid::_calc_f () 00328 { 00329 double f; 00330 Vec3_t x, v; 00331 if (_func==NULL) 00332 { 00333 _Fmin = 0.0; 00334 _Fmax = 1.0; 00335 } 00336 else 00337 { 00338 GetPoint (0, x); 00339 (*_func) (x, _Fmin, v, _udat); 00340 (*_func) (x, _Fmax, v, _udat); 00341 } 00342 for (int i=0; i<Size(); ++i) 00343 { 00344 GetPoint (i, x); 00345 if (_func==NULL) 00346 { 00347 _scalars -> SetTuple1 (i, 0); 00348 _vectors -> SetTuple3 (i, 0,0,0); 00349 } 00350 else 00351 { 00352 (*_func) (x, f, v, _udat); 00353 _scalars -> SetTuple1 (i, f); 00354 _vectors -> SetTuple3 (i, v(0), v(1), v(2)); 00355 if (f<_Fmin) _Fmin = f; 00356 if (f>_Fmax) _Fmax = f; 00357 } 00358 } 00359 } 00360 00361 }; // namespace VTK 00362 00363 #endif // MECHSYS_SGRID_H