OrderInfo.cpp 28 KB
Newer Older
1 2
#include <inmost.h>

3 4 5
#define GUARD_MPI(x) {ierr = x; if( ierr != MPI_SUCCESS ) {char str[4096]; int len; MPI_Error_string(ierr,str,&len); std::cout << #x << " not successfull: " << str << std::endl; MPI_Abort(comm,-1000);}}
#define HASH_TABLE_SIZE 2048

6 7
#if defined(USE_SOLVER)

Kirill Terekhov committed
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
namespace INMOST
{
	
	int comparator(const void *pa, const void *pb)
	{
		INMOST_DATA_ENUM_TYPE *a = (INMOST_DATA_ENUM_TYPE *) pa, *b = (INMOST_DATA_ENUM_TYPE *) pb;
		return a[0] - b[0];
	}
	
	INMOST_DATA_ENUM_TYPE binary_search_pairs(INMOST_DATA_ENUM_TYPE *link, INMOST_DATA_ENUM_TYPE size, INMOST_DATA_ENUM_TYPE find)
	{
		INMOST_DATA_ENUM_TYPE rcur = size >> 1, lcur = 0, mid, chk;
		while (rcur >= lcur)
		{
			mid = lcur + ((rcur - lcur) >> 1);
			chk = mid << 1;
			if (link[chk] < find) lcur = mid + 1;
			else if (link[chk] > find) rcur = mid - 1;
			else return chk;
		}
		return size;
	}
	
