Actual source code: ex61f.F90

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: !
  2: !        Demonstrates having each OpenMP thread manage its own PETSc objects and solves
  3: !           - each thread is ONLY allowed to access objects that IT created
  4: !                  that is, threads cannot shared objects
  5: !
  6: !        Run with "export OMP_NUM_THREADS=16 ./ex61f"
  7: !           to use 16 independent threads
  8: !
  9: !        ./configure --with-threadsafety --with-log=0 --with-openmp
 10: !
 11:          module omp_module
 12:          implicit none
 13:          contains
 14:          subroutine split_indices(total,num_pieces,ibeg,iend)
 15:            implicit none

 17:            integer :: total
 18:            integer :: num_pieces
 19:            integer :: ibeg(num_pieces), iend(num_pieces)
 20:            integer :: itmp1, itmp2, ioffset, i

 22:            itmp1 = total/num_pieces
 23:            itmp2 = mod(total,num_pieces)
 24:            ioffset = 0
 25:            do i=1,itmp2
 26:               ibeg(i) = ioffset + 1
 27:               iend(i) = ioffset + (itmp1+1)
 28:               ioffset = iend(i)
 29:            enddo
 30:            do i=itmp2+1,num_pieces
 31:               ibeg(i) = ioffset + 1
 32:               if (ibeg(i) > total) then
 33:                  iend(i) = ibeg(i) - 1
 34:               else
 35:                  iend(i) = ioffset + itmp1
 36:                  ioffset = iend(i)
 37:               endif
 38:            enddo

 40:          end subroutine split_indices
 41:        end module omp_module

 43:       module assert_mod
 44:       implicit none
 45:       contains
 46:       subroutine assert(lcond,msg,icase)
 47:       logical,intent(in) :: lcond
 48:       character(len=*), intent(in) :: msg
 49:       integer, intent(in) :: icase

 51:       if (.not.lcond) then
 52:          write(*,*) msg, icase
 53:          stop 'assertion error '
 54:       endif
 55:       return
 56:       end subroutine assert
 57:       end module assert_mod

 59:       program tpetsc
 60:       use assert_mod
 61:       use omp_module
 62:       implicit none
 63: !     ----------------------------
 64: !     test concurrent petsc solver
 65: !     ----------------------------

 67:       integer, parameter :: NCASES=100
 68:       integer, parameter :: MAXTHREADS=64
 69:       real(8), parameter :: tol = 1.0d-6

 71:       integer, dimension(MAXTHREADS) :: ibeg,iend

 73: #include <petsc/finclude/petsc.h90>

 75: !$   integer, external :: omp_get_num_threads

 77:       Mat, target :: Mcol_f_mat(MAXTHREADS)
 78:       Vec, target :: Mcol_f_vecb(MAXTHREADS)
 79:       Vec, target :: Mcol_f_vecx(MAXTHREADS)
 80:       KSP, target :: Mcol_f_ksp(MAXTHREADS)
 81:       PC,  target :: MPC(MAXTHREADS)

 83:       Mat, pointer :: col_f_mat
 84:       Vec, pointer :: col_f_vecb
 85:       Vec, pointer :: col_f_vecx
 86:       KSP, pointer :: col_f_ksp
 87:       PC, pointer :: pc

 89:       PetscInt :: ith, nthreads
 90:       PetscErrorCode ierr

 92:       integer, parameter :: nz_per_row = 9
 93:       integer, parameter :: n =100
 94:       integer :: i,j,ij,ij2,ii,jj,nz,ip, dx,dy,icase
 95:       integer, dimension(n*n*nz_per_row) :: ilist,jlist
 96:       PetscScalar :: aij
 97:       PetscScalar, dimension(n*n*nz_per_row) :: alist
 98:       logical :: isvalid_ii, isvalid_jj, is_diag

100:       PetscInt nrow
101:       PetscInt ncol
102:       PetscInt blocksize
103:       PetscScalar, dimension(0:(n*n-1)) :: x, b
104:       real(8) :: err(NCASES)

106:        integer :: t1,t2,count_rate
107:        real :: ttime



111:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
112:       call assert(ierr.eq.0,'PetscInitialize return ',ierr)

114:       nrow = n*n
115:       ncol = nrow

