Actual source code: ex54.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] = "Creates a matrix from quadrilateral finite elements in 2D, Laplacian \n\
  3:   -ne <size>       : problem size in number of elements (eg, -ne 31 gives 32^2 grid)\n\
  4:   -alpha <v>      : scaling of material coeficient in embedded circle\n\n";

  6: #include <petscksp.h>

 10: int main(int argc,char **args)
 11: {
 12:   Mat            Amat,Pmat;
 14:   PetscInt       i,m,M,its,Istart,Iend,j,Ii,ix,ne=4;
 15:   PetscReal      x,y,h;
 16:   Vec            xx,bb;
 17:   KSP            ksp;
 18:   PetscReal      soft_alpha = 1.e-3;
 19:   MPI_Comm       comm;
 20:   PetscMPIInt    npe,mype;
 21:   PetscScalar    DD[4][4],DD2[4][4];
 22: #if defined(PETSC_USE_LOG)
 23:   PetscLogStage stage;
 24: #endif
 25: #define DIAG_S 0.0
 26:   PetscScalar DD1[4][4] = { {5.0+DIAG_S, -2.0, -1.0, -2.0},
 27:                             {-2.0, 5.0+DIAG_S, -2.0, -1.0},
 28:                             {-1.0, -2.0, 5.0+DIAG_S, -2.0},
 29:                             {-2.0, -1.0, -2.0, 5.0+DIAG_S} };

 31:   PetscInitialize(&argc,&args,(char*)0,help);
 32:   comm = PETSC_COMM_WORLD;
 33:   MPI_Comm_rank(comm, &mype);
 34:   MPI_Comm_size(comm, &npe);
 35:   PetscOptionsGetInt(NULL,NULL,"-ne",&ne,NULL);
 36:   h     = 1./ne;
 37:   /* ne*ne; number of global elements */
 38:   PetscOptionsGetReal(NULL,NULL,"-alpha",&soft_alpha,NULL);
 39:   M    = (ne+1)*(ne+1); /* global number of nodes */
 40:   /* create stiffness matrix */
 41:   MatCreateAIJ(comm,PETSC_DECIDE,PETSC_DECIDE,M,M,
 42:                       18,NULL,6,NULL,&Amat);
 43:   MatCreateAIJ(comm,PETSC_DECIDE,PETSC_DECIDE,M,M,
 44:                       18,NULL,6,NULL,&Pmat);
 45:   MatGetOwnershipRange(Amat,&Istart,&Iend);
 46:   m    = Iend-Istart;
 47:   /* Generate vectors */
 48:   VecCreate(comm,&xx);
 49:   VecSetSizes(xx,m,M);
 50:   VecSetFromOptions(xx);
 51:   VecDuplicate(xx,&bb);
 52:   VecSet(bb,.0);
 53:   /* generate element matrices -- see ex56.c on how to use different data set */
 54:   {
 55:     DD1[0][0] =  0.66666666666666663;
 56:     DD1[0][1] = -0.16666666666666669;
 57:     DD1[0][2] = -0.33333333333333343;
 58:     DD1[0][3] = -0.16666666666666666;
 59:     DD1[1][0] = -0.16666666666666669;
 60:     DD1[1][1] =  0.66666666666666663;
 61:     DD1[1][2] = -0.16666666666666666;
 62:     DD1[1][3] = -0.33333333333333343;
 63:     DD1[2][0] = -0.33333333333333343;
 64:     DD1[2][1] = -0.16666666666666666;
 65:     DD1[2][2] =  0.66666666666666663;
 66:     DD1[2][3] = -0.16666666666666663;
 67:     DD1[3][0] = -0.16666666666666666;
 68:     DD1[3][1] = -0.33333333333333343;
 69:     DD1[3][2] = -0.16666666666666663;
 70:     DD1[3][3] =  0.66666666666666663;

 72:     /* BC version of element */
 73:     for (i=0;i<4;i++) {
 74:       for (j=0;j<4;j++) {
 75:         if (i<2 || j < 2) {
 76:           if (i==j) DD2[i][j] = .1*DD1[i][j];
 77:           else DD2[i][j] = 0.0;
 78:         } else DD2[i][j] = DD1[i][j];
 79:       }
 80:     }
 81:   }
 82:   {
 83:     PetscReal *coords;
 84:     PC             pc;
 85:     PetscMalloc1(2*m,&coords);
 86:     /* forms the element stiffness for the Laplacian and coordinates */
 87:     for (Ii=Istart,ix=0; Ii<Iend; Ii++,ix++) {
 88:       j = Ii/(ne+1); i = Ii%(ne+1);
 89:       /* coords */
 90:       x            = h*(Ii % (ne+1)); y = h*(Ii/(ne+1));
 91:       coords[2*ix] = x; coords[2*ix+1] = y;
 92:       if (i<ne && j<ne) {
 93:         PetscInt jj,ii,idx[4];
 94:         /* radius */
 95:         PetscReal radius = PetscSqrtReal((x-.5+h/2)*(x-.5+h/2) + (y-.5+h/2)*(y-.5+h/2));
 96:         PetscReal alpha  = 1.0;

 98:         idx[0] = Ii; idx[1] = Ii+1; idx[2] = Ii + (ne+1) + 1; idx[3] =  Ii + (ne+1);
 99:         if (radius < 0.25) alpha = soft_alpha;


102:         for (ii=0; ii<4; ii++) {
103:           for (jj=0; jj<4; jj++) DD[ii][jj] = alpha*DD1[ii][jj];
104:         }
105:         MatSetValues(Pmat,4,idx,4,idx,(const PetscScalar*)DD,ADD_VALUES);
106:         if (j>0) {
107:           MatSetValues(Amat,4,idx,4,idx,(const PetscScalar*)DD,ADD_VALUES);
108:         } else {
109:           /* a BC */
110:           for (ii=0;ii<4;ii++) {
111:             for (jj=0;jj<4;jj++) DD[ii][jj] = alpha*DD2[ii][jj];
112:           }
113:           MatSetValues(Amat,4,idx,4,idx,(const PetscScalar*)DD,ADD_VALUES);
114:         }
115:       }
116:       if (j>0) {
117:         PetscScalar v  = h*h;
118:         PetscInt    jj = Ii;
119:         VecSetValues(bb,1,&jj,&v,INSERT_VALUES);
120:       }
121:     }
122:     MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);
123:     MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);
124:     MatAssemblyBegin(Pmat,MAT_FINAL_ASSEMBLY);
125:     MatAssemblyEnd(Pmat,MAT_FINAL_ASSEMBLY);
126:     VecAssemblyBegin(bb);
127:     VecAssemblyEnd(bb);