	void Solver::OrderInfo::Integrate(INMOST_DATA_REAL_TYPE *inout, INMOST_DATA_ENUM_TYPE num) const
	{
33
#if defined(USE_MPI)
Kirill Terekhov committed
34
		if (GetSize() == 1) return;
35 36 37
#if defined(USE_OMP)
#pragma omp single
#endif
Kirill Terekhov committed
38 39 40 41 42 43
		{
			int ierr = 0;
			dynarray<INMOST_DATA_REAL_TYPE, 1024> temp(num);
			memcpy(temp.data(), inout, sizeof(INMOST_DATA_REAL_TYPE) * num);
			GUARD_MPI(MPI_Allreduce(temp.data(), inout, num, INMOST_MPI_DATA_REAL_TYPE, MPI_SUM, comm));
		}
44
#else
Kirill Terekhov committed
45 46
		(void) inout;
		(void) num;
47
#endif
Kirill Terekhov committed
48 49 50 51 52 53 54 55 56 57
	}
	
	
	void Solver::OrderInfo::PrepareMatrix(Sparse::Matrix &m, INMOST_DATA_ENUM_TYPE overlap)
	{
		have_matrix = true;
		m.isParallel() = true;
		INMOST_DATA_ENUM_TYPE two[2];
		INMOST_DATA_ENUM_TYPE mbeg, mend;
		int initial_rank;
58
#if defined(USE_MPI)
Kirill Terekhov committed
59 60 61 62 63 64 65 66 67 68 69
		int ierr = 0;
		if (comm != INMOST_MPI_COMM_WORLD)
		{
			MPI_Comm_free(&comm);
			comm = INMOST_MPI_COMM_WORLD;
		}
		if (m.GetCommunicator() == INMOST_MPI_COMM_WORLD)
			comm = INMOST_MPI_COMM_WORLD;
		else MPI_Comm_dup(m.GetCommunicator(), &comm);
		MPI_Comm_rank(comm, &rank);
		MPI_Comm_size(comm, &size);
70
#else
Kirill Terekhov committed
71 72 73
		(void) overlap;
		rank = 0;
		size = 1;
74
#endif
Kirill Terekhov committed
75 76 77 78 79 80 81 82 83 84
		initial_rank = rank;
		//std::vector<MPI_Request> requests;
		global_overlap.resize(size * 2);
		global_to_proc.resize(size + 1);
		m.GetInterval(mbeg, mend);
		global_to_proc[0] = 0;
		initial_matrix_begin = mbeg;
		initial_matrix_end = mend;
		two[0] = mbeg;
		two[1] = mend;
85
#if defined(USE_MPI)
Kirill Terekhov committed
86
		GUARD_MPI(MPI_Allgather(two, 2, INMOST_MPI_DATA_ENUM_TYPE, &global_overlap[0], 2, INMOST_MPI_DATA_ENUM_TYPE, comm));
87
#else
Kirill Terekhov committed
88 89 90
		local_vector_begin = initial_matrix_begin = local_matrix_begin = global_overlap[0] = mbeg;
		local_vector_end   = initial_matrix_end   = local_matrix_end   = global_overlap[1] = mend;
		global_to_proc[1] = mend;
91 92
#endif
#if defined(USE_MPI)
Kirill Terekhov committed
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 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331
		//reorder processors if needed
		{
			//starts of local indexes should appear in asscending order
			bool reorder = false;
			for (int k = 0; k < size - 1; k++)
				if (global_overlap[2 * k] > global_overlap[2 * (k + 1)])
				{
					reorder = true;
					break;
				}
			if (reorder)
			{
				storage_type temp(size * 2);
				//assemble array that includes rank
				for (int k = 0; k < size; ++k)
				{
					temp[2 * k + 0] = global_overlap[2 * k];
					temp[2 * k + 1] = k;
				}
				//sort array
				qsort(&temp[0], size, sizeof(INMOST_DATA_ENUM_TYPE) * 2, comparator);
				//create new group
				MPI_Group oldg, newg;
				MPI_Comm newcomm;
				std::vector<int> ranks(size);
				for (int k = 0; k < size; ++k)
					ranks[k] = temp[2 * k + 1];
				GUARD_MPI(MPI_Comm_group(comm, &oldg));
				GUARD_MPI(MPI_Group_incl(oldg, size, &ranks[0], &newg));
				GUARD_MPI(MPI_Comm_create(comm, newg, &newcomm));
				if (comm != INMOST_MPI_COMM_WORLD)
				{
					GUARD_MPI(MPI_Comm_free(&comm));
				}
				comm = newcomm;
				//compute new rank
				MPI_Comm_rank(comm, &rank);
				//sort array pairs, so we don't need to exchange them again
				qsort(&global_overlap[0], size, sizeof(INMOST_DATA_ENUM_TYPE) * 2, comparator);
			}
			//now check that there are no overlaps of local indexes
			//every mend must be equal to every mbeg
			reorder = false;
			for (int k = 0; k < size - 1; k++)
				if (global_overlap[2 * k + 1] != global_overlap[2 * (k + 1)])
				{
					//check that end is strictly greater then begin
					if (global_overlap[2 * k + 1] < global_overlap[2 * (k + 1)])
					{
						if (initial_rank == 0)
						{
							std::cout << __FILE__ << ":" << __LINE__ << " Matrix index intervals are not complete:";
							std::cout << " processor " << k + 0 << " interval " << global_overlap[2 * (k + 0)] << ":"
							<< global_overlap[2 * (k + 0) + 1];
							std::cout << " processor " << k + 1 << " interval " << global_overlap[2 * (k + 1)] << ":"
							<< global_overlap[2 * (k + 1) + 1];
							std::cout << std::endl;
							MPI_Abort(comm, -1000);
						}
					}
					reorder = true;
				}
			if (reorder)
			{
				storage_type old_overlap(global_overlap);
				//move local bounds to get non-overlapping regions
				for (int k = 0; k < size - 1; k++)
					while (global_overlap[2 * k + 1] > global_overlap[2 * (k + 1)])
					{
						//move bounds to equalize sizes
						if (global_overlap[2 * k + 1] - global_overlap[2 * k] < global_overlap[2 * (k + 1) + 1] - global_overlap[2 * (k + 1)])
							global_overlap[2 * k + 1]--; //move right bound of the current processor to left
						else
							global_overlap[2 * (k + 1)]++; //move left bound of the next processor to right
					}
				
				//TODO: if we need to merge overlapping parts of the matrices - do it here
			}
			local_matrix_begin = global_overlap[2 * rank + 0];
			local_matrix_end = global_overlap[2 * rank + 1];
			for (int k = 0; k < size; k++)
				global_to_proc[k + 1] = global_overlap[2 * k + 1];
		}
		MPI_Status stat;
		INMOST_DATA_ENUM_TYPE ext_pos = local_matrix_end;
		//may replace std::map here
		small_hash<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE, HASH_TABLE_SIZE> global_to_local;
		std::vector<std::pair<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> > current_global_to_local;
		std::vector<Sparse::Row::entry> send_row_data, recv_row_data;
		std::vector<INMOST_DATA_ENUM_TYPE> send_row_sizes, recv_row_sizes;
		std::vector<INMOST_DATA_ENUM_TYPE> incoming(4 * size);
		std::vector<MPI_Request> requests;
		INMOST_DATA_ENUM_TYPE total_send = 0, total_recv = 0;
		INMOST_DATA_ENUM_TYPE local_start = local_matrix_begin, local_end = local_matrix_end;
		for (INMOST_DATA_ENUM_TYPE it = 0; it < overlap + 1; it++)
		{
			total_send = 0, total_recv = 0;
			current_global_to_local.clear();
			for (INMOST_DATA_ENUM_TYPE k = local_start; k < local_end; ++k)
			{
				Sparse::Row &r = m[k];
				INMOST_DATA_ENUM_TYPE jend = r.Size(), ind;
				for (INMOST_DATA_ENUM_TYPE j = 0; j < jend; ++j)
				{
					ind = r.GetIndex(j);
					if (ind < local_matrix_begin || ind >= local_matrix_end)
					{
						INMOST_DATA_ENUM_TYPE &recv_pos = global_to_local[ind];
						if (recv_pos == 0) //this number was not assigned yet
						{
							recv_pos = ext_pos++;
							if (it < overlap) current_global_to_local.push_back(std::make_pair(ind, recv_pos));
						}
					}
				}
			}
			if (it == overlap)
				current_global_to_local = global_to_local.serialize();
			std::sort(current_global_to_local.begin(), current_global_to_local.end());
			//if( !current_global_to_local.empty() )
			{
				//check all the indexes that comes from other processors
				//for every processor we need arrays:
				// processor -> (array of index positions where to receive))
				// processor -> (array of index positions from where to send)
				memset(&incoming[0], 0, sizeof(INMOST_DATA_ENUM_TYPE) * size * 2);
				vector_exchange_recv.clear();
				vector_exchange_recv.push_back(0);
				if (!current_global_to_local.empty())
				{
					INMOST_DATA_ENUM_TYPE proc_beg = GetProcessor(current_global_to_local.begin()->first), proc_end =
					GetProcessor(current_global_to_local.rbegin()->first) + 1;
					INMOST_DATA_ENUM_TYPE current_ind = 0;
					for (INMOST_DATA_ENUM_TYPE proc = proc_beg; proc < proc_end; proc++)
					{
						bool first = true;
						INMOST_DATA_ENUM_TYPE numind = static_cast<INMOST_DATA_ENUM_TYPE>(vector_exchange_recv.size() + 1);
						while (current_ind < current_global_to_local.size() && current_global_to_local[current_ind].first < global_to_proc[proc + 1])
						{
							INMOST_DATA_ENUM_TYPE k = current_global_to_local[current_ind].first;
							if (first)
							{
								vector_exchange_recv.push_back(proc);
								vector_exchange_recv.push_back(1);
								vector_exchange_recv.push_back(k);
								first = false;
							}
							else
							{
								vector_exchange_recv[numind]++;
								vector_exchange_recv.push_back(k);
							}
							current_ind++;
						}
						if (!first)
						{
							incoming[proc]++;
							incoming[proc + size] += vector_exchange_recv[numind];
							vector_exchange_recv[0]++;
						}
					}
				}
				
				GUARD_MPI(MPI_Allreduce(&incoming[0], &incoming[2 * size], size * 2, INMOST_MPI_DATA_ENUM_TYPE, MPI_SUM, comm));
				//std::cout << GetRank() << " MPI_Allreduce " << __FILE__ << ":" << __LINE__ << " incoming " << incoming[size*2+rank] << " size " << incoming[size*3+rank] << std::endl;
				//prepare array that helps exchanging vector values
				requests.resize(2 * vector_exchange_recv[0] + incoming[size * 2 + rank]);
				INMOST_DATA_ENUM_TYPE j = 1;
				for (INMOST_DATA_ENUM_TYPE k = 0; k < vector_exchange_recv[0]; k++) //send rows that i want to receive
				{
					total_recv += vector_exchange_recv[j + 1];
					GUARD_MPI(MPI_Isend(&vector_exchange_recv[j + 1], 1, INMOST_MPI_DATA_ENUM_TYPE, vector_exchange_recv[j], size + vector_exchange_recv[j], comm, &requests[k])); //send number of rows
					GUARD_MPI(MPI_Isend(&vector_exchange_recv[j + 2], vector_exchange_recv[j + 1], INMOST_MPI_DATA_ENUM_TYPE, vector_exchange_recv[j], 2 * size + vector_exchange_recv[j], comm, &requests[k + vector_exchange_recv[0]])); //send row positions
					j += vector_exchange_recv[j + 1] + 2;
				}
				
				recv_row_sizes.resize(incoming[size * 3 + rank]);
				vector_exchange_send.resize(1 + incoming[size * 2 + rank] * 2 + incoming[size * 3 + rank]);
				vector_exchange_send[0] = 0;
				j = 1;
				for (INMOST_DATA_ENUM_TYPE k = 0;
					 k < incoming[size * 2 + rank]; k++) //receive rows that others want from me
				{
					INMOST_DATA_ENUM_TYPE msgsize;
					GUARD_MPI(MPI_Recv(&msgsize, 1, INMOST_MPI_DATA_ENUM_TYPE, MPI_ANY_SOURCE, size + rank, comm, &stat)); //recv number of rows
					vector_exchange_send[j++] = stat.MPI_SOURCE;
					vector_exchange_send[j++] = msgsize;
					//std::cout << GetRank() << " MPI_Irecv size " << msgsize << " rank " << stat.MPI_SOURCE << " tag " << 2*size+rank << __FILE__ << ":" << __LINE__ << std::endl;
					GUARD_MPI(MPI_Irecv(&vector_exchange_send[j], msgsize, INMOST_MPI_DATA_ENUM_TYPE, stat.MPI_SOURCE,2 * size + rank, comm, &requests[2 * vector_exchange_recv[0] + k])); //recv rows
					j += msgsize;
					total_send += msgsize;
					vector_exchange_send[0]++;
				}
				assert(total_send == incoming[size * 3 + rank]);
				assert(vector_exchange_send[0] == incoming[size * 2 + rank]);
				if (2 * vector_exchange_recv[0] + incoming[size * 2 + rank] > 0) GUARD_MPI(MPI_Waitall(2 * vector_exchange_recv[0] + incoming[size * 2 + rank], &requests[0], MPI_STATUSES_IGNORE));
			}
			/*
			 else
			 {
			 vector_exchange_recv.resize(1,0);
			 vector_exchange_send.resize(1,0);
			 }
			 */
			if (it == overlap)
			{
				//std::cout << rank << " reorder " << std::endl;
				//now we need to reorder off-diagonal parts of the matrix
				for (INMOST_DATA_ENUM_TYPE k = local_matrix_begin; k < local_end; ++k)
					for (Sparse::Row::iterator jt = m[k].Begin(); jt != m[k].End(); ++jt)
						if (global_to_local.is_present(jt->first))
							jt->first = global_to_local[jt->first];
						else
						{
							assert(jt->first >= local_matrix_begin);
							assert(jt->first < local_matrix_end);
						}
				local_vector_begin = local_matrix_begin;
				local_vector_end = ext_pos;
				{
					// change indexes for recv array
					INMOST_DATA_ENUM_TYPE i, j = 1, k;
					//for(k = 0; k < GetRank(); k++) MPI_Barrier(comm);
					//std::cout << "rank " << GetRank() << std::endl;
					//std::cout << "recv:" << std::endl;
					for (i = 0; i < vector_exchange_recv[0]; i++)
					{
						//std::cout << "proc " << vector_exchange_recv[j] << " size " << vector_exchange_recv[j+1] << std::endl;
						j++; //skip processor number
						for (k = 0; k < vector_exchange_recv[j]; ++k)
						{
							assert(global_to_local.is_present(vector_exchange_recv[j + k + 1]));
							vector_exchange_recv[j + k + 1] = global_to_local[vector_exchange_recv[j + k + 1]];
							assert(vector_exchange_recv[j + k + 1] >= local_matrix_end);
						}
						j += vector_exchange_recv[j] + 1; //add vector length + size position
					}
					//check that indexes in send array are in local matrix bounds
					//std::cout << "send:" << std::endl;
332
#ifndef NDEBUG
Kirill Terekhov committed
333 334 335 336 337 338 339 340 341 342 343 344
					j = 1;
					for (i = 0; i < vector_exchange_send[0]; i++)
					{
						//std::cout << "proc " << vector_exchange_send[j] << " size " << vector_exchange_send[j+1] << std::endl;
						j++; //skip processor number
						for (k = 0; k < vector_exchange_send[j]; ++k)
						{
							assert(vector_exchange_send[j + k + 1] >= local_matrix_begin);
							assert(vector_exchange_send[j + k + 1] < local_matrix_end);
						}
						j += vector_exchange_send[j] + 1; //add vector length + size position
					}
345
#endif
Kirill Terekhov committed
346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461
					//for(k = GetRank(); k < GetSize(); k++) MPI_Barrier(comm);
				}
				//prepare array local->global
				extended_indexes.resize(local_vector_end - local_matrix_end);
				for (std::vector<std::pair<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> >::iterator jt = current_global_to_local.begin(); jt != current_global_to_local.end(); ++jt)
					extended_indexes[jt->second - local_matrix_end] = jt->first;
				
				send_storage.resize(total_send);
				recv_storage.resize(total_recv);
				send_requests.resize(vector_exchange_send[0]);
				recv_requests.resize(vector_exchange_recv[0]);
				
			}
			else
			{
				send_row_sizes.resize(total_send);
				recv_row_sizes.resize(total_recv);
				
				INMOST_DATA_ENUM_TYPE j = 1, q = 0, f = 0, total_rows_send = 0, total_rows_recv = 0;
				for (INMOST_DATA_ENUM_TYPE k = 0; k < vector_exchange_recv[0]; k++) //recv sizes of rows
				{
					GUARD_MPI(MPI_Irecv(&recv_row_sizes[q], vector_exchange_recv[j + 1], INMOST_MPI_DATA_ENUM_TYPE, vector_exchange_recv[j], 3 * size + vector_exchange_recv[j], comm,&requests[k]));
					q += vector_exchange_recv[j + 1];
					j += vector_exchange_recv[j + 1] + 2;
				}
				j = 1;
				q = 0;
				
				for (INMOST_DATA_ENUM_TYPE k = 0; k < vector_exchange_send[0]; k++) //send sizes of rows
				{
					for (INMOST_DATA_ENUM_TYPE r = 0; r < vector_exchange_send[j + 1]; r++)
					{
						send_row_sizes[q + r] = m[vector_exchange_send[j + 2 + r]].Size();
						total_rows_send += m[vector_exchange_send[j + 2 + r]].Size();
					}
					GUARD_MPI(MPI_Isend(&send_row_sizes[q], vector_exchange_send[j + 1], INMOST_MPI_DATA_ENUM_TYPE,vector_exchange_send[j], 3 * size + rank, comm, &requests[vector_exchange_recv[0] + k])); //recv rows
					//remember processor numbers here
					q += vector_exchange_send[j + 1];
					j += vector_exchange_send[j + 1] + 2;
				}
				send_row_data.clear();
				send_row_data.reserve(total_rows_send);
				
				
				j = 1;
				for (INMOST_DATA_ENUM_TYPE k = 0; k < vector_exchange_send[0]; k++) //accumulate data in array
				{
					for (INMOST_DATA_ENUM_TYPE r = 0; r < vector_exchange_send[j + 1]; r++)
						send_row_data.insert(send_row_data.end(), m[vector_exchange_send[j + 2 + r]].Begin(), m[vector_exchange_send[j + 2 + r]].End());
					j += vector_exchange_send[j + 1] + 2;
				}
				//replace by mpi_waitsome
				if (vector_exchange_recv[0] + vector_exchange_send[0] > 0)
					GUARD_MPI(MPI_Waitall(vector_exchange_recv[0] + vector_exchange_send[0], &requests[0], MPI_STATUSES_IGNORE));
				j = 1;
				q = 0;
				for (INMOST_DATA_ENUM_TYPE k = 0;
					 k < vector_exchange_recv[0]; k++) //compute total size of data to receive
				{
					for (INMOST_DATA_ENUM_TYPE r = 0; r < vector_exchange_recv[j + 1]; r++)
						total_rows_recv += recv_row_sizes[q + r];
					q += vector_exchange_recv[j + 1];
					j += vector_exchange_recv[j + 1] + 2;
				}
				recv_row_data.resize(total_rows_recv);
				j = 1;
				q = 0;
				f = 0;
				for (INMOST_DATA_ENUM_TYPE k = 0; k < vector_exchange_recv[0]; k++) //receive row data
				{
					INMOST_DATA_ENUM_TYPE local_size = 0;
					for (INMOST_DATA_ENUM_TYPE r = 0; r < vector_exchange_recv[j + 1]; r++)
						local_size += recv_row_sizes[q + r];
					GUARD_MPI(MPI_Irecv(&recv_row_data[f], local_size, Sparse::GetRowEntryType(), vector_exchange_recv[j], 4 * size + vector_exchange_recv[j], comm, &requests[k]));
					q += vector_exchange_recv[j + 1];
					j += vector_exchange_recv[j + 1] + 2;
					f += local_size;
				}
				j = 1;
				q = 0;
				f = 0;
				for (INMOST_DATA_ENUM_TYPE k = 0; k < vector_exchange_send[0]; k++) //receive row data
				{
					INMOST_DATA_ENUM_TYPE local_size = 0;
					for (INMOST_DATA_ENUM_TYPE r = 0; r < vector_exchange_send[j + 1]; r++)
						local_size += send_row_sizes[q + r];
					GUARD_MPI( MPI_Isend(&send_row_data[f], local_size, Sparse::GetRowEntryType(), vector_exchange_send[j],4 * size + rank, comm, &requests[k + vector_exchange_recv[0]]));
					q += vector_exchange_send[j + 1];
					j += vector_exchange_send[j + 1] + 2;
					f += local_size;
				}
				local_start = local_end;
				m.SetInterval(local_matrix_begin, ext_pos);
				local_end = ext_pos;
				if (vector_exchange_recv[0] + vector_exchange_send[0] > 0)
					GUARD_MPI(MPI_Waitall(vector_exchange_recv[0] + vector_exchange_send[0], &requests[0],MPI_STATUSES_IGNORE));
				j = 1;
				q = 0;
				f = 0;
				for (INMOST_DATA_ENUM_TYPE k = 0; k < vector_exchange_recv[0]; k++) //extend matrix
				{
					for (INMOST_DATA_ENUM_TYPE r = 0; r < vector_exchange_recv[j + 1]; r++)
					{
						m[global_to_local[vector_exchange_recv[j + 2 + r]]] = Sparse::Row(&recv_row_data[f],&recv_row_data[f] + recv_row_sizes[q + r]);
						f += recv_row_sizes[q + r];
					}
					q += vector_exchange_recv[j + 1];
					j += vector_exchange_recv[j + 1] + 2;
				}
			}
			//std::cout << it << "/" << overlap << " done" << std::endl;
		}
		two[0] = local_matrix_begin;
		two[1] = local_end;
		GUARD_MPI(MPI_Allgather(two, 2, INMOST_MPI_DATA_ENUM_TYPE, &global_overlap[0], 2, INMOST_MPI_DATA_ENUM_TYPE,comm));
		//std::cout << __FUNCTION__ << " done" << std::endl;
462
#endif
Kirill Terekhov committed
463 464 465 466 467 468 469 470 471 472 473 474 475
	}
	
