Actual source code: test1f.F90

slepc-3.21.2 2024-09-25
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> ./test1f [-help]
 11: !
 12: !  Description: Simple example that tests BV interface functions.
 13: !
 14: ! ----------------------------------------------------------------------
 15: !
 16:       program main
 17: #include <slepc/finclude/slepcbv.h>
 18:       use slepcbv
 19:       implicit none

 21: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: !     Declarations
 23: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 25: #define KMAX 35

 27:       Vec            t,v
 28:       Mat            Q,M
 29:       BV             X,Y;
 30:       PetscMPIInt    rank
 31:       PetscInt       i,j,n,k,l,izero,ione
 32:       PetscScalar    z(KMAX),val
 33:       PetscScalar, pointer :: qq(:,:)
 34:       PetscScalar    one,mone,two,zero
 35:       PetscReal      nrm
 36:       PetscBool      flg
 37:       PetscErrorCode ierr

 39: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 40: !     Beginning of program
 41: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 43:       n = 10
 44:       k = 5
 45:       l = 3
 46:       one = 1.0
 47:       mone = -1.0
 48:       two = 2.0
 49:       zero = 0.0
 50:       izero = 0
 51:       ione = 1
 52:       PetscCallA(SlepcInitialize(PETSC_NULL_CHARACTER,ierr))
 53:       if (ierr .ne. 0) then
 54:         print*,'SlepcInitialize failed'
 55:         stop
 56:       endif
 57:       PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
 58:       PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))
 59:       PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-k',k,flg,ierr))
 60:       PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-l',l,flg,ierr))
 61:       if (k .gt. KMAX) then; SETERRA(PETSC_COMM_SELF,1,'Program currently limited to k=35'); endif
 62:       if (rank .eq. 0) then
 63:         write(*,110) k,n
 64:       endif
 65:  110  format (/'Test BV with',I3,' columns of length',I3,' (Fortran)')

 67: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68: !     Initialize data
 69: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

 71: !     ** Create template vector
 72:       PetscCallA(VecCreate(PETSC_COMM_WORLD,t,ierr))
 73:       PetscCallA(VecSetSizes(t,PETSC_DECIDE,n,ierr))
 74:       PetscCallA(VecSetFromOptions(t,ierr))

 76: !     ** Create BV object X
 77:       PetscCallA(BVCreate(PETSC_COMM_WORLD,X,ierr))
 78:       PetscCallA(BVSetSizesFromVec(X,t,k,ierr))
 79:       PetscCallA(BVSetFromOptions(X,ierr))

 81: !     ** Fill X entries
 82:       do j=0,k-1
 83:         PetscCallA(BVGetColumn(X,j,v,ierr))
 84:         PetscCallA(VecSet(v,zero,ierr))
 85:         do i=0,3
 86:           if (i+j<n) then
 87:             val = 3*i+j-2
 88:             PetscCallA(VecSetValue(v,i+j,val,INSERT_VALUES,ierr))
 89:           end if
 90:         end do
 91:         PetscCallA(VecAssemblyBegin(v,ierr))
 92:         PetscCallA(VecAssemblyEnd(v,ierr))
 93:         PetscCallA(BVRestoreColumn(X,j,v,ierr))
 94:       end do

 96: !     ** Create BV object Y
 97:       PetscCallA(BVCreate(PETSC_COMM_WORLD,Y,ierr))
 98:       PetscCallA(BVSetSizesFromVec(Y,t,l,ierr))
 99:       PetscCallA(BVSetFromOptions(Y,ierr))

101: !     ** Fill Y entries
102:       do j=0,l-1
103:         PetscCallA(BVGetColumn(Y,j,v,ierr))
104:         val = real(j+1)/4.0
105:         PetscCallA(VecSet(v,val,ierr))
106:         PetscCallA(BVRestoreColumn(Y,j,v,ierr))
107:       end do

109: !     ** Create Mat
110:       PetscCallA(MatCreateSeqDense(PETSC_COMM_SELF,k,l,PETSC_NULL_SCALAR,Q,ierr))
111:       PetscCallA(MatDenseGetArrayF90(Q,qq,ierr))
112:       do i=1,k
113:         do j=1,l
114:           if (i<j) then
115:             qq(i,j) = 2.0
116:           else
117:             qq(i,j) = -0.5
118:           end if
119:         end do
120:       end do
121:       PetscCallA(MatDenseRestoreArrayF90(Q,qq,ierr))

123: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: !     Test several operations
125: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

127: !     ** Test BVMult
128:       PetscCallA(BVMult(Y,two,one,X,Q,ierr))

130: !     ** Test BVMultVec
131:       PetscCallA(BVGetColumn(Y,izero,v,ierr))
132:       z(1) = 2.0
133:       do i=2,k
134:         z(i) = -0.5*z(i-1)
135:       end do
136:       PetscCallA(BVMultVec(X,mone,one,v,z,ierr))
137:       PetscCallA(BVRestoreColumn(Y,izero,v,ierr))

139: !     ** Test BVDot
140:       PetscCallA(MatCreateSeqDense(PETSC_COMM_SELF,l,k,PETSC_NULL_SCALAR,M,ierr))
141:       PetscCallA(BVDot(X,Y,M,ierr))

143: !     ** Test BVDotVec
144:       PetscCallA(BVGetColumn(Y,izero,v,ierr))
145:       PetscCallA(BVDotVec(X,v,z,ierr))
146:       PetscCallA(BVRestoreColumn(Y,izero,v,ierr))

148: !     ** Test BVMultInPlace and BVScale
149:       PetscCallA(BVMultInPlace(X,Q,ione,l,ierr))
150:       PetscCallA(BVScale(X,two,ierr))

152: !     ** Test BVNorm
153:       PetscCallA(BVNormColumn(X,izero,NORM_2,nrm,ierr))
154:       if (rank .eq. 0) then
155:         write(*,120) nrm
156:       endif
157:  120  format ('2-Norm of X[0] = ',f8.4)
158:       PetscCallA(BVNorm(X,NORM_FROBENIUS,nrm,ierr))
159:       if (rank .eq. 0) then
160:         write(*,130) nrm
161:       endif
162:  130  format ('Frobenius Norm of X = ',f8.4)

164: !     *** Clean up
165:       PetscCallA(BVDestroy(X,ierr))
166:       PetscCallA(BVDestroy(Y,ierr))
167:       PetscCallA(VecDestroy(t,ierr))
168:       PetscCallA(MatDestroy(Q,ierr))
169:       PetscCallA(MatDestroy(M,ierr))
170:       PetscCallA(SlepcFinalize(ierr))
171:       end

173: !/*TEST
174: !
175: !   test:
176: !      suffix: 1
177: !      nsize: 1
178: !      args: -bv_type {{vecs contiguous svec mat}separate output}
179: !      output_file: output/test1f_1.out
180: !
181: !   test:
182: !      suffix: 2
183: !      nsize: 2
184: !      args: -bv_type {{vecs contiguous svec mat}separate output}
185: !      output_file: output/test1f_1.out
186: !
187: !TEST*/