129:     /* Setup solver */
130:     KSPCreate(PETSC_COMM_WORLD,&ksp);
131:     KSPSetFromOptions(ksp);

133:     /* finish KSP/PC setup */
134:     KSPSetOperators(ksp, Amat, Amat);

136:     KSPGetPC(ksp,&pc);
137:     PCSetCoordinates(pc, 2, m, coords);
138:     PetscFree(coords);
139:   }

141:   if (!PETSC_TRUE) {
142:     PetscViewer viewer;
143:     PetscViewerASCIIOpen(comm, "Amat.m", &viewer);
144:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
145:     MatView(Amat,viewer);
146:     PetscViewerPopFormat(viewer);
147:     PetscViewerDestroy(&viewer);
148:   }

150:   /* solve */
151: #if defined(PETSC_USE_LOG)
152:   PetscLogStageRegister("Solve", &stage);
153:   PetscLogStagePush(stage);
154: #endif
155:   VecSet(xx,.0);

157:   KSPSetUp(ksp);

159:   KSPSolve(ksp,bb,xx);

161: #if defined(PETSC_USE_LOG)
162:   PetscLogStagePop();
163: #endif

165:   KSPGetIterationNumber(ksp,&its);

167:   if (!PETSC_TRUE) {
168:     PetscReal   norm,norm2;
169:     PetscViewer viewer;
170:     Vec         res;
171:     PetscViewerASCIIOpen(comm, "rhs.m", &viewer);
172:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
173:     VecView(bb,viewer);
174:     PetscViewerPopFormat(viewer);
175:     PetscViewerDestroy(&viewer);
176:     VecNorm(bb, NORM_2, &norm2);

178:     PetscViewerASCIIOpen(comm, "solution.m", &viewer);
179:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
180:     VecView(xx,viewer);
181:     PetscViewerPopFormat(viewer);
182:     PetscViewerDestroy(&viewer);

184:     VecDuplicate(xx, &res);
185:     MatMult(Amat, xx, res);
186:     VecAXPY(bb, -1.0, res);
187:     VecDestroy(&res);
188:     VecNorm(bb,NORM_2,&norm);
189:     PetscPrintf(PETSC_COMM_WORLD,"[%d]%s |b-Ax|/|b|=%e, |b|=%e\n",0,__FUNCT__,norm/norm2,norm2);

191:     PetscViewerASCIIOpen(comm, "residual.m", &viewer);
192:     PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
193:     VecView(bb,viewer);
194:     PetscViewerPopFormat(viewer);
195:     PetscViewerDestroy(&viewer);
196:   }

198:   /* Free work space */
199:   KSPDestroy(&ksp);
200:   VecDestroy(&xx);
201:   VecDestroy(&bb);
202:   MatDestroy(&Amat);
203:   MatDestroy(&Pmat);

205:   PetscFinalize();
206:   return 0;
207: }