/* * 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 "operator.h" #include "engine.h" #include "operator_extension.h" #include "processfields.h" #include "tools/array_ops.h" #include "fparser.hh" Operator* Operator::New() { cout << "Create FDTD operator" << endl; Operator* op = new Operator(); op->Init(); return op; } Operator::Operator() { m_MeshType = ProcessFields::CARTESIAN_MESH; Exc = 0; dT = 0; m_InvaildTimestep = false; } Operator::~Operator() { for (size_t n=0;n2)) return 0.0; if (pos[n]>=numLines[n]) return 0.0; if (dualMesh==false) return fabs(MainOp->GetIndexDelta(n,pos[n]))*gridDelta; else return fabs(MainOp->GetIndexWidth(n,pos[n]))*gridDelta; } double Operator::GetDiscLine(int n, unsigned int pos, bool dualMesh) const { if ((n<0) || (n>2)) return 0.0; if (pos>=numLines[n]) return 0.0; if (dualMesh==false) return discLines[n][pos]; else return (discLines[n][pos] + 0.5*fabs(MainOp->GetIndexDelta(n,pos))); } 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); } bool Operator::SnapToMesh(double* dcoord, unsigned int* uicoord, bool lower, bool* inside) { bool ok=true; unsigned int numLines[3]; for (int n=0;n<3;++n) { numLines[n] = GetNumberOfLines(n); if (inside) //set defaults inside[n] = true; uicoord[n]=0; if (dcoord[n]discLines[n][numLines[n]-1]) { ok=false; uicoord[n]=numLines[n]-1; if (lower) uicoord[n]=numLines[n]-2; if (inside) inside[n] = false; } else if (dcoord[n]==discLines[n][numLines[n]-1]) { uicoord[n]=numLines[n]-1; if (lower) uicoord[n]=numLines[n]-2; } else for (unsigned int i=1;i=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]Volt_Count << "\t (" << Exc->Volt_Count_Dir[0] << ", " << Exc->Volt_Count_Dir[1] << ", " << Exc->Volt_Count_Dir[2] << ")" << endl; cout << "Current excitations\t: " << Exc->Curr_Count << "\t (" << Exc->Curr_Count_Dir[0] << ", " << Exc->Curr_Count_Dir[1] << ", " << Exc->Curr_Count_Dir[2] << ")" << endl; cout << "-----------------------------------" << endl; cout << "Number of PEC edges\t: " << m_Nr_PEC[0]+m_Nr_PEC[1]+m_Nr_PEC[2] << endl; cout << "in " << GetDirName(0) << " direction\t\t: " << m_Nr_PEC[0] << endl; cout << "in " << GetDirName(1) << " direction\t\t: " << m_Nr_PEC[1] << endl; cout << "in " << GetDirName(2) << " direction\t\t: " << m_Nr_PEC[2] << endl; cout << "-----------------------------------" << endl; cout << "Timestep (s)\t\t: " << dT ; if (opt_dT) cout <<"\t(" << opt_dT << ")"; cout << endl; cout << "Timestep method name\t: " << m_Used_TS_Name << endl; cout << "Nyquist criteria (TS)\t: " << Exc->GetNyquistNum() << endl; cout << "Nyquist criteria (s)\t: " << Exc->GetNyquistNum()*dT << endl; cout << "Excitation Length (TS)\t: " << Exc->Length << endl; cout << "Excitation Length (s)\t: " << Exc->Length*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 ofstream file(filename.c_str(),ios_base::out); if (file.is_open()==false) { cerr << "Operator::DumpOperator2File: Can't open file: " << filename << endl; return; } 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]; } string names[] = {"vv", "vi", "iv" , "ii", "exc"}; FDTD_FLOAT**** array[] = {vv,vi,iv,ii,exc}; ProcessFields::DumpMultiVectorArray2VTK(file, names , array , 5, discLines, numLines, 6, "Operator dump" , (ProcessFields::MeshType)m_MeshType, discLines_scaling); Delete_N_3DArray(exc,numLines); file.close(); 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 ) { ofstream file( filename.c_str() ); if (!file.is_open()) { cerr << "Operator::DumpPEC2File: Can't open file: " << filename << endl; return; } cout << "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]( numLines); mue = Create_N_3DArray( numLines); kappa = Create_N_3DArray( numLines); sigma = Create_N_3DArray( numLines); for (pos[0]=0;pos[0]GetGrid(); for (int n=0;n<3;++n) { discLines[n] = grid->GetLines(n,discLines[n],numLines[n],true); if (n==1) if (numLines[n]<3) {cerr << "CartOperator::SetGeometryCSX: 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::SetGeometryCSX: grid delta unit must not be <=0 !!!" << endl; Reset(); return false;} else gridDelta=grid->GetDeltaUnit(); MainOp->SetGridDelta(1); MainOp->AddCellAdrOp(); return true; } 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::InitExcitation() { delete Exc; Exc = new Excitation( dT ); } 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) { GetVV(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]); GetVI(n,pos[0],pos[1],pos[2]) = (dT/EC_C[n][i])/(1+dT*EC_G[n][i]/2/EC_C[n][i]); } else { GetVV(n,pos[0],pos[1],pos[2]) = 0; GetVI(n,pos[0],pos[1],pos[2]) = 0; } if (EC_L[n][i]>0) { GetII(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]); GetIV(n,pos[0],pos[1],pos[2]) = (dT/EC_L[n][i])/(1+dT*EC_R[n][i]/2/EC_L[n][i]); } else { GetII(n,pos[0],pos[1],pos[2]) = 0; GetIV(n,pos[0],pos[1],pos[2]) = 0; } } int Operator::CalcECOperator() { Init_EC(); 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(); //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==NULL) 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]GetIndexDelta(n,pos[n]); double deltaP=MainOp->GetIndexDelta(nP,pos[nP]); double deltaPP=MainOp->GetIndexDelta(nPP,pos[nPP]); double delta_M=MainOp->GetIndexDelta(n,pos[n]-1); double deltaP_M=MainOp->GetIndexDelta(nP,pos[nP]-1); double deltaPP_M=MainOp->GetIndexDelta(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,(unsigned int*)loc_pos,true); // { // cerr << ny << " " << pos[0] << " " << pos[1] << " " << pos[2] << ": " << A_n << endl; // exit(0); // } CSProperties* prop = CSX->GetPropertyByCoordPriority(shiftCoord,CSProperties::MATERIAL); if (prop) { CSPropMaterial* mat = prop->ToMaterial(); EffMat[0] = mat->GetEpsilonWeighted(n,shiftCoord)*A_n; EffMat[1] = mat->GetKappaWeighted(n,shiftCoord)*A_n; } else { EffMat[0] = 1*A_n; EffMat[1] = 0; } 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,(unsigned int*)loc_pos,true); // cerr << A_n << endl; prop = CSX->GetPropertyByCoordPriority(shiftCoord,CSProperties::MATERIAL); if (prop) { CSPropMaterial* mat = prop->ToMaterial(); EffMat[0] += mat->GetEpsilonWeighted(n,shiftCoord)*A_n; EffMat[1] += mat->GetKappaWeighted(n,shiftCoord)*A_n; } else { EffMat[0] += 1*A_n; EffMat[1] += 0; } 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,(unsigned int*)loc_pos,true); prop = CSX->GetPropertyByCoordPriority(shiftCoord,CSProperties::MATERIAL); if (prop) { CSPropMaterial* mat = prop->ToMaterial(); EffMat[0] += mat->GetEpsilonWeighted(n,shiftCoord)*A_n; EffMat[1] += mat->GetKappaWeighted(n,shiftCoord)*A_n; } else { EffMat[0] += 1*A_n; EffMat[1] += 0; } 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,(unsigned int*)loc_pos,true); prop = CSX->GetPropertyByCoordPriority(shiftCoord,CSProperties::MATERIAL); if (prop) { CSPropMaterial* mat = prop->ToMaterial(); EffMat[0] += mat->GetEpsilonWeighted(n,shiftCoord)*A_n; EffMat[1] += mat->GetKappaWeighted(n,shiftCoord)*A_n; } else { EffMat[0] += 1*A_n; EffMat[1] += 0; } 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,(unsigned int*)loc_pos,true); prop = CSX->GetPropertyByCoordPriority(shiftCoord,CSProperties::MATERIAL); if (prop) { CSPropMaterial* mat = prop->ToMaterial(); EffMat[2] = delta_ny / mat->GetMueWeighted(n,shiftCoord); if (mat->GetSigmaWeighted(n,shiftCoord)) EffMat[3] = delta_ny / mat->GetSigmaWeighted(n,shiftCoord); else EffMat[3] = 0; } else { EffMat[2] = delta_ny; 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,(unsigned int*)loc_pos,true); prop = CSX->GetPropertyByCoordPriority(shiftCoord,CSProperties::MATERIAL); if (prop) { CSPropMaterial* mat = prop->ToMaterial(); EffMat[2] += delta_ny / mat->GetMueWeighted(n,shiftCoord); if (mat->GetSigmaWeighted(n,shiftCoord)) EffMat[3] += delta_ny/mat->GetSigmaWeighted(n,shiftCoord); else EffMat[3] = 0; } else { EffMat[2] += 1*delta_ny; 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; } 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 1 //use the new timestep-calc (1) or the old one (0) return CalcTimestep_Var3(); //the biggest one for cartesian meshes #else return CalcTimestep_Var1(); #endif } ////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 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; } } } } if (dT==0) { cerr << "Operator::CalcTimestep: Timestep is zero... this is not supposed to happen!!! exit!" << endl; exit(3); } // cerr << "Operator Timestep: " << dT << 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; } } } } if (dT==0) { cerr << "Operator::CalcTimestep: Timestep is zero... this is not supposed to happen!!! exit!" << endl; exit(3); } // cerr << "Operator Timestep: " << dT << endl; return 0; } bool Operator::CalcFieldExcitation() { if (dT==0) return false; if (Exc==0) return false; unsigned int pos[3]; double delta[3]; double amp=0; vector volt_vIndex[3]; vector volt_vExcit; vector volt_vDelay; vector volt_vDir; double volt_coord[3]; vector curr_vIndex[3]; vector curr_vExcit; vector curr_vDelay; vector curr_vDir; double curr_coord[3]; vector vec_prop = CSX->GetPropertyByType(CSProperties::ELECTRODE); if (vec_prop.size()==0) { cerr << "Operator::CalcFieldExcitation: Warning, no excitation properties found" << endl; return false; } CSPropElectrode* elec=NULL; CSProperties* prop=NULL; int priority=0; for (pos[2]=0;pos[2]GetIndexDelta(2,pos[2])); for (pos[1]=0;pos[1]GetIndexDelta(1,pos[1])); for (pos[0]=0;pos[0]GetIndexDelta(0,pos[0])); //electric field excite for (int n=0;n<3;++n) { volt_coord[0] = discLines[0][pos[0]]; volt_coord[1] = discLines[1][pos[1]]; volt_coord[2] = discLines[2][pos[2]]; volt_coord[n]+=delta[n]*0.5; for (size_t p=0;pToElectrode(); if (prop->CheckCoordInPrimitive(volt_coord,priority)==false) elec=NULL; if (elec!=NULL) { if ((elec->GetActiveDir(n)) && ( (elec->GetExcitType()==0) || (elec->GetExcitType()==1) ))//&& (pos[n]GetWeightedExcitation(n,volt_coord)*GetMeshDelta(n,pos);// delta[n]*gridDelta; if (amp!=0) { volt_vExcit.push_back(amp); volt_vDelay.push_back((unsigned int)(elec->GetDelay()/dT)); volt_vDir.push_back(n); volt_vIndex[0].push_back(pos[0]); volt_vIndex[1].push_back(pos[1]); volt_vIndex[2].push_back(pos[2]); } if (elec->GetExcitType()==1) //hard excite { GetVV(n,pos[0],pos[1],pos[2]) = 0; GetVI(n,pos[0],pos[1],pos[2]) = 0; } } } } } //magnetic field excite for (int n=0;n<3;++n) { if ((pos[0]>=numLines[0]-1) || (pos[1]>=numLines[1]-1) || (pos[2]>=numLines[2]-1)) continue; //skip the last H-Line which is outside the FDTD-domain int nP = (n+1)%3; int nPP = (n+2)%3; curr_coord[0] = discLines[0][pos[0]]; curr_coord[1] = discLines[1][pos[1]]; curr_coord[2] = discLines[2][pos[2]]; curr_coord[nP] +=delta[nP]*0.5; curr_coord[nPP] +=delta[nPP]*0.5; for (size_t p=0;pToElectrode(); if (prop->CheckCoordInPrimitive(curr_coord,priority)==false) elec=NULL; if (elec!=NULL) { if ((elec->GetActiveDir(n)) && ( (elec->GetExcitType()==2) || (elec->GetExcitType()==3) )) { amp = elec->GetWeightedExcitation(n,curr_coord)*GetMeshDelta(n,pos,true);// delta[n]*gridDelta; if (amp!=0) { curr_vExcit.push_back(amp); curr_vDelay.push_back((unsigned int)(elec->GetDelay()/dT)); curr_vDir.push_back(n); curr_vIndex[0].push_back(pos[0]); curr_vIndex[1].push_back(pos[1]); curr_vIndex[2].push_back(pos[2]); } if (elec->GetExcitType()==3) //hard excite { GetII(n,pos[0],pos[1],pos[2]) = 0; GetIV(n,pos[0],pos[1],pos[2]) = 0; } } } } } } } } //special treatment for primitives of type curve (treated as wires) see also Calc_PEC double p1[3]; double p2[3]; double deltaN=0.0; struct Grid_Path path; for (size_t p=0;pToElectrode(); for (size_t n=0;nGetQtyPrimitives();++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;tSetPos(pos[0],pos[1],pos[2]); deltaN=fabs(MainOp->GetIndexDelta(n,pos[n])); volt_coord[0] = discLines[0][pos[0]]; volt_coord[1] = discLines[1][pos[1]]; volt_coord[2] = discLines[2][pos[2]]; volt_coord[n] += 0.5*deltaN; // cerr << n << " " << coord[0] << " " << coord[1] << " " << coord[2] << endl; if (elec!=NULL) { if ((elec->GetActiveDir(n)) && (pos[n]GetExcitType()==0) || (elec->GetExcitType()==1) )) { amp = elec->GetWeightedExcitation(n,volt_coord)*deltaN*gridDelta; if (amp!=0) { volt_vExcit.push_back(amp); volt_vDelay.push_back((unsigned int)(elec->GetDelay()/dT)); volt_vDir.push_back(n); volt_vIndex[0].push_back(pos[0]); volt_vIndex[1].push_back(pos[1]); volt_vIndex[2].push_back(pos[2]); } if (elec->GetExcitType()==1) //hard excite { GetVV(n,pos[0],pos[1],pos[2]) = 0; GetVI(n,pos[0],pos[1],pos[2]) = 0; } } } } } } } } // set voltage excitations Exc->setupVoltageExcitation( volt_vIndex, volt_vExcit, volt_vDelay, volt_vDir ); // set current excitations Exc->setupCurrentExcitation( curr_vIndex, curr_vExcit, curr_vDelay, curr_vDir ); return true; } 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]; double delta; unsigned int pos[3]; for (pos[0]=startX;pos[0]<=stopX;++pos[0]) { for (pos[1]=0;pos[1]GetIndexDelta(n,pos[n]); coord[n]= discLines[n][pos[n]] + delta*0.5; CSProperties* prop = CSX->GetPropertyByCoordPriority(coord, (CSProperties::PropertyType)(CSProperties::MATERIAL | CSProperties::METAL)); if (prop) { if (prop->GetType()==CSProperties::METAL) //set to PEC { GetVV(n,pos[0],pos[1],pos[2]) = 0; GetVI(n,pos[0],pos[1],pos[2]) = 0; ++counter[n]; // cerr << "CartOperator::CalcPEC: PEC found at " << pos[0] << " ; " << pos[1] << " ; " << pos[2] << endl; } } } } } } } 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