Actual source code: ex17.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: static char help[] = "Solves a linear system with KSP.  This problem is\n\
  3: intended to test the complex numbers version of various solvers.\n\n";

  5: #include <petscksp.h>

  7: typedef enum {TEST_1,TEST_2,TEST_3,HELMHOLTZ_1,HELMHOLTZ_2} TestType;
  8: extern PetscErrorCode FormTestMatrix(Mat,PetscInt,TestType);

 12: int main(int argc,char **args)
 13: {
 14:   Vec            x,b,u;      /* approx solution, RHS, exact solution */
 15:   Mat            A;            /* linear system matrix */
 16:   KSP            ksp;         /* KSP context */
 18:   PetscInt       n    = 10,its, dim,p = 1,use_random;
 19:   PetscScalar    none = -1.0,pfive = 0.5;
 20:   PetscReal      norm;
 21:   PetscRandom    rctx;
 22:   TestType       type;
 23:   PetscBool      flg;

 25:   PetscInitialize(&argc,&args,(char*)0,help);
 26:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 27:   PetscOptionsGetInt(NULL,"-p",&p,NULL);
 28:   switch (p) {
 29:   case 1:  type = TEST_1;      dim = n;   break;
 30:   case 2:  type = TEST_2;      dim = n;   break;
 31:   case 3:  type = TEST_3;      dim = n;   break;
 32:   case 4:  type = HELMHOLTZ_1; dim = n*n; break;
 33:   case 5:  type = HELMHOLTZ_2; dim = n*n; break;
 34:   default: type = TEST_1;      dim = n;
 35:   }

 37:   /* Create vectors */
 38:   VecCreate(PETSC_COMM_WORLD,&x);
 39:   VecSetSizes(x,PETSC_DECIDE,dim);
 40:   VecSetFromOptions(x);
 41:   VecDuplicate(x,&b);
 42:   VecDuplicate(x,&u);

 44:   use_random = 1;
 45:   flg        = PETSC_FALSE;
 46:   PetscOptionsGetBool(NULL,"-norandom",&flg,NULL);
 47:   if (flg) {
 48:     use_random = 0;
 49:     VecSet(u,pfive);
 50:   } else {
 51:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 52:     PetscRandomSetFromOptions(rctx);
 53:     VecSetRandom(u,rctx);
 54:   }

 56:   /* Create and assemble matrix */
 57:   MatCreate(PETSC_COMM_WORLD,&A);
 58:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim);
 59:   MatSetFromOptions(A);
 60:   MatSetUp(A);
 61:   FormTestMatrix(A,n,type);
 62:   MatMult(A,u,b);
 63:   flg  = PETSC_FALSE;
 64:   PetscOptionsGetBool(NULL,"-printout",&flg,NULL);
 65:   if (flg) {
 66:     MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 67:     VecView(u,PETSC_VIEWER_STDOUT_WORLD);
 68:     VecView(b,PETSC_VIEWER_STDOUT_WORLD);
 69:   }

 71:   /* Create KSP context; set operators and options; solve linear system */
 72:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 73:   KSPSetOperators(ksp,A,A);
 74:   KSPSetFromOptions(ksp);
 75:   KSPSolve(ksp,b,x);
 76:   KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);

 78:   /* Check error */
 79:   VecAXPY(x,none,u);
 80:   VecNorm(x,NORM_2,&norm);
 81:   KSPGetIterationNumber(ksp,&its);
 82:   PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g,Iterations %D\n",(double)norm,its);

 84:   /* Free work space */
 85:   VecDestroy(&x); VecDestroy(&u);
 86:   VecDestroy(&b); MatDestroy(&A);
 87:   if (use_random) {PetscRandomDestroy(&rctx);}
 88:   KSPDestroy(&ksp);
 89:   PetscFinalize();
 90:   return 0;
 91: }

 95: PetscErrorCode FormTestMatrix(Mat A,PetscInt n,TestType type)
 96: {
 97: #if !defined(PETSC_USE_COMPLEX)
 98:   SETERRQ(PetscObjectComm((PetscObject)A),1,"FormTestMatrix: These problems require complex numbers.");
 99: #else

101:   PetscScalar    val[5];
103:   PetscInt       i,j,Ii,J,col[5],Istart,Iend;

105:   MatGetOwnershipRange(A,&Istart,&Iend);
106:   if (type == TEST_1) {
107:     val[0] = 1.0; val[1] = 4.0; val[2] = -2.0;
108:     for (i=1; i<n-1; i++) {
109:       col[0] = i-1; col[1] = i; col[2] = i+1;
110:       MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
111:     }
112:     i    = n-1; col[0] = n-2; col[1] = n-1;
113:     MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
114:     i    = 0; col[0] = 0; col[1] = 1; val[0] = 4.0; val[1] = -2.0;
115:     MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
116:   } else if (type == TEST_2) {
117:     val[0] = 1.0; val[1] = 0.0; val[2] = 2.0; val[3] = 1.0;
118:     for (i=2; i<n-1; i++) {
119:       col[0] = i-2; col[1] = i-1; col[2] = i; col[3] = i+1;
120:       MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
121:     }
122:     i    = n-1; col[0] = n-3; col[1] = n-2; col[2] = n-1;
123:     MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
124:     i    = 1; col[0] = 0; col[1] = 1; col[2] = 2;
125:     MatSetValues(A,1,&i,3,col,&val[1],INSERT_VALUES);
126:     i    = 0;
127:     MatSetValues(A,1,&i,2,col,&val[2],INSERT_VALUES);
128:   } else if (type == TEST_3) {
129:     val[0] = PETSC_i * 2.0;
130:     val[1] = 4.0; val[2] = 0.0; val[3] = 1.0; val[4] = 0.7;
131:     for (i=1; i<n-3; i++) {
132:       col[0] = i-1; col[1] = i; col[2] = i+1; col[3] = i+2; col[4] = i+3;
133:       MatSetValues(A,1,&i,5,col,val,INSERT_VALUES);
134:     }
135:     i    = n-3; col[0] = n-4; col[1] = n-3; col[2] = n-2; col[3] = n-1;
136:     MatSetValues(A,1,&i,4,col,val,INSERT_VALUES);
137:     i    = n-2; col[0] = n-3; col[1] = n-2; col[2] = n-1;
138:     MatSetValues(A,1,&i,3,col,val,INSERT_VALUES);
139:     i    = n-1; col[0] = n-2; col[1] = n-1;
140:     MatSetValues(A,1,&i,2,col,val,INSERT_VALUES);
141:     i    = 0; col[0] = 0; col[1] = 1; col[2] = 2; col[3] = 3;
142:     MatSetValues(A,1,&i,4,col,&val[1],INSERT_VALUES);
143:   } else if (type == HELMHOLTZ_1) {
144:     /* Problem domain: unit square: (0,1) x (0,1)
145:        Solve Helmholtz equation:
146:           -delta u - sigma1*u + i*sigma2*u = f,
147:            where delta = Laplace operator
148:        Dirichlet b.c.'s on all sides
149:      */
150:     PetscRandom rctx;
151:     PetscReal   h2,sigma1 = 5.0;
152:     PetscScalar sigma2;
153:     PetscOptionsGetReal(NULL,"-sigma1",&sigma1,NULL);
154:     PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
155:     PetscRandomSetFromOptions(rctx);
156:     PetscRandomSetInterval(rctx,0.0,PETSC_i);
157:     h2   = 1.0/((n+1)*(n+1));
158:     for (Ii=Istart; Ii<Iend; Ii++) {
159:       *val = -1.0; i = Ii/n; j = Ii - i*n;
160:       if (i>0) {
161:         J = Ii-n; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
162:       }
163:       if (i<n-1) {
164:         J = Ii+n; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
165:       }
166:       if (j>0) {
167:         J = Ii-1; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
168:       }
169:       if (j<n-1) {
170:         J = Ii+1; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
171:       }
172:       PetscRandomGetValue(rctx,&sigma2);
173:       *val = 4.0 - sigma1*h2 + sigma2*h2;
174:       MatSetValues(A,1,&Ii,1,&Ii,val,ADD_VALUES);
175:     }
176:     PetscRandomDestroy(&rctx);
177:   } else if (type == HELMHOLTZ_2) {
178:     /* Problem domain: unit square: (0,1) x (0,1)
179:        Solve Helmholtz equation:
180:           -delta u - sigma1*u = f,
181:            where delta = Laplace operator
182:        Dirichlet b.c.'s on 3 sides
183:        du/dn = i*alpha*u on (1,y), 0<y<1
184:      */
185:     PetscReal   h2,sigma1 = 200.0;
186:     PetscScalar alpha_h;
187:     PetscOptionsGetReal(NULL,"-sigma1",&sigma1,NULL);
188:     h2      = 1.0/((n+1)*(n+1));
189:     alpha_h = (PETSC_i * 10.0) / (PetscReal)(n+1);  /* alpha_h = alpha * h */
190:     for (Ii=Istart; Ii<Iend; Ii++) {
191:       *val = -1.0; i = Ii/n; j = Ii - i*n;
192:       if (i>0) {
193:         J = Ii-n; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
194:       }
195:       if (i<n-1) {
196:         J = Ii+n; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
197:       }
198:       if (j>0) {
199:         J = Ii-1; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
200:       }
201:       if (j<n-1) {
202:         J = Ii+1; MatSetValues(A,1,&Ii,1,&J,val,ADD_VALUES);
203:       }
204:       *val = 4.0 - sigma1*h2;
205:       if (!((Ii+1)%n)) *val += alpha_h;
206:       MatSetValues(A,1,&Ii,1,&Ii,val,ADD_VALUES);
207:     }
208:   } else SETERRQ(PetscObjectComm((PetscObject)A),1,"FormTestMatrix: unknown test matrix type");

210:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
211:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
212: #endif

214:   return 0;
215: }