| 1 | /*******************************************************************************
|
|---|
| 2 |
|
|---|
| 3 | BoxManager:
|
|---|
| 4 |
|
|---|
| 5 | purpose:: organize arbitrary information in a spatial way
|
|---|
| 6 |
|
|---|
| 7 | misc. notes/considerations/open questions:
|
|---|
| 8 |
|
|---|
| 9 | (1) In the struct code, we want to use Box Manager instead of
|
|---|
| 10 | current box neighbor stuff (see Struct function
|
|---|
| 11 | hypre_CreateCommInfoFromStencil. For example, to get neighbors of
|
|---|
| 12 | box b, we can call Intersect with a larger box than b).
|
|---|
| 13 |
|
|---|
| 14 | (2) will associate a Box Manager with the struct grid (implement
|
|---|
| 15 | under the struct grid)
|
|---|
| 16 |
|
|---|
| 17 | (3) will interface with the Box Manager in the struct coarsen routine
|
|---|
| 18 |
|
|---|
| 19 | the coarsen routine:
|
|---|
| 20 |
|
|---|
| 21 | (a) get all the box manager entries from the current level,
|
|---|
| 22 | coarsen them, and create a new box manager for the coarse grid,
|
|---|
| 23 | adding the boxes via AddEntry
|
|---|
| 24 |
|
|---|
| 25 | (b) check the max_distance value and see if we have
|
|---|
| 26 | all the neighbor info we need in the current box manager.
|
|---|
| 27 |
|
|---|
| 28 | (c) if (b) is no, then call GatherEntries as needed on the coarse
|
|---|
| 29 | box manager
|
|---|
| 30 |
|
|---|
| 31 |
|
|---|
| 32 | (d) call assemble for the new coarse box manager (note: if gather
|
|---|
| 33 | entries has not been called, then no communication is required
|
|---|
| 34 |
|
|---|
| 35 | (4) We will associate an assumed partition with the box manager
|
|---|
| 36 | (this will be created in the box manager assemble routine)
|
|---|
| 37 |
|
|---|
| 38 | (5) how will we use the box manager with sstruct? "on the side" as
|
|---|
| 39 | the boxmap is now, or through the struct grid (at issue is modifying
|
|---|
| 40 | the "info" associated with an entry after the box manager has
|
|---|
| 41 | already been assembled through the underlying struct grid)
|
|---|
| 42 |
|
|---|
| 43 | (6) Populating the "info" in sstruct, we need to eliminate global
|
|---|
| 44 | storage when computing offsets (can probably use MPI_SCan as in
|
|---|
| 45 | parcsr_ls/par_coarse_parms.c)
|
|---|
| 46 |
|
|---|
| 47 | (7) In SStruct we will have a separate box manager for the
|
|---|
| 48 | neighbor box information
|
|---|
| 49 |
|
|---|
| 50 | ********************************************************************************/
|
|---|
| 51 |
|
|---|
| 52 | #include "headers.h"
|
|---|
| 53 |
|
|---|
| 54 |
|
|---|
| 55 | /*--------------------------------------------------------------------------
|
|---|
| 56 | hypre_BoxManEntrySetInfo
|
|---|
| 57 | *--------------------------------------------------------------------------*/
|
|---|
| 58 |
|
|---|
| 59 |
|
|---|
| 60 | int hypre_BoxManEntrySetInfo ( hypre_BoxManEntry *entry , void *info )
|
|---|
| 61 | {
|
|---|
| 62 |
|
|---|
| 63 | hypre_BoxManEntryInfo(entry) = info;
|
|---|
| 64 |
|
|---|
| 65 | return hypre_error_flag;
|
|---|
| 66 |
|
|---|
| 67 | }
|
|---|
| 68 |
|
|---|
| 69 | /*--------------------------------------------------------------------------
|
|---|
| 70 | hypre_BoxManEntryGetInfo
|
|---|
| 71 | *--------------------------------------------------------------------------*/
|
|---|
| 72 |
|
|---|
| 73 |
|
|---|
| 74 | int hypre_BoxManEntryGetInfo ( hypre_BoxManEntry *entry , void **info_ptr )
|
|---|
| 75 | {
|
|---|
| 76 |
|
|---|
| 77 | *info_ptr = hypre_BoxManEntryInfo(entry);
|
|---|
| 78 |
|
|---|
| 79 | return hypre_error_flag;
|
|---|
| 80 |
|
|---|
| 81 | }
|
|---|
| 82 |
|
|---|
| 83 | /*--------------------------------------------------------------------------
|
|---|
| 84 | hypre_BoxManEntryGetExtents
|
|---|
| 85 | *--------------------------------------------------------------------------*/
|
|---|
| 86 |
|
|---|
| 87 |
|
|---|
| 88 | int hypre_BoxManEntryGetExtents ( hypre_BoxManEntry *entry , hypre_Index imin ,
|
|---|
| 89 | hypre_Index imax )
|
|---|
| 90 | {
|
|---|
| 91 |
|
|---|
| 92 |
|
|---|
| 93 | hypre_IndexRef entry_imin = hypre_BoxManEntryIMin(entry);
|
|---|
| 94 | hypre_IndexRef entry_imax = hypre_BoxManEntryIMax(entry);
|
|---|
| 95 |
|
|---|
| 96 | int d;
|
|---|
| 97 |
|
|---|
| 98 | for (d = 0; d < 3; d++)
|
|---|
| 99 | {
|
|---|
| 100 | hypre_IndexD(imin, d) = hypre_IndexD(entry_imin, d);
|
|---|
| 101 | hypre_IndexD(imax, d) = hypre_IndexD(entry_imax, d);
|
|---|
| 102 | }
|
|---|
| 103 |
|
|---|
| 104 |
|
|---|
| 105 | return hypre_error_flag;
|
|---|
| 106 | }
|
|---|
| 107 |
|
|---|
| 108 |
|
|---|
| 109 | /*--------------------------------------------------------------------------
|
|---|
| 110 | hypre_BoxManEntryCopy
|
|---|
| 111 | *--------------------------------------------------------------------------*/
|
|---|
| 112 |
|
|---|
| 113 | int hypre_BoxManEntryCopy( hypre_BoxManEntry *fromentry ,
|
|---|
| 114 | hypre_BoxManEntry *toentry)
|
|---|
| 115 | {
|
|---|
| 116 | int d;
|
|---|
| 117 |
|
|---|
| 118 | hypre_Index imin;
|
|---|
| 119 | hypre_Index imax;
|
|---|
| 120 |
|
|---|
| 121 | hypre_IndexRef toentry_imin;
|
|---|
| 122 | hypre_IndexRef toentry_imax;
|
|---|
| 123 |
|
|---|
| 124 |
|
|---|
| 125 | /* copy extents */
|
|---|
| 126 | hypre_BoxManEntryGetExtents( fromentry, imin, imax );
|
|---|
| 127 |
|
|---|
| 128 | toentry_imin = hypre_BoxManEntryIMin(toentry);
|
|---|
| 129 | toentry_imax = hypre_BoxManEntryIMax(toentry);
|
|---|
| 130 |
|
|---|
| 131 | for (d = 0; d < 3; d++)
|
|---|
| 132 | {
|
|---|
| 133 | hypre_IndexD(toentry_imin, d) = hypre_IndexD(imin, d);
|
|---|
| 134 | hypre_IndexD(toentry_imax, d) = hypre_IndexD(imax, d);
|
|---|
| 135 | }
|
|---|
| 136 |
|
|---|
| 137 | /* copy proc and id */
|
|---|
| 138 | hypre_BoxManEntryProc(toentry) = hypre_BoxManEntryProc(fromentry);
|
|---|
| 139 | hypre_BoxManEntryId(toentry) = hypre_BoxManEntryId(fromentry);
|
|---|
| 140 |
|
|---|
| 141 | /*copy ghost */
|
|---|
| 142 | for (d = 0; d < 6; d++)
|
|---|
| 143 | {
|
|---|
| 144 | hypre_BoxManEntryNumGhost(toentry)[d] = hypre_BoxManEntryNumGhost(fromentry)[d];
|
|---|
| 145 | }
|
|---|
| 146 |
|
|---|
| 147 | /* copy info */
|
|---|
| 148 | hypre_BoxManEntryInfo(toentry) = hypre_BoxManEntryInfo(fromentry) ;
|
|---|
| 149 |
|
|---|
| 150 | /* copy list pointer */
|
|---|
| 151 | hypre_BoxManEntryNext(toentry) = hypre_BoxManEntryNext(fromentry);
|
|---|
| 152 |
|
|---|
| 153 |
|
|---|
| 154 |
|
|---|
| 155 | return hypre_error_flag;
|
|---|
| 156 | }
|
|---|
| 157 |
|
|---|
| 158 | /*--------------------------------------------------------------------------
|
|---|
| 159 | hypre_BoxManDeleteMultipleEntries
|
|---|
| 160 |
|
|---|
| 161 | Delete multiple entries from the manager. The indices correcpond to the
|
|---|
| 162 | ordering of the entires. Assumes indices given in ascending order -
|
|---|
| 163 | this is meant for internal use inside the Assemble routime.
|
|---|
| 164 |
|
|---|
| 165 | *--------------------------------------------------------------------------*/
|
|---|
| 166 |
|
|---|
| 167 | int hypre_BoxManDeleteMultipleEntries( hypre_BoxManager *manager,
|
|---|
| 168 | int* indices , int num )
|
|---|
| 169 | {
|
|---|
| 170 |
|
|---|
| 171 | int i, j, start;
|
|---|
| 172 | int array_size = hypre_BoxManNEntries(manager);
|
|---|
| 173 |
|
|---|
| 174 | hypre_BoxManEntry *entries = hypre_BoxManEntries(manager);
|
|---|
| 175 |
|
|---|
| 176 | if (num > 0)
|
|---|
| 177 | {
|
|---|
| 178 | start = indices[0];
|
|---|
| 179 |
|
|---|
| 180 | j = 0;
|
|---|
| 181 |
|
|---|
| 182 | for (i = start; (i + j) < array_size; i++)
|
|---|
| 183 | {
|
|---|
| 184 | if (j < num)
|
|---|
| 185 | {
|
|---|
| 186 | while ((i+j) == indices[j]) /* see if deleting consecutive items */
|
|---|
| 187 | {
|
|---|
| 188 | j++; /*increase the shift*/
|
|---|
| 189 | if (j == num) break;
|
|---|
| 190 | }
|
|---|
| 191 | }
|
|---|
| 192 |
|
|---|
| 193 | if ( (i+j) < array_size) /* if deleting the last item then no moving */
|
|---|
| 194 | {
|
|---|
| 195 | hypre_BoxManEntryCopy(&entries[i+j], &entries[i]);
|
|---|
| 196 | }
|
|---|
| 197 | }
|
|---|
| 198 | hypre_BoxManNEntries(manager) = array_size - num;
|
|---|
| 199 | }
|
|---|
| 200 |
|
|---|
| 201 | return hypre_error_flag;
|
|---|
| 202 |
|
|---|
| 203 | }
|
|---|
| 204 |
|
|---|
| 205 |
|
|---|
| 206 | /*--------------------------------------------------------------------------
|
|---|
| 207 | hypre_BoxManCreate:
|
|---|
| 208 |
|
|---|
| 209 | Allocate and initialize the box manager structure.
|
|---|
| 210 |
|
|---|
| 211 | Notes:
|
|---|
| 212 |
|
|---|
| 213 | (1) max_nentries indicates how much storage you think you will need
|
|---|
| 214 | for adding entries with BoxManAddEntry
|
|---|
| 215 |
|
|---|
| 216 | (2) info_size indicates the size (in bytes) of the info object that
|
|---|
| 217 | will be attached to each entry in this box manager. (In the future, may
|
|---|
| 218 | want to let the info size be stored in the entry - so specified with
|
|---|
| 219 | AddEntry call. Then we should adjust the ExchangeData function to
|
|---|
| 220 | allow for different sizes - probably a nontrivial change.)
|
|---|
| 221 |
|
|---|
| 222 | (3) we are *not* collecting the bounding box
|
|---|
| 223 |
|
|---|
| 224 | (4) comm is needed for later calls to addentry - also used in the assemble
|
|---|
| 225 |
|
|---|
| 226 |
|
|---|
| 227 | *--------------------------------------------------------------------------*/
|
|---|
| 228 |
|
|---|
| 229 | int hypre_BoxManCreate ( int max_nentries , int info_size, int dim, MPI_Comm comm,
|
|---|
| 230 | hypre_BoxManager **manager_ptr )
|
|---|
| 231 |
|
|---|
| 232 | {
|
|---|
| 233 |
|
|---|
| 234 | hypre_BoxManager *manager;
|
|---|
| 235 |
|
|---|
| 236 | int i, d;
|
|---|
| 237 |
|
|---|
| 238 | /* allocate object */
|
|---|
| 239 | manager = hypre_CTAlloc(hypre_BoxManager, 1);
|
|---|
| 240 |
|
|---|
| 241 |
|
|---|
| 242 | /* initialize */
|
|---|
| 243 | hypre_BoxManComm(manager) = comm;
|
|---|
| 244 | hypre_BoxManMaxNEntries(manager) = max_nentries;
|
|---|
| 245 | hypre_BoxManEntryInfoSize(manager) = info_size;
|
|---|
| 246 | hypre_BoxManDim(manager) = dim;
|
|---|
| 247 |
|
|---|
| 248 | for (d = 0; d < 3; d++)
|
|---|
| 249 | {
|
|---|
| 250 | hypre_BoxManIndexesD(manager, d) = NULL;
|
|---|
| 251 | }
|
|---|
| 252 |
|
|---|
| 253 |
|
|---|
| 254 | hypre_BoxManNEntries(manager) = 0;
|
|---|
| 255 | hypre_BoxManEntries(manager) = hypre_CTAlloc(hypre_BoxManEntry, max_nentries);
|
|---|
| 256 | hypre_BoxManIndexTable(manager) = NULL;
|
|---|
| 257 | hypre_BoxManSortTable(manager) = NULL;
|
|---|
| 258 |
|
|---|
| 259 | hypre_BoxManNumProcsSort(manager) = 0;
|
|---|
| 260 | hypre_BoxManIdsSort(manager) = hypre_CTAlloc(int, max_nentries);
|
|---|
| 261 | hypre_BoxManProcsSort(manager) = hypre_CTAlloc(int, max_nentries);
|
|---|
| 262 | hypre_BoxManProcsSortOffsets(manager) = NULL;
|
|---|
| 263 |
|
|---|
| 264 | hypre_BoxManFirstLocal(manager) = 0;
|
|---|
| 265 | hypre_BoxManLocalProcOffset(manager) = 0;
|
|---|
| 266 |
|
|---|
| 267 | hypre_BoxManIsGatherCalled(manager) = 0;
|
|---|
| 268 | hypre_BoxManGatherRegions(manager) = hypre_BoxArrayCreate(0);
|
|---|
| 269 | hypre_BoxManAllGlobalKnown(manager) = 0;
|
|---|
| 270 |
|
|---|
| 271 | hypre_BoxManNumMyEntries(manager) = 0;
|
|---|
| 272 | hypre_BoxManMyIds(manager) = NULL;
|
|---|
| 273 | hypre_BoxManMyEntries(manager) = NULL;
|
|---|
| 274 |
|
|---|
| 275 | hypre_BoxManAssumedPartition(manager) = NULL;
|
|---|
| 276 |
|
|---|
| 277 | #ifdef HYPRE_NO_GLOBAL_PARTITION
|
|---|
| 278 | hypre_BoxManMyIds(manager) = hypre_CTAlloc(int, max_nentries);
|
|---|
| 279 | hypre_BoxManMyEntries(manager) = hypre_CTAlloc(hypre_BoxManEntry *, max_nentries);
|
|---|
| 280 | #endif
|
|---|
| 281 |
|
|---|
| 282 | /* ghost points: we choose a default that will give zero everywhere..*/
|
|---|
| 283 | for (i = 0; i < 6; i++)
|
|---|
| 284 | {
|
|---|
| 285 | hypre_BoxManNumGhost(manager)[i] = 0;
|
|---|
| 286 | }
|
|---|
| 287 |
|
|---|
| 288 |
|
|---|
| 289 | /* return */
|
|---|
| 290 | *manager_ptr = manager;
|
|---|
| 291 |
|
|---|
| 292 | return hypre_error_flag;
|
|---|
| 293 |
|
|---|
| 294 | }
|
|---|
| 295 |
|
|---|
| 296 | /*--------------------------------------------------------------------------
|
|---|
| 297 | hypre_BoxManIncSize:
|
|---|
| 298 |
|
|---|
| 299 | Increase storage for entries (for future calls to BoxManAddEntry).
|
|---|
| 300 |
|
|---|
| 301 |
|
|---|
| 302 | Notes:
|
|---|
| 303 |
|
|---|
| 304 | In addition, we will dynamically allocate more memory
|
|---|
| 305 | if needed when a call to BoxManAddEntry is made and there is not
|
|---|
| 306 | enough storage available.
|
|---|
| 307 |
|
|---|
| 308 | *--------------------------------------------------------------------------*/
|
|---|
| 309 |
|
|---|
| 310 | int hypre_BoxManIncSize ( hypre_BoxManager *manager , int inc_size)
|
|---|
| 311 | {
|
|---|
| 312 |
|
|---|
| 313 |
|
|---|
| 314 | int max_nentries = hypre_BoxManMaxNEntries(manager);
|
|---|
| 315 | int *ids = hypre_BoxManIdsSort(manager);
|
|---|
| 316 | int *procs = hypre_BoxManProcsSort(manager);
|
|---|
| 317 |
|
|---|
| 318 | hypre_BoxManEntry *entries = hypre_BoxManEntries(manager);
|
|---|
| 319 |
|
|---|
| 320 | /* increase size */
|
|---|
| 321 | max_nentries += inc_size;
|
|---|
| 322 | entries = hypre_TReAlloc(entries, hypre_BoxManEntry, max_nentries);
|
|---|
| 323 | ids = hypre_TReAlloc(ids, int, max_nentries);
|
|---|
| 324 | procs = hypre_TReAlloc(procs, int, max_nentries);
|
|---|
| 325 |
|
|---|
| 326 |
|
|---|
| 327 |
|
|---|
| 328 | /* update manager */
|
|---|
| 329 | hypre_BoxManMaxNEntries(manager) = max_nentries;
|
|---|
| 330 | hypre_BoxManEntries(manager) = entries;
|
|---|
| 331 | hypre_BoxManIdsSort(manager) = ids;
|
|---|
| 332 | hypre_BoxManProcsSort(manager) = procs;
|
|---|
| 333 |
|
|---|
| 334 |
|
|---|
| 335 |
|
|---|
| 336 | #ifdef HYPRE_NO_GLOBAL_PARTITION
|
|---|
| 337 | {
|
|---|
| 338 | int *my_ids = hypre_BoxManMyIds(manager);
|
|---|
| 339 | hypre_BoxManEntry **my_entries = hypre_BoxManMyEntries(manager);
|
|---|
| 340 |
|
|---|
| 341 | my_ids = hypre_TReAlloc(my_ids, int, max_nentries);
|
|---|
| 342 |
|
|---|
| 343 | my_entries = hypre_TReAlloc(my_entries, hypre_BoxManEntry *, max_nentries);
|
|---|
| 344 |
|
|---|
| 345 | hypre_BoxManMyIds(manager) = my_ids;
|
|---|
| 346 | hypre_BoxManMyEntries(manager) = my_entries;
|
|---|
| 347 | }
|
|---|
| 348 |
|
|---|
| 349 | #endif
|
|---|
| 350 |
|
|---|
| 351 |
|
|---|
| 352 | return hypre_error_flag;
|
|---|
| 353 |
|
|---|
| 354 | }
|
|---|
| 355 |
|
|---|
| 356 |
|
|---|
| 357 |
|
|---|
| 358 | /*--------------------------------------------------------------------------
|
|---|
| 359 | hypre_BoxManDestroy:
|
|---|
| 360 |
|
|---|
| 361 | De-allocate the box manager structure.
|
|---|
| 362 |
|
|---|
| 363 | *--------------------------------------------------------------------------*/
|
|---|
| 364 |
|
|---|
| 365 | int hypre_BoxManDestroy ( hypre_BoxManager *manager )
|
|---|
| 366 |
|
|---|
| 367 | {
|
|---|
| 368 | int d;
|
|---|
| 369 |
|
|---|
| 370 | if (manager)
|
|---|
| 371 | {
|
|---|
| 372 |
|
|---|
| 373 | for (d = 0; d < 3; d++)
|
|---|
| 374 | {
|
|---|
| 375 | hypre_TFree(hypre_BoxManIndexesD(manager, d));
|
|---|
| 376 | }
|
|---|
| 377 |
|
|---|
| 378 | hypre_TFree(hypre_BoxManEntries(manager));
|
|---|
| 379 | hypre_TFree(hypre_BoxManIndexTable(manager));
|
|---|
| 380 | hypre_TFree(hypre_BoxManSortTable(manager));
|
|---|
| 381 |
|
|---|
| 382 | hypre_TFree(hypre_BoxManIdsSort(manager));
|
|---|
| 383 | hypre_TFree(hypre_BoxManProcsSort(manager));
|
|---|
| 384 | hypre_TFree(hypre_BoxManProcsSortOffsets(manager));
|
|---|
| 385 |
|
|---|
| 386 | hypre_BoxArrayDestroy(hypre_BoxManGatherRegions(manager));
|
|---|
| 387 |
|
|---|
| 388 | hypre_TFree(hypre_BoxManMyIds(manager));
|
|---|
| 389 | hypre_TFree(hypre_BoxManMyEntries(manager));
|
|---|
| 390 |
|
|---|
| 391 | hypre_StructAssumedPartitionDestroy(hypre_BoxManAssumedPartition(manager));
|
|---|
| 392 |
|
|---|
| 393 | hypre_TFree(manager);
|
|---|
| 394 | }
|
|---|
| 395 |
|
|---|
| 396 |
|
|---|
| 397 | return hypre_error_flag;
|
|---|
| 398 |
|
|---|
| 399 | }
|
|---|
| 400 |
|
|---|
| 401 |
|
|---|
| 402 | /*--------------------------------------------------------------------------
|
|---|
| 403 | hypre_BoxManAddEntry:
|
|---|
| 404 |
|
|---|
| 405 | Add a box (entry) to the box manager. Each entry is given a
|
|---|
| 406 | unique id (proc_id, box_id). Need to assemble after adding entries.
|
|---|
| 407 |
|
|---|
| 408 | Notes:
|
|---|
| 409 |
|
|---|
| 410 | (1) The id assigned may be any integer - though since (proc_id,
|
|---|
| 411 | box_id) is unique, duplicates will be eliminated in the assemble.
|
|---|
| 412 |
|
|---|
| 413 |
|
|---|
| 414 | (2) if there is not enough storage available for this entry, then
|
|---|
| 415 | increase the amount automatically
|
|---|
| 416 |
|
|---|
| 417 |
|
|---|
| 418 | *--------------------------------------------------------------------------*/
|
|---|
| 419 |
|
|---|
| 420 | int hypre_BoxManAddEntry( hypre_BoxManager *manager , hypre_Index imin ,
|
|---|
| 421 | hypre_Index imax , int proc_id, int box_id,
|
|---|
| 422 | void *info )
|
|---|
| 423 |
|
|---|
| 424 | {
|
|---|
| 425 | int myid;
|
|---|
| 426 | int nentries = hypre_BoxManNEntries(manager);
|
|---|
| 427 | hypre_BoxManEntry *entries = hypre_BoxManEntries(manager);
|
|---|
| 428 | hypre_BoxManEntry *entry;
|
|---|
| 429 | hypre_IndexRef entry_imin;
|
|---|
| 430 | hypre_IndexRef entry_imax;
|
|---|
| 431 | int d;
|
|---|
| 432 | int *num_ghost = hypre_BoxManNumGhost(manager);
|
|---|
| 433 |
|
|---|
| 434 |
|
|---|
| 435 | MPI_Comm_rank(hypre_BoxManComm(manager), &myid );
|
|---|
| 436 |
|
|---|
| 437 | /* check to make sure that there is enough storage available
|
|---|
| 438 | for this new entry - if not add space for 5 more*/
|
|---|
| 439 |
|
|---|
| 440 | if (nentries + 1 > hypre_BoxManMaxNEntries(manager))
|
|---|
| 441 | {
|
|---|
| 442 | hypre_BoxManIncSize( manager, 5);
|
|---|
| 443 | }
|
|---|
| 444 |
|
|---|
| 445 | /* we add this to the end entry list - get pointer to location*/
|
|---|
| 446 | entry = &entries[nentries];
|
|---|
| 447 | entry_imin = hypre_BoxManEntryIMin(entry);
|
|---|
| 448 | entry_imax = hypre_BoxManEntryIMax(entry);
|
|---|
| 449 |
|
|---|
| 450 | /* copy information into entry */
|
|---|
| 451 | for (d = 0; d < 3; d++)
|
|---|
| 452 | {
|
|---|
| 453 | hypre_IndexD(entry_imin, d) = hypre_IndexD(imin, d);
|
|---|
| 454 | hypre_IndexD(entry_imax, d) = hypre_IndexD(imax, d);
|
|---|
| 455 | }
|
|---|
| 456 |
|
|---|
| 457 | hypre_BoxManEntryProc(entry) = proc_id;
|
|---|
| 458 | hypre_BoxManEntryId(entry) = box_id;
|
|---|
| 459 | hypre_BoxManEntryInfo(entry) = info;
|
|---|
| 460 |
|
|---|
| 461 |
|
|---|
| 462 | /* inherit and inject the numghost from mananger into the entry (as in boxmap) */
|
|---|
| 463 | for (d = 0; d < 6; d++)
|
|---|
| 464 | {
|
|---|
| 465 | hypre_BoxManEntryNumGhost(entry)[d] = num_ghost[d];
|
|---|
| 466 | }
|
|---|
| 467 | hypre_BoxManEntryNext(entry)= NULL;
|
|---|
| 468 |
|
|---|
| 469 | /* add proc and id to procs_sort and ids_sort array */
|
|---|
| 470 | hypre_BoxManProcsSort(manager)[nentries] = proc_id;
|
|---|
| 471 | hypre_BoxManIdsSort(manager)[nentries] = box_id;
|
|---|
| 472 |
|
|---|
| 473 |
|
|---|
| 474 | #ifdef HYPRE_NO_GLOBAL_PARTITION
|
|---|
| 475 | /* here we need to keep track of my entries separately just to improve
|
|---|
| 476 | speed at the beginning of the assemble - then this gets deleted when
|
|---|
| 477 | the sort_table is created. */
|
|---|
| 478 |
|
|---|
| 479 | if (proc_id == myid)
|
|---|
| 480 | {
|
|---|
| 481 | int *my_ids = hypre_BoxManMyIds(manager);
|
|---|
| 482 | hypre_BoxManEntry **my_entries = hypre_BoxManMyEntries(manager);
|
|---|
| 483 | int num_my_entries = hypre_BoxManNumMyEntries(manager);
|
|---|
| 484 |
|
|---|
| 485 | my_ids[num_my_entries] = box_id;
|
|---|
| 486 | my_entries[num_my_entries] = &entries[nentries];
|
|---|
| 487 | num_my_entries++;
|
|---|
| 488 |
|
|---|
| 489 | hypre_BoxManNumMyEntries(manager) = num_my_entries;
|
|---|
| 490 |
|
|---|
| 491 |
|
|---|
| 492 | }
|
|---|
| 493 |
|
|---|
| 494 | #endif
|
|---|
| 495 |
|
|---|
| 496 | /* increment number of entries */
|
|---|
| 497 | hypre_BoxManNEntries(manager) = nentries + 1;
|
|---|
| 498 |
|
|---|
| 499 |
|
|---|
| 500 | return hypre_error_flag;
|
|---|
| 501 |
|
|---|
| 502 | }
|
|---|
| 503 |
|
|---|
| 504 |
|
|---|
| 505 | /*--------------------------------------------------------------------------
|
|---|
| 506 | hypre_BoxManGetEntry:
|
|---|
| 507 |
|
|---|
| 508 | Given an id: (proc_id, box_id), return a pointer to the box entry.
|
|---|
| 509 |
|
|---|
| 510 | Notes:
|
|---|
| 511 |
|
|---|
| 512 | (1) Use of this is generally to get back something that has been
|
|---|
| 513 | added by the above function. If no entry is found, an error is returned.
|
|---|
| 514 |
|
|---|
| 515 | (2) This functionality will replace that previously provided by
|
|---|
| 516 | hypre_BoxManFindBoxProcEntry.
|
|---|
| 517 |
|
|---|
| 518 | (3) Need to store entry information such that this information is
|
|---|
| 519 | easily found. (During the assemble, we will sort on proc_id, then
|
|---|
| 520 | box_id, and provide a pointer to the entries. Then we can do a
|
|---|
| 521 | search into the proc_id, and then into the box_id.)
|
|---|
| 522 |
|
|---|
| 523 |
|
|---|
| 524 | *--------------------------------------------------------------------------*/
|
|---|
| 525 |
|
|---|
| 526 | int hypre_BoxManGetEntry( hypre_BoxManager *manager , int proc, int id,
|
|---|
| 527 | hypre_BoxManEntry **entry_ptr )
|
|---|
| 528 |
|
|---|
| 529 | {
|
|---|
| 530 |
|
|---|
| 531 |
|
|---|
| 532 | /* find proc_id in procs array. then find id in ids array, then grab
|
|---|
| 533 | the corresponding entry */
|
|---|
| 534 |
|
|---|
| 535 | hypre_BoxManEntry *entry;
|
|---|
| 536 |
|
|---|
| 537 | int myid;
|
|---|
| 538 | int i, offset;
|
|---|
| 539 | int start, finish;
|
|---|
| 540 | int location;
|
|---|
| 541 | int first_local = hypre_BoxManFirstLocal(manager);
|
|---|
| 542 | int *procs_sort = hypre_BoxManProcsSort(manager);
|
|---|
| 543 | int *ids_sort = hypre_BoxManIdsSort(manager);
|
|---|
| 544 | int nentries = hypre_BoxManNEntries(manager);
|
|---|
| 545 | int num_proc = hypre_BoxManNumProcsSort(manager);
|
|---|
| 546 | int *proc_offsets = hypre_BoxManProcsSortOffsets(manager);
|
|---|
| 547 |
|
|---|
| 548 | MPI_Comm_rank(hypre_BoxManComm(manager), &myid );
|
|---|
| 549 |
|
|---|
| 550 | if (nentries)
|
|---|
| 551 | {
|
|---|
| 552 |
|
|---|
| 553 |
|
|---|
| 554 | if ( hypre_BoxManSortTable(manager) == NULL)
|
|---|
| 555 | hypre_error_in_arg(1);
|
|---|
| 556 |
|
|---|
| 557 |
|
|---|
| 558 | /* check to see if it is the local id first - this will be the
|
|---|
| 559 | * case most of the time (currently it is only used in this
|
|---|
| 560 | manner)*/
|
|---|
| 561 | if (proc == myid)
|
|---|
| 562 | {
|
|---|
| 563 | start = first_local;
|
|---|
| 564 | if (start >= 0 )
|
|---|
| 565 | {
|
|---|
| 566 | finish = proc_offsets[hypre_BoxManLocalProcOffset(manager)+1];
|
|---|
| 567 | }
|
|---|
| 568 | }
|
|---|
| 569 |
|
|---|
| 570 | else /* otherwise find proc (TO DO: just have procs_sort not
|
|---|
| 571 | contain duplicates - then we could do a binary search and
|
|---|
| 572 | save on some storage - this has to be changed in assemble,
|
|---|
| 573 | then also memory management in addentry - but currently this
|
|---|
| 574 | because proc = myid)*/
|
|---|
| 575 | {
|
|---|
| 576 | start = -1;
|
|---|
| 577 | for (i = 0; i< num_proc; i++)
|
|---|
| 578 | {
|
|---|
| 579 | offset = proc_offsets[i];
|
|---|
| 580 | if (proc == procs_sort[offset])
|
|---|
| 581 | {
|
|---|
| 582 | start = offset;
|
|---|
| 583 | finish = proc_offsets[i+1];
|
|---|
| 584 | break;
|
|---|
| 585 | }
|
|---|
| 586 | }
|
|---|
| 587 | }
|
|---|
| 588 | if (start >= 0 )
|
|---|
| 589 | {
|
|---|
| 590 | /* now look for the id - returns -1 if not found*/
|
|---|
| 591 | location = hypre_BinarySearch(&ids_sort[start], id, finish-start);
|
|---|
| 592 | }
|
|---|
| 593 | else
|
|---|
| 594 | {
|
|---|
| 595 | location = -1;
|
|---|
| 596 | }
|
|---|
| 597 |
|
|---|
| 598 | }
|
|---|
| 599 | else
|
|---|
| 600 | {
|
|---|
| 601 | location = -1;
|
|---|
| 602 | }
|
|---|
| 603 |
|
|---|
| 604 |
|
|---|
| 605 | if (location >= 0 )
|
|---|
| 606 | {
|
|---|
| 607 | /* this location is relative to where we started searching - so fix if non-negative */
|
|---|
| 608 | location += start;
|
|---|
| 609 | /* now grab entry */
|
|---|
| 610 | entry = hypre_BoxManSortTable(manager)[location];
|
|---|
| 611 | }
|
|---|
| 612 | else
|
|---|
| 613 | entry = NULL;
|
|---|
| 614 |
|
|---|
| 615 | *entry_ptr = entry;
|
|---|
| 616 |
|
|---|
| 617 | return hypre_error_flag;
|
|---|
| 618 |
|
|---|
| 619 | }
|
|---|
| 620 |
|
|---|
| 621 |
|
|---|
| 622 | /*--------------------------------------------------------------------------
|
|---|
| 623 | hypre_BoxManGetAllEntries
|
|---|
| 624 |
|
|---|
| 625 | Return a list of all of the entries in the box manager (and the
|
|---|
| 626 | number of entries). These are not sorted in any particular way.
|
|---|
| 627 |
|
|---|
| 628 | We do not just return a pointer to the entries array in the box manager as
|
|---|
| 629 | we don't want the user to "mess with" the entries array - e.g., delete it :)
|
|---|
| 630 |
|
|---|
| 631 | *--------------------------------------------------------------------------*/
|
|---|
| 632 |
|
|---|
| 633 |
|
|---|
| 634 | int hypre_BoxManGetAllEntries( hypre_BoxManager *manager , int *num_entries,
|
|---|
| 635 | hypre_BoxManEntry ***entries_ptr)
|
|---|
| 636 |
|
|---|
| 637 | {
|
|---|
| 638 |
|
|---|
| 639 |
|
|---|
| 640 | hypre_BoxManEntry **entries;
|
|---|
| 641 | int i, nentries;
|
|---|
| 642 |
|
|---|
| 643 | hypre_BoxManEntry *boxman_entries = hypre_BoxManEntries(manager);
|
|---|
| 644 |
|
|---|
| 645 |
|
|---|
| 646 | nentries = hypre_BoxManNEntries(manager);
|
|---|
| 647 | entries = hypre_CTAlloc(hypre_BoxManEntry *, nentries);
|
|---|
| 648 | for (i= 0; i< nentries; i++)
|
|---|
| 649 | {
|
|---|
| 650 | entries[i]= &boxman_entries[i];;
|
|---|
| 651 | }
|
|---|
| 652 |
|
|---|
| 653 |
|
|---|
| 654 | /* return */
|
|---|
| 655 | *entries_ptr = entries;
|
|---|
| 656 | *num_entries = nentries;
|
|---|
| 657 |
|
|---|
| 658 |
|
|---|
| 659 | return hypre_error_flag;
|
|---|
| 660 |
|
|---|
| 661 | }
|
|---|
| 662 |
|
|---|
| 663 | /*--------------------------------------------------------------------------
|
|---|
| 664 | hypre_BoxManGatherEntries:
|
|---|
| 665 |
|
|---|
| 666 | All global entries that lie within the boxes supplied to this
|
|---|
| 667 | function are gathered from other processors during the assemble and
|
|---|
| 668 | stored in a processor's local box manager. Multiple calls may be
|
|---|
| 669 | made to this function. The box extents supplied here are not retained
|
|---|
| 670 | after the assemble.
|
|---|
| 671 |
|
|---|
| 672 |
|
|---|
| 673 | Note:
|
|---|
| 674 |
|
|---|
| 675 | (1) This affects whether or not calls to BoxManIntersect() can be
|
|---|
| 676 | answered correctly. In other words, the user needs to anticipate the
|
|---|
| 677 | areas of the grid where BoxManIntersect() calls will be made, and
|
|---|
| 678 | make sure that information has been collected.
|
|---|
| 679 |
|
|---|
| 680 | (2) when this is called, the boolean "is_gather_entries" is set and
|
|---|
| 681 | the box is added to gather_regions array.
|
|---|
| 682 |
|
|---|
| 683 | *--------------------------------------------------------------------------*/
|
|---|
| 684 |
|
|---|
| 685 | int hypre_BoxManGatherEntries(hypre_BoxManager *manager , hypre_Index imin ,
|
|---|
| 686 | hypre_Index imax )
|
|---|
| 687 | {
|
|---|
| 688 |
|
|---|
| 689 | hypre_Box *box;
|
|---|
| 690 |
|
|---|
| 691 | hypre_BoxArray *gather_regions;
|
|---|
| 692 |
|
|---|
| 693 |
|
|---|
| 694 | /* initialize */
|
|---|
| 695 | hypre_BoxManIsGatherCalled(manager) = 1;
|
|---|
| 696 | gather_regions = hypre_BoxManGatherRegions(manager);
|
|---|
| 697 |
|
|---|
| 698 |
|
|---|
| 699 | /* add the box to the gather region array */
|
|---|
| 700 | box = hypre_BoxCreate();
|
|---|
| 701 | hypre_BoxSetExtents( box, imin, imax );
|
|---|
| 702 | hypre_AppendBox( box, gather_regions); /* this is a copy */
|
|---|
| 703 |
|
|---|
| 704 |
|
|---|
| 705 | /* clean up */
|
|---|
| 706 | hypre_BoxDestroy(box);
|
|---|
| 707 | hypre_BoxManGatherRegions(manager) = gather_regions; /* may have been
|
|---|
| 708 | a re-alloc */
|
|---|
| 709 |
|
|---|
| 710 | return hypre_error_flag;
|
|---|
| 711 |
|
|---|
| 712 | }
|
|---|
| 713 |
|
|---|
| 714 |
|
|---|
| 715 | /*--------------------------------------------------------------------------
|
|---|
| 716 | hypre_BoxManAssemble:
|
|---|
| 717 |
|
|---|
| 718 | In the assemble, we populate the local box manager with global box
|
|---|
| 719 | information to be used by calls to BoxManIntersect(). Global box
|
|---|
| 720 | information is gathered that corresponds to the regions input by calls
|
|---|
| 721 | to hypre_BoxManGatherEntries().
|
|---|
| 722 |
|
|---|
| 723 | Notes:
|
|---|
| 724 |
|
|---|
| 725 |
|
|---|
| 726 | (1) In the assumed partition (AP) case, the boxes gathered are those
|
|---|
| 727 | that correspond to boxes living in the assumed partition regions
|
|---|
| 728 | that intersect the regions input to hypre_BoxManGatherEntries().
|
|---|
| 729 | (We will have to check for duplicates here as a box can be in more
|
|---|
| 730 | than one AP.)
|
|---|
| 731 |
|
|---|
| 732 | (2) If hypre_BoxManGatherEntries() has *not* been called, then only
|
|---|
| 733 | the box information provided via calls to hypre_BoxManAddEntry will
|
|---|
| 734 | be in the box manager. (There is a global communication to check if
|
|---|
| 735 | GatherEntires has been called on any processor). In the non-AP
|
|---|
| 736 | case, if GatherEntries is called on *any* processor, then all
|
|---|
| 737 | processors get all boxes (via allgatherv).
|
|---|
| 738 |
|
|---|
| 739 | (3) Need to check for duplicate boxes (and eliminate) - based on
|
|---|
| 740 | pair (proc_id, box_id). Also sort this identifier pair so that
|
|---|
| 741 | GetEntry calls can be made more easily.
|
|---|
| 742 |
|
|---|
| 743 | (5) ****TO DO****Particularly in the AP case, might want to think
|
|---|
| 744 | about a "smart" algorithm to decide whether point-to-point
|
|---|
| 745 | communications or an AllGather is the best way to collect the needed
|
|---|
| 746 | entries resulting from calls to GatherEntries(). If this was done
|
|---|
| 747 | well, then the AP and non-AP would not have to be treated
|
|---|
| 748 | separately at all!
|
|---|
| 749 |
|
|---|
| 750 |
|
|---|
| 751 | **Assumptions:
|
|---|
| 752 |
|
|---|
| 753 | 1. A processor has used "add entry" to put all of the boxes
|
|---|
| 754 | that it owns into its box manager
|
|---|
| 755 |
|
|---|
| 756 | 2. The assemble routine is only called once for a box manager (i.e., you
|
|---|
| 757 | don't assemble, then add more entries and then assemble again)
|
|---|
| 758 |
|
|---|
| 759 | *--------------------------------------------------------------------------*/
|
|---|
| 760 |
|
|---|
| 761 | int hypre_BoxManAssemble ( hypre_BoxManager *manager)
|
|---|
| 762 |
|
|---|
| 763 | {
|
|---|
| 764 |
|
|---|
| 765 | int myid, nprocs;
|
|---|
| 766 | int is_gather, global_is_gather;
|
|---|
| 767 | int nentries;
|
|---|
| 768 | int *procs_sort, *ids_sort;
|
|---|
| 769 | int i,j, k;
|
|---|
| 770 |
|
|---|
| 771 | hypre_BoxManEntry *entries;
|
|---|
| 772 |
|
|---|
| 773 | hypre_BoxArray *gather_regions;
|
|---|
| 774 |
|
|---|
| 775 | MPI_Comm comm = hypre_BoxManComm(manager);
|
|---|
| 776 |
|
|---|
| 777 |
|
|---|
| 778 | /* initilize */
|
|---|
| 779 | MPI_Comm_rank(comm, &myid);
|
|---|
| 780 | MPI_Comm_size(comm, &nprocs);
|
|---|
| 781 |
|
|---|
| 782 | gather_regions = hypre_BoxManGatherRegions(manager);
|
|---|
| 783 | nentries = hypre_BoxManNEntries(manager);
|
|---|
| 784 | entries = hypre_BoxManEntries(manager);
|
|---|
| 785 | procs_sort = hypre_BoxManProcsSort(manager);
|
|---|
| 786 |
|
|---|
| 787 | ids_sort = hypre_BoxManIdsSort(manager);
|
|---|
| 788 |
|
|---|
| 789 |
|
|---|
| 790 | /* do we need to gather entries -check to see if ANY processor
|
|---|
| 791 | called a gather?*/
|
|---|
| 792 | is_gather = hypre_BoxManIsGatherCalled(manager);
|
|---|
| 793 | MPI_Allreduce(&is_gather, &global_is_gather, 1, MPI_INT, MPI_LOR, comm);
|
|---|
| 794 |
|
|---|
| 795 | /* ----------------------------GATHER? ------------------------------------*/
|
|---|
| 796 |
|
|---|
| 797 | if (global_is_gather)
|
|---|
| 798 | {
|
|---|
| 799 |
|
|---|
| 800 | /* if AP, use AP to find out who owns the data we need. In the
|
|---|
| 801 | non-AP, then just gather everything. */
|
|---|
| 802 |
|
|---|
| 803 | /* Goal: to gather the entries from the relevant processor and
|
|---|
| 804 | * add to the *entries array. Also add the proc and id to the
|
|---|
| 805 | * procs_sort and ids_sort arrays */
|
|---|
| 806 |
|
|---|
| 807 | #ifdef HYPRE_NO_GLOBAL_PARTITION
|
|---|
| 808 | {
|
|---|
| 809 | int size, index;
|
|---|
| 810 | int *tmp_proc_ids;
|
|---|
| 811 | int proc_count, proc_alloc;
|
|---|
| 812 | int *proc_array;
|
|---|
| 813 | int *ap_proc_ids;
|
|---|
| 814 | int count;
|
|---|
| 815 | int global_num_boxes;
|
|---|
| 816 | int max_response_size;
|
|---|
| 817 | int non_info_size, entry_size_bytes;
|
|---|
| 818 | int *neighbor_proc_ids = NULL, *neighbor_boxnums = NULL;
|
|---|
| 819 | int *response_buf_starts;
|
|---|
| 820 | int *response_buf;
|
|---|
| 821 | int response_size, tmp_int;
|
|---|
| 822 |
|
|---|
| 823 | int *send_buf = NULL;
|
|---|
| 824 | int *send_buf_starts = NULL;
|
|---|
| 825 | int *order_index = NULL;
|
|---|
| 826 | int *delete_array = NULL;
|
|---|
| 827 | int tmp_id, start, d, proc, id;
|
|---|
| 828 | int *tmp_int_ptr;
|
|---|
| 829 | int *contact_proc_ids = NULL;
|
|---|
| 830 |
|
|---|
| 831 | int max_regions, max_refinements, ologp;
|
|---|
| 832 |
|
|---|
| 833 | int grow, grow_array[6];
|
|---|
| 834 | int sendbuf6[6], recvbuf6[6];
|
|---|
| 835 |
|
|---|
| 836 | int *local_boxnums;
|
|---|
| 837 |
|
|---|
| 838 | int dim = hypre_BoxManDim(manager);
|
|---|
| 839 |
|
|---|
| 840 | int *my_ids = hypre_BoxManMyIds(manager);
|
|---|
| 841 | int num_my_entries = hypre_BoxManNumMyEntries(manager);
|
|---|
| 842 |
|
|---|
| 843 | void *entry_response_buf;
|
|---|
| 844 | void *index_ptr;
|
|---|
| 845 |
|
|---|
| 846 | double gamma;
|
|---|
| 847 | double local_volume, global_volume;
|
|---|
| 848 | double sendbuf2[2], recvbuf2[2];
|
|---|
| 849 |
|
|---|
| 850 | hypre_BoxArray *gather_regions;
|
|---|
| 851 | hypre_BoxArray *local_boxes;
|
|---|
| 852 |
|
|---|
| 853 | hypre_Box *box, *grow_box, *bounding_box;
|
|---|
| 854 |
|
|---|
| 855 | hypre_StructAssumedPart *ap;
|
|---|
| 856 |
|
|---|
| 857 | hypre_DataExchangeResponse response_obj, response_obj2;
|
|---|
| 858 |
|
|---|
| 859 | hypre_BoxManEntry **my_entries = hypre_BoxManMyEntries(manager) ;
|
|---|
| 860 | hypre_BoxManEntry *entry_ptr;
|
|---|
| 861 |
|
|---|
| 862 | hypre_Index imin, imax;
|
|---|
| 863 |
|
|---|
| 864 | hypre_IndexRef min_ref, max_ref;
|
|---|
| 865 |
|
|---|
| 866 | /* 1. Need to create an assumed partition */
|
|---|
| 867 |
|
|---|
| 868 | /* create an array of local boxes. also get the bounding box
|
|---|
| 869 | and the global box size (as a double). we don't use the bounding box
|
|---|
| 870 | from the grid since we care about the zero volume boxes */
|
|---|
| 871 | local_boxes = hypre_BoxArrayCreate(num_my_entries);
|
|---|
| 872 | local_boxnums = hypre_CTAlloc(int, num_my_entries);
|
|---|
| 873 | bounding_box = hypre_BoxCreate();
|
|---|
| 874 | grow_box = hypre_BoxCreate();
|
|---|
| 875 | local_volume = 0.0;
|
|---|
| 876 | for (d=0; d<3; d++) /* initialize */
|
|---|
| 877 | {
|
|---|
| 878 | hypre_IndexD(imin, d) = pow(2,30) ;
|
|---|
| 879 | hypre_IndexD(imax, d) = -pow(2,30);
|
|---|
| 880 | }
|
|---|
| 881 |
|
|---|
| 882 | for (i=0; i< num_my_entries; i++)
|
|---|
| 883 | {
|
|---|
| 884 |
|
|---|
| 885 | entry_ptr = my_entries[i];
|
|---|
| 886 | min_ref = hypre_BoxManEntryIMin(entry_ptr);
|
|---|
| 887 | max_ref = hypre_BoxManEntryIMax(entry_ptr);
|
|---|
| 888 | box = hypre_BoxArrayBox(local_boxes, i);
|
|---|
| 889 | hypre_BoxSetExtents( box, min_ref, max_ref );
|
|---|
| 890 |
|
|---|
| 891 | local_boxnums[i] = hypre_BoxManEntryId(entry_ptr);
|
|---|
| 892 |
|
|---|
| 893 | local_volume += (double) hypre_BoxVolume(box);
|
|---|
| 894 |
|
|---|
| 895 | /* check for zero volume boxes */
|
|---|
| 896 |
|
|---|
| 897 | if (hypre_BoxVolume(box) == 0) /* zero volume boxes - still count */
|
|---|
| 898 | {
|
|---|
| 899 | hypre_CopyBox(box, grow_box);
|
|---|
| 900 | for (d = 0; d < 3; d++)
|
|---|
| 901 | {
|
|---|
| 902 | if(!hypre_BoxSizeD(box, d))
|
|---|
| 903 | {
|
|---|
| 904 | grow = (hypre_BoxIMinD(box, d) - hypre_BoxIMaxD(box, d) + 1)/2;
|
|---|
| 905 | grow_array[2*d] = grow;
|
|---|
| 906 | grow_array[2*d+1] = grow;
|
|---|
| 907 | }
|
|---|
| 908 | else
|
|---|
| 909 | {
|
|---|
| 910 | grow_array[2*d] = 0;
|
|---|
| 911 | grow_array[2*d+1] = 0;
|
|---|
| 912 | }
|
|---|
| 913 | }
|
|---|
| 914 | /* expand the box */
|
|---|
| 915 | hypre_BoxExpand(grow_box, grow_array);
|
|---|
| 916 | box = grow_box; /*pointer copy*/
|
|---|
| 917 | }
|
|---|
| 918 | /*now we have a vol > 0 box - look for bounding box*/
|
|---|
| 919 |
|
|---|
| 920 | for (d = 0; d < 3; d++) /* for each dimension */
|
|---|
| 921 | {
|
|---|
| 922 | hypre_IndexD(imin, d) = hypre_min( hypre_IndexD(imin, d),
|
|---|
| 923 | hypre_BoxIMinD(box, d));
|
|---|
| 924 | hypre_IndexD(imax, d) = hypre_max( hypre_IndexD(imax, d),
|
|---|
| 925 | hypre_BoxIMaxD(box, d));
|
|---|
| 926 | }
|
|---|
| 927 |
|
|---|
| 928 |
|
|---|
| 929 |
|
|---|
| 930 | }/* end of local boxes */
|
|---|
| 931 | hypre_BoxSetExtents(bounding_box, imin, imax);
|
|---|
| 932 | /* if dim < 3 then set the extra dimensions to zero */
|
|---|
| 933 | for (d = dim; d < 3; d++)
|
|---|
| 934 | {
|
|---|
| 935 | hypre_BoxIMinD(bounding_box, d) = 0;
|
|---|
| 936 | hypre_BoxIMaxD(bounding_box, d) = 0;
|
|---|
| 937 | }
|
|---|
| 938 |
|
|---|
| 939 | /* now get the bounding box */
|
|---|
| 940 | for (d = 0; d < 3; d++)
|
|---|
| 941 | {
|
|---|
| 942 | sendbuf6[d] = hypre_BoxIMinD(bounding_box, d);
|
|---|
| 943 | sendbuf6[d+3] = -hypre_BoxIMaxD(bounding_box, d);
|
|---|
| 944 | }
|
|---|
| 945 | MPI_Allreduce(sendbuf6, recvbuf6, 6, MPI_INT, MPI_MIN, comm);
|
|---|
| 946 |
|
|---|
| 947 | for (d = 0; d < 3; d++)
|
|---|
| 948 | {
|
|---|
| 949 | hypre_BoxIMinD(bounding_box, d) = recvbuf6[d];
|
|---|
| 950 | hypre_BoxIMaxD(bounding_box, d) = -recvbuf6[d+3];
|
|---|
| 951 | }
|
|---|
| 952 |
|
|---|
| 953 | if (hypre_BoxVolume(bounding_box) == 0)
|
|---|
| 954 | {
|
|---|
| 955 | /* this shouldn't happen! */
|
|---|
| 956 | hypre_error(HYPRE_ERROR_GENERIC);
|
|---|
| 957 | return hypre_error_flag;
|
|---|
| 958 | }
|
|---|
| 959 |
|
|---|
| 960 | /* get the number of global entries and the global volume */
|
|---|
| 961 |
|
|---|
| 962 | sendbuf2[0] = local_volume;
|
|---|
| 963 | sendbuf2[1] = (double) num_my_entries;
|
|---|
| 964 |
|
|---|
| 965 | MPI_Allreduce(&sendbuf2, &recvbuf2, 2, MPI_DOUBLE, MPI_SUM, comm);
|
|---|
| 966 |
|
|---|
| 967 | global_volume = recvbuf2[0];
|
|---|
| 968 | global_num_boxes = (int) recvbuf2[1];
|
|---|
| 969 |
|
|---|
| 970 | /* estimates for the assumed partition */
|
|---|
| 971 | d = nprocs/2;
|
|---|
| 972 | ologp = 0;
|
|---|
| 973 | while ( d > 0)
|
|---|
| 974 | {
|
|---|
| 975 | d = d/2; /* note - d is an int - so this is floored */
|
|---|
| 976 | ologp++;
|
|---|
| 977 | }
|
|---|
| 978 |
|
|---|
| 979 | max_regions = hypre_min(pow(2, ologp+1), 10*ologp);
|
|---|
| 980 | max_refinements = ologp;
|
|---|
| 981 | gamma = .6; /* percentage a region must be full to
|
|---|
| 982 | avoid refinement */
|
|---|
| 983 |
|
|---|
| 984 |
|
|---|
| 985 | hypre_StructAssumedPartitionCreate(dim, bounding_box,
|
|---|
| 986 | global_volume,
|
|---|
| 987 | global_num_boxes,
|
|---|
| 988 | local_boxes, local_boxnums,
|
|---|
| 989 | max_regions, max_refinements,
|
|---|
| 990 | gamma,
|
|---|
| 991 | comm,
|
|---|
| 992 | &ap);
|
|---|
| 993 |
|
|---|
| 994 |
|
|---|
| 995 |
|
|---|
| 996 |
|
|---|
| 997 | hypre_BoxManAssumedPartition(manager) = ap;
|
|---|
| 998 |
|
|---|
| 999 | hypre_BoxArrayDestroy(local_boxes);
|
|---|
| 1000 | hypre_TFree(local_boxnums);
|
|---|
| 1001 | hypre_BoxDestroy(bounding_box);
|
|---|
| 1002 | hypre_BoxDestroy(grow_box);
|
|---|
| 1003 |
|
|---|
| 1004 |
|
|---|
| 1005 |
|
|---|
| 1006 | /* 2. Need to be able to find our own entry, given the box
|
|---|
| 1007 | number - for the second data exchange do some sorting.
|
|---|
| 1008 | Then we can use my_ids to quickly find an entry. This will be
|
|---|
| 1009 | freed when the sort table is created (it's redundant).
|
|---|
| 1010 | (Note: if we are creating the AP here, then this sorting
|
|---|
| 1011 | may need to be done at the beginning of this function)*/
|
|---|
| 1012 |
|
|---|
| 1013 | hypre_entryqsort2(my_ids, my_entries, 0, num_my_entries - 1);
|
|---|
| 1014 |
|
|---|
| 1015 | /* go thru gather regions and find out which processor's AP region
|
|---|
| 1016 | they intersect */
|
|---|
| 1017 | gather_regions = hypre_BoxManGatherRegions(manager);
|
|---|
| 1018 |
|
|---|
| 1019 | /*allocate space to store info from one box */
|
|---|
| 1020 | proc_count = 0;
|
|---|
| 1021 | proc_alloc = 8;
|
|---|
| 1022 | proc_array = hypre_CTAlloc(int, proc_alloc);
|
|---|
| 1023 |
|
|---|
| 1024 |
|
|---|
| 1025 | /* probably there will mostly be one proc per box -allocate space for 2*/
|
|---|
| 1026 | size = 2*hypre_BoxArraySize(gather_regions);
|
|---|
| 1027 | tmp_proc_ids = hypre_CTAlloc(int, size);
|
|---|
| 1028 | count = 0;
|
|---|
| 1029 |
|
|---|
| 1030 | /* loop through all boxes */
|
|---|
| 1031 | hypre_ForBoxI(i, gather_regions)
|
|---|
| 1032 | {
|
|---|
| 1033 |
|
|---|
| 1034 | hypre_StructAssumedPartitionGetProcsFromBox(ap,
|
|---|
| 1035 | hypre_BoxArrayBox(gather_regions, i),
|
|---|
| 1036 | &proc_count,
|
|---|
| 1037 | &proc_alloc,
|
|---|
| 1038 | &proc_array);
|
|---|
| 1039 |
|
|---|
| 1040 | if ((count + proc_count) > size)
|
|---|
| 1041 | {
|
|---|
| 1042 | size = size + proc_count + 2*(hypre_BoxArraySize(gather_regions)-i);
|
|---|
| 1043 | tmp_proc_ids = hypre_TReAlloc(tmp_proc_ids, int, size);
|
|---|
| 1044 | }
|
|---|
| 1045 | for (j = 0; j< proc_count; j++)
|
|---|
| 1046 | {
|
|---|
| 1047 | tmp_proc_ids[count] = proc_array[j];
|
|---|
| 1048 | count++;
|
|---|
| 1049 | }
|
|---|
| 1050 | }
|
|---|
| 1051 |
|
|---|
| 1052 | hypre_TFree(proc_array);
|
|---|
| 1053 |
|
|---|
| 1054 | /* now get rid of redundencies in tmp_proc_ids (since a box can lie in more than one AP
|
|---|
| 1055 | - put in ap_proc_ids*/
|
|---|
| 1056 | qsort0(tmp_proc_ids, 0, count-1);
|
|---|
| 1057 | proc_count = 0;
|
|---|
| 1058 | ap_proc_ids = hypre_CTAlloc(int, count);
|
|---|
| 1059 |
|
|---|
| 1060 | if (count)
|
|---|
| 1061 | {
|
|---|
| 1062 | ap_proc_ids[0] = tmp_proc_ids[0];
|
|---|
| 1063 | proc_count++;
|
|---|
| 1064 | }
|
|---|
| 1065 | for (i = 1; i < count; i++)
|
|---|
| 1066 | {
|
|---|
| 1067 | if (tmp_proc_ids[i] != ap_proc_ids[proc_count-1])
|
|---|
| 1068 | {
|
|---|
| 1069 | ap_proc_ids[proc_count] = tmp_proc_ids[i];
|
|---|
| 1070 | proc_count++;
|
|---|
| 1071 | }
|
|---|
| 1072 | }
|
|---|
| 1073 | hypre_TFree(tmp_proc_ids);
|
|---|
| 1074 |
|
|---|
| 1075 | /* 3. now we have a sorted list with no duplicates in ap_proc_ids */
|
|---|
| 1076 | /* for each of these processor ids, we need to get the entries in their
|
|---|
| 1077 | box manager */
|
|---|
| 1078 |
|
|---|
| 1079 | /* note: we are not ever re-assembling, so we won't worry about whether we
|
|---|
| 1080 | * already have the info from some of these processors. */
|
|---|
| 1081 |
|
|---|
| 1082 |
|
|---|
| 1083 | /* need to do an exchange data*/
|
|---|
| 1084 |
|
|---|
| 1085 | /* if we simply returned the boxes in the
|
|---|
| 1086 | AP region, we will not have the entry information- in particular,
|
|---|
| 1087 | we will not have the "info" obj. So we have to get this info by doing a
|
|---|
| 1088 | second communication where we contact the actual owers of the boxes and
|
|---|
| 1089 | request the entry info...So:
|
|---|
| 1090 |
|
|---|
| 1091 | (1) exchange #1: contact the AP processor, get the proc ids and box numbers in that AP region
|
|---|
| 1092 |
|
|---|
| 1093 | (2) exchange #2: use this info to contact the owner processor and from them get the
|
|---|
| 1094 | entire entry infomation: box, info, etc.*/
|
|---|
| 1095 |
|
|---|
| 1096 |
|
|---|
| 1097 | /* exchange #1 - we send nothing, and the contacted proc
|
|---|
| 1098 | * returns all of its proc and box ids in its AP */
|
|---|
| 1099 |
|
|---|
| 1100 | /* build response object*/
|
|---|
| 1101 | response_obj.fill_response = hypre_FillResponseBoxMapAssemble1;
|
|---|
| 1102 | response_obj.data1 = ap; /* needed to fill responses*/
|
|---|
| 1103 | response_obj.data2 = NULL;
|
|---|
| 1104 |
|
|---|
| 1105 | send_buf = NULL;
|
|---|
| 1106 | send_buf_starts = hypre_CTAlloc(int, proc_count + 1);
|
|---|
| 1107 | for (i=0; i< proc_count+1; i++)
|
|---|
| 1108 | {
|
|---|
| 1109 | send_buf_starts[i] = 0;
|
|---|
| 1110 | }
|
|---|
| 1111 |
|
|---|
| 1112 | response_buf = NULL; /*this and the next are allocated in exchange data */
|
|---|
| 1113 | response_buf_starts = NULL;
|
|---|
| 1114 |
|
|---|
| 1115 | /*we expect back a pair of ints: proc id, box number for each box owned */
|
|---|
| 1116 | size = 2*sizeof(int);
|
|---|
| 1117 |
|
|---|
| 1118 | /* this needs to be the same for all processors */
|
|---|
| 1119 | max_response_size = (global_num_boxes/nprocs)*2;
|
|---|
| 1120 |
|
|---|
| 1121 |
|
|---|
| 1122 | hypre_DataExchangeList(proc_count, ap_proc_ids,
|
|---|
| 1123 | send_buf, send_buf_starts,
|
|---|
| 1124 | 0, size, &response_obj, max_response_size, 3,
|
|---|
| 1125 | comm, (void**) &response_buf, &response_buf_starts);
|
|---|
| 1126 |
|
|---|
| 1127 |
|
|---|
| 1128 | /*how many entries do we have?*/
|
|---|
| 1129 | size = response_buf_starts[proc_count];
|
|---|
| 1130 |
|
|---|
| 1131 | /* find a new list of procsessors to contact, eliminate duplicates */
|
|---|
| 1132 | neighbor_proc_ids = hypre_CTAlloc(int, size);
|
|---|
| 1133 | neighbor_boxnums = hypre_CTAlloc(int, size);
|
|---|
| 1134 | order_index = hypre_CTAlloc(int, size);
|
|---|
| 1135 |
|
|---|
| 1136 | index = 0;
|
|---|
| 1137 | for (i=0; i< size; i++) /* for each neighbor box */
|
|---|
| 1138 | {
|
|---|
| 1139 | neighbor_proc_ids[i] = response_buf[index++];
|
|---|
| 1140 | neighbor_boxnums[i] = response_buf[index++];
|
|---|
| 1141 |
|
|---|
| 1142 | order_index[i] = i;
|
|---|
| 1143 | }
|
|---|
| 1144 |
|
|---|
| 1145 | /*clean */
|
|---|
| 1146 | hypre_TFree(send_buf_starts);
|
|---|
| 1147 | hypre_TFree(ap_proc_ids);
|
|---|
| 1148 | hypre_TFree(response_buf_starts);
|
|---|
| 1149 | hypre_TFree(response_buf);
|
|---|
| 1150 |
|
|---|
| 1151 | /* look for duplicates and delete - also delete any with our proc_id */
|
|---|
| 1152 |
|
|---|
| 1153 | /*sort on proc_id - move boxnums and order_index*/
|
|---|
| 1154 | hypre_qsort3i(neighbor_proc_ids, neighbor_boxnums, order_index, 0, size-1);
|
|---|
| 1155 |
|
|---|
| 1156 | delete_array = hypre_CTAlloc(int, size);
|
|---|
| 1157 | index = 0;
|
|---|
| 1158 | proc_count = 0; /* keep track of # of distinct procs */
|
|---|
| 1159 |
|
|---|
| 1160 | /*now within each proc_id, we need to sort the boxnums and order index*/
|
|---|
| 1161 | if (size)
|
|---|
| 1162 | {
|
|---|
| 1163 | tmp_id = neighbor_proc_ids[0];
|
|---|
| 1164 | proc_count++;
|
|---|
| 1165 | }
|
|---|
| 1166 | start = 0;
|
|---|
| 1167 | for (i=1; i< size; i++)
|
|---|
| 1168 | {
|
|---|
| 1169 | if (neighbor_proc_ids[i] != tmp_id)
|
|---|
| 1170 | {
|
|---|
| 1171 | proc_count++;
|
|---|
| 1172 | /*now find duplicate boxnums - or delete all if myid = tmp_id*/
|
|---|
| 1173 | if (myid == tmp_id)
|
|---|
| 1174 | {
|
|---|
| 1175 | proc_count--;
|
|---|
| 1176 |
|
|---|
| 1177 | for (j=start; j< i; j++)
|
|---|
| 1178 | {
|
|---|
| 1179 | delete_array[index++] = j;
|
|---|
| 1180 | }
|
|---|
| 1181 | }
|
|---|
| 1182 | else
|
|---|
| 1183 | {
|
|---|
| 1184 | hypre_qsort2i(neighbor_boxnums, order_index, start, i-1);
|
|---|
| 1185 | for (j=start+1; j< i; j++)
|
|---|
| 1186 | {
|
|---|
| 1187 | if (neighbor_boxnums[j] == neighbor_boxnums[j-1])
|
|---|
| 1188 | {
|
|---|
| 1189 | delete_array[index++] = j;
|
|---|
| 1190 | }
|
|---|
| 1191 | }
|
|---|
| 1192 | }
|
|---|
| 1193 |
|
|---|
| 1194 | /* update start and tmp_id */
|
|---|
| 1195 | start = i;
|
|---|
| 1196 | tmp_id = neighbor_proc_ids[i];
|
|---|
| 1197 | }
|
|---|
| 1198 | }
|
|---|
| 1199 |
|
|---|
| 1200 | /* final sort and purge (the last group doesn't
|
|---|
| 1201 | get caught in the above loop) */
|
|---|
| 1202 | if (size)
|
|---|
| 1203 | {
|
|---|
| 1204 | /*now find duplicate boxnums - or delete all if myid = tmp_id*/
|
|---|
| 1205 | if (myid == tmp_id)
|
|---|
| 1206 | {
|
|---|
| 1207 | proc_count--;
|
|---|
| 1208 | for (j=start; j< size; j++)
|
|---|
| 1209 | {
|
|---|
| 1210 | delete_array[index++] = j;
|
|---|
| 1211 | }
|
|---|
| 1212 | }
|
|---|
| 1213 | else
|
|---|
| 1214 | {
|
|---|
| 1215 | hypre_qsort2i(neighbor_boxnums, order_index, start, size-1);
|
|---|
| 1216 | for (j=start+1; j< size; j++)
|
|---|
| 1217 | {
|
|---|
| 1218 | if (neighbor_boxnums[j] == neighbor_boxnums[j-1])
|
|---|
| 1219 | {
|
|---|
| 1220 | delete_array[index++] = j;
|
|---|
| 1221 | }
|
|---|
| 1222 | }
|
|---|
| 1223 | }
|
|---|
| 1224 |
|
|---|
| 1225 | }
|
|---|
| 1226 |
|
|---|
| 1227 | /*now delete (index) duplicates from ids and boxnums */
|
|---|
| 1228 | /* size = num of neighbor boxes
|
|---|
| 1229 | proc_count = # unique proc ids */
|
|---|
| 1230 | if (index)
|
|---|
| 1231 | {
|
|---|
| 1232 | start = delete_array[0];
|
|---|
| 1233 | j = 0;
|
|---|
| 1234 | for (i = start; (i + j) < size; i++)
|
|---|
| 1235 | {
|
|---|
| 1236 | if (j < index)
|
|---|
| 1237 | {
|
|---|
| 1238 | while ((i+j) == delete_array[j]) /* see if deleting consec. items */
|
|---|
| 1239 | {
|
|---|
| 1240 | j++; /*increase the shift*/
|
|---|
| 1241 | if (j == index) break;
|
|---|
| 1242 | }
|
|---|
| 1243 | }
|
|---|
| 1244 | if ((i+j) < size) /* if deleting the last item then no moving */
|
|---|
| 1245 | {
|
|---|
| 1246 | neighbor_boxnums[i] = neighbor_boxnums[i+j];
|
|---|
| 1247 | neighbor_proc_ids[i] = neighbor_proc_ids[i+j];
|
|---|
| 1248 | }
|
|---|
| 1249 |
|
|---|
| 1250 | }
|
|---|
| 1251 | }
|
|---|
| 1252 |
|
|---|
| 1253 | /* finished removing duplicate (proc, id) combinations */
|
|---|
| 1254 | /* new number of boxes is (size - index) */
|
|---|
| 1255 |
|
|---|
| 1256 | /* new contact: (send_buf is neighbor_boxnums)*/
|
|---|
| 1257 | size = size - index;
|
|---|
| 1258 | send_buf_starts = hypre_CTAlloc(int, proc_count + 1);
|
|---|
| 1259 | send_buf_starts[0] = 0;
|
|---|
| 1260 | contact_proc_ids = hypre_CTAlloc(int, proc_count);
|
|---|
| 1261 |
|
|---|
| 1262 |
|
|---|
| 1263 | /* don't check if we have an entry before getting more info -
|
|---|
| 1264 | eliminate duplicates later */
|
|---|
| 1265 |
|
|---|
| 1266 | /* create list of unique proc ids and starts buffer */
|
|---|
| 1267 | if (proc_count)
|
|---|
| 1268 | {
|
|---|
| 1269 | contact_proc_ids[0] = neighbor_proc_ids[0];
|
|---|
| 1270 | }
|
|---|
| 1271 |
|
|---|
| 1272 | index = 0; /* for indexing into unique contact_proc_ids */
|
|---|
| 1273 | for (i= 1; i< size; i++)
|
|---|
| 1274 | {
|
|---|
| 1275 | if (neighbor_proc_ids[i] != contact_proc_ids[index])
|
|---|
| 1276 | {
|
|---|
| 1277 | index++;
|
|---|
| 1278 | contact_proc_ids[index] = neighbor_proc_ids[i];
|
|---|
| 1279 | send_buf_starts[index] = i;
|
|---|
| 1280 | }
|
|---|
| 1281 | }
|
|---|
| 1282 | send_buf_starts[proc_count] = size;
|
|---|
| 1283 |
|
|---|
| 1284 | /* exchange #2 - now we contact processors with a box id and that processor needs
|
|---|
| 1285 | to send us the entry information*/
|
|---|
| 1286 |
|
|---|
| 1287 | entry_response_buf = NULL; /*this and the next are allocated in exchange data */
|
|---|
| 1288 | response_buf_starts = NULL;
|
|---|
| 1289 |
|
|---|
| 1290 | response_obj2.fill_response = hypre_FillResponseBoxMapAssemble2;
|
|---|
| 1291 | response_obj2.data1 = manager; /* needed to fill responses*/
|
|---|
| 1292 | response_obj2.data2 = NULL;
|
|---|
| 1293 |
|
|---|
| 1294 |
|
|---|
| 1295 | /*how big is an entry? extents: 6 ints, proc: 1 int, id: 1
|
|---|
| 1296 | * int , num_ghost: 6 ints, info: info_size is in bytes*/
|
|---|
| 1297 | /* note: for now, we do not need to send num_ghost - this
|
|---|
| 1298 | is just copied in addentry anyhow */
|
|---|
| 1299 | non_info_size = 8;
|
|---|
| 1300 | entry_size_bytes = non_info_size*sizeof(int) + hypre_BoxManEntryInfoSize(manager);
|
|---|
| 1301 |
|
|---|
| 1302 | /* use same max_response_size as previous exchange */
|
|---|
| 1303 |
|
|---|
| 1304 | hypre_DataExchangeList(proc_count, contact_proc_ids,
|
|---|
| 1305 | neighbor_boxnums, send_buf_starts, sizeof(int),
|
|---|
| 1306 | entry_size_bytes, &response_obj2, max_response_size, 4,
|
|---|
| 1307 | comm, &entry_response_buf, &response_buf_starts);
|
|---|
| 1308 |
|
|---|
| 1309 |
|
|---|
| 1310 | /* now we can add entries that are in response_buf
|
|---|
| 1311 | - we check for duplicates later */
|
|---|
| 1312 |
|
|---|
| 1313 | /*how many entries do we have?*/
|
|---|
| 1314 | response_size = response_buf_starts[proc_count];
|
|---|
| 1315 |
|
|---|
| 1316 | /* do we need more storage ?*/
|
|---|
| 1317 | if (nentries + response_size > hypre_BoxManMaxNEntries(manager))
|
|---|
| 1318 | {
|
|---|
| 1319 | int inc_size;
|
|---|
| 1320 |
|
|---|
| 1321 | inc_size = (response_size + nentries - hypre_BoxManMaxNEntries(manager));
|
|---|
| 1322 | hypre_BoxManIncSize ( manager, inc_size);
|
|---|
| 1323 |
|
|---|
| 1324 | entries = hypre_BoxManEntries(manager);
|
|---|
| 1325 | procs_sort = hypre_BoxManProcsSort(manager);
|
|---|
| 1326 | ids_sort = hypre_BoxManIdsSort(manager);
|
|---|
| 1327 | }
|
|---|
| 1328 |
|
|---|
| 1329 | index_ptr = entry_response_buf; /* point into response buf */
|
|---|
| 1330 | for (i = 0; i < response_size; i++)
|
|---|
| 1331 | {
|
|---|
| 1332 | size = sizeof(int);
|
|---|
| 1333 | /* imin */
|
|---|
| 1334 | for (d = 0; d < 3; d++)
|
|---|
| 1335 | {
|
|---|
| 1336 | memcpy( &tmp_int, index_ptr, size);
|
|---|
| 1337 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1338 | hypre_IndexD(imin, d) = tmp_int;
|
|---|
| 1339 | }
|
|---|
| 1340 |
|
|---|
| 1341 | /*imax */
|
|---|
| 1342 | for (d = 0; d < 3; d++)
|
|---|
| 1343 | {
|
|---|
| 1344 | memcpy( &tmp_int, index_ptr, size);
|
|---|
| 1345 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1346 | hypre_IndexD(imax, d) = tmp_int;
|
|---|
| 1347 | }
|
|---|
| 1348 |
|
|---|
| 1349 | /* proc */
|
|---|
| 1350 | tmp_int_ptr = (int *) index_ptr;
|
|---|
| 1351 | proc = *tmp_int_ptr;
|
|---|
| 1352 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1353 |
|
|---|
| 1354 | /* id */
|
|---|
| 1355 | tmp_int_ptr = (int *) index_ptr;
|
|---|
| 1356 | id = *tmp_int_ptr;
|
|---|
| 1357 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1358 |
|
|---|
| 1359 | /* info */
|
|---|
| 1360 | /* index_ptr is at info */
|
|---|
| 1361 |
|
|---|
| 1362 | hypre_BoxManAddEntry( manager , imin ,
|
|---|
| 1363 | imax , proc, id,
|
|---|
| 1364 | index_ptr );
|
|---|
| 1365 |
|
|---|
| 1366 | /* start of next entry */
|
|---|
| 1367 | index_ptr = (void *) ((char *) index_ptr + hypre_BoxManEntryInfoSize(manager));
|
|---|
| 1368 | }
|
|---|
| 1369 |
|
|---|
| 1370 | /* clean up from this section of code*/
|
|---|
| 1371 | hypre_TFree(entry_response_buf);
|
|---|
| 1372 | hypre_TFree(response_buf_starts);
|
|---|
| 1373 | hypre_TFree(send_buf_starts);
|
|---|
| 1374 | hypre_TFree(contact_proc_ids);
|
|---|
| 1375 | hypre_TFree(neighbor_boxnums);
|
|---|
| 1376 | hypre_TFree(neighbor_proc_ids);
|
|---|
| 1377 | hypre_TFree(order_index);
|
|---|
| 1378 | hypre_TFree(delete_array);
|
|---|
| 1379 |
|
|---|
| 1380 |
|
|---|
| 1381 | /*we don't need special access to my entries anymore - because we will
|
|---|
| 1382 | create the sort table */
|
|---|
| 1383 | {
|
|---|
| 1384 |
|
|---|
| 1385 | hypre_TFree(hypre_BoxManMyIds(manager));
|
|---|
| 1386 | hypre_TFree(hypre_BoxManMyEntries(manager));
|
|---|
| 1387 |
|
|---|
| 1388 | hypre_BoxManMyIds(manager) = NULL;
|
|---|
| 1389 | hypre_BoxManMyEntries(manager) = NULL;
|
|---|
| 1390 | }
|
|---|
| 1391 |
|
|---|
| 1392 | }
|
|---|
| 1393 |
|
|---|
| 1394 | /* end of gathering for the AP case */
|
|---|
| 1395 | #else
|
|---|
| 1396 |
|
|---|
| 1397 | /* beginning of gathering for the non-AP case */
|
|---|
| 1398 | {
|
|---|
| 1399 |
|
|---|
| 1400 | /* collect global data */
|
|---|
| 1401 | int entry_size_bytes;
|
|---|
| 1402 | int send_count, send_count_bytes;
|
|---|
| 1403 | int *displs, *recv_counts;
|
|---|
| 1404 | int recv_buf_size, recv_buf_size_bytes;
|
|---|
| 1405 | int d;
|
|---|
| 1406 | int size, non_info_size;
|
|---|
| 1407 | int proc, id;
|
|---|
| 1408 | int tmp_int;
|
|---|
| 1409 | int *tmp_int_ptr;
|
|---|
| 1410 | int mydispls;
|
|---|
| 1411 |
|
|---|
| 1412 | void *send_buf = NULL;
|
|---|
| 1413 | void *recv_buf = NULL;
|
|---|
| 1414 |
|
|---|
| 1415 | hypre_BoxManEntry *entry;
|
|---|
| 1416 |
|
|---|
| 1417 | hypre_IndexRef index;
|
|---|
| 1418 |
|
|---|
| 1419 | hypre_Index imin, imax;
|
|---|
| 1420 |
|
|---|
| 1421 | void *index_ptr;
|
|---|
| 1422 |
|
|---|
| 1423 | /*how big is an entry? extents: 6 ints, proc: 1 int, id: 1
|
|---|
| 1424 | * int , num_ghost: 6 ints, info: info_size is in bytes*/
|
|---|
| 1425 | /* note: for now, we do not need to send num_ghost - this
|
|---|
| 1426 | is just copied in addentry anyhow */
|
|---|
| 1427 | non_info_size = 8;
|
|---|
| 1428 | entry_size_bytes = non_info_size*sizeof(int) + hypre_BoxManEntryInfoSize(manager);
|
|---|
| 1429 |
|
|---|
| 1430 | /* figure out how many entries each proc has - let the group know */
|
|---|
| 1431 | send_count = nentries;
|
|---|
| 1432 | send_count_bytes = send_count*entry_size_bytes;
|
|---|
| 1433 | recv_counts = hypre_CTAlloc(int, nprocs);
|
|---|
| 1434 |
|
|---|
| 1435 | MPI_Allgather(&send_count_bytes, 1, MPI_INT,
|
|---|
| 1436 | recv_counts, 1, MPI_INT, comm);
|
|---|
| 1437 |
|
|---|
| 1438 | displs = hypre_CTAlloc(int, nprocs);
|
|---|
| 1439 | displs[0] = 0;
|
|---|
| 1440 | recv_buf_size_bytes = recv_counts[0];
|
|---|
| 1441 | for (i = 1; i < nprocs; i++)
|
|---|
| 1442 | {
|
|---|
| 1443 | displs[i] = displs[i-1] + recv_counts[i-1];
|
|---|
| 1444 | recv_buf_size_bytes += recv_counts[i];
|
|---|
| 1445 | }
|
|---|
| 1446 | recv_buf_size = recv_buf_size_bytes/ entry_size_bytes;
|
|---|
| 1447 | mydispls = displs[myid]/entry_size_bytes;
|
|---|
| 1448 |
|
|---|
| 1449 | /* populate the send buffer with my entries*/
|
|---|
| 1450 | send_buf = hypre_MAlloc(send_count_bytes);
|
|---|
| 1451 | recv_buf = hypre_MAlloc(recv_buf_size_bytes);
|
|---|
| 1452 |
|
|---|
| 1453 | index_ptr = send_buf; /* step through send_buf with this pointer */
|
|---|
| 1454 | /* loop over my entries */
|
|---|
| 1455 | for (i=0; i < send_count; i++)
|
|---|
| 1456 | {
|
|---|
| 1457 | entry = &entries[i];
|
|---|
| 1458 |
|
|---|
| 1459 | size = sizeof(int);
|
|---|
| 1460 |
|
|---|
| 1461 | /* imin */
|
|---|
| 1462 | index = hypre_BoxManEntryIMin(entry);
|
|---|
| 1463 | for (d = 0; d < 3; d++)
|
|---|
| 1464 | {
|
|---|
| 1465 | tmp_int = hypre_IndexD(index, d);
|
|---|
| 1466 | memcpy( index_ptr, &tmp_int, size);
|
|---|
| 1467 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1468 | }
|
|---|
| 1469 |
|
|---|
| 1470 | /* imax */
|
|---|
| 1471 | index = hypre_BoxManEntryIMax(entry);
|
|---|
| 1472 | for (d = 0; d < 3; d++)
|
|---|
| 1473 | {
|
|---|
| 1474 | tmp_int = hypre_IndexD(index, d);
|
|---|
| 1475 | memcpy( index_ptr, &tmp_int, size);
|
|---|
| 1476 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1477 | }
|
|---|
| 1478 |
|
|---|
| 1479 | /* proc */
|
|---|
| 1480 | tmp_int = hypre_BoxManEntryProc(entry);
|
|---|
| 1481 | memcpy( index_ptr, &tmp_int, size);
|
|---|
| 1482 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1483 |
|
|---|
| 1484 | /* id */
|
|---|
| 1485 | tmp_int = hypre_BoxManEntryId(entry);
|
|---|
| 1486 | memcpy( index_ptr, &tmp_int, size);
|
|---|
| 1487 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1488 |
|
|---|
| 1489 | /* num_ghost (6 integers) - Don't send */
|
|---|
| 1490 | /* size = 6*size;
|
|---|
| 1491 | memcpy(index_ptr, hypre_BoxManEntryNumGhost(entry), size);
|
|---|
| 1492 | index_ptr = (void *) ((char *) index_ptr + size);*/
|
|---|
| 1493 |
|
|---|
| 1494 | /*info*/
|
|---|
| 1495 | size = hypre_BoxManEntryInfoSize(manager);
|
|---|
| 1496 | memcpy(index_ptr, hypre_BoxManEntryInfo(entry), size);
|
|---|
| 1497 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1498 |
|
|---|
| 1499 | } /* end of loop over my entries */
|
|---|
| 1500 |
|
|---|
| 1501 | /* now send_buf is full */
|
|---|
| 1502 |
|
|---|
| 1503 |
|
|---|
| 1504 | MPI_Allgatherv(send_buf, send_count_bytes, MPI_BYTE,
|
|---|
| 1505 | recv_buf, recv_counts, displs, MPI_BYTE, comm);
|
|---|
| 1506 |
|
|---|
| 1507 | /* unpack recv_buf into entries - note that the recv_buf size
|
|---|
| 1508 | includes our own entries - which we don't need to readd (so
|
|---|
| 1509 | size needs to be recv_buf_size - nentries)
|
|---|
| 1510 | */
|
|---|
| 1511 |
|
|---|
| 1512 | if (recv_buf_size > hypre_BoxManMaxNEntries(manager))
|
|---|
| 1513 | {
|
|---|
| 1514 | int inc_size;
|
|---|
| 1515 |
|
|---|
| 1516 | inc_size = (recv_buf_size - hypre_BoxManMaxNEntries(manager));
|
|---|
| 1517 | hypre_BoxManIncSize ( manager, inc_size);
|
|---|
| 1518 |
|
|---|
| 1519 | nentries = hypre_BoxManNEntries(manager);
|
|---|
| 1520 | entries = hypre_BoxManEntries(manager);
|
|---|
| 1521 | procs_sort = hypre_BoxManProcsSort(manager);
|
|---|
| 1522 | ids_sort = hypre_BoxManIdsSort(manager);
|
|---|
| 1523 | }
|
|---|
| 1524 |
|
|---|
| 1525 | index_ptr = recv_buf; /* point into recv buf */
|
|---|
| 1526 |
|
|---|
| 1527 | for (i = 0; i < recv_buf_size; i++)
|
|---|
| 1528 | {
|
|---|
| 1529 | /* don't add my own entries */
|
|---|
| 1530 | if ( i == mydispls)
|
|---|
| 1531 | {
|
|---|
| 1532 | /*increment i to skip over my entries */
|
|---|
| 1533 | i += (send_count - 1); /* will be increased by 1 by loop */
|
|---|
| 1534 |
|
|---|
| 1535 | /* fix pointer */
|
|---|
| 1536 | size = entry_size_bytes * send_count;
|
|---|
| 1537 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1538 | continue;
|
|---|
| 1539 |
|
|---|
| 1540 | }
|
|---|
| 1541 |
|
|---|
| 1542 |
|
|---|
| 1543 | size = sizeof(int);
|
|---|
| 1544 | /* imin */
|
|---|
| 1545 | for (d = 0; d < 3; d++)
|
|---|
| 1546 | {
|
|---|
| 1547 | memcpy( &tmp_int, index_ptr, size);
|
|---|
| 1548 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1549 | hypre_IndexD(imin, d) = tmp_int;
|
|---|
| 1550 | }
|
|---|
| 1551 |
|
|---|
| 1552 | /*imax */
|
|---|
| 1553 | for (d = 0; d < 3; d++)
|
|---|
| 1554 | {
|
|---|
| 1555 | memcpy( &tmp_int, index_ptr, size);
|
|---|
| 1556 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1557 | hypre_IndexD(imax, d) = tmp_int;
|
|---|
| 1558 | }
|
|---|
| 1559 |
|
|---|
| 1560 | /* proc */
|
|---|
| 1561 | tmp_int_ptr = (int *) index_ptr;
|
|---|
| 1562 | proc = *tmp_int_ptr;
|
|---|
| 1563 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1564 |
|
|---|
| 1565 | /* id */
|
|---|
| 1566 | tmp_int_ptr = (int *) index_ptr;
|
|---|
| 1567 | id = *tmp_int_ptr;
|
|---|
| 1568 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 1569 |
|
|---|
| 1570 | /* num_ghost (6 integers) - didn't send */
|
|---|
| 1571 | /* size = 6*size;
|
|---|
| 1572 | memcpy(hypre_BoxManEntryNumGhost(entry), index_ptr, size);
|
|---|
| 1573 | index_ptr = (void *) ((char *) index_ptr + size); */
|
|---|
| 1574 |
|
|---|
| 1575 | /* info */
|
|---|
| 1576 | /* index_ptr is at info */
|
|---|
| 1577 |
|
|---|
| 1578 | hypre_BoxManAddEntry( manager , imin ,
|
|---|
| 1579 | imax , proc, id,
|
|---|
| 1580 | index_ptr );
|
|---|
| 1581 |
|
|---|
| 1582 | /* start of next entry */
|
|---|
| 1583 | index_ptr = (void *) ((char *) index_ptr + hypre_BoxManEntryInfoSize(manager));
|
|---|
| 1584 | }
|
|---|
| 1585 |
|
|---|
| 1586 | hypre_BoxManAllGlobalKnown(manager) = 1;
|
|---|
| 1587 |
|
|---|
| 1588 | } /* end of non-AP gather */
|
|---|
| 1589 | #endif
|
|---|
| 1590 |
|
|---|
| 1591 | }/* end of if (gather entries) for both AP and non-AP */
|
|---|
| 1592 |
|
|---|
| 1593 |
|
|---|
| 1594 | nentries = hypre_BoxManNEntries(manager);
|
|---|
| 1595 | entries = hypre_BoxManEntries(manager);
|
|---|
| 1596 |
|
|---|
| 1597 |
|
|---|
| 1598 | /* -----------------------SORT TABLE --------------------------------------*/
|
|---|
| 1599 |
|
|---|
| 1600 | /* now everything we need is in entries, also ids and procs have
|
|---|
| 1601 | * been added, but not sorted */
|
|---|
| 1602 |
|
|---|
| 1603 | /* check for and remove duplicate boxes - based on (proc, id) */
|
|---|
| 1604 | /* at the same time sort the procs_sort and ids_sort and
|
|---|
| 1605 | * create the sort_table */
|
|---|
| 1606 | {
|
|---|
| 1607 | int *order_index;
|
|---|
| 1608 | int *delete_array;
|
|---|
| 1609 | int *delete_array_orig_order;
|
|---|
| 1610 | int *expand_delete_array;
|
|---|
| 1611 | int tmp_id, start, index;
|
|---|
| 1612 | int first_local;
|
|---|
| 1613 | int num_procs_sort;
|
|---|
| 1614 | int *proc_offsets;
|
|---|
| 1615 | int myoffset;
|
|---|
| 1616 |
|
|---|
| 1617 | hypre_BoxManEntry **sort_table;
|
|---|
| 1618 |
|
|---|
| 1619 | order_index = hypre_CTAlloc(int, nentries);
|
|---|
| 1620 | delete_array = hypre_CTAlloc(int, nentries);
|
|---|
| 1621 | index = 0;
|
|---|
| 1622 |
|
|---|
| 1623 | /* these are negative if a proc does not have any local entries
|
|---|
| 1624 | in the manager */
|
|---|
| 1625 | first_local = -1;
|
|---|
| 1626 | myoffset = -1;
|
|---|
| 1627 |
|
|---|
| 1628 |
|
|---|
| 1629 | for (i=0; i< nentries; i++)
|
|---|
| 1630 | {
|
|---|
| 1631 | order_index[i] = i;
|
|---|
| 1632 | }
|
|---|
| 1633 | /* sort by proc_id */
|
|---|
| 1634 | hypre_qsort3i(procs_sort, ids_sort, order_index, 0, nentries-1);
|
|---|
| 1635 | num_procs_sort = 0;
|
|---|
| 1636 | /* get first id */
|
|---|
| 1637 | if (nentries)
|
|---|
| 1638 | {
|
|---|
| 1639 | tmp_id = procs_sort[0];
|
|---|
| 1640 | num_procs_sort++;
|
|---|
| 1641 | if (tmp_id == myid ) first_local = 0;
|
|---|
| 1642 | }
|
|---|
| 1643 |
|
|---|
| 1644 | /* now sort on ids within each processor number*/
|
|---|
| 1645 | start = 0;
|
|---|
| 1646 | for (i=1; i< nentries; i++)
|
|---|
| 1647 | {
|
|---|
| 1648 | if (procs_sort[i] != tmp_id)
|
|---|
| 1649 | {
|
|---|
| 1650 | hypre_qsort2i(ids_sort, order_index, start, i-1);
|
|---|
| 1651 | /*now find duplicate ids */
|
|---|
| 1652 | for (j=start+1; j< i; j++)
|
|---|
| 1653 | {
|
|---|
| 1654 | if (ids_sort[j] == ids_sort[j-1])
|
|---|
| 1655 | {
|
|---|
| 1656 | delete_array[index++] = j;
|
|---|
| 1657 | }
|
|---|
| 1658 | }
|
|---|
| 1659 | /* update start and tmp_id */
|
|---|
| 1660 | start = i;
|
|---|
| 1661 | tmp_id = procs_sort[i];
|
|---|
| 1662 | num_procs_sort++;
|
|---|
| 1663 | if(tmp_id == myid )
|
|---|
| 1664 | {
|
|---|
| 1665 | /* subtract the value of index as this is how many previous we will
|
|---|
| 1666 | delete */
|
|---|
| 1667 | first_local = i - index;
|
|---|
| 1668 | }
|
|---|
| 1669 | }
|
|---|
| 1670 | }
|
|---|
| 1671 | /* final sort and purge (the last group doesn't
|
|---|
| 1672 | get caught in the above loop) */
|
|---|
| 1673 | if (nentries)
|
|---|
| 1674 | {
|
|---|
| 1675 | hypre_qsort2i(ids_sort, order_index, start, nentries-1);
|
|---|
| 1676 | /*now find duplicate boxnums */
|
|---|
| 1677 | for (j=start+1; j<nentries; j++)
|
|---|
| 1678 | {
|
|---|
| 1679 | if (ids_sort[j] == ids_sort[j-1])
|
|---|
| 1680 | {
|
|---|
| 1681 | delete_array[index++] = j;
|
|---|
| 1682 | }
|
|---|
| 1683 | }
|
|---|
| 1684 | }
|
|---|
| 1685 | /* now index = the number to delete (in delete_array) */
|
|---|
| 1686 | /* need to delete dublicates in *entries, and build the sorted pointers to
|
|---|
| 1687 | the entries */
|
|---|
| 1688 |
|
|---|
| 1689 |
|
|---|
| 1690 | if (index)
|
|---|
| 1691 | {
|
|---|
| 1692 | /* need to get the indices of the entries to delete as they
|
|---|
| 1693 | appear in the entries array */
|
|---|
| 1694 | delete_array_orig_order = hypre_CTAlloc(int, index);
|
|---|
| 1695 | for (i=0; i< index; i++)
|
|---|
| 1696 | {
|
|---|
| 1697 | delete_array_orig_order[i] = order_index[delete_array[i]];
|
|---|
| 1698 | }
|
|---|
| 1699 | qsort0(delete_array_orig_order, 0, index-1);
|
|---|
| 1700 |
|
|---|
| 1701 | /* delete entries - this also updates nentires in the box manager */
|
|---|
| 1702 | hypre_BoxManDeleteMultipleEntries(manager, delete_array_orig_order, index );
|
|---|
| 1703 |
|
|---|
| 1704 | hypre_TFree(delete_array_orig_order);
|
|---|
| 1705 |
|
|---|
| 1706 | /* now delete from sort procs and sort ids -use delete_array because these
|
|---|
| 1707 | have already been sorted. also delete from order_index */
|
|---|
| 1708 | start = delete_array[0];
|
|---|
| 1709 | j = 0;
|
|---|
| 1710 | for (i = start; (i + j) < nentries; i++)
|
|---|
| 1711 | {
|
|---|
| 1712 | if (j < index)
|
|---|
| 1713 | {
|
|---|
| 1714 | while ((i+j) == delete_array[j]) /* see if deleting consec. items */
|
|---|
| 1715 | {
|
|---|
| 1716 | j++; /*increase the shift*/
|
|---|
| 1717 | if (j == index) break;
|
|---|
| 1718 | }
|
|---|
| 1719 | }
|
|---|
| 1720 | if ((i+j) < nentries) /* if deleting the last item then no moving */
|
|---|
| 1721 | {
|
|---|
| 1722 | ids_sort[i] = ids_sort[i+j];
|
|---|
| 1723 | procs_sort[i] = procs_sort[i+j];
|
|---|
| 1724 | order_index[i] = order_index[i+j];
|
|---|
| 1725 | }
|
|---|
| 1726 |
|
|---|
| 1727 | }
|
|---|
| 1728 |
|
|---|
| 1729 | nentries = hypre_BoxManNEntries(manager); /* updated if some entries were deleted */
|
|---|
| 1730 | entries = hypre_BoxManEntries(manager);
|
|---|
| 1731 |
|
|---|
| 1732 | /* need to adjust the indexes of the entries to be deleted by the number of
|
|---|
| 1733 | entries delete prior to each index (expand_delete_array) */
|
|---|
| 1734 | expand_delete_array = hypre_CTAlloc(int, nentries);
|
|---|
| 1735 | j=0;
|
|---|
| 1736 | start = delete_array[0];
|
|---|
| 1737 | for (i=0; i < nentries; i++)
|
|---|
| 1738 | {
|
|---|
| 1739 | if (i == start)
|
|---|
| 1740 | {
|
|---|
| 1741 | j++;
|
|---|
| 1742 | if (j < index) start = delete_array[j];
|
|---|
| 1743 | }
|
|---|
| 1744 | expand_delete_array[i] = j;
|
|---|
| 1745 | }
|
|---|
| 1746 |
|
|---|
| 1747 | }
|
|---|
| 1748 | else /* index ==0 (no deleting ) */
|
|---|
| 1749 | {
|
|---|
| 1750 | expand_delete_array = hypre_CTAlloc(int, nentries); /* all zeros */
|
|---|
| 1751 | }
|
|---|
| 1752 |
|
|---|
| 1753 |
|
|---|
| 1754 | /* populate sort table with links to entries - have to account for
|
|---|
| 1755 | deleted entries*/
|
|---|
| 1756 | sort_table = hypre_CTAlloc(hypre_BoxManEntry *, nentries);
|
|---|
| 1757 |
|
|---|
| 1758 | for (i= 0; i< nentries; i++)
|
|---|
| 1759 | {
|
|---|
| 1760 | sort_table[i] = &entries[order_index[i] - expand_delete_array[i]];
|
|---|
| 1761 | }
|
|---|
| 1762 |
|
|---|
| 1763 | /* create proc_offsets (mmyoffset corresponds to local id)*/
|
|---|
| 1764 | proc_offsets = hypre_CTAlloc(int, num_procs_sort + 1);
|
|---|
| 1765 | proc_offsets[0] = 0;
|
|---|
| 1766 |
|
|---|
| 1767 | if (nentries > 0)
|
|---|
| 1768 | {
|
|---|
| 1769 | j=1;
|
|---|
| 1770 | tmp_id = procs_sort[0];
|
|---|
| 1771 | if (myid == tmp_id) myoffset =0;
|
|---|
| 1772 | for (i=0; i < nentries; i++)
|
|---|
| 1773 | {
|
|---|
| 1774 | if (procs_sort[i] != tmp_id)
|
|---|
| 1775 | {
|
|---|
| 1776 | if (myid == procs_sort[i]) myoffset = j;
|
|---|
| 1777 | proc_offsets[j++] = i;
|
|---|
| 1778 | tmp_id = procs_sort[i];
|
|---|
| 1779 |
|
|---|
| 1780 | }
|
|---|
| 1781 | }
|
|---|
| 1782 | proc_offsets[j] = nentries; /* last one */
|
|---|
| 1783 | }
|
|---|
| 1784 |
|
|---|
| 1785 |
|
|---|
| 1786 | /* clean up from this section of code */
|
|---|
| 1787 | hypre_TFree(expand_delete_array);
|
|---|
| 1788 | hypre_TFree(delete_array);
|
|---|
| 1789 | hypre_TFree(order_index);
|
|---|
| 1790 |
|
|---|
| 1791 | hypre_BoxManProcsSortOffsets(manager) = proc_offsets;
|
|---|
| 1792 | hypre_BoxManNumProcsSort(manager) = num_procs_sort;
|
|---|
| 1793 | hypre_BoxManFirstLocal(manager) = first_local;
|
|---|
| 1794 | hypre_BoxManSortTable(manager) = sort_table;
|
|---|
| 1795 | hypre_BoxManLocalProcOffset(manager) = myoffset;
|
|---|
| 1796 |
|
|---|
| 1797 |
|
|---|
| 1798 | }/* end bracket for all or the sorting stuff */
|
|---|
| 1799 |
|
|---|
| 1800 |
|
|---|
| 1801 | /*------------------------------INDEX TABLE ---------------------------*/
|
|---|
| 1802 |
|
|---|
| 1803 | /* now build the index_table and indexes array */
|
|---|
| 1804 | /* Note: for now we are using the same scheme as in BoxMap */
|
|---|
| 1805 | {
|
|---|
| 1806 | int *indexes[3];
|
|---|
| 1807 | int size[3];
|
|---|
| 1808 | int iminmax[2];
|
|---|
| 1809 | int index_not_there;
|
|---|
| 1810 | int d, e;
|
|---|
| 1811 | int mystart, myfinish;
|
|---|
| 1812 | int imin[3];
|
|---|
| 1813 | int imax[3];
|
|---|
| 1814 | int start_loop[3];
|
|---|
| 1815 | int end_loop[3];
|
|---|
| 1816 | int loop, range, loop_num;
|
|---|
| 1817 | int *proc_offsets;
|
|---|
| 1818 |
|
|---|
| 1819 | hypre_BoxManEntry **index_table;
|
|---|
| 1820 | hypre_BoxManEntry *entry;
|
|---|
| 1821 | hypre_BoxManEntry **sort_table;
|
|---|
| 1822 |
|
|---|
| 1823 | hypre_IndexRef entry_imin;
|
|---|
| 1824 | hypre_IndexRef entry_imax;
|
|---|
| 1825 |
|
|---|
| 1826 |
|
|---|
| 1827 | /* initial */
|
|---|
| 1828 | nentries = hypre_BoxManNEntries(manager);
|
|---|
| 1829 | entries = hypre_BoxManEntries(manager);
|
|---|
| 1830 | sort_table = hypre_BoxManSortTable(manager);
|
|---|
| 1831 | proc_offsets = hypre_BoxManProcsSortOffsets(manager);
|
|---|
| 1832 |
|
|---|
| 1833 |
|
|---|
| 1834 | /*------------------------------------------------------
|
|---|
| 1835 | * Set up the indexes array and record the processor's
|
|---|
| 1836 | * entries. This will be used in ordering the link list
|
|---|
| 1837 | * of BoxManEntry- ones on this processor listed first.
|
|---|
| 1838 | *------------------------------------------------------*/
|
|---|
| 1839 | for (d = 0; d < 3; d++)
|
|---|
| 1840 | {
|
|---|
| 1841 | indexes[d] = hypre_CTAlloc(int, 2*nentries);/* room for min and max of
|
|---|
| 1842 | each entry in each
|
|---|
| 1843 | dimension */
|
|---|
| 1844 | size[d] = 0;
|
|---|
| 1845 | }
|
|---|
| 1846 | /* loop through each entry and get index */
|
|---|
| 1847 | for (e = 0; e < nentries; e++)
|
|---|
| 1848 | {
|
|---|
| 1849 | entry = &entries[e]; /* grab the entry - get min and max extents */
|
|---|
| 1850 | entry_imin = hypre_BoxManEntryIMin(entry);
|
|---|
| 1851 | entry_imax = hypre_BoxManEntryIMax(entry);
|
|---|
| 1852 |
|
|---|
| 1853 | for (d = 0; d < 3; d++) /* in each dim, check if the min and max
|
|---|
| 1854 | position is there already in the index
|
|---|
| 1855 | table */
|
|---|
| 1856 | {
|
|---|
| 1857 | iminmax[0] = hypre_IndexD(entry_imin, d);
|
|---|
| 1858 | iminmax[1] = hypre_IndexD(entry_imax, d) + 1;
|
|---|
| 1859 |
|
|---|
| 1860 | for (i = 0; i < 2; i++)
|
|---|
| 1861 | {
|
|---|
| 1862 | /* find the new index position in the indexes array */
|
|---|
| 1863 | index_not_there = 1;
|
|---|
| 1864 | for (j = 0; j < size[d]; j++) /* size indicates current
|
|---|
| 1865 | size of index in that dimension*/
|
|---|
| 1866 | {
|
|---|
| 1867 | if (iminmax[i] <= indexes[d][j])
|
|---|
| 1868 | {
|
|---|
| 1869 | if (iminmax[i] == indexes[d][j])
|
|---|
| 1870 | {
|
|---|
| 1871 | index_not_there = 0;
|
|---|
| 1872 | }
|
|---|
| 1873 | break;
|
|---|
| 1874 | }
|
|---|
| 1875 | }
|
|---|
| 1876 |
|
|---|
| 1877 | /* if the index is already there, don't add it again */
|
|---|
| 1878 | if (index_not_there)
|
|---|
| 1879 | {
|
|---|
| 1880 | for (k = size[d]; k > j; k--) /* make room for new index */
|
|---|
| 1881 | {
|
|---|
| 1882 | indexes[d][k] = indexes[d][k-1];
|
|---|
| 1883 | }
|
|---|
| 1884 | indexes[d][j] = iminmax[i];
|
|---|
| 1885 | size[d]++; /* increase the size in that dimension */
|
|---|
| 1886 | }
|
|---|
| 1887 | } /* end of for min and max */
|
|---|
| 1888 | } /* end of for each dimension of the entry */
|
|---|
| 1889 | } /* end of for each entry loop */
|
|---|
| 1890 |
|
|---|
| 1891 | /* size was increased at the end - adjust when finished */
|
|---|
| 1892 | if (nentries)
|
|---|
| 1893 | {
|
|---|
| 1894 | for (d = 0; d < 3; d++)
|
|---|
| 1895 | {
|
|---|
| 1896 | size[d]--;
|
|---|
| 1897 | }
|
|---|
| 1898 | }
|
|---|
| 1899 |
|
|---|
| 1900 | /*------------------------------------------------------
|
|---|
| 1901 | * Set up the table - do offprocessor then on-processor
|
|---|
| 1902 | *------------------------------------------------------*/
|
|---|
| 1903 |
|
|---|
| 1904 | /* allocate space for table */
|
|---|
| 1905 | index_table = hypre_CTAlloc(hypre_BoxManEntry *, (size[0] * size[1] * size[2]));
|
|---|
| 1906 |
|
|---|
| 1907 | /* which are my entries? (on-processor) */
|
|---|
| 1908 | mystart = hypre_BoxManFirstLocal(manager);
|
|---|
| 1909 | if (mystart >= 0 ) /* we have local entries) because
|
|---|
| 1910 | firstlocal = -1 if no local entries */
|
|---|
| 1911 | {
|
|---|
| 1912 | loop_num = 3;
|
|---|
| 1913 | /* basically we have need to do the same code fragment
|
|---|
| 1914 | repeated three times so that we can do off-proc
|
|---|
| 1915 | then on proc entries - this ordering is because creating
|
|---|
| 1916 | the linked list for overlapping boxes */
|
|---|
| 1917 |
|
|---|
| 1918 |
|
|---|
| 1919 | myfinish = proc_offsets[hypre_BoxManLocalProcOffset(manager)+1];
|
|---|
| 1920 | /* #1 do off proc. entries - lower range */
|
|---|
| 1921 | start_loop[0] = 0;
|
|---|
| 1922 | end_loop[0] = mystart;
|
|---|
| 1923 | /* #2 do off proc. entries - upper range */
|
|---|
| 1924 | start_loop[1]= myfinish;
|
|---|
| 1925 | end_loop[1] = nentries;
|
|---|
| 1926 | /* #3 do ON proc. entries */
|
|---|
| 1927 | start_loop[2]= mystart;
|
|---|
| 1928 | end_loop[2] = myfinish;
|
|---|
| 1929 | }
|
|---|
| 1930 | else /* no on-proc entries */
|
|---|
| 1931 | {
|
|---|
| 1932 | loop_num = 1;
|
|---|
| 1933 | start_loop[0] = 0;
|
|---|
| 1934 | end_loop[0] = nentries;
|
|---|
| 1935 | }
|
|---|
| 1936 |
|
|---|
| 1937 | for (loop = 0; loop < loop_num; loop++)
|
|---|
| 1938 | {
|
|---|
| 1939 | for (range = start_loop[loop]; range < end_loop[loop]; range++)
|
|---|
| 1940 | {
|
|---|
| 1941 | entry = sort_table[range];
|
|---|
| 1942 | entry_imin = hypre_BoxManEntryIMin(entry);
|
|---|
| 1943 | entry_imax = hypre_BoxManEntryIMax(entry);
|
|---|
| 1944 |
|
|---|
| 1945 | /* find the indexes corresponding to the current box - put in
|
|---|
| 1946 | imin and imax */
|
|---|
| 1947 | for (d = 0; d < 3; d++)
|
|---|
| 1948 | {
|
|---|
| 1949 | j = 0;
|
|---|
| 1950 | while (hypre_IndexD(entry_imin, d) != indexes[d][j])
|
|---|
| 1951 | {
|
|---|
| 1952 | j++;
|
|---|
| 1953 | }
|
|---|
| 1954 | hypre_IndexD(imin, d) = j;
|
|---|
| 1955 |
|
|---|
| 1956 | while (hypre_IndexD(entry_imax, d) + 1 != indexes[d][j])
|
|---|
| 1957 | {
|
|---|
| 1958 | j++;
|
|---|
| 1959 | }
|
|---|
| 1960 | hypre_IndexD(imax, d) = j;
|
|---|
| 1961 | } /* now have imin and imax location in index array*/
|
|---|
| 1962 |
|
|---|
| 1963 |
|
|---|
| 1964 | /* set up index table */
|
|---|
| 1965 | for (k = imin[2]; k < imax[2]; k++)
|
|---|
| 1966 | {
|
|---|
| 1967 | for (j = imin[1]; j < imax[1]; j++)
|
|---|
| 1968 | {
|
|---|
| 1969 | for (i = imin[0]; i < imax[0]; i++)
|
|---|
| 1970 | {
|
|---|
| 1971 | if (!(index_table[((k) * size[1] + j) * size[0] + i])) /* no entry -
|
|---|
| 1972 | add one */
|
|---|
| 1973 | {
|
|---|
| 1974 | index_table[((k) * size[1] + j) * size[0] + i] = entry;
|
|---|
| 1975 | }
|
|---|
| 1976 | else /* already and entry there - so add to link list for
|
|---|
| 1977 | BoxMapEntry- overlapping */
|
|---|
| 1978 | {
|
|---|
| 1979 | hypre_BoxManEntryNext(entry)= index_table[((k) * size[1] + j) * size[0] + i];
|
|---|
| 1980 | index_table[((k) * size[1] + j) * size[0] + i]= entry;
|
|---|
| 1981 | }
|
|---|
| 1982 | }
|
|---|
| 1983 | }
|
|---|
| 1984 | }
|
|---|
| 1985 | } /* end of subset of entries */
|
|---|
| 1986 | }/* end of three loops over subsets */
|
|---|
| 1987 |
|
|---|
| 1988 |
|
|---|
| 1989 | /* done with the index_table! */
|
|---|
| 1990 | hypre_TFree( hypre_BoxManIndexTable(manager)); /* in case this is a
|
|---|
| 1991 | re-assemble */
|
|---|
| 1992 | hypre_BoxManIndexTable(manager) = index_table;
|
|---|
| 1993 |
|
|---|
| 1994 | for (d = 0; d < 3; d++)
|
|---|
| 1995 | {
|
|---|
| 1996 | hypre_TFree(hypre_BoxManIndexesD(manager, d));
|
|---|
| 1997 | hypre_BoxManIndexesD(manager, d) = indexes[d];
|
|---|
| 1998 | hypre_BoxManSizeD(manager, d) = size[d];
|
|---|
| 1999 | hypre_BoxManLastIndexD(manager, d) = 0;
|
|---|
| 2000 | }
|
|---|
| 2001 |
|
|---|
| 2002 | } /* end of building index table group */
|
|---|
| 2003 |
|
|---|
| 2004 |
|
|---|
| 2005 |
|
|---|
| 2006 | /* clean up and update*/
|
|---|
| 2007 |
|
|---|
| 2008 | hypre_BoxManNEntries(manager) = nentries;
|
|---|
| 2009 | hypre_BoxManEntries(manager) = entries;
|
|---|
| 2010 |
|
|---|
| 2011 | hypre_BoxManIsGatherCalled(manager) = 0;
|
|---|
| 2012 | hypre_BoxArrayDestroy(gather_regions);
|
|---|
| 2013 | hypre_BoxManGatherRegions(manager) = hypre_BoxArrayCreate(0);
|
|---|
| 2014 |
|
|---|
| 2015 |
|
|---|
| 2016 | return hypre_error_flag;
|
|---|
| 2017 |
|
|---|
| 2018 | }
|
|---|
| 2019 |
|
|---|
| 2020 |
|
|---|
| 2021 |
|
|---|
| 2022 |
|
|---|
| 2023 | /*--------------------------------------------------------------------------
|
|---|
| 2024 |
|
|---|
| 2025 | hypre_BoxManIntersect:
|
|---|
| 2026 |
|
|---|
| 2027 | Given a box (lower and upper indices), return a list of boxes in the
|
|---|
| 2028 | global grid that are intersected by this box. The user must insure
|
|---|
| 2029 | that a processor owns the correct global information to do the
|
|---|
| 2030 | intersection. For now this is virtually the same as the box map intersect.
|
|---|
| 2031 |
|
|---|
| 2032 | Notes:
|
|---|
| 2033 |
|
|---|
| 2034 | (1) This function can also be used in the way that
|
|---|
| 2035 | hypre_BoxMapFindEntry was previously used - just pass in iupper=ilower.
|
|---|
| 2036 | (TEST THIS)
|
|---|
| 2037 |
|
|---|
| 2038 | *--------------------------------------------------------------------------*/
|
|---|
| 2039 |
|
|---|
| 2040 | int hypre_BoxManIntersect ( hypre_BoxManager *manager , hypre_Index ilower ,
|
|---|
| 2041 | hypre_Index iupper ,
|
|---|
| 2042 | hypre_BoxManEntry ***entries_ptr , int *nentries_ptr )
|
|---|
| 2043 |
|
|---|
| 2044 | {
|
|---|
| 2045 | int d, i, j, k;
|
|---|
| 2046 | int find_index_d;
|
|---|
| 2047 | int current_index_d;
|
|---|
| 2048 | int *man_indexes_d;
|
|---|
| 2049 | int man_index_size_d;
|
|---|
| 2050 | int cnt, nentries;
|
|---|
| 2051 | int *ii, *jj, *kk;
|
|---|
| 2052 | int *proc_ids, *ids, *unsort;
|
|---|
| 2053 | int tmp_id, start;
|
|---|
| 2054 |
|
|---|
| 2055 | int man_ilower[3] = {0, 0, 0};
|
|---|
| 2056 | int man_iupper[3] = {0, 0, 0};
|
|---|
| 2057 |
|
|---|
| 2058 | hypre_BoxManEntry **entries, **all_entries;
|
|---|
| 2059 | hypre_BoxManEntry *entry;
|
|---|
| 2060 |
|
|---|
| 2061 |
|
|---|
| 2062 |
|
|---|
| 2063 | /* loop through each dimension */
|
|---|
| 2064 | for (d = 0; d < 3; d++)
|
|---|
| 2065 | {
|
|---|
| 2066 | man_indexes_d = hypre_BoxManIndexesD(manager, d);
|
|---|
| 2067 | man_index_size_d = hypre_BoxManSizeD(manager, d);
|
|---|
| 2068 |
|
|---|
| 2069 | /* -----find location of ilower[d] in indexes-----*/
|
|---|
| 2070 | find_index_d = hypre_IndexD(ilower, d);
|
|---|
| 2071 |
|
|---|
| 2072 | /* Start looking in place indicated by last_index stored in map */
|
|---|
| 2073 | current_index_d = hypre_BoxManLastIndexD(manager, d);
|
|---|
| 2074 |
|
|---|
| 2075 | /* Loop downward if target index is less than current location */
|
|---|
| 2076 | while ( (current_index_d >= 0 ) &&
|
|---|
| 2077 | (find_index_d < man_indexes_d[current_index_d]) )
|
|---|
| 2078 | {
|
|---|
| 2079 | current_index_d --;
|
|---|
| 2080 | }
|
|---|
| 2081 |
|
|---|
| 2082 | /* Loop upward if target index is greater than current location */
|
|---|
| 2083 | while ( (current_index_d <= (man_index_size_d - 1)) &&
|
|---|
| 2084 | (find_index_d >= man_indexes_d[current_index_d + 1]) )
|
|---|
| 2085 | {
|
|---|
| 2086 | current_index_d ++;
|
|---|
| 2087 | }
|
|---|
| 2088 |
|
|---|
| 2089 | if( current_index_d > (man_index_size_d - 1) )
|
|---|
| 2090 | {
|
|---|
| 2091 | *entries_ptr = NULL;
|
|---|
| 2092 | *nentries_ptr = 0;
|
|---|
| 2093 | return hypre_error_flag;
|
|---|
| 2094 | }
|
|---|
| 2095 | else
|
|---|
| 2096 | {
|
|---|
| 2097 | man_ilower[d] = hypre_max(current_index_d, 0);
|
|---|
| 2098 | }
|
|---|
| 2099 |
|
|---|
| 2100 | /* -----find location of iupper[d] in indexes-----*/
|
|---|
| 2101 |
|
|---|
| 2102 | find_index_d = hypre_IndexD(iupper, d);
|
|---|
| 2103 |
|
|---|
| 2104 | /* Loop upward if target index is greater than current location */
|
|---|
| 2105 | while ( (current_index_d <= (man_index_size_d-1)) &&
|
|---|
| 2106 | (find_index_d >= man_indexes_d[current_index_d+1]) )
|
|---|
| 2107 | {
|
|---|
| 2108 | current_index_d ++;
|
|---|
| 2109 | }
|
|---|
| 2110 | if( current_index_d < 0 )
|
|---|
| 2111 | {
|
|---|
| 2112 | *entries_ptr = NULL;
|
|---|
| 2113 | *nentries_ptr = 0;
|
|---|
| 2114 | return hypre_error_flag;
|
|---|
| 2115 | }
|
|---|
| 2116 | else
|
|---|
| 2117 | {
|
|---|
| 2118 | man_iupper[d] = hypre_min(current_index_d, (man_index_size_d-1)) + 1;
|
|---|
| 2119 | }
|
|---|
| 2120 |
|
|---|
| 2121 | }
|
|---|
| 2122 |
|
|---|
| 2123 | /*-----------------------------------------------------------------
|
|---|
| 2124 | * If code reaches this point, then set up the entries array and
|
|---|
| 2125 | * eliminate duplicates. To eliminate duplicates, we need to
|
|---|
| 2126 | * compare the BoxMapEntry link lists. We accomplish this using
|
|---|
| 2127 | * the unique (proc, id) identifier.
|
|---|
| 2128 | *-----------------------------------------------------------------*/
|
|---|
| 2129 |
|
|---|
| 2130 | /* how many entries */
|
|---|
| 2131 | nentries = ((man_iupper[0] - man_ilower[0]) *
|
|---|
| 2132 | (man_iupper[1] - man_ilower[1]) *
|
|---|
| 2133 | (man_iupper[2] - man_ilower[2]));
|
|---|
| 2134 |
|
|---|
| 2135 | ii= hypre_CTAlloc(int, nentries);
|
|---|
| 2136 | jj= hypre_CTAlloc(int, nentries);
|
|---|
| 2137 | kk= hypre_CTAlloc(int, nentries);
|
|---|
| 2138 |
|
|---|
| 2139 | nentries = 0;
|
|---|
| 2140 | cnt= 0;
|
|---|
| 2141 | for (k = man_ilower[2]; k < man_iupper[2]; k++)
|
|---|
| 2142 | {
|
|---|
| 2143 | for (j = man_ilower[1]; j < man_iupper[1]; j++)
|
|---|
| 2144 | {
|
|---|
| 2145 | for (i = man_ilower[0]; i < man_iupper[0]; i++)
|
|---|
| 2146 | {
|
|---|
| 2147 | /* the next 3 `if' statements eliminate duplicates */
|
|---|
| 2148 | if (k > man_ilower[2])
|
|---|
| 2149 | {
|
|---|
| 2150 | if ( hypre_BoxManIndexTableEntry(manager, i, j, k) ==
|
|---|
| 2151 | hypre_BoxManIndexTableEntry(manager, i, j, (k-1)) )
|
|---|
| 2152 | {
|
|---|
| 2153 | continue;
|
|---|
| 2154 | }
|
|---|
| 2155 | }
|
|---|
| 2156 | if (j > man_ilower[1])
|
|---|
| 2157 | {
|
|---|
| 2158 | if ( hypre_BoxManIndexTableEntry(manager, i, j, k) ==
|
|---|
| 2159 | hypre_BoxManIndexTableEntry(manager, i, (j-1), k) )
|
|---|
| 2160 | {
|
|---|
| 2161 | continue;
|
|---|
| 2162 | }
|
|---|
| 2163 | }
|
|---|
| 2164 | if (i > man_ilower[0])
|
|---|
| 2165 | {
|
|---|
| 2166 | if ( hypre_BoxManIndexTableEntry(manager, i, j, k) ==
|
|---|
| 2167 | hypre_BoxManIndexTableEntry(manager, (i-1), j, k) )
|
|---|
| 2168 | {
|
|---|
| 2169 | continue;
|
|---|
| 2170 | }
|
|---|
| 2171 | }
|
|---|
| 2172 |
|
|---|
| 2173 | entry = hypre_BoxManIndexTableEntry(manager, i, j, k);
|
|---|
| 2174 | /* Record the indices for non-empty entries and count all ManEntries. */
|
|---|
| 2175 | if (entry != NULL)
|
|---|
| 2176 | {
|
|---|
| 2177 | ii[nentries]= i;
|
|---|
| 2178 | jj[nentries]= j;
|
|---|
| 2179 | kk[nentries++]= k;
|
|---|
| 2180 |
|
|---|
| 2181 | while (entry)
|
|---|
| 2182 | {
|
|---|
| 2183 | cnt++;
|
|---|
| 2184 | entry= hypre_BoxManEntryNext(entry);
|
|---|
| 2185 | }
|
|---|
| 2186 | }
|
|---|
| 2187 | }
|
|---|
| 2188 | }
|
|---|
| 2189 | }
|
|---|
| 2190 |
|
|---|
| 2191 |
|
|---|
| 2192 | /* if no link lists of BoxManEntries, Just point to the unique BoxManEntries */
|
|---|
| 2193 | if (nentries == cnt)
|
|---|
| 2194 | {
|
|---|
| 2195 | entries = hypre_CTAlloc(hypre_BoxManEntry *, nentries);
|
|---|
| 2196 | for (i= 0; i< nentries; i++)
|
|---|
| 2197 | {
|
|---|
| 2198 | entries[i]= hypre_BoxManIndexTableEntry(manager, ii[i], jj[i], kk[i]);
|
|---|
| 2199 | }
|
|---|
| 2200 | }
|
|---|
| 2201 |
|
|---|
| 2202 | /* otherwise: link lists of BoxManEntries. Sorting needed to eliminate duplicates. */
|
|---|
| 2203 | else
|
|---|
| 2204 | {
|
|---|
| 2205 | unsort = hypre_CTAlloc(int, cnt);
|
|---|
| 2206 | proc_ids = hypre_CTAlloc(int, cnt);
|
|---|
| 2207 | ids = hypre_CTAlloc(int, cnt);
|
|---|
| 2208 | all_entries = hypre_CTAlloc(hypre_BoxManEntry *, cnt);
|
|---|
| 2209 |
|
|---|
| 2210 | cnt = 0;
|
|---|
| 2211 | for (i= 0; i< nentries; i++)
|
|---|
| 2212 | {
|
|---|
| 2213 | entry = hypre_BoxManIndexTableEntry(manager, ii[i], jj[i], kk[i]);
|
|---|
| 2214 |
|
|---|
| 2215 | while (entry)
|
|---|
| 2216 | {
|
|---|
| 2217 | all_entries[cnt]= entry;
|
|---|
| 2218 | unsort[cnt] = cnt;
|
|---|
| 2219 |
|
|---|
| 2220 | ids[cnt] = hypre_BoxManEntryId(entry);
|
|---|
| 2221 | proc_ids[cnt++] = hypre_BoxManEntryProc(entry);
|
|---|
| 2222 |
|
|---|
| 2223 | entry= hypre_BoxManEntryNext(entry);
|
|---|
| 2224 | }
|
|---|
| 2225 | }
|
|---|
| 2226 | /* now we have a list of all of the entries - just want unique ones */
|
|---|
| 2227 | /* first sort on proc id */
|
|---|
| 2228 | hypre_qsort3i(proc_ids, ids, unsort, 0, cnt-1);
|
|---|
| 2229 | /* now sort on ids within each processor number*/
|
|---|
| 2230 | tmp_id = proc_ids[0];
|
|---|
| 2231 | start = 0;
|
|---|
| 2232 | for (i=1; i< cnt; i++)
|
|---|
| 2233 | {
|
|---|
| 2234 | if (proc_ids[i] != tmp_id)
|
|---|
| 2235 | {
|
|---|
| 2236 | hypre_qsort2i(ids, unsort, start, i-1);
|
|---|
| 2237 | /* update start and tmp_id */
|
|---|
| 2238 | start = i;
|
|---|
| 2239 | tmp_id = proc_ids[i];
|
|---|
| 2240 | }
|
|---|
| 2241 | }
|
|---|
| 2242 | if (cnt > 1) /*sort last group */
|
|---|
| 2243 | {
|
|---|
| 2244 | hypre_qsort2i(ids, unsort, start, nentries-1);
|
|---|
| 2245 | }
|
|---|
| 2246 |
|
|---|
| 2247 |
|
|---|
| 2248 | /* count the unique Entries */
|
|---|
| 2249 | nentries = 1;
|
|---|
| 2250 | for (i = 1; i< cnt; i++)
|
|---|
| 2251 | {
|
|---|
| 2252 | if (!(proc_ids[i] = proc_ids[i-1] && ids[i] == ids[i-1]))
|
|---|
| 2253 | {
|
|---|
| 2254 | nentries++;
|
|---|
| 2255 | }
|
|---|
| 2256 | }
|
|---|
| 2257 |
|
|---|
| 2258 | entries = hypre_CTAlloc(hypre_BoxManEntry *, nentries);
|
|---|
| 2259 |
|
|---|
| 2260 | /* extract the unique Entries */
|
|---|
| 2261 | entries[0] = all_entries[unsort[0]];
|
|---|
| 2262 | nentries = 1;
|
|---|
| 2263 |
|
|---|
| 2264 | for (i= 1; i< cnt; i++)
|
|---|
| 2265 | {
|
|---|
| 2266 | if (!(proc_ids[i] = proc_ids[i-1] && ids[i] == ids[i-1]))
|
|---|
| 2267 | {
|
|---|
| 2268 | entries[nentries++] = all_entries[unsort[i]];
|
|---|
| 2269 | }
|
|---|
| 2270 | }
|
|---|
| 2271 |
|
|---|
| 2272 | hypre_TFree(unsort);
|
|---|
| 2273 | hypre_TFree(ids);
|
|---|
| 2274 | hypre_TFree(proc_ids);
|
|---|
| 2275 | hypre_TFree(all_entries);
|
|---|
| 2276 | }
|
|---|
| 2277 |
|
|---|
| 2278 | hypre_TFree(ii);
|
|---|
| 2279 | hypre_TFree(jj);
|
|---|
| 2280 | hypre_TFree(kk);
|
|---|
| 2281 |
|
|---|
| 2282 | /* Reset the last index in the manager */
|
|---|
| 2283 | for (d = 0; d < 3; d++)
|
|---|
| 2284 | {
|
|---|
| 2285 | hypre_BoxManLastIndexD(manager, d) = man_ilower[d];
|
|---|
| 2286 | }
|
|---|
| 2287 |
|
|---|
| 2288 | *entries_ptr = entries;
|
|---|
| 2289 | *nentries_ptr = nentries;
|
|---|
| 2290 |
|
|---|
| 2291 | return hypre_error_flag;
|
|---|
| 2292 |
|
|---|
| 2293 | }
|
|---|
| 2294 |
|
|---|
| 2295 |
|
|---|
| 2296 |
|
|---|
| 2297 | /*------------------------------------------------------------------------------
|
|---|
| 2298 | * hypre_BoxManSetNumGhost
|
|---|
| 2299 | *-----------------------------------------------------------------------------*/
|
|---|
| 2300 |
|
|---|
| 2301 |
|
|---|
| 2302 | int
|
|---|
| 2303 | hypre_BoxManSetNumGhost( hypre_BoxManager *manager, int *num_ghost )
|
|---|
| 2304 | {
|
|---|
| 2305 |
|
|---|
| 2306 | int i;
|
|---|
| 2307 |
|
|---|
| 2308 | for (i = 0; i < 6; i++)
|
|---|
| 2309 | {
|
|---|
| 2310 | hypre_BoxManNumGhost(manager)[i] = num_ghost[i];
|
|---|
| 2311 | }
|
|---|
| 2312 |
|
|---|
| 2313 | return hypre_error_flag;
|
|---|
| 2314 |
|
|---|
| 2315 | }
|
|---|
| 2316 |
|
|---|
| 2317 |
|
|---|
| 2318 |
|
|---|
| 2319 | /******************************************************************************
|
|---|
| 2320 | hypre_fillResponseBoxManAssemble1
|
|---|
| 2321 |
|
|---|
| 2322 | contact message is null. need to return the id and boxnum of each box
|
|---|
| 2323 | in our assumed partition.
|
|---|
| 2324 |
|
|---|
| 2325 | *****************************************************************************/
|
|---|
| 2326 |
|
|---|
| 2327 | int
|
|---|
| 2328 | hypre_FillResponseBoxMapAssemble1(void *p_recv_contact_buf,
|
|---|
| 2329 | int contact_size, int contact_proc, void *ro,
|
|---|
| 2330 | MPI_Comm comm, void **p_send_response_buf,
|
|---|
| 2331 | int *response_message_size )
|
|---|
| 2332 | {
|
|---|
| 2333 |
|
|---|
| 2334 |
|
|---|
| 2335 | int myid, i, index;
|
|---|
| 2336 | int size, num_objects;
|
|---|
| 2337 | int *proc_ids;
|
|---|
| 2338 | int *boxnums;
|
|---|
| 2339 | int *send_response_buf = (int *) *p_send_response_buf;
|
|---|
| 2340 |
|
|---|
| 2341 |
|
|---|
| 2342 | hypre_DataExchangeResponse *response_obj = ro;
|
|---|
| 2343 | hypre_StructAssumedPart *ap = response_obj->data1;
|
|---|
| 2344 |
|
|---|
| 2345 | int overhead = response_obj->send_response_overhead;
|
|---|
| 2346 |
|
|---|
| 2347 | /*initialize stuff */
|
|---|
| 2348 | MPI_Comm_rank(comm, &myid );
|
|---|
| 2349 |
|
|---|
| 2350 | proc_ids = hypre_StructAssumedPartMyPartitionProcIds(ap);
|
|---|
| 2351 | boxnums = hypre_StructAssumedPartMyPartitionBoxnums(ap);
|
|---|
| 2352 |
|
|---|
| 2353 | /*we need to send back the list of all of the boxnums and
|
|---|
| 2354 | corresponding processor id */
|
|---|
| 2355 |
|
|---|
| 2356 | /*how many boxes do we have in the AP?*/
|
|---|
| 2357 | num_objects = hypre_StructAssumedPartMyPartitionIdsSize(ap);
|
|---|
| 2358 |
|
|---|
| 2359 | /* num_objects is then how much we need to send*/
|
|---|
| 2360 |
|
|---|
| 2361 |
|
|---|
| 2362 | /*check storage in send_buf for adding the information */
|
|---|
| 2363 | /* note: we are returning objects that are 2 ints in size */
|
|---|
| 2364 |
|
|---|
| 2365 | if ( response_obj->send_response_storage < num_objects )
|
|---|
| 2366 | {
|
|---|
| 2367 | response_obj->send_response_storage = hypre_max(num_objects, 10);
|
|---|
| 2368 | size = 2*(response_obj->send_response_storage + overhead);
|
|---|
| 2369 | send_response_buf = hypre_TReAlloc( send_response_buf, int,
|
|---|
| 2370 | size);
|
|---|
| 2371 | *p_send_response_buf = send_response_buf;
|
|---|
| 2372 | }
|
|---|
| 2373 |
|
|---|
| 2374 | /* populate send_response_buf */
|
|---|
| 2375 | index = 0;
|
|---|
| 2376 | for (i = 0; i< num_objects; i++)
|
|---|
| 2377 | {
|
|---|
| 2378 | /* processor id + boxnum */
|
|---|
| 2379 | send_response_buf[index++] = proc_ids[i];
|
|---|
| 2380 | send_response_buf[index++] = boxnums[i];
|
|---|
| 2381 | }
|
|---|
| 2382 |
|
|---|
| 2383 | /* return variables */
|
|---|
| 2384 | *response_message_size = num_objects;
|
|---|
| 2385 | *p_send_response_buf = send_response_buf;
|
|---|
| 2386 |
|
|---|
| 2387 | return hypre_error_flag;
|
|---|
| 2388 |
|
|---|
| 2389 |
|
|---|
| 2390 | }
|
|---|
| 2391 | /******************************************************************************
|
|---|
| 2392 |
|
|---|
| 2393 | hypre_fillResponseBoxManAssemble2
|
|---|
| 2394 |
|
|---|
| 2395 | contact message includes a box number (id). the response needs to be the entry
|
|---|
| 2396 | corresponding to the receiving processor and that id.
|
|---|
| 2397 |
|
|---|
| 2398 | *****************************************************************************/
|
|---|
| 2399 |
|
|---|
| 2400 | int
|
|---|
| 2401 | hypre_FillResponseBoxMapAssemble2(void *p_recv_contact_buf,
|
|---|
| 2402 | int contact_size, int contact_proc, void *ro,
|
|---|
| 2403 | MPI_Comm comm, void **p_send_response_buf,
|
|---|
| 2404 | int *response_message_size )
|
|---|
| 2405 | {
|
|---|
| 2406 |
|
|---|
| 2407 |
|
|---|
| 2408 | int myid, i, d, num_entries, location, size;
|
|---|
| 2409 | int proc_id, box_id, tmp_int;
|
|---|
| 2410 | int num_entries_found = 0;
|
|---|
| 2411 | int entry_size_bytes;
|
|---|
| 2412 |
|
|---|
| 2413 |
|
|---|
| 2414 | int *recv_contact_buf = (int *) p_recv_contact_buf;
|
|---|
| 2415 |
|
|---|
| 2416 |
|
|---|
| 2417 | void *send_response_buf = (void *) *p_send_response_buf;
|
|---|
| 2418 | void *index_ptr;
|
|---|
| 2419 |
|
|---|
| 2420 | hypre_BoxManEntry *entry;
|
|---|
| 2421 |
|
|---|
| 2422 | hypre_IndexRef index;
|
|---|
| 2423 |
|
|---|
| 2424 | hypre_DataExchangeResponse *response_obj = ro;
|
|---|
| 2425 |
|
|---|
| 2426 | hypre_BoxManager *manager = response_obj->data1;
|
|---|
| 2427 |
|
|---|
| 2428 | int overhead = response_obj->send_response_overhead;
|
|---|
| 2429 |
|
|---|
| 2430 | hypre_BoxManEntry **my_entries = hypre_BoxManMyEntries(manager) ;
|
|---|
| 2431 |
|
|---|
| 2432 | int *my_ids = hypre_BoxManMyIds(manager);
|
|---|
| 2433 | int num_my_entries = hypre_BoxManNumMyEntries(manager);
|
|---|
| 2434 |
|
|---|
| 2435 | /*initialize stuff */
|
|---|
| 2436 | MPI_Comm_rank(comm, &myid );
|
|---|
| 2437 |
|
|---|
| 2438 | entry_size_bytes = 8*sizeof(int) + hypre_BoxManEntryInfoSize(manager);
|
|---|
| 2439 |
|
|---|
| 2440 | /* num_entries is one for each boxnum sent (corres. to our proc. id only)*/
|
|---|
| 2441 | num_entries = contact_size;
|
|---|
| 2442 |
|
|---|
| 2443 | /*check storage in send_buf for adding the information */
|
|---|
| 2444 | if ( response_obj->send_response_storage < num_entries )
|
|---|
| 2445 | {
|
|---|
| 2446 | response_obj->send_response_storage = hypre_max(num_entries, 10);
|
|---|
| 2447 | size = entry_size_bytes*(response_obj->send_response_storage + overhead);
|
|---|
| 2448 | send_response_buf = hypre_ReAlloc( send_response_buf, size);
|
|---|
| 2449 | *p_send_response_buf = send_response_buf;
|
|---|
| 2450 | }
|
|---|
| 2451 |
|
|---|
| 2452 | index_ptr = send_response_buf; /* step through send_buf with this pointer */
|
|---|
| 2453 |
|
|---|
| 2454 | if (num_my_entries == 0)
|
|---|
| 2455 | {
|
|---|
| 2456 | num_entries_found = 0;
|
|---|
| 2457 | }
|
|---|
| 2458 | else
|
|---|
| 2459 | {
|
|---|
| 2460 | for (i=0; i < num_entries; i++)
|
|---|
| 2461 | {
|
|---|
| 2462 | proc_id = myid;
|
|---|
| 2463 | box_id = recv_contact_buf[i];
|
|---|
| 2464 |
|
|---|
| 2465 | /* find entry */
|
|---|
| 2466 | location = hypre_BinarySearch(my_ids, box_id, num_my_entries);
|
|---|
| 2467 | if (location >= 0)
|
|---|
| 2468 | {
|
|---|
| 2469 | entry = my_entries[location];
|
|---|
| 2470 | num_entries_found++;
|
|---|
| 2471 |
|
|---|
| 2472 | /*pack response buffer with information (check storage first) */
|
|---|
| 2473 |
|
|---|
| 2474 | size = sizeof(int);
|
|---|
| 2475 | /* imin */
|
|---|
| 2476 | index = hypre_BoxManEntryIMin(entry);
|
|---|
| 2477 | for (d = 0; d < 3; d++)
|
|---|
| 2478 | {
|
|---|
| 2479 | tmp_int = hypre_IndexD(index, d);
|
|---|
| 2480 | memcpy( index_ptr, &tmp_int, size);
|
|---|
| 2481 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 2482 | }
|
|---|
| 2483 | /* imax */
|
|---|
| 2484 | index = hypre_BoxManEntryIMax(entry);
|
|---|
| 2485 | for (d = 0; d < 3; d++)
|
|---|
| 2486 | {
|
|---|
| 2487 | tmp_int = hypre_IndexD(index, d);
|
|---|
| 2488 | memcpy( index_ptr, &tmp_int, size);
|
|---|
| 2489 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 2490 | }
|
|---|
| 2491 | /* proc */
|
|---|
| 2492 | memcpy( index_ptr, &proc_id, size);
|
|---|
| 2493 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 2494 |
|
|---|
| 2495 | /* id */
|
|---|
| 2496 | memcpy( index_ptr, &box_id, size);
|
|---|
| 2497 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 2498 |
|
|---|
| 2499 | /*info*/
|
|---|
| 2500 | size = hypre_BoxManEntryInfoSize(manager);
|
|---|
| 2501 | memcpy(index_ptr, hypre_BoxManEntryInfo(entry), size);
|
|---|
| 2502 | index_ptr = (void *) ((char *) index_ptr + size);
|
|---|
| 2503 |
|
|---|
| 2504 | }
|
|---|
| 2505 | }
|
|---|
| 2506 | }
|
|---|
| 2507 |
|
|---|
| 2508 | /* now send_response_buf is full */
|
|---|
| 2509 |
|
|---|
| 2510 | /* return variable */
|
|---|
| 2511 | *response_message_size = num_entries_found;
|
|---|
| 2512 | *p_send_response_buf = send_response_buf;
|
|---|
| 2513 |
|
|---|
| 2514 | return hypre_error_flag;
|
|---|
| 2515 |
|
|---|
| 2516 |
|
|---|
| 2517 | }
|
|---|
| 2518 |
|
|---|
| 2519 |
|
|---|
| 2520 |
|
|---|
| 2521 | /******************************************************************************
|
|---|
| 2522 |
|
|---|
| 2523 | some specialized sorting routines
|
|---|
| 2524 |
|
|---|
| 2525 | *****************************************************************************/
|
|---|
| 2526 |
|
|---|
| 2527 | /* sort on int i, move entry pointers ent */
|
|---|
| 2528 |
|
|---|
| 2529 | void hypre_entryqsort2( int *v,
|
|---|
| 2530 | hypre_BoxManEntry ** ent,
|
|---|
| 2531 | int left,
|
|---|
| 2532 | int right )
|
|---|
| 2533 | {
|
|---|
| 2534 | int i, last;
|
|---|
| 2535 |
|
|---|
| 2536 | if (left >= right)
|
|---|
| 2537 | {
|
|---|
| 2538 | return;
|
|---|
| 2539 | }
|
|---|
| 2540 | hypre_entryswap2( v, ent, left, (left+right)/2);
|
|---|
| 2541 | last = left;
|
|---|
| 2542 | for (i = left+1; i <= right; i++)
|
|---|
| 2543 | {
|
|---|
| 2544 | if (v[i] < v[left])
|
|---|
| 2545 | {
|
|---|
| 2546 | hypre_entryswap2(v, ent, ++last, i);
|
|---|
| 2547 | }
|
|---|
| 2548 | }
|
|---|
| 2549 | hypre_entryswap2(v, ent, left, last);
|
|---|
| 2550 | hypre_entryqsort2(v, ent, left, last-1);
|
|---|
| 2551 | hypre_entryqsort2(v, ent, last+1, right);
|
|---|
| 2552 | }
|
|---|
| 2553 |
|
|---|
| 2554 |
|
|---|
| 2555 | void hypre_entryswap2(int *v,
|
|---|
| 2556 | hypre_BoxManEntry ** ent,
|
|---|
| 2557 | int i,
|
|---|
| 2558 | int j )
|
|---|
| 2559 | {
|
|---|
| 2560 | int temp;
|
|---|
| 2561 |
|
|---|
| 2562 | hypre_BoxManEntry *temp_e;
|
|---|
| 2563 |
|
|---|
| 2564 | temp = v[i];
|
|---|
| 2565 | v[i] = v[j];
|
|---|
| 2566 | v[j] = temp;
|
|---|
| 2567 |
|
|---|
| 2568 | temp_e = ent[i];
|
|---|
| 2569 | ent[i] = ent[j];
|
|---|
| 2570 | ent[j] = temp_e;
|
|---|
| 2571 | }
|
|---|