diff --git a/FDTD/processcurrent.cpp b/FDTD/processcurrent.cpp index eef3367..eec845c 100644 --- a/FDTD/processcurrent.cpp +++ b/FDTD/processcurrent.cpp @@ -26,18 +26,6 @@ ProcessCurrent::~ProcessCurrent() { } -void ProcessCurrent::OpenFile(string outfile) -{ - if (file.is_open()) file.close(); - - file.open(outfile.c_str()); - if (file.is_open()==false) - { - cerr << "Can't open file: " << outfile << endl; - return; - } -} - void ProcessCurrent::DefineStartStopCoord(double* dstart, double* dstop) { if (Op->SnapToMesh(dstart,start,true,m_start_inside)==false) cerr << "ProcessCurrent::DefineStartStopCoord: Warning: Snapped line outside field domain!!" << endl; @@ -96,3 +84,8 @@ int ProcessCurrent::Process() file << setprecision(m_precision) << (0.5 + (double)Eng->GetNumberOfTimesteps())*Op->GetTimestep() << "\t" << current << endl; return GetNextInterval(); } + +void ProcessCurrent::DumpBox2File( string vtkfilenameprefix, bool /*dualMesh*/ ) const +{ + Processing::DumpBox2File( vtkfilenameprefix, true ); +} diff --git a/FDTD/processcurrent.h b/FDTD/processcurrent.h index 9650378..247ec7d 100644 --- a/FDTD/processcurrent.h +++ b/FDTD/processcurrent.h @@ -26,15 +26,12 @@ public: ProcessCurrent(Operator* op, Engine* eng); virtual ~ProcessCurrent(); - virtual void OpenFile(string outfile); - virtual void DefineStartStopCoord(double* dstart, double* dstop); virtual int Process(); + virtual void DumpBox2File( string vtkfilenameprefix, bool dualMesh = false ) const; //!< dump geometry to file protected: - ofstream file; - vector v_current; }; diff --git a/FDTD/processfields.cpp b/FDTD/processfields.cpp index 5bb3382..75a8257 100644 --- a/FDTD/processfields.cpp +++ b/FDTD/processfields.cpp @@ -349,5 +349,3 @@ bool ProcessFields::DumpVectorArray2HDF5(string filename, string name, FDTD_FLOA dataset.write( hdf5array, H5::PredType::NATIVE_FLOAT ); return true; } - - diff --git a/FDTD/processfields.h b/FDTD/processfields.h index 6f71635..71611f4 100644 --- a/FDTD/processfields.h +++ b/FDTD/processfields.h @@ -39,10 +39,10 @@ public: virtual void SetSubSampling(unsigned int subSampleRate, int dir=-1); //! Used file pattern e.g. pattern="tmp/efield_" --> "tmp/efield_000045.vtk" for timestep 45 or "tmp/efield_2.40000e9.vtk" for 2.4GHz E-field dump. (VTK FileType only) \sa SetFileType - void SetFilePattern(string fp) {filePattern=fp;} + void SetFilePattern(string fp) {m_filename=filePattern=fp;} //! Set the filename for a hdf5 data group file (HDF5 FileType only) \sa SetFileType - void SetFileName(string fn) {m_fileName=fn;} + void SetFileName(string fn) {m_filename=m_fileName=fn;} //! Define the Dump-Mode void SetDumpMode(DumpMode mode) {m_DumpMode=mode;} diff --git a/FDTD/processing.cpp b/FDTD/processing.cpp index c729859..8467da9 100644 --- a/FDTD/processing.cpp +++ b/FDTD/processing.cpp @@ -29,6 +29,7 @@ Processing::Processing(Operator* op, Engine* eng) Processing::~Processing() { + file.close(); } void Processing::Reset() @@ -122,6 +123,71 @@ double Processing::CalcLineIntegral(unsigned int* start, unsigned int* stop, int return result; } +void Processing::OpenFile( string outfile ) +{ + if (file.is_open()) + file.close(); + + file.open( outfile.c_str() ); + if (!file.is_open()) + cerr << "Can't open file: " << outfile << endl; + + m_filename = outfile; +} + +void Processing::DumpBox2File( string vtkfilenameprefix, bool dualMesh ) const +{ + string vtkfilename = vtkfilenameprefix + m_filename + ".vtk"; + + ofstream file( vtkfilename.c_str() ); + if (!file.is_open()) + { + cerr << "Processing::DumpBoxes2File(): Can't open file: " << vtkfilename << endl; + return; + } + + // normalize coordinates + double s1[3], s2[3]; + for (int i=0; i<3; i++) { + s1[i] = min(Op->GetDiscLine(i,start[i],dualMesh),Op->GetDiscLine(i,stop[i],dualMesh)); + s2[i] = max(Op->GetDiscLine(i,start[i],dualMesh),Op->GetDiscLine(i,stop[i],dualMesh)); + } + + // fix degenerate box/plane -> line (paraview display problem) + if (((s1[0] == s2[0]) && (s1[1] == s2[1])) || ((s1[0] == s2[0]) && (s1[2] == s2[2])) || ((s1[2] == s2[2]) && (s1[1] == s2[1]))) { + // line are not displayed correctly -> enlarge + for (int i=0; i<3; i++) { + double delta = min( Op->GetMeshDelta( i, start,dualMesh ), Op->GetMeshDelta( i, stop,dualMesh ) ) / Op->GetGridDelta() / 4.0; + s1[i] -= delta; + s2[i] += delta; + } + } + + file << "# vtk DataFile Version 2.0" << endl; + file << "" << endl; + file << "ASCII" << endl; + file << "DATASET POLYDATA" << endl; + + file << "POINTS 8 float" << endl; + file << s1[0] << " " << s1[1] << " " << s1[2] << endl; + file << s2[0] << " " << s1[1] << " " << s1[2] << endl; + file << s2[0] << " " << s2[1] << " " << s1[2] << endl; + file << s1[0] << " " << s2[1] << " " << s1[2] << endl; + file << s1[0] << " " << s1[1] << " " << s2[2] << endl; + file << s2[0] << " " << s1[1] << " " << s2[2] << endl; + file << s2[0] << " " << s2[1] << " " << s2[2] << endl; + file << s1[0] << " " << s2[1] << " " << s2[2] << endl; + + file << "POLYGONS 6 30" << endl; + file << "4 0 1 2 3" << endl; + file << "4 4 5 6 7" << endl; + file << "4 7 6 2 3" << endl; + file << "4 4 5 1 0" << endl; + file << "4 0 4 7 3" << endl; + file << "4 5 6 2 1" << endl; + + file.close(); +} void ProcessingArray::AddProcessing(Processing* proc) { @@ -157,3 +223,9 @@ int ProcessingArray::Process() } return nextProcess; } + +void ProcessingArray::DumpBoxes2File( string vtkfilenameprefix ) const +{ + for (size_t i=0;iDumpBox2File( vtkfilenameprefix ); +} diff --git a/FDTD/processing.h b/FDTD/processing.h index d42d912..5224b51 100644 --- a/FDTD/processing.h +++ b/FDTD/processing.h @@ -49,6 +49,9 @@ public: //! Set the dump precision void SetPrecision(unsigned int val) {m_precision = val;} + virtual void OpenFile( string outfile ); + + virtual void DumpBox2File( string vtkfilenameprefix, bool dualMesh = false ) const; //!< dump geometry to file protected: Operator* Op; Engine* Eng; @@ -70,6 +73,9 @@ protected: bool m_start_inside[3]; bool m_stop_inside[3]; + ofstream file; + string m_filename; + double CalcLineIntegral(unsigned int* start, unsigned int* stop, int field); }; @@ -88,6 +94,8 @@ public: int Process(); + void DumpBoxes2File( string vtkfilenameprefix ) const; + protected: unsigned int maxInterval; vector ProcessArray; diff --git a/FDTD/processvoltage.cpp b/FDTD/processvoltage.cpp index c3eba24..abb8fe8 100644 --- a/FDTD/processvoltage.cpp +++ b/FDTD/processvoltage.cpp @@ -24,19 +24,6 @@ ProcessVoltage::ProcessVoltage(Operator* op, Engine* eng) : Processing(op, eng) ProcessVoltage::~ProcessVoltage() { - file.close(); -} - -void ProcessVoltage::OpenFile(string outfile) -{ - if (file.is_open()) file.close(); - - file.open(outfile.c_str()); - if (file.is_open()==false) - { - cerr << "Can't open file: " << outfile << endl; - return; - } } int ProcessVoltage::Process() diff --git a/FDTD/processvoltage.h b/FDTD/processvoltage.h index 4c3a0df..7272343 100644 --- a/FDTD/processvoltage.h +++ b/FDTD/processvoltage.h @@ -27,13 +27,9 @@ public: ProcessVoltage(Operator* op, Engine* eng); virtual ~ProcessVoltage(); - virtual void OpenFile(string outfile); - virtual int Process(); protected: - ofstream file; - vector voltages; }; diff --git a/openems.cpp b/openems.cpp index 016e3c2..a918e6d 100644 --- a/openems.cpp +++ b/openems.cpp @@ -47,6 +47,7 @@ openEMS::openEMS() Enable_Dumps = true; DebugMat = false; DebugOp = false; + m_debugBox = false; endCrit = 1e-6; m_OverSampling = 4; @@ -93,6 +94,12 @@ bool openEMS::parseCommandLineArgument( const char *argv ) DebugOperator(); return true; } + else if (strcmp(argv,"--debug-boxes")==0) + { + cout << "openEMS - dumping boxes to 'box_dump*.vtk'" << endl; + DebugBox(); + return true; + } else if (strcmp(argv,"--engine=multithreaded")==0) { cout << "openEMS - enabled multithreading" << endl; @@ -364,6 +371,13 @@ int openEMS::SetupFDTD(const char* file) delete ProcTD; } } + + // dump all boxes (voltage, current, fields, ...) + if (m_debugBox) + { + PA->DumpBoxes2File("box_dump_"); + } + return 0; } diff --git a/openems.h b/openems.h index 3b144e2..475a5cd 100644 --- a/openems.h +++ b/openems.h @@ -42,6 +42,7 @@ public: void DebugMaterial() {DebugMat=true;} void DebugOperator() {DebugOp=true;} + void DebugBox() {m_debugBox=true;} protected: void SetupExcitation(TiXmlElement* Excite); @@ -53,6 +54,7 @@ protected: bool Enable_Dumps; bool DebugMat; bool DebugOp; + bool m_debugBox; double endCrit; int m_OverSampling; Operator* FDTD_Op;