117:       nthreads = 1
118: !$omp parallel
119: !$omp master
120: !$      nthreads = omp_get_num_threads()
121: !$omp end master
122: !$omp end parallel
123:       print*,'nthreads = ', nthreads,' NCASES = ', NCASES

125:       call split_indices(NCASES,nthreads,ibeg,iend)


128: !$omp parallel do                                                        &
129: !$omp private(ith,ierr)                                                  &
130: !$omp private(col_f_mat,col_f_vecb,col_f_vecx,col_f_ksp)
131:       do ith=1,nthreads
132:          col_f_mat => Mcol_f_mat(ith)
133:          col_f_vecb => Mcol_f_vecb(ith)
134:          col_f_vecx => Mcol_f_vecx(ith)
135:          col_f_ksp => Mcol_f_ksp(ith)


138:          call MatCreateSeqAIJ( PETSC_COMM_SELF, nrow,ncol, nz_per_row,     &
139:      &           PETSC_NULL_INTEGER, col_f_mat, ierr)
140:          call assert(ierr.eq.0,'matcreateseqaij return ',ierr)

142:          call VecCreateSeqWithArray(PETSC_COMM_SELF,1,nrow,       &
143:      &           PETSC_NULL_SCALAR, col_f_vecb, ierr)
144:          call assert(ierr.eq.0,'veccreateseqwitharray return ierr',ierr)

146:          call VecDuplicate(col_f_vecb, col_f_vecx,ierr)
147:          call assert(ierr.eq.0,'vecduplicate return ierr',ierr)

149:          call KSPCreate(PETSC_COMM_SELF, col_f_ksp,ierr)
150:          call assert(ierr.eq.0,'kspcreate return ierr',ierr)

152:        enddo

154: !      -----------------------
155: !      setup sparsity pattern
156: !      -----------------------
157:        nz = 0
158:        do j=1,n
159:        do i=1,n
160:         ij = i + (j-1)*n
161:         do dx=-1,1
162:         do dy=-1,1
163:           ii = i + dx
164:           jj = j + dy

166:           ij2 = ii + (jj-1)*n
167:           isvalid_ii = (1 <= ii).and.(ii <= n)
168:           isvalid_jj = (1 <= jj).and.(jj <= n)
169:           if (isvalid_ii.and.isvalid_jj) then
170:              is_diag = (ij .eq. ij2)
171:              nz = nz + 1
172:              ilist(nz) = ij
173:              jlist(nz) = ij2
174:              if (is_diag) then
175:                aij = 4.0d0
176:              else
177:                aij = -1.0d0
178:              endif
179:              alist(nz) = aij
180:            endif
181:           enddo
182:           enddo
183:          enddo
184:          enddo

186:        print*,'nz = ', nz

188: !      ---------------------------------
189: !      convert from fortan to c indexing
190: !      ---------------------------------
191:        ilist(1:nz) = ilist(1:nz) - 1
192:        jlist(1:nz) = jlist(1:nz) - 1


195: !      --------------
196: !      setup matrices
197: !      --------------
198:        call system_clock(t1,count_rate)

200: !$omp  parallel do                                                      &
201: !$omp& private(ith,icase,ip,i,j,ii,jj,aij,ierr,x,b)                      &
202: !$omp& private(col_f_mat,col_f_vecb,col_f_vecx,col_f_ksp,pc)
203:        do ith=1,nthreads
204:          col_f_mat => Mcol_f_mat(ith)
205:          col_f_vecb => Mcol_f_vecb(ith)
206:          col_f_vecx => Mcol_f_vecx(ith)
207:          col_f_ksp => Mcol_f_ksp(ith)
208:          pc => MPC(ith)

210:         do icase=ibeg(ith),iend(ith)

212:          do ip=1,nz
213:            ii = ilist(ip)
214:            jj = jlist(ip)
215:            aij = alist(ip)
216:            call MatSetValues(col_f_mat,1,ii,1,jj,aij,INSERT_VALUES,ierr)
217:            call assert(ierr.eq.0,'matsetvalues return ierr',ierr)
218:          enddo
219:          call MatAssemblyBegin(col_f_mat,MAT_FINAL_ASSEMBLY,ierr)
220:          call assert(ierr.eq.0,'matassemblybegin return ierr',ierr)

