Actual source code: ex54.c

petsc-3.8.4 2018-03-24
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>

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

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

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

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


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

127:     /* Setup solver */
128:     KSPCreate(PETSC_COMM_WORLD,&ksp);
129:     KSPSetFromOptions(ksp);

131:     /* finish KSP/PC setup */
132:     KSPSetOperators(ksp, Amat, Amat);

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

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

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

155:   KSPSetUp(ksp);

157:   KSPSolve(ksp,bb,xx);

159: #if defined(PETSC_USE_LOG)
160:   PetscLogStagePop();
161: #endif

163:   KSPGetIterationNumber(ksp,&its);

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

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

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

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

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

203:   PetscFinalize();
204:   return ierr;
205: }