summaryrefslogtreecommitdiff
path: root/openEMS/FDTD/engine_cylindermultigrid.cpp
blob: 8e5ce85cd20d44f57884d06d487b2b956a29aa40 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
/*
*	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 <http://www.gnu.org/licenses/>.
*/

#include "engine_cylindermultigrid.h"
#include "operator_cylindermultigrid.h"
#include "extensions/engine_ext_cylindermultigrid.h"

Engine_CylinderMultiGrid* Engine_CylinderMultiGrid::New(const Operator_CylinderMultiGrid* op, unsigned int numThreads)
{
	cout << "Create FDTD engine (cylindrical multi grid mesh using sse compression + multithreading)" << endl;
	Engine_CylinderMultiGrid* e = new Engine_CylinderMultiGrid(op);
	e->setNumThreads( numThreads );
	e->Init();
	return e;
}

Engine_CylinderMultiGrid::Engine_CylinderMultiGrid(const Operator_CylinderMultiGrid* op) : Engine_Cylinder(op)
{
	Op_CMG = op;

	m_WaitOnBase = new boost::barrier(2);
	m_WaitOnChild = new boost::barrier(2);
	m_WaitOnSync = new boost::barrier(2);

	m_Eng_Ext_MG = new Engine_Ext_CylinderMultiGrid(NULL,true);
	m_Eng_Ext_MG->SetBarrier(m_WaitOnBase, m_WaitOnChild, m_WaitOnSync);
	m_Eng_Ext_MG->SetEngine(this);

	Engine* eng = op->GetInnerOperator()->CreateEngine();
	m_InnerEngine = dynamic_cast<Engine_Multithread*>(eng);

	Engine_Ext_CylinderMultiGrid* m_InnerEng_Ext_MG = new Engine_Ext_CylinderMultiGrid(NULL,false);
	m_InnerEng_Ext_MG->SetBarrier(m_WaitOnBase, m_WaitOnChild, m_WaitOnSync);

	// if already has a base extension, switch places ... seems to be faster...
	for (size_t n=0; n<m_InnerEngine->m_Eng_exts.size(); ++n)
	{
		Engine_Ext_CylinderMultiGrid* eng_mg = dynamic_cast<Engine_Ext_CylinderMultiGrid*>(m_InnerEngine->m_Eng_exts.at(n));
		if (eng_mg)
		{
			m_InnerEngine->m_Eng_exts.at(n) = m_InnerEng_Ext_MG;
			m_InnerEng_Ext_MG = eng_mg;
			break;
		}
	}
	m_InnerEngine->m_Eng_exts.push_back(m_InnerEng_Ext_MG);
}

Engine_CylinderMultiGrid::~Engine_CylinderMultiGrid()
{
#ifdef MPI_SUPPORT
	delete m_InnerEngine->m_MPI_Barrier;
	m_InnerEngine->m_MPI_Barrier = NULL;
#endif

	m_Thread_NumTS = 0;
	m_startBarrier->wait();

	m_IteratorThread_Group.join_all();

	delete m_InnerEngine;
	m_InnerEngine = NULL;

	delete m_WaitOnBase;
	m_WaitOnBase = NULL;
	delete m_WaitOnChild;
	m_WaitOnChild = NULL;
	delete m_WaitOnSync;
	m_WaitOnSync = NULL;

	delete m_startBarrier;
	m_startBarrier = NULL;
	delete m_stopBarrier;
	m_stopBarrier = NULL;
}

void Engine_CylinderMultiGrid::Init()
{
	Engine_Multithread::Init();

	m_Eng_exts.push_back(m_Eng_Ext_MG);

	m_startBarrier = new boost::barrier(3); //both engines + organizer
	m_stopBarrier = new boost::barrier(3); //both engines + organizer

	boost::thread *t = NULL;

	t = new boost::thread( Engine_CylinderMultiGrid_Thread(this,m_startBarrier,m_stopBarrier,&m_Thread_NumTS, true) );
	m_IteratorThread_Group.add_thread( t );

	t = new boost::thread( Engine_CylinderMultiGrid_Thread(m_InnerEngine,m_startBarrier,m_stopBarrier,&m_Thread_NumTS, false) );
	m_IteratorThread_Group.add_thread( t );

	m_InnerEngine->SortExtensionByPriority();
	SortExtensionByPriority();

#ifdef MPI_SUPPORT
	//assign an MPI barrier to inner Engine
	m_InnerEngine->m_MPI_Barrier  = new boost::barrier(2);
#endif
}

bool Engine_CylinderMultiGrid::IterateTS(unsigned int iterTS)
{
	m_Thread_NumTS = iterTS;

	m_startBarrier->wait(); //start base and child iterations

	m_stopBarrier->wait();  //tell base and child to wait for another start event...

	//interpolate child data to base mesh...
	for (unsigned int n=0; n<Op_CMG->m_Split_Pos-1; ++n)
		InterpolVoltChild2Base(n);
	for (unsigned int n=0; n<Op_CMG->m_Split_Pos-2; ++n)
		InterpolCurrChild2Base(n);

	return true;
}

