Actual source code: test1.c

slepc-3.19.2 2023-09-05
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Test MatCreateTile.\n\n";

 13: #include <slepcsys.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            T,E,A;
 18:   PetscInt       i,Istart,Iend,n=10;

 20:   PetscFunctionBeginUser;
 21:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 22:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 23:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatCreateTile test, n=%" PetscInt_FMT "\n",n));

 25:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 26:      Create T=tridiag([-1 2 -1],n,n) and E=eye(n)
 27:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 29:   PetscCall(MatCreate(PETSC_COMM_WORLD,&T));
 30:   PetscCall(MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,n,n));
 31:   PetscCall(MatSetFromOptions(T));
 32:   PetscCall(MatSetUp(T));

 34:   PetscCall(MatGetOwnershipRange(T,&Istart,&Iend));
 35:   for (i=Istart;i<Iend;i++) {
 36:     if (i>0) PetscCall(MatSetValue(T,i,i-1,-1.0,INSERT_VALUES));
 37:     if (i<n-1) PetscCall(MatSetValue(T,i,i+1,-1.0,INSERT_VALUES));
 38:     PetscCall(MatSetValue(T,i,i,2.0,INSERT_VALUES));
 39:   }
 40:   PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY));
 41:   PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY));

 43:   PetscCall(MatCreate(PETSC_COMM_WORLD,&E));
 44:   PetscCall(MatSetSizes(E,PETSC_DECIDE,PETSC_DECIDE,n,n));
 45:   PetscCall(MatSetFromOptions(E));
 46:   PetscCall(MatSetUp(E));

 48:   PetscCall(MatGetOwnershipRange(E,&Istart,&Iend));
 49:   for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(E,i,i,1.0,INSERT_VALUES));
 50:   PetscCall(MatAssemblyBegin(E,MAT_FINAL_ASSEMBLY));
 51:   PetscCall(MatAssemblyEnd(E,MAT_FINAL_ASSEMBLY));

 53:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 54:      Create tiled matrix A = [ 2*T -E; 0 3*T ]
 55:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 56:   PetscCall(MatCreateTile(2.0,T,-1.0,E,0.0,E,3.0,T,&A));
 57:   PetscCall(MatView(A,NULL));

 59:   PetscCall(MatDestroy(&T));
 60:   PetscCall(MatDestroy(&E));
 61:   PetscCall(MatDestroy(&A));
 62:   PetscCall(SlepcFinalize());
 63:   return 0;
 64: }

 66: /*TEST

 68:    test:
 69:       suffix: 1
 70:       nsize: 1

 72:    test:
 73:       suffix: 2
 74:       nsize: 2

 76: TEST*/