/* * Copyright (C) 2010 Thorsten Liebig (Thorsten.Liebig@gmx.de) * * 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 . */ #include #include #include "operator.h" #include "engine.h" #include "extensions/operator_extension.h" #include "extensions/operator_ext_excitation.h" #include "Common/processfields.h" #include "tools/array_ops.h" #include "tools/vtk_file_writer.h" #include "fparser.hh" Operator* Operator::New() { cout << "Create FDTD operator" << endl; Operator* op = new Operator(); op->Init(); return op; } Operator::Operator() : Operator_Base() { Exc = 0; m_InvaildTimestep = false; m_TimeStepVar = 3; } Operator::~Operator() { for (size_t n=0; n2)) return 0.0; if (pos>=numLines[n]) return 0.0; if (dualMesh==false) return discLines[n][pos]; // return dual mesh node if (pos2)) return 0.0; if (pos>=numLines[n]) return 0.0; double delta=0; if (dualMesh==false) { if (pos0) delta = GetDiscLine(n,pos,true) - GetDiscLine(n,pos-1,true); else delta = GetDiscLine(n,1,false) - GetDiscLine(n,0,false); return delta; } } bool Operator::GetYeeCoords(int ny, unsigned int pos[3], double* coords, bool dualMesh) const { for (int n=0;n<3;++n) coords[n]=GetDiscLine(n,pos[n],dualMesh); coords[ny]=GetDiscLine(ny,pos[ny],!dualMesh); //check if position is inside the FDTD domain if (dualMesh==false) //main grid { if (pos[ny]>=numLines[ny]-1) return false; } else //dual grid { int nP = (ny+1)%3; int nPP = (ny+2)%3; if ((pos[nP]>=numLines[nP]-1) || (pos[nPP]>=numLines[nPP]-1)) return false; } return true; } double Operator::GetEdgeLength(int n, const unsigned int* pos, bool dualMesh) const { return GetDiscDelta(n,pos[n],dualMesh)*gridDelta; } double Operator::GetCellVolume(const unsigned int pos[3], bool dualMesh) const { double vol=1; for (int n=0;n<3;++n) vol*=GetEdgeLength(n,pos,dualMesh); return vol; } double Operator::GetNodeWidth(int ny, const int pos[3], bool dualMesh) const { if ( (pos[0]<0) || (pos[1]<0) || (pos[2]<0) ) return 0.0; unsigned int uiPos[]={pos[0],pos[1],pos[2]}; return GetNodeWidth(ny, uiPos, dualMesh); } double Operator::GetNodeArea(int ny, const unsigned int pos[3], bool dualMesh) const { int nyP = (ny+1)%3; int nyPP = (ny+2)%3; return GetNodeWidth(nyP,pos,dualMesh) * GetNodeWidth(nyPP,pos,dualMesh); } double Operator::GetNodeArea(int ny, const int pos[3], bool dualMesh) const { if ( (pos[0]<0) || (pos[1]<0) || (pos[2]<0) ) return 0.0; unsigned int uiPos[]={pos[0],pos[1],pos[2]}; return GetNodeArea(ny, uiPos, dualMesh); } unsigned int Operator::SnapToMeshLine(int ny, double coord, bool &inside, bool dualMesh) const { inside = false; if ((ny<0) || (ny>2)) return 0; if (coordGetDiscLine(ny,numLines-1)) return numLines-1; inside=true; if (dualMesh==false) { for (unsigned int n=0;nmax) && (l_stop[n]>max)) ) { return -2; } } SnapToMesh(l_start, uiStart, dualMesh, bStartIn); SnapToMesh(l_stop, uiStop, dualMesh, bStopIn); int iDim = 0; if (SnapMethod==0) { for (int n=0;n<3;++n) if (uiStop[n]>uiStart[n]) ++iDim; return iDim; } else if (SnapMethod==1) { for (int n=0;n<3;++n) { if (uiStop[n]>uiStart[n]) { if ((GetDiscLine( n, uiStart[n], dualMesh ) > l_start[n]) && (uiStart[n]>0)) --uiStart[n]; if ((GetDiscLine( n, uiStop[n], dualMesh ) < l_stop[n]) && (uiStop[n]uiStart[n]) ++iDim; } return iDim; } else if (SnapMethod==2) { for (int n=0;n<3;++n) { if (uiStop[n]>uiStart[n]) { if ((GetDiscLine( n, uiStart[n], dualMesh ) < l_start[n]) && (uiStart[n] l_stop[n]) && (uiStop[n]>0)) --uiStop[n]; } if (uiStop[n]>uiStart[n]) ++iDim; } return iDim; } else cerr << "Operator::SnapBox2Mesh: Unknown snapping method!" << endl; return -1; } struct Operator::Grid_Path Operator::FindPath(double start[], double stop[]) { struct Grid_Path path; // double dV[] = {stop[0]-start[0],stop[1]-start[1],stop[2]-start[2]}; unsigned int uiStart[3],uiStop[3],currPos[3]; SnapToMesh(start,uiStart); SnapToMesh(stop,uiStop); currPos[0]=uiStart[0]; currPos[1]=uiStart[1]; currPos[2]=uiStart[2]; double meshStart[] = {discLines[0][uiStart[0]], discLines[1][uiStart[1]], discLines[2][uiStart[2]]}; double meshStop[] = {discLines[0][uiStop[0]], discLines[1][uiStop[1]], discLines[2][uiStop[2]]}; bool UpDir = false; double foot=0,dist=0,minFoot=0,minDist=0; int minDir=0; unsigned int minPos[3]; double startFoot,stopFoot,currFoot; Point_Line_Distance(meshStart,start,stop,startFoot,dist); Point_Line_Distance(meshStop,start,stop,stopFoot,dist); currFoot=startFoot; minFoot=startFoot; double P[3]; while (minFoot=0) { P[n] = discLines[n][currPos[n]-1]; Point_Line_Distance(P,start,stop,foot,dist); if ((foot>currFoot) && (distcurrFoot) && (distuiStop[n]) { --currPos[n]; path.posPath[0].push_back(currPos[0]); path.posPath[1].push_back(currPos[1]); path.posPath[2].push_back(currPos[2]); path.dir.push_back(n); } else if (currPos[n]GetNyquistNum() << endl; cout << "Nyquist criteria (s)\t: " << Exc->GetNyquistNum()*dT << endl; cout << "-----------------------------------" << endl; } void Operator::ShowExtStat() const { if (m_Op_exts.size()==0) return; cout << "-----------------------------------" << endl; for (size_t n=0; nShowStat(cout); cout << "-----------------------------------" << endl; } void Operator::DumpOperator2File(string filename) { #ifdef OUTPUT_IN_DRAWINGUNITS double discLines_scaling = 1; #else double discLines_scaling = GetGridDelta(); #endif cout << "Operator: Dumping FDTD operator information to vtk file: " << filename << " ..." << flush; FDTD_FLOAT**** exc = Create_N_3DArray(numLines); if (Exc) { for (unsigned int n=0; nVolt_Count; ++n) exc[Exc->Volt_dir[n]][Exc->Volt_index[0][n]][Exc->Volt_index[1][n]][Exc->Volt_index[2][n]] = Exc->Volt_amp[n]; } FDTD_FLOAT**** vv_temp = Create_N_3DArray(numLines); FDTD_FLOAT**** vi_temp = Create_N_3DArray(numLines); FDTD_FLOAT**** iv_temp = Create_N_3DArray(numLines); FDTD_FLOAT**** ii_temp = Create_N_3DArray(numLines); unsigned int pos[3], n; for (n=0; n<3; n++) for (pos[0]=0; pos[0]SetMeshLines(discLines,numLines,discLines_scaling); vtk_Writer->SetHeader("openEMS - Operator dump"); vtk_Writer->SetNativeDump(true); vtk_Writer->AddVectorField("vv",vv_temp); Delete_N_3DArray(vv_temp,numLines); vtk_Writer->AddVectorField("vi",vi_temp); Delete_N_3DArray(vi_temp,numLines); vtk_Writer->AddVectorField("iv",iv_temp); Delete_N_3DArray(iv_temp,numLines); vtk_Writer->AddVectorField("ii",ii_temp); Delete_N_3DArray(ii_temp,numLines); vtk_Writer->AddVectorField("exc",exc); Delete_N_3DArray(exc,numLines); if (vtk_Writer->Write()==false) cerr << "Operator::DumpOperator2File: Error: Can't write file... skipping!" << endl; cout << " done!" << endl; } //! \brief dump PEC (perfect electric conductor) information (into VTK-file) //! visualization via paraview //! visualize only one component (x, y or z) void Operator::DumpPEC2File( string filename ) { cout << "Operator: Dumping PEC information to vtk file: " << filename << " ..." << flush; FDTD_FLOAT**** pec = Create_N_3DArray( numLines ); unsigned int pos[3]; #ifdef OUTPUT_IN_DRAWINGUNITS double scaling = 1.0/GetGridDelta(); #else double scaling = 1; #endif for (pos[0]=0; pos[0]SetMeshLines(discLines,numLines,scaling); vtk_Writer->SetHeader("openEMS - PEC dump"); vtk_Writer->SetNativeDump(true); vtk_Writer->AddVectorField("PEC",pec); Delete_N_3DArray(pec,numLines); if (vtk_Writer->Write()==false) cerr << "Operator::DumpPEC2File: Error: Can't write file... skipping!" << endl; cout << " done!" << endl; } void Operator::DumpMaterial2File(string filename) { #ifdef OUTPUT_IN_DRAWINGUNITS double discLines_scaling = 1; #else double discLines_scaling = GetGridDelta(); #endif cout << "Operator: Dumping material information to vtk file: " << filename << " ..." << flush; FDTD_FLOAT**** epsilon = Create_N_3DArray(numLines); FDTD_FLOAT**** mue = Create_N_3DArray(numLines); FDTD_FLOAT**** kappa = Create_N_3DArray(numLines); FDTD_FLOAT**** sigma = Create_N_3DArray(numLines); unsigned int pos[3]; for (pos[0]=0; pos[0]SetMeshLines(discLines,numLines,discLines_scaling); vtk_Writer->SetHeader("openEMS - material dump"); vtk_Writer->SetNativeDump(true); vtk_Writer->AddVectorField("epsilon",epsilon); Delete_N_3DArray(epsilon,numLines); vtk_Writer->AddVectorField("mue",mue); Delete_N_3DArray(mue,numLines); vtk_Writer->AddVectorField("kappa",kappa); Delete_N_3DArray(kappa,numLines); vtk_Writer->AddVectorField("sigma",sigma); Delete_N_3DArray(sigma,numLines); if (vtk_Writer->Write()==false) cerr << "Operator::DumpMaterial2File: Error: Can't write file... skipping!" << endl; cout << " done!" << endl; } bool Operator::SetupCSXGrid(CSRectGrid* grid) { for (int n=0; n<3; ++n) { discLines[n] = grid->GetLines(n,discLines[n],numLines[n],true); if (numLines[n]<3) { cerr << "CartOperator::SetupCSXGrid: you need at least 3 disc-lines in every direction (3D!)!!!" << endl; Reset(); return false; } } MainOp = new AdrOp(numLines[0],numLines[1],numLines[2]); MainOp->SetGrid(discLines[0],discLines[1],discLines[2]); if (grid->GetDeltaUnit()<=0) { cerr << "CartOperator::SetupCSXGrid: grid delta unit must not be <=0 !!!" << endl; Reset(); return false; } else gridDelta=grid->GetDeltaUnit(); MainOp->SetGridDelta(1); MainOp->AddCellAdrOp(); //delete the grid clone... delete grid; return true; } bool Operator::SetGeometryCSX(ContinuousStructure* geo) { if (geo==NULL) return false; CSX = geo; CSRectGrid* grid=CSX->GetGrid(); return SetupCSXGrid(CSRectGrid::Clone(grid)); } void Operator::InitOperator() { Delete_N_3DArray(vv,numLines); Delete_N_3DArray(vi,numLines); Delete_N_3DArray(iv,numLines); Delete_N_3DArray(ii,numLines); vv = Create_N_3DArray(numLines); vi = Create_N_3DArray(numLines); iv = Create_N_3DArray(numLines); ii = Create_N_3DArray(numLines); } void Operator::InitDataStorage() { if (m_StoreMaterial[0]) { if (g_settings.GetVerboseLevel()>0) cerr << "Operator::InitDataStorage(): Storing epsR material data..." << endl; Delete_N_3DArray(m_epsR,numLines); m_epsR = Create_N_3DArray(numLines); } if (m_StoreMaterial[1]) { if (g_settings.GetVerboseLevel()>0) cerr << "Operator::InitDataStorage(): Storing kappa material data..." << endl; Delete_N_3DArray(m_kappa,numLines); m_kappa = Create_N_3DArray(numLines); } if (m_StoreMaterial[2]) { if (g_settings.GetVerboseLevel()>0) cerr << "Operator::InitDataStorage(): Storing muR material data..." << endl; Delete_N_3DArray(m_mueR,numLines); m_mueR = Create_N_3DArray(numLines); } if (m_StoreMaterial[3]) { if (g_settings.GetVerboseLevel()>0) cerr << "Operator::InitDataStorage(): Storing sigma material data..." << endl; Delete_N_3DArray(m_sigma,numLines); m_sigma = Create_N_3DArray(numLines); } } void Operator::CleanupMaterialStorage() { if (!m_StoreMaterial[0] && m_epsR) { if (g_settings.GetVerboseLevel()>0) cerr << "Operator::CleanupMaterialStorage(): Delete epsR material data..." << endl; Delete_N_3DArray(m_epsR,numLines); m_epsR = NULL; } if (!m_StoreMaterial[1] && m_kappa) { if (g_settings.GetVerboseLevel()>0) cerr << "Operator::CleanupMaterialStorage(): Delete kappa material data..." << endl; Delete_N_3DArray(m_kappa,numLines); m_kappa = NULL; } if (!m_StoreMaterial[2] && m_mueR) { if (g_settings.GetVerboseLevel()>0) cerr << "Operator::CleanupMaterialStorage(): Delete mueR material data..." << endl; Delete_N_3DArray(m_mueR,numLines); m_mueR = NULL; } if (!m_StoreMaterial[3] && m_sigma) { if (g_settings.GetVerboseLevel()>0) cerr << "Operator::CleanupMaterialStorage(): Delete sigma material data..." << endl; Delete_N_3DArray(m_sigma,numLines); m_sigma = NULL; } } double Operator::GetDiscMaterial(int type, int n, const unsigned int pos[3]) const { switch (type) { case 0: if (m_epsR==0) return 0; return m_epsR[n][pos[0]][pos[1]][pos[2]]; case 1: if (m_kappa==0) return 0; return m_kappa[n][pos[0]][pos[1]][pos[2]]; case 2: if (m_mueR==0) return 0; return m_mueR[n][pos[0]][pos[1]][pos[2]]; case 3: if (m_sigma==0) return 0; return m_sigma[n][pos[0]][pos[1]][pos[2]]; } return 0; } void Operator::InitExcitation() { if (Exc!=NULL) Exc->Reset(dT); else { Exc = new Excitation( dT ); this->AddExtension(new Operator_Ext_Excitation(this,Exc)); } } void Operator::Calc_ECOperatorPos(int n, unsigned int* pos) { unsigned int i = MainOp->SetPos(pos[0],pos[1],pos[2]); if (EC_C[n][i]>0) { SetVV(n,pos[0],pos[1],pos[2], (1-dT*EC_G[n][i]/2/EC_C[n][i])/(1+dT*EC_G[n][i]/2/EC_C[n][i]) ); SetVI(n,pos[0],pos[1],pos[2], (dT/EC_C[n][i])/(1+dT*EC_G[n][i]/2/EC_C[n][i]) ); } else { SetVV(n,pos[0],pos[1],pos[2], 0 ); SetVI(n,pos[0],pos[1],pos[2], 0 ); } if (EC_L[n][i]>0) { SetII(n,pos[0],pos[1],pos[2], (1-dT*EC_R[n][i]/2/EC_L[n][i])/(1+dT*EC_R[n][i]/2/EC_L[n][i]) ); SetIV(n,pos[0],pos[1],pos[2], (dT/EC_L[n][i])/(1+dT*EC_R[n][i]/2/EC_L[n][i]) ); } else { SetII(n,pos[0],pos[1],pos[2], 0 ); SetIV(n,pos[0],pos[1],pos[2], 0 ); } } int Operator::CalcECOperator( DebugFlags debugFlags ) { Init_EC(); InitDataStorage(); if (Calc_EC()==0) return -1; m_InvaildTimestep = false; opt_dT = 0; if (dT>0) { double save_dT = dT; CalcTimestep(); opt_dT = dT; if (dTBuildExtension(); if (debugFlags & debugMaterial) DumpMaterial2File( "material_dump" ); if (debugFlags & debugOperator) DumpOperator2File( "operator_dump" ); if (debugFlags & debugPEC) DumpPEC2File( "PEC_dump" ); //cleanup for (int n=0; n<3; ++n) { delete[] EC_C[n]; EC_C[n]=NULL; delete[] EC_G[n]; EC_G[n]=NULL; delete[] EC_L[n]; EC_L[n]=NULL; delete[] EC_R[n]; EC_R[n]=NULL; } return 0; } void Operator::ApplyElectricBC(bool* dirs) { if (!dirs) return; unsigned int pos[3]; for (int n=0; n<3; ++n) { int nP = (n+1)%3; int nPP = (n+2)%3; for (pos[nP]=0; pos[nP]=(int)numLines[ny]-1) return (discLines[ny][numLines[ny]-2] - discLines[ny][numLines[ny]-1]); return (discLines[ny][pos+1] - discLines[ny][pos]); } double Operator::GetMaterial(int ny, const double* coords, int MatType, bool markAsUsed) const { CSProperties* prop = CSX->GetPropertyByCoordPriority(coords,CSProperties::MATERIAL,markAsUsed); CSPropMaterial* mat = dynamic_cast(prop); if (mat) { switch (MatType) { case 0: return mat->GetEpsilonWeighted(ny,coords); case 1: return mat->GetKappaWeighted(ny,coords); case 2: return mat->GetMueWeighted(ny,coords); case 3: return mat->GetSigmaWeighted(ny,coords); default: cerr << "Operator::GetMaterial: Error: unknown material type" << endl; return 0; } } switch (MatType) { case 0: return 1; case 1: return 0; case 2: return 1; case 3: return 0; default: cerr << "Operator::GetMaterial: Error: unknown material type" << endl; return 0; } } bool Operator::Calc_EffMatPos(int ny, const unsigned int* pos, double* EffMat) const { int n=ny; double coord[3]; double shiftCoord[3]; int nP = (n+1)%3; int nPP = (n+2)%3; coord[0] = discLines[0][pos[0]]; coord[1] = discLines[1][pos[1]]; coord[2] = discLines[2][pos[2]]; double delta=GetRawDiscDelta(n,pos[n]); double deltaP=GetRawDiscDelta(nP,pos[nP]); double deltaPP=GetRawDiscDelta(nPP,pos[nPP]); double delta_M=GetRawDiscDelta(n,pos[n]-1); double deltaP_M=GetRawDiscDelta(nP,pos[nP]-1); double deltaPP_M=GetRawDiscDelta(nPP,pos[nPP]-1); int loc_pos[3]={pos[0],pos[1],pos[2]}; double A_n; double area = 0; //******************************* epsilon,kappa averaging *****************************// //shift up-right shiftCoord[n] = coord[n]+delta*0.5; shiftCoord[nP] = coord[nP]+deltaP*0.25; shiftCoord[nPP] = coord[nPP]+deltaPP*0.25; A_n = GetNodeArea(ny,loc_pos,true); EffMat[0] = GetMaterial(n, shiftCoord, 0)*A_n; EffMat[1] = GetMaterial(n, shiftCoord, 1)*A_n; area+=A_n; //shift up-left shiftCoord[n] = coord[n]+delta*0.5; shiftCoord[nP] = coord[nP]-deltaP_M*0.25; shiftCoord[nPP] = coord[nPP]+deltaPP*0.25; --loc_pos[nP]; A_n = GetNodeArea(ny,loc_pos,true); EffMat[0] += GetMaterial(n, shiftCoord, 0)*A_n; EffMat[1] += GetMaterial(n, shiftCoord, 1)*A_n; area+=A_n; //shift down-right shiftCoord[n] = coord[n]+delta*0.5; shiftCoord[nP] = coord[nP]+deltaP*0.25; shiftCoord[nPP] = coord[nPP]-deltaPP_M*0.25; ++loc_pos[nP]; --loc_pos[nPP]; A_n = GetNodeArea(ny,loc_pos,true); EffMat[0] += GetMaterial(n, shiftCoord, 0)*A_n; EffMat[1] += GetMaterial(n, shiftCoord, 1)*A_n; area+=A_n; //shift down-left shiftCoord[n] = coord[n]+delta*0.5; shiftCoord[nP] = coord[nP]-deltaP_M*0.25; shiftCoord[nPP] = coord[nPP]-deltaPP_M*0.25; --loc_pos[nP]; A_n = GetNodeArea(ny,loc_pos,true); EffMat[0] += GetMaterial(n, shiftCoord, 0)*A_n; EffMat[1] += GetMaterial(n, shiftCoord, 1)*A_n; area+=A_n; EffMat[0]*=__EPS0__/area; EffMat[1]/=area; //******************************* mu,sigma averaging *****************************// loc_pos[0]=pos[0]; loc_pos[1]=pos[1]; loc_pos[2]=pos[2]; double length=0; //shift down shiftCoord[n] = coord[n]-delta_M*0.25; shiftCoord[nP] = coord[nP]+deltaP*0.5; shiftCoord[nPP] = coord[nPP]+deltaPP*0.5; --loc_pos[n]; double delta_ny = GetNodeWidth(n,loc_pos,true); EffMat[2] = delta_ny / GetMaterial(n, shiftCoord, 2); double sigma = GetMaterial(n, shiftCoord, 3); if (sigma) EffMat[3] = delta_ny / sigma; else EffMat[3] = 0; length=delta_ny; //shift up shiftCoord[n] = coord[n]+delta*0.25; shiftCoord[nP] = coord[nP]+deltaP*0.5; shiftCoord[nPP] = coord[nPP]+deltaPP*0.5; ++loc_pos[n]; delta_ny = GetNodeWidth(n,loc_pos,true); EffMat[2] += delta_ny / GetMaterial(n, shiftCoord, 2); sigma = GetMaterial(n, shiftCoord, 3); if (sigma) EffMat[3] += delta_ny / sigma; else EffMat[3] = 0; length+=delta_ny; EffMat[2] = length * __MUE0__ / EffMat[2]; if (EffMat[3]) EffMat[3]=length / EffMat[3]; for (int n=0; n<4; ++n) if (isnan(EffMat[n]) || isinf(EffMat[n])) { cerr << "Operator::Calc_EffMatPos: An effective material parameter is not a valid result, this should NOT have happend... exit..." << endl; exit(0); } return true; } bool Operator::Calc_LumpedElements() { vector props = CSX->GetPropertyByType(CSProperties::LUMPED_ELEMENT); for (size_t i=0;i(props.at(i)); if (PLE==NULL) return false; //sanity check: this should never happen! vector prims = PLE->GetAllPrimitives(); for (size_t bn=0;bn(prims.at(bn)); if (box) { //calculate lumped element parameter double C = PLE->GetCapacity(); if (C<=0) C = NAN; double R = PLE->GetResistance(); if (R<0) R = NAN; if ((isnan(R)) && (isnan(C))) { cerr << "Operator::Calc_LumpedElements(): Warning: Lumped Element R or C not specified! skipping. " << " ID: " << prims.at(bn)->GetID() << " @ Property: " << PLE->GetName() << endl; continue; } int ny = PLE->GetDirection(); if ((ny<0) || (ny>2)) { cerr << "Operator::Calc_LumpedElements(): Warning: Lumped Element direction is invalid! skipping. " << " ID: " << prims.at(bn)->GetID() << " @ Property: " << PLE->GetName() << endl; continue; } int nyP = (ny+1)%3; int nyPP = (ny+2)%3; unsigned int uiStart[3]; unsigned int uiStop[3]; // snap to the native coordinate system int Snap_Dimension = Operator::SnapBox2Mesh(box->GetStartCoord()->GetCoords(m_MeshType), box->GetStopCoord()->GetCoords(m_MeshType), uiStart, uiStop); if (Snap_Dimension<=0) { if (Snap_Dimension>=-1) cerr << "Operator::Calc_LumpedElements(): Warning: Lumped Element snapping failed! Dimension is: " << Snap_Dimension << " skipping. " << " ID: " << prims.at(bn)->GetID() << " @ Property: " << PLE->GetName() << endl; // Snap_Dimension == -2 means outside the simulation domain --> no special warning, but box probably marked as unused! continue; } if (uiStart[ny]==uiStop[ny]) { cerr << "Operator::Calc_LumpedElements(): Warning: Lumped Element with zero (snapped) length is invalid! skipping. " << " ID: " << prims.at(bn)->GetID() << " @ Property: " << PLE->GetName() << endl; continue; } //calculate geometric property for this lumped element unsigned int pos[3]; double unitGC=0; int ipos=0; for (pos[ny]=uiStart[ny];pos[ny]GetCaps(); double kappa = 0; double epsilon = 0; if (R>0) kappa = 1 / R / unitGC; if (C>0) { epsilon = C / unitGC; if (epsilon< __EPS0__) { cerr << "Operator::Calc_LumpedElements(): Warning: Lumped Element capacity is too small for its size! skipping. " << " ID: " << prims.at(bn)->GetID() << " @ Property: " << PLE->GetName() << endl; C = 0; } } for (pos[ny]=uiStart[ny];pos[ny]SetPos(pos[0],pos[1],pos[2]); if (C>0) EC_C[ny][ipos] = epsilon * GetEdgeArea(ny,pos)/GetEdgeLength(ny,pos); if (R>0) EC_G[ny][ipos] = kappa * GetEdgeArea(ny,pos)/GetEdgeLength(ny,pos); if (R==0) //make lumped element a PEC if resistance is zero { SetVV(ny,pos[0],pos[1],pos[2], 0 ); SetVI(ny,pos[0],pos[1],pos[2], 0 ); } else //recalculate operator inside the lumped element Calc_ECOperatorPos(ny,pos); } } } // setup metal caps if (caps) { for (pos[nyP]=uiStart[nyP];pos[nyP]<=uiStop[nyP];++pos[nyP]) { for (pos[nyPP]=uiStart[nyPP];pos[nyPP]<=uiStop[nyPP];++pos[nyPP]) { pos[ny]=uiStart[ny]; SetVV(nyP,pos[0],pos[1],pos[2], 0 ); SetVI(nyP,pos[0],pos[1],pos[2], 0 ); ++m_Nr_PEC[nyP]; SetVV(nyPP,pos[0],pos[1],pos[2], 0 ); SetVI(nyPP,pos[0],pos[1],pos[2], 0 ); ++m_Nr_PEC[nyPP]; pos[ny]=uiStop[ny]; SetVV(nyP,pos[0],pos[1],pos[2], 0 ); SetVI(nyP,pos[0],pos[1],pos[2], 0 ); ++m_Nr_PEC[nyP]; SetVV(nyPP,pos[0],pos[1],pos[2], 0 ); SetVI(nyPP,pos[0],pos[1],pos[2], 0 ); ++m_Nr_PEC[nyPP]; } } } box->SetPrimitiveUsed(true); } else cerr << "Operator::Calc_LumpedElements(): Warning: Primitves other than boxes are not supported for lumped elements! skipping " << prims.at(bn)->GetTypeName() << " ID: " << prims.at(bn)->GetID() << " @ Property: " << PLE->GetName() << endl; } } return true; } void Operator::DumpExciationSignals() { Exc->DumpVoltageExcite("et"); Exc->DumpCurrentExcite("ht"); } void Operator::Init_EC() { for (int n=0; n<3; ++n) { //init x-cell-array delete[] EC_C[n]; delete[] EC_G[n]; delete[] EC_L[n]; delete[] EC_R[n]; EC_C[n] = new double[MainOp->GetSize()]; EC_G[n] = new double[MainOp->GetSize()]; EC_L[n] = new double[MainOp->GetSize()]; EC_R[n] = new double[MainOp->GetSize()]; for (unsigned int i=0; iGetSize(); i++) //init all { EC_C[n][i]=0; EC_G[n][i]=0; EC_L[n][i]=0; EC_R[n][i]=0; } } } bool Operator::Calc_EC() { if (CSX==NULL) { cerr << "CartOperator::Calc_EC: CSX not given or invalid!!!" << endl; return false; } unsigned int ipos; unsigned int pos[3]; double inEC[4]; for (int n=0; n<3; ++n) { for (pos[2]=0; pos[2]SetPos(pos[0],pos[1],pos[2]); EC_C[n][ipos]=inEC[0]; EC_G[n][ipos]=inEC[1]; EC_L[n][ipos]=inEC[2]; EC_R[n][ipos]=inEC[3]; } } } } return true; } double Operator::CalcTimestep() { if (m_TimeStepVar==3) return CalcTimestep_Var3(); //the biggest one for cartesian meshes //variant 1 is default return CalcTimestep_Var1(); } ////Berechnung nach Andreas Rennings Dissertation 2008, Seite 66, Formel 4.52 double Operator::CalcTimestep_Var1() { m_Used_TS_Name = string("Rennings_1"); // cout << "Operator::CalcTimestep(): Using timestep algorithm by Andreas Rennings, Dissertation @ University Duisburg-Essen, 2008, pp. 66, eq. 4.52" << endl; dT=1e200; double newT; unsigned int pos[3]; unsigned int smallest_pos[3] = {0, 0, 0}; unsigned int smallest_n = 0; unsigned int ipos; unsigned int ipos_PM; unsigned int ipos_PPM; MainOp->SetReflection2Cell(); for (int n=0; n<3; ++n) { int nP = (n+1)%3; int nPP = (n+2)%3; for (pos[2]=0; pos[2]SetPos(pos[0],pos[1],pos[2]); ipos_PM = MainOp->Shift(nP,-1); MainOp->ResetShift(); ipos_PPM= MainOp->Shift(nPP,-1); MainOp->ResetShift(); newT = 2/sqrt( ( 4/EC_L[nP][ipos] + 4/EC_L[nP][ipos_PPM] + 4/EC_L[nPP][ipos] + 4/EC_L[nPP][ipos_PM]) / EC_C[n][ipos] ); if ((newT0.0)) { dT=newT; smallest_pos[0]=pos[0];smallest_pos[1]=pos[1];smallest_pos[2]=pos[2]; smallest_n = n; } } } } } if (dT==0) { cerr << "Operator::CalcTimestep: Timestep is zero... this is not supposed to happen!!! exit!" << endl; exit(3); } if (g_settings.GetVerboseLevel()>1) { cout << "Operator::CalcTimestep_Var1: Smallest timestep (" << dT << "s) found at position: " << smallest_n << " : " << smallest_pos[0] << ";" << smallest_pos[1] << ";" << smallest_pos[2] << endl; } return 0; } double min(double* val, unsigned int count) { if (count==0) return 0.0; double min = val[0]; for (unsigned int n=1; nSetReflection2Cell(); for (int n=0; n<3; ++n) { int nP = (n+1)%3; int nPP = (n+2)%3; for (pos[2]=0; pos[2]ResetShift(); ipos = MainOp->SetPos(pos[0],pos[1],pos[2]); wqp = 1/(EC_L[nPP][ipos]*EC_C[n][MainOp->GetShiftedPos(nP ,1)]) + 1/(EC_L[nPP][ipos]*EC_C[n][ipos]); wqp += 1/(EC_L[nP ][ipos]*EC_C[n][MainOp->GetShiftedPos(nPP,1)]) + 1/(EC_L[nP ][ipos]*EC_C[n][ipos]); ipos = MainOp->Shift(nP,-1); wqp += 1/(EC_L[nPP][ipos]*EC_C[n][MainOp->GetShiftedPos(nP ,1)]) + 1/(EC_L[nPP][ipos]*EC_C[n][ipos]); ipos = MainOp->Shift(nPP,-1); wqp += 1/(EC_L[nP ][ipos]*EC_C[n][MainOp->GetShiftedPos(nPP,1)]) + 1/(EC_L[nP ][ipos]*EC_C[n][ipos]); MainOp->ResetShift(); ipos = MainOp->SetPos(pos[0],pos[1],pos[2]); wt_4[0] = 1/(EC_L[nPP][ipos] *EC_C[nP ][ipos]); wt_4[1] = 1/(EC_L[nPP][MainOp->GetShiftedPos(nP ,-1)] *EC_C[nP ][ipos]); wt_4[2] = 1/(EC_L[nP ][ipos] *EC_C[nPP][ipos]); wt_4[3] = 1/(EC_L[nP ][MainOp->GetShiftedPos(nPP,-1)] *EC_C[nPP][ipos]); wt1 = wt_4[0]+wt_4[1]+wt_4[2]+wt_4[3] - 2*min(wt_4,4); MainOp->ResetShift(); ipos = MainOp->SetPos(pos[0],pos[1],pos[2]); wt_4[0] = 1/(EC_L[nPP][ipos] *EC_C[nP ][MainOp->GetShiftedPos(n,1)]); wt_4[1] = 1/(EC_L[nPP][MainOp->GetShiftedPos(nP ,-1)] *EC_C[nP ][MainOp->GetShiftedPos(n,1)]); wt_4[2] = 1/(EC_L[nP ][ipos] *EC_C[nPP][MainOp->GetShiftedPos(n,1)]); wt_4[3] = 1/(EC_L[nP ][MainOp->GetShiftedPos(nPP,-1)] *EC_C[nPP][MainOp->GetShiftedPos(n,1)]); wt2 = wt_4[0]+wt_4[1]+wt_4[2]+wt_4[3] - 2*min(wt_4,4); w_total = wqp + wt1 + wt2; newT = 2/sqrt( w_total ); if ((newT0.0)) { dT=newT; smallest_pos[0]=pos[0];smallest_pos[1]=pos[1];smallest_pos[2]=pos[2]; smallest_n = n; } } } } } if (dT==0) { cerr << "Operator::CalcTimestep: Timestep is zero... this is not supposed to happen!!! exit!" << endl; exit(3); } if (g_settings.GetVerboseLevel()>1) { cout << "Operator::CalcTimestep_Var3: Smallest timestep (" << dT << "s) found at position: " << smallest_n << " : " << smallest_pos[0] << ";" << smallest_pos[1] << ";" << smallest_pos[2] << endl; } return 0; } bool Operator::CalcPEC() { m_Nr_PEC[0]=0; m_Nr_PEC[1]=0; m_Nr_PEC[2]=0; CalcPEC_Range(0,numLines[0]-1,m_Nr_PEC); CalcPEC_Curves(); return true; } void Operator::CalcPEC_Range(unsigned int startX, unsigned int stopX, unsigned int* counter) { double coord[3]; unsigned int pos[3]; for (pos[0]=startX; pos[0]<=stopX; ++pos[0]) { for (pos[1]=0; pos[1]GetPropertyByCoordPriority(coord, (CSProperties::PropertyType)(CSProperties::MATERIAL | CSProperties::METAL), true); if (prop) { if (prop->GetType()==CSProperties::METAL) //set to PEC { SetVV(n,pos[0],pos[1],pos[2], 0 ); SetVI(n,pos[0],pos[1],pos[2], 0 ); ++counter[n]; } } } } } } } void Operator::CalcPEC_Curves() { //special treatment for primitives of type curve (treated as wires) double p1[3]; double p2[3]; struct Grid_Path path; vector vec_prop = CSX->GetPropertyByType(CSProperties::METAL); for (size_t p=0; pGetQtyPrimitives(); ++n) { CSPrimitives* prim = prop->GetPrimitive(n); CSPrimCurve* curv = prim->ToCurve(); if (curv) { for (size_t i=1; iGetNumberOfPoints(); ++i) { curv->GetPoint(i-1,p1); curv->GetPoint(i,p2); path = FindPath(p1,p2); if (path.dir.size()>0) prim->SetPrimitiveUsed(true); for (size_t t=0; t