xref: /petsc/src/mat/tutorials/ex4.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown /*
2*c4762a1bSJed Brown  *
3*c4762a1bSJed Brown  *  Created on: Sep 25, 2017
4*c4762a1bSJed Brown  *      Author: Fande Kong
5*c4762a1bSJed Brown  */
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown static char help[] = "Illustrate the use of MatResetPreallocation.\n";
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown #include "petscmat.h"
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown int main(int argc,char **argv)
12*c4762a1bSJed Brown {
13*c4762a1bSJed Brown   Mat             A;
14*c4762a1bSJed Brown   MPI_Comm        comm;
15*c4762a1bSJed Brown   PetscInt        n=5,m=5,*dnnz,*onnz,i,rstart,rend,M,N;
16*c4762a1bSJed Brown   PetscErrorCode  ierr;
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
19*c4762a1bSJed Brown   comm = MPI_COMM_WORLD;
20*c4762a1bSJed Brown   ierr = PetscMalloc2(m,&dnnz,m,&onnz);CHKERRQ(ierr);
21*c4762a1bSJed Brown   for (i=0; i<m; i++) {
22*c4762a1bSJed Brown     dnnz[i] = 1;
23*c4762a1bSJed Brown     onnz[i] = 1;
24*c4762a1bSJed Brown   }
25*c4762a1bSJed Brown   ierr = MatCreateAIJ(comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,PETSC_DECIDE,dnnz,PETSC_DECIDE,onnz,&A);CHKERRQ(ierr);
26*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown   /* This assembly shrinks memory because we do not insert enough number of values */
31*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown   /* MatResetPreallocation restores the memory required by users */
35*c4762a1bSJed Brown   ierr = MatResetPreallocation(A);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
39*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
40*c4762a1bSJed Brown     ierr = MatSetValue(A,i,i,2.0,INSERT_VALUES);CHKERRQ(ierr);
41*c4762a1bSJed Brown     if (rend<N) {
42*c4762a1bSJed Brown       ierr = MatSetValue(A,i,rend,1.0,INSERT_VALUES);CHKERRQ(ierr);
43*c4762a1bSJed Brown     }
44*c4762a1bSJed Brown   }
45*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = PetscFinalize();
50*c4762a1bSJed Brown   return ierr;
51*c4762a1bSJed Brown }
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown /*TEST
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown    test:
56*c4762a1bSJed Brown       suffix: 1
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown    test:
59*c4762a1bSJed Brown       suffix: 2
60*c4762a1bSJed Brown       nsize: 2
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown TEST*/
63