	void Solver::OrderInfo::RestoreMatrix(Sparse::Matrix &m)
	{
		//restore matrix size
		m.SetInterval(initial_matrix_begin, initial_matrix_end);
		//restore indexes
		for (Sparse::Matrix::iterator it = m.Begin(); it != m.End(); ++it)
			for (Sparse::Row::iterator jt = it->Begin(); jt != it->End(); ++jt)
				if (jt->first >= initial_matrix_end)
					jt->first = extended_indexes[jt->first - initial_matrix_end];
		m.isParallel() = false;
		have_matrix = false;
476
#if defined(USE_MPI)
Kirill Terekhov committed
477 478 479 480 481
		if (comm != INMOST_MPI_COMM_WORLD)
		{
			MPI_Comm_free(&comm);
			comm = INMOST_MPI_COMM_WORLD;
		}
482
#endif
Kirill Terekhov committed
483 484 485 486 487
		//std::cout << __FUNCTION__ << std::endl;
	}
	
	Solver::OrderInfo::~OrderInfo()
	{
488
#if defined(USE_MPI)
Kirill Terekhov committed
489 490
		if (comm != INMOST_MPI_COMM_WORLD)
			MPI_Comm_free(&comm);
491
#endif
Kirill Terekhov committed
492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547
	}
	