void Engine_CylinderMultiGrid::InterpolVoltChild2Base(unsigned int rPos)
{
	//interpolate voltages from child engine to the base engine...
	unsigned int pos[3];
	pos[0] = rPos;
	for (pos[1]=0; pos[1]<numLines[1]; ++pos[1])
	{
		for (pos[2]=0; pos[2]<numVectors; ++pos[2])
		{
			//r - direction
			f4_volt[0][pos[0]][pos[1]][pos[2]].v  = Op_CMG->f4_interpol_v_2p[0][pos[1]].v * m_InnerEngine->f4_volt[0][pos[0]][Op_CMG->m_interpol_pos_v_2p[0][pos[1]]][pos[2]].v
													+ Op_CMG->f4_interpol_v_2pp[0][pos[1]].v * m_InnerEngine->f4_volt[0][pos[0]][Op_CMG->m_interpol_pos_v_2pp[0][pos[1]]][pos[2]].v;

			//z - direction
			f4_volt[2][pos[0]][pos[1]][pos[2]].v  = Op_CMG->f4_interpol_v_2p[0][pos[1]].v * m_InnerEngine->f4_volt[2][pos[0]][Op_CMG->m_interpol_pos_v_2p[0][pos[1]]][pos[2]].v
													+ Op_CMG->f4_interpol_v_2pp[0][pos[1]].v * m_InnerEngine->f4_volt[2][pos[0]][Op_CMG->m_interpol_pos_v_2pp[0][pos[1]]][pos[2]].v;

			//alpha - direction
			f4_volt[1][pos[0]][pos[1]][pos[2]].v  = Op_CMG->f4_interpol_v_2p[1][pos[1]].v * m_InnerEngine->f4_volt[1][pos[0]][Op_CMG->m_interpol_pos_v_2p[1][pos[1]]][pos[2]].v
													+ Op_CMG->f4_interpol_v_2pp[1][pos[1]].v * m_InnerEngine->f4_volt[1][pos[0]][Op_CMG->m_interpol_pos_v_2pp[1][pos[1]]][pos[2]].v;
		}
	}
}

void Engine_CylinderMultiGrid::InterpolCurrChild2Base(unsigned int rPos)
{
	//interpolate voltages from child engine to the base engine...
	unsigned int pos[3];
	pos[0] = rPos;
	for (pos[1]=0; pos[1]<numLines[1]; ++pos[1])
	{
		for (pos[2]=0; pos[2]<numVectors; ++pos[2])
		{
			//r - direction
			f4_curr[0][pos[0]][pos[1]][pos[2]].v  = Op_CMG->f4_interpol_i_2p[0][pos[1]].v * m_InnerEngine->f4_curr[0][pos[0]][Op_CMG->m_interpol_pos_i_2p[0][pos[1]]][pos[2]].v
													+ Op_CMG->f4_interpol_i_2pp[0][pos[1]].v * m_InnerEngine->f4_curr[0][pos[0]][Op_CMG->m_interpol_pos_i_2pp[0][pos[1]]][pos[2]].v;

			//z - direction
			f4_curr[2][pos[0]][pos[1]][pos[2]].v  = Op_CMG->f4_interpol_i_2p[0][pos[1]].v * m_InnerEngine->f4_curr[2][pos[0]][Op_CMG->m_interpol_pos_i_2p[0][pos[1]]][pos[2]].v
													+ Op_CMG->f4_interpol_i_2pp[0][pos[1]].v * m_InnerEngine->f4_curr[2][pos[0]][Op_CMG->m_interpol_pos_i_2pp[0][pos[1]]][pos[2]].v;

			//alpha - direction
			f4_curr[1][pos[0]][pos[1]][pos[2]].v  = Op_CMG->f4_interpol_i_2p[1][pos[1]].v * m_InnerEngine->f4_curr[1][pos[0]][Op_CMG->m_interpol_pos_i_2p[1][pos[1]]][pos[2]].v
													+ Op_CMG->f4_interpol_i_2pp[1][pos[1]].v * m_InnerEngine->f4_curr[1][pos[0]][Op_CMG->m_interpol_pos_i_2pp[1][pos[1]]][pos[2]].v;
		}
	}
}

#ifdef MPI_SUPPORT
	void Engine_CylinderMultiGrid::SendReceiveVoltages()
	{
		//do the local voltage sync, child is waiting...
		Engine_Multithread::SendReceiveVoltages();

		//run inner voltage sync
		m_InnerEngine->m_MPI_Barrier->wait();
	}

	void Engine_CylinderMultiGrid::SendReceiveCurrents()
	{
		//do the local current sync, child is waiting...
		Engine_Multithread::SendReceiveCurrents();

		//run inner voltage sync
		m_InnerEngine->m_MPI_Barrier->wait();
	}

#endif

/****************************************************************************************/
Engine_CylinderMultiGrid_Thread::Engine_CylinderMultiGrid_Thread( Engine_Multithread* engine, boost::barrier *start, boost::barrier *stop, volatile unsigned int* numTS, bool isBase)
{
	m_startBarrier = start;
	m_stopBarrier = stop;
	m_Eng=engine;
	m_isBase=isBase;
	m_numTS = numTS;
}

void Engine_CylinderMultiGrid_Thread::operator()()
{
	m_startBarrier->wait(); //wait for Base engine to start the iterations...

	while (*m_numTS>0)	//m_numTS==0 request to terminate this thread...
	{
		if (m_isBase)
			m_Eng->Engine_Multithread::IterateTS(*m_numTS);
		else
			m_Eng->IterateTS(*m_numTS);
		m_stopBarrier->wait();		//sync all workers after iterations are performed
		m_startBarrier->wait();		//wait for Base engine to start the iterations again ...
	}
}