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