	void Solver::OrderInfo::Clear()
	{
		global_to_proc.clear();
		global_overlap.clear();
		vector_exchange_recv.clear();
		vector_exchange_send.clear();
		send_storage.clear();
		recv_storage.clear();
		send_requests.clear();
		recv_requests.clear();
		extended_indexes.clear();
		local_vector_begin = local_vector_end = 0;
		initial_matrix_begin = initial_matrix_end = 0;
		local_matrix_begin = local_matrix_end = 0;
		have_matrix = false;
	}
	
	void Solver::OrderInfo::PrepareVector(Sparse::Vector &v) const
	{
		if (!have_matrix) throw PrepareMatrixFirst;
		v.SetInterval(local_vector_begin, local_vector_end);
		v.isParallel() = true;
	}
	
	void Solver::OrderInfo::RestoreVector(Sparse::Vector &v) const
	{
		assert(have_matrix);
		if (v.isParallel())
		{
			v.SetInterval(initial_matrix_begin, initial_matrix_end);
			v.isParallel() = false;
		}
	}
	
	Solver::OrderInfo::OrderInfo()
	: global_to_proc(), global_overlap(), vector_exchange_recv(), vector_exchange_send(), send_storage(), recv_storage(), send_requests(), recv_requests(),
	extended_indexes()
	{
		comm = INMOST_MPI_COMM_WORLD;
		rank = 0;
		size = 1;
		initial_matrix_begin = 0;
		initial_matrix_end = 0;
		local_matrix_begin = 0;
		local_matrix_end = 0;
		local_vector_begin = 0;
		local_vector_end = 0;
		have_matrix = false;
	}
	
