Actual source code: test14f.F90

slepc-3.22.1 2024-10-28
Report Typos and Errors
  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: !  Program usage: mpiexec -n <np> ./test14f [-help] [-n <n>] [all SLEPc options]
 11: !
 12: !  Description: Simple example that tests solving a DSNHEP problem.
 13: !
 14: !  The command line options are:
 15: !    -n <n>, where <n> = matrix size
 16: !
 17: ! ----------------------------------------------------------------------
 18: !
 19:       program main
 20: #include <slepc/finclude/slepcds.h>
 21:       use slepcds
 22:       implicit none

 24: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 25: !     Declarations
 26: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 27: !
 28: !  Variables:
 29: !     A     problem matrix
 30: !     ds    dense solver context

 32:       Mat            A
 33:       DS             ds
 34:       PetscInt       n, i, ld, zero
 35:       PetscMPIInt    rank
 36:       PetscErrorCode ierr
 37:       PetscBool      flg
 38:       PetscScalar    wr(100), wi(100)
 39:       PetscReal      re, im
 40:       PetscScalar, pointer :: aa(:,:)

 42: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 43: !     Beginning of program
 44: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 46:       zero = 0
 47:       PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER,ierr))
 48:       if (ierr .ne. 0) then
 49:         print*,'SlepcInitialize failed'
 50:         stop
 51:       endif
 52:       PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
 53:       n = 10
 54:       PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))
 55:       if (n .gt. 100) then; SETERRA(PETSC_COMM_SELF,1,'Program currently limited to n=100'); endif

 57:       if (rank .eq. 0) then
 58:         write(*,110) n
 59:       endif
 60:  110  format (/'Solve a Dense System of type NHEP, n =',I3,' (Fortran)')

 62: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63: !     Create DS object
 64: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 66:       PetscCallA(DSCreate(PETSC_COMM_WORLD,ds,ierr))
 67:       PetscCallA(DSSetType(ds,DSNHEP,ierr))
 68:       PetscCallA(DSSetFromOptions(ds,ierr))
 69:       ld = n
 70:       PetscCallA(DSAllocate(ds,ld,ierr))
 71:       PetscCallA(DSSetDimensions(ds,n,zero,zero,ierr))

 73: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 74: !     Fill with Grcar matrix
 75: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 77:       PetscCallA(DSGetMat(ds,DS_MAT_A,A,ierr))
 78:       PetscCallA(MatDenseGetArrayF90(A,aa,ierr))
 79:       call FillUpMatrix(n,aa)
 80:       PetscCallA(MatDenseRestoreArrayF90(A,aa,ierr))
 81:       PetscCallA(DSRestoreMat(ds,DS_MAT_A,A,ierr))
 82:       PetscCallA(DSSetState(ds,DS_STATE_INTERMEDIATE,ierr))

 84: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85: !     Solve the problem and show eigenvalues
 86: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 88:       PetscCallA(DSSolve(ds,wr,wi,ierr))
 89: !     PetscCallA(DSSort(ds,wr,wi,PETSC_NULL_SCALAR,PETSC_NULL_SCALAR,PETSC_NULL_INTEGER,ierr))

 91:       if (rank .eq. 0) then
 92:         write(*,*) 'Computed eigenvalues ='
 93:         do i=1,n
 94: #if defined(PETSC_USE_COMPLEX)
 95:           re = PetscRealPart(wr(i))
 96:           im = PetscImaginaryPart(wr(i))
 97: #else
 98:           re = wr(i)
 99:           im = wi(i)
100: #endif
101:           if (abs(im).lt.1.d-10) then
102:             write(*,120) re
103:           else
104:             write(*,130) re, im
105:           endif
106:         end do
107:       endif
108:  120  format ('  ',F8.5)
109:  130  format ('  ',F8.5,SP,F8.5,'i')

111: !     *** Clean up
112:       PetscCallA(DSDestroy(ds,ierr))
113:       PetscCallA(SlepcFinalize(ierr))
114:       end

116: ! -----------------------------------------------------------------

118:       subroutine FillUpMatrix(n,X)
119:       PetscInt    n,i,j
120:       PetscScalar X(n,n)

122:       do i=2,n
123:         X(i,i-1) = -1.d0
124:       end do
125:       do j=0,3
126:         do i=1,n-j
127:           X(i,i+j) = 1.d0
128:         end do
129:       end do

131:       end

133: !/*TEST
134: !
135: !   test:
136: !      suffix: 1
137: !      requires: !complex
138: !
139: !TEST*/