Actual source code: ex39.c

petsc-dev 2014-02-02
Report Typos and Errors
  2: static char help[] = "Lattice Gauge 2D model.\n"
  3:                      "Parameters:\n"
  4:                      "-size n          to use a grid size of n, i.e n space and n time steps\n"
  5:                      "-beta b          controls the randomness of the gauge field\n"
  6:                      "-rho r           the quark mass (?)";

  8: #include <petscksp.h>
  9: #include <petscpcasa.h>
 10: #include <petscdmda.h>

 12: PetscErrorCode computeMaxEigVal(Mat A, PetscInt its, PetscScalar *eig);

 16: int main(int Argc,char **Args)
 17: {
 18:   PetscBool      flg;
 19:   PetscInt       n   = -6;
 20:   PetscScalar    rho = 1.0;
 21:   PetscReal      h;
 22:   PetscReal      beta = 1.0;
 23:   DM             da;
 24:   PetscRandom    rctx;
 25:   PetscMPIInt    comm_size;
 26:   Mat            H,HtH;
 27:   PetscInt       x, y, xs, ys, xm, ym;
 28:   PetscReal      r1, r2;
 29:   PetscScalar    uxy1, uxy2;
 30:   MatStencil     sxy, sxy_m;
 31:   PetscScalar    val, valconj;
 32:   Vec            b, Htb,xvec;
 33:   KSP            kspmg;
 34:   PC             pcmg;
 36:   PetscInt       ix[1]   = {0};
 37:   PetscScalar    vals[1] = {1.0};

 39:   PetscInitialize(&Argc,&Args,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,"-size",&n,&flg);
 41:   PetscOptionsGetReal(NULL,"-beta",&beta,&flg);
 42:   PetscOptionsGetScalar(NULL,"-rho",&rho,&flg);

 44:   /* Set the fudge parameters, we scale the whole thing by 1/(2*h) later */
 45:   h    = 1.;
 46:   rho *= 1./(2.*h);

 48:   /* Geometry info */
 49:   DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_PERIODIC,DMDA_BOUNDARY_PERIODIC, DMDA_STENCIL_STAR, n, n,
 50:                       PETSC_DECIDE, PETSC_DECIDE, 2 /* this is the # of dof's */,
 51:                       1, NULL, NULL, &da);

 53:   /* Random numbers */
 54:   PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 55:   PetscRandomSetFromOptions(rctx);

 57:   /* Single or multi processor ? */
 58:   MPI_Comm_size(PETSC_COMM_WORLD,&comm_size);

 60:   /* construct matrix */
 61:   DMSetMatType(da,MATAIJ);
 62:   DMCreateMatrix(da, &H);

 64:   /* get local corners for this processor */
 65:   DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0);

 67:   /* Assemble the matrix */
 68:   for (x=xs; x<xs+xm; x++) {
 69:     for (y=ys; y<ys+ym; y++) {
 70:       /* each lattice point sets only the *forward* pointing parameters (right, down),
 71:          i.e. Nabla_1^+ and Nabla_2^+.
 72:          In this way we can use only local random number creation. That means
 73:          we also have to set the corresponding backward pointing entries. */
 74:       /* Compute some normally distributed random numbers via Box-Muller */
 75:       PetscRandomGetValueReal(rctx, &r1);
 76:       r1   = 1.-r1; /* to change from [0,1) to (0,1], which we need for the log */
 77:       PetscRandomGetValueReal(rctx, &r2);
 78:       PetscReal R = PetscSqrtReal(-2.*PetscLogReal(r1));
 79:       PetscReal c = PetscCosReal(2.*PETSC_PI*r2);
 80:       PetscReal s = PetscSinReal(2.*PETSC_PI*r2);

 82:       /* use those to set the field */
 83:       uxy1 = PetscExpScalar(((PetscScalar) (R*c/beta))*PETSC_i);
 84:       uxy2 = PetscExpScalar(((PetscScalar) (R*s/beta))*PETSC_i);

 86:       sxy.i = x; sxy.j = y; /* the point where we are */

 88:       /* center action */
 89:       sxy.c = 0; /* spin 0, 0 */
 90:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy, &rho, ADD_VALUES);
 91:       sxy.c = 1; /* spin 1, 1 */
 92:       val   = -rho;
 93:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy, &val, ADD_VALUES);

 95:       sxy_m.i = x+1; sxy_m.j = y; /* right action */
 96:       sxy.c   = 0; sxy_m.c = 0; /* spin 0, 0 */
 97:       val     = -uxy1; valconj = PetscConj(val);
 98:       MatSetValuesStencil(H, 1, &sxy_m, 1, &sxy, &val, ADD_VALUES);
 99:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy_m, &valconj, ADD_VALUES);
