xref: /petsc/src/mat/tests/ex241f.F90 (revision ccfd86f17c20321558100f6af55b03dc7cd752d2)
1c4762a1bSJed Brown!     Test code contributed by Thibaut Appel <t.appel17@imperial.ac.uk>
2c4762a1bSJed Brown
3c4762a1bSJed Brown  program test_assembly
4c4762a1bSJed Brown
5c4762a1bSJed Brown#include <petsc/finclude/petscmat.h>
6c4762a1bSJed Brown
7c4762a1bSJed Brown  use PetscMat
8d33816bfSBarry Smith  use ISO_Fortran_Env, only : real64
9c4762a1bSJed Brown
10c4762a1bSJed Brown  implicit none
11c4762a1bSJed Brown  PetscInt,    parameter :: wp = real64, n = 10
12c4762a1bSJed Brown  PetscScalar, parameter :: zero = 0.0, one = 1.0
13c4762a1bSJed Brown  Mat      :: L
14c4762a1bSJed Brown  PetscInt :: istart, iend, row, i1, i0
15c4762a1bSJed Brown  PetscErrorCode :: ierr
16c4762a1bSJed Brown
17c4762a1bSJed Brown  PetscInt    cols(1),rows(1)
18c4762a1bSJed Brown  PetscScalar  vals(1)
19c4762a1bSJed Brown
20d8606c27SBarry Smith  PetscCallA(PetscInitialize(ierr))
21c4762a1bSJed Brown
22c4762a1bSJed Brown  i0 = 0
23c4762a1bSJed Brown  i1 = 1
24c4762a1bSJed Brown
25d8606c27SBarry Smith  PetscCallA(MatCreate(PETSC_COMM_WORLD,L,ierr))
26d8606c27SBarry Smith  PetscCallA(MatSetType(L,MATAIJ,ierr))
27d8606c27SBarry Smith  PetscCallA(MatSetSizes(L,PETSC_DECIDE,PETSC_DECIDE,n,n,ierr))
28c4762a1bSJed Brown
295d83a8b1SBarry Smith  PetscCallA(MatSeqAIJSetPreallocation(L,i1,PETSC_NULL_INTEGER_ARRAY,ierr))
305d83a8b1SBarry Smith  PetscCallA(MatMPIAIJSetPreallocation(L,i1,PETSC_NULL_INTEGER_ARRAY,i0,PETSC_NULL_INTEGER_ARRAY,ierr)) ! No allocated non-zero in off-diagonal part
31d8606c27SBarry Smith  PetscCallA(MatSetOption(L,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE,ierr))
32d8606c27SBarry Smith  PetscCallA(MatSetOption(L,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE,ierr))
33d8606c27SBarry Smith  PetscCallA(MatSetOption(L,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE,ierr))
34c4762a1bSJed Brown
35d8606c27SBarry Smith  PetscCallA(MatGetOwnershipRange(L,istart,iend,ierr))
36c4762a1bSJed Brown
37c4762a1bSJed Brown  ! assembling a diagonal matrix
38c4762a1bSJed Brown  do row = istart,iend-1
39c4762a1bSJed Brown
40*ccfd86f1SBarry Smith    cols = [row]; vals = [one]; rows = [row]
41d8606c27SBarry Smith    PetscCallA(MatSetValues(L,i1,rows,i1,cols,vals,ADD_VALUES,ierr))
42c4762a1bSJed Brown
43c4762a1bSJed Brown  end do
44c4762a1bSJed Brown
45d8606c27SBarry Smith  PetscCallA(MatAssemblyBegin(L,MAT_FINAL_ASSEMBLY,ierr))
46d8606c27SBarry Smith  PetscCallA(MatAssemblyEnd(L,MAT_FINAL_ASSEMBLY,ierr))
47c4762a1bSJed Brown
48d8606c27SBarry Smith  PetscCallA(MatSetOption(L,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE,ierr))
49c4762a1bSJed Brown
50d8606c27SBarry Smith  !PetscCallA(MatZeroEntries(L,ierr))
51c4762a1bSJed Brown
52c4762a1bSJed Brown  ! assembling a diagonal matrix, adding a zero value to non-diagonal part
53c4762a1bSJed Brown  do row = istart,iend-1
54c4762a1bSJed Brown
55c4762a1bSJed Brown    if (row == 0) then
56c4762a1bSJed Brown      cols = [n-1]
57c4762a1bSJed Brown      vals = [zero]
58c4762a1bSJed Brown      rows = [row]
59d8606c27SBarry Smith      PetscCallA(MatSetValues(L,i1,rows,i1,cols,vals,ADD_VALUES,ierr))
60c4762a1bSJed Brown    end if
61*ccfd86f1SBarry Smith    cols = [row]; vals = [one] ; rows = [ row]
62d8606c27SBarry Smith    PetscCallA(MatSetValues(L,i1,rows,i1,cols,vals,ADD_VALUES,ierr))
63c4762a1bSJed Brown
64c4762a1bSJed Brown  end do
65c4762a1bSJed Brown
66d8606c27SBarry Smith  PetscCallA(MatAssemblyBegin(L,MAT_FINAL_ASSEMBLY,ierr))
67d8606c27SBarry Smith  PetscCallA(MatAssemblyEnd(L,MAT_FINAL_ASSEMBLY,ierr))
68d8606c27SBarry Smith  PetscCallA(MatDestroy(L,ierr))
69c4762a1bSJed Brown
70d8606c27SBarry Smith  PetscCallA(PetscFinalize(ierr))
71c4762a1bSJed Brown
72c4762a1bSJed Brownend program test_assembly
73c4762a1bSJed Brown
74c4762a1bSJed Brown!/*TEST
75c4762a1bSJed Brown!
76c4762a1bSJed Brown!   build:
77c4762a1bSJed Brown!      requires: complex
78c4762a1bSJed Brown!
79c4762a1bSJed Brown!   test:
80c4762a1bSJed Brown!      nsize: 2
81c4762a1bSJed Brown!
82c4762a1bSJed Brown!TEST*/
83