1174b6946SBarry Smith #define PETSCKSP_DLL 2174b6946SBarry Smith 3174b6946SBarry Smith 4174b6946SBarry Smith #include "petscpc.h" /*I "petscpc.h" I*/ 5a9f2baa8SLisandro Dalcin #include "petscmg.h" /*I "petscmg.h" I*/ 6174b6946SBarry Smith #include "petscda.h" /*I "petscda.h" I*/ 77233f9f0SBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h" 87233f9f0SBarry Smith 9*8e722e36SBarry Smith typedef struct { 10*8e722e36SBarry Smith DA da; 11*8e722e36SBarry Smith PCExoticType type; 12*8e722e36SBarry Smith Mat P; /* the constructed interpolation matrix */ 13*8e722e36SBarry Smith PetscTruth directSolve; /* use direct LU factorization to construct interpolation */ 14*8e722e36SBarry Smith KSP ksp; 15*8e722e36SBarry Smith } PC_Exotic; 16174b6946SBarry Smith 17*8e722e36SBarry Smith const char *PCExoticTypes[] = {"face","wirebasket","PCExoticType","PC_Exotic",0}; 18064c8009SBarry Smith 19064c8009SBarry Smith 20064c8009SBarry Smith #undef __FUNCT__ 21064c8009SBarry Smith #define __FUNCT__ "DAGetWireBasketInterpolation" 22064c8009SBarry Smith /* 23064c8009SBarry Smith DAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space 24064c8009SBarry Smith 25064c8009SBarry Smith */ 26*8e722e36SBarry Smith PetscErrorCode DAGetWireBasketInterpolation(PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P) 27064c8009SBarry Smith { 28*8e722e36SBarry Smith DA da = exotic->da; 29064c8009SBarry Smith PetscErrorCode ierr; 30064c8009SBarry Smith 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; 31064c8009SBarry Smith PetscInt mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf; 32064c8009SBarry Smith Mat Xint, Xsurf,Xint_tmp; 33064c8009SBarry Smith IS isint,issurf,is,row,col; 34064c8009SBarry Smith ISLocalToGlobalMapping ltg; 35064c8009SBarry Smith MPI_Comm comm; 36064c8009SBarry Smith Mat A,Aii,Ais,Asi,*Aholder,iAii; 37064c8009SBarry Smith MatFactorInfo info; 38064c8009SBarry Smith PetscScalar *xsurf,*xint; 39*8e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 40064c8009SBarry Smith PetscScalar tmp; 41064c8009SBarry Smith #endif 42064c8009SBarry Smith PetscTable ht; 43064c8009SBarry Smith 44064c8009SBarry Smith PetscFunctionBegin; 45064c8009SBarry Smith ierr = DAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);CHKERRQ(ierr); 46064c8009SBarry Smith if (dof != 1) SETERRQ(PETSC_ERR_SUP,"Only for single field problems"); 47064c8009SBarry Smith if (dim != 3) SETERRQ(PETSC_ERR_SUP,"Only coded for 3d problems"); 48064c8009SBarry Smith ierr = DAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr); 49064c8009SBarry Smith ierr = DAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr); 50064c8009SBarry Smith istart = istart ? -1 : 0; 51064c8009SBarry Smith jstart = jstart ? -1 : 0; 52064c8009SBarry Smith kstart = kstart ? -1 : 0; 53064c8009SBarry Smith 54064c8009SBarry Smith /* 55064c8009SBarry Smith the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 56064c8009SBarry Smith to all the local degrees of freedom (this includes the vertices, edges and faces). 57064c8009SBarry Smith 58064c8009SBarry Smith Xint are the subset of the interpolation into the interior 59064c8009SBarry Smith 60064c8009SBarry Smith Xface are the interpolation onto faces but not into the interior 61064c8009SBarry Smith 62064c8009SBarry Smith Xsurf are the interpolation onto the vertices and edges (the surfbasket) 63064c8009SBarry Smith Xint 64064c8009SBarry Smith Symbolically one could write P = ( Xface ) after interchanging the rows to match the natural ordering on the domain 65064c8009SBarry Smith Xsurf 66064c8009SBarry Smith */ 67064c8009SBarry Smith N = (m - istart)*(n - jstart)*(p - kstart); 68064c8009SBarry Smith Nint = (m-2-istart)*(n-2-jstart)*(p-2-kstart); 69064c8009SBarry Smith Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) ); 70064c8009SBarry Smith Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8; 71064c8009SBarry Smith Nsurf = Nface + Nwire; 72064c8009SBarry Smith ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,26,PETSC_NULL,&Xint);CHKERRQ(ierr); 73064c8009SBarry Smith ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,PETSC_NULL,&Xsurf);CHKERRQ(ierr); 74064c8009SBarry Smith ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 75064c8009SBarry Smith 76064c8009SBarry Smith /* 77064c8009SBarry Smith Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 78064c8009SBarry Smith Xsurf will be all zero (thus making the coarse matrix singular). 79064c8009SBarry Smith */ 80064c8009SBarry Smith if (m-istart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3"); 81064c8009SBarry Smith if (n-jstart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3"); 82064c8009SBarry Smith if (p-kstart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3"); 83064c8009SBarry Smith 84064c8009SBarry Smith cnt = 0; 85064c8009SBarry Smith xsurf[cnt++] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + Nsurf] = 1;} xsurf[cnt++ + 2*Nsurf] = 1; 86064c8009SBarry Smith 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;} 87064c8009SBarry Smith xsurf[cnt++ + 6*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 7*Nsurf] = 1;} xsurf[cnt++ + 8*Nsurf] = 1; 88064c8009SBarry Smith for (k=1;k<p-1-kstart;k++) { 89064c8009SBarry Smith xsurf[cnt++ + 9*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 10*Nsurf] = 1;} xsurf[cnt++ + 11*Nsurf] = 1; 90064c8009SBarry Smith for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 12*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 13*Nsurf] = 1;} 91064c8009SBarry Smith xsurf[cnt++ + 14*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 15*Nsurf] = 1;} xsurf[cnt++ + 16*Nsurf] = 1; 92064c8009SBarry Smith } 93064c8009SBarry Smith xsurf[cnt++ + 17*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 18*Nsurf] = 1;} xsurf[cnt++ + 19*Nsurf] = 1; 94064c8009SBarry Smith 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;} 95064c8009SBarry Smith xsurf[cnt++ + 23*Nsurf] = 1; for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 24*Nsurf] = 1;} xsurf[cnt++ + 25*Nsurf] = 1; 96064c8009SBarry Smith 97*8e722e36SBarry Smith /* interpolations only sum to 1 when using direct solver */ 98*8e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 99064c8009SBarry Smith for (i=0; i<Nsurf; i++) { 100064c8009SBarry Smith tmp = 0.0; 101064c8009SBarry Smith for (j=0; j<26; j++) { 102064c8009SBarry Smith tmp += xsurf[i+j*Nsurf]; 103064c8009SBarry Smith } 104064c8009SBarry Smith if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %G",i,PetscAbsScalar(tmp)); 105064c8009SBarry Smith } 106064c8009SBarry Smith #endif 107064c8009SBarry Smith ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 108064c8009SBarry Smith /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/ 109064c8009SBarry Smith 110064c8009SBarry Smith 111064c8009SBarry Smith /* 112064c8009SBarry Smith I are the indices for all the needed vertices (in global numbering) 113064c8009SBarry Smith Iint are the indices for the interior values, I surf for the surface values 114064c8009SBarry Smith (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it 115064c8009SBarry Smith is NOT the local DA ordering.) 116064c8009SBarry Smith IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering 117064c8009SBarry Smith */ 118064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start)) 119064c8009SBarry Smith ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr); 120064c8009SBarry Smith ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr); 121064c8009SBarry Smith for (k=0; k<p-kstart; k++) { 122064c8009SBarry Smith for (j=0; j<n-jstart; j++) { 123064c8009SBarry Smith for (i=0; i<m-istart; i++) { 124064c8009SBarry Smith II[c++] = i + j*mwidth + k*mwidth*nwidth; 125064c8009SBarry Smith 126064c8009SBarry Smith if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) { 127064c8009SBarry Smith IIint[cint] = i + j*mwidth + k*mwidth*nwidth; 128064c8009SBarry Smith Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 129064c8009SBarry Smith } else { 130064c8009SBarry Smith IIsurf[csurf] = i + j*mwidth + k*mwidth*nwidth; 131064c8009SBarry Smith Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 132064c8009SBarry Smith } 133064c8009SBarry Smith } 134064c8009SBarry Smith } 135064c8009SBarry Smith } 136064c8009SBarry Smith if (c != N) SETERRQ(PETSC_ERR_PLIB,"c != N"); 137064c8009SBarry Smith if (cint != Nint) SETERRQ(PETSC_ERR_PLIB,"cint != Nint"); 138064c8009SBarry Smith if (csurf != Nsurf) SETERRQ(PETSC_ERR_PLIB,"csurf != Nsurf"); 139064c8009SBarry Smith ierr = DAGetISLocalToGlobalMapping(da,<g);CHKERRQ(ierr); 140064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr); 141064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr); 142064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr); 143064c8009SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 144064c8009SBarry Smith ierr = ISCreateGeneral(comm,N,II,&is);CHKERRQ(ierr); 145064c8009SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,&isint);CHKERRQ(ierr); 146064c8009SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,&issurf);CHKERRQ(ierr); 147064c8009SBarry Smith ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr); 148064c8009SBarry Smith 149064c8009SBarry Smith ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr); 150064c8009SBarry Smith A = *Aholder; 151064c8009SBarry Smith ierr = PetscFree(Aholder);CHKERRQ(ierr); 152064c8009SBarry Smith 153064c8009SBarry Smith ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr); 154064c8009SBarry Smith ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr); 155064c8009SBarry Smith ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr); 156064c8009SBarry Smith 157064c8009SBarry Smith /* 158064c8009SBarry Smith Solve for the interpolation onto the interior Xint 159064c8009SBarry Smith */ 160*8e722e36SBarry Smith ierr = MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);CHKERRQ(ierr); 161*8e722e36SBarry Smith ierr = MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr); 162*8e722e36SBarry Smith ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr); 163*8e722e36SBarry Smith if (exotic->directSolve) { 164064c8009SBarry Smith ierr = MatGetFactor(Aii,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr); 165064c8009SBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 166064c8009SBarry Smith ierr = MatGetOrdering(Aii,MATORDERING_ND,&row,&col);CHKERRQ(ierr); 167064c8009SBarry Smith ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr); 168064c8009SBarry Smith ierr = ISDestroy(row);CHKERRQ(ierr); 169064c8009SBarry Smith ierr = ISDestroy(col);CHKERRQ(ierr); 170064c8009SBarry Smith ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr); 171064c8009SBarry Smith ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr); 172064c8009SBarry Smith ierr = MatDestroy(iAii);CHKERRQ(ierr); 173*8e722e36SBarry Smith } else { 174*8e722e36SBarry Smith Vec b,x; 175*8e722e36SBarry Smith PetscScalar *xint_tmp; 176064c8009SBarry Smith 177*8e722e36SBarry Smith ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 178*8e722e36SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);CHKERRQ(ierr); 179*8e722e36SBarry Smith ierr = MatGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 180*8e722e36SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);CHKERRQ(ierr); 181*8e722e36SBarry Smith ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 182*8e722e36SBarry Smith for (i=0; i<26; i++) { 183*8e722e36SBarry Smith ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr); 184*8e722e36SBarry Smith ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr); 185*8e722e36SBarry Smith ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr); 186*8e722e36SBarry Smith ierr = VecResetArray(x);CHKERRQ(ierr); 187*8e722e36SBarry Smith ierr = VecResetArray(b);CHKERRQ(ierr); 188*8e722e36SBarry Smith } 189*8e722e36SBarry Smith ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 190*8e722e36SBarry Smith ierr = MatRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 191*8e722e36SBarry Smith ierr = VecDestroy(x);CHKERRQ(ierr); 192*8e722e36SBarry Smith ierr = VecDestroy(b);CHKERRQ(ierr); 193*8e722e36SBarry Smith } 194*8e722e36SBarry Smith ierr = MatDestroy(Xint_tmp);CHKERRQ(ierr); 195*8e722e36SBarry Smith 196*8e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 197064c8009SBarry Smith ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 198064c8009SBarry Smith for (i=0; i<Nint; i++) { 199064c8009SBarry Smith tmp = 0.0; 200064c8009SBarry Smith for (j=0; j<26; j++) { 201064c8009SBarry Smith tmp += xint[i+j*Nint]; 202064c8009SBarry Smith } 203064c8009SBarry Smith if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %G",i,PetscAbsScalar(tmp)); 204064c8009SBarry Smith } 205064c8009SBarry Smith ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 206064c8009SBarry Smith /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 207064c8009SBarry Smith #endif 208064c8009SBarry Smith 209064c8009SBarry Smith 210064c8009SBarry Smith /* total vertices total faces total edges */ 211064c8009SBarry Smith 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); 212064c8009SBarry Smith 213064c8009SBarry Smith /* 214064c8009SBarry Smith For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 215064c8009SBarry Smith */ 216064c8009SBarry Smith cnt = 0; 217064c8009SBarry Smith gl[cnt++] = 0; { gl[cnt++] = 1;} gl[cnt++] = m-istart-1; 218064c8009SBarry Smith { gl[cnt++] = mwidth; { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;} 219064c8009SBarry Smith gl[cnt++] = mwidth*(n-jstart-1); { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1; 220064c8009SBarry Smith { 221064c8009SBarry Smith gl[cnt++] = mwidth*nwidth; { gl[cnt++] = mwidth*nwidth+1;} gl[cnt++] = mwidth*nwidth+ m-istart-1; 222064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;} 223064c8009SBarry Smith 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; 224064c8009SBarry Smith } 225064c8009SBarry Smith 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; 226064c8009SBarry Smith { 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;} 227064c8009SBarry Smith 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; 228064c8009SBarry Smith 229064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 230064c8009SBarry Smith /* convert that to global numbering and get them on all processes */ 231064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,26,gl,gl);CHKERRQ(ierr); 232064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 233064c8009SBarry Smith ierr = PetscMalloc(26*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr); 234064c8009SBarry Smith ierr = MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 235064c8009SBarry Smith 236064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */ 237064c8009SBarry Smith ierr = PetscTableCreate(Ntotal/3,&ht);CHKERRQ(ierr); 238064c8009SBarry Smith for (i=0; i<26*mp*np*pp; i++){ 239064c8009SBarry Smith ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr); 240064c8009SBarry Smith } 241064c8009SBarry Smith ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr); 242064c8009SBarry Smith if (cnt != Ntotal) SETERRQ2(PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal); 243064c8009SBarry Smith ierr = PetscFree(globals);CHKERRQ(ierr); 244064c8009SBarry Smith for (i=0; i<26; i++) { 245064c8009SBarry Smith ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr); 246064c8009SBarry Smith gl[i]--; 247064c8009SBarry Smith } 248064c8009SBarry Smith ierr = PetscTableDestroy(ht);CHKERRQ(ierr); 249064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 250064c8009SBarry Smith 251064c8009SBarry Smith /* construct global interpolation matrix */ 252064c8009SBarry Smith ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr); 253064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 254064c8009SBarry Smith ierr = MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint+Nsurf,PETSC_NULL,P);CHKERRQ(ierr); 255064c8009SBarry Smith } else { 256064c8009SBarry Smith ierr = MatZeroEntries(*P);CHKERRQ(ierr); 257064c8009SBarry Smith } 258064c8009SBarry Smith ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 259064c8009SBarry Smith ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 260064c8009SBarry Smith ierr = MatSetValues(*P,Nint,IIint,26,gl,xint,INSERT_VALUES);CHKERRQ(ierr); 261064c8009SBarry Smith ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 262064c8009SBarry Smith ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 263064c8009SBarry Smith ierr = MatSetValues(*P,Nsurf,IIsurf,26,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr); 264064c8009SBarry Smith ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 265064c8009SBarry Smith ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 266064c8009SBarry Smith ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 267064c8009SBarry Smith ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr); 268064c8009SBarry Smith 269*8e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 270064c8009SBarry Smith { 271064c8009SBarry Smith Vec x,y; 272064c8009SBarry Smith PetscScalar *yy; 273064c8009SBarry Smith ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr); 274064c8009SBarry Smith ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr); 275064c8009SBarry Smith ierr = VecSet(x,1.0);CHKERRQ(ierr); 276064c8009SBarry Smith ierr = MatMult(*P,x,y);CHKERRQ(ierr); 277064c8009SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 278064c8009SBarry Smith for (i=0; i<Ng; i++) { 279064c8009SBarry Smith if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %G",i,PetscAbsScalar(yy[i])); 280064c8009SBarry Smith } 281064c8009SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 282064c8009SBarry Smith ierr = VecDestroy(x);CHKERRQ(ierr); 283064c8009SBarry Smith ierr = VecDestroy(y);CHKERRQ(ierr); 284064c8009SBarry Smith } 285064c8009SBarry Smith #endif 286064c8009SBarry Smith 287064c8009SBarry Smith ierr = MatDestroy(Aii);CHKERRQ(ierr); 288064c8009SBarry Smith ierr = MatDestroy(Ais);CHKERRQ(ierr); 289064c8009SBarry Smith ierr = MatDestroy(Asi);CHKERRQ(ierr); 290064c8009SBarry Smith ierr = MatDestroy(A);CHKERRQ(ierr); 291064c8009SBarry Smith ierr = ISDestroy(is);CHKERRQ(ierr); 292064c8009SBarry Smith ierr = ISDestroy(isint);CHKERRQ(ierr); 293064c8009SBarry Smith ierr = ISDestroy(issurf);CHKERRQ(ierr); 294064c8009SBarry Smith ierr = MatDestroy(Xint);CHKERRQ(ierr); 295064c8009SBarry Smith ierr = MatDestroy(Xsurf);CHKERRQ(ierr); 296064c8009SBarry Smith PetscFunctionReturn(0); 297064c8009SBarry Smith } 298064c8009SBarry Smith 299064c8009SBarry Smith #undef __FUNCT__ 300064c8009SBarry Smith #define __FUNCT__ "DAGetFaceInterpolation" 301064c8009SBarry Smith /* 302064c8009SBarry Smith DAGetFaceInterpolation - Gets the interpolation for a face based coarse space 303064c8009SBarry Smith 304064c8009SBarry Smith */ 305*8e722e36SBarry Smith PetscErrorCode DAGetFaceInterpolation(PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P) 306064c8009SBarry Smith { 307*8e722e36SBarry Smith DA da = exotic->da; 308064c8009SBarry Smith PetscErrorCode ierr; 309064c8009SBarry Smith 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; 310064c8009SBarry Smith PetscInt mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf; 311064c8009SBarry Smith Mat Xint, Xsurf,Xint_tmp; 312064c8009SBarry Smith IS isint,issurf,is,row,col; 313064c8009SBarry Smith ISLocalToGlobalMapping ltg; 314064c8009SBarry Smith MPI_Comm comm; 315064c8009SBarry Smith Mat A,Aii,Ais,Asi,*Aholder,iAii; 316064c8009SBarry Smith MatFactorInfo info; 317064c8009SBarry Smith PetscScalar *xsurf,*xint; 318064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 319064c8009SBarry Smith PetscScalar tmp; 320064c8009SBarry Smith #endif 321064c8009SBarry Smith PetscTable ht; 322064c8009SBarry Smith 323064c8009SBarry Smith PetscFunctionBegin; 324064c8009SBarry Smith ierr = DAGetInfo(da,&dim,0,0,0,&mp,&np,&pp,&dof,0,0,0);CHKERRQ(ierr); 325064c8009SBarry Smith if (dof != 1) SETERRQ(PETSC_ERR_SUP,"Only for single field problems"); 326064c8009SBarry Smith if (dim != 3) SETERRQ(PETSC_ERR_SUP,"Only coded for 3d problems"); 327064c8009SBarry Smith ierr = DAGetCorners(da,0,0,0,&m,&n,&p);CHKERRQ(ierr); 328064c8009SBarry Smith ierr = DAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr); 329064c8009SBarry Smith istart = istart ? -1 : 0; 330064c8009SBarry Smith jstart = jstart ? -1 : 0; 331064c8009SBarry Smith kstart = kstart ? -1 : 0; 332064c8009SBarry Smith 333064c8009SBarry Smith /* 334064c8009SBarry Smith the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 335064c8009SBarry Smith to all the local degrees of freedom (this includes the vertices, edges and faces). 336064c8009SBarry Smith 337064c8009SBarry Smith Xint are the subset of the interpolation into the interior 338064c8009SBarry Smith 339064c8009SBarry Smith Xface are the interpolation onto faces but not into the interior 340064c8009SBarry Smith 341064c8009SBarry Smith Xsurf are the interpolation onto the vertices and edges (the surfbasket) 342064c8009SBarry Smith Xint 343064c8009SBarry Smith Symbolically one could write P = ( Xface ) after interchanging the rows to match the natural ordering on the domain 344064c8009SBarry Smith Xsurf 345064c8009SBarry Smith */ 346064c8009SBarry Smith N = (m - istart)*(n - jstart)*(p - kstart); 347064c8009SBarry Smith Nint = (m-2-istart)*(n-2-jstart)*(p-2-kstart); 348064c8009SBarry Smith Nface = 2*( (m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart) ); 349064c8009SBarry Smith Nwire = 4*( (m-2-istart) + (n-2-jstart) + (p-2-kstart) ) + 8; 350064c8009SBarry Smith Nsurf = Nface + Nwire; 351064c8009SBarry Smith ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,6,PETSC_NULL,&Xint);CHKERRQ(ierr); 352064c8009SBarry Smith ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,PETSC_NULL,&Xsurf);CHKERRQ(ierr); 353064c8009SBarry Smith ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 354064c8009SBarry Smith 355064c8009SBarry Smith /* 356064c8009SBarry Smith Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 357064c8009SBarry Smith Xsurf will be all zero (thus making the coarse matrix singular). 358064c8009SBarry Smith */ 359064c8009SBarry Smith if (m-istart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in X direction must be at least 3"); 360064c8009SBarry Smith if (n-jstart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in Y direction must be at least 3"); 361064c8009SBarry Smith if (p-kstart < 3) SETERRQ(PETSC_ERR_SUP,"Number of grid points per process in Z direction must be at least 3"); 362064c8009SBarry Smith 363064c8009SBarry Smith cnt = 0; 364064c8009SBarry Smith for (j=1;j<n-1-jstart;j++) { for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 0*Nsurf] = 1;} } 365064c8009SBarry Smith for (k=1;k<p-1-kstart;k++) { 366064c8009SBarry Smith for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 1*Nsurf] = 1;} 367064c8009SBarry Smith for (j=1;j<n-1-jstart;j++) { xsurf[cnt++ + 2*Nsurf] = 1; /* these are the interior nodes */ xsurf[cnt++ + 3*Nsurf] = 1;} 368064c8009SBarry Smith for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 4*Nsurf] = 1;} 369064c8009SBarry Smith } 370064c8009SBarry Smith for (j=1;j<n-1-jstart;j++) {for (i=1; i<m-istart-1; i++) { xsurf[cnt++ + 5*Nsurf] = 1;} } 371064c8009SBarry Smith 372064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 373064c8009SBarry Smith for (i=0; i<Nsurf; i++) { 374064c8009SBarry Smith tmp = 0.0; 375064c8009SBarry Smith for (j=0; j<6; j++) { 376064c8009SBarry Smith tmp += xsurf[i+j*Nsurf]; 377064c8009SBarry Smith } 378064c8009SBarry Smith if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %G",i,PetscAbsScalar(tmp)); 379064c8009SBarry Smith } 380064c8009SBarry Smith #endif 381064c8009SBarry Smith ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 382064c8009SBarry Smith /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/ 383064c8009SBarry Smith 384064c8009SBarry Smith 385064c8009SBarry Smith /* 386064c8009SBarry Smith I are the indices for all the needed vertices (in global numbering) 387064c8009SBarry Smith Iint are the indices for the interior values, I surf for the surface values 388064c8009SBarry Smith (This is just for the part of the global matrix obtained with MatGetSubMatrix(), it 389064c8009SBarry Smith is NOT the local DA ordering.) 390064c8009SBarry Smith IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering 391064c8009SBarry Smith */ 392064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start)) 393064c8009SBarry Smith ierr = PetscMalloc3(N,PetscInt,&II,Nint,PetscInt,&Iint,Nsurf,PetscInt,&Isurf);CHKERRQ(ierr); 394064c8009SBarry Smith ierr = PetscMalloc2(Nint,PetscInt,&IIint,Nsurf,PetscInt,&IIsurf);CHKERRQ(ierr); 395064c8009SBarry Smith for (k=0; k<p-kstart; k++) { 396064c8009SBarry Smith for (j=0; j<n-jstart; j++) { 397064c8009SBarry Smith for (i=0; i<m-istart; i++) { 398064c8009SBarry Smith II[c++] = i + j*mwidth + k*mwidth*nwidth; 399064c8009SBarry Smith 400064c8009SBarry Smith if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) { 401064c8009SBarry Smith IIint[cint] = i + j*mwidth + k*mwidth*nwidth; 402064c8009SBarry Smith Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 403064c8009SBarry Smith } else { 404064c8009SBarry Smith IIsurf[csurf] = i + j*mwidth + k*mwidth*nwidth; 405064c8009SBarry Smith Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 406064c8009SBarry Smith } 407064c8009SBarry Smith } 408064c8009SBarry Smith } 409064c8009SBarry Smith } 410064c8009SBarry Smith if (c != N) SETERRQ(PETSC_ERR_PLIB,"c != N"); 411064c8009SBarry Smith if (cint != Nint) SETERRQ(PETSC_ERR_PLIB,"cint != Nint"); 412064c8009SBarry Smith if (csurf != Nsurf) SETERRQ(PETSC_ERR_PLIB,"csurf != Nsurf"); 413064c8009SBarry Smith ierr = DAGetISLocalToGlobalMapping(da,<g);CHKERRQ(ierr); 414064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr); 415064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr); 416064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr); 417064c8009SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 418064c8009SBarry Smith ierr = ISCreateGeneral(comm,N,II,&is);CHKERRQ(ierr); 419064c8009SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,&isint);CHKERRQ(ierr); 420064c8009SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,&issurf);CHKERRQ(ierr); 421064c8009SBarry Smith ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr); 422064c8009SBarry Smith 423064c8009SBarry Smith ierr = MatGetSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr); 424064c8009SBarry Smith A = *Aholder; 425064c8009SBarry Smith ierr = PetscFree(Aholder);CHKERRQ(ierr); 426064c8009SBarry Smith 427064c8009SBarry Smith ierr = MatGetSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr); 428064c8009SBarry Smith ierr = MatGetSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr); 429064c8009SBarry Smith ierr = MatGetSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr); 430064c8009SBarry Smith 431064c8009SBarry Smith /* 432064c8009SBarry Smith Solve for the interpolation onto the interior Xint 433064c8009SBarry Smith */ 434064c8009SBarry Smith ierr = MatDuplicate(Xint,MAT_DO_NOT_COPY_VALUES,&Xint_tmp);CHKERRQ(ierr); 435064c8009SBarry Smith ierr = MatMatMult(Ais,Xsurf,MAT_REUSE_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr); 436064c8009SBarry Smith ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr); 437064c8009SBarry Smith 438*8e722e36SBarry Smith if (exotic->directSolve) { 439064c8009SBarry Smith ierr = MatGetFactor(Aii,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr); 440064c8009SBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 441064c8009SBarry Smith ierr = MatGetOrdering(Aii,MATORDERING_ND,&row,&col);CHKERRQ(ierr); 442064c8009SBarry Smith ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr); 443064c8009SBarry Smith ierr = ISDestroy(row);CHKERRQ(ierr); 444064c8009SBarry Smith ierr = ISDestroy(col);CHKERRQ(ierr); 445064c8009SBarry Smith ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr); 446064c8009SBarry Smith ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr); 447064c8009SBarry Smith ierr = MatDestroy(iAii);CHKERRQ(ierr); 448064c8009SBarry Smith } else { 449064c8009SBarry Smith Vec b,x; 450064c8009SBarry Smith PetscScalar *xint_tmp; 451064c8009SBarry Smith 452064c8009SBarry Smith ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 453064c8009SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&x);CHKERRQ(ierr); 454064c8009SBarry Smith ierr = MatGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 455064c8009SBarry Smith ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Nint,0,&b);CHKERRQ(ierr); 456*8e722e36SBarry Smith ierr = KSPSetOperators(exotic->ksp,Aii,Aii,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 457064c8009SBarry Smith for (i=0; i<6; i++) { 458064c8009SBarry Smith ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr); 459064c8009SBarry Smith ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr); 460*8e722e36SBarry Smith ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr); 461064c8009SBarry Smith ierr = VecResetArray(x);CHKERRQ(ierr); 462064c8009SBarry Smith ierr = VecResetArray(b);CHKERRQ(ierr); 463064c8009SBarry Smith } 464064c8009SBarry Smith ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 465064c8009SBarry Smith ierr = MatRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 466064c8009SBarry Smith ierr = VecDestroy(x);CHKERRQ(ierr); 467064c8009SBarry Smith ierr = VecDestroy(b);CHKERRQ(ierr); 468064c8009SBarry Smith } 469064c8009SBarry Smith ierr = MatDestroy(Xint_tmp);CHKERRQ(ierr); 470064c8009SBarry Smith 471064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 472064c8009SBarry Smith ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 473064c8009SBarry Smith for (i=0; i<Nint; i++) { 474064c8009SBarry Smith tmp = 0.0; 475064c8009SBarry Smith for (j=0; j<6; j++) { 476064c8009SBarry Smith tmp += xint[i+j*Nint]; 477064c8009SBarry Smith } 478064c8009SBarry Smith if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %G",i,PetscAbsScalar(tmp)); 479064c8009SBarry Smith } 480064c8009SBarry Smith ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 481064c8009SBarry Smith /* ierr =MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 482064c8009SBarry Smith #endif 483064c8009SBarry Smith 484064c8009SBarry Smith 485064c8009SBarry Smith /* total faces */ 486064c8009SBarry Smith Ntotal = mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1); 487064c8009SBarry Smith 488064c8009SBarry Smith /* 489064c8009SBarry Smith For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 490064c8009SBarry Smith */ 491064c8009SBarry Smith cnt = 0; 492064c8009SBarry Smith { gl[cnt++] = mwidth+1;} 493064c8009SBarry Smith { 494064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth+1;} 495064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;} 496064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;} 497064c8009SBarry Smith } 498064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;} 499064c8009SBarry Smith 500064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 501064c8009SBarry Smith /* convert that to global numbering and get them on all processes */ 502064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,6,gl,gl);CHKERRQ(ierr); 503064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 504064c8009SBarry Smith ierr = PetscMalloc(6*mp*np*pp*sizeof(PetscInt),&globals);CHKERRQ(ierr); 505064c8009SBarry Smith ierr = MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,((PetscObject)da)->comm);CHKERRQ(ierr); 506064c8009SBarry Smith 507064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */ 508064c8009SBarry Smith ierr = PetscTableCreate(Ntotal/3,&ht);CHKERRQ(ierr); 509064c8009SBarry Smith for (i=0; i<6*mp*np*pp; i++){ 510064c8009SBarry Smith ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr); 511064c8009SBarry Smith } 512064c8009SBarry Smith ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr); 513064c8009SBarry Smith if (cnt != Ntotal) SETERRQ2(PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal); 514064c8009SBarry Smith ierr = PetscFree(globals);CHKERRQ(ierr); 515064c8009SBarry Smith for (i=0; i<6; i++) { 516064c8009SBarry Smith ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr); 517064c8009SBarry Smith gl[i]--; 518064c8009SBarry Smith } 519064c8009SBarry Smith ierr = PetscTableDestroy(ht);CHKERRQ(ierr); 520064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 521064c8009SBarry Smith 522064c8009SBarry Smith /* construct global interpolation matrix */ 523064c8009SBarry Smith ierr = MatGetLocalSize(Aglobal,&Ng,PETSC_NULL);CHKERRQ(ierr); 524064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 525064c8009SBarry Smith ierr = MatCreateMPIAIJ(((PetscObject)da)->comm,Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,PETSC_NULL,Nint,PETSC_NULL,P);CHKERRQ(ierr); 526064c8009SBarry Smith } else { 527064c8009SBarry Smith ierr = MatZeroEntries(*P);CHKERRQ(ierr); 528064c8009SBarry Smith } 529064c8009SBarry Smith ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 530064c8009SBarry Smith ierr = MatGetArray(Xint,&xint);CHKERRQ(ierr); 531064c8009SBarry Smith ierr = MatSetValues(*P,Nint,IIint,6,gl,xint,INSERT_VALUES);CHKERRQ(ierr); 532064c8009SBarry Smith ierr = MatRestoreArray(Xint,&xint);CHKERRQ(ierr); 533064c8009SBarry Smith ierr = MatGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 534064c8009SBarry Smith ierr = MatSetValues(*P,Nsurf,IIsurf,6,gl,xsurf,INSERT_VALUES);CHKERRQ(ierr); 535064c8009SBarry Smith ierr = MatRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 536064c8009SBarry Smith ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 537064c8009SBarry Smith ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 538064c8009SBarry Smith ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr); 539064c8009SBarry Smith 540064c8009SBarry Smith 541064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 542064c8009SBarry Smith { 543064c8009SBarry Smith Vec x,y; 544064c8009SBarry Smith PetscScalar *yy; 545064c8009SBarry Smith ierr = VecCreateMPI(((PetscObject)da)->comm,Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr); 546064c8009SBarry Smith ierr = VecCreateMPI(((PetscObject)da)->comm,PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr); 547064c8009SBarry Smith ierr = VecSet(x,1.0);CHKERRQ(ierr); 548064c8009SBarry Smith ierr = MatMult(*P,x,y);CHKERRQ(ierr); 549064c8009SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 550064c8009SBarry Smith for (i=0; i<Ng; i++) { 551064c8009SBarry Smith if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ2(PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %G",i,PetscAbsScalar(yy[i])); 552064c8009SBarry Smith } 553064c8009SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 554064c8009SBarry Smith ierr = VecDestroy(x);CHKERRQ(ierr); 555064c8009SBarry Smith ierr = VecDestroy(y);CHKERRQ(ierr); 556064c8009SBarry Smith } 557064c8009SBarry Smith #endif 558064c8009SBarry Smith 559064c8009SBarry Smith ierr = MatDestroy(Aii);CHKERRQ(ierr); 560064c8009SBarry Smith ierr = MatDestroy(Ais);CHKERRQ(ierr); 561064c8009SBarry Smith ierr = MatDestroy(Asi);CHKERRQ(ierr); 562064c8009SBarry Smith ierr = MatDestroy(A);CHKERRQ(ierr); 563064c8009SBarry Smith ierr = ISDestroy(is);CHKERRQ(ierr); 564064c8009SBarry Smith ierr = ISDestroy(isint);CHKERRQ(ierr); 565064c8009SBarry Smith ierr = ISDestroy(issurf);CHKERRQ(ierr); 566064c8009SBarry Smith ierr = MatDestroy(Xint);CHKERRQ(ierr); 567064c8009SBarry Smith ierr = MatDestroy(Xsurf);CHKERRQ(ierr); 568064c8009SBarry Smith PetscFunctionReturn(0); 569064c8009SBarry Smith } 570174b6946SBarry Smith 5717233f9f0SBarry Smith 572174b6946SBarry Smith #undef __FUNCT__ 5737233f9f0SBarry Smith #define __FUNCT__ "PCExoticSetType" 5747233f9f0SBarry Smith /*@ 5757233f9f0SBarry Smith PCExoticSetType - Sets the type of coarse grid interpolation to use 5767233f9f0SBarry Smith 5777233f9f0SBarry Smith Collective on PC 5787233f9f0SBarry Smith 5797233f9f0SBarry Smith Input Parameters: 5807233f9f0SBarry Smith + pc - the preconditioner context 5817233f9f0SBarry Smith - type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face) 5827233f9f0SBarry Smith 583563e08c6SBarry Smith Notes: The face based interpolation has 1 degree of freedom per face and ignores the 584563e08c6SBarry Smith edge and vertex values completely in the coarse problem. For any seven point 585563e08c6SBarry Smith stencil the interpolation of a constant on all faces into the interior is that constant. 586563e08c6SBarry Smith 587563e08c6SBarry Smith The wirebasket interpolation has 1 degree of freedom per vertex, per edge and 588563e08c6SBarry Smith per face. A constant on the subdomain boundary is interpolated as that constant 589563e08c6SBarry Smith in the interior of the domain. 590563e08c6SBarry Smith 591563e08c6SBarry Smith The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence 592563e08c6SBarry Smith if A is nonsingular A_c is also nonsingular. 593563e08c6SBarry Smith 594563e08c6SBarry Smith Both interpolations are suitable for only scalar problems. 595563e08c6SBarry Smith 5967233f9f0SBarry Smith Level: intermediate 5977233f9f0SBarry Smith 5987233f9f0SBarry Smith 5997233f9f0SBarry Smith .seealso: PCEXOTIC, PCExoticType() 6007233f9f0SBarry Smith @*/ 6017233f9f0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCExoticSetType(PC pc,PCExoticType type) 6027233f9f0SBarry Smith { 6037233f9f0SBarry Smith PetscErrorCode ierr,(*f)(PC,PCExoticType); 6047233f9f0SBarry Smith 6057233f9f0SBarry Smith PetscFunctionBegin; 6067233f9f0SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6077233f9f0SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCExoticSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 6087233f9f0SBarry Smith if (f) { 6097233f9f0SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 6107233f9f0SBarry Smith } 6117233f9f0SBarry Smith PetscFunctionReturn(0); 6127233f9f0SBarry Smith } 6137233f9f0SBarry Smith 6147233f9f0SBarry Smith #undef __FUNCT__ 6157233f9f0SBarry Smith #define __FUNCT__ "PCExoticSetType_Exotic" 6167233f9f0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCExoticSetType_Exotic(PC pc,PCExoticType type) 6177233f9f0SBarry Smith { 6187233f9f0SBarry Smith PC_MG **mg = (PC_MG**)pc->data; 619bedee52eSLisandro Dalcin PC_Exotic *ctx = (PC_Exotic*) mg[0]->innerctx; 6207233f9f0SBarry Smith 6217233f9f0SBarry Smith PetscFunctionBegin; 6227233f9f0SBarry Smith ctx->type = type; 6237233f9f0SBarry Smith PetscFunctionReturn(0); 6247233f9f0SBarry Smith } 6257233f9f0SBarry Smith 6267233f9f0SBarry Smith #undef __FUNCT__ 6277233f9f0SBarry Smith #define __FUNCT__ "PCSetUp_Exotic" 6287233f9f0SBarry Smith PetscErrorCode PCSetUp_Exotic(PC pc) 629174b6946SBarry Smith { 630174b6946SBarry Smith PetscErrorCode ierr; 63196bdf778SBarry Smith Mat A; 6327233f9f0SBarry Smith PC_MG **mg = (PC_MG**)pc->data; 6337233f9f0SBarry Smith PC_Exotic *ex = (PC_Exotic*) mg[0]->innerctx; 63496bdf778SBarry Smith MatReuse reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX; 635174b6946SBarry Smith 636174b6946SBarry Smith PetscFunctionBegin; 637174b6946SBarry Smith ierr = PCGetOperators(pc,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr); 6387233f9f0SBarry Smith if (ex->type == PC_EXOTIC_FACE) { 639*8e722e36SBarry Smith ierr = DAGetFaceInterpolation(ex,A,reuse,&ex->P);CHKERRQ(ierr); 6407233f9f0SBarry Smith } else if (ex->type == PC_EXOTIC_WIREBASKET) { 641*8e722e36SBarry Smith ierr = DAGetWireBasketInterpolation(ex,A,reuse,&ex->P);CHKERRQ(ierr); 6427233f9f0SBarry Smith } else SETERRQ1(PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type); 64396bdf778SBarry Smith ierr = PCMGSetInterpolation(pc,1,ex->P);CHKERRQ(ierr); 6447233f9f0SBarry Smith ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 645174b6946SBarry Smith PetscFunctionReturn(0); 646174b6946SBarry Smith } 647174b6946SBarry Smith 648174b6946SBarry Smith #undef __FUNCT__ 6497233f9f0SBarry Smith #define __FUNCT__ "PCDestroy_Exotic" 6507233f9f0SBarry Smith PetscErrorCode PCDestroy_Exotic(PC pc) 651174b6946SBarry Smith { 652174b6946SBarry Smith PetscErrorCode ierr; 6537233f9f0SBarry Smith PC_MG **mg = (PC_MG**)pc->data; 654bedee52eSLisandro Dalcin PC_Exotic *ctx = (PC_Exotic*) mg[0]->innerctx; 655174b6946SBarry Smith 656174b6946SBarry Smith PetscFunctionBegin; 65796bdf778SBarry Smith if (ctx->da) {ierr = DADestroy(ctx->da);CHKERRQ(ierr);} 65896bdf778SBarry Smith if (ctx->P) {ierr = MatDestroy(ctx->P);CHKERRQ(ierr);} 659*8e722e36SBarry Smith if (ctx->ksp) {ierr = KSPDestroy(ctx->ksp);CHKERRQ(ierr);} 6607233f9f0SBarry Smith ierr = PetscFree(ctx);CHKERRQ(ierr); 6617233f9f0SBarry Smith ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 662174b6946SBarry Smith PetscFunctionReturn(0); 663174b6946SBarry Smith } 664174b6946SBarry Smith 665174b6946SBarry Smith #undef __FUNCT__ 6667233f9f0SBarry Smith #define __FUNCT__ "PCSetUp_Exotic_Error" 6677233f9f0SBarry Smith PetscErrorCode PCSetUp_Exotic_Error(PC pc) 668174b6946SBarry Smith { 669174b6946SBarry Smith PetscFunctionBegin; 670f91d8e95SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You are using the Exotic preconditioner but never called PCSetDA()"); 671174b6946SBarry Smith PetscFunctionReturn(0); 672174b6946SBarry Smith } 673174b6946SBarry Smith 674174b6946SBarry Smith #undef __FUNCT__ 675f91d8e95SBarry Smith #define __FUNCT__ "PCSetDA" 6767233f9f0SBarry Smith /*@ 677f91d8e95SBarry Smith PCSetDA - Sets the DA that is to be used by the PCEXOTIC or certain other preconditioners 6787233f9f0SBarry Smith 6797233f9f0SBarry Smith Collective on PC 6807233f9f0SBarry Smith 6817233f9f0SBarry Smith Input Parameters: 6827233f9f0SBarry Smith + pc - the preconditioner context 6837233f9f0SBarry Smith - da - the da 6847233f9f0SBarry Smith 6857233f9f0SBarry Smith Level: intermediate 6867233f9f0SBarry Smith 6877233f9f0SBarry Smith 6887233f9f0SBarry Smith .seealso: PCEXOTIC, PCExoticType() 6897233f9f0SBarry Smith @*/ 690f91d8e95SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSetDA(PC pc,DA da) 691174b6946SBarry Smith { 6927233f9f0SBarry Smith PetscErrorCode ierr,(*f)(PC,DA); 693174b6946SBarry Smith 694174b6946SBarry Smith PetscFunctionBegin; 695174b6946SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 696174b6946SBarry Smith PetscValidHeaderSpecific(da,DM_COOKIE,1); 697f91d8e95SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSetDA_C",(void (**)(void))&f);CHKERRQ(ierr); 6987233f9f0SBarry Smith if (f) { 6997233f9f0SBarry Smith ierr = (*f)(pc,da);CHKERRQ(ierr); 7007233f9f0SBarry Smith } 7017233f9f0SBarry Smith PetscFunctionReturn(0); 7027233f9f0SBarry Smith } 703174b6946SBarry Smith 7047233f9f0SBarry Smith 7057233f9f0SBarry Smith #undef __FUNCT__ 706f91d8e95SBarry Smith #define __FUNCT__ "PCSetDA_Exotic" 707f91d8e95SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSetDA_Exotic(PC pc,DA da) 7087233f9f0SBarry Smith { 7097233f9f0SBarry Smith PetscErrorCode ierr; 7107233f9f0SBarry Smith PC_MG **mg = (PC_MG**)pc->data; 711bedee52eSLisandro Dalcin PC_Exotic *ctx = (PC_Exotic*) mg[0]->innerctx; 7127233f9f0SBarry Smith 7137233f9f0SBarry Smith PetscFunctionBegin; 7147233f9f0SBarry Smith ctx->da = da; 7157233f9f0SBarry Smith pc->ops->setup = PCSetUp_Exotic; 716174b6946SBarry Smith ierr = PetscObjectReference((PetscObject)da);CHKERRQ(ierr); 717174b6946SBarry Smith PetscFunctionReturn(0); 718174b6946SBarry Smith } 719174b6946SBarry Smith 7207233f9f0SBarry Smith #undef __FUNCT__ 7217233f9f0SBarry Smith #define __FUNCT__ "PCView_Exotic" 7227233f9f0SBarry Smith PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer) 7237233f9f0SBarry Smith { 7247233f9f0SBarry Smith PC_MG **mg = (PC_MG**)pc->data; 7257233f9f0SBarry Smith PetscErrorCode ierr; 7267233f9f0SBarry Smith PetscTruth iascii; 727bedee52eSLisandro Dalcin PC_Exotic *ctx = (PC_Exotic*) mg[0]->innerctx; 7287233f9f0SBarry Smith 7297233f9f0SBarry Smith PetscFunctionBegin; 7307233f9f0SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 7317233f9f0SBarry Smith if (iascii) { 7327233f9f0SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Exotic type = %s\n",PCExoticTypes[ctx->type]);CHKERRQ(ierr); 733*8e722e36SBarry Smith if (ctx->directSolve) { 734*8e722e36SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using direct solver to construct interpolation\n");CHKERRQ(ierr); 735*8e722e36SBarry Smith } else { 736*8e722e36SBarry Smith PetscViewer sviewer; 737*8e722e36SBarry Smith PetscMPIInt rank; 738*8e722e36SBarry Smith 739*8e722e36SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using iterative solver to construct interpolation\n");CHKERRQ(ierr); 740*8e722e36SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 741*8e722e36SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* should not need to push this twice? */ 742*8e722e36SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 743*8e722e36SBarry Smith ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 744*8e722e36SBarry Smith if (!rank) { 745*8e722e36SBarry Smith ierr = KSPView(ctx->ksp,sviewer);CHKERRQ(ierr); 746*8e722e36SBarry Smith } 747*8e722e36SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 748*8e722e36SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 749*8e722e36SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 750*8e722e36SBarry Smith } 7517233f9f0SBarry Smith } 7527233f9f0SBarry Smith ierr = PCView_MG(pc,viewer);CHKERRQ(ierr); 7537233f9f0SBarry Smith PetscFunctionReturn(0); 7547233f9f0SBarry Smith } 7557233f9f0SBarry Smith 7567233f9f0SBarry Smith #undef __FUNCT__ 7577233f9f0SBarry Smith #define __FUNCT__ "PCSetFromOptions_Exotic" 7587233f9f0SBarry Smith PetscErrorCode PCSetFromOptions_Exotic(PC pc) 7597233f9f0SBarry Smith { 7607233f9f0SBarry Smith PetscErrorCode ierr; 7617233f9f0SBarry Smith PetscTruth flg; 7627233f9f0SBarry Smith PC_MG **mg = (PC_MG**)pc->data; 7637233f9f0SBarry Smith PCExoticType mgctype; 764bedee52eSLisandro Dalcin PC_Exotic *ctx = (PC_Exotic*) mg[0]->innerctx; 7657233f9f0SBarry Smith 7667233f9f0SBarry Smith PetscFunctionBegin; 7677233f9f0SBarry Smith ierr = PetscOptionsHead("Exotic coarse space options");CHKERRQ(ierr); 7687233f9f0SBarry Smith ierr = PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 7697233f9f0SBarry Smith if (flg) { 7707233f9f0SBarry Smith ierr = PCExoticSetType(pc,mgctype);CHKERRQ(ierr); 7717233f9f0SBarry Smith } 772*8e722e36SBarry Smith ierr = PetscOptionsTruth("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,PETSC_NULL);CHKERRQ(ierr); 773*8e722e36SBarry Smith if (!ctx->directSolve) { 774*8e722e36SBarry Smith if (!ctx->ksp) { 775*8e722e36SBarry Smith const char *prefix; 776*8e722e36SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&ctx->ksp);CHKERRQ(ierr); 777*8e722e36SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 778*8e722e36SBarry Smith ierr = PetscLogObjectParent(pc,ctx->ksp);CHKERRQ(ierr); 779*8e722e36SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 780*8e722e36SBarry Smith ierr = KSPSetOptionsPrefix(ctx->ksp,prefix);CHKERRQ(ierr); 781*8e722e36SBarry Smith ierr = KSPAppendOptionsPrefix(ctx->ksp,"exotic_");CHKERRQ(ierr); 782*8e722e36SBarry Smith } 783*8e722e36SBarry Smith ierr = KSPSetFromOptions(ctx->ksp);CHKERRQ(ierr); 784*8e722e36SBarry Smith } 7857233f9f0SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7867233f9f0SBarry Smith PetscFunctionReturn(0); 7877233f9f0SBarry Smith } 7887233f9f0SBarry Smith 789174b6946SBarry Smith 790174b6946SBarry Smith /*MC 7917233f9f0SBarry Smith PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces 792174b6946SBarry Smith 7937233f9f0SBarry Smith This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse 7943b65e785SBarry Smith grid spaces. These coarse grid spaces originate in the work of Bramble, Pasciak and Schatz, "The Construction 7953b65e785SBarry Smith of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53 pages 1--24, 1989. 7963b65e785SBarry Smith They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith, 7973b65e785SBarry Smith New York University, 1990. They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis 7983b65e785SBarry Smith of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical 7993b65e785SBarry Smith Analysis, volume 31. pages 1662-1694, 1994. These were developed in the context of iterative substructuring preconditioners. 8003b65e785SBarry Smith They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund. 8013b65e785SBarry Smith They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example, 8023b65e785SBarry Smith Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco 8033b65e785SBarry Smith Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings 8043b65e785SBarry Smith of the 17th International Conference on Domain Decomposition Methods in 8053b65e785SBarry Smith Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in 8063b65e785SBarry Smith Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 255-261, 2007. 8073b65e785SBarry Smith Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy min- 8083b65e785SBarry Smith imizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer, 8093b65e785SBarry Smith Marco Discacciati, David Keyes, OlofWidlund, andWalter Zulehner, editors, Proceedings 8103b65e785SBarry Smith of the 17th International Conference on Domain Decomposition Methods 8113b65e785SBarry Smith in Science and Engineering, held in Strobl, Austria, July 3-7, 2006, number 60 in 8123b65e785SBarry Smith Springer-Verlag, Lecture Notes in Computational Science and Engineering, pages 247-254, 2007 8133b65e785SBarry Smith Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition 8143b65e785SBarry Smith for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J. 8153b65e785SBarry Smith Numer. Anal., 46(4):2153-2168, 2008. 8163b65e785SBarry Smith Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz 8173b65e785SBarry Smith algorithm for almost incompressible elasticity. Technical Report 8183b65e785SBarry Smith TR2008-912, Department of Computer Science, Courant Institute 8193b65e785SBarry Smith of Mathematical Sciences, New York University, May 2008. URL: 8203b65e785SBarry Smith http://cs.nyu.edu/csweb/Research/TechReports/TR2008-912/TR2008-912.pdf 8217233f9f0SBarry Smith 8227233f9f0SBarry Smith Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type> 8237233f9f0SBarry Smith -pc_mg_type <type> 8247233f9f0SBarry Smith 82525a35f6fSSatish Balay Level: advanced 826174b6946SBarry Smith 827f91d8e95SBarry Smith .seealso: PCMG, PCSetDA(), PCExoticType, PCExoticSetType() 828174b6946SBarry Smith M*/ 829174b6946SBarry Smith 830174b6946SBarry Smith EXTERN_C_BEGIN 831174b6946SBarry Smith #undef __FUNCT__ 8327233f9f0SBarry Smith #define __FUNCT__ "PCCreate_Exotic" 8337233f9f0SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Exotic(PC pc) 834174b6946SBarry Smith { 835174b6946SBarry Smith PetscErrorCode ierr; 8367233f9f0SBarry Smith PC_Exotic *ex; 8377233f9f0SBarry Smith PC_MG **mg; 838174b6946SBarry Smith 839174b6946SBarry Smith PetscFunctionBegin; 840f91d8e95SBarry Smith /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 841f91d8e95SBarry Smith if (pc->ops->destroy) { ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); pc->data = 0;} 842f91d8e95SBarry Smith ierr = PetscStrfree(((PetscObject)pc)->type_name);CHKERRQ(ierr); 843f91d8e95SBarry Smith ((PetscObject)pc)->type_name = 0; 844f91d8e95SBarry Smith 845174b6946SBarry Smith ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); 846174b6946SBarry Smith ierr = PCMGSetLevels(pc,2,PETSC_NULL);CHKERRQ(ierr); 847174b6946SBarry Smith ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 84896bdf778SBarry Smith ierr = PetscNew(PC_Exotic,&ex);CHKERRQ(ierr);\ 8497233f9f0SBarry Smith ex->type = PC_EXOTIC_FACE; 850a36ce94aSMatthew Knepley mg = (PC_MG**) pc->data; 8517233f9f0SBarry Smith mg[0]->innerctx = ex; 8527233f9f0SBarry Smith 8537233f9f0SBarry Smith 8547233f9f0SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Exotic; 8557233f9f0SBarry Smith pc->ops->view = PCView_Exotic; 8567233f9f0SBarry Smith pc->ops->destroy = PCDestroy_Exotic; 8577233f9f0SBarry Smith pc->ops->setup = PCSetUp_Exotic_Error; 8587233f9f0SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCExoticSetType_C","PCExoticSetType_Exotic",PCExoticSetType_Exotic);CHKERRQ(ierr); 859f91d8e95SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSetDA_C","PCSetDA_Exotic",PCSetDA_Exotic);CHKERRQ(ierr); 860174b6946SBarry Smith PetscFunctionReturn(0); 861174b6946SBarry Smith } 862174b6946SBarry Smith EXTERN_C_END 863