100:       sxy.c   = 0; sxy_m.c = 1; /* spin 0, 1 */
101:       val     = -uxy1; valconj = PetscConj(val);
102:       MatSetValuesStencil(H, 1, &sxy_m, 1, &sxy, &val, ADD_VALUES);
103:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy_m, &valconj, ADD_VALUES);
104:       sxy.c   = 1; sxy_m.c = 0; /* spin 1, 0 */
105:       val     = uxy1; valconj = PetscConj(val);
106:       MatSetValuesStencil(H, 1, &sxy_m, 1, &sxy, &val, ADD_VALUES);
107:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy_m, &valconj, ADD_VALUES);
108:       sxy.c   = 1; sxy_m.c = 1; /* spin 1, 1 */
109:       val     = uxy1; valconj = PetscConj(val);
110:       MatSetValuesStencil(H, 1, &sxy_m, 1, &sxy, &val, ADD_VALUES);
111:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy_m, &valconj, ADD_VALUES);

113:       sxy_m.i = x; sxy_m.j = y+1; /* down action */
114:       sxy.c   = 0; sxy_m.c = 0; /* spin 0, 0 */
115:       val     = -uxy2; valconj = PetscConj(val);
116:       MatSetValuesStencil(H, 1, &sxy_m, 1, &sxy, &val, ADD_VALUES);
117:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy_m, &valconj, ADD_VALUES);
118:       sxy.c   = 0; sxy_m.c = 1; /* spin 0, 1 */
119:       val     = -PETSC_i*uxy2; valconj = PetscConj(val);
120:       MatSetValuesStencil(H, 1, &sxy_m, 1, &sxy, &val, ADD_VALUES);
121:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy_m, &valconj, ADD_VALUES);
122:       sxy.c   = 1; sxy_m.c = 0; /* spin 1, 0 */
123:       val     = -PETSC_i*uxy2; valconj = PetscConj(val);
124:       MatSetValuesStencil(H, 1, &sxy_m, 1, &sxy, &val, ADD_VALUES);
125:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy_m, &valconj, ADD_VALUES);
126:       sxy.c   = 1; sxy_m.c = 1; /* spin 1, 1 */
127:       val     = PetscConj(uxy2); valconj = PetscConj(val);
128:       MatSetValuesStencil(H, 1, &sxy_m, 1, &sxy, &val, ADD_VALUES);
129:       MatSetValuesStencil(H, 1, &sxy, 1, &sxy_m, &valconj, ADD_VALUES);
130:     }
131:   }

133:   MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY);
134:   MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY);

136:   /* scale H */
137:   MatScale(H, 1./(2.*h));

139:   /* it looks like H is Hermetian */
140:   /* construct normal equations */
141:   MatMatMult(H, H, MAT_INITIAL_MATRIX, 1., &HtH);

143:   /* permutation matrix to check whether H and HtH are identical to the ones in the paper */
144: /*   Mat perm; */
145: /*   DMCreateMatrix(da, &perm); */
146: /*   PetscInt row, col; */
147: /*   PetscScalar one = 1.0; */
148: /*   for (PetscInt i=0; i<n; i++) { */
149: /*     for (PetscInt j=0; j<n; j++) { */
150: /*       row = (i*n+j)*2; col = i*n+j; */
151: /*       MatSetValues(perm, 1, &row, 1, &col, &one, INSERT_VALUES); */
152: /*       row = (i*n+j)*2+1; col = i*n+j + n*n; */
153: /*       MatSetValues(perm, 1, &row, 1, &col, &one, INSERT_VALUES); */
154: /*     } */
155: /*   } */
156: /*   MatAssemblyBegin(perm, MAT_FINAL_ASSEMBLY); */
157: /*   MatAssemblyEnd(perm, MAT_FINAL_ASSEMBLY); */

