Bug Summary

File:vec/is/sf/impls/basic/neighbor/sfneighbor.c
Warning:line 285, column 29
Array access (via field 'leaf') results in a null pointer dereference

Annotated Source Code

[?] Use j/k keys for keyboard navigation

1#include <../src/vec/is/sf/impls/basic/sfpack.h>
2#include <../src/vec/is/sf/impls/basic/sfbasic.h>
3
4#if defined(PETSC_HAVE_MPI_NEIGHBORHOOD_COLLECTIVES1)
5
6/* SF Neighbor completely relies on the two sided info built by SF Basic. Therefore we build Neighbor as a subclass of Basic */
7
8typedef struct {
9 SPPACKBASICHEADERPetscErrorCode (*Pack) (PetscInt,PetscInt,const PetscInt*,PetscInt
,PetscSFPackOpt,const void*,void*); PetscErrorCode (*UnpackAndInsert
)(PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void
*,const void*); PetscErrorCode (*UnpackAndAdd) (PetscInt,PetscInt
,const PetscInt*,PetscInt,PetscSFPackOpt,void*,const void*); PetscErrorCode
(*UnpackAndMin) (PetscInt,PetscInt,const PetscInt*,PetscInt,
PetscSFPackOpt,void*,const void*); PetscErrorCode (*UnpackAndMax
) (PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,
void*,const void*); PetscErrorCode (*UnpackAndMinloc)(PetscInt
,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void*,const
void*); PetscErrorCode (*UnpackAndMaxloc)(PetscInt,PetscInt,
const PetscInt*,PetscInt,PetscSFPackOpt,void*,const void*); PetscErrorCode
(*UnpackAndMult) (PetscInt,PetscInt,const PetscInt*,PetscInt
,PetscSFPackOpt,void*,const void*); PetscErrorCode (*UnpackAndLAND
) (PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,
void*,const void*); PetscErrorCode (*UnpackAndBAND) (PetscInt
,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void*,const
void*); PetscErrorCode (*UnpackAndLOR) (PetscInt,PetscInt,const
PetscInt*,PetscInt,PetscSFPackOpt,void*,const void*); PetscErrorCode
(*UnpackAndBOR) (PetscInt,PetscInt,const PetscInt*,PetscInt,
PetscSFPackOpt,void*,const void*); PetscErrorCode (*UnpackAndLXOR
) (PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,
void*,const void*); PetscErrorCode (*UnpackAndBXOR) (PetscInt
,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void*,const
void*); PetscErrorCode (*FetchAndInsert) (PetscInt,PetscInt,
const PetscInt*,PetscInt,PetscSFPackOpt,void*,void*); PetscErrorCode
(*FetchAndAdd) (PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt
,void*,void*); PetscErrorCode (*FetchAndMin) (PetscInt,PetscInt
,const PetscInt*,PetscInt,PetscSFPackOpt,void*,void*); PetscErrorCode
(*FetchAndMax) (PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt
,void*,void*); PetscErrorCode (*FetchAndMinloc) (PetscInt,PetscInt
,const PetscInt*,PetscInt,PetscSFPackOpt,void*,void*); PetscErrorCode
(*FetchAndMaxloc) (PetscInt,PetscInt,const PetscInt*,PetscInt
,PetscSFPackOpt,void*,void*); PetscErrorCode (*FetchAndMult) (
PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void
*,void*); PetscErrorCode (*FetchAndLAND) (PetscInt,PetscInt,const
PetscInt*,PetscInt,PetscSFPackOpt,void*,void*); PetscErrorCode
(*FetchAndBAND) (PetscInt,PetscInt,const PetscInt*,PetscInt,
PetscSFPackOpt,void*,void*); PetscErrorCode (*FetchAndLOR) (PetscInt
,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void*,void*
); PetscErrorCode (*FetchAndBOR) (PetscInt,PetscInt,const PetscInt
*,PetscInt,PetscSFPackOpt,void*,void*); PetscErrorCode (*FetchAndLXOR
) (PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,
void*,void*); PetscErrorCode (*FetchAndBXOR) (PetscInt,PetscInt
,const PetscInt*,PetscInt,PetscSFPackOpt,void*,void*); PetscMPIInt
tag; MPI_Datatype unit; MPI_Datatype basicunit; PetscBool isbuiltin
; size_t unitbytes; PetscInt bs; const void *key; PetscSFPack
next; char **root; char **leaf; PetscMPIInt half; MPI_Request
*requests
;
10 char *rootbuf; /* contiguous buffer for all root ranks */
11 char *leafbuf; /* contiguous buffer for all non-distinguished leaf ranks. Distiguished ones share root buffers. */
12} *PetscSFPack_Neighbor;
13
14typedef struct {
15 SFBASICHEADERPetscMPIInt niranks; PetscMPIInt ndiranks; PetscMPIInt *iranks
; PetscInt itotal; PetscInt *ioffset; PetscInt *irootloc; PetscSFPackOpt
rootpackopt; PetscSFPack avail; PetscSFPack inuse
;
16 MPI_Comm comms[2]; /* Communicators with distributed topology in both directions */
17 PetscBool initialized[2]; /* Are the two communicators initialized? */
18 PetscMPIInt *rootdispls,*rootcounts,*leafdispls,*leafcounts; /* displs/counts for non-distinguished ranks */
19} PetscSF_Neighbor;
20
21/*===================================================================================*/
22/* Internal utility routines */
23/*===================================================================================*/
24
25/* Get the communicator with distributed graph topology, which is not cheap to build so we do it on demand (instead of at PetscSFSetUp time) */
26static PetscErrorCode PetscSFGetDistComm_Neighbor(PetscSF sf,PetscSFDirection direction,MPI_Comm *distcomm)
27{
28 PetscErrorCode ierr;
29 PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
30 PetscInt nrootranks,ndrootranks,nleafranks,ndleafranks;
31 const PetscMPIInt *rootranks,*leafranks;
32 MPI_Comm comm;
33
34 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 34; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
35 ierr = PetscSFGetRootInfo_Basic(sf,&nrootranks,&ndrootranks,&rootranks,NULL((void*)0),NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),35,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
; /* Which ranks will access my roots (I am a destination) */
36 ierr = PetscSFGetLeafInfo_Basic(sf,&nleafranks,&ndleafranks,&leafranks,NULL((void*)0),NULL((void*)0),NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),36,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
; /* My leaves will access whose roots (I am a source) */
37
38 if (!dat->initialized[direction]) {
39 const PetscMPIInt indegree = nrootranks-ndrootranks,*sources = rootranks+ndrootranks;
40 const PetscMPIInt outdegree = nleafranks-ndleafranks,*destinations = leafranks+ndleafranks;
41 MPI_Comm *mycomm = &dat->comms[direction];
42 ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),42,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
43 if (direction == PETSCSF_LEAF2ROOT_REDUCE) {
44 ierr = MPI_Dist_graph_create_adjacent(comm,indegree,sources,dat->rootcounts/*src weights*/,outdegree,destinations,dat->leafcounts/*dest weights*/,MPI_INFO_NULL((MPI_Info)0x1c000000),1/*reorder*/,mycomm);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),44,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
45 } else { /* PETSCSF_ROOT2LEAF_BCAST, reverse src & dest */
46 ierr = MPI_Dist_graph_create_adjacent(comm,outdegree,destinations,dat->leafcounts/*src weights*/,indegree,sources,dat->rootcounts/*dest weights*/,MPI_INFO_NULL((MPI_Info)0x1c000000),1/*reorder*/,mycomm);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),46,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
47 }
48 dat->initialized[direction] = PETSC_TRUE;
49 }
50 *distcomm = dat->comms[direction];
51 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
52}
53
54static PetscErrorCode PetscSFPackGet_Neighbor(PetscSF sf,MPI_Datatype unit,const void *key,PetscSFDirection direction,PetscSFPack_Neighbor *mylink)
55{
56 PetscErrorCode ierr;
57 PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
58 PetscInt i,nrootranks,ndrootranks,nleafranks,ndleafranks;
59 const PetscInt *rootoffset,*leafoffset;
60 PetscSFPack *p;
61 PetscSFPack_Neighbor link;
62
63 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 63; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
64 ierr = PetscSFGetRootInfo_Basic(sf,&nrootranks,&ndrootranks,NULL((void*)0),&rootoffset,NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),64,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
65 ierr = PetscSFGetLeafInfo_Basic(sf,&nleafranks,&ndleafranks,NULL((void*)0),&leafoffset,NULL((void*)0),NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),65,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
66
67 /* Look for types in cache */
68 for (p=&dat->avail; (link=(PetscSFPack_Neighbor)*p); p=&link->next) {
69 PetscBool match;
70 ierr = MPIPetsc_Type_compare(unit,link->unit,&match);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),70,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
71 if (match) {
72 *p = link->next; /* Remove from available list */
73 goto found;
74 }
75 }
76
77 ierr = PetscNew(&link)PetscMallocA(1,PETSC_TRUE,77,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,(size_t)(1)*sizeof(**((&link))),((&link)))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),77,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
78 ierr = PetscSFPackSetupType((PetscSFPack)link,unit);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),78,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
79 ierr = PetscMalloc2(nrootranks,&link->root,nleafranks,&link->leaf)PetscMallocA(2,PETSC_FALSE,79,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,(size_t)(nrootranks)*sizeof(**(&link->root)),(&link
->root),(size_t)(nleafranks)*sizeof(**(&link->leaf)
),(&link->leaf))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),79,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
80 /* Double the requests. First half are used for reduce (leaf2root) communication, second half for bcast (root2leaf) communication */
81 link->half = 1;
82 ierr = PetscMalloc1(link->half*2,&link->requests)PetscMallocA(1,PETSC_FALSE,82,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,(size_t)(link->half*2)*sizeof(**(&link->requests))
,(&link->requests))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),82,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
83 ierr = PetscCommGetNewTag(PetscObjectComm((PetscObject)sf),&link->tag);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),83,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
; /* Actually, tag is not need for neighborhood collectives */
84
85 /* Allocate root and leaf buffers */
86 ierr = PetscMalloc2(rootoffset[nrootranks]*link->unitbytes,&link->rootbuf,(leafoffset[nleafranks]-leafoffset[ndleafranks])*link->unitbytes,&link->leafbuf)PetscMallocA(2,PETSC_FALSE,86,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,(size_t)(rootoffset[nrootranks]*link->unitbytes)*sizeof(*
*(&link->rootbuf)),(&link->rootbuf),(size_t)((leafoffset
[nleafranks]-leafoffset[ndleafranks])*link->unitbytes)*sizeof
(**(&link->leafbuf)),(&link->leafbuf))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),86,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
87 for (i=0; i<nrootranks; i++) link->root[i] = link->rootbuf + rootoffset[i]*link->unitbytes;
88 for (i=0; i<nleafranks; i++) {
89 if (i < ndleafranks) { /* Leaf buffers for distinguished ranks are pointers directly into root buffers */
90 if (ndrootranks != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cannot match distinguished ranks")return PetscError(((MPI_Comm)0x44000001),90,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,77,PETSC_ERROR_INITIAL,"Cannot match distinguished ranks")
;
91 link->leaf[i] = link->root[0];
92 continue;
93 }
94 link->leaf[i] = link->leafbuf + (leafoffset[i] - leafoffset[ndleafranks])*link->unitbytes;
95 }
96
97found:
98 link->key = key;
99 link->next = dat->inuse;
100 dat->inuse = (PetscSFPack)link;
101 *mylink = link;
102 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
103}
104
105/*===================================================================================*/
106/* Implementations of SF public APIs */
107/*===================================================================================*/
108static PetscErrorCode PetscSFSetUp_Neighbor(PetscSF sf)
109{
110 PetscErrorCode ierr;
111 PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
112 PetscInt i,j,nrootranks,ndrootranks,nleafranks,ndleafranks;
113 const PetscInt *rootoffset,*leafoffset;
114 PetscMPIInt m,n;
115
116 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 116; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
117 ierr = PetscSFSetUp_Basic(sf);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),117,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
118 ierr = PetscSFGetRootInfo_Basic(sf,&nrootranks,&ndrootranks,NULL((void*)0),&rootoffset,NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),118,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
119 ierr = PetscSFGetLeafInfo_Basic(sf,&nleafranks,&ndleafranks,NULL((void*)0),&leafoffset,NULL((void*)0),NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),119,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
120
121 /* Only setup MPI displs/counts for non-distinguished ranks. Distinguished ranks use shared memory */
122 ierr = PetscMalloc4(nrootranks-ndrootranks,&dat->rootdispls,nrootranks-ndrootranks,&dat->rootcounts,nleafranks-ndleafranks,&dat->leafdispls,nleafranks-ndleafranks,&dat->leafcounts)PetscMallocA(4,PETSC_FALSE,122,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,(size_t)(nrootranks-ndrootranks)*sizeof(**(&dat->rootdispls
)),(&dat->rootdispls),(size_t)(nrootranks-ndrootranks)
*sizeof(**(&dat->rootcounts)),(&dat->rootcounts
),(size_t)(nleafranks-ndleafranks)*sizeof(**(&dat->leafdispls
)),(&dat->leafdispls),(size_t)(nleafranks-ndleafranks)
*sizeof(**(&dat->leafcounts)),(&dat->leafcounts
))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),122,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
123 for (i=ndrootranks,j=0; i<nrootranks; i++,j++) {
124 ierr = PetscMPIIntCast(rootoffset[i]-rootoffset[ndrootranks],&m);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),124,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
; dat->rootdispls[j] = m;
125 ierr = PetscMPIIntCast(rootoffset[i+1]-rootoffset[i], &n);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),125,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
; dat->rootcounts[j] = n;
126 }
127
128 for (i=ndleafranks,j=0; i<nleafranks; i++,j++) {
129 ierr = PetscMPIIntCast(leafoffset[i]-leafoffset[ndleafranks],&m);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),129,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
; dat->leafdispls[j] = m;
130 ierr = PetscMPIIntCast(leafoffset[i+1]-leafoffset[i], &n);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),130,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
; dat->leafcounts[j] = n;
131 }
132 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
133}
134
135static PetscErrorCode PetscSFReset_Neighbor(PetscSF sf)
136{
137 PetscErrorCode ierr;
138 PetscInt i;
139 PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
140 PetscSFPack_Neighbor link,next;
141
142 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 142; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
143 if (dat->inuse) SETERRQ(PetscObjectComm((PetscObject)sf),PETSC_ERR_ARG_WRONGSTATE,"Outstanding operation has not been completed")return PetscError(PetscObjectComm((PetscObject)sf),143,__func__
,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,73,PETSC_ERROR_INITIAL,"Outstanding operation has not been completed"
)
;
144 ierr = PetscFree4(dat->rootdispls,dat->rootcounts,dat->leafdispls,dat->leafcounts)PetscFreeA(4,144,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,&(dat->rootdispls),&(dat->rootcounts),&(dat
->leafdispls),&(dat->leafcounts))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),144,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
145 for (i=0; i<2; i++) {
146 if (dat->initialized[i]) {
147 ierr = MPI_Comm_free(&dat->comms[i]);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),147,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
148 dat->initialized[i] = PETSC_FALSE;
149 }
150 }
151
152 ierr = PetscFree2(dat->iranks,dat->ioffset)PetscFreeA(2,152,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,&(dat->iranks),&(dat->ioffset))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),152,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
153 ierr = PetscFree(dat->irootloc)((*PetscTrFree)((void*)(dat->irootloc),153,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
) || ((dat->irootloc) = 0,0))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),153,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
154 for (link=(PetscSFPack_Neighbor)dat->avail; link; link=next) {
155 next = (PetscSFPack_Neighbor)link->next;
156 if (!link->isbuiltin) {ierr = MPI_Type_free(&link->unit);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),156,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;}
157 ierr = PetscFree2(link->root,link->leaf)PetscFreeA(2,157,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,&(link->root),&(link->leaf))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),157,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
158 ierr = PetscFree2(link->rootbuf,link->leafbuf)PetscFreeA(2,158,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,&(link->rootbuf),&(link->leafbuf))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),158,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
159 ierr = PetscFree(link->requests)((*PetscTrFree)((void*)(link->requests),159,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
) || ((link->requests) = 0,0))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),159,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
160 ierr = PetscFree(link)((*PetscTrFree)((void*)(link),160,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
) || ((link) = 0,0))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),160,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
161 }
162 dat->avail = NULL((void*)0);
163 ierr = PetscSFPackDestoryOptimization(&sf->leafpackopt);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),163,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
164 ierr = PetscSFPackDestoryOptimization(&dat->rootpackopt);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),164,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
165 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
166}
167
168static PetscErrorCode PetscSFDestroy_Neighbor(PetscSF sf)
169{
170 PetscErrorCode ierr;
171
172 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 172; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
173 ierr = PetscSFReset_Neighbor(sf);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),173,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
174 ierr = PetscFree(sf->data)((*PetscTrFree)((void*)(sf->data),174,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
) || ((sf->data) = 0,0))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),174,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
175 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
176}
177
178static PetscErrorCode PetscSFBcastAndOpBegin_Neighbor(PetscSF sf,MPI_Datatype unit,const void *rootdata,void *leafdata,MPI_Op op)
179{
180 PetscErrorCode ierr;
181 PetscSFPack_Neighbor link;
182 PetscInt i,nrootranks,ndrootranks,nleafranks,ndleafranks;
183 const PetscInt *rootoffset,*rootloc;
184 PetscMPIInt n,ind,outd;
185 PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
186 MPI_Comm distcomm;
187 MPI_Request *req;
188 void *sbuf,*rbuf;
189
190 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 190; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
191 ierr = PetscSFGetRootInfo_Basic(sf,&nrootranks,&ndrootranks,NULL((void*)0),&rootoffset,&rootloc);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),191,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
192 ierr = PetscSFGetLeafInfo_Basic(sf,&nleafranks,&ndleafranks,NULL((void*)0),NULL((void*)0),NULL((void*)0),NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),192,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
193 ierr = PetscSFPackGet_Neighbor(sf,unit,rootdata,PETSCSF_ROOT2LEAF_BCAST,&link);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),193,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
194
195 /* Pack root data */
196 for (i=0; i<nrootranks; i++) {
197 void *packstart = link->root[i];
198 ierr = PetscMPIIntCast(rootoffset[i+1]-rootoffset[i],&n);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),198,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
199 (*link->Pack)(n,link->bs,rootloc+rootoffset[i],i,dat->rootpackopt,rootdata,packstart);
200 }
201
202 /* Do neighborhood alltoallv for non-distinguished ranks */
203 req = &link->requests[PETSCSF_ROOT2LEAF_BCAST];
204 ierr = PetscSFGetDistComm_Neighbor(sf,PETSCSF_ROOT2LEAF_BCAST,&distcomm);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),204,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
205 outd = nrootranks - ndrootranks;
206 ind = nleafranks - ndleafranks;
207 sbuf = link->root ? link->root[ndrootranks] : NULL((void*)0);
208 rbuf = link->leaf ? link->leaf[ndleafranks] : NULL((void*)0);
209 ierr = MPI_Start_ineighbor_alltoallv(outd,ind,sbuf,dat->rootcounts,dat->rootdispls,unit,rbuf,dat->leafcounts,dat->leafdispls,unit,distcomm,req)((petsc_isend_ct += (PetscLogDouble)(outd),0) || (petsc_irecv_ct
+= (PetscLogDouble)(ind),0) || PetscMPITypeSizeCount((outd),
(dat->rootcounts),(unit),(&petsc_isend_len)) || PetscMPITypeSizeCount
((ind),(dat->leafcounts),(unit),(&petsc_irecv_len)) ||
(((outd) || (ind)) && MPI_Ineighbor_alltoallv((sbuf)
,(dat->rootcounts),(dat->rootdispls),(unit),(rbuf),(dat
->leafcounts),(dat->leafdispls),(unit),(distcomm),(req)
)))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),209,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
210 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
211}
212
213static PetscErrorCode PetscSFReduceBegin_Neighbor(PetscSF sf,MPI_Datatype unit,const void *leafdata,void *rootdata,MPI_Op op)
214{
215 PetscErrorCode ierr;
216 PetscInt i,nrootranks,ndrootranks,nleafranks,ndleafranks;
217 const PetscInt *leafoffset,*leafloc;
218 PetscMPIInt n,ind,outd;
219 PetscSFPack_Neighbor link;
220 PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
221 MPI_Comm distcomm;
222 MPI_Request *req;
223 void *sbuf,*rbuf;
224
225 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 225; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
226 ierr = PetscSFGetRootInfo_Basic(sf,&nrootranks,&ndrootranks,NULL((void*)0),NULL((void*)0),NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),226,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
227 ierr = PetscSFGetLeafInfo_Basic(sf,&nleafranks,&ndleafranks,NULL((void*)0),&leafoffset,&leafloc,NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),227,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
228 ierr = PetscSFPackGet_Neighbor(sf,unit,leafdata,PETSCSF_LEAF2ROOT_REDUCE,&link);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),228,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
229
230 /* Pack leaf data */
231 for (i=0; i<nleafranks; i++) {
232 void *packstart = link->leaf[i];
233 ierr = PetscMPIIntCast(leafoffset[i+1]-leafoffset[i],&n);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),233,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
234 (*link->Pack)(n,link->bs,leafloc+leafoffset[i],i,sf->leafpackopt,leafdata,packstart);
235 }
236
237 /* Do neighborhood alltoallv for non-distinguished ranks */
238 req = &link->requests[PETSCSF_LEAF2ROOT_REDUCE];
239 ierr = PetscSFGetDistComm_Neighbor(sf,PETSCSF_LEAF2ROOT_REDUCE,&distcomm);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),239,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
240 ind = nrootranks - ndrootranks;
241 outd = nleafranks - ndleafranks;
242 sbuf = link->leaf ? link->leaf[ndleafranks] : NULL((void*)0);
243 rbuf = link->root ? link->root[ndrootranks] : NULL((void*)0);
244 ierr = MPI_Start_ineighbor_alltoallv(outd,ind,sbuf,dat->leafcounts,dat->leafdispls,unit,rbuf,dat->rootcounts,dat->rootdispls,unit,distcomm,req)((petsc_isend_ct += (PetscLogDouble)(outd),0) || (petsc_irecv_ct
+= (PetscLogDouble)(ind),0) || PetscMPITypeSizeCount((outd),
(dat->leafcounts),(unit),(&petsc_isend_len)) || PetscMPITypeSizeCount
((ind),(dat->rootcounts),(unit),(&petsc_irecv_len)) ||
(((outd) || (ind)) && MPI_Ineighbor_alltoallv((sbuf)
,(dat->leafcounts),(dat->leafdispls),(unit),(rbuf),(dat
->rootcounts),(dat->rootdispls),(unit),(distcomm),(req)
)))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),244,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
245 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
246}
247
248static PetscErrorCode PetscSFFetchAndOpEnd_Neighbor(PetscSF sf,MPI_Datatype unit,void *rootdata,const void *leafdata,void *leafupdate,MPI_Op op)
249{
250 PetscErrorCode ierr;
251 PetscErrorCode (*FetchAndOp)(PetscInt,PetscInt,const PetscInt*,PetscInt,PetscSFPackOpt,void*,void*);
252 PetscSFPack_Basic link;
253 PetscInt i,nrootranks,ndrootranks,nleafranks,ndleafranks;
254 const PetscInt *rootoffset,*leafoffset,*rootloc,*leafloc;
255 const PetscMPIInt *rootranks,*leafranks;
256 MPI_Comm comm;
257 PetscMPIInt n,ind,outd;
258 PetscSF_Neighbor *dat = (PetscSF_Neighbor*)sf->data;
259 void *sbuf,*rbuf;
260
261 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 261; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
262 ierr = PetscSFPackGetInUse(sf,unit,leafdata,PETSC_OWN_POINTER,(PetscSFPack*)&link);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),262,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
263 ierr = PetscSFPackWaitall_Basic(link,PETSCSF_LEAF2ROOT_REDUCE);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),263,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
264
265 ierr = PetscSFGetRootInfo_Basic(sf,&nrootranks,&ndrootranks,&rootranks,&rootoffset,&rootloc);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),265,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
266 ierr = PetscSFGetLeafInfo_Basic(sf,&nleafranks,&ndleafranks,&leafranks,&leafoffset,&leafloc,NULL((void*)0));CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),266,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
267
268 /* Process local fetch-and-op */
269 ierr = PetscSFPackGetFetchAndOp(sf,(PetscSFPack)link,op,&FetchAndOp);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),269,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
1
Value assigned to field 'leaf'
270 for (i=0; i<nrootranks; i++) {
2
Assuming 'i' is >= 'nrootranks'
3
Loop condition is false. Execution continues on line 277
271 void *packstart = link->root[i];
272 ierr = PetscMPIIntCast(rootoffset[i+1]-rootoffset[i],&n);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),272,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
273 (*FetchAndOp)(n,link->bs,rootloc+rootoffset[i],i,dat->rootpackopt,rootdata,packstart);
274 }
275
276 /* Bcast the updated root buffer back to leaves */
277 ierr = PetscSFGetDistComm_Neighbor(sf,PETSCSF_ROOT2LEAF_BCAST,&comm);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),277,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
278 outd = nrootranks - ndrootranks;
279 ind = nleafranks - ndleafranks;
280 sbuf = link->root ? link->root[ndrootranks] : NULL((void*)0);
4
Assuming the condition is false
5
'?' condition is false
281 rbuf = link->leaf ? link->leaf[ndleafranks] : NULL((void*)0);
6
Assuming pointer value is null
7
'?' condition is false
282 ierr = MPI_Start_neighbor_alltoallv(outd,ind,sbuf,dat->rootcounts,dat->rootdispls,unit,rbuf,dat->leafcounts,dat->leafdispls,unit,comm)((petsc_isend_ct += (PetscLogDouble)(outd),0) || (petsc_irecv_ct
+= (PetscLogDouble)(ind),0) || PetscMPITypeSizeCount((outd),
(dat->rootcounts),(unit),(&petsc_isend_len)) || PetscMPITypeSizeCount
((ind),(dat->leafcounts),(unit),(&petsc_irecv_len)) ||
(((outd) || (ind)) && MPI_Neighbor_alltoallv((sbuf),
(dat->rootcounts),(dat->rootdispls),(unit),(rbuf),(dat->
leafcounts),(dat->leafdispls),(unit),(comm))))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),282,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
283
284 for (i=0; i<nleafranks; i++) {
8
Assuming 'i' is < 'nleafranks'
9
Loop condition is true. Entering loop body
285 const void *packstart = link->leaf[i];
10
Array access (via field 'leaf') results in a null pointer dereference
286 ierr = PetscMPIIntCast(leafoffset[i+1]-leafoffset[i],&n);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),286,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
287 (*link->UnpackAndInsert)(n,link->bs,leafloc+leafoffset[i],i,sf->leafpackopt,leafupdate,packstart);
288 }
289
290 ierr = PetscSFPackReclaim(sf,(PetscSFPack*)&link);CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),290,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
291 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
292}
293
294PETSC_INTERNextern __attribute__((visibility ("hidden"))) PetscErrorCode PetscSFCreate_Neighbor(PetscSF sf)
295{
296 PetscErrorCode ierr;
297 PetscSF_Neighbor *dat;
298
299 PetscFunctionBegindo { do { ; if (petscstack && (petscstack->currentsize
< 64)) { petscstack->function[petscstack->currentsize
] = __func__; petscstack->file[petscstack->currentsize]
= "/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
; petscstack->line[petscstack->currentsize] = 299; petscstack
->petscroutine[petscstack->currentsize] = PETSC_TRUE; petscstack
->currentsize++; } if (petscstack) { petscstack->hotdepth
+= (PETSC_FALSE || petscstack->hotdepth); } ; } while (0)
; ; } while (0)
;
300 sf->ops->CreateEmbeddedSF = PetscSFCreateEmbeddedSF_Basic;
301 sf->ops->CreateEmbeddedLeafSF = PetscSFCreateEmbeddedLeafSF_Basic;
302 sf->ops->BcastAndOpEnd = PetscSFBcastAndOpEnd_Basic;
303 sf->ops->ReduceEnd = PetscSFReduceEnd_Basic;
304 sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Basic;
305 sf->ops->GetLeafRanks = PetscSFGetLeafRanks_Basic;
306 sf->ops->View = PetscSFView_Basic;
307
308 sf->ops->SetUp = PetscSFSetUp_Neighbor;
309 sf->ops->Reset = PetscSFReset_Neighbor;
310 sf->ops->Destroy = PetscSFDestroy_Neighbor;
311 sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Neighbor;
312 sf->ops->ReduceBegin = PetscSFReduceBegin_Neighbor;
313 sf->ops->FetchAndOpEnd = PetscSFFetchAndOpEnd_Neighbor;
314
315 ierr = PetscNewLog(sf,&dat)(PetscMallocA(1,PETSC_TRUE,315,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,(size_t)(1)*sizeof(**(((&dat)))),(((&dat)))) || PetscLogObjectMemory
((PetscObject)sf,sizeof(**(&dat))))
;CHKERRQ(ierr)do {if (__builtin_expect(!!(ierr),0)) return PetscError(((MPI_Comm
)0x44000001),315,__func__,"/sandbox/petsc/petsc.next-tmp/src/vec/is/sf/impls/basic/neighbor/sfneighbor.c"
,ierr,PETSC_ERROR_REPEAT," ");} while (0)
;
316 sf->data = (void*)dat;
317 PetscFunctionReturn(0)do { do { ; if (petscstack && petscstack->currentsize
> 0) { petscstack->currentsize--; petscstack->function
[petscstack->currentsize] = 0; petscstack->file[petscstack
->currentsize] = 0; petscstack->line[petscstack->currentsize
] = 0; petscstack->petscroutine[petscstack->currentsize
] = PETSC_FALSE; } if (petscstack) { petscstack->hotdepth =
(((petscstack->hotdepth-1)<(0)) ? (0) : (petscstack->
hotdepth-1)); } ; } while (0); return(0);} while (0)
;
318}
319#endif