2010-03-11 15:47:40 +00:00
|
|
|
/*
|
2015-03-06 20:38:27 +00:00
|
|
|
* Copyright (C) 2010-2015 Thorsten Liebig (Thorsten.Liebig@gmx.de)
|
2010-03-11 15:47:40 +00:00
|
|
|
*
|
|
|
|
* This program is free software: you can redistribute it and/or modify
|
|
|
|
* it under the terms of the GNU General Public License as published by
|
|
|
|
* the Free Software Foundation, either version 3 of the License, or
|
|
|
|
* (at your option) any later version.
|
|
|
|
*
|
|
|
|
* This program is distributed in the hope that it will be useful,
|
|
|
|
* but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
* GNU General Public License for more details.
|
|
|
|
*
|
|
|
|
* You should have received a copy of the GNU General Public License
|
|
|
|
* along with this program. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
*/
|
|
|
|
|
2010-03-11 09:56:19 +00:00
|
|
|
#include "openems.h"
|
2010-03-15 19:50:49 +00:00
|
|
|
#include <iomanip>
|
2012-09-24 13:16:47 +00:00
|
|
|
#include <iostream>
|
|
|
|
#include <fstream>
|
2010-03-11 09:56:19 +00:00
|
|
|
#include "tools/array_ops.h"
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
#include "tools/signal.h"
|
2011-11-08 10:49:14 +00:00
|
|
|
#include "tools/useful.h"
|
2010-05-06 20:55:59 +00:00
|
|
|
#include "FDTD/operator_cylinder.h"
|
2010-09-08 05:36:32 +00:00
|
|
|
#include "FDTD/operator_cylindermultigrid.h"
|
2010-03-26 11:57:52 +00:00
|
|
|
#include "FDTD/engine_multithread.h"
|
2010-05-20 20:02:06 +00:00
|
|
|
#include "FDTD/operator_multithread.h"
|
2012-07-17 11:23:00 +00:00
|
|
|
#include "FDTD/extensions/operator_ext_excitation.h"
|
2012-07-18 11:12:25 +00:00
|
|
|
#include "FDTD/extensions/operator_ext_tfsf.h"
|
2010-12-06 14:30:47 +00:00
|
|
|
#include "FDTD/extensions/operator_ext_mur_abc.h"
|
|
|
|
#include "FDTD/extensions/operator_ext_upml.h"
|
|
|
|
#include "FDTD/extensions/operator_ext_lorentzmaterial.h"
|
2023-11-18 11:23:15 +00:00
|
|
|
#include "FDTD/extensions/operator_ext_lumpedRLC.h"
|
2012-05-08 11:58:20 +00:00
|
|
|
#include "FDTD/extensions/operator_ext_conductingsheet.h"
|
2015-05-04 18:47:19 +00:00
|
|
|
#include "FDTD/extensions/operator_ext_steadystate.h"
|
|
|
|
#include "FDTD/extensions/engine_ext_steadystate.h"
|
2010-12-02 12:51:34 +00:00
|
|
|
#include "FDTD/engine_interface_fdtd.h"
|
2011-04-13 10:16:54 +00:00
|
|
|
#include "FDTD/engine_interface_cylindrical_fdtd.h"
|
2010-12-06 09:44:25 +00:00
|
|
|
#include "Common/processvoltage.h"
|
|
|
|
#include "Common/processcurrent.h"
|
2011-01-18 09:45:03 +00:00
|
|
|
#include "Common/processfieldprobe.h"
|
2010-12-06 09:44:25 +00:00
|
|
|
#include "Common/processmodematch.h"
|
|
|
|
#include "Common/processfields_td.h"
|
2010-12-19 19:41:08 +00:00
|
|
|
#include "Common/processfields_fd.h"
|
2011-01-31 11:25:55 +00:00
|
|
|
#include "Common/processfields_sar.h"
|
2012-11-11 13:54:51 +00:00
|
|
|
#include <hdf5.h> // only for H5get_libversion()
|
2010-10-04 09:35:20 +00:00
|
|
|
#include <boost/version.hpp> // only for BOOST_LIB_VERSION
|
2011-04-08 10:27:42 +00:00
|
|
|
#include <vtkVersion.h>
|
2010-03-11 09:56:19 +00:00
|
|
|
|
|
|
|
//external libs
|
|
|
|
#include "tinyxml.h"
|
|
|
|
#include "ContinuousStructure.h"
|
2012-12-03 12:59:39 +00:00
|
|
|
#include "CSPropProbeBox.h"
|
2021-08-25 17:05:38 +00:00
|
|
|
#include "CSPrimBox.h"
|
2012-12-03 12:59:39 +00:00
|
|
|
#include "CSPropDumpBox.h"
|
2010-03-11 09:56:19 +00:00
|
|
|
|
2015-06-18 19:45:22 +00:00
|
|
|
using namespace std;
|
|
|
|
|
2010-03-25 14:08:54 +00:00
|
|
|
double CalcDiffTime(timeval t1, timeval t2)
|
|
|
|
{
|
2010-03-27 14:54:44 +00:00
|
|
|
double s_diff = t1.tv_sec - t2.tv_sec;
|
|
|
|
s_diff += (t1.tv_usec-t2.tv_usec)*1e-6;
|
2010-03-25 14:08:54 +00:00
|
|
|
return s_diff;
|
|
|
|
}
|
|
|
|
|
2010-03-11 09:56:19 +00:00
|
|
|
openEMS::openEMS()
|
|
|
|
{
|
2022-12-29 12:08:44 +00:00
|
|
|
setlocale(LC_NUMERIC, "en_US.UTF-8");
|
2010-03-11 09:56:19 +00:00
|
|
|
FDTD_Op=NULL;
|
|
|
|
FDTD_Eng=NULL;
|
2015-05-04 18:47:19 +00:00
|
|
|
Eng_Ext_SSD=NULL;
|
2011-01-31 11:00:00 +00:00
|
|
|
m_CSX=NULL;
|
2010-03-11 09:56:19 +00:00
|
|
|
PA=NULL;
|
2010-04-09 13:51:37 +00:00
|
|
|
CylinderCoords = false;
|
2010-03-11 14:48:55 +00:00
|
|
|
Enable_Dumps = true;
|
2010-03-12 19:39:04 +00:00
|
|
|
DebugMat = false;
|
2010-03-17 22:16:41 +00:00
|
|
|
DebugOp = false;
|
2010-07-08 09:28:11 +00:00
|
|
|
m_debugCSX = false;
|
2010-08-16 21:17:19 +00:00
|
|
|
m_debugBox = m_debugPEC = m_no_simulation = false;
|
2012-09-24 13:16:47 +00:00
|
|
|
m_DumpStats = false;
|
2010-03-15 15:59:37 +00:00
|
|
|
endCrit = 1e-6;
|
2010-04-04 17:48:36 +00:00
|
|
|
m_OverSampling = 4;
|
2013-08-23 15:29:10 +00:00
|
|
|
m_CellConstantMaterial=false;
|
2010-03-30 11:13:00 +00:00
|
|
|
|
2010-12-06 14:27:48 +00:00
|
|
|
m_engine = EngineType_Multithreaded; //default engine type
|
2010-03-30 11:13:00 +00:00
|
|
|
m_engine_numThreads = 0;
|
2010-10-18 11:26:25 +00:00
|
|
|
|
|
|
|
m_Abort = false;
|
2012-07-17 11:23:00 +00:00
|
|
|
m_Exc = 0;
|
2015-12-13 21:33:26 +00:00
|
|
|
|
|
|
|
m_TS_method=3;
|
|
|
|
m_TS=0;
|
|
|
|
m_TS_fac=1.0;
|
|
|
|
m_maxTime=0.0;
|
|
|
|
|
|
|
|
for (int n=0;n<6;++n)
|
|
|
|
{
|
|
|
|
m_BC_type[n] = 0;
|
|
|
|
m_PML_size[n] = 8;
|
|
|
|
m_Mur_v_ph[n] = 0;
|
|
|
|
}
|
2010-03-11 09:56:19 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
openEMS::~openEMS()
|
|
|
|
{
|
2010-04-01 14:11:25 +00:00
|
|
|
Reset();
|
2010-03-11 09:56:19 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
void openEMS::Reset()
|
|
|
|
{
|
|
|
|
if (PA) PA->DeleteAll();
|
2010-12-06 12:04:37 +00:00
|
|
|
delete PA;
|
|
|
|
PA=0;
|
|
|
|
delete FDTD_Eng;
|
|
|
|
FDTD_Eng=0;
|
|
|
|
delete FDTD_Op;
|
|
|
|
FDTD_Op=0;
|
2011-01-31 11:00:00 +00:00
|
|
|
delete m_CSX;
|
|
|
|
m_CSX=0;
|
2012-07-17 11:23:00 +00:00
|
|
|
delete m_Exc;
|
|
|
|
m_Exc=0;
|
2010-03-11 09:56:19 +00:00
|
|
|
}
|
|
|
|
|
2015-12-19 14:38:21 +00:00
|
|
|
void openEMS::showUsage()
|
|
|
|
{
|
|
|
|
cout << " Usage: openEMS <FDTD_XML_FILE> [<options>...]" << endl << endl;
|
|
|
|
cout << " <options>" << endl;
|
|
|
|
cout << "\t--disable-dumps\t\tDisable all field dumps for faster simulation" << endl;
|
|
|
|
cout << "\t--debug-material\tDump material distribution to a vtk file for debugging" << endl;
|
|
|
|
cout << "\t--debug-PEC\t\tDump metal distribution to a vtk file for debugging" << endl;
|
|
|
|
cout << "\t--debug-operator\tDump operator to vtk file for debugging" << endl;
|
|
|
|
cout << "\t--debug-boxes\t\tDump e.g. probe boxes to vtk file for debugging" << endl;
|
|
|
|
cout << "\t--debug-CSX\t\tWrite CSX geometry file to debugCSX.xml" << endl;
|
|
|
|
cout << "\t--engine=<type>\t\tChoose engine type" << endl;
|
|
|
|
cout << "\t\t--engine=fastest\t\tfastest available engine (default)" << endl;
|
|
|
|
cout << "\t\t--engine=basic\t\t\tbasic FDTD engine" << endl;
|
|
|
|
cout << "\t\t--engine=sse\t\t\tengine using sse vector extensions" << endl;
|
|
|
|
cout << "\t\t--engine=sse-compressed\t\tengine using compressed operator + sse vector extensions" << endl;
|
|
|
|
#ifdef MPI_SUPPORT
|
|
|
|
cout << "\t\t--engine=MPI\t\t\tengine using compressed operator + sse vector extensions + MPI parallel processing" << endl;
|
|
|
|
cout << "\t\t--engine=multithreaded\t\tengine using compressed operator + sse vector extensions + MPI + multithreading" << endl;
|
|
|
|
#else
|
|
|
|
cout << "\t\t--engine=multithreaded\t\tengine using compressed operator + sse vector extensions + multithreading" << endl;
|
|
|
|
#endif
|
|
|
|
cout << "\t--numThreads=<n>\tForce use n threads for multithreaded engine (needs: --engine=multithreaded)" << endl;
|
|
|
|
cout << "\t--no-simulation\t\tonly run preprocessing; do not simulate" << endl;
|
|
|
|
cout << "\t--dump-statistics\tdump simulation statistics to '" << __OPENEMS_RUN_STAT_FILE__ << "' and '" << __OPENEMS_STAT_FILE__ << "'" << endl;
|
|
|
|
cout << "\n\t Additional global arguments " << endl;
|
|
|
|
g_settings.ShowArguments(cout,"\t");
|
|
|
|
cout << endl;
|
|
|
|
}
|
|
|
|
|
2010-03-26 10:57:53 +00:00
|
|
|
//! \brief processes a command line argument
|
2010-07-16 08:41:12 +00:00
|
|
|
//! \return true if argument is known
|
|
|
|
//! \return false if argument is unknown
|
2010-03-26 10:57:53 +00:00
|
|
|
bool openEMS::parseCommandLineArgument( const char *argv )
|
|
|
|
{
|
|
|
|
if (!argv)
|
|
|
|
return false;
|
|
|
|
|
|
|
|
if (strcmp(argv,"--disable-dumps")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - disabling all field dumps" << endl;
|
|
|
|
SetEnableDumps(false);
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
else if (strcmp(argv,"--debug-material")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - dumping material to 'material_dump.vtk'" << endl;
|
|
|
|
DebugMaterial();
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
else if (strcmp(argv,"--debug-operator")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - dumping operator to 'operator_dump.vtk'" << endl;
|
|
|
|
DebugOperator();
|
|
|
|
return true;
|
|
|
|
}
|
2010-04-19 14:09:41 +00:00
|
|
|
else if (strcmp(argv,"--debug-boxes")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - dumping boxes to 'box_dump*.vtk'" << endl;
|
|
|
|
DebugBox();
|
|
|
|
return true;
|
|
|
|
}
|
2010-06-02 14:37:21 +00:00
|
|
|
else if (strcmp(argv,"--debug-PEC")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - dumping PEC info to 'PEC_dump.vtk'" << endl;
|
2017-05-01 10:34:16 +00:00
|
|
|
DebugPEC();
|
2010-06-02 14:37:21 +00:00
|
|
|
return true;
|
|
|
|
}
|
2010-07-08 09:28:11 +00:00
|
|
|
else if (strcmp(argv,"--debug-CSX")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - dumping CSX geometry to 'debugCSX.xml'" << endl;
|
2017-05-01 10:34:16 +00:00
|
|
|
DebugCSX();
|
2010-07-08 09:28:11 +00:00
|
|
|
return true;
|
|
|
|
}
|
2010-12-06 14:27:48 +00:00
|
|
|
else if (strcmp(argv,"--engine=basic")==0)
|
2010-03-26 11:57:52 +00:00
|
|
|
{
|
2010-12-06 14:27:48 +00:00
|
|
|
cout << "openEMS - enabled basic engine" << endl;
|
|
|
|
m_engine = EngineType_Basic;
|
2010-03-30 11:13:00 +00:00
|
|
|
return true;
|
|
|
|
}
|
2010-04-21 09:18:22 +00:00
|
|
|
else if (strcmp(argv,"--engine=sse")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - enabled sse engine" << endl;
|
|
|
|
m_engine = EngineType_SSE;
|
|
|
|
return true;
|
|
|
|
}
|
2010-05-19 09:41:35 +00:00
|
|
|
else if (strcmp(argv,"--engine=sse-compressed")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - enabled compressed sse engine" << endl;
|
|
|
|
m_engine = EngineType_SSE_Compressed;
|
|
|
|
return true;
|
|
|
|
}
|
2010-12-06 14:27:48 +00:00
|
|
|
else if (strcmp(argv,"--engine=multithreaded")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - enabled multithreading" << endl;
|
|
|
|
m_engine = EngineType_Multithreaded;
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
else if (strncmp(argv,"--numThreads=",13)==0)
|
|
|
|
{
|
2018-01-04 14:51:18 +00:00
|
|
|
this->SetNumberOfThreads(atoi(argv+13));
|
2010-12-06 14:27:48 +00:00
|
|
|
cout << "openEMS - fixed number of threads: " << m_engine_numThreads << endl;
|
|
|
|
return true;
|
|
|
|
}
|
2010-05-21 06:16:24 +00:00
|
|
|
else if (strcmp(argv,"--engine=fastest")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - enabled multithreading engine" << endl;
|
|
|
|
m_engine = EngineType_Multithreaded;
|
|
|
|
return true;
|
|
|
|
}
|
2010-08-16 21:17:19 +00:00
|
|
|
else if (strcmp(argv,"--no-simulation")==0)
|
|
|
|
{
|
|
|
|
cout << "openEMS - disabling simulation => preprocessing only" << endl;
|
|
|
|
m_no_simulation = true;
|
|
|
|
return true;
|
|
|
|
}
|
2012-09-24 13:16:47 +00:00
|
|
|
else if (strcmp(argv,"--dump-statistics")==0)
|
|
|
|
{
|
2012-09-25 09:24:25 +00:00
|
|
|
cout << "openEMS - dump simulation statistics to '" << __OPENEMS_RUN_STAT_FILE__ << "' and '" << __OPENEMS_STAT_FILE__ << "'" << endl;
|
2012-09-24 13:16:47 +00:00
|
|
|
m_DumpStats = true;
|
|
|
|
return true;
|
|
|
|
}
|
2010-03-26 10:57:53 +00:00
|
|
|
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
|
2023-01-06 19:01:07 +00:00
|
|
|
void openEMS::SetNumberOfThreads(int val)
|
|
|
|
{
|
|
|
|
if ((val<0) || (val>boost::thread::hardware_concurrency()))
|
|
|
|
val = boost::thread::hardware_concurrency();
|
|
|
|
m_engine_numThreads = val;
|
|
|
|
}
|
|
|
|
|
2015-12-19 14:38:21 +00:00
|
|
|
string openEMS::GetExtLibsInfo(string prefix)
|
2010-08-25 06:17:45 +00:00
|
|
|
{
|
|
|
|
stringstream str;
|
|
|
|
|
2015-12-19 14:38:21 +00:00
|
|
|
str << prefix << "Used external libraries:" << endl;
|
|
|
|
str << prefix << "\t" << ContinuousStructure::GetInfoLine(true) << endl;
|
2010-08-25 06:17:45 +00:00
|
|
|
|
2010-09-17 13:32:11 +00:00
|
|
|
// libhdf5
|
|
|
|
unsigned int major, minor, release;
|
|
|
|
if (H5get_libversion( &major, &minor, &release ) >= 0)
|
|
|
|
{
|
2015-12-19 14:38:21 +00:00
|
|
|
str << prefix << "\t" << "hdf5 -- Version: " << major << '.' << minor << '.' << release << endl;
|
|
|
|
str << prefix << "\t" << " compiled against: " H5_VERS_INFO << endl;
|
2010-09-17 13:32:11 +00:00
|
|
|
}
|
|
|
|
|
2010-10-04 09:35:20 +00:00
|
|
|
// tinyxml
|
2015-12-19 14:38:21 +00:00
|
|
|
str << prefix << "\t" << "tinyxml -- compiled against: " << TIXML_MAJOR_VERSION << '.' << TIXML_MINOR_VERSION << '.' << TIXML_PATCH_VERSION << endl;
|
2010-10-04 09:35:20 +00:00
|
|
|
|
2011-11-28 10:58:41 +00:00
|
|
|
// fparser
|
2015-12-19 14:38:21 +00:00
|
|
|
str << prefix << "\t" << "fparser" << endl;
|
2011-11-28 10:58:41 +00:00
|
|
|
|
2010-10-04 09:35:20 +00:00
|
|
|
// boost
|
2015-12-19 14:38:21 +00:00
|
|
|
str << prefix << "\t" << "boost -- compiled against: " << BOOST_LIB_VERSION << endl;
|
2011-04-08 10:27:42 +00:00
|
|
|
|
|
|
|
//vtk
|
2015-12-19 14:38:21 +00:00
|
|
|
str << prefix << "\t" << "vtk -- Version: " << vtkVersion::GetVTKMajorVersion() << "." << vtkVersion::GetVTKMinorVersion() << "." << vtkVersion::GetVTKBuildVersion() << endl;
|
|
|
|
str << prefix << "\t" << " compiled against: " << VTK_VERSION << endl;
|
2010-10-04 09:35:20 +00:00
|
|
|
|
2010-08-25 06:17:45 +00:00
|
|
|
return str.str();
|
|
|
|
}
|
|
|
|
|
2015-12-19 14:38:21 +00:00
|
|
|
void openEMS::WelcomeScreen()
|
|
|
|
{
|
|
|
|
#if defined(_LP64) || defined(_WIN64)
|
|
|
|
string bits = "64bit";
|
|
|
|
#else
|
|
|
|
string bits = "32bit";
|
|
|
|
#endif
|
|
|
|
|
|
|
|
cout << " ---------------------------------------------------------------------- " << endl;
|
2023-11-18 11:23:15 +00:00
|
|
|
cout << " | openEMS " << bits << " -- version " << GIT_VERSION << endl;
|
2022-12-30 16:18:53 +00:00
|
|
|
cout << " | (C) 2010-2023 Thorsten Liebig <thorsten.liebig@gmx.de> GPL license" << endl;
|
2015-12-19 14:38:21 +00:00
|
|
|
cout << " ---------------------------------------------------------------------- " << endl;
|
|
|
|
cout << openEMS::GetExtLibsInfo("\t") << endl;
|
|
|
|
}
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
bool openEMS::SetupBoundaryConditions()
|
2010-07-30 13:28:15 +00:00
|
|
|
{
|
2015-12-13 21:33:26 +00:00
|
|
|
FDTD_Op->SetBoundaryCondition(m_BC_type); //operator only knows about PEC and PMC, everything else is defined by extensions (see below)
|
2010-07-30 13:28:15 +00:00
|
|
|
|
|
|
|
/**************************** create all operator/engine extensions here !!!! **********************************/
|
2010-12-06 12:04:37 +00:00
|
|
|
for (int n=0; n<6; ++n)
|
2010-07-30 13:28:15 +00:00
|
|
|
{
|
2012-05-08 11:58:20 +00:00
|
|
|
FDTD_Op->SetBCSize(n, 0);
|
2015-12-13 21:33:26 +00:00
|
|
|
if (m_BC_type[n]==2) //Mur-ABC
|
2010-07-30 13:28:15 +00:00
|
|
|
{
|
2012-05-08 11:58:20 +00:00
|
|
|
FDTD_Op->SetBCSize(n, 1);
|
2010-07-30 13:28:15 +00:00
|
|
|
Operator_Ext_Mur_ABC* op_ext_mur = new Operator_Ext_Mur_ABC(FDTD_Op);
|
|
|
|
op_ext_mur->SetDirection(n/2,n%2);
|
2015-12-15 21:13:54 +00:00
|
|
|
if (m_Mur_v_ph[n]>0)
|
|
|
|
op_ext_mur->SetPhaseVelocity(m_Mur_v_ph[n]);
|
2010-07-30 13:28:15 +00:00
|
|
|
FDTD_Op->AddExtension(op_ext_mur);
|
|
|
|
}
|
2015-12-13 21:33:26 +00:00
|
|
|
if (m_BC_type[n]==3)
|
|
|
|
FDTD_Op->SetBCSize(n, m_PML_size[n]);
|
2010-07-30 13:28:15 +00:00
|
|
|
}
|
2010-10-04 09:52:39 +00:00
|
|
|
|
2015-12-15 21:13:54 +00:00
|
|
|
|
2010-10-04 09:52:39 +00:00
|
|
|
//create the upml
|
2015-12-13 21:33:26 +00:00
|
|
|
Operator_Ext_UPML::Create_UPML(FDTD_Op, m_BC_type, m_PML_size, string());
|
2010-07-30 13:28:15 +00:00
|
|
|
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2013-12-19 14:15:36 +00:00
|
|
|
Engine_Interface_FDTD* openEMS::NewEngineInterface(int multigridlevel)
|
2011-04-13 10:16:54 +00:00
|
|
|
{
|
2013-12-19 14:15:36 +00:00
|
|
|
Operator_CylinderMultiGrid* op_cyl_mg = dynamic_cast<Operator_CylinderMultiGrid*>(FDTD_Op);
|
2014-01-06 15:08:13 +00:00
|
|
|
while (op_cyl_mg && multigridlevel>0)
|
2013-12-19 14:15:36 +00:00
|
|
|
{
|
|
|
|
int mgl = op_cyl_mg->GetMultiGridLevel();
|
|
|
|
if (mgl==multigridlevel)
|
|
|
|
{
|
|
|
|
if (g_settings.GetVerboseLevel()>0)
|
2015-05-04 18:47:19 +00:00
|
|
|
cout << __func__ << ": Operator with requested multi-grid level found." << endl;
|
2014-01-06 15:08:13 +00:00
|
|
|
return new Engine_Interface_Cylindrical_FDTD(op_cyl_mg);
|
2013-12-19 14:15:36 +00:00
|
|
|
}
|
|
|
|
Operator_Cylinder* op_cyl_inner = op_cyl_mg->GetInnerOperator();
|
|
|
|
op_cyl_mg = dynamic_cast<Operator_CylinderMultiGrid*>(op_cyl_inner);
|
|
|
|
if (op_cyl_mg==NULL) //inner most operator reached
|
|
|
|
{
|
|
|
|
if (g_settings.GetVerboseLevel()>0)
|
2015-05-04 18:47:19 +00:00
|
|
|
cout << __func__ << ": Operator with highest multi-grid level chosen." << endl;
|
2014-01-06 15:08:13 +00:00
|
|
|
return new Engine_Interface_Cylindrical_FDTD(op_cyl_inner);
|
2013-12-19 14:15:36 +00:00
|
|
|
}
|
|
|
|
// try next level
|
|
|
|
}
|
|
|
|
Operator_Cylinder* op_cyl = dynamic_cast<Operator_Cylinder*>(FDTD_Op);
|
2014-01-06 15:08:13 +00:00
|
|
|
if (op_cyl)
|
|
|
|
return new Engine_Interface_Cylindrical_FDTD(op_cyl);
|
2011-11-07 11:07:55 +00:00
|
|
|
Operator_sse* op_sse = dynamic_cast<Operator_sse*>(FDTD_Op);
|
2014-01-06 15:08:13 +00:00
|
|
|
if (op_sse)
|
|
|
|
return new Engine_Interface_SSE_FDTD(op_sse);
|
|
|
|
return new Engine_Interface_FDTD(FDTD_Op);
|
2011-04-13 10:16:54 +00:00
|
|
|
}
|
|
|
|
|
2015-12-19 14:03:17 +00:00
|
|
|
void openEMS::SetVerboseLevel(int level)
|
|
|
|
{
|
|
|
|
g_settings.SetVerboseLevel(level);
|
|
|
|
}
|
|
|
|
|
2011-01-31 11:00:00 +00:00
|
|
|
bool openEMS::SetupProcessing()
|
2011-01-10 07:27:50 +00:00
|
|
|
{
|
|
|
|
//*************** setup processing ************//
|
2011-11-08 10:49:14 +00:00
|
|
|
if (g_settings.GetVerboseLevel()>0)
|
|
|
|
cout << "Setting up processing..." << endl;
|
2011-01-10 07:27:50 +00:00
|
|
|
|
2012-07-16 15:15:10 +00:00
|
|
|
unsigned int Nyquist = FDTD_Op->GetExcitationSignal()->GetNyquistNum();
|
2011-01-10 07:27:50 +00:00
|
|
|
PA = new ProcessingArray(Nyquist);
|
|
|
|
|
|
|
|
double start[3];
|
|
|
|
double stop[3];
|
2011-08-31 09:05:04 +00:00
|
|
|
bool l_MultiBox = false;
|
2011-01-31 11:00:00 +00:00
|
|
|
vector<CSProperties*> Probes = m_CSX->GetPropertyByType(CSProperties::PROBEBOX);
|
2011-01-10 07:27:50 +00:00
|
|
|
for (size_t i=0; i<Probes.size(); ++i)
|
|
|
|
{
|
2021-08-25 17:05:38 +00:00
|
|
|
CSPropProbeBox* pb = Probes.at(i)->ToProbeBox();
|
|
|
|
if (!pb)
|
|
|
|
continue;
|
2011-08-31 09:05:04 +00:00
|
|
|
//check whether one or more probe boxes are defined
|
2021-08-25 17:05:38 +00:00
|
|
|
l_MultiBox = (pb->GetQtyPrimitives()>1);
|
2011-08-31 09:05:04 +00:00
|
|
|
|
2021-08-25 17:05:38 +00:00
|
|
|
for (size_t nb=0; nb<pb->GetQtyPrimitives(); ++nb)
|
2011-01-10 07:27:50 +00:00
|
|
|
{
|
2021-08-25 17:05:38 +00:00
|
|
|
CSPrimitives* prim = pb->GetPrimitive(nb);
|
2011-08-31 09:05:04 +00:00
|
|
|
if (prim!=NULL)
|
2011-01-10 07:27:50 +00:00
|
|
|
{
|
2011-08-31 09:05:04 +00:00
|
|
|
double bnd[6] = {0,0,0,0,0,0};
|
2011-11-07 11:04:34 +00:00
|
|
|
prim->GetBoundBox(bnd,true);
|
2011-08-31 09:05:04 +00:00
|
|
|
start[0]= bnd[0];
|
|
|
|
start[1]=bnd[2];
|
|
|
|
start[2]=bnd[4];
|
|
|
|
stop[0] = bnd[1];
|
|
|
|
stop[1] =bnd[3];
|
|
|
|
stop[2] =bnd[5];
|
2021-08-25 17:05:38 +00:00
|
|
|
|
2013-01-22 08:22:01 +00:00
|
|
|
ProcessIntegral* proc = NULL;
|
2021-08-25 17:05:38 +00:00
|
|
|
if (pb->GetProbeType()==0)
|
2011-01-10 07:27:50 +00:00
|
|
|
{
|
2021-08-25 17:05:38 +00:00
|
|
|
CSPrimBox* box = prim->ToBox();
|
2022-12-29 12:08:44 +00:00
|
|
|
if (!(box) || box->GetDimension()!=1)
|
2011-08-31 09:05:04 +00:00
|
|
|
{
|
2021-08-25 17:05:38 +00:00
|
|
|
cerr << "openEMS::SetupProcessing: Error: Probe primitive type or dimension not suitable ... skipping probe " << pb->GetName() << endl;
|
2011-08-31 09:05:04 +00:00
|
|
|
continue;
|
|
|
|
}
|
2021-08-25 17:05:38 +00:00
|
|
|
// use the direction and coordinates of the box
|
|
|
|
for (int n=0;n<3;++n)
|
2011-08-31 09:05:04 +00:00
|
|
|
{
|
2021-08-25 17:05:38 +00:00
|
|
|
start[n] = box->GetCoord(2*n);
|
|
|
|
stop[n] = box->GetCoord(2*n+1);
|
2011-08-31 09:05:04 +00:00
|
|
|
}
|
2021-08-25 17:05:38 +00:00
|
|
|
ProcessVoltage* procVolt = new ProcessVoltage(NewEngineInterface());
|
|
|
|
proc=procVolt;
|
|
|
|
}
|
|
|
|
else if (pb->GetProbeType()==1)
|
|
|
|
{
|
|
|
|
ProcessCurrent* procCurr = new ProcessCurrent(NewEngineInterface());
|
|
|
|
proc=procCurr;
|
|
|
|
}
|
|
|
|
else if (pb->GetProbeType()==2)
|
|
|
|
proc = new ProcessFieldProbe(NewEngineInterface(),0);
|
|
|
|
else if (pb->GetProbeType()==3)
|
|
|
|
proc = new ProcessFieldProbe(NewEngineInterface(),1);
|
|
|
|
else if ((pb->GetProbeType()==10) || (pb->GetProbeType()==11))
|
|
|
|
{
|
|
|
|
ProcessModeMatch* pmm = new ProcessModeMatch(NewEngineInterface());
|
|
|
|
pmm->SetFieldType(pb->GetProbeType()-10);
|
|
|
|
pmm->SetModeFunction(0,pb->GetAttributeValue("ModeFunctionX"));
|
|
|
|
pmm->SetModeFunction(1,pb->GetAttributeValue("ModeFunctionY"));
|
|
|
|
pmm->SetModeFunction(2,pb->GetAttributeValue("ModeFunctionZ"));
|
|
|
|
proc = pmm;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
cerr << "openEMS::SetupFDTD: Warning: Probe type " << pb->GetProbeType() << " of property '" << pb->GetName() << "' is unknown..." << endl;
|
|
|
|
continue;
|
|
|
|
}
|
|
|
|
if (CylinderCoords)
|
|
|
|
proc->SetMeshType(Processing::CYLINDRICAL_MESH);
|
|
|
|
if ((pb->GetProbeType()==1) || (pb->GetProbeType()==3))
|
|
|
|
{
|
|
|
|
proc->SetDualTime(true);
|
|
|
|
proc->SetDualMesh(true);
|
2011-01-10 07:27:50 +00:00
|
|
|
}
|
2021-08-25 17:05:38 +00:00
|
|
|
if (pb->GetProbeType()==11)
|
|
|
|
proc->SetDualTime(true);
|
|
|
|
proc->SetProcessInterval(Nyquist/m_OverSampling);
|
|
|
|
if (pb->GetStartTime()>0 || pb->GetStopTime()>0)
|
|
|
|
proc->SetProcessStartStopTime(pb->GetStartTime(), pb->GetStopTime());
|
|
|
|
proc->AddFrequency(pb->GetFDSamples());
|
|
|
|
proc->GetNormalDir(pb->GetNormalDir());
|
|
|
|
if (l_MultiBox==false)
|
|
|
|
proc->SetName(pb->GetName());
|
2011-01-10 07:27:50 +00:00
|
|
|
else
|
2021-08-25 17:05:38 +00:00
|
|
|
proc->SetName(pb->GetName(),nb);
|
|
|
|
proc->DefineStartStopCoord(start,stop);
|
|
|
|
if (g_settings.showProbeDiscretization())
|
|
|
|
proc->ShowSnappedCoords();
|
|
|
|
proc->SetWeight(pb->GetWeighting());
|
|
|
|
PA->AddProcessing(proc);
|
|
|
|
prim->SetPrimitiveUsed(true);
|
2011-01-10 07:27:50 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2011-01-31 11:00:00 +00:00
|
|
|
vector<CSProperties*> DumpProps = m_CSX->GetPropertyByType(CSProperties::DUMPBOX);
|
2011-01-10 07:27:50 +00:00
|
|
|
for (size_t i=0; i<DumpProps.size(); ++i)
|
|
|
|
{
|
|
|
|
ProcessFields* ProcField=NULL;
|
|
|
|
|
2011-08-31 09:05:04 +00:00
|
|
|
//check whether one or more probe boxes are defined
|
|
|
|
l_MultiBox = (DumpProps.at(i)->GetQtyPrimitives()>1);
|
|
|
|
|
|
|
|
for (size_t nb=0; nb<DumpProps.at(i)->GetQtyPrimitives(); ++nb)
|
2011-01-10 07:27:50 +00:00
|
|
|
{
|
2011-08-31 09:05:04 +00:00
|
|
|
|
|
|
|
CSPrimitives* prim = DumpProps.at(i)->GetPrimitive(nb);
|
|
|
|
if (prim!=NULL)
|
2011-01-10 07:27:50 +00:00
|
|
|
{
|
2011-08-31 09:05:04 +00:00
|
|
|
double bnd[6] = {0,0,0,0,0,0};
|
2011-11-07 11:04:34 +00:00
|
|
|
prim->GetBoundBox(bnd,true);
|
2011-08-31 09:05:04 +00:00
|
|
|
start[0]= bnd[0];
|
|
|
|
start[1]=bnd[2];
|
|
|
|
start[2]=bnd[4];
|
|
|
|
stop[0] = bnd[1];
|
|
|
|
stop[1] =bnd[3];
|
|
|
|
stop[2] =bnd[5];
|
|
|
|
CSPropDumpBox* db = DumpProps.at(i)->ToDumpBox();
|
|
|
|
if (db)
|
2011-01-10 07:27:50 +00:00
|
|
|
{
|
2017-05-28 10:01:04 +00:00
|
|
|
if ((db->GetDumpType()>=0) && (db->GetDumpType()<=5))
|
2013-12-19 14:15:36 +00:00
|
|
|
ProcField = new ProcessFieldsTD(NewEngineInterface(db->GetMultiGridLevel()));
|
2017-05-28 10:01:04 +00:00
|
|
|
else if ((db->GetDumpType()>=10) && (db->GetDumpType()<=15))
|
2013-12-19 14:15:36 +00:00
|
|
|
ProcField = new ProcessFieldsFD(NewEngineInterface(db->GetMultiGridLevel()));
|
2012-11-29 15:45:48 +00:00
|
|
|
else if ( ((db->GetDumpType()>=20) && (db->GetDumpType()<=22)) || (db->GetDumpType()==29) )
|
2013-05-15 14:02:30 +00:00
|
|
|
{
|
2013-12-19 14:15:36 +00:00
|
|
|
ProcessFieldsSAR* procSAR = new ProcessFieldsSAR(NewEngineInterface(db->GetMultiGridLevel()));
|
2013-05-15 14:02:30 +00:00
|
|
|
ProcField = procSAR;
|
|
|
|
string method = db->GetAttributeValue("SAR_Method");
|
|
|
|
if (!method.empty())
|
|
|
|
procSAR->SetSARAveragingMethod(method);
|
2013-08-23 15:28:18 +00:00
|
|
|
// use (center)-cell based conductivity only
|
|
|
|
procSAR->SetUseCellConductivity(true);
|
2013-05-15 14:02:30 +00:00
|
|
|
}
|
2011-01-10 07:27:50 +00:00
|
|
|
else
|
2011-08-31 09:05:04 +00:00
|
|
|
cerr << "openEMS::SetupFDTD: unknown dump box type... skipping!" << endl;
|
|
|
|
if (ProcField)
|
2011-01-31 11:25:55 +00:00
|
|
|
{
|
2011-08-31 09:05:04 +00:00
|
|
|
ProcField->SetEnable(Enable_Dumps);
|
|
|
|
ProcField->SetProcessInterval(Nyquist/m_OverSampling);
|
2015-03-06 20:38:27 +00:00
|
|
|
if (db->GetStopTime()>0 || db->GetStartTime()>0)
|
|
|
|
ProcField->SetProcessStartStopTime(db->GetStartTime(), db->GetStopTime());
|
2011-08-31 09:05:04 +00:00
|
|
|
if ((db->GetDumpType()==1) || (db->GetDumpType()==11))
|
|
|
|
{
|
|
|
|
ProcField->SetDualTime(true);
|
|
|
|
//make dualMesh the default mesh for h-field dumps, maybe overwritten by interpolation type (node-interpolation)
|
|
|
|
ProcField->SetDualMesh(true);
|
|
|
|
}
|
|
|
|
if (db->GetDumpType()>=10)
|
|
|
|
{
|
|
|
|
ProcField->AddFrequency(db->GetFDSamples());
|
|
|
|
ProcField->SetDumpType((ProcessFields::DumpType)(db->GetDumpType()-10));
|
|
|
|
}
|
|
|
|
else
|
|
|
|
ProcField->SetDumpType((ProcessFields::DumpType)db->GetDumpType());
|
|
|
|
|
|
|
|
if (db->GetDumpType()==20)
|
|
|
|
ProcField->SetDumpType(ProcessFields::SAR_LOCAL_DUMP);
|
2012-11-29 15:45:48 +00:00
|
|
|
if (db->GetDumpType()==21)
|
|
|
|
ProcField->SetDumpType(ProcessFields::SAR_1G_DUMP);
|
|
|
|
if (db->GetDumpType()==22)
|
|
|
|
ProcField->SetDumpType(ProcessFields::SAR_10G_DUMP);
|
|
|
|
if (db->GetDumpType()==29)
|
|
|
|
ProcField->SetDumpType(ProcessFields::SAR_RAW_DATA);
|
2011-08-31 09:05:04 +00:00
|
|
|
|
|
|
|
//SetupMaterialStorages() has previewed storage needs... refresh here to prevent cleanup!!!
|
2017-05-28 10:01:04 +00:00
|
|
|
if ( ProcField->NeedPermittivity() && Enable_Dumps)
|
|
|
|
FDTD_Op->SetMaterialStoreFlags(0,true);
|
|
|
|
if ( ProcField->NeedConductivity() && Enable_Dumps)
|
2011-08-31 09:05:04 +00:00
|
|
|
FDTD_Op->SetMaterialStoreFlags(1,true);
|
2017-05-28 10:01:04 +00:00
|
|
|
if ( ProcField->NeedPermeability() && Enable_Dumps)
|
|
|
|
FDTD_Op->SetMaterialStoreFlags(2,true);
|
2011-08-31 09:05:04 +00:00
|
|
|
|
|
|
|
ProcField->SetDumpMode((Engine_Interface_Base::InterpolationType)db->GetDumpMode());
|
|
|
|
ProcField->SetFileType((ProcessFields::FileType)db->GetFileType());
|
|
|
|
if (CylinderCoords)
|
|
|
|
ProcField->SetMeshType(Processing::CYLINDRICAL_MESH);
|
|
|
|
if (db->GetSubSampling())
|
|
|
|
for (int n=0; n<3; ++n)
|
|
|
|
ProcField->SetSubSampling(db->GetSubSampling(n),n);
|
|
|
|
if (db->GetOptResolution())
|
|
|
|
for (int n=0; n<3; ++n)
|
|
|
|
ProcField->SetOptResolution(db->GetOptResolution(n),n);
|
|
|
|
|
|
|
|
if (l_MultiBox==false)
|
|
|
|
ProcField->SetName(db->GetName());
|
|
|
|
else
|
|
|
|
ProcField->SetName(db->GetName(),nb);
|
|
|
|
|
|
|
|
ProcField->SetFileName(ProcField->GetName());
|
|
|
|
ProcField->DefineStartStopCoord(start,stop);
|
|
|
|
if (g_settings.showProbeDiscretization())
|
|
|
|
ProcField->ShowSnappedCoords();
|
|
|
|
PA->AddProcessing(ProcField);
|
|
|
|
prim->SetPrimitiveUsed(true);
|
2011-01-31 11:25:55 +00:00
|
|
|
}
|
2011-01-10 07:27:50 +00:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2011-01-31 11:00:00 +00:00
|
|
|
bool openEMS::SetupMaterialStorages()
|
2011-01-10 07:27:50 +00:00
|
|
|
{
|
2011-01-31 11:00:00 +00:00
|
|
|
vector<CSProperties*> DumpProps = m_CSX->GetPropertyByType(CSProperties::DUMPBOX);
|
2011-01-10 07:27:50 +00:00
|
|
|
for (size_t i=0; i<DumpProps.size(); ++i)
|
|
|
|
{
|
|
|
|
CSPropDumpBox* db = DumpProps.at(i)->ToDumpBox();
|
|
|
|
if (!db)
|
|
|
|
continue;
|
|
|
|
if (db->GetQtyPrimitives()==0)
|
|
|
|
continue;
|
|
|
|
//check for current density dump types
|
2013-01-30 13:15:35 +00:00
|
|
|
if ( ((db->GetDumpType()==2) || (db->GetDumpType()==12) || // current density storage
|
2013-02-06 15:47:18 +00:00
|
|
|
(db->GetDumpType()==20) || (db->GetDumpType()==21) || (db->GetDumpType()==22)) && // SAR dump types
|
2013-01-30 13:15:35 +00:00
|
|
|
Enable_Dumps )
|
2011-01-10 07:27:50 +00:00
|
|
|
FDTD_Op->SetMaterialStoreFlags(1,true); //tell operator to store kappa material data
|
2017-05-28 10:01:04 +00:00
|
|
|
if ( ((db->GetDumpType()==4) || (db->GetDumpType()==14)) || Enable_Dumps) // electric flux density storage
|
|
|
|
FDTD_Op->SetMaterialStoreFlags(0,true); //tell operator to store epsR material data
|
|
|
|
if ( ((db->GetDumpType()==5) || (db->GetDumpType()==15)) || Enable_Dumps) // magnetic flux density storage
|
|
|
|
FDTD_Op->SetMaterialStoreFlags(2,true); //tell operator to store mueR material data
|
2011-01-10 07:27:50 +00:00
|
|
|
}
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
void openEMS::SetupCylinderMultiGrid(std::string val)
|
|
|
|
{
|
|
|
|
m_CC_MultiGrid.clear();
|
|
|
|
m_CC_MultiGrid = SplitString2Double(val,',');
|
|
|
|
}
|
|
|
|
|
|
|
|
bool openEMS::SetupOperator()
|
2011-02-11 12:52:48 +00:00
|
|
|
{
|
|
|
|
if (CylinderCoords)
|
|
|
|
{
|
2015-12-13 21:33:26 +00:00
|
|
|
if (m_CC_MultiGrid.size()>0)
|
2011-02-11 12:52:48 +00:00
|
|
|
{
|
2015-12-13 21:33:26 +00:00
|
|
|
FDTD_Op = Operator_CylinderMultiGrid::New(m_CC_MultiGrid, m_engine_numThreads);
|
2011-03-14 09:37:12 +00:00
|
|
|
if (FDTD_Op==NULL)
|
|
|
|
FDTD_Op = Operator_Cylinder::New(m_engine_numThreads);
|
2011-02-11 12:52:48 +00:00
|
|
|
}
|
|
|
|
else
|
|
|
|
FDTD_Op = Operator_Cylinder::New(m_engine_numThreads);
|
|
|
|
}
|
|
|
|
else if (m_engine == EngineType_SSE)
|
|
|
|
{
|
|
|
|
FDTD_Op = Operator_sse::New();
|
|
|
|
}
|
|
|
|
else if (m_engine == EngineType_SSE_Compressed)
|
|
|
|
{
|
|
|
|
FDTD_Op = Operator_SSE_Compressed::New();
|
|
|
|
}
|
|
|
|
else if (m_engine == EngineType_Multithreaded)
|
|
|
|
{
|
|
|
|
FDTD_Op = Operator_Multithread::New(m_engine_numThreads);
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
FDTD_Op = Operator::New();
|
|
|
|
}
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
void openEMS::Set_BC_Type(int idx, int type)
|
|
|
|
{
|
|
|
|
if ((idx<0) || (idx>5))
|
|
|
|
return;
|
|
|
|
m_BC_type[idx] = type;
|
|
|
|
}
|
|
|
|
|
2015-12-30 13:17:40 +00:00
|
|
|
int openEMS::Get_BC_Type(int idx)
|
|
|
|
{
|
|
|
|
if ((idx<0) || (idx>5))
|
|
|
|
return -1;
|
|
|
|
return m_BC_type[idx];
|
|
|
|
}
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
void openEMS::Set_BC_PML(int idx, unsigned int size)
|
|
|
|
{
|
|
|
|
if ((idx<0) || (idx>5))
|
|
|
|
return;
|
|
|
|
m_BC_type[idx] = 3;
|
|
|
|
m_PML_size[idx] = size;
|
|
|
|
}
|
2011-01-10 07:27:50 +00:00
|
|
|
|
2015-12-30 13:17:40 +00:00
|
|
|
int openEMS::Get_PML_Size(int idx)
|
|
|
|
{
|
|
|
|
if ((idx<0) || (idx>5))
|
|
|
|
return -1;
|
|
|
|
if (m_BC_type[idx]!=3)
|
|
|
|
return -1; // return -1 if BC was *not* a PML
|
|
|
|
return m_PML_size[idx];
|
|
|
|
}
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
void openEMS::Set_Mur_PhaseVel(int idx, double val)
|
|
|
|
{
|
|
|
|
if ((idx<0) || (idx>5))
|
|
|
|
return;
|
|
|
|
m_Mur_v_ph[idx] = val;
|
|
|
|
}
|
|
|
|
|
|
|
|
bool openEMS::ParseFDTDSetup(std::string file)
|
2010-03-11 09:56:19 +00:00
|
|
|
{
|
|
|
|
Reset();
|
|
|
|
|
2011-11-08 10:49:14 +00:00
|
|
|
if (g_settings.GetVerboseLevel()>0)
|
|
|
|
cout << "Read openEMS xml file: " << file << " ..." << endl;
|
2010-10-06 14:27:31 +00:00
|
|
|
|
2010-03-11 09:56:19 +00:00
|
|
|
TiXmlDocument doc(file);
|
|
|
|
if (!doc.LoadFile())
|
|
|
|
{
|
|
|
|
cerr << "openEMS: Error File-Loading failed!!! File: " << file << endl;
|
|
|
|
exit(-1);
|
|
|
|
}
|
|
|
|
|
2011-11-08 10:49:14 +00:00
|
|
|
if (g_settings.GetVerboseLevel()>0)
|
|
|
|
cout << "Read openEMS Settings..." << endl;
|
2010-03-12 20:14:17 +00:00
|
|
|
TiXmlElement* openEMSxml = doc.FirstChildElement("openEMS");
|
|
|
|
if (openEMSxml==NULL)
|
|
|
|
{
|
|
|
|
cerr << "Can't read openEMS ... " << endl;
|
|
|
|
exit(-1);
|
|
|
|
}
|
|
|
|
TiXmlElement* FDTD_Opts = openEMSxml->FirstChildElement("FDTD");
|
2010-04-07 14:31:23 +00:00
|
|
|
|
2010-03-11 09:56:19 +00:00
|
|
|
if (FDTD_Opts==NULL)
|
|
|
|
{
|
2010-03-12 20:14:17 +00:00
|
|
|
cerr << "Can't read openEMS FDTD Settings... " << endl;
|
2010-03-11 09:56:19 +00:00
|
|
|
exit(-1);
|
|
|
|
}
|
2016-01-27 17:37:07 +00:00
|
|
|
|
|
|
|
if (g_settings.GetVerboseLevel()>0)
|
|
|
|
cout << "Read Geometry..." << endl;
|
|
|
|
ContinuousStructure* csx = new ContinuousStructure();
|
|
|
|
string EC(csx->ReadFromXML(openEMSxml));
|
|
|
|
if (EC.empty()==false)
|
|
|
|
cerr << EC << endl;
|
|
|
|
this->SetCSX(csx);
|
|
|
|
|
|
|
|
|
|
|
|
return this->Parse_XML_FDTDSetup(FDTD_Opts);
|
|
|
|
}
|
|
|
|
|
|
|
|
bool openEMS::Parse_XML_FDTDSetup(TiXmlElement* FDTD_Opts)
|
|
|
|
{
|
2012-01-17 10:35:00 +00:00
|
|
|
double dhelp=0;
|
|
|
|
FDTD_Opts->QueryDoubleAttribute("NumberOfTimesteps",&dhelp);
|
|
|
|
if (dhelp<0)
|
2015-12-13 21:33:26 +00:00
|
|
|
this->SetNumberOfTimeSteps(0);
|
2010-03-29 08:12:38 +00:00
|
|
|
else
|
2015-12-13 21:33:26 +00:00
|
|
|
this->SetNumberOfTimeSteps((unsigned int)dhelp);
|
2010-04-07 14:31:23 +00:00
|
|
|
|
2012-01-17 10:35:00 +00:00
|
|
|
int ihelp = 0;
|
|
|
|
FDTD_Opts->QueryIntAttribute("CylinderCoords",&ihelp);
|
|
|
|
if (ihelp==1)
|
2015-12-15 21:13:54 +00:00
|
|
|
{
|
2015-12-13 21:33:26 +00:00
|
|
|
this->SetCylinderCoords(true);
|
|
|
|
const char* cchelp = FDTD_Opts->Attribute("MultiGrid");
|
|
|
|
if (cchelp!=NULL)
|
|
|
|
this->SetupCylinderMultiGrid(string(cchelp));
|
2015-12-15 21:13:54 +00:00
|
|
|
}
|
2015-12-13 21:33:26 +00:00
|
|
|
|
2019-02-13 18:26:27 +00:00
|
|
|
dhelp = 0;
|
|
|
|
FDTD_Opts->QueryDoubleAttribute("MaxTime",&dhelp);
|
|
|
|
if (dhelp>0)
|
|
|
|
this->SetMaxTime(dhelp);
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
dhelp = 0;
|
|
|
|
FDTD_Opts->QueryDoubleAttribute("endCriteria",&dhelp);
|
|
|
|
if (dhelp==0)
|
|
|
|
this->SetEndCriteria(1e-6);
|
|
|
|
else
|
|
|
|
this->SetEndCriteria(dhelp);
|
2010-04-09 13:51:37 +00:00
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
ihelp = 0;
|
|
|
|
FDTD_Opts->QueryIntAttribute("OverSampling",&ihelp);
|
2016-07-30 20:04:15 +00:00
|
|
|
if (ihelp>1)
|
2015-12-13 21:33:26 +00:00
|
|
|
this->SetOverSampling(ihelp);
|
2010-03-11 09:56:19 +00:00
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
// check for cell constant material averaging
|
|
|
|
if (FDTD_Opts->QueryIntAttribute("CellConstantMaterial",&ihelp)==TIXML_SUCCESS)
|
|
|
|
this->SetCellConstantMaterial(ihelp==1);
|
2010-04-04 17:48:36 +00:00
|
|
|
|
2010-03-11 09:56:19 +00:00
|
|
|
TiXmlElement* BC = FDTD_Opts->FirstChildElement("BoundaryCond");
|
|
|
|
if (BC==NULL)
|
|
|
|
{
|
|
|
|
cerr << "Can't read openEMS boundary cond Settings... " << endl;
|
|
|
|
exit(-3);
|
|
|
|
}
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
// const char* tmp = BC->Attribute("PML_Grading");
|
|
|
|
// string pml_gradFunc;
|
|
|
|
// if (tmp)
|
|
|
|
// pml_gradFunc = string(tmp);
|
|
|
|
|
|
|
|
string bound_names[] = {"xmin","xmax","ymin","ymax","zmin","zmax"};
|
2015-12-15 21:13:54 +00:00
|
|
|
string s_bc;
|
2015-12-13 21:33:26 +00:00
|
|
|
for (int n=0; n<6; ++n)
|
|
|
|
{
|
|
|
|
int EC = BC->QueryIntAttribute(bound_names[n].c_str(),&ihelp);
|
|
|
|
if (EC==TIXML_SUCCESS)
|
2015-12-15 21:13:54 +00:00
|
|
|
{
|
2015-12-13 21:33:26 +00:00
|
|
|
this->Set_BC_Type(n, ihelp);
|
|
|
|
continue;
|
2015-12-15 21:13:54 +00:00
|
|
|
}
|
2015-12-13 21:33:26 +00:00
|
|
|
if (EC==TIXML_WRONG_TYPE)
|
|
|
|
{
|
|
|
|
const char* tmp = BC->Attribute(bound_names[n].c_str());
|
|
|
|
if (tmp)
|
|
|
|
s_bc = string(tmp);
|
2015-12-15 21:13:54 +00:00
|
|
|
else
|
|
|
|
cerr << "openEMS::SetupBoundaryConditions: Warning, boundary condition for \"" << bound_names[n] << "\" unknown... set to PEC " << endl;
|
2015-12-13 21:33:26 +00:00
|
|
|
if (s_bc=="PEC")
|
|
|
|
this->Set_BC_Type(n, 0);
|
|
|
|
else if (s_bc=="PMC")
|
|
|
|
this->Set_BC_Type(n, 1);
|
|
|
|
else if (s_bc=="MUR")
|
|
|
|
this->Set_BC_Type(n, 2);
|
|
|
|
else if (strncmp(s_bc.c_str(),"PML_=",4)==0)
|
|
|
|
this->Set_BC_PML(n, atoi(s_bc.c_str()+4));
|
|
|
|
else
|
|
|
|
cerr << "openEMS::SetupBoundaryConditions: Warning, boundary condition for \"" << bound_names[n] << "\" unknown... set to PEC " << endl;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
cerr << "openEMS::SetupBoundaryConditions: Warning, boundary condition for \"" << bound_names[n] << "\" not found... set to PEC " << endl;
|
|
|
|
}
|
|
|
|
|
|
|
|
//read general mur phase velocity
|
2015-12-15 21:13:54 +00:00
|
|
|
if (BC->QueryDoubleAttribute("MUR_PhaseVelocity",&dhelp) == TIXML_SUCCESS)
|
2015-12-13 21:33:26 +00:00
|
|
|
for (int n=0;n<6;++n)
|
|
|
|
this->Set_Mur_PhaseVel(n, dhelp);
|
|
|
|
|
|
|
|
string mur_v_ph_names[6] = {"MUR_PhaseVelocity_xmin", "MUR_PhaseVelocity_xmax", "MUR_PhaseVelocity_ymin", "MUR_PhaseVelocity_ymax", "MUR_PhaseVelocity_zmin", "MUR_PhaseVelocity_zmax"};
|
|
|
|
for (int n=0; n<6; ++n)
|
|
|
|
if (BC->QueryDoubleAttribute(mur_v_ph_names[n].c_str(),&dhelp) == TIXML_SUCCESS)
|
|
|
|
this->Set_Mur_PhaseVel(n, dhelp);
|
|
|
|
|
|
|
|
TiXmlElement* m_Excite_Elem = FDTD_Opts->FirstChildElement("Excitation");
|
|
|
|
if (!m_Excite_Elem)
|
|
|
|
{
|
|
|
|
cerr << "Excitation::setupExcitation: Error, can't read openEMS excitation settings... " << endl;
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
|
|
|
|
Excitation* exc = this->InitExcitation();
|
|
|
|
double f0=0, fc=0, f_max=0;
|
|
|
|
ihelp = -1;
|
|
|
|
m_Excite_Elem->QueryIntAttribute("Type",&ihelp);
|
|
|
|
switch (ihelp)
|
|
|
|
{
|
|
|
|
case Excitation::GaissianPulse:
|
|
|
|
m_Excite_Elem->QueryDoubleAttribute("f0",&f0);
|
|
|
|
m_Excite_Elem->QueryDoubleAttribute("fc",&fc);
|
|
|
|
exc->SetupGaussianPulse(f0, fc);
|
|
|
|
break;
|
|
|
|
case Excitation::Sinusoidal: // sinusoidal excite
|
|
|
|
m_Excite_Elem->QueryDoubleAttribute("f0",&f0);
|
|
|
|
exc->SetupSinusoidal(f0);
|
|
|
|
break;
|
|
|
|
case Excitation::DiracPulse:
|
|
|
|
FDTD_Opts->QueryDoubleAttribute("f_max",&f_max);
|
|
|
|
exc->SetupDiracPulse(f_max);
|
|
|
|
break;
|
|
|
|
case Excitation::Step:
|
|
|
|
FDTD_Opts->QueryDoubleAttribute("f_max",&f_max);
|
|
|
|
exc->SetupStepExcite(f_max);
|
|
|
|
break;
|
|
|
|
case Excitation::CustomExcite:
|
|
|
|
m_Excite_Elem->QueryDoubleAttribute("f0",&f0);
|
|
|
|
FDTD_Opts->QueryDoubleAttribute("f_max",&f_max);
|
|
|
|
exc->SetupCustomExcite(m_Excite_Elem->Attribute("Function"), f0, f_max);
|
|
|
|
break;
|
2010-03-11 09:56:19 +00:00
|
|
|
}
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
if (FDTD_Opts->QueryIntAttribute("TimeStepMethod",&ihelp)==TIXML_SUCCESS)
|
|
|
|
this->SetTimeStepMethod(ihelp);
|
|
|
|
if (FDTD_Opts->QueryDoubleAttribute("TimeStep",&dhelp)==TIXML_SUCCESS)
|
|
|
|
this->SetTimeStep(dhelp);
|
|
|
|
if (FDTD_Opts->QueryDoubleAttribute("TimeStepFactor",&dhelp)==TIXML_SUCCESS)
|
|
|
|
this->SetTimeStepFactor(dhelp);
|
2016-12-02 18:03:35 +00:00
|
|
|
return true;
|
2015-12-13 21:33:26 +00:00
|
|
|
}
|
|
|
|
|
2015-12-19 14:03:02 +00:00
|
|
|
void openEMS::SetGaussExcite(double f0, double fc)
|
|
|
|
{
|
|
|
|
this->InitExcitation();
|
|
|
|
m_Exc->SetupGaussianPulse(f0, fc);
|
|
|
|
}
|
|
|
|
|
2021-12-21 12:41:16 +00:00
|
|
|
void openEMS::SetSinusExcite(double f0)
|
|
|
|
{
|
|
|
|
this->InitExcitation();
|
|
|
|
m_Exc->SetupSinusoidal(f0);
|
|
|
|
}
|
|
|
|
|
|
|
|
void openEMS::SetDiracExcite(double f_max)
|
|
|
|
{
|
|
|
|
this->InitExcitation();
|
|
|
|
m_Exc->SetupDiracPulse(f_max);
|
|
|
|
}
|
|
|
|
|
|
|
|
void openEMS::SetStepExcite(double f_max)
|
|
|
|
{
|
|
|
|
this->InitExcitation();
|
|
|
|
m_Exc->SetupStepExcite(f_max);
|
|
|
|
}
|
|
|
|
|
2023-12-16 11:13:56 +00:00
|
|
|
void openEMS::SetCustomExcite(std::string str, double f0, double fmax)
|
|
|
|
{
|
|
|
|
this->InitExcitation();
|
|
|
|
m_Exc->SetupCustomExcite(str, f0, fmax);
|
|
|
|
}
|
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
Excitation* openEMS::InitExcitation()
|
|
|
|
{
|
|
|
|
delete m_Exc;
|
|
|
|
m_Exc = new Excitation();
|
|
|
|
return m_Exc;
|
|
|
|
}
|
|
|
|
|
|
|
|
void openEMS::SetCSX(ContinuousStructure* csx)
|
|
|
|
{
|
|
|
|
delete m_CSX;
|
|
|
|
m_CSX = csx;
|
|
|
|
}
|
|
|
|
|
|
|
|
int openEMS::SetupFDTD()
|
|
|
|
{
|
|
|
|
timeval startTime;
|
|
|
|
gettimeofday(&startTime,NULL);
|
|
|
|
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_EXIT_FORCE);
|
|
|
|
|
2015-12-30 13:17:56 +00:00
|
|
|
if (m_CSX==NULL)
|
|
|
|
{
|
|
|
|
cerr << "openEMS::SetupFDTD: Error: CSXCAD is not set!" << endl;
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
2015-12-30 13:17:56 +00:00
|
|
|
return 3;
|
|
|
|
}
|
|
|
|
if (m_CSX==NULL)
|
|
|
|
{
|
|
|
|
cerr << "openEMS::SetupFDTD: Error: CSXCAD is not set!" << endl;
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
2015-12-30 13:17:56 +00:00
|
|
|
return 3;
|
|
|
|
}
|
|
|
|
std::string ec = m_CSX->Update();
|
|
|
|
if (!ec.empty())
|
|
|
|
cerr << ec << endl;
|
2011-11-16 10:24:25 +00:00
|
|
|
if (g_settings.GetVerboseLevel()>2)
|
2011-03-14 07:50:28 +00:00
|
|
|
m_CSX->ShowPropertyStatus(cerr);
|
|
|
|
|
2010-11-02 15:32:00 +00:00
|
|
|
if (CylinderCoords)
|
2011-01-31 11:00:00 +00:00
|
|
|
if (m_CSX->GetCoordInputType()!=CYLINDRICAL)
|
2010-11-02 15:32:00 +00:00
|
|
|
{
|
2010-12-06 12:04:37 +00:00
|
|
|
cerr << "openEMS::SetupFDTD: Warning: Coordinate system found in the CSX file is not a cylindrical. Forcing to cylindrical coordinate system!" << endl;
|
2011-01-31 11:00:00 +00:00
|
|
|
m_CSX->SetCoordInputType(CYLINDRICAL); //tell CSX to use cylinder-coords
|
2010-11-02 15:32:00 +00:00
|
|
|
}
|
|
|
|
|
2010-07-08 09:28:11 +00:00
|
|
|
if (m_debugCSX)
|
2011-03-28 12:27:28 +00:00
|
|
|
m_CSX->Write2XML("debugCSX.xml");
|
2010-07-08 09:28:11 +00:00
|
|
|
|
2010-03-11 09:56:19 +00:00
|
|
|
//*************** setup operator ************//
|
2015-12-13 21:33:26 +00:00
|
|
|
if (SetupOperator()==false)
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
{
|
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
2011-03-21 14:10:09 +00:00
|
|
|
return 2;
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
}
|
2010-04-09 13:51:37 +00:00
|
|
|
|
2013-12-28 20:02:49 +00:00
|
|
|
// default material averaging is quarter cell averaging
|
|
|
|
FDTD_Op->SetQuarterCellMaterialAvg();
|
|
|
|
|
2013-08-23 15:29:10 +00:00
|
|
|
if (m_CellConstantMaterial)
|
2013-02-18 09:38:55 +00:00
|
|
|
{
|
|
|
|
FDTD_Op->SetCellConstantMaterial();
|
|
|
|
if (g_settings.GetVerboseLevel()>0)
|
2015-05-04 18:47:19 +00:00
|
|
|
cout << "Enabling constant cell material assumption." << endl;
|
2013-02-18 09:38:55 +00:00
|
|
|
}
|
|
|
|
|
2015-12-19 14:03:30 +00:00
|
|
|
if (m_Exc==NULL)
|
|
|
|
{
|
|
|
|
cerr << "openEMS::SetupFDTD: Error, excitation is not defined! Abort!" << endl;
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
2015-12-19 14:03:30 +00:00
|
|
|
return 3;
|
|
|
|
}
|
|
|
|
|
2012-07-17 11:23:00 +00:00
|
|
|
FDTD_Op->SetExcitationSignal(m_Exc);
|
|
|
|
FDTD_Op->AddExtension(new Operator_Ext_Excitation(FDTD_Op));
|
2012-07-25 08:50:56 +00:00
|
|
|
if (!CylinderCoords)
|
2013-04-12 14:07:39 +00:00
|
|
|
FDTD_Op->AddExtension(new Operator_Ext_TFSF(FDTD_Op));
|
2012-07-17 11:23:00 +00:00
|
|
|
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
if (FDTD_Op->SetGeometryCSX(m_CSX)==false)
|
|
|
|
{
|
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
|
|
|
return(2);
|
|
|
|
}
|
2010-03-11 09:56:19 +00:00
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
SetupBoundaryConditions();
|
2010-04-27 21:06:42 +00:00
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
FDTD_Op->SetTimeStepMethod(m_TS_method);
|
2013-03-27 10:58:24 +00:00
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
if (m_TS>0)
|
|
|
|
FDTD_Op->SetTimestep(m_TS);
|
|
|
|
if (m_TS_fac<1)
|
|
|
|
FDTD_Op->SetTimestepFactor(m_TS_fac);
|
2010-09-03 09:36:59 +00:00
|
|
|
|
2015-05-04 18:47:19 +00:00
|
|
|
// Is a steady state detection requested
|
|
|
|
Operator_Ext_SteadyState* Op_Ext_SSD = NULL;
|
|
|
|
if (m_Exc->GetSignalPeriod()>0)
|
|
|
|
{
|
|
|
|
cout << "Create a steady state detection using a period of " << m_Exc->GetSignalPeriod() << " s" << endl;
|
|
|
|
Op_Ext_SSD = new Operator_Ext_SteadyState(FDTD_Op, m_Exc->GetSignalPeriod());
|
|
|
|
unsigned int pos[3];
|
|
|
|
for (int p=0;p<3;++p)
|
|
|
|
pos[p] = FDTD_Op->GetNumberOfLines(p)/2;
|
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 0);
|
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 1);
|
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 2);
|
|
|
|
|
|
|
|
for (int n=0;n<3;++n)
|
|
|
|
{
|
|
|
|
for (int p=0;p<3;++p)
|
|
|
|
pos[p] = FDTD_Op->GetNumberOfLines(p)/2;
|
2015-09-29 17:57:21 +00:00
|
|
|
|
|
|
|
pos[n] *= 1/4;
|
2015-05-04 18:47:19 +00:00
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 0);
|
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 1);
|
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 2);
|
|
|
|
|
2015-09-29 17:57:21 +00:00
|
|
|
pos[n] *= 3/4;
|
2015-05-04 18:47:19 +00:00
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 0);
|
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 1);
|
|
|
|
Op_Ext_SSD->Add_E_Probe(pos, 2);
|
|
|
|
}
|
|
|
|
FDTD_Op->AddExtension(Op_Ext_SSD);
|
|
|
|
}
|
|
|
|
|
2013-03-19 13:02:06 +00:00
|
|
|
if ((m_CSX->GetQtyPropertyType(CSProperties::LORENTZMATERIAL)>0) || (m_CSX->GetQtyPropertyType(CSProperties::DEBYEMATERIAL)>0))
|
2012-04-27 14:47:27 +00:00
|
|
|
FDTD_Op->AddExtension(new Operator_Ext_LorentzMaterial(FDTD_Op));
|
2012-05-08 11:58:20 +00:00
|
|
|
if (m_CSX->GetQtyPropertyType(CSProperties::CONDUCTINGSHEET)>0)
|
2015-12-13 21:33:26 +00:00
|
|
|
FDTD_Op->AddExtension(new Operator_Ext_ConductingSheet(FDTD_Op, m_Exc->GetMaxFreq()));
|
2023-11-18 11:23:15 +00:00
|
|
|
if (m_CSX->GetQtyPropertyType(CSProperties::LUMPED_ELEMENT)>0)
|
|
|
|
FDTD_Op->AddExtension(new Operator_Ext_LumpedRLC(FDTD_Op));
|
|
|
|
|
2012-04-27 14:47:27 +00:00
|
|
|
|
2011-01-10 07:27:50 +00:00
|
|
|
//check all properties to request material storage during operator creation...
|
2011-01-31 11:00:00 +00:00
|
|
|
SetupMaterialStorages();
|
2011-01-07 15:12:07 +00:00
|
|
|
|
2011-01-10 07:27:50 +00:00
|
|
|
/******************* create the EC-FDTD operator *****************************/
|
2010-10-27 09:17:58 +00:00
|
|
|
Operator::DebugFlags debugFlags = Operator::None;
|
|
|
|
if (DebugMat)
|
|
|
|
debugFlags |= Operator::debugMaterial;
|
|
|
|
if (DebugOp)
|
|
|
|
debugFlags |= Operator::debugOperator;
|
|
|
|
if (m_debugPEC)
|
|
|
|
debugFlags |= Operator::debugPEC;
|
2011-01-10 07:27:50 +00:00
|
|
|
|
2010-10-27 09:17:58 +00:00
|
|
|
FDTD_Op->CalcECOperator( debugFlags );
|
2011-01-10 07:27:50 +00:00
|
|
|
/*******************************************************************************/
|
2010-12-06 12:04:37 +00:00
|
|
|
|
2011-01-10 07:27:50 +00:00
|
|
|
//reset flags for material storage, if no dump-box resets it to true, it will be cleaned up...
|
|
|
|
FDTD_Op->SetMaterialStoreFlags(0,false);
|
2011-01-07 15:12:07 +00:00
|
|
|
FDTD_Op->SetMaterialStoreFlags(1,false);
|
2011-01-10 07:27:50 +00:00
|
|
|
FDTD_Op->SetMaterialStoreFlags(2,false);
|
|
|
|
FDTD_Op->SetMaterialStoreFlags(3,false);
|
2011-01-07 15:12:07 +00:00
|
|
|
|
2015-12-13 21:33:26 +00:00
|
|
|
unsigned int maxTime_TS = (unsigned int)(m_maxTime/FDTD_Op->GetTimestep());
|
|
|
|
if ((m_maxTime>0) && (maxTime_TS<NrTS))
|
2010-05-21 14:55:04 +00:00
|
|
|
NrTS = maxTime_TS;
|
2010-03-11 09:56:19 +00:00
|
|
|
|
2015-05-04 18:47:19 +00:00
|
|
|
if (!m_Exc->buildExcitationSignal(NrTS))
|
2010-05-03 16:33:14 +00:00
|
|
|
exit(2);
|
2012-07-17 11:23:00 +00:00
|
|
|
m_Exc->DumpVoltageExcite("et");
|
|
|
|
m_Exc->DumpCurrentExcite("ht");
|
2011-03-08 12:35:38 +00:00
|
|
|
|
2010-06-05 22:47:56 +00:00
|
|
|
timeval OpDoneTime;
|
|
|
|
gettimeofday(&OpDoneTime,NULL);
|
2010-03-11 09:56:19 +00:00
|
|
|
|
2011-11-08 10:49:14 +00:00
|
|
|
if (g_settings.GetVerboseLevel()>0)
|
|
|
|
{
|
|
|
|
FDTD_Op->ShowStat();
|
|
|
|
FDTD_Op->ShowExtStat();
|
|
|
|
cout << "Creation time for operator: " << CalcDiffTime(OpDoneTime,startTime) << " s" << endl;
|
|
|
|
}
|
2011-11-14 13:12:36 +00:00
|
|
|
cout << "FDTD simulation size: " << FDTD_Op->GetNumberOfLines(0) << "x" << FDTD_Op->GetNumberOfLines(1) << "x" << FDTD_Op->GetNumberOfLines(2) << " --> " << FDTD_Op->GetNumberCells() << " FDTD cells " << endl;
|
2012-07-17 11:23:00 +00:00
|
|
|
cout << "FDTD timestep is: " <<FDTD_Op->GetTimestep() << " s; Nyquist rate: " << m_Exc->GetNyquistNum() << " timesteps @" << CalcNyquistFrequency(m_Exc->GetNyquistNum(),FDTD_Op->GetTimestep()) << " Hz" << endl;
|
|
|
|
if (m_Exc->GetNyquistNum()>1000)
|
2011-11-08 10:49:14 +00:00
|
|
|
cerr << "openEMS::SetupFDTD: Warning, the timestep seems to be very small --> long simulation. Check your mesh!?" << endl;
|
2010-03-11 09:56:19 +00:00
|
|
|
|
2015-09-03 20:53:31 +00:00
|
|
|
if (m_Exc->GetSignalPeriod()==0)
|
|
|
|
{
|
|
|
|
cout << "Excitation signal length is: " << m_Exc->GetLength() << " timesteps (" << m_Exc->GetLength()*FDTD_Op->GetTimestep() << "s)" << endl;
|
|
|
|
cout << "Max. number of timesteps: " << NrTS << " ( --> " << (double)NrTS/(double)(m_Exc->GetLength()) << " * Excitation signal length)" << endl;
|
|
|
|
if ( ((double)NrTS/(double)m_Exc->GetLength() < 3) && (m_Exc->GetExciteType()==0))
|
|
|
|
cerr << "openEMS::SetupFDTD: Warning, max. number of timesteps is smaller than three times the excitation. " << endl << \
|
|
|
|
"\tYou may want to choose a higher number of max. timesteps... " << endl;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
int p = int(m_Exc->GetSignalPeriod()/FDTD_Op->GetTimestep());
|
|
|
|
cout << "Excitation signal period is: " << p << " timesteps (" << m_Exc->GetSignalPeriod() << "s)" << endl;
|
|
|
|
cout << "Max. number of timesteps: " << NrTS << " ( --> " << (double)NrTS/(double)(m_Exc->GetLength()) << " * Excitation signal period)" << endl;
|
|
|
|
if (NrTS/p < 3)
|
|
|
|
cerr << "openEMS::SetupFDTD: Warning, max. number of timesteps is smaller than three times the excitation signal period. " << endl << \
|
|
|
|
"\tYou may want to choose a higher number of max. timesteps... " << endl;
|
|
|
|
}
|
2012-01-17 09:22:16 +00:00
|
|
|
|
2010-12-06 12:04:37 +00:00
|
|
|
if (m_no_simulation)
|
|
|
|
{
|
2010-10-19 07:05:51 +00:00
|
|
|
// simulation was disabled (to generate debug output only)
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
2010-10-19 07:05:51 +00:00
|
|
|
return 1;
|
|
|
|
}
|
|
|
|
|
2010-03-11 09:56:19 +00:00
|
|
|
//create FDTD engine
|
2010-05-20 20:02:06 +00:00
|
|
|
FDTD_Eng = FDTD_Op->CreateEngine();
|
2010-03-26 11:57:52 +00:00
|
|
|
|
2015-05-04 18:47:19 +00:00
|
|
|
if (Op_Ext_SSD)
|
|
|
|
{
|
|
|
|
Eng_Ext_SSD = dynamic_cast<Engine_Ext_SteadyState*>(Op_Ext_SSD->GetEngineExtention());
|
|
|
|
Eng_Ext_SSD->SetEngineInterface(this->NewEngineInterface());
|
|
|
|
}
|
|
|
|
|
2011-01-10 07:27:50 +00:00
|
|
|
//setup all processing classes
|
2011-01-31 11:00:00 +00:00
|
|
|
if (SetupProcessing()==false)
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
{
|
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
2011-01-10 07:27:50 +00:00
|
|
|
return 2;
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
}
|
2010-04-19 14:09:41 +00:00
|
|
|
|
2011-01-10 07:27:50 +00:00
|
|
|
// Cleanup all unused material storages...
|
2011-01-07 15:12:07 +00:00
|
|
|
FDTD_Op->CleanupMaterialStorage();
|
|
|
|
|
2011-01-10 07:27:50 +00:00
|
|
|
//check and warn for unused properties and primitives
|
2011-01-31 11:00:00 +00:00
|
|
|
m_CSX->WarnUnusedPrimitves(cerr);
|
2010-05-29 15:40:18 +00:00
|
|
|
|
2010-04-19 14:09:41 +00:00
|
|
|
// dump all boxes (voltage, current, fields, ...)
|
|
|
|
if (m_debugBox)
|
|
|
|
{
|
|
|
|
PA->DumpBoxes2File("box_dump_");
|
|
|
|
}
|
|
|
|
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
2010-03-11 09:56:19 +00:00
|
|
|
return 0;
|
|
|
|
}
|
|
|
|
|
2010-06-21 14:14:41 +00:00
|
|
|
string FormatTime(int sec)
|
|
|
|
{
|
|
|
|
stringstream ss;
|
|
|
|
if (sec<60)
|
|
|
|
{
|
2010-06-25 13:22:01 +00:00
|
|
|
ss << setw(9) << sec << "s";
|
2010-06-21 14:14:41 +00:00
|
|
|
return ss.str();
|
|
|
|
}
|
|
|
|
if (sec<3600)
|
|
|
|
{
|
2010-06-25 13:22:01 +00:00
|
|
|
ss << setw(6) << sec/60 << "m" << setw(2) << setfill('0') << sec%60 << "s";
|
2010-06-21 14:14:41 +00:00
|
|
|
return ss.str();
|
|
|
|
}
|
|
|
|
ss << setw(3) << sec/3600 << "h" << setw(2) << setfill('0') << (sec%3600)/60 << "m" << setw(2) << setfill('0') << sec%60 << "s";
|
|
|
|
return ss.str();
|
|
|
|
}
|
|
|
|
|
2010-10-18 11:26:25 +00:00
|
|
|
bool openEMS::CheckAbortCond()
|
|
|
|
{
|
|
|
|
if (m_Abort) //abort was set externally
|
|
|
|
return true;
|
|
|
|
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
//check whether SIGINT is received
|
|
|
|
if (Signal::ReceivedSIGINT())
|
|
|
|
{
|
|
|
|
cerr << "openEMS::CheckAbortCond(): Received SIGINT, aborting simulation gracefully..." << endl;
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
2010-10-18 11:26:25 +00:00
|
|
|
//check whether the file "ABORT" exist in current working directory
|
|
|
|
ifstream ifile("ABORT");
|
|
|
|
if (ifile)
|
|
|
|
{
|
|
|
|
ifile.close();
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
cerr << "openEMS::CheckAbortCond(): Found file \"ABORT\", aborting simulation gracefully..." << endl;
|
2010-10-18 11:26:25 +00:00
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
|
2010-03-11 09:56:19 +00:00
|
|
|
void openEMS::RunFDTD()
|
|
|
|
{
|
2010-03-29 20:11:24 +00:00
|
|
|
cout << "Running FDTD engine... this may take a while... grab a cup of coffee?!?" << endl;
|
2010-03-25 14:08:54 +00:00
|
|
|
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_EXIT_GRACEFUL);
|
|
|
|
|
2010-04-03 15:36:50 +00:00
|
|
|
//special handling of a field processing, needed to realize the end criteria...
|
2011-04-13 10:16:54 +00:00
|
|
|
ProcessFields* ProcField = new ProcessFields(NewEngineInterface());
|
2010-04-03 15:36:50 +00:00
|
|
|
PA->AddProcessing(ProcField);
|
2010-03-15 15:59:37 +00:00
|
|
|
double maxE=0,currE=0;
|
2010-04-03 15:36:50 +00:00
|
|
|
|
2011-02-16 09:11:39 +00:00
|
|
|
//init processings
|
|
|
|
PA->InitAll();
|
|
|
|
|
2010-04-03 15:36:50 +00:00
|
|
|
//add all timesteps to end-crit field processing with max excite amplitude
|
2012-07-16 15:15:10 +00:00
|
|
|
unsigned int maxExcite = FDTD_Op->GetExcitationSignal()->GetMaxExcitationTimestep();
|
|
|
|
// for (unsigned int n=0; n<FDTD_Op->Exc->Volt_Count; ++n)
|
|
|
|
// ProcField->AddStep(FDTD_Op->Exc->Volt_delay[n]+maxExcite);
|
|
|
|
ProcField->AddStep(maxExcite);
|
2010-04-03 15:36:50 +00:00
|
|
|
|
2010-03-15 15:59:37 +00:00
|
|
|
double change=1;
|
|
|
|
int prevTS=0,currTS=0;
|
2012-09-25 09:24:25 +00:00
|
|
|
double numCells = FDTD_Op->GetNumberCells();
|
|
|
|
double speed = 0;
|
2010-03-15 15:59:37 +00:00
|
|
|
double t_diff;
|
2012-09-25 09:24:25 +00:00
|
|
|
double t_run;
|
2010-03-27 14:54:44 +00:00
|
|
|
|
|
|
|
timeval currTime;
|
|
|
|
gettimeofday(&currTime,NULL);
|
|
|
|
timeval startTime = currTime;
|
|
|
|
timeval prevTime= currTime;
|
|
|
|
|
2012-09-25 09:24:25 +00:00
|
|
|
if (m_DumpStats)
|
|
|
|
InitRunStatistics(__OPENEMS_RUN_STAT_FILE__);
|
2010-03-11 09:56:19 +00:00
|
|
|
//*************** simulate ************//
|
2010-03-27 14:54:44 +00:00
|
|
|
|
2010-12-17 14:13:43 +00:00
|
|
|
PA->PreProcess();
|
2010-03-11 09:56:19 +00:00
|
|
|
int step=PA->Process();
|
2010-03-29 08:12:38 +00:00
|
|
|
if ((step<0) || (step>(int)NrTS)) step=NrTS;
|
2010-10-18 11:26:25 +00:00
|
|
|
while ((FDTD_Eng->GetNumberOfTimesteps()<NrTS) && (change>endCrit) && !CheckAbortCond())
|
2010-03-11 09:56:19 +00:00
|
|
|
{
|
|
|
|
FDTD_Eng->IterateTS(step);
|
|
|
|
step=PA->Process();
|
2010-04-03 15:36:50 +00:00
|
|
|
|
2015-05-04 18:47:19 +00:00
|
|
|
if ((Eng_Ext_SSD==NULL) && ProcField->CheckTimestep())
|
2010-04-03 15:36:50 +00:00
|
|
|
{
|
2011-11-07 11:07:55 +00:00
|
|
|
currE = ProcField->CalcTotalEnergyEstimate();
|
2010-04-03 15:36:50 +00:00
|
|
|
if (currE>maxE)
|
|
|
|
maxE=currE;
|
|
|
|
}
|
|
|
|
|
2010-03-12 19:39:04 +00:00
|
|
|
// cout << " do " << step << " steps; current: " << eng.GetNumberOfTimesteps() << endl;
|
2010-03-15 15:59:37 +00:00
|
|
|
currTS = FDTD_Eng->GetNumberOfTimesteps();
|
2010-03-29 08:12:38 +00:00
|
|
|
if ((step<0) || (step>(int)(NrTS - currTS))) step=NrTS - currTS;
|
2010-03-15 15:59:37 +00:00
|
|
|
|
2010-03-25 14:08:54 +00:00
|
|
|
gettimeofday(&currTime,NULL);
|
|
|
|
|
|
|
|
t_diff = CalcDiffTime(currTime,prevTime);
|
2012-09-25 09:24:25 +00:00
|
|
|
|
2010-03-15 15:59:37 +00:00
|
|
|
if (t_diff>4)
|
|
|
|
{
|
2012-09-25 09:24:25 +00:00
|
|
|
t_run = CalcDiffTime(currTime,startTime);
|
|
|
|
speed = numCells*(currTS-prevTS)/t_diff;
|
|
|
|
cout << "[@" << FormatTime(t_run) << "] Timestep: " << setw(12) << currTS ;
|
|
|
|
cout << " || Speed: " << setw(6) << setprecision(1) << std::fixed << speed*1e-6 << " MC/s (" << setw(4) << setprecision(3) << std::scientific << t_diff/(currTS-prevTS) << " s/TS)" ;
|
2015-05-04 18:47:19 +00:00
|
|
|
if (Eng_Ext_SSD==NULL)
|
|
|
|
{
|
|
|
|
currE = ProcField->CalcTotalEnergyEstimate();
|
|
|
|
if (currE>maxE)
|
|
|
|
maxE=currE;
|
|
|
|
if (maxE)
|
|
|
|
change = currE/maxE;
|
|
|
|
cout << " || Energy: ~" << setw(6) << setprecision(2) << std::scientific << currE << " (-" << setw(5) << setprecision(2) << std::fixed << fabs(10.0*log10(change)) << "dB)" << endl;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
change = Eng_Ext_SSD->GetLastDiff();
|
|
|
|
cout << " || SteadyState: " << setw(6) << setprecision(2) << std::fixed << 10.0*log10(change) << " dB" << endl;
|
|
|
|
}
|
2010-03-15 15:59:37 +00:00
|
|
|
prevTime=currTime;
|
|
|
|
prevTS=currTS;
|
2010-06-28 16:05:03 +00:00
|
|
|
|
|
|
|
PA->FlushNext();
|
2012-09-25 09:24:25 +00:00
|
|
|
|
|
|
|
if (m_DumpStats)
|
|
|
|
DumpRunStatistics(__OPENEMS_RUN_STAT_FILE__, t_run, currTS, speed, currE);
|
2023-01-06 19:01:07 +00:00
|
|
|
FDTD_Eng->NextInterval(speed);
|
2010-03-15 15:59:37 +00:00
|
|
|
}
|
2010-03-11 09:56:19 +00:00
|
|
|
}
|
2012-07-16 15:15:10 +00:00
|
|
|
if ((change>endCrit) && (FDTD_Op->GetExcitationSignal()->GetExciteType()==0))
|
2012-01-17 09:22:16 +00:00
|
|
|
cerr << "RunFDTD: Warning: Max. number of timesteps was reached before the end-criteria of -" << fabs(10.0*log10(endCrit)) << "dB was reached... " << endl << \
|
|
|
|
"\tYou may want to choose a higher number of max. timesteps... " << endl;
|
|
|
|
|
2010-03-25 14:08:54 +00:00
|
|
|
gettimeofday(&currTime,NULL);
|
|
|
|
t_diff = CalcDiffTime(currTime,startTime);
|
2010-03-11 09:56:19 +00:00
|
|
|
|
2010-03-12 19:39:04 +00:00
|
|
|
cout << "Time for " << FDTD_Eng->GetNumberOfTimesteps() << " iterations with " << FDTD_Op->GetNumberCells() << " cells : " << t_diff << " sec" << endl;
|
2012-09-25 09:24:25 +00:00
|
|
|
cout << "Speed: " << numCells*(double)FDTD_Eng->GetNumberOfTimesteps()/t_diff*1e-6 << " MCells/s " << endl;
|
2012-09-24 13:16:47 +00:00
|
|
|
|
2012-09-25 09:24:25 +00:00
|
|
|
if (m_DumpStats)
|
|
|
|
DumpStatistics(__OPENEMS_STAT_FILE__, t_diff);
|
2012-11-06 11:46:54 +00:00
|
|
|
|
|
|
|
//*************** postproc ************//
|
|
|
|
PA->PostProcess();
|
Handle SIGINT for openEMS and Python, with graceful exit support.
Currently, openEMS doesn't have any special code to handle SIGINT (which
is raised by pressing Control-C). By default, the program is terminated
without saving data. This worked okay in the past, but now its
limitations are becoming obvious.
1. When openEMS is used as a Python module, Control-C stops working
because SIGINT is now managed by Python in order to generate
KeyboardInterrupt exceptions, normally this isn't a problem, but if
we are running an external C++ (Cython) function such as openEMS, the
Python interpreter mainloop has no control until we return. As a
result, SIGINT is received but never handled. In Cython, programs are
expected to call PyErr_CheckSignals() in its blocking loop periodically
to temporally transfer control back to Python to handle signals. But
this introduces a dependency of Cython in the FDTD mainloop.
2. During a simulation, it's not possible to abort it gracefully by
pressing Control-C, this is a limitation of openEMS itself, it's
always a force exit. Currently the only supported method for graceful
exit is creating a file called "ABORT" in the simulation directory.
If we already need to implement a signal handler, adding a graceful
exit at the same time would be a good idea.
This commit installs SIGINT handlers during SetupFDTD() and RunFDTD().
1. In RunFDTD(), if SIGINT is received once, a status flag is set, which
is then checked in CheckAbortCond(), allowing a graceful exit with the
same effect of an "ABORT" file. If SIGINT is received twice, openEMS
force exit without saving data (just like the old default behavior).
2. In SetupFDTD(), if SIGINT is received, openEMS immediately force
exit without saving data, identical to the old behavior. In a huge
simulation, initializing and compressing operators may have a long
time. so we want an early exit before RunFDTD().
3. Before RunFDTD() and SetupFDTD() return, the original signal handler
for SIGINT is restored. This is important since when we're acting as
a shared library. When a program (such as the Python interpreter) calls
us, changing the SIGINT handler unilaterally may overwrite the original
handler and affect the functionality of the original program. For
example, Python would never be able to raise KeyboardInterrupt again.
Thus, we save the original handler and restore it later.
Signed-off-by: Yifeng Li <tomli@tomli.me>
2023-05-07 09:59:13 +00:00
|
|
|
|
|
|
|
Signal::SetupHandlerForSIGINT(SIGNAL_ORIGINAL);
|
2012-09-24 13:16:47 +00:00
|
|
|
}
|
|
|
|
|
|
|
|
bool openEMS::DumpStatistics(const string& filename, double time)
|
|
|
|
{
|
|
|
|
ofstream stat_file;
|
|
|
|
stat_file.open(filename.c_str());
|
|
|
|
|
|
|
|
if (!stat_file.is_open())
|
|
|
|
{
|
|
|
|
cerr << "openEMS::DumpStatistics: Error, opening file failed..." << endl;
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
stat_file << std::setprecision( 16 );
|
|
|
|
stat_file << FDTD_Op->GetNumberCells() << "\t% number of cells" << endl;
|
|
|
|
stat_file << FDTD_Op->GetTimestep() << "\t% timestep (s)" << endl;
|
|
|
|
stat_file << FDTD_Eng->GetNumberOfTimesteps() << "\t% number of iterations" << endl;
|
2022-03-11 21:17:27 +00:00
|
|
|
stat_file << FDTD_Eng->GetNumberOfTimesteps()*FDTD_Op->GetTimestep() << "\t% total numerical time (s)" << endl;
|
2012-09-24 13:16:47 +00:00
|
|
|
stat_file << time << "\t% simulation time (s)" << endl;
|
|
|
|
stat_file << (double)FDTD_Op->GetNumberCells()*(double)FDTD_Eng->GetNumberOfTimesteps()/time << "\t% speed (cells/s)" << endl;
|
|
|
|
|
|
|
|
stat_file.close();
|
|
|
|
return true;
|
2010-03-11 09:56:19 +00:00
|
|
|
}
|
2012-09-25 09:24:25 +00:00
|
|
|
|
|
|
|
bool openEMS::InitRunStatistics(const string& filename)
|
|
|
|
{
|
|
|
|
ofstream stat_file;
|
|
|
|
stat_file.open(filename.c_str(), ios_base::out);
|
|
|
|
|
|
|
|
if (!stat_file.is_open())
|
|
|
|
{
|
|
|
|
cerr << "openEMS::InitRunStatistics: Error, opening file failed..." << endl;
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
stat_file << "%time\ttimestep\tspeed\tenergy" << endl;
|
|
|
|
stat_file.close();
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
|
|
|
bool openEMS::DumpRunStatistics(const string& filename, double time, unsigned int ts, double speed, double energy)
|
|
|
|
{
|
|
|
|
ofstream stat_file;
|
|
|
|
stat_file.open(filename.c_str(), ios_base::app);
|
|
|
|
|
|
|
|
if (!stat_file.is_open())
|
|
|
|
{
|
|
|
|
cerr << "openEMS::DumpRunStatistics: Error, opening file failed..." << endl;
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
stat_file << time << "\t" << ts << "\t" << speed << "\t" << energy << endl;
|
|
|
|
stat_file.close();
|
|
|
|
return true;
|
|
|
|
}
|