Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test MatCreateTile.\n\n";
12 :
13 : #include <slepcsys.h>
14 :
15 3 : int main(int argc,char **argv)
16 : {
17 3 : Mat T,E,A;
18 3 : PetscInt i,Istart,Iend,n=10;
19 :
20 3 : PetscFunctionBeginUser;
21 3 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
22 3 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
23 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatCreateTile test, n=%" PetscInt_FMT "\n",n));
24 :
25 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26 : Create T=tridiag([-1 2 -1],n,n) and E=eye(n)
27 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
28 :
29 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&T));
30 3 : PetscCall(MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,n,n));
31 3 : PetscCall(MatSetFromOptions(T));
32 :
33 3 : PetscCall(MatGetOwnershipRange(T,&Istart,&Iend));
34 23 : for (i=Istart;i<Iend;i++) {
35 20 : if (i>0) PetscCall(MatSetValue(T,i,i-1,-1.0,INSERT_VALUES));
36 20 : if (i<n-1) PetscCall(MatSetValue(T,i,i+1,-1.0,INSERT_VALUES));
37 20 : PetscCall(MatSetValue(T,i,i,2.0,INSERT_VALUES));
38 : }
39 3 : PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY));
40 3 : PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY));
41 :
42 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&E));
43 3 : PetscCall(MatSetSizes(E,PETSC_DECIDE,PETSC_DECIDE,n,n));
44 3 : PetscCall(MatSetFromOptions(E));
45 :
46 3 : PetscCall(MatGetOwnershipRange(E,&Istart,&Iend));
47 23 : for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(E,i,i,1.0,INSERT_VALUES));
48 3 : PetscCall(MatAssemblyBegin(E,MAT_FINAL_ASSEMBLY));
49 3 : PetscCall(MatAssemblyEnd(E,MAT_FINAL_ASSEMBLY));
50 :
51 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52 : Create tiled matrix A = [ 2*T -E; 0 3*T ]
53 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54 3 : PetscCall(MatCreateTile(2.0,T,-1.0,E,0.0,E,3.0,T,&A));
55 3 : PetscCall(MatView(A,NULL));
56 :
57 3 : PetscCall(MatDestroy(&T));
58 3 : PetscCall(MatDestroy(&E));
59 3 : PetscCall(MatDestroy(&A));
60 3 : PetscCall(SlepcFinalize());
61 : return 0;
62 : }
63 :
64 : /*TEST
65 :
66 : test:
67 : suffix: 1
68 : nsize: 1
69 :
70 : test:
71 : suffix: 2
72 : nsize: 2
73 :
74 : TEST*/
|