Actual source code: ex40.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 <petscdm.h>
 11: #include <petscdmadda.h>

 13: PetscErrorCode computeMinEigVal(Mat A, PetscInt its, PetscScalar *eig);

 17: int main(int Argc,char **Args)
 18: {
 19:   PetscBool      flg;
 20:   PetscInt       n   = 6,i;
 21:   PetscScalar    rho = 1.0;
 22:   PetscReal      h;
 23:   PetscReal      beta = 1.0;
 24:   DM             adda;
 25:   PetscInt       nodes[2];
 26:   PetscBool      periodic[2];
 27:   PetscInt       refine[2];
 28:   PetscRandom    rctx;
 29:   PetscMPIInt    comm_size;
 30:   Mat            H;
 31:   PetscInt       *lcs, *lce;
 32:   PetscInt       x, y;
 33:   PetscReal      r1, r2;
 34:   PetscScalar    uxy1, uxy2;
 35:   ADDAIdx        sxy, sxy_m;
 36:   PetscScalar    val, valconj;
 37:   Mat            HtH;
 38:   Vec            b, Htb;
 39:   Vec            xvec;
 40:   KSP            kspmg;
 41:   PC             pcmg;

 44:   PetscInitialize(&Argc,&Args,(char*)0,help);
 45:   PetscOptionsGetInt(NULL,"-size",&n,&flg);
 46:   PetscOptionsGetReal(NULL,"-beta",&beta,&flg);
 47:   PetscOptionsGetScalar(NULL,"-rho",&rho,&flg);

 49:   /* Set the fudge parameters, we scale the whole thing by 1/(2*h) later */
 50:   h    = 1.;
 51:   rho *= 1./(2.*h);

 53:   /* Geometry info */
 54:   for (i=0; i<2; i++) {
 55:     nodes[i]    = n;
 56:     periodic[i] = PETSC_TRUE;
 57:     refine[i]   = 3;
 58:   }
 59:   DMADDACreate(PETSC_COMM_WORLD, 2, nodes, NULL, 2,periodic, &adda);
 60:   DMADDASetRefinement(adda, refine, 2);

 62:   /* Random numbers */
 63:   PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
 64:   PetscRandomSetFromOptions(rctx);

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

 69:   /* construct matrix */
 70:   DMSetMatType(adda,MATAIJ);
 71:   DMCreateMatrix(adda, &H);

 73:   /* get local corners for this processor, user is responsible for freeing lcs,lce */
 74:   DMADDAGetCorners(adda, &lcs, &lce);

 76:   /* Allocate space for the indices that we use to construct the matrix */
 77:   PetscMalloc1(2, &(sxy.x));
 78:   PetscMalloc1(2, &(sxy_m.x));

 80:   /* Assemble the matrix */
 81:   for (x=lcs[0]; x<lce[0]; x++) {
 82:     for (y=lcs[1]; y<lce[1]; y++) {
 83:       /* each lattice point sets only the *forward* pointing parameters (right, down),
 84:          i.e. Nabla_1^+ and Nabla_2^+.
 85:          In this way we can use only local random number creation. That means
 86:          we also have to set the corresponding backward pointing entries. */
 87:       /* Compute some normally distributed random numbers via Box-Muller */
 88:       PetscRandomGetValueReal(rctx, &r1);
 89:       r1   = 1.-r1; /* to change from [0,1) to (0,1], which we need for the log */
 90:       PetscRandomGetValueReal(rctx, &r2);
 91:       PetscReal R = PetscSqrtReal(-2.*log(r1));
 92:       PetscReal c = PetscCosReal(2.*PETSC_PI*r2);
 93:       PetscReal s = PetscSinReal(2.*PETSC_PI*r2);

 95:       /* use those to set the field */
 96:       uxy1 = PetscExpScalar(((PetscScalar) (R*c/beta))*PETSC_i);
 97:       uxy2 = PetscExpScalar(((PetscScalar) (R*s/beta))*PETSC_i);

 99:       sxy.x[0] = x; sxy.x[1] = y; /* the point where we are */

101:       /* center action */
102:       sxy.d = 0; /* spin 0, 0 */
103:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy, &rho, ADD_VALUES);
104:       sxy.d = 1; /* spin 1, 1 */
105:       val   = -rho;
106:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy, &val, ADD_VALUES);

108:       sxy_m.x[0] = x+1; sxy_m.x[1] = y; /* right action */
109:       sxy.d      = 0; sxy_m.d = 0; /* spin 0, 0 */
110:       val        = -uxy1; valconj = PetscConj(val);

112:       DMADDAMatSetValues(H, adda, 1, &sxy_m, adda, 1, &sxy, &val, ADD_VALUES);
113:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy_m, &valconj, ADD_VALUES);

