Actual source code: ex2.c
1: /*
2: * Example code testing SeqDense matrices with an LDA (leading dimension
3: * of the user-allocated arrray) larger than M.
4: * This example tests the functionality of MatInsertValues, MatMult,
5: * and MatMultTranspose.
6: */
7: #include <stdlib.h>
8: #include petscmat.h
12: int main(int argc,char **argv)
13: {
14: Mat A,A11,A12,A21,A22;
15: Vec X,X1,X2,Y,Z,Z1,Z2;
16: PetscScalar *a,*b,*x,*y,*z,v,one=1,mone=-1;
17: PetscReal nrm;
19: PetscInt size=8,size1=6,size2=2, i,j;
21: PetscInitialize(&argc,&argv,0,0);
23: /*
24: * Create matrix and three vectors: these are all normal
25: */
26: PetscMalloc(size*size*sizeof(PetscScalar),&a);
27: PetscMalloc(size*size*sizeof(PetscScalar),&b);
28: for (i=0; i<size; i++) {
29: for (j=0; j<size; j++) {
30: a[i+j*size] = rand(); b[i+j*size] = a[i+j*size];
31: }
32: }
33: MatCreate(MPI_COMM_SELF,size,size,size,size,&A);
34: MatSetType(A,MATSEQDENSE);
35: MatSeqDenseSetPreallocation(A,a);
36: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
37: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
39: PetscMalloc(size*sizeof(PetscScalar),&x);
40: for (i=0; i<size; i++) {
41: x[i] = one;
42: }
43: VecCreateSeqWithArray(MPI_COMM_SELF,size,x,&X);
44: VecAssemblyBegin(X);
45: VecAssemblyEnd(X);
47: PetscMalloc(size*sizeof(PetscScalar),&y);
48: VecCreateSeqWithArray(MPI_COMM_SELF,size,y,&Y);
49: VecAssemblyBegin(Y);
50: VecAssemblyEnd(Y);
52: PetscMalloc(size*sizeof(PetscScalar),&z);
53: VecCreateSeqWithArray(MPI_COMM_SELF,size,z,&Z);
54: VecAssemblyBegin(Z);
55: VecAssemblyEnd(Z);
57: /*
58: * Now create submatrices and subvectors
59: */
60: MatCreate(MPI_COMM_SELF,size1,size1,size1,size1,&A11);
61: MatSetType(A11,MATSEQDENSE);
62: MatSeqDenseSetPreallocation(A11,b);
63: MatSeqDenseSetLDA(A11,size);
64: MatAssemblyBegin(A11,MAT_FINAL_ASSEMBLY);
65: MatAssemblyEnd(A11,MAT_FINAL_ASSEMBLY);
67: MatCreate(MPI_COMM_SELF,size1,size2,size1,size2,&A12);
68: MatSetType(A12,MATSEQDENSE);
69: MatSeqDenseSetPreallocation(A12,b+size1*size);
70: MatSeqDenseSetLDA(A12,size);
71: MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY);
72: MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY);
74: MatCreate(MPI_COMM_SELF,size2,size1,size2,size1,&A21);
75: MatSetType(A21,MATSEQDENSE);
76: MatSeqDenseSetPreallocation(A21,b+size1);
77: MatSeqDenseSetLDA(A21,size);
78: MatAssemblyBegin(A21,MAT_FINAL_ASSEMBLY);
79: MatAssemblyEnd(A21,MAT_FINAL_ASSEMBLY);
81: MatCreate(MPI_COMM_SELF,size2,size2,size2,size2,&A22);
82: MatSetType(A22,MATSEQDENSE);
83: MatSeqDenseSetPreallocation(A22,b+size1*size+size1);
84: MatSeqDenseSetLDA(A22,size);
85: MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY);
86: MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY);
88: VecCreateSeqWithArray(MPI_COMM_SELF,size1,x,&X1);
89: VecCreateSeqWithArray(MPI_COMM_SELF,size2,x+size1,&X2);
90: VecCreateSeqWithArray(MPI_COMM_SELF,size1,z,&Z1);
91: VecCreateSeqWithArray(MPI_COMM_SELF,size2,z+size1,&Z2);
93: /*
94: * Now multiple matrix times input in two ways;
95: * compare the results
96: */
97: MatMult(A,X,Y);
98: MatMult(A11,X1,Z1);
99: MatMultAdd(A12,X2,Z1,Z1);
100: MatMult(A22,X2,Z2);
101: MatMultAdd(A21,X1,Z2,Z2);
102: VecAXPY(&mone,Y,Z);
103: VecNorm(Z,NORM_2,&nrm);
104: printf("Test1; error norm=%e\n",nrm);
105: /*
106: printf("MatMult the usual way:\n"); VecView(Y,0);
107: printf("MatMult by subblock:\n"); VecView(Z,0);
108: */
110: /*
111: * Next test: change both matrices
112: */
113: v = rand(); i=1; j=size-2;
114: MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);
115: j -= size1;
116: MatSetValues(A12,1,&i,1,&j,&v,INSERT_VALUES);
117: v = rand(); i=j=size1+1;
118: MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);
119: i=j=1;
120: MatSetValues(A22,1,&i,1,&j,&v,INSERT_VALUES);
121: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
122: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
123: MatAssemblyBegin(A12,MAT_FINAL_ASSEMBLY);
124: MatAssemblyEnd(A12,MAT_FINAL_ASSEMBLY);
125: MatAssemblyBegin(A22,MAT_FINAL_ASSEMBLY);
126: MatAssemblyEnd(A22,MAT_FINAL_ASSEMBLY);
128: MatMult(A,X,Y);
129: MatMult(A11,X1,Z1);
130: MatMultAdd(A12,X2,Z1,Z1);
131: MatMult(A22,X2,Z2);
132: MatMultAdd(A21,X1,Z2,Z2);
133: VecAXPY(&mone,Y,Z);
134: VecNorm(Z,NORM_2,&nrm);
135: printf("Test2; error norm=%e\n",nrm);
137: /*
138: * Transpose product
139: */
140: MatMultTranspose(A,X,Y);
141: MatMultTranspose(A11,X1,Z1);
142: MatMultTransposeAdd(A21,X2,Z1,Z1);
143: MatMultTranspose(A22,X2,Z2);
144: MatMultTransposeAdd(A12,X1,Z2,Z2);
145: VecAXPY(&mone,Y,Z);
146: VecNorm(Z,NORM_2,&nrm);
147: printf("Test3; error norm=%e\n",nrm);
149: PetscFinalize();
150: return 0;
151: }