159: /*   Mat Hperm; */
160: /*   MatPtAP(H, perm, MAT_INITIAL_MATRIX, 1.0, &Hperm); */
161: /*   PetscPrintf(PETSC_COMM_WORLD, "Matrix H after construction\n"); */
162: /*   MatView(Hperm, PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD)); */

164: /*   Mat HtHperm; */
165: /*   MatPtAP(HtH, perm, MAT_INITIAL_MATRIX, 1.0, &HtHperm); */
166: /*   PetscPrintf(PETSC_COMM_WORLD, "Matrix HtH:\n"); */
167: /*   MatView(HtHperm, PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD)); */

169:   /* right hand side */
170:   DMCreateGlobalVector(da, &b);
171:   VecSet(b,0.0);
172:   VecSetValues(b, 1, ix, vals, INSERT_VALUES);
173:   VecAssemblyBegin(b);
174:   VecAssemblyEnd(b);
175: /*   VecSetRandom(b, rctx); */
176:   VecDuplicate(b, &Htb);
177:   MatMultTranspose(H, b, Htb);

179:   /* construct solver */
180:   KSPCreate(PETSC_COMM_WORLD,&kspmg);
181:   KSPSetType(kspmg, KSPCG);

183:   KSPGetPC(kspmg,&pcmg);
184:   PCSetType(pcmg,PCASA);

186:   /* maybe user wants to override some of the choices */
187:   KSPSetFromOptions(kspmg);

189:   KSPSetOperators(kspmg, HtH, HtH, DIFFERENT_NONZERO_PATTERN);

191:   DMDASetRefinementFactor(da, 3, 3, 3);
192:   PCSetDM(pcmg,da);

194:   PCASASetTolerances(pcmg, 1.e-6, 1.e-10,PETSC_DEFAULT,PETSC_DEFAULT);

196:   VecDuplicate(b, &xvec);
197:   VecSet(xvec, 0.0);

199:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200:                       Solve the linear system
201:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

203:   KSPSolve(kspmg, Htb, xvec);

205: /*   VecView(xvec, PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD)); */

207:   KSPDestroy(&kspmg);
208:   VecDestroy(&xvec);

210:   /*   seems to be destroyed by KSPDestroy */
211:   VecDestroy(&b);
212:   VecDestroy(&Htb);
213:   MatDestroy(&HtH);
214:   MatDestroy(&H);

216:   DMDestroy(&da);
217:   PetscRandomDestroy(&rctx);
218:   PetscFinalize();
219:   return 0;
220: }

222: /* --------------------------------------------------------------------- */
225: PetscErrorCode computeMaxEigVal(Mat A, PetscInt its, PetscScalar *eig)
226: {
228:   PetscRandom    rctx;      /* random number generator context */
229:   Vec            x0, x, x_1, tmp;
230:   PetscScalar    lambda_its, lambda_its_1;
231:   PetscInt       i;

234:   PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
235:   PetscRandomSetFromOptions(rctx);
236:   MatGetVecs(A, &x_1, &x);
237:   VecSetRandom(x, rctx);
238:   VecDuplicate(x, &x0);
239:   VecCopy(x, x0);

241:   MatMult(A, x, x_1);
242:   for (i=0; i<its; i++) {
243:     tmp  = x; x = x_1; x_1 = tmp;
244:     MatMult(A, x, x_1);
245:   }
246:   VecDot(x0, x, &lambda_its);
247:   VecDot(x0, x_1, &lambda_its_1);

249:   *eig = lambda_its_1/lambda_its;

251:   VecDestroy(&x0);
252:   VecDestroy(&x);
253:   VecDestroy(&x_1);
254:   return(0);
255: }