/*@ PetscSplitOwnership - Given a global (or local) length determines a local (or global) length via a simple formula Collective on MPI_Comm (if n or N is PETSC_DECIDE) Input Parameters: + comm - MPI communicator that shares the object being divided . n - local length (or PETSC_DECIDE to have it set) - N - global length (or PETSC_DECIDE) Level: developer Notes: n and N cannot be both PETSC_DECIDE If one processor calls this with n or N of PETSC_DECIDE then all processors must. Otherwise, an error is thrown in debug mode while the program will hang in optimized (i.e. configured --with-debugging=0) mode. .seealso: PetscSplitOwnershipBlock() @*/ PetscErrorCode PetscSplitOwnership(MPI_Comm comm,PetscInt *n,PetscInt *N) { PetscErrorCode ierr; PetscMPIInt size,rank; PetscFunctionBegin; if (*N == PETSC_DECIDE && *n == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Both n and N cannot be PETSC_DECIDE\n likely a call to VecSetSizes() or MatSetSizes() is wrong.\nSee https://www.mcs.anl.gov/petsc/documentation/faq.html#split"); #if defined(PETSC_USE_DEBUG) { PetscMPIInt l[2],g[2]; l[0] = (*n == PETSC_DECIDE) ? 1 : 0; l[1] = (*N == PETSC_DECIDE) ? 1 : 0; ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); ierr = MPIU_Allreduce(l,g,2,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); if (g[0] && g[0] != size) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"All processes must supply PETSC_DECIDE for local size"); if (g[1] && g[1] != size) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"All processes must supply PETSC_DECIDE for global size"); } #endif if (*N == PETSC_DECIDE) { ierr = MPIU_Allreduce(n,N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); } else if (*n == PETSC_DECIDE) { ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); *n = *N/size + ((*N % size) > rank); #if defined(PETSC_USE_DEBUG) } else { PetscInt tmp; ierr = MPIU_Allreduce(n,&tmp,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); if (tmp != *N) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Sum of local lengths %D does not equal global length %D, my local length %D\n likely a call to VecSetSizes() or MatSetSizes() is wrong.\nSee https://www.mcs.anl.gov/petsc/documentation/faq.html#split",tmp,*N,*n); #endif } PetscFunctionReturn(0); }
PetscErrorCode VecMin_MPI(Vec xin,PetscInt *idx,PetscReal *z) { PetscErrorCode ierr; PetscReal work; PetscFunctionBegin; /* Find the local Min */ ierr = VecMin_Seq(xin,idx,&work);CHKERRQ(ierr); /* Find the global Min */ if (!idx) { ierr = MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); } else { PetscReal work2[2],z2[2]; PetscInt rstart; ierr = VecGetOwnershipRange(xin,&rstart,NULL);CHKERRQ(ierr); work2[0] = work; work2[1] = *idx + rstart; ierr = MPIU_Allreduce(work2,z2,2,MPIU_REAL,MPIU_MININDEX_OP,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); *z = z2[0]; *idx = (PetscInt)z2[1]; } PetscFunctionReturn(0); }
PetscErrorCode KSPAGMRESRoddecInitNeighboor(KSP ksp) { MPI_Comm comm; KSP_AGMRES *agmres = (KSP_AGMRES*)(ksp->data); PetscErrorCode ierr; PetscMPIInt First, Last, rank, size; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject)agmres->vecs[0], &comm);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); ierr = MPIU_Allreduce(&rank, &First, 1, MPIU_INT, MPI_MIN, comm);CHKERRQ(ierr); ierr = MPIU_Allreduce(&rank, &Last, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr); if ((rank != Last) && (!rank)) { agmres->Ileft = rank - 1; agmres->Iright = rank + 1; } else { if (rank == Last) { agmres->Ileft = rank - 1; agmres->Iright = First; } else { { agmres->Ileft = Last; agmres->Iright = rank + 1; } } } agmres->rank = rank; agmres->size = size; agmres->First = First; agmres->Last = Last; PetscFunctionReturn(0); }
/*@ DMDAGetBoundingBox - Returns the global bounding box for the DMDA. Collective on DMDA Input Parameter: . dm - the DM Output Parameters: + gmin - global minimum coordinates (length dim, optional) - gmax - global maximim coordinates (length dim, optional) Level: beginner .keywords: distributed array, get, coordinates .seealso: DMDAGetCoordinateDA(), DMGetCoordinates(), DMDAGetLocalBoundingBox() @*/ PetscErrorCode DMDAGetBoundingBox(DM dm,PetscReal gmin[],PetscReal gmax[]) { PetscErrorCode ierr; PetscMPIInt count; PetscReal lmin[3],lmax[3]; PetscFunctionBegin; PetscValidHeaderSpecific(dm,DM_CLASSID,1); ierr = PetscMPIIntCast(dm->dim,&count);CHKERRQ(ierr); ierr = DMDAGetLocalBoundingBox(dm,lmin,lmax);CHKERRQ(ierr); if (gmin) {ierr = MPIU_Allreduce(lmin,gmin,count,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);} if (gmax) {ierr = MPIU_Allreduce(lmax,gmax,count,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);} PetscFunctionReturn(0); }
PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm,Vec *vec) { PetscSection gSection; PetscInt localSize, bs, blockSize = -1, pStart, pEnd, p; PetscErrorCode ierr; PetscFunctionBegin; ierr = DMGetDefaultGlobalSection(dm, &gSection);CHKERRQ(ierr); ierr = PetscSectionGetChart(gSection, &pStart, &pEnd);CHKERRQ(ierr); for (p = pStart; p < pEnd; ++p) { PetscInt dof, cdof; ierr = PetscSectionGetDof(gSection, p, &dof);CHKERRQ(ierr); ierr = PetscSectionGetConstraintDof(gSection, p, &cdof);CHKERRQ(ierr); if ((blockSize < 0) && (dof > 0) && (dof-cdof > 0)) blockSize = dof-cdof; if ((dof > 0) && (dof-cdof != blockSize)) { blockSize = 1; break; } } if (blockSize < 0) blockSize = PETSC_MAX_INT; ierr = MPIU_Allreduce(&blockSize, &bs, 1, MPIU_INT, MPI_MIN, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); if (blockSize == PETSC_MAX_INT) blockSize = 1; /* Everyone was empty */ ierr = PetscSectionGetConstrainedStorageSize(gSection, &localSize);CHKERRQ(ierr); if (localSize%blockSize) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %d and local storage size %d", blockSize, localSize); ierr = VecCreate(PetscObjectComm((PetscObject)dm), vec);CHKERRQ(ierr); ierr = VecSetSizes(*vec, localSize, PETSC_DETERMINE);CHKERRQ(ierr); ierr = VecSetBlockSize(*vec, bs);CHKERRQ(ierr); ierr = VecSetType(*vec,dm->vectype);CHKERRQ(ierr); ierr = VecSetDM(*vec, dm);CHKERRQ(ierr); /* ierr = VecSetLocalToGlobalMapping(*vec, dm->ltogmap);CHKERRQ(ierr); */ PetscFunctionReturn(0); }
/*@ MatComputeBandwidth - Calculate the full bandwidth of the matrix, meaning the width 2k+1 where k diagonals on either side are sufficient to contain all the matrix nonzeros. Collective on Mat Input Parameters: + A - The Mat - fraction - An optional percentage of the Frobenius norm of the matrix that the bandwidth should enclose Output Parameter: . bw - The matrix bandwidth Level: beginner .seealso: DMPlexCreate(), DMPlexSetConeSize(), DMPlexSetChart() @*/ PetscErrorCode MatComputeBandwidth(Mat A, PetscReal fraction, PetscInt *bw) { PetscInt lbw[2] = {0, 0}, gbw[2]; PetscInt rStart, rEnd, r; PetscErrorCode ierr; PetscFunctionBegin; PetscValidHeaderSpecific(A, MAT_CLASSID, 1); PetscValidLogicalCollectiveReal(A,fraction,2); PetscValidPointer(bw, 3); if ((fraction > 0.0) && (fraction < 1.0)) SETERRQ(PetscObjectComm((PetscObject) A), PETSC_ERR_SUP, "We do not yet support a fractional bandwidth"); ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); for (r = rStart; r < rEnd; ++r) { const PetscInt *cols; PetscInt ncols; ierr = MatGetRow(A, r, &ncols, &cols, NULL);CHKERRQ(ierr); if (ncols) { lbw[0] = PetscMax(lbw[0], r - cols[0]); lbw[1] = PetscMax(lbw[1], cols[ncols-1] - r); } ierr = MatRestoreRow(A, r, &ncols, &cols, NULL);CHKERRQ(ierr); } ierr = MPIU_Allreduce(lbw, gbw, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject) A));CHKERRQ(ierr); *bw = 2*PetscMax(gbw[0], gbw[1]) + 1; PetscFunctionReturn(0); }
PetscErrorCode KSPMonitorRange_Private(KSP ksp,PetscInt it,PetscReal *per) { PetscErrorCode ierr; Vec resid; PetscReal rmax,pwork; PetscInt i,n,N; const PetscScalar *r; PetscFunctionBegin; ierr = KSPBuildResidual(ksp,NULL,NULL,&resid); CHKERRQ(ierr); ierr = VecNorm(resid,NORM_INFINITY,&rmax); CHKERRQ(ierr); ierr = VecGetLocalSize(resid,&n); CHKERRQ(ierr); ierr = VecGetSize(resid,&N); CHKERRQ(ierr); ierr = VecGetArrayRead(resid,&r); CHKERRQ(ierr); pwork = 0.0; for (i=0; i<n; i++) pwork += (PetscAbsScalar(r[i]) > .20*rmax); ierr = VecRestoreArrayRead(resid,&r); CHKERRQ(ierr); ierr = VecDestroy(&resid); CHKERRQ(ierr); ierr = MPIU_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ksp)); CHKERRQ(ierr); *per = *per/N; PetscFunctionReturn(0); }
/*@C PetscGatherNumberOfMessages - Computes the number of messages a node expects to receive Collective on MPI_Comm Input Parameters: + comm - Communicator . iflags - an array of integers of length sizeof(comm). A '1' in ilengths[i] represent a message from current node to ith node. Optionally NULL - ilengths - Non zero ilengths[i] represent a message to i of length ilengths[i]. Optionally NULL. Output Parameters: . nrecvs - number of messages received Level: developer Concepts: mpi utility Notes: With this info, the correct message lengths can be determined using PetscGatherMessageLengths() Either iflags or ilengths should be provided. If iflags is not provided (NULL) it can be computed from ilengths. If iflags is provided, ilengths is not required. .seealso: PetscGatherMessageLengths() @*/ PetscErrorCode PetscGatherNumberOfMessages(MPI_Comm comm,const PetscMPIInt iflags[],const PetscMPIInt ilengths[],PetscMPIInt *nrecvs) { PetscMPIInt size,rank,*recv_buf,i,*iflags_local = NULL,*iflags_localm = NULL; PetscErrorCode ierr; PetscFunctionBegin; ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); ierr = PetscMalloc2(size,&recv_buf,size,&iflags_localm);CHKERRQ(ierr); /* If iflags not provided, compute iflags from ilengths */ if (!iflags) { if (!ilengths) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Either iflags or ilengths should be provided"); iflags_local = iflags_localm; for (i=0; i<size; i++) { if (ilengths[i]) iflags_local[i] = 1; else iflags_local[i] = 0; } } else iflags_local = (PetscMPIInt*) iflags; /* Post an allreduce to determine the numer of messages the current node will receive */ ierr = MPIU_Allreduce(iflags_local,recv_buf,size,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); *nrecvs = recv_buf[rank]; ierr = PetscFree2(recv_buf,iflags_localm);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode MatGetDiagonalHermitian_Normal(Mat N,Vec v) { Mat_Normal *Na = (Mat_Normal*)N->data; Mat A = Na->A; PetscErrorCode ierr; PetscInt i,j,rstart,rend,nnz; const PetscInt *cols; PetscScalar *diag,*work,*values; const PetscScalar *mvalues; PetscFunctionBegin; ierr = PetscMalloc2(A->cmap->N,&diag,A->cmap->N,&work);CHKERRQ(ierr); ierr = PetscMemzero(work,A->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); for (i=rstart; i<rend; i++) { ierr = MatGetRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); for (j=0; j<nnz; j++) { work[cols[j]] += mvalues[j]*PetscConj(mvalues[j]); } ierr = MatRestoreRow(A,i,&nnz,&cols,&mvalues);CHKERRQ(ierr); } ierr = MPIU_Allreduce(work,diag,A->cmap->N,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)N));CHKERRQ(ierr); rstart = N->cmap->rstart; rend = N->cmap->rend; ierr = VecGetArray(v,&values);CHKERRQ(ierr); ierr = PetscMemcpy(values,diag+rstart,(rend-rstart)*sizeof(PetscScalar));CHKERRQ(ierr); ierr = VecRestoreArray(v,&values);CHKERRQ(ierr); ierr = PetscFree2(diag,work);CHKERRQ(ierr); ierr = VecScale(v,Na->scale);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@ MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or zero, decreases h until this is satisfied. Logically Collective on Vec Input Parameters: + U - base vector that is added to . a - vector that is added . h - scaling factor on a - dummy - context variable (unused) Options Database Keys: . -mat_mffd_check_positivity Level: advanced Notes: This is rarely used directly, rather it is passed as an argument to MatMFFDSetCheckh() .seealso: MatMFFDSetCheckh() @*/ PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) { PetscReal val, minval; PetscScalar *u_vec, *a_vec; PetscErrorCode ierr; PetscInt i,n; MPI_Comm comm; PetscFunctionBegin; PetscValidHeaderSpecific(U,VEC_CLASSID,2); PetscValidHeaderSpecific(a,VEC_CLASSID,3); PetscValidPointer(h,4); ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); minval = PetscAbsScalar(*h)*PetscRealConstant(1.01); for (i=0; i<n; i++) { if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { val = PetscAbsScalar(u_vec[i]/a_vec[i]); if (val < minval) minval = val; } } ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); if (val <= PetscAbsScalar(*h)) { ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr); if (PetscRealPart(*h) > 0.0) *h = 0.99*val; else *h = -0.99*val; } PetscFunctionReturn(0); }
/*@ DMDANaturalAllToGlobalCreate - Creates a scatter context that maps from a copy of the entire vector on each processor to its local part in the global vector. Collective on DMDA Input Parameter: . da - the distributed array context Output Parameter: . scatter - the scatter context Level: advanced .keywords: distributed array, global to local, begin, coarse problem .seealso: DMDAGlobalToNaturalEnd(), DMLocalToGlobalBegin(), DMDACreate2d(), DMGlobalToLocalBegin(), DMGlobalToLocalEnd(), DMDACreateNaturalVector() @*/ PetscErrorCode DMDANaturalAllToGlobalCreate(DM da,VecScatter *scatter) { PetscErrorCode ierr; DM_DA *dd = (DM_DA*)da->data; PetscInt M,m = dd->Nlocal,start; IS from,to; Vec tmplocal,global; AO ao; PetscFunctionBegin; PetscValidHeaderSpecific(da,DM_CLASSID,1); PetscValidPointer(scatter,2); ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr); /* create the scatter context */ ierr = MPIU_Allreduce(&m,&M,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)da));CHKERRQ(ierr); ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)da),dd->w,m,PETSC_DETERMINE,0,&global);CHKERRQ(ierr); ierr = VecGetOwnershipRange(global,&start,NULL);CHKERRQ(ierr); ierr = ISCreateStride(PetscObjectComm((PetscObject)da),m,start,1,&from);CHKERRQ(ierr); ierr = AOPetscToApplicationIS(ao,from);CHKERRQ(ierr); ierr = ISCreateStride(PetscObjectComm((PetscObject)da),m,start,1,&to);CHKERRQ(ierr); ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,dd->w,M,0,&tmplocal);CHKERRQ(ierr); ierr = VecScatterCreate(tmplocal,from,global,to,scatter);CHKERRQ(ierr); ierr = VecDestroy(&tmplocal);CHKERRQ(ierr); ierr = VecDestroy(&global);CHKERRQ(ierr); ierr = ISDestroy(&from);CHKERRQ(ierr); ierr = ISDestroy(&to);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* Helper rutine to handle user postenvents and recording */ static PetscErrorCode TSPostEvent(TS ts,PetscReal t,Vec U) { PetscErrorCode ierr; TSEvent event = ts->event; PetscBool terminate = PETSC_FALSE; PetscBool restart = PETSC_FALSE; PetscInt i,ctr,stepnum; PetscBool ts_terminate,ts_restart; PetscBool forwardsolve = PETSC_TRUE; /* Flag indicating that TS is doing a forward solve */ PetscFunctionBegin; if (event->postevent) { PetscObjectState state_prev,state_post; ierr = PetscObjectStateGet((PetscObject)U,&state_prev);CHKERRQ(ierr); ierr = (*event->postevent)(ts,event->nevents_zero,event->events_zero,t,U,forwardsolve,event->ctx);CHKERRQ(ierr); ierr = PetscObjectStateGet((PetscObject)U,&state_post);CHKERRQ(ierr); if (state_prev != state_post) restart = PETSC_TRUE; } /* Handle termination events and step restart */ for (i=0; i<event->nevents_zero; i++) if (event->terminate[event->events_zero[i]]) terminate = PETSC_TRUE; ierr = MPIU_Allreduce(&terminate,&ts_terminate,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); if (ts_terminate) {ierr = TSSetConvergedReason(ts,TS_CONVERGED_EVENT);CHKERRQ(ierr);} event->status = ts_terminate ? TSEVENT_NONE : TSEVENT_RESET_NEXTSTEP; ierr = MPIU_Allreduce(&restart,&ts_restart,1,MPIU_BOOL,MPI_LOR,((PetscObject)ts)->comm);CHKERRQ(ierr); if (ts_restart) ts->steprestart = PETSC_TRUE; event->ptime_prev = t; /* Reset event residual functions as states might get changed by the postevent callback */ if (event->postevent) {ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr);} /* Cache current event residual functions */ for (i=0; i < event->nevents; i++) event->fvalue_prev[i] = event->fvalue[i]; /* Record the event in the event recorder */ ierr = TSGetTimeStepNumber(ts,&stepnum);CHKERRQ(ierr); ctr = event->recorder.ctr; if (ctr == event->recsize) { ierr = TSEventRecorderResize(event);CHKERRQ(ierr); } event->recorder.time[ctr] = t; event->recorder.stepnum[ctr] = stepnum; event->recorder.nevents[ctr] = event->nevents_zero; for (i=0; i<event->nevents_zero; i++) event->recorder.eventidx[ctr][i] = event->events_zero[i]; event->recorder.ctr++; PetscFunctionReturn(0); }
/*@ ISPartitioningToNumbering - Takes an ISPartitioning and on each processor generates an IS that contains a new global node number for each index based on the partitioing. Collective on IS Input Parameters . partitioning - a partitioning as generated by MatPartitioningApply() Output Parameter: . is - on each processor the index set that defines the global numbers (in the new numbering) for all the nodes currently (before the partitioning) on that processor Level: advanced .seealso: MatPartitioningCreate(), AOCreateBasic(), ISPartitioningCount() @*/ PetscErrorCode ISPartitioningToNumbering(IS part,IS *is) { MPI_Comm comm; PetscInt i,np,npt,n,*starts = NULL,*sums = NULL,*lsizes = NULL,*newi = NULL; const PetscInt *indices = NULL; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject)part,&comm);CHKERRQ(ierr); /* count the number of partitions, i.e., virtual processors */ ierr = ISGetLocalSize(part,&n);CHKERRQ(ierr); ierr = ISGetIndices(part,&indices);CHKERRQ(ierr); np = 0; for (i=0; i<n; i++) np = PetscMax(np,indices[i]); ierr = MPIU_Allreduce(&np,&npt,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); np = npt+1; /* so that it looks like a MPI_Comm_size output */ /* lsizes - number of elements of each partition on this particular processor sums - total number of "previous" nodes for any particular partition starts - global number of first element in each partition on this processor */ ierr = PetscMalloc3(np,&lsizes,np,&starts,np,&sums);CHKERRQ(ierr); ierr = PetscMemzero(lsizes,np*sizeof(PetscInt));CHKERRQ(ierr); for (i=0; i<n; i++) lsizes[indices[i]]++; ierr = MPIU_Allreduce(lsizes,sums,np,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); ierr = MPI_Scan(lsizes,starts,np,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); for (i=0; i<np; i++) starts[i] -= lsizes[i]; for (i=1; i<np; i++) { sums[i] += sums[i-1]; starts[i] += sums[i-1]; } /* For each local index give it the new global number */ ierr = PetscMalloc1(n,&newi);CHKERRQ(ierr); for (i=0; i<n; i++) newi[i] = starts[indices[i]]++; ierr = PetscFree3(lsizes,starts,sums);CHKERRQ(ierr); ierr = ISRestoreIndices(part,&indices);CHKERRQ(ierr); ierr = ISCreateGeneral(comm,n,newi,PETSC_OWN_POINTER,is);CHKERRQ(ierr); ierr = ISSetPermutation(*is);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* PetscSplitReductionApply - Actually do the communication required for a split phase reduction */ static PetscErrorCode PetscSplitReductionApply(PetscSplitReduction *sr) { PetscErrorCode ierr; PetscInt i,numops = sr->numopsbegin,*reducetype = sr->reducetype; PetscScalar *lvalues = sr->lvalues,*gvalues = sr->gvalues; PetscInt sum_flg = 0,max_flg = 0, min_flg = 0; MPI_Comm comm = sr->comm; PetscMPIInt size,cmul = sizeof(PetscScalar)/sizeof(PetscReal); PetscFunctionBegin; if (sr->numopsend > 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Cannot call this after VecxxxEnd() has been called"); ierr = PetscLogEventBarrierBegin(VEC_ReduceBarrier,0,0,0,0,comm);CHKERRQ(ierr); ierr = MPI_Comm_size(sr->comm,&size);CHKERRQ(ierr); if (size == 1) { ierr = PetscMemcpy(gvalues,lvalues,numops*sizeof(PetscScalar));CHKERRQ(ierr); } else { /* determine if all reductions are sum, max, or min */ for (i=0; i<numops; i++) { if (reducetype[i] == REDUCE_MAX) max_flg = 1; else if (reducetype[i] == REDUCE_SUM) sum_flg = 1; else if (reducetype[i] == REDUCE_MIN) min_flg = 1; else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in PetscSplitReduction() data structure, probably memory corruption"); } if (sum_flg + max_flg + min_flg > 1) { /* after all the entires in lvalues we store the reducetype flags to indicate to the reduction operations what are sums and what are max */ for (i=0; i<numops; i++) lvalues[numops+i] = reducetype[i]; ierr = MPIU_Allreduce(lvalues,gvalues,2*numops,MPIU_SCALAR,PetscSplitReduction_Op,comm);CHKERRQ(ierr); } else if (max_flg) { /* Compute max of real and imag parts separately, presumably only the real part is used */ ierr = MPIU_Allreduce((PetscReal*)lvalues,(PetscReal*)gvalues,cmul*numops,MPIU_REAL,MPIU_MAX,comm);CHKERRQ(ierr); } else if (min_flg) { ierr = MPIU_Allreduce((PetscReal*)lvalues,(PetscReal*)gvalues,cmul*numops,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); } else { ierr = MPIU_Allreduce(lvalues,gvalues,numops,MPIU_SCALAR,MPIU_SUM,comm);CHKERRQ(ierr); } } sr->state = STATE_END; sr->numopsend = 0; ierr = PetscLogEventBarrierEnd(VEC_ReduceBarrier,0,0,0,0,comm);CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@ ISPartitioningCount - Takes a ISPartitioning and determines the number of resulting elements on each (partition) process Collective on IS Input Parameters: + partitioning - a partitioning as generated by MatPartitioningApply() - len - length of the array count, this is the total number of partitions Output Parameter: . count - array of length size, to contain the number of elements assigned to each partition, where size is the number of partitions generated (see notes below). Level: advanced Notes: By default the number of partitions generated (and thus the length of count) is the size of the communicator associated with IS, but it can be set by MatPartitioningSetNParts. The resulting array of lengths can for instance serve as input of PCBJacobiSetTotalBlocks. .seealso: MatPartitioningCreate(), AOCreateBasic(), ISPartitioningToNumbering(), MatPartitioningSetNParts(), MatPartitioningApply() @*/ PetscErrorCode ISPartitioningCount(IS part,PetscInt len,PetscInt count[]) { MPI_Comm comm; PetscInt i,n,*lsizes; const PetscInt *indices; PetscErrorCode ierr; PetscMPIInt npp; PetscFunctionBegin; ierr = PetscObjectGetComm((PetscObject)part,&comm);CHKERRQ(ierr); if (len == PETSC_DEFAULT) { PetscMPIInt size; ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); len = (PetscInt) size; } /* count the number of partitions */ ierr = ISGetLocalSize(part,&n);CHKERRQ(ierr); ierr = ISGetIndices(part,&indices);CHKERRQ(ierr); #if defined(PETSC_USE_DEBUG) { PetscInt np = 0,npt; for (i=0; i<n; i++) np = PetscMax(np,indices[i]); ierr = MPIU_Allreduce(&np,&npt,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); np = npt+1; /* so that it looks like a MPI_Comm_size output */ if (np > len) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Length of count array %D is less than number of partitions %D",len,np); } #endif /* lsizes - number of elements of each partition on this particular processor sums - total number of "previous" nodes for any particular partition starts - global number of first element in each partition on this processor */ ierr = PetscCalloc1(len,&lsizes);CHKERRQ(ierr); for (i=0; i<n; i++) lsizes[indices[i]]++; ierr = ISRestoreIndices(part,&indices);CHKERRQ(ierr); ierr = PetscMPIIntCast(len,&npp);CHKERRQ(ierr); ierr = MPIU_Allreduce(lsizes,count,npp,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); ierr = PetscFree(lsizes);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode MatIsSymmetric_IS(Mat A,PetscReal tol,PetscBool *flg) { PetscErrorCode ierr; Mat_IS *matis = (Mat_IS*)A->data; PetscBool local_sym; PetscFunctionBegin; ierr = MatIsSymmetric(matis->A,tol,&local_sym);CHKERRQ(ierr); ierr = MPIU_Allreduce(&local_sym,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); PetscFunctionReturn(0); }
static PetscErrorCode VecDot_MPI(Vec xin,Vec yin,PetscScalar *z) { PetscScalar sum,work; PetscErrorCode ierr; PetscFunctionBegin; ierr = VecDot_Seq(xin,yin,&work);CHKERRQ(ierr); ierr = MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); *z = sum; PetscFunctionReturn(0); }
PetscErrorCode VecNorm_MPI(Vec xin,NormType type,PetscReal *z) { PetscReal sum,work = 0.0; const PetscScalar *xx; PetscErrorCode ierr; PetscInt n = xin->map->n; PetscBLASInt one = 1,bn; PetscFunctionBegin; ierr = PetscBLASIntCast(n,&bn);CHKERRQ(ierr); if (type == NORM_2 || type == NORM_FROBENIUS) { ierr = VecGetArrayRead(xin,&xx);CHKERRQ(ierr); work = PetscRealPart(BLASdot_(&bn,xx,&one,xx,&one)); ierr = VecRestoreArrayRead(xin,&xx);CHKERRQ(ierr); ierr = MPIU_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); *z = PetscSqrtReal(sum); ierr = PetscLogFlops(2.0*xin->map->n);CHKERRQ(ierr); } else if (type == NORM_1) { /* Find the local part */ ierr = VecNorm_Seq(xin,NORM_1,&work);CHKERRQ(ierr); /* Find the global max */ ierr = MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); } else if (type == NORM_INFINITY) { /* Find the local max */ ierr = VecNorm_Seq(xin,NORM_INFINITY,&work);CHKERRQ(ierr); /* Find the global max */ ierr = MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); } else if (type == NORM_1_AND_2) { PetscReal temp[2]; ierr = VecNorm_Seq(xin,NORM_1,temp);CHKERRQ(ierr); ierr = VecNorm_Seq(xin,NORM_2,temp+1);CHKERRQ(ierr); temp[1] = temp[1]*temp[1]; ierr = MPIU_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); z[1] = PetscSqrtReal(z[1]); } PetscFunctionReturn(0); }
static PetscErrorCode MPIPetsc_Iallreduce(void *sendbuf,void *recvbuf,PetscMPIInt count,MPI_Datatype datatype,MPI_Op op,MPI_Comm comm,MPI_Request *request) { PETSC_UNUSED PetscErrorCode ierr; PetscFunctionBegin; #if defined(PETSC_HAVE_MPI_IALLREDUCE) ierr = MPI_Iallreduce(sendbuf,recvbuf,count,datatype,op,comm,request);CHKERRQ(ierr); #elif defined(PETSC_HAVE_MPIX_IALLREDUCE) ierr = MPIX_Iallreduce(sendbuf,recvbuf,count,datatype,op,comm,request);CHKERRQ(ierr); #else ierr = MPIU_Allreduce(sendbuf,recvbuf,count,datatype,op,comm);CHKERRQ(ierr); *request = MPI_REQUEST_NULL; #endif PetscFunctionReturn(0); }
PetscErrorCode VecMax_MPI(Vec xin,PetscInt *idx,PetscReal *z) { PetscErrorCode ierr; PetscReal work; PetscFunctionBegin; /* Find the local max */ ierr = VecMax_Seq(xin,idx,&work);CHKERRQ(ierr); /* Find the global max */ if (!idx) { ierr = MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); } else { PetscReal work2[2],z2[2]; PetscInt rstart; rstart = xin->map->rstart; work2[0] = work; work2[1] = *idx + rstart; ierr = MPIU_Allreduce(work2,z2,2,MPIU_REAL,MPIU_MAXINDEX_OP,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); *z = z2[0]; *idx = (PetscInt)z2[1]; } PetscFunctionReturn(0); }
/*@ PetscCommBuildTwoSidedSetType - set algorithm to use when building two-sided communication Logically Collective Input Arguments: + comm - PETSC_COMM_WORLD - twosided - algorithm to use in subsequent calls to PetscCommBuildTwoSided() Level: developer Note: This option is currently global, but could be made per-communicator. .seealso: PetscCommBuildTwoSided(), PetscCommBuildTwoSidedGetType() @*/ PetscErrorCode PetscCommBuildTwoSidedSetType(MPI_Comm comm,PetscBuildTwoSidedType twosided) { PetscFunctionBegin; #if defined(PETSC_USE_DEBUG) { /* We don't have a PetscObject so can't use PetscValidLogicalCollectiveEnum */ PetscMPIInt ierr; PetscMPIInt b1[2],b2[2]; b1[0] = -(PetscMPIInt)twosided; b1[1] = (PetscMPIInt)twosided; ierr = MPIU_Allreduce(b1,b2,2,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); if (-b2[0] != b2[1]) SETERRQ(comm,PETSC_ERR_ARG_WRONG,"Enum value must be same on all processes"); } #endif _twosided_type = twosided; PetscFunctionReturn(0); }
PetscErrorCode VecMDot_MPI(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z) { PetscScalar awork[128],*work = awork; PetscErrorCode ierr; PetscFunctionBegin; if (nv > 128) { ierr = PetscMalloc1(nv,&work);CHKERRQ(ierr); } ierr = VecMDot_Seq(xin,nv,y,work);CHKERRQ(ierr); ierr = MPIU_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));CHKERRQ(ierr); if (nv > 128) { ierr = PetscFree(work);CHKERRQ(ierr); } PetscFunctionReturn(0); }
/*@ PetscSplitOwnershipBlock - Given a global (or local) length determines a local (or global) length via a simple formula. Splits so each processors local size is divisible by the block size. Collective on MPI_Comm (if N is PETSC_DECIDE) Input Parameters: + comm - MPI communicator that shares the object being divided . bs - block size . n - local length (or PETSC_DECIDE to have it set) - N - global length (or PETSC_DECIDE) Level: developer Notes: n and N cannot be both PETSC_DECIDE If one processor calls this with N of PETSC_DECIDE then all processors must, otherwise the program will hang. .seealso: PetscSplitOwnership() @*/ PetscErrorCode PetscSplitOwnershipBlock(MPI_Comm comm,PetscInt bs,PetscInt *n,PetscInt *N) { PetscErrorCode ierr; PetscMPIInt size,rank; PetscFunctionBegin; if (*N == PETSC_DECIDE && *n == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Both n and N cannot be PETSC_DECIDE"); if (*N == PETSC_DECIDE) { if (*n % bs != 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"local size %D not divisible by block size %D",*n,bs); ierr = MPIU_Allreduce(n,N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); } else if (*n == PETSC_DECIDE) { PetscInt Nbs = *N/bs; ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); *n = bs*(Nbs/size + ((Nbs % size) > rank)); } PetscFunctionReturn(0); }
PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) { Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; PetscErrorCode ierr; PetscBool flag; PetscFunctionBegin; /* If the matrix dimensions are not equal,or no of nonzeros */ if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { flag = PETSC_FALSE; } /* if the a->i are the same */ ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); /* if a->j are the same */ ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); ierr = MPIU_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); PetscFunctionReturn(0); }
/*@ MatChop - Set all values in the matrix less than the tolerance to zero Input Parameters: + A - The matrix - tol - The zero tolerance Output Parameters: . A - The chopped matrix Level: intermediate .seealso: MatCreate(), MatZeroEntries() @*/ PetscErrorCode MatChop(Mat A, PetscReal tol) { PetscScalar *newVals; PetscInt *newCols; PetscInt rStart, rEnd, numRows, maxRows, r, colMax = 0; PetscErrorCode ierr; PetscFunctionBegin; ierr = MatGetOwnershipRange(A, &rStart, &rEnd);CHKERRQ(ierr); for (r = rStart; r < rEnd; ++r) { PetscInt ncols; ierr = MatGetRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); colMax = PetscMax(colMax, ncols);CHKERRQ(ierr); ierr = MatRestoreRow(A, r, &ncols, NULL, NULL);CHKERRQ(ierr); } numRows = rEnd - rStart; ierr = MPIU_Allreduce(&numRows, &maxRows, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)A));CHKERRQ(ierr); ierr = PetscMalloc2(colMax,&newCols,colMax,&newVals);CHKERRQ(ierr); for (r = rStart; r < rStart+maxRows; ++r) { const PetscScalar *vals; const PetscInt *cols; PetscInt ncols, newcols, c; if (r < rEnd) { ierr = MatGetRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); for (c = 0; c < ncols; ++c) { newCols[c] = cols[c]; newVals[c] = PetscAbsScalar(vals[c]) < tol ? 0.0 : vals[c]; } newcols = ncols; ierr = MatRestoreRow(A, r, &ncols, &cols, &vals);CHKERRQ(ierr); ierr = MatSetValues(A, 1, &r, newcols, newCols, newVals, INSERT_VALUES);CHKERRQ(ierr); } ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); } ierr = PetscFree2(newCols,newVals);CHKERRQ(ierr); PetscFunctionReturn(0); }
static PetscErrorCode MatColoringApply_Greedy(MatColoring mc,ISColoring *iscoloring) { PetscErrorCode ierr; PetscInt finalcolor,finalcolor_global; ISColoringValue *colors; PetscInt ncolstotal,ncols; PetscReal *wts; PetscInt i,*lperm; PetscFunctionBegin; ierr = MatGetSize(mc->mat,NULL,&ncolstotal);CHKERRQ(ierr); ierr = MatGetLocalSize(mc->mat,NULL,&ncols);CHKERRQ(ierr); if (!mc->user_weights) { ierr = MatColoringCreateWeights(mc,&wts,&lperm);CHKERRQ(ierr); } else { wts = mc->user_weights; lperm = mc->user_lperm; } ierr = PetscMalloc1(ncols,&colors);CHKERRQ(ierr); if (mc->dist == 1) { ierr = GreedyColoringLocalDistanceOne_Private(mc,wts,lperm,colors);CHKERRQ(ierr); } else if (mc->dist == 2) { ierr = GreedyColoringLocalDistanceTwo_Private(mc,wts,lperm,colors);CHKERRQ(ierr); } else SETERRQ(PetscObjectComm((PetscObject)mc),PETSC_ERR_ARG_OUTOFRANGE,"Only distance 1 and distance 2 supported by MatColoringGreedy"); finalcolor=0; for (i=0;i<ncols;i++) { if (colors[i] > finalcolor) finalcolor=colors[i]; } finalcolor_global=0; ierr = MPIU_Allreduce(&finalcolor,&finalcolor_global,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)mc));CHKERRQ(ierr); ierr = PetscLogEventBegin(MATCOLORING_ISCreate,mc,0,0,0);CHKERRQ(ierr); ierr = ISColoringCreate(PetscObjectComm((PetscObject)mc),finalcolor_global+1,ncols,colors,PETSC_OWN_POINTER,iscoloring);CHKERRQ(ierr); ierr = PetscLogEventEnd(MATCOLORING_ISCreate,mc,0,0,0);CHKERRQ(ierr); if (!mc->user_weights) { ierr = PetscFree(wts);CHKERRQ(ierr); ierr = PetscFree(lperm);CHKERRQ(ierr); } PetscFunctionReturn(0); }
static PetscErrorCode MatStashScatterBegin_Ref(Mat mat,MatStash *stash,PetscInt *owners) { PetscInt *owner,*startv,*starti,tag1=stash->tag1,tag2=stash->tag2,bs2; PetscInt size=stash->size,nsends; PetscErrorCode ierr; PetscInt count,*sindices,**rindices,i,j,idx,lastidx,l; PetscScalar **rvalues,*svalues; MPI_Comm comm = stash->comm; MPI_Request *send_waits,*recv_waits,*recv_waits1,*recv_waits2; PetscMPIInt *sizes,*nlengths,nreceives; PetscInt *sp_idx,*sp_idy; PetscScalar *sp_val; PetscMatStashSpace space,space_next; PetscFunctionBegin; { /* make sure all processors are either in INSERTMODE or ADDMODE */ InsertMode addv; ierr = MPIU_Allreduce((PetscEnum*)&mat->insertmode,(PetscEnum*)&addv,1,MPIU_ENUM,MPI_BOR,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); if (addv == (ADD_VALUES|INSERT_VALUES)) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"Some processors inserted others added"); mat->insertmode = addv; /* in case this processor had no cache */ } bs2 = stash->bs*stash->bs; /* first count number of contributors to each processor */ ierr = PetscCalloc1(size,&sizes);CHKERRQ(ierr); ierr = PetscCalloc1(size,&nlengths);CHKERRQ(ierr); ierr = PetscMalloc1(stash->n+1,&owner);CHKERRQ(ierr); i = j = 0; lastidx = -1; space = stash->space_head; while (space) { space_next = space->next; sp_idx = space->idx; for (l=0; l<space->local_used; l++) { /* if indices are NOT locally sorted, need to start search at the beginning */ if (lastidx > (idx = sp_idx[l])) j = 0; lastidx = idx; for (; j<size; j++) { if (idx >= owners[j] && idx < owners[j+1]) { nlengths[j]++; owner[i] = j; break; } } i++; } space = space_next; } /* Now check what procs get messages - and compute nsends. */ for (i=0, nsends=0; i<size; i++) { if (nlengths[i]) { sizes[i] = 1; nsends++; } } {PetscMPIInt *onodes,*olengths; /* Determine the number of messages to expect, their lengths, from from-ids */ ierr = PetscGatherNumberOfMessages(comm,sizes,nlengths,&nreceives);CHKERRQ(ierr); ierr = PetscGatherMessageLengths(comm,nsends,nreceives,nlengths,&onodes,&olengths);CHKERRQ(ierr); /* since clubbing row,col - lengths are multiplied by 2 */ for (i=0; i<nreceives; i++) olengths[i] *=2; ierr = PetscPostIrecvInt(comm,tag1,nreceives,onodes,olengths,&rindices,&recv_waits1);CHKERRQ(ierr); /* values are size 'bs2' lengths (and remove earlier factor 2 */ for (i=0; i<nreceives; i++) olengths[i] = olengths[i]*bs2/2; ierr = PetscPostIrecvScalar(comm,tag2,nreceives,onodes,olengths,&rvalues,&recv_waits2);CHKERRQ(ierr); ierr = PetscFree(onodes);CHKERRQ(ierr); ierr = PetscFree(olengths);CHKERRQ(ierr);} /* do sends: 1) starts[i] gives the starting index in svalues for stuff going to the ith processor */ ierr = PetscMalloc2(bs2*stash->n,&svalues,2*(stash->n+1),&sindices);CHKERRQ(ierr); ierr = PetscMalloc1(2*nsends,&send_waits);CHKERRQ(ierr); ierr = PetscMalloc2(size,&startv,size,&starti);CHKERRQ(ierr); /* use 2 sends the first with all_a, the next with all_i and all_j */ startv[0] = 0; starti[0] = 0; for (i=1; i<size; i++) { startv[i] = startv[i-1] + nlengths[i-1]; starti[i] = starti[i-1] + 2*nlengths[i-1]; } i = 0; space = stash->space_head; while (space) { space_next = space->next; sp_idx = space->idx; sp_idy = space->idy; sp_val = space->val; for (l=0; l<space->local_used; l++) { j = owner[i]; if (bs2 == 1) { svalues[startv[j]] = sp_val[l]; } else { PetscInt k; PetscScalar *buf1,*buf2; buf1 = svalues+bs2*startv[j]; buf2 = space->val + bs2*l; for (k=0; k<bs2; k++) buf1[k] = buf2[k]; } sindices[starti[j]] = sp_idx[l]; sindices[starti[j]+nlengths[j]] = sp_idy[l]; startv[j]++; starti[j]++; i++; } space = space_next; } startv[0] = 0; for (i=1; i<size; i++) startv[i] = startv[i-1] + nlengths[i-1]; for (i=0,count=0; i<size; i++) { if (sizes[i]) { ierr = MPI_Isend(sindices+2*startv[i],2*nlengths[i],MPIU_INT,i,tag1,comm,send_waits+count++);CHKERRQ(ierr); ierr = MPI_Isend(svalues+bs2*startv[i],bs2*nlengths[i],MPIU_SCALAR,i,tag2,comm,send_waits+count++);CHKERRQ(ierr); } } #if defined(PETSC_USE_INFO) ierr = PetscInfo1(NULL,"No of messages: %d \n",nsends);CHKERRQ(ierr); for (i=0; i<size; i++) { if (sizes[i]) { ierr = PetscInfo2(NULL,"Mesg_to: %d: size: %d bytes\n",i,nlengths[i]*(bs2*sizeof(PetscScalar)+2*sizeof(PetscInt)));CHKERRQ(ierr); } } #endif ierr = PetscFree(nlengths);CHKERRQ(ierr); ierr = PetscFree(owner);CHKERRQ(ierr); ierr = PetscFree2(startv,starti);CHKERRQ(ierr); ierr = PetscFree(sizes);CHKERRQ(ierr); /* recv_waits need to be contiguous for MatStashScatterGetMesg_Private() */ ierr = PetscMalloc1(2*nreceives,&recv_waits);CHKERRQ(ierr); for (i=0; i<nreceives; i++) { recv_waits[2*i] = recv_waits1[i]; recv_waits[2*i+1] = recv_waits2[i]; } stash->recv_waits = recv_waits; ierr = PetscFree(recv_waits1);CHKERRQ(ierr); ierr = PetscFree(recv_waits2);CHKERRQ(ierr); stash->svalues = svalues; stash->sindices = sindices; stash->rvalues = rvalues; stash->rindices = rindices; stash->send_waits = send_waits; stash->nsends = nsends; stash->nrecvs = nreceives; stash->reproduce_count = 0; PetscFunctionReturn(0); }
PetscErrorCode TSEventHandler(TS ts) { PetscErrorCode ierr; TSEvent event; PetscReal t; Vec U; PetscInt i; PetscReal dt,dt_min; PetscInt rollback=0,in[2],out[2]; PetscInt fvalue_sign,fvalueprev_sign; PetscFunctionBegin; PetscValidHeaderSpecific(ts,TS_CLASSID,1); if (!ts->event) PetscFunctionReturn(0); event = ts->event; ierr = TSGetTime(ts,&t);CHKERRQ(ierr); ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); if (event->status == TSEVENT_NONE) { if (ts->steps == 1) /* After first successful step */ event->timestep_orig = ts->ptime - ts->ptime_prev; event->timestep_prev = dt; } if (event->status == TSEVENT_RESET_NEXTSTEP) { /* Restore time step */ dt = PetscMin(event->timestep_orig,event->timestep_prev); ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); event->status = TSEVENT_NONE; } if (event->status == TSEVENT_NONE) { event->ptime_end = t; } ierr = VecLockPush(U);CHKERRQ(ierr); ierr = (*event->eventhandler)(ts,t,U,event->fvalue,event->ctx);CHKERRQ(ierr); ierr = VecLockPop(U);CHKERRQ(ierr); for (i=0; i < event->nevents; i++) { if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { event->status = TSEVENT_ZERO; event->fvalue_right[i] = event->fvalue[i]; continue; } fvalue_sign = PetscSign(PetscRealPart(event->fvalue[i])); fvalueprev_sign = PetscSign(PetscRealPart(event->fvalue_prev[i])); if (fvalueprev_sign != 0 && (fvalue_sign != fvalueprev_sign) && (PetscAbsScalar(event->fvalue_prev[i]) > event->vtol[i])) { switch (event->direction[i]) { case -1: if (fvalue_sign < 0) { rollback = 1; /* Compute new time step */ dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); if (event->monitor) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); } event->fvalue_right[i] = event->fvalue[i]; event->side[i] = 1; if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; event->status = TSEVENT_LOCATED_INTERVAL; } break; case 1: if (fvalue_sign > 0) { rollback = 1; /* Compute new time step */ dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); if (event->monitor) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); } event->fvalue_right[i] = event->fvalue[i]; event->side[i] = 1; if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; event->status = TSEVENT_LOCATED_INTERVAL; } break; case 0: rollback = 1; /* Compute new time step */ dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); if (event->monitor) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Event %D interval detected [%g - %g]\n",event->iterctr,i,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); } event->fvalue_right[i] = event->fvalue[i]; event->side[i] = 1; if (!event->iterctr) event->zerocrossing[i] = PETSC_TRUE; event->status = TSEVENT_LOCATED_INTERVAL; break; } } } in[0] = event->status; in[1] = rollback; ierr = MPIU_Allreduce(in,out,2,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); event->status = (TSEventStatus)out[0]; rollback = out[1]; if (rollback) event->status = TSEVENT_LOCATED_INTERVAL; event->nevents_zero = 0; if (event->status == TSEVENT_ZERO) { for (i=0; i < event->nevents; i++) { if (PetscAbsScalar(event->fvalue[i]) < event->vtol[i]) { event->events_zero[event->nevents_zero++] = i; if (event->monitor) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: Event %D zero crossing at time %g located in %D iterations\n",i,(double)t,event->iterctr);CHKERRQ(ierr); } event->zerocrossing[i] = PETSC_FALSE; } event->side[i] = 0; } ierr = TSPostEvent(ts,t,U);CHKERRQ(ierr); dt = event->ptime_end - t; if (PetscAbsReal(dt) < PETSC_SMALL) { /* we hit the event, continue with the candidate time step */ dt = event->timestep_prev; event->status = TSEVENT_NONE; } ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr); event->iterctr = 0; PetscFunctionReturn(0); } if (event->status == TSEVENT_LOCATED_INTERVAL) { ierr = TSRollBack(ts);CHKERRQ(ierr); ierr = TSSetConvergedReason(ts,TS_CONVERGED_ITERATING);CHKERRQ(ierr); event->status = TSEVENT_PROCESSING; event->ptime_right = t; } else { for (i=0; i < event->nevents; i++) { if (event->status == TSEVENT_PROCESSING && event->zerocrossing[i]) { /* Compute new time step */ dt = TSEventComputeStepSize(event->ptime_prev,t,event->ptime_right,event->fvalue_prev[i],event->fvalue[i],event->fvalue_right[i],event->side[i],dt); event->side[i] = -1; } event->fvalue_prev[i] = event->fvalue[i]; } if (event->monitor && event->status == TSEVENT_PROCESSING) { ierr = PetscViewerASCIIPrintf(event->monitor,"TSEvent: iter %D - Stepping forward as no event detected in interval [%g - %g]\n",event->iterctr,(double)event->ptime_prev,(double)t);CHKERRQ(ierr); } event->ptime_prev = t; } if (event->status == TSEVENT_PROCESSING) event->iterctr++; ierr = MPIU_Allreduce(&dt,&dt_min,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); ierr = TSSetTimeStep(ts,dt_min);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode MatISGetMPIXAIJ_IS(Mat mat, MatReuse reuse, Mat *M) { Mat_IS *matis = (Mat_IS*)(mat->data); Mat local_mat; /* info on mat */ PetscInt bs,rows,cols,lrows,lcols; PetscInt local_rows,local_cols; PetscBool isdense,issbaij,isseqaij; PetscMPIInt nsubdomains; /* values insertion */ PetscScalar *array; /* work */ PetscErrorCode ierr; PetscFunctionBegin; /* get info from mat */ ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&nsubdomains);CHKERRQ(ierr); if (nsubdomains == 1) { if (reuse == MAT_INITIAL_MATRIX) { ierr = MatDuplicate(matis->A,MAT_COPY_VALUES,&(*M));CHKERRQ(ierr); } else { ierr = MatCopy(matis->A,*M,SAME_NONZERO_PATTERN);CHKERRQ(ierr); } PetscFunctionReturn(0); } ierr = MatGetSize(mat,&rows,&cols);CHKERRQ(ierr); ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); ierr = MatGetLocalSize(mat,&lrows,&lcols);CHKERRQ(ierr); ierr = MatGetSize(matis->A,&local_rows,&local_cols);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQDENSE,&isdense);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)matis->A,MATSEQSBAIJ,&issbaij);CHKERRQ(ierr); if (reuse == MAT_INITIAL_MATRIX) { MatType new_mat_type; PetscBool issbaij_red; /* determining new matrix type */ ierr = MPIU_Allreduce(&issbaij,&issbaij_red,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); if (issbaij_red) { new_mat_type = MATSBAIJ; } else { if (bs>1) { new_mat_type = MATBAIJ; } else { new_mat_type = MATAIJ; } } ierr = MatCreate(PetscObjectComm((PetscObject)mat),M);CHKERRQ(ierr); ierr = MatSetSizes(*M,lrows,lcols,rows,cols);CHKERRQ(ierr); ierr = MatSetBlockSize(*M,bs);CHKERRQ(ierr); ierr = MatSetType(*M,new_mat_type);CHKERRQ(ierr); ierr = MatISSetMPIXAIJPreallocation_Private(mat,*M,PETSC_FALSE);CHKERRQ(ierr); } else { PetscInt mbs,mrows,mcols,mlrows,mlcols; /* some checks */ ierr = MatGetBlockSize(*M,&mbs);CHKERRQ(ierr); ierr = MatGetSize(*M,&mrows,&mcols);CHKERRQ(ierr); ierr = MatGetLocalSize(*M,&mlrows,&mlcols);CHKERRQ(ierr); if (mrows != rows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of rows (%d != %d)",rows,mrows); if (mcols != cols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of cols (%d != %d)",cols,mcols); if (mlrows != lrows) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local rows (%d != %d)",lrows,mlrows); if (mlcols != lcols) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong number of local cols (%d != %d)",lcols,mlcols); if (mbs != bs) SETERRQ2(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot reuse matrix. Wrong block size (%d != %d)",bs,mbs); ierr = MatZeroEntries(*M);CHKERRQ(ierr); } if (issbaij) { ierr = MatConvert(matis->A,MATSEQBAIJ,MAT_INITIAL_MATRIX,&local_mat);CHKERRQ(ierr); } else { ierr = PetscObjectReference((PetscObject)matis->A);CHKERRQ(ierr); local_mat = matis->A; } /* Set values */ ierr = MatSetLocalToGlobalMapping(*M,mat->rmap->mapping,mat->cmap->mapping);CHKERRQ(ierr); if (isdense) { /* special case for dense local matrices */ PetscInt i,*dummy_rows,*dummy_cols; if (local_rows != local_cols) { ierr = PetscMalloc2(local_rows,&dummy_rows,local_cols,&dummy_cols);CHKERRQ(ierr); for (i=0;i<local_rows;i++) dummy_rows[i] = i; for (i=0;i<local_cols;i++) dummy_cols[i] = i; } else { ierr = PetscMalloc1(local_rows,&dummy_rows);CHKERRQ(ierr); for (i=0;i<local_rows;i++) dummy_rows[i] = i; dummy_cols = dummy_rows; } ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); ierr = MatDenseGetArray(local_mat,&array);CHKERRQ(ierr); ierr = MatSetValuesLocal(*M,local_rows,dummy_rows,local_cols,dummy_cols,array,ADD_VALUES);CHKERRQ(ierr); ierr = MatDenseRestoreArray(local_mat,&array);CHKERRQ(ierr); if (dummy_rows != dummy_cols) { ierr = PetscFree2(dummy_rows,dummy_cols);CHKERRQ(ierr); } else { ierr = PetscFree(dummy_rows);CHKERRQ(ierr); } } else if (isseqaij) { PetscInt i,nvtxs,*xadj,*adjncy; PetscBool done; ierr = MatGetRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); ierr = MatSeqAIJGetArray(local_mat,&array);CHKERRQ(ierr); for (i=0;i<nvtxs;i++) { ierr = MatSetValuesLocal(*M,1,&i,xadj[i+1]-xadj[i],adjncy+xadj[i],array+xadj[i],ADD_VALUES);CHKERRQ(ierr); } ierr = MatRestoreRowIJ(local_mat,0,PETSC_FALSE,PETSC_FALSE,&nvtxs,(const PetscInt**)&xadj,(const PetscInt**)&adjncy,&done);CHKERRQ(ierr); if (!done) SETERRQ1(PetscObjectComm((PetscObject)local_mat),PETSC_ERR_PLIB,"Error in MatRestoreRowIJ called in %s\n",__FUNCT__); ierr = MatSeqAIJRestoreArray(local_mat,&array);CHKERRQ(ierr); } else { /* very basic values insertion for all other matrix types */ PetscInt i; for (i=0;i<local_rows;i++) { PetscInt j; const PetscInt *local_indices_cols; ierr = MatGetRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); ierr = MatSetValuesLocal(*M,1,&i,j,local_indices_cols,array,ADD_VALUES);CHKERRQ(ierr); ierr = MatRestoreRow(local_mat,i,&j,&local_indices_cols,(const PetscScalar**)&array);CHKERRQ(ierr); } } ierr = MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatDestroy(&local_mat);CHKERRQ(ierr); ierr = MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); if (isdense) { ierr = MatSetOption(*M,MAT_ROW_ORIENTED,PETSC_TRUE);CHKERRQ(ierr); } PetscFunctionReturn(0); }
/*@ ISPairToList - convert an IS pair encoding an integer map to a list of ISs. Each IS on the output list contains the preimage for each index on the second input IS. The ISs on the output list are constructed on the subcommunicators of the input IS pair. Each subcommunicator corresponds to the preimage of some index j -- this subcomm contains exactly the ranks that assign some indices i to j. This is essentially the inverse of ISListToPair(). Collective on indis. Input arguments: + xis - domain IS - yis - range IS Output arguments: + listlen - length of islist - islist - list of ISs breaking up indis by color Note: + xis and yis must be of the same length and have congruent communicators. - The resulting ISs have subcommunicators in a "deadlock-free" order (see ISListToPair()). Level: advanced .seealso ISListToPair() @*/ PetscErrorCode ISPairToList(IS xis, IS yis, PetscInt *listlen, IS **islist) { PetscErrorCode ierr; IS indis = xis, coloris = yis; PetscInt *inds, *colors, llen, ilen, lstart, lend, lcount,l; PetscMPIInt rank, size, llow, lhigh, low, high,color,subsize; const PetscInt *ccolors, *cinds; MPI_Comm comm, subcomm; PetscFunctionBegin; PetscValidHeaderSpecific(xis, IS_CLASSID, 1); PetscValidHeaderSpecific(yis, IS_CLASSID, 2); PetscCheckSameComm(xis,1,yis,2); PetscValidIntPointer(listlen,3); PetscValidPointer(islist,4); ierr = PetscObjectGetComm((PetscObject)xis,&comm);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); ierr = MPI_Comm_rank(comm, &size);CHKERRQ(ierr); /* Extract, copy and sort the local indices and colors on the color. */ ierr = ISGetLocalSize(coloris, &llen);CHKERRQ(ierr); ierr = ISGetLocalSize(indis, &ilen);CHKERRQ(ierr); if (llen != ilen) SETERRQ2(comm, PETSC_ERR_ARG_SIZ, "Incompatible IS sizes: %D and %D", ilen, llen); ierr = ISGetIndices(coloris, &ccolors);CHKERRQ(ierr); ierr = ISGetIndices(indis, &cinds);CHKERRQ(ierr); ierr = PetscMalloc2(ilen,&inds,llen,&colors);CHKERRQ(ierr); ierr = PetscMemcpy(inds,cinds,ilen*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscMemcpy(colors,ccolors,llen*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscSortIntWithArray(llen, colors, inds);CHKERRQ(ierr); /* Determine the global extent of colors. */ llow = 0; lhigh = -1; lstart = 0; lcount = 0; while (lstart < llen) { lend = lstart+1; while (lend < llen && colors[lend] == colors[lstart]) ++lend; llow = PetscMin(llow,colors[lstart]); lhigh = PetscMax(lhigh,colors[lstart]); ++lcount; } ierr = MPIU_Allreduce(&llow,&low,1,MPI_INT,MPI_MIN,comm);CHKERRQ(ierr); ierr = MPIU_Allreduce(&lhigh,&high,1,MPI_INT,MPI_MAX,comm);CHKERRQ(ierr); *listlen = 0; if (low <= high) { if (lcount > 0) { *listlen = lcount; if (!*islist) { ierr = PetscMalloc1(lcount, islist);CHKERRQ(ierr); } } /* Traverse all possible global colors, and participate in the subcommunicators for the locally-supported colors. */ lcount = 0; lstart = 0; lend = 0; for (l = low; l <= high; ++l) { /* Find the range of indices with the same color, which is not smaller than l. Observe that, since colors is sorted, and is a subsequence of [low,high], as soon as we find a new color, it is >= l. */ if (lstart < llen) { /* The start of the next locally-owned color is identified. Now look for the end. */ if (lstart == lend) { lend = lstart+1; while (lend < llen && colors[lend] == colors[lstart]) ++lend; } /* Now check whether the identified color segment matches l. */ if (colors[lstart] < l) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Locally owned color %D at location %D is < than the next global color %D", colors[lstart], lcount, l); } color = (PetscMPIInt)(colors[lstart] == l); /* Check whether a proper subcommunicator exists. */ ierr = MPIU_Allreduce(&color,&subsize,1,MPI_INT,MPI_SUM,comm);CHKERRQ(ierr); if (subsize == 1) subcomm = PETSC_COMM_SELF; else if (subsize == size) subcomm = comm; else { /* a proper communicator is necessary, so we create it. */ ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); } if (colors[lstart] == l) { /* If we have l among the local colors, we create an IS to hold the corresponding indices. */ ierr = ISCreateGeneral(subcomm, lend-lstart,inds+lstart,PETSC_COPY_VALUES,*islist+lcount);CHKERRQ(ierr); /* Position lstart at the beginning of the next local color. */ lstart = lend; /* Increment the counter of the local colors split off into an IS. */ ++lcount; } if (subsize > 0 && subsize < size) { /* Irrespective of color, destroy the split off subcomm: a subcomm used in the IS creation above is duplicated into a proper PETSc comm. */ ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); } } /* for (l = low; l < high; ++l) */ } /* if (low <= high) */ ierr = PetscFree2(inds,colors);CHKERRQ(ierr); PetscFunctionReturn(0); }