/*@ DMLabelHasPoint - Determine whether a label assigns a value to a point Input Parameters: + label - the DMLabel - point - the point Output Parameter: . contains - Flag indicating whether the label maps this point to a value Note: The user must call DMLabelCreateIndex() before this function. Level: developer .seealso: DMLabelCreateIndex(), DMLabelGetValue(), DMLabelSetValue() @*/ PetscErrorCode DMLabelHasPoint(DMLabel label, PetscInt point, PetscBool *contains) { PetscErrorCode ierr; PetscFunctionBeginHot; PetscValidPointer(contains, 3); ierr = DMLabelMakeAllValid_Private(label);CHKERRQ(ierr); #if defined(PETSC_USE_DEBUG) if (!label->bt) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Must call DMLabelCreateIndex() before DMLabelHasPoint()"); if ((point < label->pStart) || (point >= label->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Label point %D is not in [%D, %D)", point, label->pStart, label->pEnd); #endif *contains = PetscBTLookup(label->bt, point - label->pStart) ? PETSC_TRUE : PETSC_FALSE; PetscFunctionReturn(0); }
static PetscErrorCode PCBDDCAdjGetNextLayer_Private(PetscInt* queue_tip,PetscInt n_prev,PetscBT touched,PetscInt* xadj,PetscInt* adjncy,PetscInt* n_added) { PetscInt i,j,n; PetscErrorCode ierr; PetscFunctionBegin; n = 0; for (i=-n_prev;i<0;i++) { PetscInt start_dof = queue_tip[i]; for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { PetscInt dof = adjncy[j]; if (!PetscBTLookup(touched,dof)) { ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); queue_tip[n] = dof; n++; } } } *n_added = n; PetscFunctionReturn(0); }
static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Once(Mat C,PetscInt is_max,IS is[]) { Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; PetscErrorCode ierr; PetscMPIInt size,rank,tag1,tag2,*len_s,nrqr,nrqs,*id_r1,*len_r1,flag,len; const PetscInt *idx_i; PetscInt idx,isz,col,*n,*data1,**data1_start,*data2,*data2_i,*data,*data_i, Mbs,i,j,k,*odata1,*odata2, proc_id,**odata2_ptr,*ctable=0,*btable,len_max,len_est; PetscInt proc_end=0,*iwork,len_unused,nodata2; PetscInt ois_max; /* max no of is[] in each of processor */ char *t_p; MPI_Comm comm; MPI_Request *s_waits1,*s_waits2,r_req; MPI_Status *s_status,r_status; PetscBT *table; /* mark indices of this processor's is[] */ PetscBT table_i; PetscBT otable; /* mark indices of other processors' is[] */ PetscInt bs=C->rmap->bs,Bn = c->B->cmap->n,Bnbs = Bn/bs,*Bowners; IS garray_local,garray_gl; PetscFunctionBegin; comm = ((PetscObject)C)->comm; size = c->size; rank = c->rank; Mbs = c->Mbs; ierr = PetscObjectGetNewTag((PetscObject)C,&tag1);CHKERRQ(ierr); ierr = PetscObjectGetNewTag((PetscObject)C,&tag2);CHKERRQ(ierr); /* create tables used in step 1: table[i] - mark c->garray of proc [i] step 3: table[i] - mark indices of is[i] when whose=MINE table[0] - mark incideces of is[] when whose=OTHER */ len = PetscMax(is_max, size);CHKERRQ(ierr); ierr = PetscMalloc2(len,PetscBT,&table,(Mbs/PETSC_BITS_PER_BYTE+1)*len,char,&t_p);CHKERRQ(ierr); for (i=0; i<len; i++) { table[i] = t_p + (Mbs/PETSC_BITS_PER_BYTE+1)*i; } ierr = MPI_Allreduce(&is_max,&ois_max,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); /* 1. Send this processor's is[] to other processors */ /*---------------------------------------------------*/ /* allocate spaces */ ierr = PetscMalloc(is_max*sizeof(PetscInt),&n);CHKERRQ(ierr); len = 0; for (i=0; i<is_max; i++) { ierr = ISGetLocalSize(is[i],&n[i]);CHKERRQ(ierr); len += n[i]; } if (!len) { is_max = 0; } else { len += 1 + is_max; /* max length of data1 for one processor */ } ierr = PetscMalloc((size*len+1)*sizeof(PetscInt),&data1);CHKERRQ(ierr); ierr = PetscMalloc(size*sizeof(PetscInt*),&data1_start);CHKERRQ(ierr); for (i=0; i<size; i++) data1_start[i] = data1 + i*len; ierr = PetscMalloc4(size,PetscInt,&len_s,size,PetscInt,&btable,size,PetscInt,&iwork,size+1,PetscInt,&Bowners);CHKERRQ(ierr); /* gather c->garray from all processors */ ierr = ISCreateGeneral(comm,Bnbs,c->garray,&garray_local);CHKERRQ(ierr); ierr = ISAllGather(garray_local, &garray_gl);CHKERRQ(ierr); ierr = ISDestroy(garray_local);CHKERRQ(ierr); ierr = MPI_Allgather(&Bnbs,1,MPIU_INT,Bowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); Bowners[0] = 0; for (i=0; i<size; i++) Bowners[i+1] += Bowners[i]; if (is_max){ /* hash table ctable which maps c->row to proc_id) */ ierr = PetscMalloc(Mbs*sizeof(PetscInt),&ctable);CHKERRQ(ierr); for (proc_id=0,j=0; proc_id<size; proc_id++) { for (; j<C->rmap->range[proc_id+1]/bs; j++) { ctable[j] = proc_id; } } /* hash tables marking c->garray */ ierr = ISGetIndices(garray_gl,&idx_i); for (i=0; i<size; i++){ table_i = table[i]; ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); for (j = Bowners[i]; j<Bowners[i+1]; j++){ /* go through B cols of proc[i]*/ ierr = PetscBTSet(table_i,idx_i[j]);CHKERRQ(ierr); } } ierr = ISRestoreIndices(garray_gl,&idx_i);CHKERRQ(ierr); } /* if (is_max) */ ierr = ISDestroy(garray_gl);CHKERRQ(ierr); /* evaluate communication - mesg to who, length, and buffer space */ for (i=0; i<size; i++) len_s[i] = 0; /* header of data1 */ for (proc_id=0; proc_id<size; proc_id++){ iwork[proc_id] = 0; *data1_start[proc_id] = is_max; data1_start[proc_id]++; for (j=0; j<is_max; j++) { if (proc_id == rank){ *data1_start[proc_id] = n[j]; } else { *data1_start[proc_id] = 0; } data1_start[proc_id]++; } } for (i=0; i<is_max; i++) { ierr = ISGetIndices(is[i],&idx_i);CHKERRQ(ierr); for (j=0; j<n[i]; j++){ idx = idx_i[j]; *data1_start[rank] = idx; data1_start[rank]++; /* for local proccessing */ proc_end = ctable[idx]; for (proc_id=0; proc_id<=proc_end; proc_id++){ /* for others to process */ if (proc_id == rank ) continue; /* done before this loop */ if (proc_id < proc_end && !PetscBTLookup(table[proc_id],idx)) continue; /* no need for sending idx to [proc_id] */ *data1_start[proc_id] = idx; data1_start[proc_id]++; len_s[proc_id]++; } } /* update header data */ for (proc_id=0; proc_id<size; proc_id++){ if (proc_id== rank) continue; *(data1 + proc_id*len + 1 + i) = len_s[proc_id] - iwork[proc_id]; iwork[proc_id] = len_s[proc_id] ; } ierr = ISRestoreIndices(is[i],&idx_i);CHKERRQ(ierr); } nrqs = 0; nrqr = 0; for (i=0; i<size; i++){ data1_start[i] = data1 + i*len; if (len_s[i]){ nrqs++; len_s[i] += 1 + is_max; /* add no. of header msg */ } } for (i=0; i<is_max; i++) { ierr = ISDestroy(is[i]);CHKERRQ(ierr); } ierr = PetscFree(n);CHKERRQ(ierr); ierr = PetscFree(ctable);CHKERRQ(ierr); /* Determine the number of messages to expect, their lengths, from from-ids */ ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&nrqr);CHKERRQ(ierr); ierr = PetscGatherMessageLengths(comm,nrqs,nrqr,len_s,&id_r1,&len_r1);CHKERRQ(ierr); /* Now post the sends */ ierr = PetscMalloc2(size,MPI_Request,&s_waits1,size,MPI_Request,&s_waits2);CHKERRQ(ierr); k = 0; for (proc_id=0; proc_id<size; proc_id++){ /* send data1 to processor [proc_id] */ if (len_s[proc_id]){ ierr = MPI_Isend(data1_start[proc_id],len_s[proc_id],MPIU_INT,proc_id,tag1,comm,s_waits1+k);CHKERRQ(ierr); k++; } } /* 2. Receive other's is[] and process. Then send back */ /*-----------------------------------------------------*/ len = 0; for (i=0; i<nrqr; i++){ if (len_r1[i] > len)len = len_r1[i]; } ierr = PetscFree(len_r1);CHKERRQ(ierr); ierr = PetscFree(id_r1);CHKERRQ(ierr); for (proc_id=0; proc_id<size; proc_id++) len_s[proc_id] = iwork[proc_id] = 0; ierr = PetscMalloc((len+1)*sizeof(PetscInt),&odata1);CHKERRQ(ierr); ierr = PetscMalloc(size*sizeof(PetscInt**),&odata2_ptr);CHKERRQ(ierr); ierr = PetscBTCreate(Mbs,otable);CHKERRQ(ierr); len_max = ois_max*(Mbs+1); /* max space storing all is[] for each receive */ len_est = 2*len_max; /* estimated space of storing is[] for all receiving messages */ ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr); nodata2 = 0; /* nodata2+1: num of PetscMalloc(,&odata2_ptr[]) called */ odata2_ptr[nodata2] = odata2; len_unused = len_est; /* unused space in the array odata2_ptr[nodata2]-- needs to be >= len_max */ k = 0; while (k < nrqr){ /* Receive messages */ ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag1,comm,&flag,&r_status);CHKERRQ(ierr); if (flag){ ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRQ(ierr); proc_id = r_status.MPI_SOURCE; ierr = MPI_Irecv(odata1,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr); ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr); /* Process messages */ /* make sure there is enough unused space in odata2 array */ if (len_unused < len_max){ /* allocate more space for odata2 */ ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr); odata2_ptr[++nodata2] = odata2; len_unused = len_est; } ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,odata1,OTHER,odata2,&otable);CHKERRQ(ierr); len = 1 + odata2[0]; for (i=0; i<odata2[0]; i++){ len += odata2[1 + i]; } /* Send messages back */ ierr = MPI_Isend(odata2,len,MPIU_INT,proc_id,tag2,comm,s_waits2+k);CHKERRQ(ierr); k++; odata2 += len; len_unused -= len; len_s[proc_id] = len; /* num of messages sending back to [proc_id] by this proc */ } } ierr = PetscFree(odata1);CHKERRQ(ierr); ierr = PetscBTDestroy(otable);CHKERRQ(ierr); /* 3. Do local work on this processor's is[] */ /*-------------------------------------------*/ /* make sure there is enough unused space in odata2(=data) array */ len_max = is_max*(Mbs+1); /* max space storing all is[] for this processor */ if (len_unused < len_max){ /* allocate more space for odata2 */ ierr = PetscMalloc((len_est+1)*sizeof(PetscInt),&odata2);CHKERRQ(ierr); odata2_ptr[++nodata2] = odata2; len_unused = len_est; } data = odata2; ierr = MatIncreaseOverlap_MPISBAIJ_Local(C,data1_start[rank],MINE,data,table);CHKERRQ(ierr); ierr = PetscFree(data1_start);CHKERRQ(ierr); /* 4. Receive work done on other processors, then merge */ /*------------------------------------------------------*/ /* get max number of messages that this processor expects to recv */ ierr = MPI_Allreduce(len_s,iwork,size,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); ierr = PetscMalloc((iwork[rank]+1)*sizeof(PetscInt),&data2);CHKERRQ(ierr); ierr = PetscFree4(len_s,btable,iwork,Bowners);CHKERRQ(ierr); k = 0; while (k < nrqs){ /* Receive messages */ ierr = MPI_Iprobe(MPI_ANY_SOURCE,tag2,comm,&flag,&r_status); if (flag){ ierr = MPI_Get_count(&r_status,MPIU_INT,&len);CHKERRQ(ierr); proc_id = r_status.MPI_SOURCE; ierr = MPI_Irecv(data2,len,MPIU_INT,proc_id,r_status.MPI_TAG,comm,&r_req);CHKERRQ(ierr); ierr = MPI_Wait(&r_req,&r_status);CHKERRQ(ierr); if (len > 1+is_max){ /* Add data2 into data */ data2_i = data2 + 1 + is_max; for (i=0; i<is_max; i++){ table_i = table[i]; data_i = data + 1 + is_max + Mbs*i; isz = data[1+i]; for (j=0; j<data2[1+i]; j++){ col = data2_i[j]; if (!PetscBTLookupSet(table_i,col)) {data_i[isz++] = col;} } data[1+i] = isz; if (i < is_max - 1) data2_i += data2[1+i]; } } k++; } } ierr = PetscFree(data2);CHKERRQ(ierr); ierr = PetscFree2(table,t_p);CHKERRQ(ierr); /* phase 1 sends are complete */ ierr = PetscMalloc(size*sizeof(MPI_Status),&s_status);CHKERRQ(ierr); if (nrqs) {ierr = MPI_Waitall(nrqs,s_waits1,s_status);CHKERRQ(ierr);} ierr = PetscFree(data1);CHKERRQ(ierr); /* phase 2 sends are complete */ if (nrqr){ierr = MPI_Waitall(nrqr,s_waits2,s_status);CHKERRQ(ierr);} ierr = PetscFree2(s_waits1,s_waits2);CHKERRQ(ierr); ierr = PetscFree(s_status);CHKERRQ(ierr); /* 5. Create new is[] */ /*--------------------*/ for (i=0; i<is_max; i++) { data_i = data + 1 + is_max + Mbs*i; ierr = ISCreateGeneral(PETSC_COMM_SELF,data[1+i],data_i,is+i);CHKERRQ(ierr); } for (k=0; k<=nodata2; k++){ ierr = PetscFree(odata2_ptr[k]);CHKERRQ(ierr); } ierr = PetscFree(odata2_ptr);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* Would computation be reduced by swapping the loop 'for each is' and 'for each row'? */ static PetscErrorCode MatIncreaseOverlap_MPISBAIJ_Local(Mat C,PetscInt *data,PetscInt whose,PetscInt *nidx,PetscBT *table) { Mat_MPISBAIJ *c = (Mat_MPISBAIJ*)C->data; Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)(c->A)->data; Mat_SeqBAIJ *b = (Mat_SeqBAIJ*)(c->B)->data; PetscErrorCode ierr; PetscInt row,mbs,Mbs,*nidx_i,col,col_max,isz,isz0,*ai,*aj,*bi,*bj,*garray,rstart,l; PetscInt a_start,a_end,b_start,b_end,i,j,k,is_max,*idx_i,n; PetscBT table0; /* mark the indices of input is[] for look up */ PetscBT table_i; /* poits to i-th table. When whose=OTHER, a single table is used for all is[] */ PetscFunctionBegin; Mbs = c->Mbs; mbs = a->mbs; ai = a->i; aj = a->j; bi = b->i; bj = b->j; garray = c->garray; rstart = c->rstartbs; is_max = data[0]; ierr = PetscBTCreate(Mbs,table0);CHKERRQ(ierr); nidx[0] = is_max; idx_i = data + is_max + 1; /* ptr to input is[0] array */ nidx_i = nidx + is_max + 1; /* ptr to output is[0] array */ for (i=0; i<is_max; i++) { /* for each is */ isz = 0; n = data[1+i]; /* size of input is[i] */ /* initialize and set table_i(mark idx and nidx) and table0(only mark idx) */ if (whose == MINE){ /* process this processor's is[] */ table_i = table[i]; nidx_i = nidx + 1+ is_max + Mbs*i; } else { /* process other processor's is[] - only use one temp table */ table_i = table[0]; } ierr = PetscBTMemzero(Mbs,table_i);CHKERRQ(ierr); ierr = PetscBTMemzero(Mbs,table0);CHKERRQ(ierr); if (n==0) { nidx[1+i] = 0; /* size of new is[i] */ continue; } isz0 = 0; col_max = 0; for (j=0; j<n; j++){ col = idx_i[j]; if (col >= Mbs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"index col %D >= Mbs %D",col,Mbs); if(!PetscBTLookupSet(table_i,col)) { ierr = PetscBTSet(table0,col);CHKERRQ(ierr); if (whose == MINE) {nidx_i[isz0] = col;} if (col_max < col) col_max = col; isz0++; } } if (whose == MINE) {isz = isz0;} k = 0; /* no. of indices from input is[i] that have been examined */ for (row=0; row<mbs; row++){ a_start = ai[row]; a_end = ai[row+1]; b_start = bi[row]; b_end = bi[row+1]; if (PetscBTLookup(table0,row+rstart)){ /* row is on input is[i]: do row search: collect all col in this row */ for (l = a_start; l<a_end ; l++){ /* Amat */ col = aj[l] + rstart; if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;} } for (l = b_start; l<b_end ; l++){ /* Bmat */ col = garray[bj[l]]; if (!PetscBTLookupSet(table_i,col)) {nidx_i[isz++] = col;} } k++; if (k >= isz0) break; /* for (row=0; row<mbs; row++) */ } else { /* row is not on input is[i]: do col serach: add row onto nidx_i if there is a col in nidx_i */ for (l = a_start; l<a_end ; l++){ /* Amat */ col = aj[l] + rstart; if (col > col_max) break; if (PetscBTLookup(table0,col)){ if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;} break; /* for l = start; l<end ; l++) */ } } for (l = b_start; l<b_end ; l++){ /* Bmat */ col = garray[bj[l]]; if (col > col_max) break; if (PetscBTLookup(table0,col)){ if (!PetscBTLookupSet(table_i,row+rstart)) {nidx_i[isz++] = row+rstart;} break; /* for l = start; l<end ; l++) */ } } } } if (i < is_max - 1){ idx_i += n; /* ptr to input is[i+1] array */ nidx_i += isz; /* ptr to output is[i+1] array */ } nidx[1+i] = isz; /* size of new is[i] */ } /* for each is */ ierr = PetscBTDestroy(table0);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* PCISSetUp - */ PetscErrorCode PCISSetUp(PC pc, PetscBool computesolvers) { PC_IS *pcis = (PC_IS*)(pc->data); Mat_IS *matis; MatReuse reuse; PetscErrorCode ierr; PetscBool flg,issbaij; Vec counter; PetscFunctionBegin; ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATIS,&flg);CHKERRQ(ierr); if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Preconditioner type of Neumann Neumman requires matrix of type MATIS"); matis = (Mat_IS*)pc->pmat->data; /* first time creation, get info on substructuring */ if (!pc->setupcalled) { PetscInt n_I; PetscInt *idx_I_local,*idx_B_local,*idx_I_global,*idx_B_global; PetscBT bt; PetscInt i,j; /* get info on mapping */ ierr = PetscObjectReference((PetscObject)pc->pmat->rmap->mapping);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingDestroy(&pcis->mapping);CHKERRQ(ierr); pcis->mapping = pc->pmat->rmap->mapping; ierr = ISLocalToGlobalMappingGetSize(pcis->mapping,&pcis->n);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingGetInfo(pcis->mapping,&(pcis->n_neigh),&(pcis->neigh),&(pcis->n_shared),&(pcis->shared));CHKERRQ(ierr); /* Identifying interior and interface nodes, in local numbering */ ierr = PetscBTCreate(pcis->n,&bt);CHKERRQ(ierr); for (i=0;i<pcis->n_neigh;i++) for (j=0;j<pcis->n_shared[i];j++) { ierr = PetscBTSet(bt,pcis->shared[i][j]);CHKERRQ(ierr); } /* Creating local and global index sets for interior and inteface nodes. */ ierr = PetscMalloc1(pcis->n,&idx_I_local);CHKERRQ(ierr); ierr = PetscMalloc1(pcis->n,&idx_B_local);CHKERRQ(ierr); for (i=0, pcis->n_B=0, n_I=0; i<pcis->n; i++) { if (!PetscBTLookup(bt,i)) { idx_I_local[n_I] = i; n_I++; } else { idx_B_local[pcis->n_B] = i; pcis->n_B++; } } /* Getting the global numbering */ idx_B_global = idx_I_local + n_I; /* Just avoiding allocating extra memory, since we have vacant space */ idx_I_global = idx_B_local + pcis->n_B; ierr = ISLocalToGlobalMappingApply(pcis->mapping,pcis->n_B,idx_B_local,idx_B_global);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingApply(pcis->mapping,n_I,idx_I_local,idx_I_global);CHKERRQ(ierr); /* Creating the index sets */ ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_local,PETSC_COPY_VALUES, &pcis->is_B_local);CHKERRQ(ierr); ierr = ISCreateGeneral(PETSC_COMM_SELF,pcis->n_B,idx_B_global,PETSC_COPY_VALUES,&pcis->is_B_global);CHKERRQ(ierr); ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_local,PETSC_COPY_VALUES, &pcis->is_I_local);CHKERRQ(ierr); ierr = ISCreateGeneral(PETSC_COMM_SELF,n_I,idx_I_global,PETSC_COPY_VALUES,&pcis->is_I_global);CHKERRQ(ierr); /* Freeing memory */ ierr = PetscFree(idx_B_local);CHKERRQ(ierr); ierr = PetscFree(idx_I_local);CHKERRQ(ierr); ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); /* Creating work vectors and arrays */ ierr = VecDuplicate(matis->x,&pcis->vec1_N);CHKERRQ(ierr); ierr = VecDuplicate(pcis->vec1_N,&pcis->vec2_N);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n-pcis->n_B,&pcis->vec1_D);CHKERRQ(ierr); ierr = VecDuplicate(pcis->vec1_D,&pcis->vec2_D);CHKERRQ(ierr); ierr = VecDuplicate(pcis->vec1_D,&pcis->vec3_D);CHKERRQ(ierr); ierr = VecDuplicate(pcis->vec1_D,&pcis->vec4_D);CHKERRQ(ierr); ierr = VecCreateSeq(PETSC_COMM_SELF,pcis->n_B,&pcis->vec1_B);CHKERRQ(ierr); ierr = VecDuplicate(pcis->vec1_B,&pcis->vec2_B);CHKERRQ(ierr); ierr = VecDuplicate(pcis->vec1_B,&pcis->vec3_B);CHKERRQ(ierr); ierr = MatCreateVecs(pc->pmat,&pcis->vec1_global,0);CHKERRQ(ierr); ierr = PetscMalloc1(pcis->n,&pcis->work_N);CHKERRQ(ierr); /* scaling vector */ if (!pcis->D) { /* it can happen that the user passed in a scaling vector via PCISSetSubdomainDiagonalScaling */ ierr = VecDuplicate(pcis->vec1_B,&pcis->D);CHKERRQ(ierr); ierr = VecSet(pcis->D,pcis->scaling_factor);CHKERRQ(ierr); } /* Creating the scatter contexts */ ierr = VecScatterCreate(pcis->vec1_N,pcis->is_I_local,pcis->vec1_D,(IS)0,&pcis->N_to_D);CHKERRQ(ierr); ierr = VecScatterCreate(pcis->vec1_global,pcis->is_I_global,pcis->vec1_D,(IS)0,&pcis->global_to_D);CHKERRQ(ierr); ierr = VecScatterCreate(pcis->vec1_N,pcis->is_B_local,pcis->vec1_B,(IS)0,&pcis->N_to_B);CHKERRQ(ierr); ierr = VecScatterCreate(pcis->vec1_global,pcis->is_B_global,pcis->vec1_B,(IS)0,&pcis->global_to_B);CHKERRQ(ierr); /* map from boundary to local */ ierr = ISLocalToGlobalMappingCreateIS(pcis->is_B_local,&pcis->BtoNmap);CHKERRQ(ierr); } /* Extracting the blocks A_II, A_BI, A_IB and A_BB from A. If the numbering is such that interior nodes come first than the interface ones, we have [ A_II | A_IB ] A = [------+------] [ A_BI | A_BB ] */ reuse = MAT_INITIAL_MATRIX; if (pcis->reusesubmatrices && pc->setupcalled) { if (pc->flag == SAME_NONZERO_PATTERN) { reuse = MAT_REUSE_MATRIX; } else { reuse = MAT_INITIAL_MATRIX; } } if (reuse == MAT_INITIAL_MATRIX) { ierr = MatDestroy(&pcis->A_II);CHKERRQ(ierr); ierr = MatDestroy(&pcis->A_IB);CHKERRQ(ierr); ierr = MatDestroy(&pcis->A_BI);CHKERRQ(ierr); ierr = MatDestroy(&pcis->A_BB);CHKERRQ(ierr); } ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_I_local,reuse,&pcis->A_II);CHKERRQ(ierr); ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_B_local,reuse,&pcis->A_BB);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); if (!issbaij) { ierr = MatGetSubMatrix(matis->A,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr); ierr = MatGetSubMatrix(matis->A,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr); } else { Mat newmat; ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&newmat);CHKERRQ(ierr); ierr = MatGetSubMatrix(newmat,pcis->is_I_local,pcis->is_B_local,reuse,&pcis->A_IB);CHKERRQ(ierr); ierr = MatGetSubMatrix(newmat,pcis->is_B_local,pcis->is_I_local,reuse,&pcis->A_BI);CHKERRQ(ierr); ierr = MatDestroy(&newmat);CHKERRQ(ierr); } /* Creating scaling vector D */ ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_is_use_stiffness_scaling",&pcis->use_stiffness_scaling,NULL);CHKERRQ(ierr); if (pcis->use_stiffness_scaling) { ierr = MatGetDiagonal(pcis->A_BB,pcis->D);CHKERRQ(ierr); } ierr = MatCreateVecs(pc->pmat,&counter,0);CHKERRQ(ierr); /* temporary auxiliar vector */ ierr = VecSet(counter,0.0);CHKERRQ(ierr); ierr = VecScatterBegin(pcis->global_to_B,pcis->D,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecScatterEnd(pcis->global_to_B,pcis->D,counter,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); ierr = VecScatterBegin(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(pcis->global_to_B,counter,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecPointwiseDivide(pcis->D,pcis->D,pcis->vec1_B);CHKERRQ(ierr); ierr = VecDestroy(&counter);CHKERRQ(ierr); /* See historical note 01, at the bottom of this file. */ /* Creating the KSP contexts for the local Dirichlet and Neumann problems */ if (computesolvers) { PC pc_ctx; pcis->pure_neumann = matis->pure_neumann; /* Dirichlet */ ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_D);CHKERRQ(ierr); ierr = KSPSetErrorIfNotConverged(pcis->ksp_D,pc->erroriffailure);CHKERRQ(ierr); ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_D,(PetscObject)pc,1);CHKERRQ(ierr); ierr = KSPSetOperators(pcis->ksp_D,pcis->A_II,pcis->A_II);CHKERRQ(ierr); ierr = KSPSetOptionsPrefix(pcis->ksp_D,"is_localD_");CHKERRQ(ierr); ierr = KSPGetPC(pcis->ksp_D,&pc_ctx);CHKERRQ(ierr); ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); ierr = KSPSetType(pcis->ksp_D,KSPPREONLY);CHKERRQ(ierr); ierr = KSPSetFromOptions(pcis->ksp_D);CHKERRQ(ierr); /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ ierr = KSPSetUp(pcis->ksp_D);CHKERRQ(ierr); /* Neumann */ ierr = KSPCreate(PETSC_COMM_SELF,&pcis->ksp_N);CHKERRQ(ierr); ierr = KSPSetErrorIfNotConverged(pcis->ksp_N,pc->erroriffailure);CHKERRQ(ierr); ierr = PetscObjectIncrementTabLevel((PetscObject)pcis->ksp_N,(PetscObject)pc,1);CHKERRQ(ierr); ierr = KSPSetOperators(pcis->ksp_N,matis->A,matis->A);CHKERRQ(ierr); ierr = KSPSetOptionsPrefix(pcis->ksp_N,"is_localN_");CHKERRQ(ierr); ierr = KSPGetPC(pcis->ksp_N,&pc_ctx);CHKERRQ(ierr); ierr = PCSetType(pc_ctx,PCLU);CHKERRQ(ierr); ierr = KSPSetType(pcis->ksp_N,KSPPREONLY);CHKERRQ(ierr); ierr = KSPSetFromOptions(pcis->ksp_N);CHKERRQ(ierr); { PetscBool damp_fixed = PETSC_FALSE, remove_nullspace_fixed = PETSC_FALSE, set_damping_factor_floating = PETSC_FALSE, not_damp_floating = PETSC_FALSE, not_remove_nullspace_floating = PETSC_FALSE; PetscReal fixed_factor, floating_factor; ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&fixed_factor,&damp_fixed);CHKERRQ(ierr); if (!damp_fixed) fixed_factor = 0.0; ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_damp_fixed",&damp_fixed,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_remove_nullspace_fixed",&remove_nullspace_fixed,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetReal(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating", &floating_factor,&set_damping_factor_floating);CHKERRQ(ierr); if (!set_damping_factor_floating) floating_factor = 0.0; ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_set_damping_factor_floating",&set_damping_factor_floating,NULL);CHKERRQ(ierr); if (!set_damping_factor_floating) floating_factor = 1.e-12; ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_damp_floating",¬_damp_floating,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(((PetscObject)pc_ctx)->options,((PetscObject)pc_ctx)->prefix,"-pc_is_not_remove_nullspace_floating",¬_remove_nullspace_floating,NULL);CHKERRQ(ierr); if (pcis->pure_neumann) { /* floating subdomain */ if (!(not_damp_floating)) { ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); } if (!(not_remove_nullspace_floating)) { MatNullSpace nullsp; ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr); ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); } } else { /* fixed subdomain */ if (damp_fixed) { ierr = PCFactorSetShiftType(pc_ctx,MAT_SHIFT_NONZERO);CHKERRQ(ierr); ierr = PCFactorSetShiftAmount(pc_ctx,floating_factor);CHKERRQ(ierr); } if (remove_nullspace_fixed) { MatNullSpace nullsp; ierr = MatNullSpaceCreate(PETSC_COMM_SELF,PETSC_TRUE,0,NULL,&nullsp);CHKERRQ(ierr); ierr = MatSetNullSpace(matis->A,nullsp);CHKERRQ(ierr); ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); } } } /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */ ierr = KSPSetUp(pcis->ksp_N);CHKERRQ(ierr); } PetscFunctionReturn(0); }
/*@ ISDifference - Computes the difference between two index sets. Collective on IS Input Parameter: + is1 - first index, to have items removed from it - is2 - index values to be removed Output Parameters: . isout - is1 - is2 Notes: Negative values are removed from the lists. is2 may have values that are not in is1. This requires O(imax-imin) memory and O(imax-imin) work, where imin and imax are the bounds on the indices in is1. Level: intermediate Concepts: index sets^difference Concepts: IS^difference .seealso: ISDestroy(), ISView(), ISSum(), ISExpand() @*/ PetscErrorCode ISDifference(IS is1,IS is2,IS *isout) { PetscErrorCode ierr; PetscInt i,n1,n2,imin,imax,nout,*iout; const PetscInt *i1,*i2; PetscBT mask; MPI_Comm comm; PetscFunctionBegin; PetscValidHeaderSpecific(is1,IS_CLASSID,1); PetscValidHeaderSpecific(is2,IS_CLASSID,2); PetscValidPointer(isout,3); ierr = ISGetIndices(is1,&i1);CHKERRQ(ierr); ierr = ISGetLocalSize(is1,&n1);CHKERRQ(ierr); /* Create a bit mask array to contain required values */ if (n1) { imin = PETSC_MAX_INT; imax = 0; for (i=0; i<n1; i++) { if (i1[i] < 0) continue; imin = PetscMin(imin,i1[i]); imax = PetscMax(imax,i1[i]); } } else imin = imax = 0; ierr = PetscBTCreate(imax-imin,&mask);CHKERRQ(ierr); /* Put the values from is1 */ for (i=0; i<n1; i++) { if (i1[i] < 0) continue; ierr = PetscBTSet(mask,i1[i] - imin);CHKERRQ(ierr); } ierr = ISRestoreIndices(is1,&i1);CHKERRQ(ierr); /* Remove the values from is2 */ ierr = ISGetIndices(is2,&i2);CHKERRQ(ierr); ierr = ISGetLocalSize(is2,&n2);CHKERRQ(ierr); for (i=0; i<n2; i++) { if (i2[i] < imin || i2[i] > imax) continue; ierr = PetscBTClear(mask,i2[i] - imin);CHKERRQ(ierr); } ierr = ISRestoreIndices(is2,&i2);CHKERRQ(ierr); /* Count the number in the difference */ nout = 0; for (i=0; i<imax-imin+1; i++) { if (PetscBTLookup(mask,i)) nout++; } /* create the new IS containing the difference */ ierr = PetscMalloc1(nout,&iout);CHKERRQ(ierr); nout = 0; for (i=0; i<imax-imin+1; i++) { if (PetscBTLookup(mask,i)) iout[nout++] = i + imin; } ierr = PetscObjectGetComm((PetscObject)is1,&comm);CHKERRQ(ierr); ierr = ISCreateGeneral(comm,nout,iout,PETSC_OWN_POINTER,isout);CHKERRQ(ierr); ierr = PetscBTDestroy(&mask);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* - Checks face match - Flips non-matching - Inserts faces of support cells in FIFO */ static PetscErrorCode DMPlexCheckFace_Internal(DM dm, PetscInt *faceFIFO, PetscInt *fTop, PetscInt *fBottom, PetscInt cStart, PetscInt fStart, PetscInt fEnd, PetscBT seenCells, PetscBT flippedCells, PetscBT seenFaces) { const PetscInt *support, *coneA, *coneB, *coneOA, *coneOB; PetscInt supportSize, coneSizeA, coneSizeB, posA = -1, posB = -1; PetscInt face, dim, seenA, flippedA, seenB, flippedB, mismatch, c; PetscErrorCode ierr; PetscFunctionBegin; face = faceFIFO[(*fTop)++]; ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); if (supportSize < 2) PetscFunctionReturn(0); if (supportSize != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Faces should separate only two cells, not %d", supportSize); seenA = PetscBTLookup(seenCells, support[0]-cStart); flippedA = PetscBTLookup(flippedCells, support[0]-cStart) ? 1 : 0; seenB = PetscBTLookup(seenCells, support[1]-cStart); flippedB = PetscBTLookup(flippedCells, support[1]-cStart) ? 1 : 0; ierr = DMPlexGetConeSize(dm, support[0], &coneSizeA);CHKERRQ(ierr); ierr = DMPlexGetConeSize(dm, support[1], &coneSizeB);CHKERRQ(ierr); ierr = DMPlexGetCone(dm, support[0], &coneA);CHKERRQ(ierr); ierr = DMPlexGetCone(dm, support[1], &coneB);CHKERRQ(ierr); ierr = DMPlexGetConeOrientation(dm, support[0], &coneOA);CHKERRQ(ierr); ierr = DMPlexGetConeOrientation(dm, support[1], &coneOB);CHKERRQ(ierr); for (c = 0; c < coneSizeA; ++c) { if (!PetscBTLookup(seenFaces, coneA[c]-fStart)) { faceFIFO[(*fBottom)++] = coneA[c]; ierr = PetscBTSet(seenFaces, coneA[c]-fStart);CHKERRQ(ierr); } if (coneA[c] == face) posA = c; if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart); } if (posA < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[0]); for (c = 0; c < coneSizeB; ++c) { if (!PetscBTLookup(seenFaces, coneB[c]-fStart)) { faceFIFO[(*fBottom)++] = coneB[c]; ierr = PetscBTSet(seenFaces, coneB[c]-fStart);CHKERRQ(ierr); } if (coneB[c] == face) posB = c; if (*fBottom > fEnd-fStart) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face %d was pushed exceeding capacity %d > %d", coneA[c], *fBottom, fEnd-fStart); } if (posB < 0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %d could not be located in cell %d", face, support[1]); if (dim == 1) { mismatch = posA == posB; } else { mismatch = coneOA[posA] == coneOB[posB]; } if (mismatch ^ (flippedA ^ flippedB)) { if (seenA && seenB) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen cells %d and %d do not match: Fault mesh is non-orientable", support[0], support[1]); if (!seenA && !flippedA) { ierr = PetscBTSet(flippedCells, support[0]-cStart);CHKERRQ(ierr); } else if (!seenB && !flippedB) { ierr = PetscBTSet(flippedCells, support[1]-cStart);CHKERRQ(ierr); } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); ierr = PetscBTSet(seenCells, support[0]-cStart);CHKERRQ(ierr); ierr = PetscBTSet(seenCells, support[1]-cStart);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@ DMPlexOrient - Give a consistent orientation to the input mesh Input Parameters: . dm - The DM Note: The orientation data for the DM are change in-place. $ This routine will fail for non-orientable surfaces, such as the Moebius strip. Level: advanced .seealso: DMCreate(), DMPLEX @*/ PetscErrorCode DMPlexOrient(DM dm) { MPI_Comm comm; PetscSF sf; const PetscInt *lpoints; const PetscSFNode *rpoints; PetscSFNode *rorntComp = NULL, *lorntComp = NULL; PetscInt *numNeighbors, **neighbors; PetscSFNode *nrankComp; PetscBool *match, *flipped; PetscBT seenCells, flippedCells, seenFaces; PetscInt *faceFIFO, fTop, fBottom, *cellComp, *faceComp; PetscInt numLeaves, numRoots, dim, h, cStart, cEnd, c, cell, fStart, fEnd, face, off, totNeighbors = 0; PetscMPIInt rank, size, numComponents, comp = 0; PetscBool flg, flg2; PetscViewer viewer = NULL, selfviewer = NULL; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view", &flg);CHKERRQ(ierr); ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-orientation_view_synchronized", &flg2);CHKERRQ(ierr); ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); ierr = PetscSFGetGraph(sf, &numRoots, &numLeaves, &lpoints, &rpoints);CHKERRQ(ierr); /* Truth Table mismatch flips do action mismatch flipA ^ flipB action F 0 flips no F F F F 1 flip yes F T T F 2 flips no T F T T 0 flips yes T T F T 1 flip no T 2 flips yes */ ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); ierr = DMPlexGetVTKCellHeight(dm, &h);CHKERRQ(ierr); ierr = DMPlexGetHeightStratum(dm, h, &cStart, &cEnd);CHKERRQ(ierr); ierr = DMPlexGetHeightStratum(dm, h+1, &fStart, &fEnd);CHKERRQ(ierr); ierr = PetscBTCreate(cEnd - cStart, &seenCells);CHKERRQ(ierr); ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); ierr = PetscBTCreate(cEnd - cStart, &flippedCells);CHKERRQ(ierr); ierr = PetscBTMemzero(cEnd - cStart, flippedCells);CHKERRQ(ierr); ierr = PetscBTCreate(fEnd - fStart, &seenFaces);CHKERRQ(ierr); ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); ierr = PetscCalloc3(fEnd - fStart, &faceFIFO, cEnd-cStart, &cellComp, fEnd-fStart, &faceComp);CHKERRQ(ierr); /* OLD STYLE - Add an integer array over cells and faces (component) for connected component number Foreach component - Mark the initial cell as seen - Process component as usual - Set component for all seenCells - Wipe seenCells and seenFaces (flippedCells can stay) - Generate parallel adjacency for component using SF and seenFaces - Collect numComponents adj data from each proc to 0 - Build same serial graph - Use same solver - Use Scatterv to to send back flipped flags for each component - Negate flippedCells by component NEW STYLE - Create the adj on each process - Bootstrap to complete graph on proc 0 */ /* Loop over components */ for (cell = cStart; cell < cEnd; ++cell) cellComp[cell-cStart] = -1; do { /* Look for first unmarked cell */ for (cell = cStart; cell < cEnd; ++cell) if (cellComp[cell-cStart] < 0) break; if (cell >= cEnd) break; /* Initialize FIFO with first cell in component */ { const PetscInt *cone; PetscInt coneSize; fTop = fBottom = 0; ierr = DMPlexGetConeSize(dm, cell, &coneSize);CHKERRQ(ierr); ierr = DMPlexGetCone(dm, cell, &cone);CHKERRQ(ierr); for (c = 0; c < coneSize; ++c) { faceFIFO[fBottom++] = cone[c]; ierr = PetscBTSet(seenFaces, cone[c]-fStart);CHKERRQ(ierr); } ierr = PetscBTSet(seenCells, cell-cStart);CHKERRQ(ierr); } /* Consider each face in FIFO */ while (fTop < fBottom) { ierr = DMPlexCheckFace_Internal(dm, faceFIFO, &fTop, &fBottom, cStart, fStart, fEnd, seenCells, flippedCells, seenFaces);CHKERRQ(ierr); } /* Set component for cells and faces */ for (cell = 0; cell < cEnd-cStart; ++cell) { if (PetscBTLookup(seenCells, cell)) cellComp[cell] = comp; } for (face = 0; face < fEnd-fStart; ++face) { if (PetscBTLookup(seenFaces, face)) faceComp[face] = comp; } /* Wipe seenCells and seenFaces for next component */ ierr = PetscBTMemzero(fEnd - fStart, seenFaces);CHKERRQ(ierr); ierr = PetscBTMemzero(cEnd - cStart, seenCells);CHKERRQ(ierr); ++comp; } while (1); numComponents = comp; if (flg) { PetscViewer v; ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr); ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for serial flipped cells:\n", rank);CHKERRQ(ierr); ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); ierr = PetscViewerFlush(v);CHKERRQ(ierr); ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr); } /* Now all subdomains are oriented, but we need a consistent parallel orientation */ if (numLeaves >= 0) { /* Store orientations of boundary faces*/ ierr = PetscCalloc2(numRoots,&rorntComp,numRoots,&lorntComp);CHKERRQ(ierr); for (face = fStart; face < fEnd; ++face) { const PetscInt *cone, *support, *ornt; PetscInt coneSize, supportSize; ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); if (supportSize != 1) continue; ierr = DMPlexGetSupport(dm, face, &support);CHKERRQ(ierr); ierr = DMPlexGetCone(dm, support[0], &cone);CHKERRQ(ierr); ierr = DMPlexGetConeSize(dm, support[0], &coneSize);CHKERRQ(ierr); ierr = DMPlexGetConeOrientation(dm, support[0], &ornt);CHKERRQ(ierr); for (c = 0; c < coneSize; ++c) if (cone[c] == face) break; if (dim == 1) { /* Use cone position instead, shifted to -1 or 1 */ if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = 1-c*2; else rorntComp[face].rank = c*2-1; } else { if (PetscBTLookup(flippedCells, support[0]-cStart)) rorntComp[face].rank = ornt[c] < 0 ? -1 : 1; else rorntComp[face].rank = ornt[c] < 0 ? 1 : -1; } rorntComp[face].index = faceComp[face-fStart]; } /* Communicate boundary edge orientations */ ierr = PetscSFBcastBegin(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr); ierr = PetscSFBcastEnd(sf, MPIU_2INT, rorntComp, lorntComp);CHKERRQ(ierr); } /* Get process adjacency */ ierr = PetscMalloc2(numComponents, &numNeighbors, numComponents, &neighbors);CHKERRQ(ierr); viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)dm)); if (flg2) {ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);} ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); for (comp = 0; comp < numComponents; ++comp) { PetscInt l, n; numNeighbors[comp] = 0; ierr = PetscMalloc1(PetscMax(numLeaves, 0), &neighbors[comp]);CHKERRQ(ierr); /* I know this is p^2 time in general, but for bounded degree its alright */ for (l = 0; l < numLeaves; ++l) { const PetscInt face = lpoints[l]; /* Find a representative face (edge) separating pairs of procs */ if ((face >= fStart) && (face < fEnd) && (faceComp[face-fStart] == comp)) { const PetscInt rrank = rpoints[l].rank; const PetscInt rcomp = lorntComp[face].index; for (n = 0; n < numNeighbors[comp]; ++n) if ((rrank == rpoints[neighbors[comp][n]].rank) && (rcomp == lorntComp[lpoints[neighbors[comp][n]]].index)) break; if (n >= numNeighbors[comp]) { PetscInt supportSize; ierr = DMPlexGetSupportSize(dm, face, &supportSize);CHKERRQ(ierr); if (supportSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary faces should see one cell, not %d", supportSize); if (flg) {ierr = PetscViewerASCIIPrintf(selfviewer, "[%d]: component %d, Found representative leaf %d (face %d) connecting to face %d on (%d, %d) with orientation %d\n", rank, comp, l, face, rpoints[l].index, rrank, rcomp, lorntComp[face].rank);CHKERRQ(ierr);} neighbors[comp][numNeighbors[comp]++] = l; } } } totNeighbors += numNeighbors[comp]; } ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&selfviewer);CHKERRQ(ierr); ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); if (flg2) {ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);} ierr = PetscMalloc2(totNeighbors, &nrankComp, totNeighbors, &match);CHKERRQ(ierr); for (comp = 0, off = 0; comp < numComponents; ++comp) { PetscInt n; for (n = 0; n < numNeighbors[comp]; ++n, ++off) { const PetscInt face = lpoints[neighbors[comp][n]]; const PetscInt o = rorntComp[face].rank*lorntComp[face].rank; if (o < 0) match[off] = PETSC_TRUE; else if (o > 0) match[off] = PETSC_FALSE; else SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid face %d (%d, %d) neighbor: %d comp: %d", face, rorntComp[face], lorntComp[face], neighbors[comp][n], comp); nrankComp[off].rank = rpoints[neighbors[comp][n]].rank; nrankComp[off].index = lorntComp[lpoints[neighbors[comp][n]]].index; } ierr = PetscFree(neighbors[comp]);CHKERRQ(ierr); } /* Collect the graph on 0 */ if (numLeaves >= 0) { Mat G; PetscBT seenProcs, flippedProcs; PetscInt *procFIFO, pTop, pBottom; PetscInt *N = NULL, *Noff; PetscSFNode *adj = NULL; PetscBool *val = NULL; PetscMPIInt *recvcounts = NULL, *displs = NULL, *Nc, p, o; PetscMPIInt size = 0; ierr = PetscCalloc1(numComponents, &flipped);CHKERRQ(ierr); if (!rank) {ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);} ierr = PetscCalloc4(size, &recvcounts, size+1, &displs, size, &Nc, size+1, &Noff);CHKERRQ(ierr); ierr = MPI_Gather(&numComponents, 1, MPI_INT, Nc, 1, MPI_INT, 0, comm);CHKERRQ(ierr); for (p = 0; p < size; ++p) { displs[p+1] = displs[p] + Nc[p]; } if (!rank) {ierr = PetscMalloc1(displs[size],&N);CHKERRQ(ierr);} ierr = MPI_Gatherv(numNeighbors, numComponents, MPIU_INT, N, Nc, displs, MPIU_INT, 0, comm);CHKERRQ(ierr); for (p = 0, o = 0; p < size; ++p) { recvcounts[p] = 0; for (c = 0; c < Nc[p]; ++c, ++o) recvcounts[p] += N[o]; displs[p+1] = displs[p] + recvcounts[p]; } if (!rank) {ierr = PetscMalloc2(displs[size], &adj, displs[size], &val);CHKERRQ(ierr);} ierr = MPI_Gatherv(nrankComp, totNeighbors, MPIU_2INT, adj, recvcounts, displs, MPIU_2INT, 0, comm);CHKERRQ(ierr); ierr = MPI_Gatherv(match, totNeighbors, MPIU_BOOL, val, recvcounts, displs, MPIU_BOOL, 0, comm);CHKERRQ(ierr); ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr); if (!rank) { for (p = 1; p <= size; ++p) {Noff[p] = Noff[p-1] + Nc[p-1];} if (flg) { PetscInt n; for (p = 0, off = 0; p < size; ++p) { for (c = 0; c < Nc[p]; ++c) { ierr = PetscPrintf(PETSC_COMM_SELF, "Proc %d Comp %d:\n", p, c);CHKERRQ(ierr); for (n = 0; n < N[Noff[p]+c]; ++n, ++off) { ierr = PetscPrintf(PETSC_COMM_SELF, " edge (%d, %d) (%d):\n", adj[off].rank, adj[off].index, val[off]);CHKERRQ(ierr); } } } } /* Symmetrize the graph */ ierr = MatCreate(PETSC_COMM_SELF, &G);CHKERRQ(ierr); ierr = MatSetSizes(G, Noff[size], Noff[size], Noff[size], Noff[size]);CHKERRQ(ierr); ierr = MatSetUp(G);CHKERRQ(ierr); for (p = 0, off = 0; p < size; ++p) { for (c = 0; c < Nc[p]; ++c) { const PetscInt r = Noff[p]+c; PetscInt n; for (n = 0; n < N[r]; ++n, ++off) { const PetscInt q = Noff[adj[off].rank] + adj[off].index; const PetscScalar o = val[off] ? 1.0 : 0.0; ierr = MatSetValues(G, 1, &r, 1, &q, &o, INSERT_VALUES);CHKERRQ(ierr); ierr = MatSetValues(G, 1, &q, 1, &r, &o, INSERT_VALUES);CHKERRQ(ierr); } } } ierr = MatAssemblyBegin(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(G, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = PetscBTCreate(Noff[size], &seenProcs);CHKERRQ(ierr); ierr = PetscBTMemzero(Noff[size], seenProcs);CHKERRQ(ierr); ierr = PetscBTCreate(Noff[size], &flippedProcs);CHKERRQ(ierr); ierr = PetscBTMemzero(Noff[size], flippedProcs);CHKERRQ(ierr); ierr = PetscMalloc1(Noff[size], &procFIFO);CHKERRQ(ierr); pTop = pBottom = 0; for (p = 0; p < Noff[size]; ++p) { if (PetscBTLookup(seenProcs, p)) continue; /* Initialize FIFO with next proc */ procFIFO[pBottom++] = p; ierr = PetscBTSet(seenProcs, p);CHKERRQ(ierr); /* Consider each proc in FIFO */ while (pTop < pBottom) { const PetscScalar *ornt; const PetscInt *neighbors; PetscInt proc, nproc, seen, flippedA, flippedB, mismatch, numNeighbors, n; proc = procFIFO[pTop++]; flippedA = PetscBTLookup(flippedProcs, proc) ? 1 : 0; ierr = MatGetRow(G, proc, &numNeighbors, &neighbors, &ornt);CHKERRQ(ierr); /* Loop over neighboring procs */ for (n = 0; n < numNeighbors; ++n) { nproc = neighbors[n]; mismatch = PetscRealPart(ornt[n]) > 0.5 ? 0 : 1; seen = PetscBTLookup(seenProcs, nproc); flippedB = PetscBTLookup(flippedProcs, nproc) ? 1 : 0; if (mismatch ^ (flippedA ^ flippedB)) { if (seen) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Previously seen procs %d and %d do not match: Fault mesh is non-orientable", proc, nproc); if (!flippedB) { ierr = PetscBTSet(flippedProcs, nproc);CHKERRQ(ierr); } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent mesh orientation: Fault mesh is non-orientable"); } else if (mismatch && flippedA && flippedB) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Attempt to flip already flipped cell: Fault mesh is non-orientable"); if (!seen) { procFIFO[pBottom++] = nproc; ierr = PetscBTSet(seenProcs, nproc);CHKERRQ(ierr); } } } } ierr = PetscFree(procFIFO);CHKERRQ(ierr); ierr = MatDestroy(&G);CHKERRQ(ierr); ierr = PetscFree2(adj, val);CHKERRQ(ierr); ierr = PetscBTDestroy(&seenProcs);CHKERRQ(ierr); } /* Scatter flip flags */ { PetscBool *flips = NULL; if (!rank) { ierr = PetscMalloc1(Noff[size], &flips);CHKERRQ(ierr); for (p = 0; p < Noff[size]; ++p) { flips[p] = PetscBTLookup(flippedProcs, p) ? PETSC_TRUE : PETSC_FALSE; if (flg && flips[p]) {ierr = PetscPrintf(comm, "Flipping Proc+Comp %d:\n", p);CHKERRQ(ierr);} } for (p = 0; p < size; ++p) { displs[p+1] = displs[p] + Nc[p]; } } ierr = MPI_Scatterv(flips, Nc, displs, MPIU_BOOL, flipped, numComponents, MPIU_BOOL, 0, comm);CHKERRQ(ierr); ierr = PetscFree(flips);CHKERRQ(ierr); } if (!rank) {ierr = PetscBTDestroy(&flippedProcs);CHKERRQ(ierr);} ierr = PetscFree(N);CHKERRQ(ierr); ierr = PetscFree4(recvcounts, displs, Nc, Noff);CHKERRQ(ierr); ierr = PetscFree2(nrankComp, match);CHKERRQ(ierr); /* Decide whether to flip cells in each component */ for (c = 0; c < cEnd-cStart; ++c) {if (flipped[cellComp[c]]) {ierr = PetscBTNegate(flippedCells, c);CHKERRQ(ierr);}} ierr = PetscFree(flipped);CHKERRQ(ierr); } if (flg) { PetscViewer v; ierr = PetscViewerASCIIGetStdout(comm, &v);CHKERRQ(ierr); ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr); ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d]BT for parallel flipped cells:\n", rank);CHKERRQ(ierr); ierr = PetscBTView(cEnd-cStart, flippedCells, v);CHKERRQ(ierr); ierr = PetscViewerFlush(v);CHKERRQ(ierr); ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr); } /* Reverse flipped cells in the mesh */ for (c = cStart; c < cEnd; ++c) { if (PetscBTLookup(flippedCells, c-cStart)) { ierr = DMPlexReverseCell(dm, c);CHKERRQ(ierr); } } ierr = PetscBTDestroy(&seenCells);CHKERRQ(ierr); ierr = PetscBTDestroy(&flippedCells);CHKERRQ(ierr); ierr = PetscBTDestroy(&seenFaces);CHKERRQ(ierr); ierr = PetscFree2(numNeighbors, neighbors);CHKERRQ(ierr); ierr = PetscFree2(rorntComp, lorntComp);CHKERRQ(ierr); ierr = PetscFree3(faceFIFO, cellComp, faceComp);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode PCBDDCGraphSetUp(PCBDDCGraph graph, PetscInt custom_minimal_size, IS neumann_is, IS dirichlet_is, PetscInt n_ISForDofs, IS ISForDofs[], IS custom_primal_vertices) { IS subset,subset_n; MPI_Comm comm; const PetscInt *is_indices; PetscInt n_neigh,*neigh,*n_shared,**shared,*queue_global; PetscInt i,j,k,s,total_counts,nodes_touched,is_size; PetscMPIInt commsize; PetscBool same_set,mirrors_found; PetscErrorCode ierr; PetscFunctionBegin; PetscValidLogicalCollectiveInt(graph->l2gmap,custom_minimal_size,2); if (neumann_is) { PetscValidHeaderSpecific(neumann_is,IS_CLASSID,3); PetscCheckSameComm(graph->l2gmap,1,neumann_is,3); } graph->has_dirichlet = PETSC_FALSE; if (dirichlet_is) { PetscValidHeaderSpecific(dirichlet_is,IS_CLASSID,4); PetscCheckSameComm(graph->l2gmap,1,dirichlet_is,4); graph->has_dirichlet = PETSC_TRUE; } PetscValidLogicalCollectiveInt(graph->l2gmap,n_ISForDofs,5); for (i=0;i<n_ISForDofs;i++) { PetscValidHeaderSpecific(ISForDofs[i],IS_CLASSID,6); PetscCheckSameComm(graph->l2gmap,1,ISForDofs[i],6); } if (custom_primal_vertices) { PetscValidHeaderSpecific(custom_primal_vertices,IS_CLASSID,6); PetscCheckSameComm(graph->l2gmap,1,custom_primal_vertices,7); } ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&comm);CHKERRQ(ierr); ierr = MPI_Comm_size(comm,&commsize);CHKERRQ(ierr); /* custom_minimal_size */ graph->custom_minimal_size = custom_minimal_size; /* get info l2gmap and allocate work vectors */ ierr = ISLocalToGlobalMappingGetInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); /* check if we have any local periodic nodes (periodic BCs) */ mirrors_found = PETSC_FALSE; if (graph->nvtxs && n_neigh) { for (i=0; i<n_shared[0]; i++) graph->count[shared[0][i]] += 1; for (i=0; i<n_shared[0]; i++) { if (graph->count[shared[0][i]] > 1) { mirrors_found = PETSC_TRUE; break; } } } /* compute local mirrors (if any) */ if (mirrors_found) { IS to,from; PetscInt *local_indices,*global_indices; ierr = ISCreateStride(PETSC_COMM_SELF,graph->nvtxs,0,1,&to);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,to,&from);CHKERRQ(ierr); /* get arrays of local and global indices */ ierr = PetscMalloc1(graph->nvtxs,&local_indices);CHKERRQ(ierr); ierr = ISGetIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr); ierr = PetscMemcpy(local_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); ierr = ISRestoreIndices(to,(const PetscInt**)&is_indices);CHKERRQ(ierr); ierr = PetscMalloc1(graph->nvtxs,&global_indices);CHKERRQ(ierr); ierr = ISGetIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr); ierr = PetscMemcpy(global_indices,is_indices,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); ierr = ISRestoreIndices(from,(const PetscInt**)&is_indices);CHKERRQ(ierr); /* allocate space for mirrors */ ierr = PetscMalloc2(graph->nvtxs,&graph->mirrors,graph->nvtxs,&graph->mirrors_set);CHKERRQ(ierr); ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); graph->mirrors_set[0] = 0; k=0; for (i=0;i<n_shared[0];i++) { j=shared[0][i]; if (graph->count[j] > 1) { graph->mirrors[j]++; k++; } } /* allocate space for set of mirrors */ ierr = PetscMalloc1(k,&graph->mirrors_set[0]);CHKERRQ(ierr); for (i=1;i<graph->nvtxs;i++) graph->mirrors_set[i]=graph->mirrors_set[i-1]+graph->mirrors[i-1]; /* fill arrays */ ierr = PetscMemzero(graph->mirrors,graph->nvtxs*sizeof(PetscInt));CHKERRQ(ierr); for (j=0;j<n_shared[0];j++) { i=shared[0][j]; if (graph->count[i] > 1) graph->mirrors_set[i][graph->mirrors[i]++]=global_indices[i]; } ierr = PetscSortIntWithArray(graph->nvtxs,global_indices,local_indices);CHKERRQ(ierr); for (i=0;i<graph->nvtxs;i++) { if (graph->mirrors[i] > 0) { ierr = PetscFindInt(graph->mirrors_set[i][0],graph->nvtxs,global_indices,&k);CHKERRQ(ierr); j = global_indices[k]; while ( k > 0 && global_indices[k-1] == j) k--; for (j=0;j<graph->mirrors[i];j++) { graph->mirrors_set[i][j]=local_indices[k+j]; } ierr = PetscSortInt(graph->mirrors[i],graph->mirrors_set[i]);CHKERRQ(ierr); } } ierr = PetscFree(local_indices);CHKERRQ(ierr); ierr = PetscFree(global_indices);CHKERRQ(ierr); ierr = ISDestroy(&to);CHKERRQ(ierr); ierr = ISDestroy(&from);CHKERRQ(ierr); } ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr); /* Count total number of neigh per node */ k = 0; for (i=1;i<n_neigh;i++) { k += n_shared[i]; for (j=0;j<n_shared[i];j++) { graph->count[shared[i][j]] += 1; } } /* Allocate space for storing the set of neighbours for each node */ if (graph->nvtxs) { ierr = PetscMalloc1(k,&graph->neighbours_set[0]);CHKERRQ(ierr); } for (i=1;i<graph->nvtxs;i++) { /* dont count myself */ graph->neighbours_set[i]=graph->neighbours_set[i-1]+graph->count[i-1]; } /* Get information for sharing subdomains */ ierr = PetscMemzero(graph->count,graph->nvtxs*sizeof(*graph->count));CHKERRQ(ierr); for (i=1;i<n_neigh;i++) { /* dont count myself */ s = n_shared[i]; for (j=0;j<s;j++) { k = shared[i][j]; graph->neighbours_set[k][graph->count[k]] = neigh[i]; graph->count[k] += 1; } } /* sort set of sharing subdomains */ for (i=0;i<graph->nvtxs;i++) { ierr = PetscSortRemoveDupsInt(&graph->count[i],graph->neighbours_set[i]);CHKERRQ(ierr); } /* free memory allocated by ISLocalToGlobalMappingGetInfo */ ierr = ISLocalToGlobalMappingRestoreInfo(graph->l2gmap,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); /* Get info for dofs splitting User can specify just a subset; an additional field is considered as a complementary field */ for (i=0;i<graph->nvtxs;i++) graph->which_dof[i] = n_ISForDofs; /* by default a dof belongs to the complement set */ for (i=0;i<n_ISForDofs;i++) { ierr = ISGetLocalSize(ISForDofs[i],&is_size);CHKERRQ(ierr); ierr = ISGetIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr); for (j=0;j<is_size;j++) { if (is_indices[j] > -1 && is_indices[j] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ graph->which_dof[is_indices[j]] = i; } } ierr = ISRestoreIndices(ISForDofs[i],(const PetscInt**)&is_indices);CHKERRQ(ierr); } /* Take into account Neumann nodes */ if (neumann_is) { ierr = ISGetLocalSize(neumann_is,&is_size);CHKERRQ(ierr); ierr = ISGetIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); for (i=0;i<is_size;i++) { if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ graph->special_dof[is_indices[i]] = PCBDDCGRAPH_NEUMANN_MARK; } } ierr = ISRestoreIndices(neumann_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); } /* Take into account Dirichlet nodes (they overwrite any neumann boundary mark previously set) */ if (dirichlet_is) { ierr = ISGetLocalSize(dirichlet_is,&is_size);CHKERRQ(ierr); ierr = ISGetIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); for (i=0;i<is_size;i++){ if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs) { /* out of bounds indices (if any) are skipped */ if (commsize > graph->commsizelimit) { /* dirichlet nodes treated as internal */ ierr = PetscBTSet(graph->touched,is_indices[i]);CHKERRQ(ierr); graph->subset[is_indices[i]] = 0; } graph->special_dof[is_indices[i]] = PCBDDCGRAPH_DIRICHLET_MARK; } } ierr = ISRestoreIndices(dirichlet_is,(const PetscInt**)&is_indices);CHKERRQ(ierr); } /* mark local periodic nodes (if any) and adapt CSR graph (if any) */ if (graph->mirrors) { for (i=0;i<graph->nvtxs;i++) if (graph->mirrors[i]) graph->special_dof[i] = PCBDDCGRAPH_LOCAL_PERIODIC_MARK; if (graph->xadj) { PetscInt *new_xadj,*new_adjncy; /* sort CSR graph */ for (i=0;i<graph->nvtxs;i++) ierr = PetscSortInt(graph->xadj[i+1]-graph->xadj[i],&graph->adjncy[graph->xadj[i]]);CHKERRQ(ierr); /* adapt local CSR graph in case of local periodicity */ k = 0; for (i=0;i<graph->nvtxs;i++) for (j=graph->xadj[i];j<graph->xadj[i+1];j++) k += graph->mirrors[graph->adjncy[j]]; ierr = PetscMalloc1(graph->nvtxs+1,&new_xadj);CHKERRQ(ierr); ierr = PetscMalloc1(k+graph->xadj[graph->nvtxs],&new_adjncy);CHKERRQ(ierr); new_xadj[0] = 0; for (i=0;i<graph->nvtxs;i++) { k = graph->xadj[i+1]-graph->xadj[i]; ierr = PetscMemcpy(&new_adjncy[new_xadj[i]],&graph->adjncy[graph->xadj[i]],k*sizeof(PetscInt));CHKERRQ(ierr); new_xadj[i+1] = new_xadj[i]+k; for (j=graph->xadj[i];j<graph->xadj[i+1];j++) { k = graph->mirrors[graph->adjncy[j]]; ierr = PetscMemcpy(&new_adjncy[new_xadj[i+1]],graph->mirrors_set[graph->adjncy[j]],k*sizeof(PetscInt));CHKERRQ(ierr); new_xadj[i+1] += k; } k = new_xadj[i+1]-new_xadj[i]; ierr = PetscSortRemoveDupsInt(&k,&new_adjncy[new_xadj[i]]);CHKERRQ(ierr); new_xadj[i+1] = new_xadj[i]+k; } /* set new CSR into graph */ ierr = PetscFree(graph->xadj);CHKERRQ(ierr); ierr = PetscFree(graph->adjncy);CHKERRQ(ierr); graph->xadj = new_xadj; graph->adjncy = new_adjncy; } } /* mark special nodes (if any) -> each will become a single node equivalence class */ if (custom_primal_vertices) { ierr = ISGetLocalSize(custom_primal_vertices,&is_size);CHKERRQ(ierr); ierr = ISGetIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); for (i=0,j=0;i<is_size;i++){ if (is_indices[i] > -1 && is_indices[i] < graph->nvtxs && graph->special_dof[is_indices[i]] != PCBDDCGRAPH_DIRICHLET_MARK) { /* out of bounds indices (if any) are skipped */ graph->special_dof[is_indices[i]] = PCBDDCGRAPH_SPECIAL_MARK-j; j++; } } ierr = ISRestoreIndices(custom_primal_vertices,(const PetscInt**)&is_indices);CHKERRQ(ierr); } /* mark interior nodes (if commsize > graph->commsizelimit) as touched and belonging to partition number 0 */ if (commsize > graph->commsizelimit) { for (i=0;i<graph->nvtxs;i++) { if (!graph->count[i]) { ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); graph->subset[i] = 0; } } } /* init graph structure and compute default subsets */ nodes_touched = 0; for (i=0;i<graph->nvtxs;i++) { if (PetscBTLookup(graph->touched,i)) { nodes_touched++; } } i = 0; graph->ncc = 0; total_counts = 0; /* allocated space for queues */ if (commsize == graph->commsizelimit) { ierr = PetscMalloc2(graph->nvtxs+1,&graph->cptr,graph->nvtxs,&graph->queue);CHKERRQ(ierr); } else { PetscInt nused = graph->nvtxs - nodes_touched; ierr = PetscMalloc2(nused+1,&graph->cptr,nused,&graph->queue);CHKERRQ(ierr); } while (nodes_touched<graph->nvtxs) { /* find first untouched node in local ordering */ while (PetscBTLookup(graph->touched,i)) i++; ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); graph->subset[i] = graph->ncc+1; graph->cptr[graph->ncc] = total_counts; graph->queue[total_counts] = i; total_counts++; nodes_touched++; /* now find all other nodes having the same set of sharing subdomains */ for (j=i+1;j<graph->nvtxs;j++) { /* check for same number of sharing subdomains, dof number and same special mark */ if (!PetscBTLookup(graph->touched,j) && graph->count[i] == graph->count[j] && graph->which_dof[i] == graph->which_dof[j] && graph->special_dof[i] == graph->special_dof[j]) { /* check for same set of sharing subdomains */ same_set = PETSC_TRUE; for (k=0;k<graph->count[j];k++){ if (graph->neighbours_set[i][k] != graph->neighbours_set[j][k]) { same_set = PETSC_FALSE; } } /* I found a friend of mine */ if (same_set) { ierr = PetscBTSet(graph->touched,j);CHKERRQ(ierr); graph->subset[j] = graph->ncc+1; nodes_touched++; graph->queue[total_counts] = j; total_counts++; } } } graph->ncc++; } /* set default number of subsets (at this point no info on csr and/or local_subs has been taken into account, so n_subsets = ncc */ graph->n_subsets = graph->ncc; ierr = PetscMalloc1(graph->n_subsets,&graph->subset_ncc);CHKERRQ(ierr); for (i=0;i<graph->n_subsets;i++) { graph->subset_ncc[i] = 1; } /* final pointer */ graph->cptr[graph->ncc] = total_counts; /* For consistency reasons (among neighbours), I need to sort (by global ordering) each connected component */ /* Get a reference node (min index in global ordering) for each subset for tagging messages */ ierr = PetscMalloc1(graph->ncc,&graph->subset_ref_node);CHKERRQ(ierr); ierr = PetscMalloc1(graph->cptr[graph->ncc],&queue_global);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingApply(graph->l2gmap,graph->cptr[graph->ncc],graph->queue,queue_global);CHKERRQ(ierr); for (j=0;j<graph->ncc;j++) { ierr = PetscSortIntWithArray(graph->cptr[j+1]-graph->cptr[j],&queue_global[graph->cptr[j]],&graph->queue[graph->cptr[j]]);CHKERRQ(ierr); graph->subset_ref_node[j] = graph->queue[graph->cptr[j]]; } ierr = PetscFree(queue_global);CHKERRQ(ierr); graph->queue_sorted = PETSC_TRUE; /* save information on subsets (needed when analyzing the connected components) */ if (graph->ncc) { ierr = PetscMalloc2(graph->ncc,&graph->subset_size,graph->ncc,&graph->subset_idxs);CHKERRQ(ierr); ierr = PetscMalloc1(graph->cptr[graph->ncc],&graph->subset_idxs[0]);CHKERRQ(ierr); ierr = PetscMemzero(graph->subset_idxs[0],graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr); for (j=1;j<graph->ncc;j++) { graph->subset_size[j-1] = graph->cptr[j] - graph->cptr[j-1]; graph->subset_idxs[j] = graph->subset_idxs[j-1] + graph->subset_size[j-1]; } graph->subset_size[graph->ncc-1] = graph->cptr[graph->ncc] - graph->cptr[graph->ncc-1]; ierr = PetscMemcpy(graph->subset_idxs[0],graph->queue,graph->cptr[graph->ncc]*sizeof(PetscInt));CHKERRQ(ierr); } /* renumber reference nodes */ ierr = ISCreateGeneral(PetscObjectComm((PetscObject)(graph->l2gmap)),graph->ncc,graph->subset_ref_node,PETSC_COPY_VALUES,&subset_n);CHKERRQ(ierr); ierr = ISLocalToGlobalMappingApplyIS(graph->l2gmap,subset_n,&subset);CHKERRQ(ierr); ierr = ISDestroy(&subset_n);CHKERRQ(ierr); ierr = ISRenumber(subset,NULL,NULL,&subset_n);CHKERRQ(ierr); ierr = ISDestroy(&subset);CHKERRQ(ierr); ierr = ISGetLocalSize(subset_n,&k);CHKERRQ(ierr); if (k != graph->ncc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Invalid size of new subset! %D != %D",k,graph->ncc); ierr = ISGetIndices(subset_n,&is_indices);CHKERRQ(ierr); ierr = PetscMemcpy(graph->subset_ref_node,is_indices,graph->ncc*sizeof(PetscInt));CHKERRQ(ierr); ierr = ISRestoreIndices(subset_n,&is_indices);CHKERRQ(ierr); ierr = ISDestroy(&subset_n);CHKERRQ(ierr); /* free workspace */ graph->setupcalled = PETSC_TRUE; PetscFunctionReturn(0); }
PetscErrorCode PCBDDCGraphComputeConnectedComponentsLocal(PCBDDCGraph graph) { PetscInt ncc,cum_queue,n; PetscMPIInt commsize; PetscErrorCode ierr; PetscFunctionBegin; if (!graph->setupcalled) SETERRQ(PetscObjectComm((PetscObject)graph->l2gmap),PETSC_ERR_ORDER,"PCBDDCGraphSetUp should be called first"); /* quiet return if there isn't any local info */ if (!graph->xadj && !graph->n_local_subs) { PetscFunctionReturn(0); } /* reset any previous search of connected components */ ierr = PetscBTMemzero(graph->nvtxs,graph->touched);CHKERRQ(ierr); ierr = MPI_Comm_size(PetscObjectComm((PetscObject)graph->l2gmap),&commsize);CHKERRQ(ierr); if (commsize > graph->commsizelimit) { PetscInt i; for (i=0;i<graph->nvtxs;i++) { if (graph->special_dof[i] == PCBDDCGRAPH_DIRICHLET_MARK || !graph->count[i]) { ierr = PetscBTSet(graph->touched,i);CHKERRQ(ierr); } } } /* begin search for connected components */ cum_queue = 0; ncc = 0; for (n=0;n<graph->n_subsets;n++) { PetscInt pid = n+1; /* partition labeled by 0 is discarded */ PetscInt found = 0,prev = 0,first = 0,ncc_pid = 0; while (found != graph->subset_size[n]) { PetscInt added = 0; if (!prev) { /* search for new starting dof */ while (PetscBTLookup(graph->touched,graph->subset_idxs[n][first])) first++; ierr = PetscBTSet(graph->touched,graph->subset_idxs[n][first]);CHKERRQ(ierr); graph->queue[cum_queue] = graph->subset_idxs[n][first]; graph->cptr[ncc] = cum_queue; prev = 1; cum_queue++; found++; ncc_pid++; ncc++; } ierr = PCBDDCGraphComputeCC_Private(graph,pid,graph->queue + cum_queue,prev,&added);CHKERRQ(ierr); if (!added) { graph->subset_ncc[n] = ncc_pid; graph->cptr[ncc] = cum_queue; } prev = added; found += added; cum_queue += added; if (added && found == graph->subset_size[n]) { graph->subset_ncc[n] = ncc_pid; graph->cptr[ncc] = cum_queue; } } } graph->ncc = ncc; graph->queue_sorted = PETSC_FALSE; PetscFunctionReturn(0); }
PETSC_STATIC_INLINE PetscErrorCode PCBDDCGraphComputeCC_Private(PCBDDCGraph graph,PetscInt pid,PetscInt* queue_tip,PetscInt n_prev,PetscInt* n_added) { PetscInt i,j,n; PetscInt *xadj = graph->xadj,*adjncy = graph->adjncy; PetscBT touched = graph->touched; PetscBool havecsr = (PetscBool)(!!xadj); PetscBool havesubs = (PetscBool)(!!graph->n_local_subs); PetscErrorCode ierr; PetscFunctionBegin; n = 0; if (havecsr && !havesubs) { for (i=-n_prev;i<0;i++) { PetscInt start_dof = queue_tip[i]; /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs */ if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) { for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */ PetscInt dof = graph->subset_idxs[pid-1][j]; if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) { ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); queue_tip[n] = dof; n++; } } } else { for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { PetscInt dof = adjncy[j]; if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid) { ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); queue_tip[n] = dof; n++; } } } } } else if (havecsr && havesubs) { PetscInt sid = graph->local_subs[queue_tip[-n_prev]]; for (i=-n_prev;i<0;i++) { PetscInt start_dof = queue_tip[i]; /* we assume that if a dof has a size 1 adjacency list and the corresponding entry is negative, it is connected to all dofs belonging to the local sub */ if (xadj[start_dof+1]-xadj[start_dof] == 1 && adjncy[xadj[start_dof]] < 0) { for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */ PetscInt dof = graph->subset_idxs[pid-1][j]; if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) { ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); queue_tip[n] = dof; n++; } } } else { for (j=xadj[start_dof];j<xadj[start_dof+1];j++) { PetscInt dof = adjncy[j]; if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) { ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); queue_tip[n] = dof; n++; } } } } } else { /* sub info only */ PetscInt sid = graph->local_subs[queue_tip[-n_prev]]; for (j=0;j<graph->subset_size[pid-1];j++) { /* pid \in [1,graph->n_subsets] */ PetscInt dof = graph->subset_idxs[pid-1][j]; if (!PetscBTLookup(touched,dof) && graph->subset[dof] == pid && graph->local_subs[dof] == sid) { ierr = PetscBTSet(touched,dof);CHKERRQ(ierr); queue_tip[n] = dof; n++; } } } *n_added = n; PetscFunctionReturn(0); }
PetscErrorCode PCBDDCGraphComputeConnectedComponents(PCBDDCGraph graph) { PetscBool adapt_interface_reduced; MPI_Comm interface_comm; PetscMPIInt size; PetscInt i; PetscErrorCode ierr; PetscFunctionBegin; /* compute connected components locally */ ierr = PetscObjectGetComm((PetscObject)(graph->l2gmap),&interface_comm);CHKERRQ(ierr); ierr = PCBDDCGraphComputeConnectedComponentsLocal(graph);CHKERRQ(ierr); /* check consistency of connected components among neighbouring subdomains -> it adapt them in case it is needed */ ierr = MPI_Comm_size(interface_comm,&size);CHKERRQ(ierr); adapt_interface_reduced = PETSC_FALSE; if (size > 1) { PetscInt i; PetscBool adapt_interface = PETSC_FALSE; for (i=0;i<graph->n_subsets;i++) { /* We are not sure that on a given subset of the local interface, with two connected components, the latters be the same among sharing subdomains */ if (graph->subset_ncc[i] > 1) { adapt_interface = PETSC_TRUE; break; } } ierr = MPIU_Allreduce(&adapt_interface,&adapt_interface_reduced,1,MPIU_BOOL,MPI_LOR,interface_comm);CHKERRQ(ierr); } if (graph->n_subsets && adapt_interface_reduced) { PetscBT subset_cc_adapt; MPI_Request *send_requests,*recv_requests; PetscInt *send_buffer,*recv_buffer; PetscInt sum_requests,start_of_recv,start_of_send; PetscInt *cum_recv_counts; PetscInt *labels; PetscInt ncc,cum_queue,mss,mns,j,k,s; PetscInt **refine_buffer=NULL,*private_labels = NULL; ierr = PetscMalloc1(graph->nvtxs,&labels);CHKERRQ(ierr); ierr = PetscMemzero(labels,graph->nvtxs*sizeof(*labels));CHKERRQ(ierr); for (i=0;i<graph->ncc;i++) for (j=graph->cptr[i];j<graph->cptr[i+1];j++) labels[graph->queue[j]] = i; /* allocate some space */ ierr = PetscMalloc1(graph->n_subsets+1,&cum_recv_counts);CHKERRQ(ierr); ierr = PetscMemzero(cum_recv_counts,(graph->n_subsets+1)*sizeof(*cum_recv_counts));CHKERRQ(ierr); /* first count how many neighbours per connected component I will receive from */ cum_recv_counts[0] = 0; for (i=0;i<graph->n_subsets;i++) cum_recv_counts[i+1] = cum_recv_counts[i]+graph->count[graph->subset_idxs[i][0]]; ierr = PetscMalloc1(cum_recv_counts[graph->n_subsets],&recv_buffer);CHKERRQ(ierr); ierr = PetscMalloc2(cum_recv_counts[graph->n_subsets],&send_requests,cum_recv_counts[graph->n_subsets],&recv_requests);CHKERRQ(ierr); for (i=0;i<cum_recv_counts[graph->n_subsets];i++) { send_requests[i] = MPI_REQUEST_NULL; recv_requests[i] = MPI_REQUEST_NULL; } /* exchange with my neighbours the number of my connected components on the subset of interface */ sum_requests = 0; for (i=0;i<graph->n_subsets;i++) { PetscMPIInt neigh,tag; PetscInt count,*neighs; count = graph->count[graph->subset_idxs[i][0]]; neighs = graph->neighbours_set[graph->subset_idxs[i][0]]; ierr = PetscMPIIntCast(2*graph->subset_ref_node[i],&tag);CHKERRQ(ierr); for (k=0;k<count;k++) { ierr = PetscMPIIntCast(neighs[k],&neigh);CHKERRQ(ierr); ierr = MPI_Isend(&graph->subset_ncc[i],1,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr); ierr = MPI_Irecv(&recv_buffer[sum_requests],1,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr); sum_requests++; } } ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); /* determine the subsets I have to adapt (those having more than 1 cc) */ ierr = PetscBTCreate(graph->n_subsets,&subset_cc_adapt);CHKERRQ(ierr); ierr = PetscBTMemzero(graph->n_subsets,subset_cc_adapt);CHKERRQ(ierr); for (i=0;i<graph->n_subsets;i++) { if (graph->subset_ncc[i] > 1) { ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr); continue; } for (j=cum_recv_counts[i];j<cum_recv_counts[i+1];j++){ if (recv_buffer[j] > 1) { ierr = PetscBTSet(subset_cc_adapt,i);CHKERRQ(ierr); break; } } } ierr = PetscFree(recv_buffer);CHKERRQ(ierr); /* determine send/recv buffers sizes */ j = 0; mss = 0; for (i=0;i<graph->n_subsets;i++) { if (PetscBTLookup(subset_cc_adapt,i)) { j += graph->subset_size[i]; mss = PetscMax(graph->subset_size[i],mss); } } k = 0; mns = 0; for (i=0;i<graph->n_subsets;i++) { if (PetscBTLookup(subset_cc_adapt,i)) { k += (cum_recv_counts[i+1]-cum_recv_counts[i])*graph->subset_size[i]; mns = PetscMax(cum_recv_counts[i+1]-cum_recv_counts[i],mns); } } ierr = PetscMalloc2(j,&send_buffer,k,&recv_buffer);CHKERRQ(ierr); /* fill send buffer (order matters: subset_idxs ordered by global ordering) */ j = 0; for (i=0;i<graph->n_subsets;i++) if (PetscBTLookup(subset_cc_adapt,i)) for (k=0;k<graph->subset_size[i];k++) send_buffer[j++] = labels[graph->subset_idxs[i][k]]; /* now exchange the data */ start_of_recv = 0; start_of_send = 0; sum_requests = 0; for (i=0;i<graph->n_subsets;i++) { if (PetscBTLookup(subset_cc_adapt,i)) { PetscMPIInt neigh,tag; PetscInt size_of_send = graph->subset_size[i]; j = graph->subset_idxs[i][0]; ierr = PetscMPIIntCast(2*graph->subset_ref_node[i]+1,&tag);CHKERRQ(ierr); for (k=0;k<graph->count[j];k++) { ierr = PetscMPIIntCast(graph->neighbours_set[j][k],&neigh);CHKERRQ(ierr); ierr = MPI_Isend(&send_buffer[start_of_send],size_of_send,MPIU_INT,neigh,tag,interface_comm,&send_requests[sum_requests]);CHKERRQ(ierr); ierr = MPI_Irecv(&recv_buffer[start_of_recv],size_of_send,MPIU_INT,neigh,tag,interface_comm,&recv_requests[sum_requests]);CHKERRQ(ierr); start_of_recv += size_of_send; sum_requests++; } start_of_send += size_of_send; } } ierr = MPI_Waitall(sum_requests,recv_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); /* refine connected components */ start_of_recv = 0; /* allocate some temporary space */ if (mss) { ierr = PetscMalloc1(mss,&refine_buffer);CHKERRQ(ierr); ierr = PetscMalloc2(mss*(mns+1),&refine_buffer[0],mss,&private_labels);CHKERRQ(ierr); } ncc = 0; cum_queue = 0; graph->cptr[0] = 0; for (i=0;i<graph->n_subsets;i++) { if (PetscBTLookup(subset_cc_adapt,i)) { PetscInt subset_counter = 0; PetscInt sharingprocs = cum_recv_counts[i+1]-cum_recv_counts[i]+1; /* count myself */ PetscInt buffer_size = graph->subset_size[i]; /* compute pointers */ for (j=1;j<buffer_size;j++) refine_buffer[j] = refine_buffer[j-1] + sharingprocs; /* analyze contributions from subdomains that share the i-th subset The stricture of refine_buffer is suitable to find intersections of ccs among sharingprocs. supposing the current subset is shared by 3 processes and has dimension 5 with global dofs 0,1,2,3,4 (local 0,4,3,1,2) sharing procs connected components: neigh 0: [0 1 4], [2 3], labels [4,7] (2 connected components) neigh 1: [0 1], [2 3 4], labels [3 2] (2 connected components) neigh 2: [0 4], [1], [2 3], labels [1 5 6] (3 connected components) refine_buffer will be filled as: [ 4, 3, 1; 4, 2, 1; 7, 2, 6; 4, 3, 5; 7, 2, 6; ]; The connected components in local ordering are [0], [1], [2 3], [4] */ /* fill temp_buffer */ for (k=0;k<buffer_size;k++) refine_buffer[k][0] = labels[graph->subset_idxs[i][k]]; for (j=0;j<sharingprocs-1;j++) { for (k=0;k<buffer_size;k++) refine_buffer[k][j+1] = recv_buffer[start_of_recv+k]; start_of_recv += buffer_size; } ierr = PetscMemzero(private_labels,buffer_size*sizeof(PetscInt));CHKERRQ(ierr); for (j=0;j<buffer_size;j++) { if (!private_labels[j]) { /* found a new cc */ PetscBool same_set; graph->cptr[ncc] = cum_queue; ncc++; subset_counter++; private_labels[j] = subset_counter; graph->queue[cum_queue++] = graph->subset_idxs[i][j]; for (k=j+1;k<buffer_size;k++) { /* check for other nodes in new cc */ same_set = PETSC_TRUE; for (s=0;s<sharingprocs;s++) { if (refine_buffer[j][s] != refine_buffer[k][s]) { same_set = PETSC_FALSE; break; } } if (same_set) { private_labels[k] = subset_counter; graph->queue[cum_queue++] = graph->subset_idxs[i][k]; } } } } graph->cptr[ncc] = cum_queue; graph->subset_ncc[i] = subset_counter; graph->queue_sorted = PETSC_FALSE; } else { /* this subset does not need to be adapted */ ierr = PetscMemcpy(graph->queue+cum_queue,graph->subset_idxs[i],graph->subset_size[i]*sizeof(PetscInt));CHKERRQ(ierr); ncc++; cum_queue += graph->subset_size[i]; graph->cptr[ncc] = cum_queue; } } graph->cptr[ncc] = cum_queue; graph->ncc = ncc; if (mss) { ierr = PetscFree2(refine_buffer[0],private_labels);CHKERRQ(ierr); ierr = PetscFree(refine_buffer);CHKERRQ(ierr); } ierr = PetscFree(labels);CHKERRQ(ierr); ierr = MPI_Waitall(sum_requests,send_requests,MPI_STATUSES_IGNORE);CHKERRQ(ierr); ierr = PetscFree2(send_requests,recv_requests);CHKERRQ(ierr); ierr = PetscFree2(send_buffer,recv_buffer);CHKERRQ(ierr); ierr = PetscFree(cum_recv_counts);CHKERRQ(ierr); ierr = PetscBTDestroy(&subset_cc_adapt);CHKERRQ(ierr); } /* Determine if we are in 2D or 3D */ if (!graph->twodimset) { PetscBool twodim = PETSC_TRUE; for (i=0;i<graph->ncc;i++) { PetscInt repdof = graph->queue[graph->cptr[i]]; PetscInt ccsize = graph->cptr[i+1]-graph->cptr[i]; if (graph->count[repdof] > 1 && ccsize > graph->custom_minimal_size) { twodim = PETSC_FALSE; break; } } ierr = MPIU_Allreduce(&twodim,&graph->twodim,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)graph->l2gmap));CHKERRQ(ierr); graph->twodimset = PETSC_TRUE; } PetscFunctionReturn(0); }