source: CIVL/examples/mpi-omp/AMG2013/sstruct_mv/HYPRE_sstruct_graph.c@ beab7f2

main test-branch
Last change on this file since beab7f2 was ea777aa, checked in by Alex Wilton <awilton@…>, 3 years ago

Moved examples, include, build_default.properties, common.xml, and README out from dev.civl.com into the root of the repo.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@5704 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 23.3 KB
Line 
1/*BHEADER**********************************************************************
2 * Copyright (c) 2008, Lawrence Livermore National Security, LLC.
3 * Produced at the Lawrence Livermore National Laboratory.
4 * This file is part of HYPRE. See file COPYRIGHT for details.
5 *
6 * HYPRE is free software; you can redistribute it and/or modify it under the
7 * terms of the GNU Lesser General Public License (as published by the Free
8 * Software Foundation) version 2.1 dated February 1999.
9 *
10 * $Revision: 2.4 $
11 ***********************************************************************EHEADER*/
12
13
14
15/******************************************************************************
16 *
17 * HYPRE_SStructGraph interface
18 *
19 *****************************************************************************/
20
21#include "headers.h"
22
23/*--------------------------------------------------------------------------
24 * HYPRE_SStructGraphCreate
25 *--------------------------------------------------------------------------*/
26
27int
28HYPRE_SStructGraphCreate( MPI_Comm comm,
29 HYPRE_SStructGrid grid,
30 HYPRE_SStructGraph *graph_ptr )
31{
32 int ierr = 0;
33
34 hypre_SStructGraph *graph;
35 int nparts;
36 hypre_SStructStencil ***stencils;
37 hypre_SStructPGrid **pgrids;
38 int nvars;
39 int part, var;
40
41 graph = hypre_TAlloc(hypre_SStructGraph, 1);
42
43 hypre_SStructGraphComm(graph) = comm;
44 hypre_SStructGraphNDim(graph) = hypre_SStructGridNDim(grid);
45 hypre_SStructGridRef(grid, &hypre_SStructGraphGrid(graph));
46 nparts = hypre_SStructGridNParts(grid);
47 hypre_SStructGraphNParts(graph) = nparts;
48 pgrids = hypre_SStructGridPGrids(grid);
49 hypre_SStructGraphPGrids(graph) = pgrids;
50 stencils = hypre_TAlloc(hypre_SStructStencil **, nparts);
51 for (part = 0; part < nparts; part++)
52 {
53 nvars = hypre_SStructPGridNVars(pgrids[part]);
54 stencils[part] = hypre_TAlloc(hypre_SStructStencil *, nvars);
55 for (var = 0; var < nvars; var++)
56 {
57 stencils[part][var] = NULL;
58 }
59 }
60 hypre_SStructGraphStencils(graph) = stencils;
61
62 hypre_SStructGraphNUVEntries(graph) = 0;
63 hypre_SStructGraphAUVEntries(graph) = 0;
64 hypre_SStructGraphIUVEntries(graph) = NULL;
65
66 hypre_SStructGraphUVEntries(graph) = NULL;
67 hypre_SStructGraphTotUEntries(graph) = 0;
68 hypre_SStructGraphRefCount(graph) = 1;
69 hypre_SStructGraphObjectType(graph) = HYPRE_SSTRUCT;
70
71 *graph_ptr = graph;
72
73 return ierr;
74}
75
76/*--------------------------------------------------------------------------
77 * HYPRE_SStructGraphDestroy
78 *--------------------------------------------------------------------------*/
79
80int
81HYPRE_SStructGraphDestroy( HYPRE_SStructGraph graph )
82{
83 int ierr = 0;
84
85 int nparts;
86 hypre_SStructPGrid **pgrids;
87 hypre_SStructStencil ***stencils;
88 int nUventries;
89 int *iUventries;
90
91 hypre_SStructUVEntry **Uventries;
92 hypre_SStructUVEntry *Uventry;
93 int nvars;
94 int part, var, i;
95
96 if (graph)
97 {
98 hypre_SStructGraphRefCount(graph) --;
99 if (hypre_SStructGraphRefCount(graph) == 0)
100 {
101 nparts = hypre_SStructGraphNParts(graph);
102 pgrids = hypre_SStructGraphPGrids(graph);
103 stencils = hypre_SStructGraphStencils(graph);
104 nUventries = hypre_SStructGraphNUVEntries(graph);
105 iUventries = hypre_SStructGraphIUVEntries(graph);
106
107 Uventries = hypre_SStructGraphUVEntries(graph);
108 for (part = 0; part < nparts; part++)
109 {
110 nvars = hypre_SStructPGridNVars(pgrids[part]);
111 for (var = 0; var < nvars; var++)
112 {
113 HYPRE_SStructStencilDestroy(stencils[part][var]);
114 }
115 hypre_TFree(stencils[part]);
116 }
117 HYPRE_SStructGridDestroy(hypre_SStructGraphGrid(graph));
118 hypre_TFree(stencils);
119 for (i = 0; i < nUventries; i++)
120 {
121 Uventry = Uventries[iUventries[i]];
122 if (Uventry)
123 {
124 hypre_TFree(hypre_SStructUVEntryUEntries(Uventry));
125 hypre_TFree(Uventry);
126 }
127 Uventries[iUventries[i]] = NULL;
128 }
129 hypre_TFree(iUventries);
130 hypre_TFree(Uventries);
131 hypre_TFree(graph);
132 }
133 }
134
135 return ierr;
136}
137
138/*--------------------------------------------------------------------------
139 * HYPRE_SStructGraphSetStencil
140 *--------------------------------------------------------------------------*/
141
142int
143HYPRE_SStructGraphSetStencil( HYPRE_SStructGraph graph,
144 int part,
145 int var,
146 HYPRE_SStructStencil stencil )
147{
148 int ierr = 0;
149
150 hypre_SStructStencilRef(stencil,
151 &hypre_SStructGraphStencil(graph, part, var));
152
153 return ierr;
154}
155
156/*--------------------------------------------------------------------------
157 * HYPRE_SStructGraphAddEntries-
158 * THIS IS FOR A NON-OVERLAPPING GRID GRAPH.
159 *--------------------------------------------------------------------------*/
160
161int
162HYPRE_SStructGraphAddEntries( HYPRE_SStructGraph graph,
163 int part,
164 int *index,
165 int var,
166 int to_part,
167 int *to_index,
168 int to_var )
169{
170 int ierr = 0;
171
172 hypre_SStructGrid *grid = hypre_SStructGraphGrid(graph);
173 int ndim = hypre_SStructGridNDim(grid);
174 int nUventries = hypre_SStructGraphNUVEntries(graph);
175 int aUventries = hypre_SStructGraphAUVEntries(graph);
176 int *iUventries = hypre_SStructGraphIUVEntries(graph);
177
178 int type = hypre_SStructGraphObjectType(graph);
179
180 hypre_SStructUVEntry **Uventries = hypre_SStructGraphUVEntries(graph);
181 hypre_SStructUVEntry *Uventry;
182 int nUentries;
183 hypre_SStructUEntry *Uentries;
184
185 hypre_BoxMapEntry *map_entry;
186 hypre_Index cindex;
187 HYPRE_BigInt big_rank, startrank;
188 int rank, i;
189 int box, to_box, to_proc;
190
191 if (!nUventries)
192 {
193 /* allocate space for non-stencil entries GEC1102
194 * the size equal to the ghost local size of grid */
195 aUventries = hypre_SStructGridGhlocalSize(grid);
196 iUventries = hypre_TAlloc(int, aUventries);
197
198 Uventries = hypre_CTAlloc(hypre_SStructUVEntry *, aUventries);
199 hypre_SStructGraphAUVEntries(graph) = aUventries;
200 hypre_SStructGraphIUVEntries(graph) = iUventries;
201
202 hypre_SStructGraphUVEntries(graph) = Uventries;
203 }
204 else if (nUventries >= aUventries)
205 {
206 /* need more space in iUventries array */
207 aUventries += 1000;
208 iUventries = hypre_TReAlloc(iUventries, int, aUventries);
209 hypre_SStructGraphAUVEntries(graph) = aUventries;
210 hypre_SStructGraphIUVEntries(graph) = iUventries;
211
212 }
213
214 /* compute location (rank) for Uventry */
215 hypre_CopyToCleanIndex(index, ndim, cindex);
216 hypre_SStructGridFindMapEntry(grid, part, cindex, var, &map_entry);
217
218 /* GEC0203 getting the rank */
219 hypre_SStructMapEntryGetGlobalRank(map_entry, cindex, &big_rank, type);
220
221 /* GEC 0902 filling up the iUventries with local ghrank
222 * since HYPRE_SSTRUCT is chosen */
223
224 if (type == HYPRE_SSTRUCT || type == HYPRE_STRUCT)
225 {
226 startrank = hypre_SStructGridGhstartRank(grid);
227 }
228 if (type == HYPRE_PARCSR)
229 {
230 startrank = hypre_SStructGridStartRank(grid);
231 }
232
233 rank = (int)(big_rank-startrank);
234
235 iUventries[nUventries] = rank;
236
237 if (Uventries[rank] == NULL)
238 {
239 Uventry = hypre_TAlloc(hypre_SStructUVEntry, 1);
240 hypre_SStructUVEntryPart(Uventry) = part;
241 hypre_CopyToCleanIndex(index, ndim, hypre_SStructUVEntryIndex(Uventry));
242 hypre_SStructUVEntryVar(Uventry) = var;
243 hypre_SStructMapEntryGetBox(map_entry, &box);
244 hypre_SStructUVEntryBox(Uventry)= box;
245 nUentries = 1;
246 Uentries = hypre_TAlloc(hypre_SStructUEntry, nUentries);
247 }
248 else
249 {
250 Uventry = Uventries[rank];
251 nUentries = hypre_SStructUVEntryNUEntries(Uventry) + 1;
252 Uentries = hypre_SStructUVEntryUEntries(Uventry);
253 Uentries = hypre_TReAlloc(Uentries, hypre_SStructUEntry, nUentries);
254 }
255 hypre_SStructUVEntryNUEntries(Uventry) = nUentries;
256 hypre_SStructUVEntryUEntries(Uventry) = Uentries;
257
258 i = nUentries - 1;
259 hypre_SStructUVEntryToPart(Uventry, i) = to_part;
260 hypre_CopyToCleanIndex(to_index, ndim,
261 hypre_SStructUVEntryToIndex(Uventry, i));
262 hypre_SStructUVEntryToVar(Uventry, i) = to_var;
263
264 hypre_CopyToCleanIndex(to_index, ndim, cindex);
265 hypre_SStructGridFindMapEntry(grid, to_part, cindex, to_var, &map_entry);
266 hypre_SStructMapEntryGetBox(map_entry, &to_box);
267 hypre_SStructUVEntryToBox(Uventry, i)= to_box;
268 hypre_SStructMapEntryGetProcess(map_entry, &to_proc);
269 hypre_SStructUVEntryToProc(Uventry, i)= to_proc;
270
271 Uventries[rank] = Uventry; /* GEC1102 where rank labels Uventries */
272
273 hypre_SStructGraphNUVEntries(graph) ++;
274 hypre_SStructGraphUVEntries(graph) = Uventries;
275 hypre_SStructGraphTotUEntries(graph) ++;
276
277 return ierr;
278}
279
280/*--------------------------------------------------------------------------
281 * HYPRE_SStructGraphAssemble
282 *
283 * This routine mainly computes the column numbers for the non-stencil
284 * graph entries (i.e., those created by GraphAddEntries calls). The
285 * routine will compute as many of these numbers on-process, but if
286 * the information needed to compute a column is not stored locally,
287 * it will be computed off-process instead.
288 *
289 * To do this, column info is first requested from other processes
290 * (tag=1 communications). While waiting for this info, requests from
291 * other processes are filled (tag=2). Simultaneously, a fanin/fanout
292 * procedure (tag=0) is invoked to determine when to stop: each
293 * process participates in the send portion of the fanin once it has
294 * received all of its requested column data and once it has completed
295 * its receive portion of the fanin; each process then participates in
296 * the fanout.
297 *
298 *--------------------------------------------------------------------------*/
299
300int
301HYPRE_SStructGraphAssemble( HYPRE_SStructGraph graph )
302{
303 int ierr = 0;
304
305 MPI_Comm comm = hypre_SStructGraphComm(graph);
306 hypre_SStructGrid *grid = hypre_SStructGraphGrid(graph);
307 int nUventries = hypre_SStructGraphNUVEntries(graph);
308 int *iUventries = hypre_SStructGraphIUVEntries(graph);
309 hypre_SStructUVEntry **Uventries = hypre_SStructGraphUVEntries(graph);
310 int totUentries = hypre_SStructGraphTotUEntries(graph);
311
312 int type = hypre_SStructGraphObjectType(graph);
313
314
315 hypre_SStructUVEntry *Uventry;
316 hypre_SStructUEntry *Uentry;
317 int nUentries;
318 int to_part;
319 hypre_IndexRef to_index;
320 int to_var;
321 int to_box;
322 int to_proc;
323 int proc;
324 HYPRE_BigInt rank;
325 hypre_BoxMapEntry *map_entry;
326
327 /* type 0 communications used to determine completion (NULL messages) */
328 MPI_Request *t0requests;
329 MPI_Status *t0statuses;
330 int *t0hiprocs;
331 int t0nhi;
332 int t0loproc;
333
334 /* type 1 communications used to get column data from other processes */
335 MPI_Request *t1requests, tmprequest;
336 MPI_Status *t1statuses;
337 int **t1sendbufs;
338 int **t1recvbufs;
339 int *t1bufsizes;
340 int *t1bufprocs;
341 hypre_SStructUEntry ***t1Uentries;
342 int t1totsize;
343 int t1ncomms;
344 int t1complete;
345 struct UentryLink
346 {
347 hypre_SStructUEntry *uentry;
348 struct UentryLink *next;
349 };
350 struct UentryLink *t1links;
351 struct UentryLink *t1linksptr;
352 struct UentryLink **t1lists;
353
354 /* type 2 communications used to send column data to other processes */
355 MPI_Status t2status;
356 int *t2commbuf;
357 int t2bufsize, t2maxbufsize;
358 int t2flag;
359
360 int nprocs, myproc;
361 int fanin_complete, complete;
362 int i, j;
363
364 MPI_Comm_size(comm, &nprocs);
365 MPI_Comm_rank(comm, &myproc);
366
367 /*---------------------------------------------------------
368 * Sort the iUventries array and eliminate duplicates.
369 *---------------------------------------------------------*/
370
371 if (nUventries > 1)
372 {
373 qsort0(iUventries, 0, nUventries - 1);
374
375 j = 1;
376 for (i = 1; i < nUventries; i++)
377 {
378 if (iUventries[i] > iUventries[i-1])
379 {
380 iUventries[j] = iUventries[i];
381 j++;
382 }
383 }
384 nUventries = j;
385 hypre_SStructGraphNUVEntries(graph) = nUventries;
386 }
387 /*---------------------------------------------------------
388 * Compute non-stencil column numbers (if possible), and
389 * start building requests for needed off-process info.
390 *---------------------------------------------------------*/
391
392 t1totsize = 0;
393 t1ncomms = 0;
394 for (i = 0; i < nUventries; i++)
395 {
396 Uventry = Uventries[iUventries[i]];
397 nUentries = hypre_SStructUVEntryNUEntries(Uventry);
398 for (j = 0; j < nUentries; j++)
399 {
400 to_part = hypre_SStructUVEntryToPart(Uventry, j);
401 to_index = hypre_SStructUVEntryToIndex(Uventry, j);
402 to_var = hypre_SStructUVEntryToVar(Uventry, j);
403 /*
404 hypre_SStructGridFindMapEntry(grid, to_part, to_index, to_var,
405 &map_entry);
406 */
407
408 /*---------------------------------------------------------
409 * used in future? The to_box corresponds to the first
410 * map_entry on the map_entry link list.
411 *---------------------------------------------------------*/
412 to_box = hypre_SStructUVEntryToBox(Uventry, j);
413 to_proc = hypre_SStructUVEntryToProc(Uventry, j);
414 hypre_SStructGridBoxProcFindMapEntry(grid, to_part, to_var, to_box,
415 to_proc, &map_entry);
416 if (map_entry != NULL)
417 {
418 /* compute ranks locally */
419
420 hypre_SStructMapEntryGetGlobalRank(map_entry, to_index, &rank, type);
421 hypre_SStructUVEntryRank(Uventry, j) = rank;
422
423 }
424 else
425 {
426 printf("off process computation\n");
427
428 /* compute ranks off-process: start building type 1 requests */
429 hypre_SStructMapEntryGetProcess(map_entry, &proc);
430
431 /* initialize some things */
432 if (t1totsize == 0)
433 {
434 t1links = hypre_TAlloc(struct UentryLink, totUentries);
435 t1linksptr = t1links;
436 t1lists = hypre_CTAlloc(struct UentryLink *, nprocs);
437 t1bufprocs = hypre_TAlloc(int, 10);
438 }
439
440 /* add a new process to the requests list */
441 if (t1lists[proc] == NULL)
442 {
443 if (t1ncomms%10 == 0)
444 {
445 t1bufprocs = hypre_TReAlloc(t1bufprocs, int, (t1ncomms+10));
446 }
447 t1bufprocs[t1ncomms] = proc;
448 t1ncomms++;
449 }
450
451 /* set up the new link */
452 (t1linksptr -> uentry) = hypre_SStructUVEntryUEntry(Uventry, j);
453 (t1linksptr -> next) = t1lists[proc];
454
455 /* insert link into the appropriate list */
456 t1lists[proc] = t1linksptr;
457 t1linksptr++;
458
459 t1totsize++;
460 }
461 }
462 }
463
464 /* set up remaining type 1 request info */
465
466 if (t1ncomms > 0)
467 {
468 t1bufsizes = hypre_TAlloc(int, t1ncomms);
469 t1Uentries = hypre_TAlloc(hypre_SStructUEntry **, t1ncomms);
470 t1Uentries[0] = hypre_TAlloc(hypre_SStructUEntry *, t1totsize);
471 j = 0;
472 for (i = 0; i < t1ncomms; i++)
473 {
474 proc = t1bufprocs[i];
475 t1linksptr = t1lists[proc];
476 t1bufsizes[i] = 0;
477 while (t1linksptr != NULL)
478 {
479 t1Uentries[0][j++] = (t1linksptr -> uentry);
480 t1linksptr = (t1linksptr -> next);
481 t1bufsizes[i]++;
482 }
483 }
484 hypre_TFree(t1links);
485 hypre_TFree(t1lists);
486
487 t1sendbufs = hypre_TAlloc(int *, t1ncomms);
488 t1sendbufs[0] = hypre_TAlloc(int, t1totsize*6);
489 t1recvbufs = hypre_TAlloc(int *, t1ncomms);
490 t1recvbufs[0] = hypre_TAlloc(int, t1totsize);
491
492 for (j = 0; j < t1totsize; j++)
493 {
494 Uentry = t1Uentries[0][j];
495 to_part = hypre_SStructUEntryToPart(Uentry);
496 to_index = hypre_SStructUEntryToIndex(Uentry);
497 to_var = hypre_SStructUEntryToVar(Uentry);
498 to_box = hypre_SStructUEntryToBox(Uentry);
499 t1sendbufs[0][6*j ] = to_part;
500 t1sendbufs[0][6*j+1] = hypre_IndexD(to_index, 0);
501 t1sendbufs[0][6*j+2] = hypre_IndexD(to_index, 1);
502 t1sendbufs[0][6*j+3] = hypre_IndexD(to_index, 2);
503 t1sendbufs[0][6*j+4] = to_var;
504 t1sendbufs[0][6*j+5] = to_box;
505 }
506
507 /* GEC1002 commenting this out to replace it by something else
508 * since I think that a 6 is missing in sending buffer
509 *
510 * for (i = 1; i < t1ncomms; i++)
511 * {
512 * t1sendbufs[i] = t1sendbufs[i-1] + t1bufsizes[i-1];
513 * t1recvbufs[i] = t1sendbufs[i-1] + t1bufsizes[i-1];
514 * t1Uentries[i] = t1Uentries[i-1] + t1bufsizes[i-1];
515 * } */
516
517 for (i = 1; i < t1ncomms; i++)
518 {
519 t1sendbufs[i] = t1sendbufs[i-1] + 6*t1bufsizes[i-1];
520 t1recvbufs[i] = t1sendbufs[i-1] + t1bufsizes[i-1];
521 t1Uentries[i] = t1Uentries[i-1] + t1bufsizes[i-1];
522 }
523 }
524
525 /*---------------------------------------------------------
526 * Exchange column data with other processes
527 *---------------------------------------------------------*/
528
529 t1complete = 1;
530 if (t1ncomms > 0)
531 {
532 t1requests = hypre_CTAlloc(MPI_Request, t1ncomms);
533 t1statuses = hypre_CTAlloc(MPI_Status, t1ncomms);
534
535 /* post type 1 requests (note: tag=1 on send; tag=2 on receive) */
536 for (i = 0; i < t1ncomms; i++)
537 {
538
539 MPI_Irecv(t1recvbufs[i], t1bufsizes[i], MPI_INT, t1bufprocs[i],
540 2, comm, &t1requests[i]);
541 }
542 for (i = 0; i < t1ncomms; i++)
543 {
544 /* Note: No need to check below that these have completed */
545 MPI_Isend(t1sendbufs[i], t1bufsizes[i]*6, MPI_INT, t1bufprocs[i],
546 1, comm, &tmprequest);
547 }
548
549 t1complete = 0;
550 }
551
552 complete = 1;
553 if (nprocs > 1)
554 {
555 /* compute type 0 communication info */
556 for (i = 1, t0nhi = 0; i < nprocs; i *= 2, t0nhi++);
557 t0hiprocs = hypre_TAlloc(int, t0nhi);
558 t0nhi = 0;
559 proc = myproc;
560 for (i = 1; i < nprocs; i *= 2)
561 {
562 if ((proc % 2) == 0)
563 {
564 /* added if test to avoid t0hi ranks exceeding nprocs */
565 if ((myproc + i) < nprocs)
566 {
567 t0hiprocs[t0nhi] = myproc + i;
568 t0nhi++;
569 }
570 proc /= 2;
571 }
572 else
573 {
574 t0loproc = myproc - i;
575 break;
576 }
577 }
578
579 t0requests = hypre_CTAlloc(MPI_Request, (t0nhi + 1));
580 t0statuses = hypre_CTAlloc(MPI_Status, (t0nhi + 1));
581
582 /* post type 0 fanin receives */
583 for (i = 0; i < t0nhi; i++)
584 {
585 MPI_Irecv(NULL, 0, MPI_INT, t0hiprocs[i], 0, comm, &t0requests[i]);
586 }
587
588 t2commbuf = NULL;
589 t2maxbufsize = 0;
590
591 fanin_complete = 0;
592 complete = 0;
593 }
594
595 while (!complete)
596 {
597 /* fill as many column data requests as possible */
598 /* (note: tag=1 on receive; tag=2 on send) */
599 MPI_Iprobe(MPI_ANY_SOURCE, 1, comm, &t2flag, &t2status);
600 while (t2flag)
601 {
602 proc = t2status.MPI_SOURCE;
603
604 MPI_Get_count(&t2status, MPI_INT, &t2bufsize);
605 if (t2bufsize > t2maxbufsize)
606 {
607 t2commbuf = hypre_TReAlloc(t2commbuf, int, t2bufsize);
608 t2maxbufsize = t2bufsize;
609 }
610 MPI_Recv(t2commbuf, t2bufsize, MPI_INT, proc, 1, comm, &t2status);
611
612 t2bufsize /= 6;
613 for (j = 0; j < t2bufsize; j++)
614 {
615 to_part = t2commbuf[6*j];
616 hypre_IndexD(to_index, 0) = t2commbuf[6*j+1];
617 hypre_IndexD(to_index, 1) = t2commbuf[6*j+2];
618 hypre_IndexD(to_index, 2) = t2commbuf[6*j+3];
619 to_var = t2commbuf[6*j+4];
620 to_box = t2commbuf[6*j+5]; /* future use? */
621 hypre_SStructGridFindMapEntry(grid, to_part, to_index, to_var,
622 &map_entry);
623
624 hypre_SStructMapEntryGetGlobalRank(map_entry, to_index, &rank, type);
625
626 t2commbuf[j] = rank;
627 }
628 MPI_Send(t2commbuf, t2bufsize, MPI_INT, proc, 2, comm);
629
630 MPI_Iprobe(MPI_ANY_SOURCE, 1, comm, &t2flag, &t2status);
631 }
632
633 /* attempt to complete type 0 and type 1 communications */
634 if (!t1complete)
635 {
636 MPI_Testall(t1ncomms, t1requests, &t1complete, t1statuses);
637 }
638 else if (!fanin_complete)
639 {
640 MPI_Testall(t0nhi, t0requests, &fanin_complete, t0statuses);
641 if (fanin_complete)
642 {
643 if (myproc > 0)
644 {
645 /* post type 0 fanin send */
646 MPI_Send(NULL, 0, MPI_INT, t0loproc, 0, comm);
647
648 /* post type 0 fanout receive */
649 MPI_Irecv(NULL, 0, MPI_INT, t0loproc, 0, comm, t0requests);
650 }
651 }
652 }
653 else
654 {
655 MPI_Test(t0requests, &complete, t0statuses);
656 if (complete)
657 {
658 /* post type 0 fanout sends */
659 for (i = 0; i < t0nhi; i++)
660 {
661 MPI_Send(NULL, 0, MPI_INT, t0hiprocs[i], 0, comm);
662 }
663 }
664 }
665 }
666
667 /* copy column data computed off-process into local data structure */
668 if (t1ncomms > 0)
669 {
670 for (j = 0; j < t1totsize; j++)
671 {
672 Uentry = t1Uentries[0][j];
673 hypre_SStructUEntryRank(Uentry) = t1recvbufs[0][j];
674 }
675 }
676
677 if (nprocs > 1)
678 {
679 hypre_TFree(t0requests);
680 hypre_TFree(t0statuses);
681 hypre_TFree(t0hiprocs);
682
683 hypre_TFree(t2commbuf);
684 }
685
686 if (t1ncomms > 0)
687 {
688 hypre_TFree(t1requests);
689 hypre_TFree(t1statuses);
690 hypre_TFree(t1sendbufs[0]);
691 hypre_TFree(t1sendbufs);
692 hypre_TFree(t1recvbufs[0]);
693 hypre_TFree(t1recvbufs);
694 hypre_TFree(t1bufsizes);
695 hypre_TFree(t1bufprocs);
696 hypre_TFree(t1Uentries[0]);
697 hypre_TFree(t1Uentries);
698 }
699
700 return ierr;
701}
702/*****************************************************************
703 *
704 *
705 *****************************************************************/
706int HYPRE_SStructGraphSetObjectType(HYPRE_SStructGraph graph,
707 int type)
708{
709 int ierr = 0;
710 hypre_SStructGraphObjectType(graph) = type;
711 return ierr;
712}
Note: See TracBrowser for help on using the repository browser.