#include <../src/mat/impls/aij/mpi/mpiaij.h> #include <../src/mat/impls/baij/mpi/mpibaij.h> #undef __FUNCT__ #define __FUNCT__ "MatFDColoringApply_BAIJ" PetscErrorCode MatFDColoringApply_BAIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) { PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; PetscErrorCode ierr; PetscInt k,cstart,cend,l,row,col,nz,spidx,i,j; PetscScalar dx=0.0,*y,*xx,*w3_array,*dy_i,*dy=coloring->dy; PetscScalar *vscale_array; PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; void *fctx=coloring->fctx; PetscBool flg=PETSC_FALSE; PetscInt ctype=coloring->ctype,nxloc,nrows_k; PetscScalar *valaddr; MatEntry *Jentry=coloring->matentry; const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; PetscInt bs=J->rmap->bs; PetscFunctionBegin; ierr = MatSetUnfactored(J);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); if (flg) { ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); } else { PetscBool assembled; ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); if (assembled) { ierr = MatZeroEntries(J);CHKERRQ(ierr); } } /* create vscale for storing dx */ if (!vscale) { if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { /* contain the "diagonal" on processor scalings followed by the off processor* - garray must be non-bloced! */ Mat_MPIBAIJ *baij=(Mat_MPIBAIJ*)J->data; PetscInt *garray; ierr = PetscMalloc(baij->B->cmap->n*sizeof(PetscInt),&garray);CHKERRQ(ierr); for (i=0; iB->cmap->n/bs; i++) { for (j=0; jgarray[i]+j; } } ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,baij->B->cmap->n,garray,&vscale);CHKERRQ(ierr); ierr = PetscFree(garray);CHKERRQ(ierr); } else if (ctype == IS_COLORING_GHOSTED) { ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); } coloring->vscale = vscale; } /* (1) Set w1 = F(x1) */ if (!coloring->fset) { ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); } else { coloring->fset = PETSC_FALSE; } /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); if (coloring->htype[0] == 'w') { /* vscale = dx is a constant scalar */ ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); } else { ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); for (col=0; col= 0.0) dx = umin; else if (PetscRealPart(dx) < 0.0 ) dx = -umin; } dx *= epsilon; vscale_array[col] = 1.0/dx; } ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); } if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] != 'w') { ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); } /* (3) Loop over each color */ if (!coloring->w3) { ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); } w3 = coloring->w3; ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ if (vscale) { ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); } nz = 0; for (k=0; kcurrentcolor = k; /* (3-1) Loop over each column associated with color adding the perturbation to the vector w3 = x1 + dx. */ ierr = VecCopy(x1,w3);CHKERRQ(ierr); dy_i = dy; for (i=0; ihtype[0] == 'w') { for (l=0; lcolumns[k][l]; /* local column (in global index!) of the matrix we are probing for */ w3_array[col] += 1.0/dx; if (i) { w3_array[col-1] -= 1.0/dx; /* resume original w3[col-1] */ } } } else { /* htype == 'ds' */ vscale_array -= cstart; /* shift pointer so global index can be used */ for (l=0; lcolumns[k][l]; /* local column (in global index!) of the matrix we are probing for */ w3_array[col] += 1.0/vscale_array[col]; if (i) { w3_array[col-1] -= 1.0/vscale_array[col-1]; /* resume original w3[col-1] */ } } vscale_array += cstart; } if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); /* (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) w2 = F(x1 + dx) - F(x1) */ ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = VecPlaceArray(w2,dy_i);CHKERRQ(ierr); /* place w2 to the array dy_i */ ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); //--------------------------------------------- ierr = VecResetArray(w2);CHKERRQ(ierr); dy_i += nxloc; /* points to dy+i*nxloc */ } /* (3-3) Loop over rows of vector, putting results into Jacobian matrix */ nrows_k = nrows[k]; ierr = VecGetArray(w2,&y);CHKERRQ(ierr); if (coloring->htype[0] == 'w') { for (l=0; lcurrentcolor = -1; PetscFunctionReturn(0); } extern PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat); #undef __FUNCT__ #define __FUNCT__ "MatFDColoringApply_AIJ" PetscErrorCode MatFDColoringApply_AIJ(Mat J,MatFDColoring coloring,Vec x1,MatStructure *flag,void *sctx) { PetscErrorCode (*f)(void*,Vec,Vec,void*) = (PetscErrorCode (*)(void*,Vec,Vec,void*))coloring->f; PetscErrorCode ierr; PetscInt k,cstart,cend,l,row,col,nz; PetscScalar dx=0.0,*y,*xx,*w3_array; PetscScalar *vscale_array; PetscReal epsilon=coloring->error_rel,umin=coloring->umin,unorm; Vec w1=coloring->w1,w2=coloring->w2,w3,vscale=coloring->vscale; void *fctx=coloring->fctx; PetscBool flg=PETSC_FALSE; PetscInt ctype=coloring->ctype,nxloc,nrows_k; MatEntry *Jentry=coloring->matentry; const PetscInt ncolors=coloring->ncolors,*ncolumns=coloring->ncolumns,*nrows=coloring->nrows; PetscFunctionBegin; ierr = MatSetUnfactored(J);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,"-mat_fd_coloring_dont_rezero",&flg,NULL);CHKERRQ(ierr); if (flg) { ierr = PetscInfo(coloring,"Not calling MatZeroEntries()\n");CHKERRQ(ierr); } else { PetscBool assembled; ierr = MatAssembled(J,&assembled);CHKERRQ(ierr); if (assembled) { ierr = MatZeroEntries(J);CHKERRQ(ierr); } } /* create vscale for storing dx */ if (!vscale) { if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] == 'd') { Mat_MPIAIJ *aij=(Mat_MPIAIJ*)J->data; ierr = VecCreateGhost(PetscObjectComm((PetscObject)J),J->cmap->n,PETSC_DETERMINE,aij->B->cmap->n,aij->garray,&vscale);CHKERRQ(ierr); } else if (ctype == IS_COLORING_GHOSTED) { ierr = VecDuplicate(x1,&vscale);CHKERRQ(ierr); } coloring->vscale = vscale; } /* (1) Set w1 = F(x1) */ if (!coloring->fset) { ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = (*f)(sctx,x1,w1,fctx);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); } else { coloring->fset = PETSC_FALSE; } /* (2) Compute vscale = 1./dx - the local scale factors, including ghost points */ if (coloring->htype[0] == 'w') { /* vscale = dx is a constant scalar */ ierr = VecNorm(x1,NORM_2,&unorm);CHKERRQ(ierr); dx = 1.0/(PetscSqrtReal(1.0 + unorm)*epsilon); } else { ierr = VecGetLocalSize(x1,&nxloc);CHKERRQ(ierr); ierr = VecGetArray(x1,&xx);CHKERRQ(ierr); ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); for (col=0; col= 0.0) dx = umin; else if (PetscRealPart(dx) < 0.0 ) dx = -umin; } dx *= epsilon; vscale_array[col] = 1.0/dx; } ierr = VecRestoreArray(x1,&xx);CHKERRQ(ierr); ierr = VecRestoreArray(vscale,&vscale_array);CHKERRQ(ierr); } if (ctype == IS_COLORING_GLOBAL && coloring->htype[0] != 'w') { ierr = VecGhostUpdateBegin(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecGhostUpdateEnd(vscale,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); } /* (3) Loop over each color */ if (!coloring->w3) { ierr = VecDuplicate(x1,&coloring->w3);CHKERRQ(ierr); ierr = PetscLogObjectParent((PetscObject)coloring,(PetscObject)coloring->w3);CHKERRQ(ierr); } w3 = coloring->w3; ierr = VecGetOwnershipRange(x1,&cstart,&cend);CHKERRQ(ierr); /* used by ghosted vscale */ if (vscale) { ierr = VecGetArray(vscale,&vscale_array);CHKERRQ(ierr); } nz = 0; for (k=0; kcurrentcolor = k; /* (3-1) Loop over each column associated with color adding the perturbation to the vector w3 = x1 + dx. */ ierr = VecCopy(x1,w3);CHKERRQ(ierr); ierr = VecGetArray(w3,&w3_array);CHKERRQ(ierr); if (ctype == IS_COLORING_GLOBAL) w3_array -= cstart; /* shift pointer so global index can be used */ if (coloring->htype[0] == 'w') { for (l=0; lcolumns[k][l]; /* local column (in global index!) of the matrix we are probing for */ w3_array[col] += 1.0/dx; } } else { /* htype == 'ds' */ vscale_array -= cstart; /* shift pointer so global index can be used */ for (l=0; lcolumns[k][l]; /* local column (in global index!) of the matrix we are probing for */ w3_array[col] += 1.0/vscale_array[col]; } vscale_array += cstart; } if (ctype == IS_COLORING_GLOBAL) w3_array += cstart; ierr = VecRestoreArray(w3,&w3_array);CHKERRQ(ierr); /* (3-2) Evaluate function at w3 = x1 + dx (here dx is a vector of perturbations) w2 = F(x1 + dx) - F(x1) */ ierr = PetscLogEventBegin(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = (*f)(sctx,w3,w2,fctx);CHKERRQ(ierr); ierr = PetscLogEventEnd(MAT_FDColoringFunction,0,0,0,0);CHKERRQ(ierr); ierr = VecAXPY(w2,-1.0,w1);CHKERRQ(ierr); /* (3-3) Loop over rows of vector, putting results into Jacobian matrix */ nrows_k = nrows[k]; ierr = VecGetArray(w2,&y);CHKERRQ(ierr); if (coloring->htype[0] == 'w') { for (l=0; lcurrentcolor = -1; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatFDColoringCreate_MPIXAIJ" PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat mat,ISColoring iscoloring,MatFDColoring c) { PetscErrorCode ierr; PetscMPIInt size,*ncolsonproc,*disp,nn; PetscInt i,n,nrows,nrows_i,j,k,m,ncols,col; const PetscInt *is,*A_ci,*A_cj,*B_ci,*B_cj,*row = NULL,*ltog=NULL; PetscInt nis=iscoloring->n,nctot,*cols; PetscInt *rowhit,cstart,cend,colb; IS *isa; ISLocalToGlobalMapping map=mat->cmap->mapping; PetscInt ctype=c->ctype; Mat A,B; PetscScalar *A_val,*B_val; PetscInt spidx; PetscInt *spidxA,*spidxB,nz,bs,bs2; PetscScalar **valaddrhit; MatEntry *Jentry; PetscBool isBAIJ; #if defined(PETSC_USE_CTABLE) PetscTable colmap=NULL; #else PetscInt *colmap=NULL; /* local col number of off-diag col */ #endif PetscFunctionBegin; if (ctype == IS_COLORING_GHOSTED) { if (!map) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"When using ghosted differencing matrix must have local to global mapping provided with MatSetLocalToGlobalMapping"); ierr = ISLocalToGlobalMappingGetIndices(map,<og);CHKERRQ(ierr); } ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr); if (!isBAIJ) { Mat_MPIAIJ *aij=(Mat_MPIAIJ*)mat->data; Mat_SeqAIJ *spA,*spB; bs = 1; /* only bs=1 is supported for non MPIBAIJ matrix */ A=aij->A; spA=(Mat_SeqAIJ*)A->data; A_val=spA->a; B=aij->B; spB=(Mat_SeqAIJ*)B->data; B_val=spB->a; nz = spA->nz + spB->nz; /* total nonzero entries of mat */ if (!aij->colmap) { /* Allow access to data structures of local part of matrix - creates aij->colmap which maps global column number to local number in part B */ ierr = MatCreateColmap_MPIAIJ_Private(mat);CHKERRQ(ierr); colmap = aij->colmap; } ierr = MatGetColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); ierr = MatGetColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); } else { Mat_MPIBAIJ *aij=(Mat_MPIBAIJ*)mat->data; Mat_SeqBAIJ *spA,*spB; A=aij->A; spA=(Mat_SeqBAIJ*)A->data; A_val=spA->a; B=aij->B; spB=(Mat_SeqBAIJ*)B->data; B_val=spB->a; nz = spA->nz + spB->nz; /* total nonzero entries of mat */ if (!aij->colmap) { ierr = MatCreateColmap_MPIBAIJ_Private(mat);CHKERRQ(ierr); colmap = aij->colmap; } ierr = MatGetColumnIJ_SeqBAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); ierr = MatGetColumnIJ_SeqBAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); } m = mat->rmap->n/bs; cstart = mat->cmap->rstart/bs; cend = mat->cmap->rend/bs; c->M = mat->rmap->N/bs; /* set the global rows and columns and local rows */ c->N = mat->cmap->N/bs; c->m = m; c->rstart = mat->rmap->rstart/bs; c->ncolors = nis; ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); ierr = PetscMalloc(nis*sizeof(PetscInt*),&c->columns);CHKERRQ(ierr); ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); ierr = PetscLogObjectMemory((PetscObject)c,3*nis*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscMalloc(nz*sizeof(MatEntry),&Jentry);CHKERRQ(ierr); ierr = PetscLogObjectMemory((PetscObject)c,nz*sizeof(MatEntry));CHKERRQ(ierr); c->matentry = Jentry; ierr = PetscMalloc2(m+1,PetscInt,&rowhit,m+1,PetscScalar*,&valaddrhit);CHKERRQ(ierr); nz = 0; ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); for (i=0; incolumns[i] = n; /* local number of columns of this color on this process */ if (n) { ierr = PetscMalloc(n*sizeof(PetscInt),&c->columns[i]);CHKERRQ(ierr); ierr = PetscLogObjectMemory((PetscObject)c,n*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscMemcpy(c->columns[i],is,n*sizeof(PetscInt));CHKERRQ(ierr); /* convert global column indices to local ones! */ } else { c->columns[i] = 0; } if (ctype == IS_COLORING_GLOBAL) { /* Determine nctot, the total (parallel) number of columns of this color */ ierr = MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);CHKERRQ(ierr); ierr = PetscMalloc2(size,PetscMPIInt,&ncolsonproc,size,PetscMPIInt,&disp);CHKERRQ(ierr); /* ncolsonproc[j]: local ncolumns on proc[j] of this color */ ierr = PetscMPIIntCast(n,&nn);CHKERRQ(ierr); ierr = MPI_Allgather(&nn,1,MPI_INT,ncolsonproc,1,MPI_INT,PetscObjectComm((PetscObject)mat));CHKERRQ(ierr); nctot = 0; for (j=0; j= cstart && col < cend) { /* column is in diagonal block of matrix A */ row = A_cj + A_ci[col-cstart]; nrows = A_ci[col-cstart+1] - A_ci[col-cstart]; nrows_i += nrows; /* loop over columns of A marking them in rowhit */ for (k=0; knrows[i] = nrows_i; for (j=0; jcmap->n*sizeof(PetscScalar),&c->dy);CHKERRQ(ierr); } else { ierr = MatRestoreColumnIJ_SeqAIJ_Color(A,0,PETSC_FALSE,PETSC_FALSE,&ncols,&A_ci,&A_cj,&spidxA,NULL);CHKERRQ(ierr); ierr = MatRestoreColumnIJ_SeqAIJ_Color(B,0,PETSC_FALSE,PETSC_FALSE,&ncols,&B_ci,&B_cj,&spidxB,NULL);CHKERRQ(ierr); } if (ctype == IS_COLORING_GHOSTED) { ierr = ISLocalToGlobalMappingRestoreIndices(map,<og);CHKERRQ(ierr); } PetscFunctionReturn(0); }