	Solver::OrderInfo::OrderInfo(const OrderInfo &other)
	: global_to_proc(other.global_to_proc), global_overlap(other.global_overlap), 	vector_exchange_recv(other.vector_exchange_recv), vector_exchange_send(other.vector_exchange_send),
	extended_indexes(other.extended_indexes)
	{
548
#if defined(USE_MPI)
Kirill Terekhov committed
549 550 551
		if (other.comm == INMOST_MPI_COMM_WORLD)
			comm = INMOST_MPI_COMM_WORLD;
		else MPI_Comm_dup(other.comm, &comm);
552
#else
Kirill Terekhov committed
553
		comm = other.comm;
554
#endif
Kirill Terekhov committed
555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571
		rank = other.rank;
		size = other.size;
		initial_matrix_begin = other.initial_matrix_begin;
		initial_matrix_end = other.initial_matrix_end;
		local_vector_begin = other.local_vector_begin;
		local_vector_end = other.local_vector_end;
		local_matrix_begin = other.local_matrix_begin;
		local_matrix_end = other.local_matrix_end;
		have_matrix = other.have_matrix;
		send_storage.resize(other.send_storage.size());
		recv_storage.resize(other.recv_storage.size());
		send_requests.resize(other.send_requests.size());
		recv_requests.resize(other.recv_requests.size());
	}
	
