Actual source code: ex30.c
petsc-3.7.7 2017-09-25
1: static char help[] = " * Example code testing SeqDense matrices with an LDA (leading dimension of the user-allocated arrray) larger than M.\n";
3: #include <petscksp.h>
7: int main(int argc,char **argv)
8: {
9: KSP solver;
10: PC pc;
11: Mat A,B;
12: Vec X,Y,Z;
13: MatScalar *a;
14: PetscScalar *b,*x,*y,*z;
15: PetscReal nrm;
16: PetscErrorCode ierr,size=8,lda=10, i,j;
18: PetscInitialize(&argc,&argv,0,help);
19: /* Create matrix and three vectors: these are all normal */
20: PetscMalloc1(lda*size,&b);
21: for (i=0; i<size; i++) {
22: for (j=0; j<size; j++) {
23: b[i+j*lda] = rand();
24: }
25: }
26: MatCreate(MPI_COMM_SELF,&A);
27: MatSetSizes(A,size,size,size,size);
28: MatSetType(A,MATSEQDENSE);
29: MatSeqDenseSetPreallocation(A,NULL);
31: MatDenseGetArray(A,&a);
32: for (i=0; i<size; i++) {
33: for (j=0; j<size; j++) {
34: a[i+j*size] = b[i+j*lda];
35: }
36: }
37: MatDenseRestoreArray(A,&a);
38: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
39: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
41: MatCreate(MPI_COMM_SELF,&B);
42: MatSetSizes(B,size,size,size,size);
43: MatSetType(B,MATSEQDENSE);
44: MatSeqDenseSetPreallocation(B,b);
45: MatSeqDenseSetLDA(B,lda);
46: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
47: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
49: PetscMalloc1(size,&x);
50: for (i=0; i<size; i++) x[i] = 1.0;
51: VecCreateSeqWithArray(MPI_COMM_SELF,1,size,x,&X);
52: VecAssemblyBegin(X);
53: VecAssemblyEnd(X);
55: PetscMalloc1(size,&y);
56: VecCreateSeqWithArray(MPI_COMM_SELF,1,size,y,&Y);
57: VecAssemblyBegin(Y);
58: VecAssemblyEnd(Y);
60: PetscMalloc1(size,&z);
61: VecCreateSeqWithArray(MPI_COMM_SELF,1,size,z,&Z);
62: VecAssemblyBegin(Z);
63: VecAssemblyEnd(Z);
65: /*
66: * Solve with A and B
67: */
68: KSPCreate(MPI_COMM_SELF,&solver);
69: KSPSetType(solver,KSPPREONLY);
70: KSPGetPC(solver,&pc);
71: PCSetType(pc,PCLU);
72: KSPSetOperators(solver,A,A);
73: KSPSolve(solver,X,Y);
74: KSPSetOperators(solver,B,B);
75: KSPSolve(solver,X,Z);
76: VecAXPY(Z,-1.0,Y);
77: VecNorm(Z,NORM_2,&nrm);
78: printf("Test1; error norm=%e\n",nrm);
80: /* Free spaces */
81: PetscFree(b);
82: PetscFree(x);
83: PetscFree(y);
84: PetscFree(z);
85: VecDestroy(&X);
86: VecDestroy(&Y);
87: VecDestroy(&Z);
88: MatDestroy(&A);
89: MatDestroy(&B);
90: KSPDestroy(&solver);
92: PetscFinalize();
93: return 0;
94: }