115:       sxy.d = 0; sxy_m.d = 1; /* spin 0, 1 */
116:       val   = -uxy1; valconj = PetscConj(val);

118:       DMADDAMatSetValues(H, adda, 1, &sxy_m, adda, 1, &sxy, &val, ADD_VALUES);
119:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy_m, &valconj, ADD_VALUES);

121:       sxy.d = 1; sxy_m.d = 0; /* spin 1, 0 */
122:       val   = uxy1; valconj = PetscConj(val);

124:       DMADDAMatSetValues(H, adda, 1, &sxy_m, adda, 1, &sxy, &val, ADD_VALUES);
125:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy_m, &valconj, ADD_VALUES);

127:       sxy.d = 1; sxy_m.d = 1; /* spin 1, 1 */
128:       val   = uxy1; valconj = PetscConj(val);

130:       DMADDAMatSetValues(H, adda, 1, &sxy_m, adda, 1, &sxy, &val, ADD_VALUES);
131:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy_m, &valconj, ADD_VALUES);

133:       sxy_m.x[0] = x; sxy_m.x[1] = y+1; /* down action */
134:       sxy.d = 0; sxy_m.d = 0; /* spin 0, 0 */
135:       val   = -uxy2; valconj = PetscConj(val);

137:       DMADDAMatSetValues(H, adda, 1, &sxy_m, adda, 1, &sxy, &val, ADD_VALUES);
138:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy_m, &valconj, ADD_VALUES);

140:       sxy.d = 0; sxy_m.d = 1; /* spin 0, 1 */
141:       val   = -PETSC_i*uxy2; valconj = PetscConj(val);

143:       DMADDAMatSetValues(H, adda, 1, &sxy_m, adda, 1, &sxy, &val, ADD_VALUES);
144:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy_m, &valconj, ADD_VALUES);

146:       sxy.d = 1; sxy_m.d = 0; /* spin 1, 0 */
147:       val   = -PETSC_i*uxy2; valconj = PetscConj(val);

149:       DMADDAMatSetValues(H, adda, 1, &sxy_m, adda, 1, &sxy, &val, ADD_VALUES);
150:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy_m, &valconj, ADD_VALUES);

152:       sxy.d = 1; sxy_m.d = 1; /* spin 1, 1 */
153:       val   = PetscConj(uxy2); valconj = PetscConj(val);

155:       DMADDAMatSetValues(H, adda, 1, &sxy_m, adda, 1, &sxy, &val, ADD_VALUES);
156:       DMADDAMatSetValues(H, adda, 1, &sxy, adda, 1, &sxy_m, &valconj, ADD_VALUES);
157:     }
158:   }

160:   PetscFree(sxy.x);
161:   PetscFree(sxy_m.x);

163:   PetscFree(lcs);
164:   PetscFree(lce);

166:   MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY);
167:   MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY);

169:   /* scale H */
170:   MatScale(H, 1./(2.*h));

172:   /* construct normal equations */
173:   MatMatMult(H, H, MAT_INITIAL_MATRIX, 1., &HtH);

175:   PetscScalar mineval;
176:   computeMinEigVal(HtH, 1000, &mineval);
177:   PetscPrintf(PETSC_COMM_WORLD, "Minimum eigenvalue of H^{dag} H is %f\n", (double)PetscAbsScalar(mineval));