	Solver::OrderInfo &Solver::OrderInfo::operator=(OrderInfo const &other)
	{
572
#if defined(USE_MPI)
Kirill Terekhov committed
573 574 575
		if (other.comm == INMOST_MPI_COMM_WORLD)
			comm = INMOST_MPI_COMM_WORLD;
		else MPI_Comm_dup(other.comm, &comm);
576
#else
Kirill Terekhov committed
577
		comm = other.comm;
578
#endif
Kirill Terekhov committed
579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630
		global_to_proc = other.global_to_proc;
		global_overlap = other.global_overlap;
		vector_exchange_recv = other.vector_exchange_recv;
		vector_exchange_send = other.vector_exchange_send;
		extended_indexes = other.extended_indexes;
		rank = other.rank;
		size = other.size;
		initial_matrix_begin = other.initial_matrix_begin;
		initial_matrix_end = other.initial_matrix_end;
		local_vector_begin = other.local_vector_begin;
		local_vector_end = other.local_vector_end;
		local_matrix_begin = other.local_matrix_begin;
		local_matrix_end = other.local_matrix_end;
		have_matrix = other.have_matrix;
		send_storage.resize(other.send_storage.size());
		recv_storage.resize(other.recv_storage.size());
		send_requests.resize(other.send_requests.size());
		recv_requests.resize(other.recv_requests.size());
		return *this;
	}
	
	
	INMOST_DATA_ENUM_TYPE Solver::OrderInfo::GetProcessor(INMOST_DATA_ENUM_TYPE gind) const
	{
		assert(have_matrix);
		storage_type::const_iterator find = std::lower_bound(global_to_proc.begin(), global_to_proc.end(), gind);
		assert(find != global_to_proc.end());
		if ((*find) == gind && find + 1 != global_to_proc.end())
			return static_cast<INMOST_DATA_ENUM_TYPE>(find - global_to_proc.begin());
		else return static_cast<INMOST_DATA_ENUM_TYPE>(find - global_to_proc.begin()) - 1;
	}
	
