xref: /petsc/src/mat/tutorials/ex9.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests MatCreateComposite()\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*T
5*c4762a1bSJed Brown    Concepts: Mat^composite matrices
6*c4762a1bSJed Brown    Processors: n
7*c4762a1bSJed Brown T*/
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown /*
10*c4762a1bSJed Brown   Include "petscmat.h" so that we can use matrices.
11*c4762a1bSJed Brown   automatically includes:
12*c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h    - vectors
13*c4762a1bSJed Brown      petscmat.h    - matrices
14*c4762a1bSJed Brown      petscis.h     - index sets            petscviewer.h - viewers
15*c4762a1bSJed Brown */
16*c4762a1bSJed Brown #include <petscmat.h>
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown int main(int argc,char **args)
19*c4762a1bSJed Brown {
20*c4762a1bSJed Brown   Mat              *A,B;           /* matrix */
21*c4762a1bSJed Brown   PetscErrorCode   ierr;
22*c4762a1bSJed Brown   Vec              x,y,v,v2,z,z2;
23*c4762a1bSJed Brown   PetscReal        rnorm;
24*c4762a1bSJed Brown   PetscInt         n = 20;         /* size of the matrix */
25*c4762a1bSJed Brown   PetscInt         nmat = 3;       /* number of matrices */
26*c4762a1bSJed Brown   PetscInt         i;
27*c4762a1bSJed Brown   PetscRandom      rctx;
28*c4762a1bSJed Brown   MatCompositeType type;
29*c4762a1bSJed Brown   PetscScalar      scalings[5]={2,3,4,5,6};
30*c4762a1bSJed Brown 
31*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
32*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nmat",&nmat,NULL);CHKERRQ(ierr);
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   /*
36*c4762a1bSJed Brown      Create random matrices
37*c4762a1bSJed Brown   */
38*c4762a1bSJed Brown   ierr = PetscMalloc1(nmat+3,&A);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rctx);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n/2,3,NULL,3,NULL,&A[0]);CHKERRQ(ierr);
41*c4762a1bSJed Brown   for (i = 1; i < nmat+1; i++) {
42*c4762a1bSJed Brown     ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,3,NULL,3,NULL,&A[i]);CHKERRQ(ierr);
43*c4762a1bSJed Brown   }
44*c4762a1bSJed Brown   ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n/2,n,3,NULL,3,NULL,&A[nmat+1]);CHKERRQ(ierr);
45*c4762a1bSJed Brown   for (i = 0; i < nmat+2; i++) {
46*c4762a1bSJed Brown     ierr = MatSetRandom(A[i],rctx);CHKERRQ(ierr);
47*c4762a1bSJed Brown   }
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown   ierr = MatCreateVecs(A[1],&x,&y);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = VecDuplicate(y,&z);CHKERRQ(ierr);
51*c4762a1bSJed Brown   ierr = VecDuplicate(z,&z2);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = MatCreateVecs(A[0],&v,NULL);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = VecDuplicate(v,&v2);CHKERRQ(ierr);
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown   /* Test MatMult of an ADDITIVE MatComposite B made up of A[1],A[2],A[3] with separate scalings */
56*c4762a1bSJed Brown 
57*c4762a1bSJed Brown   /* Do MatMult with A[1],A[2],A[3] by hand and store the result in z */
58*c4762a1bSJed Brown   ierr = VecSet(x,1.0);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = MatMult(A[1],x,z);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = VecScale(z,scalings[1]);CHKERRQ(ierr);
61*c4762a1bSJed Brown   for (i = 2; i < nmat+1; i++) {
62*c4762a1bSJed Brown     ierr = MatMult(A[i],x,z2);CHKERRQ(ierr);
63*c4762a1bSJed Brown     ierr = VecAXPY(z,scalings[i],z2);CHKERRQ(ierr);
64*c4762a1bSJed Brown   }
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown   /* Do MatMult using MatComposite and store the result in y */
67*c4762a1bSJed Brown   ierr = VecSet(y,0.0);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = MatCreateComposite(PETSC_COMM_WORLD,nmat,A+1,&B);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = MatCompositeSetScalings(B,&scalings[1]);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = MatMultAdd(B,x,y,y);CHKERRQ(ierr);
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   /* Diff y and z */
74*c4762a1bSJed Brown   ierr = VecAXPY(y,-1.0,z);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = VecNorm(y,NORM_2,&rnorm);CHKERRQ(ierr);
76*c4762a1bSJed Brown   if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) {
77*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with composite add %g\n",(double)rnorm);CHKERRQ(ierr);
78*c4762a1bSJed Brown   }
79*c4762a1bSJed Brown 
80*c4762a1bSJed Brown   /* Test MatCompositeMerge on ADDITIVE MatComposite */
81*c4762a1bSJed Brown   ierr = MatCompositeSetMatStructure(B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); /* default */
82*c4762a1bSJed Brown   ierr = MatCompositeMerge(B);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = MatMult(B,x,y);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = VecAXPY(y,-1.0,z);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = VecNorm(y,NORM_2,&rnorm);CHKERRQ(ierr);
87*c4762a1bSJed Brown   if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) {
88*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with composite add after merge %g\n",(double)rnorm);CHKERRQ(ierr);
89*c4762a1bSJed Brown   }
90*c4762a1bSJed Brown 
91*c4762a1bSJed Brown   /*
92*c4762a1bSJed Brown      Test n x n/2 multiplicative composite B made up of A[0],A[1],A[2] with separate scalings
93*c4762a1bSJed Brown   */
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   /* Do MatMult with A[0],A[1],A[2] by hand and store the result in z */
96*c4762a1bSJed Brown   ierr = VecSet(v,1.0);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = MatMult(A[0],v,z);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = VecScale(z,scalings[0]);CHKERRQ(ierr);
99*c4762a1bSJed Brown   for (i = 1; i < nmat; i++) {
100*c4762a1bSJed Brown     ierr = MatMult(A[i],z,y);CHKERRQ(ierr);
101*c4762a1bSJed Brown     ierr = VecScale(y,scalings[i]);CHKERRQ(ierr);
102*c4762a1bSJed Brown     ierr = VecCopy(y,z);CHKERRQ(ierr);
103*c4762a1bSJed Brown   }
104*c4762a1bSJed Brown 
105*c4762a1bSJed Brown   /* Do MatMult using MatComposite and store the result in y */
106*c4762a1bSJed Brown   ierr = MatCreateComposite(PETSC_COMM_WORLD,nmat,A,&B);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = MatCompositeSetType(B,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = MatCompositeSetMergeType(B,MAT_COMPOSITE_MERGE_LEFT);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = MatCompositeSetScalings(B,&scalings[0]);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); /* do MatCompositeMerge() if -mat_composite_merge 1 */
113*c4762a1bSJed Brown   ierr = MatMult(B,v,y);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   /* Diff y and z */
117*c4762a1bSJed Brown   ierr = VecAXPY(y,-1.0,z);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = VecNorm(y,NORM_2,&rnorm);CHKERRQ(ierr);
119*c4762a1bSJed Brown   if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) {
120*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with composite multiplicative %g\n",(double)rnorm);CHKERRQ(ierr);
121*c4762a1bSJed Brown   }
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown   /*
124*c4762a1bSJed Brown      Test n/2 x n multiplicative composite B made up of A[2], A[3], A[4] without separate scalings
125*c4762a1bSJed Brown   */
126*c4762a1bSJed Brown   ierr = VecSet(x,1.0);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = MatMult(A[2],x,z);CHKERRQ(ierr);
128*c4762a1bSJed Brown   for (i = 3; i < nmat+1; i++) {
129*c4762a1bSJed Brown     ierr = MatMult(A[i],z,y);CHKERRQ(ierr);
130*c4762a1bSJed Brown     ierr = VecCopy(y,z);CHKERRQ(ierr);
131*c4762a1bSJed Brown   }
132*c4762a1bSJed Brown   ierr = MatMult(A[nmat+1],z,v);CHKERRQ(ierr);
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown   ierr = MatCreateComposite(PETSC_COMM_WORLD,nmat,A+2,&B);CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr = MatCompositeSetType(B,MAT_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
136*c4762a1bSJed Brown   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); /* do MatCompositeMerge() if -mat_composite_merge 1 */
139*c4762a1bSJed Brown   ierr = MatMult(B,x,v2);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown   ierr = VecAXPY(v2,-1.0,v);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = VecNorm(v2,NORM_2,&rnorm);CHKERRQ(ierr);
144*c4762a1bSJed Brown   if (rnorm > 10000.0*PETSC_MACHINE_EPSILON) {
145*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with composite multiplicative %g\n",(double)rnorm);CHKERRQ(ierr);
146*c4762a1bSJed Brown   }
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown   /*
149*c4762a1bSJed Brown      Test get functions
150*c4762a1bSJed Brown   */
151*c4762a1bSJed Brown   ierr = MatCreateComposite(PETSC_COMM_WORLD,nmat,A,&B);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = MatCompositeGetNumberMat(B,&n);CHKERRQ(ierr);
153*c4762a1bSJed Brown   if (nmat != n) {
154*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with GetNumberMat %D != %D\n",nmat,n);CHKERRQ(ierr);
155*c4762a1bSJed Brown   }
156*c4762a1bSJed Brown   ierr = MatCompositeGetMat(B,0,&A[nmat+2]);CHKERRQ(ierr);
157*c4762a1bSJed Brown   if (A[0] != A[nmat+2]) {
158*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with GetMat\n");CHKERRQ(ierr);
159*c4762a1bSJed Brown   }
160*c4762a1bSJed Brown   ierr = MatCompositeGetType(B,&type);CHKERRQ(ierr);
161*c4762a1bSJed Brown   if (type != MAT_COMPOSITE_ADDITIVE) {
162*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Error with GetType\n");CHKERRQ(ierr);
163*c4762a1bSJed Brown   }
164*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown   /*
167*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
168*c4762a1bSJed Brown      are no longer needed.
169*c4762a1bSJed Brown   */
170*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
171*c4762a1bSJed Brown   ierr = VecDestroy(&y);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = VecDestroy(&v);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = VecDestroy(&v2);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = VecDestroy(&z);CHKERRQ(ierr);
175*c4762a1bSJed Brown   ierr = VecDestroy(&z2);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rctx);CHKERRQ(ierr);
177*c4762a1bSJed Brown   for (i = 0; i < nmat+2; i++) {
178*c4762a1bSJed Brown     ierr = MatDestroy(&A[i]);CHKERRQ(ierr);
179*c4762a1bSJed Brown   }
180*c4762a1bSJed Brown   ierr = PetscFree(A);CHKERRQ(ierr);
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown   ierr = PetscFinalize();
183*c4762a1bSJed Brown   return ierr;
184*c4762a1bSJed Brown }
185*c4762a1bSJed Brown 
186*c4762a1bSJed Brown /*TEST
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown    test:
189*c4762a1bSJed Brown       nsize: 2
190*c4762a1bSJed Brown       requires: double
191*c4762a1bSJed Brown       args: -mat_composite_merge {{0 1}shared output} -mat_composite_merge_mvctx {{0 1}shared output}
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown TEST*/
194