1 #define PETSCKSP_DLL 2 3 4 #include "petscpc.h" /*I "petscpc.h" I*/ 5 #include "petscmg.h" /*I "petscmg.h" I*/ 6 #include "petscdm.h" /*I "petscdm.h" I*/ 7 #include "../src/ksp/pc/impls/mg/mgimpl.h" 8 9 typedef struct { 10 PCExoticType type; 11 Mat P; /* the constructed interpolation matrix */ 12 PetscBool directSolve; /* use direct LU factorization to construct interpolation */ 13 KSP ksp; 14 } PC_Exotic; 15 16 const char *PCExoticTypes[] = {"face","wirebasket","PCExoticType","PC_Exotic",0}; 17 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "DMDAGetWireBasketInterpolation" 21 /* 22 DMDAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space 23 24 */ 25 PetscErrorCode DMDAGetWireBasketInterpolation(DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P) 26 { 27 PetscErrorCode ierr; 28 PetscInt dim,i,j,k,m,n,p,dof,Nint,Nface,Nwire,Nsurf,*Iint,*Isurf,cint = 0,csurf = 0,istart,jstart,kstart,*II,N,c = 0; 29 PetscInt mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf; 30 Mat Xint, Xsurf,Xint_tmp; 31 IS isint,issurf,is,row,col; 32 ISLocalToGlobalMapping ltg; 33 MPI_Comm comm; 34 Mat A,Aii,Ais,Asi,*Aholder,iAii; 35 MatFactorInfo info; 36 PetscScalar *xsurf,*xint; 37 #if defined(PETSC_USE_DEBUG_foo) 38 PetscScalar tmp; 39 #endif 40 PetscTable ht; 41 42 PetscFunctionBegin; 43 ierr = DMDAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);CHKERRQ(ierr); 44 if (dof != 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only for single field problems"); 45 if (dim != 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only coded for 3d problems"); 46 ierr = DMDAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr); 47 ierr = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr); 48 istart = istart ? -1 : 0; 49 jstart = jstart ? -1 : 0; 50 kstart = kstart ? -1 : 0; 51 52 /* 53 the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 54 to all the local degrees of freedom (this includes the vertices, edges and faces). 55 56 Xint are the subset of the interpolation into the interior 57 58 Xface are the interpolation onto faces but not into the interior 59 60 Xsurf are the interpolation onto the vertices and edges (the surfbasket) 61 Xint 62 Symbolically one could write P = ( Xface ) after interchanging the rows to match the natural ordering on the domain 63 Xsurf 64 */ 65 N = (m - istart)*(n - jstart)*(p - kstart); 66 Nint = (m-2-istart)*(n-2-jstart)*(p-2-kstart); 67 Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) ); 68 Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8; 69 Nsurf = Nface + Nwire; 70 ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,26,PETSC_NULL,&Xint);CHKERRQ(ierr); 71 ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,PETSC_NULL,&Xsurf);CHKERRQ(ierr); 72 ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 73 74 /* 75 Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 76 Xsurf will be all zero (thus making the coarse matrix singular). 77 */ 78 if (m-istart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3"); 79 if (n-jstart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3"); 80 if (p-kstart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3"); 81 82 cnt = 0; 83 xsurf[cnt++] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + Nsurf] = 1;} xsurf[cnt++ + 2*Nsurf] = 1; 84 for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 3*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;} xsurf[cnt++ + 5*Nsurf] = 1;} 85 xsurf[cnt++ + 6*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 7*Nsurf] = 1;} xsurf[cnt++ + 8*Nsurf] = 1; 86 for (k=1;k<p-1-kstart;k++) { 87 xsurf[cnt++ + 9*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 10*Nsurf] = 1;} xsurf[cnt++ + 11*Nsurf] = 1; 88 for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 12*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 13*Nsurf] = 1;} 89 xsurf[cnt++ + 14*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 15*Nsurf] = 1;} xsurf[cnt++ + 16*Nsurf] = 1; 90 } 91 xsurf[cnt++ + 17*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 18*Nsurf] = 1;} xsurf[cnt++ + 19*Nsurf] = 1; 92 for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 20*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 21*Nsurf] = 1;} xsurf[cnt++ + 22*Nsurf] = 1;} 93 xsurf[cnt++ + 23*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 24*Nsurf] = 1;} xsurf[cnt++ + 25*Nsurf] = 1; 94 95 /* interpolations only sum to 1 when using direct solver */ 96 #if defined(PETSC_USE_DEBUG_foo) 97 for (i=0; i<Nsurf; i++) { 98 tmp = 0.0; 99 for (j=0; j<26; j++) { 100 tmp += xsurf[i+j*Nsurf]; 101 } 102 if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %G",i,PetscAbsScalar(tmp)); 103 } 104 #endif 105 ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 106 /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/ 107 108 109 /* 110 I are the indices for all the needed vertices (in global numbering) 111 Iint are the indices for the interior values, I surf for the surface values 112 (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it 113 is NOT the local DMDA ordering.) 114 IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering 115 */ 116 #define Endpoint(a,start,b) (a == 0 || a == (b-1-start)) 117 ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr); 118 ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr); 119 for (k=0; k<p-kstart; k++) { 120 for (j=0; j<n-jstart; j++) { 121 for (i=0; i<m-istart; i++) { 122 II[c++] = i + j*mwidth + k*mwidth*nwidth; 123 124 if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) { 125 IIint[cint] = i + j*mwidth + k*mwidth*nwidth; 126 Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 127 } else { 128 IIsurf[csurf] = i + j*mwidth + k*mwidth*nwidth; 129 Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 130 } 131 } 132 } 133 } 134 if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N"); 135 if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint"); 136 if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf"); 137 ierr = DMGetLocalToGlobalMapping(da,<g);CHKERRQ(ierr); 138 ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr); 139 ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr); 140 ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr); 141 ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 142 ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 143 ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr); 144 ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr); 145 ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr); 146 147 ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr); 148 A = *Aholder; 149 ierr = PetscFree(Aholder);CHKERRQ(ierr); 150 151 ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr); 152 ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr); 153 ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr); 154 155 /* 156 Solve for the interpolation onto the interior Xint 157 */ 158 ierr = MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);CHKERRQ(ierr); 159 ierr = MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr); 160 ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr); 161 if (exotic->directSolve) { 162 ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr); 163 ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 164 ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr); 165 ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr); 166 ierr = ISDestroy(row);CHKERRQ(ierr); 167 ierr = ISDestroy(col);CHKERRQ(ierr); 168 ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr); 169 ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr); 170 ierr = MatDestroy(iAii);CHKERRQ(ierr); 171 } else { 172 Vec b,x; 173 PetscScalar *xint_tmp; 174 175 ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 176 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);CHKERRQ(ierr); 177 ierr = MatGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 178 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);CHKERRQ(ierr); 179 ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 180 for (i=0; i<26; i++) { 181 ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr); 182 ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr); 183 ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr); 184 ierr = VecResetArray(x);CHKERRQ(ierr); 185 ierr = VecResetArray(b);CHKERRQ(ierr); 186 } 187 ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 188 ierr = MatRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 189 ierr = VecDestroy(x);CHKERRQ(ierr); 190 ierr = VecDestroy(b);CHKERRQ(ierr); 191 } 192 ierr = MatDestroy(Xint_tmp);CHKERRQ(ierr); 193 194 #if defined(PETSC_USE_DEBUG_foo) 195 ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 196 for (i=0; i<Nint; i++) { 197 tmp = 0.0; 198 for (j=0; j<26; j++) { 199 tmp += xint[i+j*Nint]; 200 } 201 if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %G",i,PetscAbsScalar(tmp)); 202 } 203 ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 204 /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 205 #endif 206 207 208 /* total vertices total faces total edges */ 209 Ntotal = (mp + 1)*(np + 1)*(pp + 1) + mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1) + mp*(np+1)*(pp+1) + np*(mp+1)*(pp+1) + pp*(mp+1)*(np+1); 210 211 /* 212 For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 213 */ 214 cnt = 0; 215 gl[cnt++] = 0; { gl[cnt++] = 1;} gl[cnt++] = m-istart-1; 216 { gl[cnt++] = mwidth; { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;} 217 gl[cnt++] = mwidth*(n-jstart-1); { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1; 218 { 219 gl[cnt++] = mwidth*nwidth; { gl[cnt++] = mwidth*nwidth+1;} gl[cnt++] = mwidth*nwidth+ m-istart-1; 220 { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;} 221 gl[cnt++] = mwidth*nwidth+ mwidth*(n-jstart-1); { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1) + m-istart-1; 222 } 223 gl[cnt++] = mwidth*nwidth*(p-kstart-1); { gl[cnt++] = mwidth*nwidth*(p-kstart-1)+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1) + m-istart-1; 224 { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth; { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1)+mwidth+m-istart-1;} 225 gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth*(n-jstart-1); { gl[cnt++] = mwidth*nwidth*(p-kstart-1)+ mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth*(n-jstart-1) + m-istart-1; 226 227 /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 228 /* convert that to global numbering and get them on all processes */ 229 ierr = ISLocalToGlobalMappingApply(ltg,26,gl,gl);CHKERRQ(ierr); 230 /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 231 ierr = PetscMalloc(26*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr); 232 ierr = MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 233 234 /* Number the coarse grid points from 0 to Ntotal */ 235 ierr = PetscTableCreate(Ntotal/3,&ht);CHKERRQ(ierr); 236 for (i=0; i<26*mp*np*pp; i++){ 237 ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr); 238 } 239 ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr); 240 if (cnt != Ntotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal); 241 ierr = PetscFree(globals);CHKERRQ(ierr); 242 for (i=0; i<26; i++) { 243 ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr); 244 gl[i]--; 245 } 246 ierr = PetscTableDestroy(ht);CHKERRQ(ierr); 247 /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 248 249 /* construct global interpolation matrix */ 250 ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr); 251 if (reuse == MAT_INITIAL_MATRIX) { 252 ierr = MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint+Nsurf,PETSC_NULL,P);CHKERRQ(ierr); 253 } else { 254 ierr = MatZeroEntries(*P);CHKERRQ(ierr); 255 } 256 ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 257 ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 258 ierr = MatSetValues(*P,Nint,IIint,26,gl,xint,INSERT_VALUES);CHKERRQ(ierr); 259 ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 260 ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 261 ierr = MatSetValues(*P,Nsurf,IIsurf,26,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr); 262 ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 263 ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 264 ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 265 ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr); 266 267 #if defined(PETSC_USE_DEBUG_foo) 268 { 269 Vec x,y; 270 PetscScalar *yy; 271 ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr); 272 ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr); 273 ierr = VecSet(x,1.0);CHKERRQ(ierr); 274 ierr = MatMult(*P,x,y);CHKERRQ(ierr); 275 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 276 for (i=0; i<Ng; i++) { 277 if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %G",i,PetscAbsScalar(yy[i])); 278 } 279 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 280 ierr = VecDestroy(x);CHKERRQ(ierr); 281 ierr = VecDestroy(y);CHKERRQ(ierr); 282 } 283 #endif 284 285 ierr = MatDestroy(Aii);CHKERRQ(ierr); 286 ierr = MatDestroy(Ais);CHKERRQ(ierr); 287 ierr = MatDestroy(Asi);CHKERRQ(ierr); 288 ierr = MatDestroy(A);CHKERRQ(ierr); 289 ierr = ISDestroy(is);CHKERRQ(ierr); 290 ierr = ISDestroy(isint);CHKERRQ(ierr); 291 ierr = ISDestroy(issurf);CHKERRQ(ierr); 292 ierr = MatDestroy(Xint);CHKERRQ(ierr); 293 ierr = MatDestroy(Xsurf);CHKERRQ(ierr); 294 PetscFunctionReturn(0); 295 } 296 297 #undef __FUNCT__ 298 #define __FUNCT__ "DMDAGetFaceInterpolation" 299 /* 300 DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space 301 302 */ 303 PetscErrorCode DMDAGetFaceInterpolation(DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P) 304 { 305 PetscErrorCode ierr; 306 PetscInt dim,i,j,k,m,n,p,dof,Nint,Nface,Nwire,Nsurf,*Iint,*Isurf,cint = 0,csurf = 0,istart,jstart,kstart,*II,N,c = 0; 307 PetscInt mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf; 308 Mat Xint, Xsurf,Xint_tmp; 309 IS isint,issurf,is,row,col; 310 ISLocalToGlobalMapping ltg; 311 MPI_Comm comm; 312 Mat A,Aii,Ais,Asi,*Aholder,iAii; 313 MatFactorInfo info; 314 PetscScalar *xsurf,*xint; 315 #if defined(PETSC_USE_DEBUG_foo) 316 PetscScalar tmp; 317 #endif 318 PetscTable ht; 319 320 PetscFunctionBegin; 321 ierr = DMDAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);CHKERRQ(ierr); 322 if (dof != 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only for single field problems"); 323 if (dim != 3) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only coded for 3d problems"); 324 ierr = DMDAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr); 325 ierr = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr); 326 istart = istart ? -1 : 0; 327 jstart = jstart ? -1 : 0; 328 kstart = kstart ? -1 : 0; 329 330 /* 331 the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 332 to all the local degrees of freedom (this includes the vertices, edges and faces). 333 334 Xint are the subset of the interpolation into the interior 335 336 Xface are the interpolation onto faces but not into the interior 337 338 Xsurf are the interpolation onto the vertices and edges (the surfbasket) 339 Xint 340 Symbolically one could write P = ( Xface ) after interchanging the rows to match the natural ordering on the domain 341 Xsurf 342 */ 343 N = (m - istart)*(n - jstart)*(p - kstart); 344 Nint = (m-2-istart)*(n-2-jstart)*(p-2-kstart); 345 Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) ); 346 Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8; 347 Nsurf = Nface + Nwire; 348 ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,6,PETSC_NULL,&Xint);CHKERRQ(ierr); 349 ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,PETSC_NULL,&Xsurf);CHKERRQ(ierr); 350 ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 351 352 /* 353 Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 354 Xsurf will be all zero (thus making the coarse matrix singular). 355 */ 356 if (m-istart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3"); 357 if (n-jstart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3"); 358 if (p-kstart < 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3"); 359 360 cnt = 0; 361 for (j=1;j<n-1-jstart;j++) { for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 0*Nsurf] = 1;} } 362 for (k=1;k<p-1-kstart;k++) { 363 for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 1*Nsurf] = 1;} 364 for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 2*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 3*Nsurf] = 1;} 365 for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;} 366 } 367 for (j=1;j<n-1-jstart;j++) {for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 5*Nsurf] = 1;} } 368 369 #if defined(PETSC_USE_DEBUG_foo) 370 for (i=0; i<Nsurf; i++) { 371 tmp = 0.0; 372 for (j=0; j<6; j++) { 373 tmp += xsurf[i+j*Nsurf]; 374 } 375 if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %G",i,PetscAbsScalar(tmp)); 376 } 377 #endif 378 ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 379 /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/ 380 381 382 /* 383 I are the indices for all the needed vertices (in global numbering) 384 Iint are the indices for the interior values, I surf for the surface values 385 (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it 386 is NOT the local DMDA ordering.) 387 IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering 388 */ 389 #define Endpoint(a,start,b) (a == 0 || a == (b-1-start)) 390 ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr); 391 ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr); 392 for (k=0; k<p-kstart; k++) { 393 for (j=0; j<n-jstart; j++) { 394 for (i=0; i<m-istart; i++) { 395 II[c++] = i + j*mwidth + k*mwidth*nwidth; 396 397 if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) { 398 IIint[cint] = i + j*mwidth + k*mwidth*nwidth; 399 Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 400 } else { 401 IIsurf[csurf] = i + j*mwidth + k*mwidth*nwidth; 402 Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 403 } 404 } 405 } 406 } 407 if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N"); 408 if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint"); 409 if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf"); 410 ierr = DMGetLocalToGlobalMapping(da,<g);CHKERRQ(ierr); 411 ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr); 412 ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr); 413 ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr); 414 ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 415 ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 416 ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr); 417 ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr); 418 ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr); 419 420 ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr); 421 A = *Aholder; 422 ierr = PetscFree(Aholder);CHKERRQ(ierr); 423 424 ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr); 425 ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr); 426 ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr); 427 428 /* 429 Solve for the interpolation onto the interior Xint 430 */ 431 ierr = MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);CHKERRQ(ierr); 432 ierr = MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr); 433 ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr); 434 435 if (exotic->directSolve) { 436 ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr); 437 ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 438 ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr); 439 ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr); 440 ierr = ISDestroy(row);CHKERRQ(ierr); 441 ierr = ISDestroy(col);CHKERRQ(ierr); 442 ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr); 443 ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr); 444 ierr = MatDestroy(iAii);CHKERRQ(ierr); 445 } else { 446 Vec b,x; 447 PetscScalar *xint_tmp; 448 449 ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 450 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);CHKERRQ(ierr); 451 ierr = MatGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 452 ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);CHKERRQ(ierr); 453 ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 454 for (i=0; i<6; i++) { 455 ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr); 456 ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr); 457 ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr); 458 ierr = VecResetArray(x);CHKERRQ(ierr); 459 ierr = VecResetArray(b);CHKERRQ(ierr); 460 } 461 ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 462 ierr = MatRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 463 ierr = VecDestroy(x);CHKERRQ(ierr); 464 ierr = VecDestroy(b);CHKERRQ(ierr); 465 } 466 ierr = MatDestroy(Xint_tmp);CHKERRQ(ierr); 467 468 #if defined(PETSC_USE_DEBUG_foo) 469 ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 470 for (i=0; i<Nint; i++) { 471 tmp = 0.0; 472 for (j=0; j<6; j++) { 473 tmp += xint[i+j*Nint]; 474 } 475 if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %G",i,PetscAbsScalar(tmp)); 476 } 477 ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 478 /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 479 #endif 480 481 482 /* total faces */ 483 Ntotal = mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1); 484 485 /* 486 For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 487 */ 488 cnt = 0; 489 { gl[cnt++] = mwidth+1;} 490 { 491 { gl[cnt++] = mwidth*nwidth+1;} 492 { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;} 493 { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;} 494 } 495 { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;} 496 497 /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 498 /* convert that to global numbering and get them on all processes */ 499 ierr = ISLocalToGlobalMappingApply(ltg,6,gl,gl);CHKERRQ(ierr); 500 /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 501 ierr = PetscMalloc(6*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr); 502 ierr = MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 503 504 /* Number the coarse grid points from 0 to Ntotal */ 505 ierr = PetscTableCreate(Ntotal/3,&ht);CHKERRQ(ierr); 506 for (i=0; i<6*mp*np*pp; i++){ 507 ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr); 508 } 509 ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr); 510 if (cnt != Ntotal) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal); 511 ierr = PetscFree(globals);CHKERRQ(ierr); 512 for (i=0; i<6; i++) { 513 ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr); 514 gl[i]--; 515 } 516 ierr = PetscTableDestroy(ht);CHKERRQ(ierr); 517 /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 518 519 /* construct global interpolation matrix */ 520 ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr); 521 if (reuse == MAT_INITIAL_MATRIX) { 522 ierr = MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint,PETSC_NULL,P);CHKERRQ(ierr); 523 } else { 524 ierr = MatZeroEntries(*P);CHKERRQ(ierr); 525 } 526 ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 527 ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 528 ierr = MatSetValues(*P,Nint,IIint,6,gl,xint,INSERT_VALUES);CHKERRQ(ierr); 529 ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 530 ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 531 ierr = MatSetValues(*P,Nsurf,IIsurf,6,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr); 532 ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 533 ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 534 ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 535 ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr); 536 537 538 #if defined(PETSC_USE_DEBUG_foo) 539 { 540 Vec x,y; 541 PetscScalar *yy; 542 ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr); 543 ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr); 544 ierr = VecSet(x,1.0);CHKERRQ(ierr); 545 ierr = MatMult(*P,x,y);CHKERRQ(ierr); 546 ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 547 for (i=0; i<Ng; i++) { 548 if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %G",i,PetscAbsScalar(yy[i])); 549 } 550 ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 551 ierr = VecDestroy(x);CHKERRQ(ierr); 552 ierr = VecDestroy(y);CHKERRQ(ierr); 553 } 554 #endif 555 556 ierr = MatDestroy(Aii);CHKERRQ(ierr); 557 ierr = MatDestroy(Ais);CHKERRQ(ierr); 558 ierr = MatDestroy(Asi);CHKERRQ(ierr); 559 ierr = MatDestroy(A);CHKERRQ(ierr); 560 ierr = ISDestroy(is);CHKERRQ(ierr); 561 ierr = ISDestroy(isint);CHKERRQ(ierr); 562 ierr = ISDestroy(issurf);CHKERRQ(ierr); 563 ierr = MatDestroy(Xint);CHKERRQ(ierr); 564 ierr = MatDestroy(Xsurf);CHKERRQ(ierr); 565 PetscFunctionReturn(0); 566 } 567 568 569 #undef __FUNCT__ 570 #define __FUNCT__ "PCExoticSetType" 571 /*@ 572 PCExoticSetType - Sets the type of coarse grid interpolation to use 573 574 Logically Collective on PC 575 576 Input Parameters: 577 + pc - the preconditioner context 578 - type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face) 579 580 Notes: The face based interpolation has 1 degree of freedom per face and ignores the 581 edge and vertex values completely in the coarse problem. For any seven point 582 stencil the interpolation of a constant on all faces into the interior is that constant. 583 584 The wirebasket interpolation has 1 degree of freedom per vertex, per edge and 585 per face. A constant on the subdomain boundary is interpolated as that constant 586 in the interior of the domain. 587 588 The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence 589 if A is nonsingular A_c is also nonsingular. 590 591 Both interpolations are suitable for only scalar problems. 592 593 Level: intermediate 594 595 596 .seealso: PCEXOTIC, PCExoticType() 597 @*/ 598 PetscErrorCode PCExoticSetType(PC pc,PCExoticType type) 599 { 600 PetscErrorCode ierr; 601 602 PetscFunctionBegin; 603 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 604 PetscValidLogicalCollectiveEnum(pc,type,2); 605 ierr = PetscTryMethod(pc,"PCExoticSetType_C",(PC,PCExoticType),(pc,type));CHKERRQ(ierr); 606 PetscFunctionReturn(0); 607 } 608 609 #undef __FUNCT__ 610 #define __FUNCT__ "PCExoticSetType_Exotic" 611 PetscErrorCode PCExoticSetType_Exotic(PC pc,PCExoticType type) 612 { 613 PC_MG *mg = (PC_MG*)pc->data; 614 PC_Exotic *ctx = (PC_Exotic*) mg->innerctx; 615 616 PetscFunctionBegin; 617 ctx->type = type; 618 PetscFunctionReturn(0); 619 } 620 621 #undef __FUNCT__ 622 #define __FUNCT__ "PCSetUp_Exotic" 623 PetscErrorCode PCSetUp_Exotic(PC pc) 624 { 625 PetscErrorCode ierr; 626 Mat A; 627 PC_MG *mg = (PC_MG*)pc->data; 628 PC_Exotic *ex = (PC_Exotic*) mg->innerctx; 629 MatReuse reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX; 630 631 PetscFunctionBegin; 632 if (!pc->dm) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to call PCSetDM() before using this PC"); 633 ierr = PCGetOperators(pc,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr); 634 if (ex->type == PC_EXOTIC_FACE) { 635 ierr = DMDAGetFaceInterpolation(pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr); 636 } else if (ex->type == PC_EXOTIC_WIREBASKET) { 637 ierr = DMDAGetWireBasketInterpolation(pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr); 638 } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type); 639 ierr = PCMGSetInterpolation(pc,1,ex->P);CHKERRQ(ierr); 640 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 641 PetscFunctionReturn(0); 642 } 643 644 #undef __FUNCT__ 645 #define __FUNCT__ "PCDestroy_Exotic" 646 PetscErrorCode PCDestroy_Exotic(PC pc) 647 { 648 PetscErrorCode ierr; 649 PC_MG *mg = (PC_MG*)pc->data; 650 PC_Exotic *ctx = (PC_Exotic*) mg->innerctx; 651 652 PetscFunctionBegin; 653 if (ctx->P) {ierr = MatDestroy(ctx->P);CHKERRQ(ierr);} 654 if (ctx->ksp) {ierr = KSPDestroy(ctx->ksp);CHKERRQ(ierr);} 655 ierr = PetscFree(ctx);CHKERRQ(ierr); 656 ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 657 PetscFunctionReturn(0); 658 } 659 660 #undef __FUNCT__ 661 #define __FUNCT__ "PCView_Exotic" 662 PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer) 663 { 664 PC_MG *mg = (PC_MG*)pc->data; 665 PetscErrorCode ierr; 666 PetscBool iascii; 667 PC_Exotic *ctx = (PC_Exotic*) mg->innerctx; 668 669 PetscFunctionBegin; 670 ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 671 if (iascii) { 672 ierr = PetscViewerASCIIPrintf(viewer," Exotic type = %s\n",PCExoticTypes[ctx->type]);CHKERRQ(ierr); 673 if (ctx->directSolve) { 674 ierr = PetscViewerASCIIPrintf(viewer," Using direct solver to construct interpolation\n");CHKERRQ(ierr); 675 } else { 676 PetscViewer sviewer; 677 PetscMPIInt rank; 678 679 ierr = PetscViewerASCIIPrintf(viewer," Using iterative solver to construct interpolation\n");CHKERRQ(ierr); 680 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 681 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* should not need to push this twice? */ 682 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 683 ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 684 if (!rank) { 685 ierr = KSPView(ctx->ksp,sviewer);CHKERRQ(ierr); 686 } 687 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 688 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 689 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 690 } 691 } 692 ierr = PCView_MG(pc,viewer);CHKERRQ(ierr); 693 PetscFunctionReturn(0); 694 } 695 696 #undef __FUNCT__ 697 #define __FUNCT__ "PCSetFromOptions_Exotic" 698 PetscErrorCode PCSetFromOptions_Exotic(PC pc) 699 { 700 PetscErrorCode ierr; 701 PetscBool flg; 702 PC_MG *mg = (PC_MG*)pc->data; 703 PCExoticType mgctype; 704 PC_Exotic *ctx = (PC_Exotic*) mg->innerctx; 705 706 PetscFunctionBegin; 707 ierr = PetscOptionsHead("Exotic coarse space options");CHKERRQ(ierr); 708 ierr = PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 709 if (flg) { 710 ierr = PCExoticSetType(pc,mgctype);CHKERRQ(ierr); 711 } 712 ierr = PetscOptionsBool("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,PETSC_NULL);CHKERRQ(ierr); 713 if (!ctx->directSolve) { 714 if (!ctx->ksp) { 715 const char *prefix; 716 ierr = KSPCreate(PETSC_COMM_SELF,&ctx->ksp);CHKERRQ(ierr); 717 ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 718 ierr = PetscLogObjectParent(pc,ctx->ksp);CHKERRQ(ierr); 719 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 720 ierr = KSPSetOptionsPrefix(ctx->ksp,prefix);CHKERRQ(ierr); 721 ierr = KSPAppendOptionsPrefix(ctx->ksp,"exotic_");CHKERRQ(ierr); 722 } 723 ierr = KSPSetFromOptions(ctx->ksp);CHKERRQ(ierr); 724 } 725 ierr = PetscOptionsTail();CHKERRQ(ierr); 726 PetscFunctionReturn(0); 727 } 728 729 730 /*MC 731 PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces 732 733 This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse 734 grid spaces. 735 736 Notes: By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES 737 738 References: These coarse grid spaces originate in the work of Bramble, Pasciak and Schatz, "The Construction 739 of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53 pages 1--24, 1989. 740 They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith, 741 New York University, 1990. They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis 742 of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical 743 Analysis, volume 31. pages 1662-1694, 1994. These were developed in the context of iterative substructuring preconditioners. 744 They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund. 745 They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example, 746 Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco 747 Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings 748 of the 17th International Conference on Domain Decomposition Methods in 749 Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in 750 Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 255-261, 2007. 751 Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy min- 752 imizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer, 753 Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings 754 of the 17th International Conference on Domain Decomposition Methods 755 in Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in 756 Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 247-254, 2007 757 Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition 758 for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J. 759 Numer. Anal., 46(4):2153-2168, 2008. 760 Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz 761 algorithm for almost incompressible elasticity. Technical Report 762 TR2008-912, Department of Computer Science, Courant Institute 763 of Mathematical Sciences, New York University, May 2008. URL: 764 http://cs.nyu.edu/csweb/Research/TechReports/TR2008-912/TR2008-912.pdf 765 766 Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type> 767 -pc_mg_type <type> 768 769 Level: advanced 770 771 .seealso: PCMG, PCSetDM(), PCExoticType, PCExoticSetType() 772 M*/ 773 774 EXTERN_C_BEGIN 775 #undef __FUNCT__ 776 #define __FUNCT__ "PCCreate_Exotic" 777 PetscErrorCode PCCreate_Exotic(PC pc) 778 { 779 PetscErrorCode ierr; 780 PC_Exotic *ex; 781 PC_MG *mg; 782 783 PetscFunctionBegin; 784 /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 785 if (pc->ops->destroy) { ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); pc->data = 0;} 786 ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr); 787 ((PetscObject)pc)->type_name = 0; 788 789 ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); 790 ierr = PCMGSetLevels(pc,2,PETSC_NULL);CHKERRQ(ierr); 791 ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 792 ierr = PetscNew(PC_Exotic,&ex);CHKERRQ(ierr);\ 793 ex->type = PC_EXOTIC_FACE; 794 mg = (PC_MG*) pc->data; 795 mg->innerctx = ex; 796 797 798 pc->ops->setfromoptions = PCSetFromOptions_Exotic; 799 pc->ops->view = PCView_Exotic; 800 pc->ops->destroy = PCDestroy_Exotic; 801 pc->ops->setup = PCSetUp_Exotic; 802 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCExoticSetType_C","PCExoticSetType_Exotic",PCExoticSetType_Exotic);CHKERRQ(ierr); 803 PetscFunctionReturn(0); 804 } 805 EXTERN_C_END 806