	void Solver::OrderInfo::GetOverlapRegion(INMOST_DATA_ENUM_TYPE proc, INMOST_DATA_ENUM_TYPE &mbeg,
											 INMOST_DATA_ENUM_TYPE &mend) const
	{
		assert(have_matrix);
		mbeg = global_overlap[proc * 2 + 0];
		mend = global_overlap[proc * 2 + 1];
	}
	
	void Solver::OrderInfo::GetLocalRegion(INMOST_DATA_ENUM_TYPE proc, INMOST_DATA_ENUM_TYPE &mbeg,
										   INMOST_DATA_ENUM_TYPE &mend) const
	{
		assert(have_matrix);
		mbeg = global_to_proc[proc + 0];
		mend = global_to_proc[proc + 1];
	}
	
	
	void Solver::OrderInfo::Update(Sparse::Vector &x)
	{
		//std::cout << __FUNCTION__ << " start" << std::endl;
631
#if defined(USE_MPI)
Kirill Terekhov committed
632
		if (GetSize() == 1) return;
633 634 635
#if defined(USE_OMP)
#pragma omp single
#endif
Kirill Terekhov committed
636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674
		{
			//use MPI_Put/MPI_Get to update vector
			assert(x.isParallel()); //the vector was prepared
			INMOST_DATA_ENUM_TYPE i, j = 1, k, l = 0;
			int ierr;
			for (i = 0; i < vector_exchange_recv[0]; i++)
			{
				//std::cout << GetRank() << " MPI_Irecv size " << vector_exchange_recv[j+1] << " dest " << vector_exchange_recv[j] << " tag " << vector_exchange_recv[j]*size+rank << std::endl;
				GUARD_MPI(MPI_Irecv(&recv_storage[l], vector_exchange_recv[j + 1], INMOST_MPI_DATA_REAL_TYPE, vector_exchange_recv[j], vector_exchange_recv[j] * size + rank, comm, &recv_requests[i]));
				l += vector_exchange_recv[j + 1];
				j += vector_exchange_recv[j + 1] + 2;
			}
			j = 1, l = 0;
			for (i = 0; i < vector_exchange_send[0]; i++)
			{
				//std::cout << GetRank() << " MPI_Isend size " << vector_exchange_send[j+1] << " dest " << vector_exchange_send[j] << " tag " << rank*size+vector_exchange_send[j] << std::endl;
				for (k = 0; k < vector_exchange_send[j + 1]; k++)
					send_storage[l + k] = x[vector_exchange_send[k + j + 2]];
				GUARD_MPI(MPI_Isend(&send_storage[l], vector_exchange_send[j + 1], INMOST_MPI_DATA_REAL_TYPE, vector_exchange_send[j], rank * size + vector_exchange_send[j], comm, &send_requests[i]));
				l += vector_exchange_send[j + 1];
				j += vector_exchange_send[j + 1] + 2;
			}
			if (vector_exchange_recv[0] > 0)
			{
				GUARD_MPI(MPI_Waitall(static_cast<int>(recv_requests.size()), &recv_requests[0], MPI_STATUSES_IGNORE));
				j = 1, l = 0;
				for (i = 0; i < vector_exchange_recv[0]; i++)
				{
					for (k = 0; k < vector_exchange_recv[j + 1]; k++)
						x[vector_exchange_recv[k + j + 2]] = recv_storage[l + k];
					l += vector_exchange_recv[j + 1];
					j += vector_exchange_recv[j + 1] + 2;
				}
			}
			if (vector_exchange_send[0] > 0)
			{
				GUARD_MPI(MPI_Waitall(static_cast<int>(send_requests.size()), &send_requests[0], MPI_STATUSES_IGNORE));
			}
		}
675
#else
Kirill Terekhov committed
676
		(void) x;
677
#endif
Kirill Terekhov committed
678 679 680 681 682 683
		//std::cout << __FUNCTION__ << " end" << std::endl;
	}
	
