Actual source code: ex44f.F90

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1:       program main              !   Solves the linear system  J x = f
  2:  #include <petsc/finclude/petsc.h>
  3:       use petscksp
  4:       implicit none
  5:       Vec x,f
  6:       Mat J
  7:       DM da
  8:       KSP ksp
  9:       PetscErrorCode ierr
 10:       PetscInt eight,one

 12:       eight = 8
 13:       one = 1
 14:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 15:       if (ierr .ne. 0) then
 16:         print*,'Unable to initialize PETSc'
 17:         stop
 18:       endif
 19:       call DMDACreate1d(MPI_COMM_WORLD,DM_BOUNDARY_NONE,eight,one,one,PETSC_NULL_INTEGER,da,ierr);CHKERRA(ierr)
 20:       call DMSetFromOptions(da,ierr)
 21:       call DMSetUp(da,ierr)
 22:       call DMCreateGlobalVector(da,x,ierr);CHKERRA(ierr)
 23:       call VecDuplicate(x,f,ierr);CHKERRA(ierr)
 24:       call DMSetMatType(da,MATAIJ,ierr);CHKERRA(ierr)
 25:       call DMCreateMatrix(da,J,ierr);CHKERRA(ierr)

 27:       call ComputeRHS(da,f,ierr);CHKERRA(ierr)
 28:       call ComputeMatrix(da,J,ierr);CHKERRA(ierr)

 30:       call KSPCreate(MPI_COMM_WORLD,ksp,ierr);CHKERRA(ierr)
 31:       call KSPSetOperators(ksp,J,J,ierr);CHKERRA(ierr)
 32:       call KSPSetFromOptions(ksp,ierr);CHKERRA(ierr)
 33:       call KSPSolve(ksp,f,x,ierr);CHKERRA(ierr)

 35:       call MatDestroy(J,ierr);CHKERRA(ierr)
 36:       call VecDestroy(x,ierr);CHKERRA(ierr)
 37:       call VecDestroy(f,ierr);CHKERRA(ierr)
 38:       call KSPDestroy(ksp,ierr);CHKERRA(ierr)
 39:       call DMDestroy(da,ierr);CHKERRA(ierr)
 40:       call PetscFinalize(ierr);CHKERRA(ierr)
 41:       end

 43: ! AVX512 crashes without this..
 44:       subroutine knl_workarround(xx)
 45:       PetscScalar xx,sd
 46:       common /cb/ sd
 47:       data sd /0/
 48:       sd = sd+xx
 49:       end

 51:       subroutine  ComputeRHS(da,x,ierr)
 52:       use petscdmda
 53:       implicit none
 54:       DM da
 55:       Vec x
 56:       PetscErrorCode ierr
 57:       PetscInt xs,xm,i,mx
 58:       PetscScalar hx
 59:       PetscScalar, pointer :: xx(:)
 60:       call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
 61:      &     PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
 62:      &     PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr)
 63:       call DMDAGetCorners(da,xs,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,xm,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr)
 64:       hx     = 1.0_PETSC_REAL_KIND/(mx-1)
 65:       call DMDAVecGetArrayF90(da,x,xx,ierr);CHKERRQ(ierr)
 66:       do i=xs,xs+xm-1
 67:         call knl_workarround(xx(i))
 68:         xx(i) = i*hx
 69:       enddo
 70:       call DMDAVecRestoreArrayF90(da,x,xx,ierr);CHKERRQ(ierr)
 71:       return
 72:       end

 74:       subroutine ComputeMatrix(da,J,ierr)
 75:       use petscdm
 76:       use petscmat
 77:       implicit none
 78:       Mat J
 79:       DM da
 80:       PetscErrorCode ierr
 81:       PetscInt xs,xm,i,mx
 82:       PetscScalar hx,one

 84:       one = 1.0
 85:       call DMDAGetInfo(da,PETSC_NULL_INTEGER,mx,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER, &
 86:      &  PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,   &
 87:      &  PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr)
 88:       call DMDAGetCorners(da,xs,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,xm,PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,ierr);CHKERRQ(ierr)
 89:       hx     = 1.0_PETSC_REAL_KIND/(mx-1)
 90:       do i=xs,xs+xm-1
 91:         if ((i .eq. 0) .or. (i .eq. mx-1)) then
 92:           call MatSetValue(J,i,i,one,INSERT_VALUES,ierr);CHKERRQ(ierr)
 93:         else
 94:           call MatSetValue(J,i,i-1,-hx,INSERT_VALUES,ierr);CHKERRQ(ierr)
 95:           call MatSetValue(J,i,i+1,-hx,INSERT_VALUES,ierr);CHKERRQ(ierr)
 96:           call MatSetValue(J,i,i,2*hx,INSERT_VALUES,ierr);CHKERRQ(ierr)
 97:         endif
 98:       enddo
 99:       call MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
100:       call MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY,ierr);CHKERRQ(ierr)
101:       return
102:       end