Actual source code: ex30.c
petsc-3.11.4 2019-09-28
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>
5: int main(int argc,char **argv)
6: {
7: KSP solver;
8: PC pc;
9: Mat A,B;
10: Vec X,Y,Z;
11: MatScalar *a;
12: PetscScalar *b,*x,*y,*z;
13: PetscReal nrm;
14: PetscErrorCode ierr,size=8,lda=10, i,j;
16: PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
17: /* Create matrix and three vectors: these are all normal */
18: PetscMalloc1(lda*size,&b);
19: for (i=0; i<size; i++) {
20: for (j=0; j<size; j++) {
21: b[i+j*lda] = rand();
22: }
23: }
24: MatCreate(MPI_COMM_SELF,&A);
25: MatSetSizes(A,size,size,size,size);
26: MatSetType(A,MATSEQDENSE);
27: MatSeqDenseSetPreallocation(A,NULL);
29: MatDenseGetArray(A,&a);
30: for (i=0; i<size; i++) {
31: for (j=0; j<size; j++) {
32: a[i+j*size] = b[i+j*lda];
33: }
34: }
35: MatDenseRestoreArray(A,&a);
36: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
37: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
39: MatCreate(MPI_COMM_SELF,&B);
40: MatSetSizes(B,size,size,size,size);
41: MatSetType(B,MATSEQDENSE);
42: MatSeqDenseSetPreallocation(B,b);
43: MatSeqDenseSetLDA(B,lda);
44: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
45: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
47: PetscMalloc1(size,&x);
48: for (i=0; i<size; i++) x[i] = 1.0;
49: VecCreateSeqWithArray(MPI_COMM_SELF,1,size,x,&X);
50: VecAssemblyBegin(X);
51: VecAssemblyEnd(X);
53: PetscMalloc1(size,&y);
54: VecCreateSeqWithArray(MPI_COMM_SELF,1,size,y,&Y);
55: VecAssemblyBegin(Y);
56: VecAssemblyEnd(Y);
58: PetscMalloc1(size,&z);
59: VecCreateSeqWithArray(MPI_COMM_SELF,1,size,z,&Z);
60: VecAssemblyBegin(Z);
61: VecAssemblyEnd(Z);
63: /*
64: * Solve with A and B
65: */
66: KSPCreate(MPI_COMM_SELF,&solver);
67: KSPSetType(solver,KSPPREONLY);
68: KSPGetPC(solver,&pc);
69: PCSetType(pc,PCLU);
70: KSPSetOperators(solver,A,A);
71: KSPSolve(solver,X,Y);
72: KSPSetOperators(solver,B,B);
73: KSPSolve(solver,X,Z);
74: VecAXPY(Z,-1.0,Y);
75: VecNorm(Z,NORM_2,&nrm);
76: PetscPrintf(PETSC_COMM_SELF,"Test1; error norm=%e\n",nrm);
78: /* Free spaces */
79: PetscFree(b);
80: PetscFree(x);
81: PetscFree(y);
82: PetscFree(z);
83: VecDestroy(&X);
84: VecDestroy(&Y);
85: VecDestroy(&Z);
86: MatDestroy(&A);
87: MatDestroy(&B);
88: KSPDestroy(&solver);
90: PetscFinalize();
91: return ierr;
92: }
95: /*TEST
97: test:
99: TEST*/