1174b6946SBarry Smith 2c6db04a5SJed Brown #include <petscdmda.h> /*I "petscdmda.h" I*/ 3af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 482c86c8fSBarry Smith #include <petscctable.h> 57233f9f0SBarry Smith 68e722e36SBarry Smith typedef struct { 78e722e36SBarry Smith PCExoticType type; 88e722e36SBarry Smith Mat P; /* the constructed interpolation matrix */ 9ace3abfcSBarry Smith PetscBool directSolve; /* use direct LU factorization to construct interpolation */ 108e722e36SBarry Smith KSP ksp; 118e722e36SBarry Smith } PC_Exotic; 12174b6946SBarry Smith 130a545947SLisandro Dalcin const char *const PCExoticTypes[] = {"face","wirebasket","PCExoticType","PC_Exotic",NULL}; 14064c8009SBarry Smith 15064c8009SBarry Smith /* 16aa219208SBarry Smith DMDAGetWireBasketInterpolation - Gets the interpolation for a wirebasket based coarse space 17064c8009SBarry Smith 18064c8009SBarry Smith */ 19c0decd05SBarry Smith PetscErrorCode DMDAGetWireBasketInterpolation(PC pc,DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P) 20064c8009SBarry Smith { 21064c8009SBarry Smith PetscErrorCode ierr; 22064c8009SBarry 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; 2328d20b34SBarry Smith PetscInt mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[26],*globals,Ng,*IIint,*IIsurf,Nt; 24064c8009SBarry Smith Mat Xint, Xsurf,Xint_tmp; 25064c8009SBarry Smith IS isint,issurf,is,row,col; 26064c8009SBarry Smith ISLocalToGlobalMapping ltg; 27064c8009SBarry Smith MPI_Comm comm; 28064c8009SBarry Smith Mat A,Aii,Ais,Asi,*Aholder,iAii; 29064c8009SBarry Smith MatFactorInfo info; 30064c8009SBarry Smith PetscScalar *xsurf,*xint; 311683a169SBarry Smith const PetscScalar *rxint; 328e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 33064c8009SBarry Smith PetscScalar tmp; 34064c8009SBarry Smith #endif 35064c8009SBarry Smith PetscTable ht; 36064c8009SBarry Smith 37064c8009SBarry Smith PetscFunctionBegin; 380a545947SLisandro Dalcin ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,&mp,&np,&pp,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 39ce94432eSBarry Smith if (dof != 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only for single field problems"); 40ce94432eSBarry Smith if (dim != 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only coded for 3d problems"); 410a545947SLisandro Dalcin ierr = DMDAGetCorners(da,NULL,NULL,NULL,&m,&n,&p);CHKERRQ(ierr); 42aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr); 43064c8009SBarry Smith istart = istart ? -1 : 0; 44064c8009SBarry Smith jstart = jstart ? -1 : 0; 45064c8009SBarry Smith kstart = kstart ? -1 : 0; 46064c8009SBarry Smith 47064c8009SBarry Smith /* 48064c8009SBarry Smith the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 49064c8009SBarry Smith to all the local degrees of freedom (this includes the vertices, edges and faces). 50064c8009SBarry Smith 51064c8009SBarry Smith Xint are the subset of the interpolation into the interior 52064c8009SBarry Smith 53064c8009SBarry Smith Xface are the interpolation onto faces but not into the interior 54064c8009SBarry Smith 55064c8009SBarry Smith Xsurf are the interpolation onto the vertices and edges (the surfbasket) 56064c8009SBarry Smith Xint 57064c8009SBarry Smith Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain 58064c8009SBarry Smith Xsurf 59064c8009SBarry Smith */ 60064c8009SBarry Smith N = (m - istart)*(n - jstart)*(p - kstart); 61064c8009SBarry Smith Nint = (m-2-istart)*(n-2-jstart)*(p-2-kstart); 62064c8009SBarry Smith Nface = 2*((m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart)); 63064c8009SBarry Smith Nwire = 4*((m-2-istart) + (n-2-jstart) + (p-2-kstart)) + 8; 64064c8009SBarry Smith Nsurf = Nface + Nwire; 650298fd71SBarry Smith ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,26,NULL,&Xint);CHKERRQ(ierr); 660298fd71SBarry Smith ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,26,NULL,&Xsurf);CHKERRQ(ierr); 678c778c55SBarry Smith ierr = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 68064c8009SBarry Smith 69064c8009SBarry Smith /* 70064c8009SBarry Smith Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 71064c8009SBarry Smith Xsurf will be all zero (thus making the coarse matrix singular). 72064c8009SBarry Smith */ 73e32f2f54SBarry Smith 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"); 74e32f2f54SBarry Smith 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"); 75e32f2f54SBarry Smith 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"); 76064c8009SBarry Smith 77064c8009SBarry Smith cnt = 0; 782fa5cd67SKarl Rupp 792fa5cd67SKarl Rupp xsurf[cnt++] = 1; 802fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + Nsurf] = 1; 812fa5cd67SKarl Rupp xsurf[cnt++ + 2*Nsurf] = 1; 822fa5cd67SKarl Rupp 832fa5cd67SKarl Rupp for (j=1; j<n-1-jstart; j++) { 842fa5cd67SKarl Rupp xsurf[cnt++ + 3*Nsurf] = 1; 852fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 4*Nsurf] = 1; 862fa5cd67SKarl Rupp xsurf[cnt++ + 5*Nsurf] = 1; 87064c8009SBarry Smith } 882fa5cd67SKarl Rupp 892fa5cd67SKarl Rupp xsurf[cnt++ + 6*Nsurf] = 1; 902fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 7*Nsurf] = 1; 912fa5cd67SKarl Rupp xsurf[cnt++ + 8*Nsurf] = 1; 922fa5cd67SKarl Rupp 932fa5cd67SKarl Rupp for (k=1; k<p-1-kstart; k++) { 942fa5cd67SKarl Rupp xsurf[cnt++ + 9*Nsurf] = 1; 952fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 10*Nsurf] = 1; 962fa5cd67SKarl Rupp xsurf[cnt++ + 11*Nsurf] = 1; 972fa5cd67SKarl Rupp 982fa5cd67SKarl Rupp for (j=1; j<n-1-jstart; j++) { 992fa5cd67SKarl Rupp xsurf[cnt++ + 12*Nsurf] = 1; 1002fa5cd67SKarl Rupp /* these are the interior nodes */ 1012fa5cd67SKarl Rupp xsurf[cnt++ + 13*Nsurf] = 1; 1022fa5cd67SKarl Rupp } 1032fa5cd67SKarl Rupp 1042fa5cd67SKarl Rupp xsurf[cnt++ + 14*Nsurf] = 1; 1052fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 15*Nsurf] = 1; 1062fa5cd67SKarl Rupp xsurf[cnt++ + 16*Nsurf] = 1; 1072fa5cd67SKarl Rupp } 1082fa5cd67SKarl Rupp 1092fa5cd67SKarl Rupp xsurf[cnt++ + 17*Nsurf] = 1; 1102fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 18*Nsurf] = 1; 1112fa5cd67SKarl Rupp xsurf[cnt++ + 19*Nsurf] = 1; 1122fa5cd67SKarl Rupp 1132fa5cd67SKarl Rupp for (j=1;j<n-1-jstart;j++) { 1142fa5cd67SKarl Rupp xsurf[cnt++ + 20*Nsurf] = 1; 1152fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 21*Nsurf] = 1; 1162fa5cd67SKarl Rupp xsurf[cnt++ + 22*Nsurf] = 1; 1172fa5cd67SKarl Rupp } 1182fa5cd67SKarl Rupp 1192fa5cd67SKarl Rupp xsurf[cnt++ + 23*Nsurf] = 1; 1202fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 24*Nsurf] = 1; 1212fa5cd67SKarl Rupp xsurf[cnt++ + 25*Nsurf] = 1; 1222fa5cd67SKarl Rupp 1238e722e36SBarry Smith /* interpolations only sum to 1 when using direct solver */ 1248e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 125064c8009SBarry Smith for (i=0; i<Nsurf; i++) { 126064c8009SBarry Smith tmp = 0.0; 1272fa5cd67SKarl Rupp for (j=0; j<26; j++) tmp += xsurf[i+j*Nsurf]; 128*98921bdaSJacob Faibussowitsch if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %g",i,(double)PetscAbsScalar(tmp)); 129064c8009SBarry Smith } 130064c8009SBarry Smith #endif 1318c778c55SBarry Smith ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 132064c8009SBarry Smith /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/ 133064c8009SBarry Smith 134064c8009SBarry Smith /* 135064c8009SBarry Smith I are the indices for all the needed vertices (in global numbering) 136064c8009SBarry Smith Iint are the indices for the interior values, I surf for the surface values 1377dae84e0SHong Zhang (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it 138aa219208SBarry Smith is NOT the local DMDA ordering.) 139064c8009SBarry Smith IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering 140064c8009SBarry Smith */ 141064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start)) 142dcca6d9dSJed Brown ierr = PetscMalloc3(N,&II,Nint,&Iint,Nsurf,&Isurf);CHKERRQ(ierr); 143dcca6d9dSJed Brown ierr = PetscMalloc2(Nint,&IIint,Nsurf,&IIsurf);CHKERRQ(ierr); 144064c8009SBarry Smith for (k=0; k<p-kstart; k++) { 145064c8009SBarry Smith for (j=0; j<n-jstart; j++) { 146064c8009SBarry Smith for (i=0; i<m-istart; i++) { 147064c8009SBarry Smith II[c++] = i + j*mwidth + k*mwidth*nwidth; 148064c8009SBarry Smith 149064c8009SBarry Smith if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) { 150064c8009SBarry Smith IIint[cint] = i + j*mwidth + k*mwidth*nwidth; 151064c8009SBarry Smith Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 152064c8009SBarry Smith } else { 153064c8009SBarry Smith IIsurf[csurf] = i + j*mwidth + k*mwidth*nwidth; 154064c8009SBarry Smith Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 155064c8009SBarry Smith } 156064c8009SBarry Smith } 157064c8009SBarry Smith } 158064c8009SBarry Smith } 159e32f2f54SBarry Smith if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N"); 160e32f2f54SBarry Smith if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint"); 161e32f2f54SBarry Smith if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf"); 1621411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<g);CHKERRQ(ierr); 163064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr); 164064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr); 165064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr); 166064c8009SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 16770b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 16870b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr); 16970b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr); 170064c8009SBarry Smith ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr); 171064c8009SBarry Smith 1727dae84e0SHong Zhang ierr = MatCreateSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr); 173064c8009SBarry Smith A = *Aholder; 174064c8009SBarry Smith ierr = PetscFree(Aholder);CHKERRQ(ierr); 175064c8009SBarry Smith 1767dae84e0SHong Zhang ierr = MatCreateSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr); 1777dae84e0SHong Zhang ierr = MatCreateSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr); 1787dae84e0SHong Zhang ierr = MatCreateSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr); 179064c8009SBarry Smith 180064c8009SBarry Smith /* 181064c8009SBarry Smith Solve for the interpolation onto the interior Xint 182064c8009SBarry Smith */ 18392e4c2edSBarry Smith ierr = MatMatMult(Ais,Xsurf,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr); 1848e722e36SBarry Smith ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr); 1858e722e36SBarry Smith if (exotic->directSolve) { 1862692d6eeSBarry Smith ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr); 187064c8009SBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 1882692d6eeSBarry Smith ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr); 189064c8009SBarry Smith ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr); 190fcfd50ebSBarry Smith ierr = ISDestroy(&row);CHKERRQ(ierr); 191fcfd50ebSBarry Smith ierr = ISDestroy(&col);CHKERRQ(ierr); 192064c8009SBarry Smith ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr); 193064c8009SBarry Smith ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr); 1946bf464f9SBarry Smith ierr = MatDestroy(&iAii);CHKERRQ(ierr); 1958e722e36SBarry Smith } else { 1968e722e36SBarry Smith Vec b,x; 1978e722e36SBarry Smith PetscScalar *xint_tmp; 198064c8009SBarry Smith 1998c778c55SBarry Smith ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr); 2000a545947SLisandro Dalcin ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,NULL,&x);CHKERRQ(ierr); 2018c778c55SBarry Smith ierr = MatDenseGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 2020a545947SLisandro Dalcin ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,NULL,&b);CHKERRQ(ierr); 20323ee1639SBarry Smith ierr = KSPSetOperators(exotic->ksp,Aii,Aii);CHKERRQ(ierr); 2048e722e36SBarry Smith for (i=0; i<26; i++) { 2058e722e36SBarry Smith ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr); 2068e722e36SBarry Smith ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr); 2078e722e36SBarry Smith ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr); 208c0decd05SBarry Smith ierr = KSPCheckSolve(exotic->ksp,pc,x);CHKERRQ(ierr); 2098e722e36SBarry Smith ierr = VecResetArray(x);CHKERRQ(ierr); 2108e722e36SBarry Smith ierr = VecResetArray(b);CHKERRQ(ierr); 2118e722e36SBarry Smith } 2128c778c55SBarry Smith ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr); 2138c778c55SBarry Smith ierr = MatDenseRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 2146bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 2156bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 2168e722e36SBarry Smith } 2176bf464f9SBarry Smith ierr = MatDestroy(&Xint_tmp);CHKERRQ(ierr); 2188e722e36SBarry Smith 2198e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 2201683a169SBarry Smith ierr = MatDenseGetArrayRead(Xint,&rxint);CHKERRQ(ierr); 221064c8009SBarry Smith for (i=0; i<Nint; i++) { 222064c8009SBarry Smith tmp = 0.0; 2231683a169SBarry Smith for (j=0; j<26; j++) tmp += rxint[i+j*Nint]; 2242fa5cd67SKarl Rupp 225*98921bdaSJacob Faibussowitsch if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %g",i,(double)PetscAbsScalar(tmp)); 226064c8009SBarry Smith } 2271683a169SBarry Smith ierr = MatDenseRestoreArrayRead(Xint,&rxint);CHKERRQ(ierr); 228064c8009SBarry Smith /* ierr = MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 229064c8009SBarry Smith #endif 230064c8009SBarry Smith 231064c8009SBarry Smith /* total vertices total faces total edges */ 232064c8009SBarry 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); 233064c8009SBarry Smith 234064c8009SBarry Smith /* 235064c8009SBarry Smith For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 236064c8009SBarry Smith */ 237064c8009SBarry Smith cnt = 0; 2382fa5cd67SKarl Rupp 239064c8009SBarry Smith gl[cnt++] = 0; { gl[cnt++] = 1;} gl[cnt++] = m-istart-1; 240064c8009SBarry Smith { gl[cnt++] = mwidth; { gl[cnt++] = mwidth+1;} gl[cnt++] = mwidth + m-istart-1;} 241064c8009SBarry Smith gl[cnt++] = mwidth*(n-jstart-1); { gl[cnt++] = mwidth*(n-jstart-1)+1;} gl[cnt++] = mwidth*(n-jstart-1) + m-istart-1; 242064c8009SBarry Smith { 243064c8009SBarry Smith gl[cnt++] = mwidth*nwidth; { gl[cnt++] = mwidth*nwidth+1;} gl[cnt++] = mwidth*nwidth+ m-istart-1; 244064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;} 245064c8009SBarry 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; 246064c8009SBarry Smith } 247064c8009SBarry 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; 248064c8009SBarry 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;} 249064c8009SBarry 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; 250064c8009SBarry Smith 251064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 252064c8009SBarry Smith /* convert that to global numbering and get them on all processes */ 253064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,26,gl,gl);CHKERRQ(ierr); 254064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 255785e854fSJed Brown ierr = PetscMalloc1(26*mp*np*pp,&globals);CHKERRQ(ierr); 256ffc4695bSBarry Smith ierr = MPI_Allgather(gl,26,MPIU_INT,globals,26,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr); 257064c8009SBarry Smith 258064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */ 2590298fd71SBarry Smith ierr = MatGetSize(Aglobal,&Nt,NULL);CHKERRQ(ierr); 26012dd4568SBarry Smith ierr = PetscTableCreate(Ntotal/3,Nt+1,&ht);CHKERRQ(ierr); 261064c8009SBarry Smith for (i=0; i<26*mp*np*pp; i++) { 262064c8009SBarry Smith ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr); 263064c8009SBarry Smith } 264064c8009SBarry Smith ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr); 265*98921bdaSJacob Faibussowitsch if (cnt != Ntotal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal); 266064c8009SBarry Smith ierr = PetscFree(globals);CHKERRQ(ierr); 267064c8009SBarry Smith for (i=0; i<26; i++) { 268064c8009SBarry Smith ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr); 269064c8009SBarry Smith gl[i]--; 270064c8009SBarry Smith } 2716bc0bbbfSBarry Smith ierr = PetscTableDestroy(&ht);CHKERRQ(ierr); 272064c8009SBarry Smith /* PetscIntView(26,gl,PETSC_VIEWER_STDOUT_WORLD); */ 273064c8009SBarry Smith 274064c8009SBarry Smith /* construct global interpolation matrix */ 2750298fd71SBarry Smith ierr = MatGetLocalSize(Aglobal,&Ng,NULL);CHKERRQ(ierr); 276064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 277ce94432eSBarry Smith ierr = MatCreateAIJ(PetscObjectComm((PetscObject)da),Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,NULL,Nint+Nsurf,NULL,P);CHKERRQ(ierr); 278064c8009SBarry Smith } else { 279064c8009SBarry Smith ierr = MatZeroEntries(*P);CHKERRQ(ierr); 280064c8009SBarry Smith } 281064c8009SBarry Smith ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 2821683a169SBarry Smith ierr = MatDenseGetArrayRead(Xint,&rxint);CHKERRQ(ierr); 2831683a169SBarry Smith ierr = MatSetValues(*P,Nint,IIint,26,gl,rxint,INSERT_VALUES);CHKERRQ(ierr); 2841683a169SBarry Smith ierr = MatDenseRestoreArrayRead(Xint,&rxint);CHKERRQ(ierr); 2851683a169SBarry Smith ierr = MatDenseGetArrayRead(Xsurf,&rxint);CHKERRQ(ierr); 2861683a169SBarry Smith ierr = MatSetValues(*P,Nsurf,IIsurf,26,gl,rxint,INSERT_VALUES);CHKERRQ(ierr); 2871683a169SBarry Smith ierr = MatDenseRestoreArrayRead(Xsurf,&rxint);CHKERRQ(ierr); 288064c8009SBarry Smith ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 289064c8009SBarry Smith ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 290064c8009SBarry Smith ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr); 291064c8009SBarry Smith 2928e722e36SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 293064c8009SBarry Smith { 294064c8009SBarry Smith Vec x,y; 295064c8009SBarry Smith PetscScalar *yy; 296ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr); 297ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr); 298064c8009SBarry Smith ierr = VecSet(x,1.0);CHKERRQ(ierr); 299064c8009SBarry Smith ierr = MatMult(*P,x,y);CHKERRQ(ierr); 300064c8009SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 301064c8009SBarry Smith for (i=0; i<Ng; i++) { 302*98921bdaSJacob Faibussowitsch if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %g",i,(double)PetscAbsScalar(yy[i])); 303064c8009SBarry Smith } 304064c8009SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 305064c8009SBarry Smith ierr = VecDestroy(x);CHKERRQ(ierr); 306064c8009SBarry Smith ierr = VecDestroy(y);CHKERRQ(ierr); 307064c8009SBarry Smith } 308064c8009SBarry Smith #endif 309064c8009SBarry Smith 310fcfd50ebSBarry Smith ierr = MatDestroy(&Aii);CHKERRQ(ierr); 311fcfd50ebSBarry Smith ierr = MatDestroy(&Ais);CHKERRQ(ierr); 312fcfd50ebSBarry Smith ierr = MatDestroy(&Asi);CHKERRQ(ierr); 313fcfd50ebSBarry Smith ierr = MatDestroy(&A);CHKERRQ(ierr); 314fcfd50ebSBarry Smith ierr = ISDestroy(&is);CHKERRQ(ierr); 315fcfd50ebSBarry Smith ierr = ISDestroy(&isint);CHKERRQ(ierr); 316fcfd50ebSBarry Smith ierr = ISDestroy(&issurf);CHKERRQ(ierr); 317fcfd50ebSBarry Smith ierr = MatDestroy(&Xint);CHKERRQ(ierr); 318fcfd50ebSBarry Smith ierr = MatDestroy(&Xsurf);CHKERRQ(ierr); 319064c8009SBarry Smith PetscFunctionReturn(0); 320064c8009SBarry Smith } 321064c8009SBarry Smith 322064c8009SBarry Smith /* 323aa219208SBarry Smith DMDAGetFaceInterpolation - Gets the interpolation for a face based coarse space 324064c8009SBarry Smith 325064c8009SBarry Smith */ 326c0decd05SBarry Smith PetscErrorCode DMDAGetFaceInterpolation(PC pc,DM da,PC_Exotic *exotic,Mat Aglobal,MatReuse reuse,Mat *P) 327064c8009SBarry Smith { 328064c8009SBarry Smith PetscErrorCode ierr; 329064c8009SBarry 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; 33028d20b34SBarry Smith PetscInt mwidth,nwidth,pwidth,cnt,mp,np,pp,Ntotal,gl[6],*globals,Ng,*IIint,*IIsurf,Nt; 331064c8009SBarry Smith Mat Xint, Xsurf,Xint_tmp; 332064c8009SBarry Smith IS isint,issurf,is,row,col; 333064c8009SBarry Smith ISLocalToGlobalMapping ltg; 334064c8009SBarry Smith MPI_Comm comm; 335064c8009SBarry Smith Mat A,Aii,Ais,Asi,*Aholder,iAii; 336064c8009SBarry Smith MatFactorInfo info; 337064c8009SBarry Smith PetscScalar *xsurf,*xint; 3381683a169SBarry Smith const PetscScalar *rxint; 339064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 340064c8009SBarry Smith PetscScalar tmp; 341064c8009SBarry Smith #endif 342064c8009SBarry Smith PetscTable ht; 343064c8009SBarry Smith 344064c8009SBarry Smith PetscFunctionBegin; 3450a545947SLisandro Dalcin ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,&mp,&np,&pp,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 346ce94432eSBarry Smith if (dof != 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only for single field problems"); 347ce94432eSBarry Smith if (dim != 3) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only coded for 3d problems"); 3480a545947SLisandro Dalcin ierr = DMDAGetCorners(da,NULL,NULL,NULL,&m,&n,&p);CHKERRQ(ierr); 349aa219208SBarry Smith ierr = DMDAGetGhostCorners(da,&istart,&jstart,&kstart,&mwidth,&nwidth,&pwidth);CHKERRQ(ierr); 350064c8009SBarry Smith istart = istart ? -1 : 0; 351064c8009SBarry Smith jstart = jstart ? -1 : 0; 352064c8009SBarry Smith kstart = kstart ? -1 : 0; 353064c8009SBarry Smith 354064c8009SBarry Smith /* 355064c8009SBarry Smith the columns of P are the interpolation of each coarse grid point (one for each vertex and edge) 356064c8009SBarry Smith to all the local degrees of freedom (this includes the vertices, edges and faces). 357064c8009SBarry Smith 358064c8009SBarry Smith Xint are the subset of the interpolation into the interior 359064c8009SBarry Smith 360064c8009SBarry Smith Xface are the interpolation onto faces but not into the interior 361064c8009SBarry Smith 362064c8009SBarry Smith Xsurf are the interpolation onto the vertices and edges (the surfbasket) 363064c8009SBarry Smith Xint 364064c8009SBarry Smith Symbolically one could write P = (Xface) after interchanging the rows to match the natural ordering on the domain 365064c8009SBarry Smith Xsurf 366064c8009SBarry Smith */ 367064c8009SBarry Smith N = (m - istart)*(n - jstart)*(p - kstart); 368064c8009SBarry Smith Nint = (m-2-istart)*(n-2-jstart)*(p-2-kstart); 369064c8009SBarry Smith Nface = 2*((m-2-istart)*(n-2-jstart) + (m-2-istart)*(p-2-kstart) + (n-2-jstart)*(p-2-kstart)); 370064c8009SBarry Smith Nwire = 4*((m-2-istart) + (n-2-jstart) + (p-2-kstart)) + 8; 371064c8009SBarry Smith Nsurf = Nface + Nwire; 3720298fd71SBarry Smith ierr = MatCreateSeqDense(MPI_COMM_SELF,Nint,6,NULL,&Xint);CHKERRQ(ierr); 3730298fd71SBarry Smith ierr = MatCreateSeqDense(MPI_COMM_SELF,Nsurf,6,NULL,&Xsurf);CHKERRQ(ierr); 3748c778c55SBarry Smith ierr = MatDenseGetArray(Xsurf,&xsurf);CHKERRQ(ierr); 375064c8009SBarry Smith 376064c8009SBarry Smith /* 377064c8009SBarry Smith Require that all 12 edges and 6 faces have at least one grid point. Otherwise some of the columns of 378064c8009SBarry Smith Xsurf will be all zero (thus making the coarse matrix singular). 379064c8009SBarry Smith */ 380e32f2f54SBarry Smith 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"); 381e32f2f54SBarry Smith 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"); 382e32f2f54SBarry Smith 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"); 383064c8009SBarry Smith 384064c8009SBarry Smith cnt = 0; 3852fa5cd67SKarl Rupp for (j=1; j<n-1-jstart; j++) { 3862fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 0*Nsurf] = 1; 387064c8009SBarry Smith } 3882fa5cd67SKarl Rupp 3892fa5cd67SKarl Rupp for (k=1; k<p-1-kstart; k++) { 3902fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 1*Nsurf] = 1; 3912fa5cd67SKarl Rupp for (j=1; j<n-1-jstart; j++) { 3922fa5cd67SKarl Rupp xsurf[cnt++ + 2*Nsurf] = 1; 3932fa5cd67SKarl Rupp /* these are the interior nodes */ 3942fa5cd67SKarl Rupp xsurf[cnt++ + 3*Nsurf] = 1; 3952fa5cd67SKarl Rupp } 3962fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 4*Nsurf] = 1; 3972fa5cd67SKarl Rupp } 3982fa5cd67SKarl Rupp for (j=1;j<n-1-jstart;j++) { 3992fa5cd67SKarl Rupp for (i=1; i<m-istart-1; i++) xsurf[cnt++ + 5*Nsurf] = 1; 4002fa5cd67SKarl Rupp } 401064c8009SBarry Smith 402064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 403064c8009SBarry Smith for (i=0; i<Nsurf; i++) { 404064c8009SBarry Smith tmp = 0.0; 4052fa5cd67SKarl Rupp for (j=0; j<6; j++) tmp += xsurf[i+j*Nsurf]; 4062fa5cd67SKarl Rupp 407*98921bdaSJacob Faibussowitsch if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xsurf interpolation at i %D value %g",i,(double)PetscAbsScalar(tmp)); 408064c8009SBarry Smith } 409064c8009SBarry Smith #endif 4108c778c55SBarry Smith ierr = MatDenseRestoreArray(Xsurf,&xsurf);CHKERRQ(ierr); 411064c8009SBarry Smith /* ierr = MatView(Xsurf,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);*/ 412064c8009SBarry Smith 413064c8009SBarry Smith /* 414064c8009SBarry Smith I are the indices for all the needed vertices (in global numbering) 415064c8009SBarry Smith Iint are the indices for the interior values, I surf for the surface values 4167dae84e0SHong Zhang (This is just for the part of the global matrix obtained with MatCreateSubMatrix(), it 417aa219208SBarry Smith is NOT the local DMDA ordering.) 418064c8009SBarry Smith IIint and IIsurf are the same as the Iint, Isurf except they are in the global numbering 419064c8009SBarry Smith */ 420064c8009SBarry Smith #define Endpoint(a,start,b) (a == 0 || a == (b-1-start)) 421dcca6d9dSJed Brown ierr = PetscMalloc3(N,&II,Nint,&Iint,Nsurf,&Isurf);CHKERRQ(ierr); 422dcca6d9dSJed Brown ierr = PetscMalloc2(Nint,&IIint,Nsurf,&IIsurf);CHKERRQ(ierr); 423064c8009SBarry Smith for (k=0; k<p-kstart; k++) { 424064c8009SBarry Smith for (j=0; j<n-jstart; j++) { 425064c8009SBarry Smith for (i=0; i<m-istart; i++) { 426064c8009SBarry Smith II[c++] = i + j*mwidth + k*mwidth*nwidth; 427064c8009SBarry Smith 428064c8009SBarry Smith if (!Endpoint(i,istart,m) && !Endpoint(j,jstart,n) && !Endpoint(k,kstart,p)) { 429064c8009SBarry Smith IIint[cint] = i + j*mwidth + k*mwidth*nwidth; 430064c8009SBarry Smith Iint[cint++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 431064c8009SBarry Smith } else { 432064c8009SBarry Smith IIsurf[csurf] = i + j*mwidth + k*mwidth*nwidth; 433064c8009SBarry Smith Isurf[csurf++] = i + j*(m-istart) + k*(m-istart)*(n-jstart); 434064c8009SBarry Smith } 435064c8009SBarry Smith } 436064c8009SBarry Smith } 437064c8009SBarry Smith } 438e32f2f54SBarry Smith if (c != N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"c != N"); 439e32f2f54SBarry Smith if (cint != Nint) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"cint != Nint"); 440e32f2f54SBarry Smith if (csurf != Nsurf) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csurf != Nsurf"); 4411411c6eeSJed Brown ierr = DMGetLocalToGlobalMapping(da,<g);CHKERRQ(ierr); 442064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,N,II,II);CHKERRQ(ierr); 443064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,Nint,IIint,IIint);CHKERRQ(ierr); 444064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,Nsurf,IIsurf,IIsurf);CHKERRQ(ierr); 445064c8009SBarry Smith ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 44670b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,N,II,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 44770b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,Nint,Iint,PETSC_COPY_VALUES,&isint);CHKERRQ(ierr); 44870b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,Nsurf,Isurf,PETSC_COPY_VALUES,&issurf);CHKERRQ(ierr); 449064c8009SBarry Smith ierr = PetscFree3(II,Iint,Isurf);CHKERRQ(ierr); 450064c8009SBarry Smith 451671e8369SHong Zhang ierr = ISSort(is);CHKERRQ(ierr); 4527dae84e0SHong Zhang ierr = MatCreateSubMatrices(Aglobal,1,&is,&is,MAT_INITIAL_MATRIX,&Aholder);CHKERRQ(ierr); 453064c8009SBarry Smith A = *Aholder; 454064c8009SBarry Smith ierr = PetscFree(Aholder);CHKERRQ(ierr); 455064c8009SBarry Smith 4567dae84e0SHong Zhang ierr = MatCreateSubMatrix(A,isint,isint,MAT_INITIAL_MATRIX,&Aii);CHKERRQ(ierr); 4577dae84e0SHong Zhang ierr = MatCreateSubMatrix(A,isint,issurf,MAT_INITIAL_MATRIX,&Ais);CHKERRQ(ierr); 4587dae84e0SHong Zhang ierr = MatCreateSubMatrix(A,issurf,isint,MAT_INITIAL_MATRIX,&Asi);CHKERRQ(ierr); 459064c8009SBarry Smith 460064c8009SBarry Smith /* 461064c8009SBarry Smith Solve for the interpolation onto the interior Xint 462064c8009SBarry Smith */ 46392e4c2edSBarry Smith ierr = MatMatMult(Ais,Xsurf,MAT_INITIAL_MATRIX,PETSC_DETERMINE,&Xint_tmp);CHKERRQ(ierr); 464064c8009SBarry Smith ierr = MatScale(Xint_tmp,-1.0);CHKERRQ(ierr); 465064c8009SBarry Smith 4668e722e36SBarry Smith if (exotic->directSolve) { 4672692d6eeSBarry Smith ierr = MatGetFactor(Aii,MATSOLVERPETSC,MAT_FACTOR_LU,&iAii);CHKERRQ(ierr); 468064c8009SBarry Smith ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 4692692d6eeSBarry Smith ierr = MatGetOrdering(Aii,MATORDERINGND,&row,&col);CHKERRQ(ierr); 470064c8009SBarry Smith ierr = MatLUFactorSymbolic(iAii,Aii,row,col,&info);CHKERRQ(ierr); 471fcfd50ebSBarry Smith ierr = ISDestroy(&row);CHKERRQ(ierr); 472fcfd50ebSBarry Smith ierr = ISDestroy(&col);CHKERRQ(ierr); 473064c8009SBarry Smith ierr = MatLUFactorNumeric(iAii,Aii,&info);CHKERRQ(ierr); 474064c8009SBarry Smith ierr = MatMatSolve(iAii,Xint_tmp,Xint);CHKERRQ(ierr); 475fcfd50ebSBarry Smith ierr = MatDestroy(&iAii);CHKERRQ(ierr); 476064c8009SBarry Smith } else { 477064c8009SBarry Smith Vec b,x; 478064c8009SBarry Smith PetscScalar *xint_tmp; 479064c8009SBarry Smith 4808c778c55SBarry Smith ierr = MatDenseGetArray(Xint,&xint);CHKERRQ(ierr); 4810a545947SLisandro Dalcin ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,NULL,&x);CHKERRQ(ierr); 4828c778c55SBarry Smith ierr = MatDenseGetArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 4830a545947SLisandro Dalcin ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,1,Nint,NULL,&b);CHKERRQ(ierr); 48423ee1639SBarry Smith ierr = KSPSetOperators(exotic->ksp,Aii,Aii);CHKERRQ(ierr); 485064c8009SBarry Smith for (i=0; i<6; i++) { 486064c8009SBarry Smith ierr = VecPlaceArray(x,xint+i*Nint);CHKERRQ(ierr); 487064c8009SBarry Smith ierr = VecPlaceArray(b,xint_tmp+i*Nint);CHKERRQ(ierr); 4888e722e36SBarry Smith ierr = KSPSolve(exotic->ksp,b,x);CHKERRQ(ierr); 489c0decd05SBarry Smith ierr = KSPCheckSolve(exotic->ksp,pc,x);CHKERRQ(ierr); 490064c8009SBarry Smith ierr = VecResetArray(x);CHKERRQ(ierr); 491064c8009SBarry Smith ierr = VecResetArray(b);CHKERRQ(ierr); 492064c8009SBarry Smith } 4938c778c55SBarry Smith ierr = MatDenseRestoreArray(Xint,&xint);CHKERRQ(ierr); 4948c778c55SBarry Smith ierr = MatDenseRestoreArray(Xint_tmp,&xint_tmp);CHKERRQ(ierr); 4956bf464f9SBarry Smith ierr = VecDestroy(&x);CHKERRQ(ierr); 4966bf464f9SBarry Smith ierr = VecDestroy(&b);CHKERRQ(ierr); 497064c8009SBarry Smith } 4986bf464f9SBarry Smith ierr = MatDestroy(&Xint_tmp);CHKERRQ(ierr); 499064c8009SBarry Smith 500064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 5011683a169SBarry Smith ierr = MatDenseGetArrayRead(Xint,&rxint);CHKERRQ(ierr); 502064c8009SBarry Smith for (i=0; i<Nint; i++) { 503064c8009SBarry Smith tmp = 0.0; 5041683a169SBarry Smith for (j=0; j<6; j++) tmp += rxint[i+j*Nint]; 5052fa5cd67SKarl Rupp 506*98921bdaSJacob Faibussowitsch if (PetscAbsScalar(tmp-1.0) > 1.e-10) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong Xint interpolation at i %D value %g",i,(double)PetscAbsScalar(tmp)); 507064c8009SBarry Smith } 5081683a169SBarry Smith ierr = MatDenseRestoreArrayRead(Xint,&rxint);CHKERRQ(ierr); 509064c8009SBarry Smith /* ierr = MatView(Xint,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 510064c8009SBarry Smith #endif 511064c8009SBarry Smith 512064c8009SBarry Smith /* total faces */ 513064c8009SBarry Smith Ntotal = mp*np*(pp+1) + mp*pp*(np+1) + np*pp*(mp+1); 514064c8009SBarry Smith 515064c8009SBarry Smith /* 516064c8009SBarry Smith For each vertex, edge, face on process (in the same orderings as used above) determine its local number including ghost points 517064c8009SBarry Smith */ 518064c8009SBarry Smith cnt = 0; 519064c8009SBarry Smith { gl[cnt++] = mwidth+1;} 520064c8009SBarry Smith { 521064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth+1;} 522064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth + mwidth; /* these are the interior nodes */ gl[cnt++] = mwidth*nwidth + mwidth+m-istart-1;} 523064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth+mwidth*(n-jstart-1)+1;} 524064c8009SBarry Smith } 525064c8009SBarry Smith { gl[cnt++] = mwidth*nwidth*(p-kstart-1) + mwidth+1;} 526064c8009SBarry Smith 527064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 528064c8009SBarry Smith /* convert that to global numbering and get them on all processes */ 529064c8009SBarry Smith ierr = ISLocalToGlobalMappingApply(ltg,6,gl,gl);CHKERRQ(ierr); 530064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 531785e854fSJed Brown ierr = PetscMalloc1(6*mp*np*pp,&globals);CHKERRQ(ierr); 532ffc4695bSBarry Smith ierr = MPI_Allgather(gl,6,MPIU_INT,globals,6,MPIU_INT,PetscObjectComm((PetscObject)da));CHKERRMPI(ierr); 533064c8009SBarry Smith 534064c8009SBarry Smith /* Number the coarse grid points from 0 to Ntotal */ 5350298fd71SBarry Smith ierr = MatGetSize(Aglobal,&Nt,NULL);CHKERRQ(ierr); 53628d20b34SBarry Smith ierr = PetscTableCreate(Ntotal/3,Nt+1,&ht);CHKERRQ(ierr); 537064c8009SBarry Smith for (i=0; i<6*mp*np*pp; i++) { 538064c8009SBarry Smith ierr = PetscTableAddCount(ht,globals[i]+1);CHKERRQ(ierr); 539064c8009SBarry Smith } 540064c8009SBarry Smith ierr = PetscTableGetCount(ht,&cnt);CHKERRQ(ierr); 541*98921bdaSJacob Faibussowitsch if (cnt != Ntotal) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Hash table size %D not equal to total number coarse grid points %D",cnt,Ntotal); 542064c8009SBarry Smith ierr = PetscFree(globals);CHKERRQ(ierr); 543064c8009SBarry Smith for (i=0; i<6; i++) { 544064c8009SBarry Smith ierr = PetscTableFind(ht,gl[i]+1,&gl[i]);CHKERRQ(ierr); 545064c8009SBarry Smith gl[i]--; 546064c8009SBarry Smith } 5476bc0bbbfSBarry Smith ierr = PetscTableDestroy(&ht);CHKERRQ(ierr); 548064c8009SBarry Smith /* PetscIntView(6,gl,PETSC_VIEWER_STDOUT_WORLD); */ 549064c8009SBarry Smith 550064c8009SBarry Smith /* construct global interpolation matrix */ 5510298fd71SBarry Smith ierr = MatGetLocalSize(Aglobal,&Ng,NULL);CHKERRQ(ierr); 552064c8009SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 553ce94432eSBarry Smith ierr = MatCreateAIJ(PetscObjectComm((PetscObject)da),Ng,PETSC_DECIDE,PETSC_DECIDE,Ntotal,Nint+Nsurf,NULL,Nint,NULL,P);CHKERRQ(ierr); 554064c8009SBarry Smith } else { 555064c8009SBarry Smith ierr = MatZeroEntries(*P);CHKERRQ(ierr); 556064c8009SBarry Smith } 557064c8009SBarry Smith ierr = MatSetOption(*P,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 5581683a169SBarry Smith ierr = MatDenseGetArrayRead(Xint,&rxint);CHKERRQ(ierr); 5591683a169SBarry Smith ierr = MatSetValues(*P,Nint,IIint,6,gl,rxint,INSERT_VALUES);CHKERRQ(ierr); 5601683a169SBarry Smith ierr = MatDenseRestoreArrayRead(Xint,&rxint);CHKERRQ(ierr); 5611683a169SBarry Smith ierr = MatDenseGetArrayRead(Xsurf,&rxint);CHKERRQ(ierr); 5621683a169SBarry Smith ierr = MatSetValues(*P,Nsurf,IIsurf,6,gl,rxint,INSERT_VALUES);CHKERRQ(ierr); 5631683a169SBarry Smith ierr = MatDenseRestoreArrayRead(Xsurf,&rxint);CHKERRQ(ierr); 564064c8009SBarry Smith ierr = MatAssemblyBegin(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 565064c8009SBarry Smith ierr = MatAssemblyEnd(*P,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 566064c8009SBarry Smith ierr = PetscFree2(IIint,IIsurf);CHKERRQ(ierr); 567064c8009SBarry Smith 568064c8009SBarry Smith #if defined(PETSC_USE_DEBUG_foo) 569064c8009SBarry Smith { 570064c8009SBarry Smith Vec x,y; 571064c8009SBarry Smith PetscScalar *yy; 572ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),Ng,PETSC_DETERMINE,&y);CHKERRQ(ierr); 573ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)da),PETSC_DETERMINE,Ntotal,&x);CHKERRQ(ierr); 574064c8009SBarry Smith ierr = VecSet(x,1.0);CHKERRQ(ierr); 575064c8009SBarry Smith ierr = MatMult(*P,x,y);CHKERRQ(ierr); 576064c8009SBarry Smith ierr = VecGetArray(y,&yy);CHKERRQ(ierr); 577064c8009SBarry Smith for (i=0; i<Ng; i++) { 578*98921bdaSJacob Faibussowitsch if (PetscAbsScalar(yy[i]-1.0) > 1.e-10) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Wrong p interpolation at i %D value %g",i,(double)PetscAbsScalar(yy[i])); 579064c8009SBarry Smith } 580064c8009SBarry Smith ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr); 581064c8009SBarry Smith ierr = VecDestroy(x);CHKERRQ(ierr); 582064c8009SBarry Smith ierr = VecDestroy(y);CHKERRQ(ierr); 583064c8009SBarry Smith } 584064c8009SBarry Smith #endif 585064c8009SBarry Smith 586fcfd50ebSBarry Smith ierr = MatDestroy(&Aii);CHKERRQ(ierr); 587fcfd50ebSBarry Smith ierr = MatDestroy(&Ais);CHKERRQ(ierr); 588fcfd50ebSBarry Smith ierr = MatDestroy(&Asi);CHKERRQ(ierr); 589fcfd50ebSBarry Smith ierr = MatDestroy(&A);CHKERRQ(ierr); 590fcfd50ebSBarry Smith ierr = ISDestroy(&is);CHKERRQ(ierr); 591fcfd50ebSBarry Smith ierr = ISDestroy(&isint);CHKERRQ(ierr); 592fcfd50ebSBarry Smith ierr = ISDestroy(&issurf);CHKERRQ(ierr); 593fcfd50ebSBarry Smith ierr = MatDestroy(&Xint);CHKERRQ(ierr); 594fcfd50ebSBarry Smith ierr = MatDestroy(&Xsurf);CHKERRQ(ierr); 595064c8009SBarry Smith PetscFunctionReturn(0); 596064c8009SBarry Smith } 597174b6946SBarry Smith 5987233f9f0SBarry Smith /*@ 5997233f9f0SBarry Smith PCExoticSetType - Sets the type of coarse grid interpolation to use 6007233f9f0SBarry Smith 6013f9fe445SBarry Smith Logically Collective on PC 6027233f9f0SBarry Smith 6037233f9f0SBarry Smith Input Parameters: 6047233f9f0SBarry Smith + pc - the preconditioner context 6057233f9f0SBarry Smith - type - either PC_EXOTIC_FACE or PC_EXOTIC_WIREBASKET (defaults to face) 6067233f9f0SBarry Smith 60795452b02SPatrick Sanan Notes: 60895452b02SPatrick Sanan The face based interpolation has 1 degree of freedom per face and ignores the 609563e08c6SBarry Smith edge and vertex values completely in the coarse problem. For any seven point 610563e08c6SBarry Smith stencil the interpolation of a constant on all faces into the interior is that constant. 611563e08c6SBarry Smith 612563e08c6SBarry Smith The wirebasket interpolation has 1 degree of freedom per vertex, per edge and 613563e08c6SBarry Smith per face. A constant on the subdomain boundary is interpolated as that constant 614563e08c6SBarry Smith in the interior of the domain. 615563e08c6SBarry Smith 616563e08c6SBarry Smith The coarse grid matrix is obtained via the Galerkin computation A_c = R A R^T, hence 617563e08c6SBarry Smith if A is nonsingular A_c is also nonsingular. 618563e08c6SBarry Smith 619563e08c6SBarry Smith Both interpolations are suitable for only scalar problems. 620563e08c6SBarry Smith 6217233f9f0SBarry Smith Level: intermediate 6227233f9f0SBarry Smith 6237233f9f0SBarry Smith .seealso: PCEXOTIC, PCExoticType() 6247233f9f0SBarry Smith @*/ 6257087cfbeSBarry Smith PetscErrorCode PCExoticSetType(PC pc,PCExoticType type) 6267233f9f0SBarry Smith { 6274ac538c5SBarry Smith PetscErrorCode ierr; 6287233f9f0SBarry Smith 6297233f9f0SBarry Smith PetscFunctionBegin; 6300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 631c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 6324ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCExoticSetType_C",(PC,PCExoticType),(pc,type));CHKERRQ(ierr); 6337233f9f0SBarry Smith PetscFunctionReturn(0); 6347233f9f0SBarry Smith } 6357233f9f0SBarry Smith 636f7a08781SBarry Smith static PetscErrorCode PCExoticSetType_Exotic(PC pc,PCExoticType type) 6377233f9f0SBarry Smith { 638f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 63931567311SBarry Smith PC_Exotic *ctx = (PC_Exotic*) mg->innerctx; 6407233f9f0SBarry Smith 6417233f9f0SBarry Smith PetscFunctionBegin; 6427233f9f0SBarry Smith ctx->type = type; 6437233f9f0SBarry Smith PetscFunctionReturn(0); 6447233f9f0SBarry Smith } 6457233f9f0SBarry Smith 6467233f9f0SBarry Smith PetscErrorCode PCSetUp_Exotic(PC pc) 647174b6946SBarry Smith { 648174b6946SBarry Smith PetscErrorCode ierr; 64996bdf778SBarry Smith Mat A; 650f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 65131567311SBarry Smith PC_Exotic *ex = (PC_Exotic*) mg->innerctx; 65296bdf778SBarry Smith MatReuse reuse = (ex->P) ? MAT_REUSE_MATRIX : MAT_INITIAL_MATRIX; 653174b6946SBarry Smith 654174b6946SBarry Smith PetscFunctionBegin; 655ce94432eSBarry Smith if (!pc->dm) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Need to call PCSetDM() before using this PC"); 65623ee1639SBarry Smith ierr = PCGetOperators(pc,NULL,&A);CHKERRQ(ierr); 6577233f9f0SBarry Smith if (ex->type == PC_EXOTIC_FACE) { 658c0decd05SBarry Smith ierr = DMDAGetFaceInterpolation(pc,pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr); 6597233f9f0SBarry Smith } else if (ex->type == PC_EXOTIC_WIREBASKET) { 660c0decd05SBarry Smith ierr = DMDAGetWireBasketInterpolation(pc,pc->dm,ex,A,reuse,&ex->P);CHKERRQ(ierr); 661*98921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_PLIB,"Unknown exotic coarse space %d",ex->type); 66296bdf778SBarry Smith ierr = PCMGSetInterpolation(pc,1,ex->P);CHKERRQ(ierr); 663d0660788SBarry Smith /* if PC has attached DM we must remove it or the PCMG will use it to compute incorrect sized vectors and interpolations */ 6640298fd71SBarry Smith ierr = PCSetDM(pc,NULL);CHKERRQ(ierr); 6657233f9f0SBarry Smith ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 666174b6946SBarry Smith PetscFunctionReturn(0); 667174b6946SBarry Smith } 668174b6946SBarry Smith 6697233f9f0SBarry Smith PetscErrorCode PCDestroy_Exotic(PC pc) 670174b6946SBarry Smith { 671174b6946SBarry Smith PetscErrorCode ierr; 672f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 67331567311SBarry Smith PC_Exotic *ctx = (PC_Exotic*) mg->innerctx; 674174b6946SBarry Smith 675174b6946SBarry Smith PetscFunctionBegin; 676fcfd50ebSBarry Smith ierr = MatDestroy(&ctx->P);CHKERRQ(ierr); 677fcfd50ebSBarry Smith ierr = KSPDestroy(&ctx->ksp);CHKERRQ(ierr); 6787233f9f0SBarry Smith ierr = PetscFree(ctx);CHKERRQ(ierr); 6797233f9f0SBarry Smith ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 680174b6946SBarry Smith PetscFunctionReturn(0); 681174b6946SBarry Smith } 682174b6946SBarry Smith 6837233f9f0SBarry Smith PetscErrorCode PCView_Exotic(PC pc,PetscViewer viewer) 6847233f9f0SBarry Smith { 685f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 6867233f9f0SBarry Smith PetscErrorCode ierr; 687ace3abfcSBarry Smith PetscBool iascii; 68831567311SBarry Smith PC_Exotic *ctx = (PC_Exotic*) mg->innerctx; 6897233f9f0SBarry Smith 6907233f9f0SBarry Smith PetscFunctionBegin; 691251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 6927233f9f0SBarry Smith if (iascii) { 6937233f9f0SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Exotic type = %s\n",PCExoticTypes[ctx->type]);CHKERRQ(ierr); 6948e722e36SBarry Smith if (ctx->directSolve) { 6958e722e36SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using direct solver to construct interpolation\n");CHKERRQ(ierr); 6968e722e36SBarry Smith } else { 6978e722e36SBarry Smith PetscViewer sviewer; 6988e722e36SBarry Smith PetscMPIInt rank; 6998e722e36SBarry Smith 7008e722e36SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using iterative solver to construct interpolation\n");CHKERRQ(ierr); 7018e722e36SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 7028e722e36SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); /* should not need to push this twice? */ 7033f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 704ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRMPI(ierr); 705dd400576SPatrick Sanan if (rank == 0) { 7068e722e36SBarry Smith ierr = KSPView(ctx->ksp,sviewer);CHKERRQ(ierr); 7078e722e36SBarry Smith } 7083f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 7098e722e36SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7108e722e36SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 7118e722e36SBarry Smith } 7127233f9f0SBarry Smith } 7137233f9f0SBarry Smith ierr = PCView_MG(pc,viewer);CHKERRQ(ierr); 7147233f9f0SBarry Smith PetscFunctionReturn(0); 7157233f9f0SBarry Smith } 7167233f9f0SBarry Smith 7174416b707SBarry Smith PetscErrorCode PCSetFromOptions_Exotic(PetscOptionItems *PetscOptionsObject,PC pc) 7187233f9f0SBarry Smith { 7197233f9f0SBarry Smith PetscErrorCode ierr; 720ace3abfcSBarry Smith PetscBool flg; 721f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 7227233f9f0SBarry Smith PCExoticType mgctype; 72331567311SBarry Smith PC_Exotic *ctx = (PC_Exotic*) mg->innerctx; 7247233f9f0SBarry Smith 7257233f9f0SBarry Smith PetscFunctionBegin; 726e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Exotic coarse space options");CHKERRQ(ierr); 7277233f9f0SBarry Smith ierr = PetscOptionsEnum("-pc_exotic_type","face or wirebasket","PCExoticSetType",PCExoticTypes,(PetscEnum)ctx->type,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 7287233f9f0SBarry Smith if (flg) { 7297233f9f0SBarry Smith ierr = PCExoticSetType(pc,mgctype);CHKERRQ(ierr); 7307233f9f0SBarry Smith } 7310298fd71SBarry Smith ierr = PetscOptionsBool("-pc_exotic_direct_solver","use direct solver to construct interpolation","None",ctx->directSolve,&ctx->directSolve,NULL);CHKERRQ(ierr); 7328e722e36SBarry Smith if (!ctx->directSolve) { 7338e722e36SBarry Smith if (!ctx->ksp) { 7348e722e36SBarry Smith const char *prefix; 7358e722e36SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&ctx->ksp);CHKERRQ(ierr); 736422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(ctx->ksp,pc->erroriffailure);CHKERRQ(ierr); 7378e722e36SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 7383bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ctx->ksp);CHKERRQ(ierr); 7398e722e36SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 7408e722e36SBarry Smith ierr = KSPSetOptionsPrefix(ctx->ksp,prefix);CHKERRQ(ierr); 7418e722e36SBarry Smith ierr = KSPAppendOptionsPrefix(ctx->ksp,"exotic_");CHKERRQ(ierr); 7428e722e36SBarry Smith } 7438e722e36SBarry Smith ierr = KSPSetFromOptions(ctx->ksp);CHKERRQ(ierr); 7448e722e36SBarry Smith } 7457233f9f0SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7467233f9f0SBarry Smith PetscFunctionReturn(0); 7477233f9f0SBarry Smith } 7487233f9f0SBarry Smith 749174b6946SBarry Smith /*MC 7507233f9f0SBarry Smith PCEXOTIC - Two level overlapping Schwarz preconditioner with exotic (non-standard) coarse grid spaces 751174b6946SBarry Smith 7527233f9f0SBarry Smith This uses the PCMG infrastructure restricted to two levels and the face and wirebasket based coarse 75324c3aa18SBarry Smith grid spaces. 75424c3aa18SBarry Smith 75595452b02SPatrick Sanan Notes: 75695452b02SPatrick Sanan 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 75724c3aa18SBarry Smith 75896a0c994SBarry Smith References: 75996a0c994SBarry Smith + 1. - These coarse grid spaces originate in the work of Bramble, Pasciak and Schatz, "The Construction 76096a0c994SBarry Smith of Preconditioners for Elliptic Problems by Substructing IV", Mathematics of Computation, volume 53, 1989. 76196a0c994SBarry Smith . 2. - They were generalized slightly in "Domain Decomposition Method for Linear Elasticity", Ph. D. thesis, Barry Smith, 7624f02bc6aSBarry Smith New York University, 1990. 76396a0c994SBarry Smith . 3. - They were then explored in great detail in Dryja, Smith, Widlund, "Schwarz Analysis 7643b65e785SBarry Smith of Iterative Substructuring Methods for Elliptic Problems in Three Dimensions, SIAM Journal on Numerical 76596a0c994SBarry Smith Analysis, volume 31. 1994. These were developed in the context of iterative substructuring preconditioners. 76696a0c994SBarry Smith . 4. - They were then ingeniously applied as coarse grid spaces for overlapping Schwarz methods by Dohrmann and Widlund. 7673b65e785SBarry Smith They refer to them as GDSW (generalized Dryja, Smith, Widlund preconditioners). See, for example, 7683b65e785SBarry Smith Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Extending theory for domain decomposition algorithms to irregular subdomains. In Ulrich Langer, Marco 7693b65e785SBarry Smith Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings 7703b65e785SBarry Smith of the 17th International Conference on Domain Decomposition Methods in 77196a0c994SBarry Smith Science and Engineering, held in Strobl, Austria, 2006, number 60 in 77296a0c994SBarry Smith Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007. 77396a0c994SBarry Smith . 5. - Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. A family of energy minimizing coarse spaces for overlapping Schwarz preconditioners. In Ulrich Langer, 7743b65e785SBarry Smith Marco Discacciati, David Keyes, Olof Widlund, and Walter Zulehner, editors, Proceedings 7753b65e785SBarry Smith of the 17th International Conference on Domain Decomposition Methods 77696a0c994SBarry Smith in Science and Engineering, held in Strobl, Austria, 2006, number 60 in 77796a0c994SBarry Smith Springer Verlag, Lecture Notes in Computational Science and Engineering, 2007 77896a0c994SBarry Smith . 6. - Clark R. Dohrmann, Axel Klawonn, and Olof B. Widlund. Domain decomposition 7793b65e785SBarry Smith for less regular subdomains: Overlapping Schwarz in two dimensions. SIAM J. 78096a0c994SBarry Smith Numer. Anal., 46(4), 2008. 78196a0c994SBarry Smith - 7. - Clark R. Dohrmann and Olof B. Widlund. An overlapping Schwarz 7823b65e785SBarry Smith algorithm for almost incompressible elasticity. Technical Report 78396a0c994SBarry Smith TR2008 912, Department of Computer Science, Courant Institute 7843b65e785SBarry Smith of Mathematical Sciences, New York University, May 2008. URL: 7857233f9f0SBarry Smith 7867233f9f0SBarry Smith Options Database: The usual PCMG options are supported, such as -mg_levels_pc_type <type> -mg_coarse_pc_type <type> 7877233f9f0SBarry Smith -pc_mg_type <type> 7887233f9f0SBarry Smith 78925a35f6fSSatish Balay Level: advanced 790174b6946SBarry Smith 7916c699258SBarry Smith .seealso: PCMG, PCSetDM(), PCExoticType, PCExoticSetType() 792174b6946SBarry Smith M*/ 793174b6946SBarry Smith 7948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Exotic(PC pc) 795174b6946SBarry Smith { 796174b6946SBarry Smith PetscErrorCode ierr; 7977233f9f0SBarry Smith PC_Exotic *ex; 798f3fbd535SBarry Smith PC_MG *mg; 799174b6946SBarry Smith 800174b6946SBarry Smith PetscFunctionBegin; 801f91d8e95SBarry Smith /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 8022fa5cd67SKarl Rupp if (pc->ops->destroy) { 8032fa5cd67SKarl Rupp ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 8040a545947SLisandro Dalcin pc->data = NULL; 8052fa5cd67SKarl Rupp } 806503cfb0cSBarry Smith ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr); 8070a545947SLisandro Dalcin ((PetscObject)pc)->type_name = NULL; 808f91d8e95SBarry Smith 809174b6946SBarry Smith ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); 8100298fd71SBarry Smith ierr = PCMGSetLevels(pc,2,NULL);CHKERRQ(ierr); 8112134b1e4SBarry Smith ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);CHKERRQ(ierr); 812b00a9115SJed Brown ierr = PetscNew(&ex);CHKERRQ(ierr); \ 8137233f9f0SBarry Smith ex->type = PC_EXOTIC_FACE; 814f3fbd535SBarry Smith mg = (PC_MG*) pc->data; 81531567311SBarry Smith mg->innerctx = ex; 8167233f9f0SBarry Smith 8177233f9f0SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Exotic; 8187233f9f0SBarry Smith pc->ops->view = PCView_Exotic; 8197233f9f0SBarry Smith pc->ops->destroy = PCDestroy_Exotic; 8206c699258SBarry Smith pc->ops->setup = PCSetUp_Exotic; 8212fa5cd67SKarl Rupp 822bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCExoticSetType_C",PCExoticSetType_Exotic);CHKERRQ(ierr); 823174b6946SBarry Smith PetscFunctionReturn(0); 824174b6946SBarry Smith } 825