222:          call MatAssemblyEnd(col_f_mat,MAT_FINAL_ASSEMBLY,ierr)
223:          call assert(ierr.eq.0,'matassemblyend return ierr',ierr)

225:          call KSPSetOperators(col_f_ksp,col_f_mat,col_f_mat,ierr)
226:          call assert(ierr.eq.0,'KSPSetOperators return ierr',ierr)

228:          ! set linear solver
229:          call KSPGetPC(col_f_ksp,pc,ierr)
230:          call assert(ierr.eq.0,'KSPGetPC return ierr ', ierr)

232:          call PCSetType(pc,PCLU,ierr)
233:          call assert(ierr.eq.0,'PCSetType return ierr ',ierr)

235:          ! associate petsc vector with allocated array
236:          x(:) = 1
237: !$omp    critical
238:          call VecPlaceArray(col_f_vecx,x,ierr)
239: !$omp    end critical
240:          call assert(ierr.eq.0,'VecPlaceArray col_f_vecx return ',ierr)

242:          b(:) = 0
243:          do ip=1,nz
244:            i = ilist(ip)
245:            j = jlist(ip)
246:            aij = alist(ip)
247:            b(i) = b(i) + aij*x(j)
248:          enddo
249:          x = 0
250: !$omp    critical
251:          call VecPlaceArray(col_f_vecb,b,ierr)
252: !$omp    end critical
253:          call assert(ierr.eq.0,'VecPlaceArray col_f_vecb return ',ierr)



257: !  -----------------------------------------------------------
258: !  key test, need to solve multiple linear systems in parallel
259: !  -----------------------------------------------------------
260:          call KSPSetFromOptions(col_f_ksp,ierr)
261:          call assert(ierr.eq.0,'KSPSetFromOptions return ierr ',ierr)

263:          call KSPSetUp(col_f_ksp,ierr)
264:          call assert(ierr.eq.0,'KSPSetup return ierr ',ierr)


267:          call KSPSolve(col_f_ksp,col_f_vecb,col_f_vecx,ierr)
268:          call assert(ierr.eq.0,'KSPSolve return ierr ',ierr)


271: !        ------------
272: !        check answer
273: !        ------------
274:          err(icase) = maxval(abs(x(:)-1.0d0))

276: !$omp    critical
277:          call VecResetArray(col_f_vecx,ierr)
278: !$omp    end critical
279:          call assert(ierr.eq.0,'VecResetArray return ierr ',ierr)

281: !$omp    critical
282:          call VecResetArray(col_f_vecb,ierr)
283: !$omp    end critical
284:          call assert(ierr.eq.0,'VecResetArray return ierr ',ierr)

286:        enddo
287:        enddo

289: !$omp parallel do                                                        &
290: !$omp private(ith,ierr)                                                  &
291: !$omp private(col_f_mat,col_f_vecb,col_f_vecx,col_f_ksp)
292:       do ith=1,nthreads
293:          col_f_mat => Mcol_f_mat(ith)
294:          col_f_vecb => Mcol_f_vecb(ith)
295:          col_f_vecx => Mcol_f_vecx(ith)
296:          col_f_ksp => Mcol_f_ksp(ith)


299:          call MatDestroy(col_f_mat, ierr)
300:          call assert(ierr.eq.0,'matdestroy return ',ierr)

302:          call VecDestroy(col_f_vecb, ierr)
303:          call assert(ierr.eq.0,'vecdestroy return ierr',ierr)

305:          call VecDestroy(col_f_vecx,ierr)
306:          call assert(ierr.eq.0,'vecdestroy return ierr',ierr)

308:          call KSPDestroy(col_f_ksp,ierr)
309:          call assert(ierr.eq.0,'kspdestroy return ierr',ierr)

311:        enddo

313:        call system_clock(t2,count_rate)
314:        ttime = real(t2-t1)/real(count_rate)
315:        write(*,*) 'total time is ',ttime

317:        write(*,*) 'maxval(err) ', maxval(err)
318:        do icase=1,NCASES
319:         if (err(icase) > tol) then
320:          write(*,*) 'icase,err(icase) ',icase,err(icase)
321:         endif
322:        enddo

324:        call PetscFinalize(ierr)
325:        call assert(ierr.eq.0,'petscfinalize return ierr',ierr)

327:        stop 'all done'
328:        end program tpetsc