/* * Copyright (C) 2012-2014 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 "nf2ff_calc.h" #include "../tools/array_ops.h" #include "../tools/useful.h" #include #include #include #include #include #include #include using namespace std; nf2ff_calc_thread::nf2ff_calc_thread(nf2ff_calc* nfc, unsigned int start, unsigned int stop, unsigned int threadID, nf2ff_data &data) { m_nf_calc = nfc; m_start = start; m_stop = stop; m_threadID = threadID; m_data = data; } void nf2ff_calc_thread::operator()() { m_nf_calc->m_Barrier->wait(); // start int ny = m_data.ny; int nP = (ny+1)%3; int nPP = (ny+2)%3; unsigned int* numLines = m_data.numLines; float* normDir = m_data.normDir; float **lines = m_data.lines; float* edge_length_P = m_data.edge_length_P; float* edge_length_PP = m_data.edge_length_PP; unsigned int pos[3]; unsigned int pos_t=0; unsigned int num_t=m_stop-m_start+1; complex**** Js=m_data.Js; complex**** Ms=m_data.Ms; complex**** E_field=m_data.E_field; complex**** H_field=m_data.H_field; int mesh_type = m_data.mesh_type; // calc Js and Ms (eq. 8.15a/b) pos[ny]=0; for (pos_t=0; pos_t** m_Nt=m_data.m_Nt; complex** m_Np=m_data.m_Np; complex** m_Lt=m_data.m_Lt; complex** m_Lp=m_data.m_Lp; float center[3] = {m_nf_calc->m_centerCoord[0],m_nf_calc->m_centerCoord[1],m_nf_calc->m_centerCoord[2]}; if (mesh_type==1) { center[0] = m_nf_calc->m_centerCoord[0]*cos(m_nf_calc->m_centerCoord[1]); center[1] = m_nf_calc->m_centerCoord[0]*sin(m_nf_calc->m_centerCoord[1]); } // calc local Nt,Np,Lt and Lp float area; float cosT_cosP,cosP_sinT; float cosT_sinP,sinT_sinP; float sinT,sinP; float cosP,cosT; float r_cos_psi; float k = 2*M_PI*m_nf_calc->m_freq/__C0__*sqrt(m_nf_calc->m_permittivity*m_nf_calc->m_permeability); complex exp_jkr; complex _I_(0,1); for (unsigned int tn=0;tnm_numTheta;++tn) for (unsigned int pn=0;pnm_numPhi;++pn) { sinT = sin(m_nf_calc->m_theta[tn]); sinP = sin(m_nf_calc->m_phi[pn]); cosT = cos(m_nf_calc->m_theta[tn]); cosP = cos(m_nf_calc->m_phi[pn]); cosT_cosP = cosT*cosP; cosT_sinP = cosT*sinP; cosP_sinT = cosP*sinT; sinT_sinP = sinP*sinT; for (pos_t=0; pos_tm_Barrier->wait(); //combine all thread local Nt,Np,Lt and Lp m_nf_calc->m_Barrier->wait(); //wait for termination } /***********************************************************************/ nf2ff_calc::nf2ff_calc(float freq, vector theta, vector phi, vector center) { m_freq = freq; m_permittivity = 1; m_permeability = 1; m_numTheta = theta.size(); m_theta = new float[m_numTheta]; for (size_t n=0;n >(numLines); m_E_phi = Create2DArray >(numLines); m_H_theta = Create2DArray >(numLines); m_H_phi = Create2DArray >(numLines); m_P_rad = Create2DArray(numLines); if (center.size()==3) { m_centerCoord[0]=center.at(0); m_centerCoord[1]=center.at(1); m_centerCoord[2]=center.at(2); } else if (center.size()>0) { cerr << "nf2ff_calc::nf2ff_calc: Warning: Center coordinates error, ignoring!" << endl; m_centerCoord[0]=m_centerCoord[1]=m_centerCoord[2]=0.0; } else m_centerCoord[0]=m_centerCoord[1]=m_centerCoord[2]=0.0; m_radPower = 0; m_maxDir = 0; m_radius = 1; for (int n=0;n<3;++n) { m_MirrorType[n] = MIRROR_OFF; m_MirrorPos[n] = 0.0; } m_Barrier = NULL; m_numThreads = boost::thread::hardware_concurrency(); } nf2ff_calc::~nf2ff_calc() { delete[] m_phi; m_phi = NULL; delete[] m_theta; m_theta = NULL; unsigned int numLines[2] = {m_numTheta, m_numPhi}; Delete2DArray(m_E_theta,numLines); m_E_theta = NULL; Delete2DArray(m_E_phi,numLines); m_E_phi = NULL; Delete2DArray(m_H_theta,numLines); m_H_theta = NULL; Delete2DArray(m_H_phi,numLines); m_H_phi = NULL; Delete2DArray(m_P_rad,numLines); m_P_rad = NULL; delete m_Barrier; m_Barrier = NULL; } int nf2ff_calc::GetNormalDir(unsigned int* numLines) { int ny = -1; int nP,nPP; for (int n=0;n<3;++n) { nP = (n+1)%3; nPP = (n+2)%3; if ((numLines[n]==1) && (numLines[nP]>2) && (numLines[nPP]>2)) ny=n; } return ny; } void nf2ff_calc::SetMirror(int type, int dir, float pos) { if ((dir<0) || (dir>3)) { cerr << "nf2ff_calc::SetMirror: Error, invalid direction!" << endl; return; } if ((type!=MIRROR_PEC) && (type!=MIRROR_PMC)) { cerr << "nf2ff_calc::SetMirror: Error, invalid type!" << endl; return; } m_MirrorType[dir] = type; m_MirrorPos[dir] = pos; } bool nf2ff_calc::AddMirrorPlane(int n, float **lines, unsigned int* numLines, complex**** E_field, complex**** H_field, int MeshType) { float E_factor[3] = {1,1,1}; float H_factor[3] = {1,1,1}; int nP = (n+1)%3; int nPP = (n+2)%3; // mirror in ny direction for (unsigned int i=0;iAddSinglePlane(lines, numLines, E_field, H_field, MeshType); } bool nf2ff_calc::AddPlane(float **lines, unsigned int* numLines, complex**** E_field, complex**** H_field, int MeshType) { this->AddSinglePlane(lines, numLines, E_field, H_field, MeshType); for (int n=0;n<3;++n) { int nP = (n+1)%3; int nPP = (n+2)%3; // check if a single mirror plane is on if ((m_MirrorType[n]!=MIRROR_OFF) && (m_MirrorType[nP]==MIRROR_OFF) && (m_MirrorType[nPP]==MIRROR_OFF)) { this->AddMirrorPlane(n, lines, numLines, E_field, H_field, MeshType); for (unsigned int i=0;iAddMirrorPlane(nP, lines, numLines, E_field, H_field, MeshType); this->AddMirrorPlane(nPP, lines, numLines, E_field, H_field, MeshType); this->AddMirrorPlane(nP, lines, numLines, E_field, H_field, MeshType); for (unsigned int i=0;iAddMirrorPlane(0, lines, numLines, E_field, H_field, MeshType); this->AddMirrorPlane(1, lines, numLines, E_field, H_field, MeshType); this->AddMirrorPlane(0, lines, numLines, E_field, H_field, MeshType); this->AddMirrorPlane(2, lines, numLines, E_field, H_field, MeshType); this->AddMirrorPlane(0, lines, numLines, E_field, H_field, MeshType); this->AddMirrorPlane(1, lines, numLines, E_field, H_field, MeshType); this->AddMirrorPlane(0, lines, numLines, E_field, H_field, MeshType); for (unsigned int i=0;i**** E_field, complex**** H_field, int MeshType) { //find normal direction int ny = this->GetNormalDir(numLines); if (ny<0) { cerr << "nf2ff_calc::AddPlane: Error can't determine normal direction..." << endl; return false; } int nP = (ny+1)%3; int nPP = (ny+2)%3; complex**** Js = Create_N_3DArray >(numLines); complex**** Ms = Create_N_3DArray >(numLines); float normDir[3]= {0,0,0}; if (lines[ny][0]>=m_centerCoord[ny]) normDir[ny]=1; else normDir[ny]=-1; unsigned int pos[3]; float edge_length_P[numLines[nP]]; for (unsigned int n=1;n power = 0; double area; for (pos[0]=0; pos[0] jpt = AssignJobs2Threads(numLines[nP], m_numThreads, true); m_numThreads = jpt.size(); nf2ff_data thread_data[m_numThreads]; m_Barrier = new boost::barrier(m_numThreads+1); // numThread workers + 1 controller unsigned int start=0; unsigned int stop=jpt.at(0)-1; for (unsigned int n=0; n >(numAngles); thread_data[n].m_Np=Create2DArray >(numAngles); thread_data[n].m_Lt=Create2DArray >(numAngles); thread_data[n].m_Lp=Create2DArray >(numAngles); boost::thread *t = new boost::thread( nf2ff_calc_thread(this,start,stop,n,thread_data[n]) ); m_thread_group.add_thread( t ); start = stop+1; if (nwait(); //start // threads: calc Js and Ms (eq. 8.15a/b) // threads calc their local Nt,Np,Lt and Lp m_Barrier->wait(); //combine all thread local Nt,Np,Lt and Lp complex** Nt = Create2DArray >(numAngles); complex** Np = Create2DArray >(numAngles); complex** Lt = Create2DArray >(numAngles); complex** Lp = Create2DArray >(numAngles); for (unsigned int n=0; nwait(); //wait for termination m_thread_group.join_all(); // wait for termination delete m_Barrier; m_Barrier = NULL; //cleanup Js & Ms Delete_N_3DArray(Js,numLines); Delete_N_3DArray(Ms,numLines); // calc equations 8.23a/b and 8.24a/b float k = 2*M_PI*m_freq/__C0__*sqrt(m_permittivity*m_permeability); complex factor(0,k/4.0/M_PI/m_radius); complex f_exp(0,-1*k*m_radius); factor *= exp(f_exp); float fZ0 = __Z0__ * sqrt(m_permeability/m_permittivity); complex Z0 = fZ0; float P_max = 0; for (unsigned int tn=0;tnP_max) P_max = m_P_rad[tn][pn]; } //cleanup Nx and Lx Delete2DArray(Nt,numAngles); Delete2DArray(Np,numAngles); Delete2DArray(Lt,numAngles); Delete2DArray(Lp,numAngles); m_maxDir = 4*M_PI*P_max / m_radPower; return true; }