179:   /* permutation matrix to check whether H and HtH are identical to the ones in the paper */
180: /*   Mat perm; */
181: /*   ADDACreatematrix(adda, MATSEQAIJ, &perm); */
182: /*   PetscInt row, col; */
183: /*   PetscScalar one = 1.0; */
184: /*   for (PetscInt i=0; i<n; i++) { */
185: /*     for (PetscInt j=0; j<n; j++) { */
186: /*       row = (i*n+j)*2; col = i*n+j; */
187: /*       MatSetValues(perm, 1, &row, 1, &col, &one, INSERT_VALUES); */
188: /*       row = (i*n+j)*2+1; col = i*n+j + n*n; */
189: /*       MatSetValues(perm, 1, &row, 1, &col, &one, INSERT_VALUES); */
190: /*     } */
191: /*   } */
192: /*   MatAssemblyBegin(perm, MAT_FINAL_ASSEMBLY); */
193: /*   MatAssemblyEnd(perm, MAT_FINAL_ASSEMBLY); */

195: /*   Mat Hperm; */
196: /*   MatPtAP(H, perm, MAT_INITIAL_MATRIX, 1.0, &Hperm); */
197: /*   PetscPrintf(PETSC_COMM_WORLD, "Matrix H after construction\n"); */
198: /*   MatView(Hperm, PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD)); */

200: /*   Mat HtHperm; */
201: /*   MatPtAP(HtH, perm, MAT_INITIAL_MATRIX, 1.0, &HtHperm); */
202: /*   PetscPrintf(PETSC_COMM_WORLD, "Matrix HtH:\n"); */
203: /*   MatView(HtHperm, PETSC_VIEWER_STDOUT_(PETSC_COMM_WORLD)); */

205:   /* right hand side */
206:   DMCreateGlobalVector(adda, &b);
207:   VecSet(b,0.0);
208:   PetscInt    ix[1]   = {0};
209:   PetscScalar vals[1] = {1.0};
210:   VecSetValues(b, 1, ix, vals, INSERT_VALUES);
211:   VecAssemblyBegin(b);
212:   VecAssemblyEnd(b);
213: /*   VecSetRandom(b, rctx); */
214:   VecDuplicate(b, &Htb);
215:   MatMultTranspose(H, b, Htb);

217:   /* construct solver */
218:   KSPCreate(PETSC_COMM_WORLD,&kspmg);
219:   KSPSetType(kspmg, KSPCG);

221:   KSPGetPC(kspmg,&pcmg);
222:   PCSetType(pcmg,PCASA);

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

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

229:   PCSetDM(pcmg,adda);
230:   DMDestroy(&adda);

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

234:   VecDuplicate(b, &xvec);
235:   VecSet(xvec, 0.0);

237:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238:                       Solve the linear system
239:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

241:   KSPSolve(kspmg, Htb, xvec);

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

245:   KSPDestroy(&kspmg);

247:   VecDestroy(&xvec);
248:   VecDestroy(&b);
249:   VecDestroy(&Htb);
250:   MatDestroy(&H);
251:   MatDestroy(&HtH);

253:   PetscRandomDestroy(&rctx);
254:   PetscFinalize();
255:   return 0;
256: }

258: /* --------------------------------------------------------------------- */
261: PetscErrorCode computeMinEigVal(Mat A, PetscInt its, PetscScalar *eig)
262: {
264:   PetscRandom    rctx;      /* random number generator context */
265:   Vec            x0, x, x_1, tmp;
266:   PetscScalar    lambda_its, lambda_its_1;
267:   PetscReal      norm;
268:   Mat            G;
269:   PetscInt       i;

272:   PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
273:   PetscRandomSetFromOptions(rctx);

275:   /* compute G = I-1/norm(A)*A */
276:   MatNorm(A, NORM_1, &norm);
277:   MatConvert(A, MATSAME, MAT_INITIAL_MATRIX, &G);
278:   MatShift(G, -norm);
279:   MatScale(G, -1./norm);

281:   MatGetVecs(G, &x_1, &x);
282:   VecSetRandom(x, rctx);
283:   VecDuplicate(x, &x0);
284:   VecCopy(x, x0);

286:   MatMult(G, x, x_1);
287:   for (i=0; i<its; i++) {
288:     tmp  = x; x = x_1; x_1 = tmp;
289:     MatMult(G, x, x_1);
290:   }
291:   VecDot(x0, x, &lambda_its);
292:   VecDot(x0, x_1, &lambda_its_1);

294:   *eig = norm*(1.-lambda_its_1/lambda_its);

296:   VecDestroy(&x0);
297:   VecDestroy(&x);
298:   VecDestroy(&x_1);
299:   PetscRandomDestroy(&rctx);
300:   MatDestroy(&G);
301:   return(0);
302: }