	void Solver::OrderInfo::Accumulate(Sparse::Vector &x)
	{
		//std::cout << __FUNCTION__ << " start" << std::endl;
684
#if defined(USE_MPI)
Kirill Terekhov committed
685
		if (GetSize() == 1) return;
686 687 688
#if defined(USE_OMP)
#pragma omp single
#endif
Kirill Terekhov committed
689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733
		{
			//use MPI_Put/MPI_Get to update vector
			assert(x.isParallel()); //the vector was prepared
			INMOST_DATA_ENUM_TYPE i, j = 1, k, l = 0;
			int ierr;
			for (i = 0; i < vector_exchange_send[0]; i++)
			{
				//std::cout << GetRank() << " MPI_Irecv size " << vector_exchange_send[j+1] << " dest " << vector_exchange_send[j] << " tag " << vector_exchange_send[j]*size+rank << std::endl;
				GUARD_MPI(MPI_Irecv(&send_storage[l], vector_exchange_send[j + 1], INMOST_MPI_DATA_REAL_TYPE,
									vector_exchange_send[j], vector_exchange_send[j] * size + rank, comm,
									&send_requests[i]));
				l += vector_exchange_send[j + 1];
				j += vector_exchange_send[j + 1] + 2;
			}
			j = 1, l = 0;
			for (i = 0; i < vector_exchange_recv[0]; i++)
			{
				for (k = 0; k < vector_exchange_recv[j + 1]; k++)
					recv_storage[l + k] = x[vector_exchange_recv[k + j + 2]];
				//std::cout << GetRank() << " MPI_Isend size " << vector_exchange_recv[j+1] << " dest " << vector_exchange_recv[j] << " tag " << rank*size+vector_exchange_recv[j] << std::endl;
				GUARD_MPI(MPI_Isend(&recv_storage[l], vector_exchange_recv[j + 1], INMOST_MPI_DATA_REAL_TYPE,
									vector_exchange_recv[j], rank * size + vector_exchange_recv[j], comm,
									&recv_requests[i]));
				l += vector_exchange_recv[j + 1];
				j += vector_exchange_recv[j + 1] + 2;
			}
			if (vector_exchange_send[0] > 0)
			{
				//std::cout << GetRank() << " Waitall send " << send_requests.size() << std::endl;
				GUARD_MPI(MPI_Waitall(static_cast<int>(send_requests.size()), &send_requests[0], MPI_STATUSES_IGNORE));
				j = 1, l = 0;
				for (i = 0; i < vector_exchange_send[0]; i++)
				{
					for (k = 0; k < vector_exchange_send[j + 1]; k++)
						x[vector_exchange_send[k + j + 2]] += send_storage[l + k];
					l += vector_exchange_send[j + 1];
					j += vector_exchange_send[j + 1] + 2;
				}
			}
			if (vector_exchange_recv[0] > 0)
			{
				//std::cout << GetRank() << " Waitall recv " << recv_requests.size() << std::endl;
				GUARD_MPI(MPI_Waitall(static_cast<int>(recv_requests.size()), &recv_requests[0], MPI_STATUSES_IGNORE));
			}
		}
734
#else
Kirill Terekhov committed
735
		(void) x;
736
#endif
Kirill Terekhov committed
737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753
		//std::cout << __FUNCTION__ << " end" << std::endl;
	}
	/*
	 void Solver::OrderInfo::ScalarProd(Vector const & left, Vector const & right, INMOST_DATA_ENUM_TYPE index_begin, INMOST_DATA_ENUM_TYPE index_end, INMOST_DATA_REAL_TYPE & sum) const
	 {
	 INMOST_DATA_INTEGER_TYPE ibeg = index_begin, iend = index_end;
	 #if defined(USE_OMP)
	 #pragma omp for reduction(+:sum)
	 #endif
	 for(INMOST_DATA_INTEGER_TYPE i = ibeg; i < iend; ++i)
	 {
	 sum += left[i]*right[i];
	 }
	 Integrate(&sum,1);
	 }
	 */
	
754 755
}

Kirill Terekhov committed
756
#endif//USE_SOLVER